1*ba92ff59SBarry Smith #include <petsctao.h> /*I "petsctao.h" I*/ 2aaa7dc30SBarry Smith #include <petsc-private/matimpl.h> 3aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h> 4a7e14dcfSSatish Balay 5a7e14dcfSSatish Balay #undef __FUNCT__ 6a7e14dcfSSatish Balay #define __FUNCT__ "VecGetSubVec" 7a7e14dcfSSatish Balay /*@ 8a7e14dcfSSatish Balay VecGetSubVec - Gets a subvector using the IS 9a7e14dcfSSatish Balay 10a7e14dcfSSatish Balay Input Parameters: 11a7e14dcfSSatish Balay + vfull - the full matrix 12a7e14dcfSSatish Balay . is - the index set for the subvector 13a7e14dcfSSatish Balay . reduced_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, TAO_SUBSET_MATRIXFREE) 14a7e14dcfSSatish Balay - maskvalue - the value to set the unused vector elements to (for TAO_SUBSET_MASK or TAO_SUBSET_MATRIXFREE) 15a7e14dcfSSatish Balay 16a7e14dcfSSatish Balay Output Parameters: 17a7e14dcfSSatish Balay . vreduced - the subvector 18a7e14dcfSSatish Balay 19a7e14dcfSSatish Balay Note: 20a7e14dcfSSatish Balay maskvalue should usually be 0.0, unless a pointwise divide will be used. 21a7e14dcfSSatish Balay @*/ 22a7e14dcfSSatish Balay PetscErrorCode VecGetSubVec(Vec vfull, IS is, PetscInt reduced_type, PetscReal maskvalue, Vec *vreduced) 23a7e14dcfSSatish Balay { 24a7e14dcfSSatish Balay PetscErrorCode ierr; 25a7e14dcfSSatish Balay PetscInt nfull,nreduced,nreduced_local,rlow,rhigh,flow,fhigh; 26a7e14dcfSSatish Balay PetscInt i,nlocal; 27a7e14dcfSSatish Balay PetscReal *fv,*rv; 28a7e14dcfSSatish Balay const PetscInt *s; 29a7e14dcfSSatish Balay IS ident; 30a7e14dcfSSatish Balay VecType vtype; 31a7e14dcfSSatish Balay VecScatter scatter; 32a7e14dcfSSatish Balay MPI_Comm comm; 33a7e14dcfSSatish Balay 34a7e14dcfSSatish Balay PetscFunctionBegin; 35a7e14dcfSSatish Balay PetscValidHeaderSpecific(vfull,VEC_CLASSID,1); 36a7e14dcfSSatish Balay PetscValidHeaderSpecific(is,IS_CLASSID,2); 37a7e14dcfSSatish Balay 38a7e14dcfSSatish Balay ierr = VecGetSize(vfull, &nfull);CHKERRQ(ierr); 39a7e14dcfSSatish Balay ierr = ISGetSize(is, &nreduced);CHKERRQ(ierr); 40a7e14dcfSSatish Balay 41a7e14dcfSSatish Balay if (nreduced == nfull) { 42a7e14dcfSSatish Balay ierr = VecDestroy(vreduced);CHKERRQ(ierr); 43a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr); 44a7e14dcfSSatish Balay ierr = VecCopy(vfull,*vreduced);CHKERRQ(ierr); 45a7e14dcfSSatish Balay } else { 46a7e14dcfSSatish Balay switch (reduced_type) { 47a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 48a7e14dcfSSatish Balay ierr = VecGetType(vfull,&vtype);CHKERRQ(ierr); 49a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr); 50a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nreduced_local);CHKERRQ(ierr); 51a7e14dcfSSatish Balay ierr = PetscObjectGetComm((PetscObject)vfull,&comm);CHKERRQ(ierr); 52a7e14dcfSSatish Balay if (*vreduced) { 53a7e14dcfSSatish Balay ierr = VecDestroy(vreduced);CHKERRQ(ierr); 54a7e14dcfSSatish Balay } 55a7e14dcfSSatish Balay ierr = VecCreate(comm,vreduced);CHKERRQ(ierr); 56a7e14dcfSSatish Balay ierr = VecSetType(*vreduced,vtype);CHKERRQ(ierr); 57a7e14dcfSSatish Balay 58a7e14dcfSSatish Balay ierr = VecSetSizes(*vreduced,nreduced_local,nreduced);CHKERRQ(ierr); 59a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(*vreduced,&rlow,&rhigh);CHKERRQ(ierr); 60a7e14dcfSSatish Balay ierr = ISCreateStride(comm,nreduced_local,rlow,1,&ident);CHKERRQ(ierr); 61a7e14dcfSSatish Balay ierr = VecScatterCreate(vfull,is,*vreduced,ident,&scatter);CHKERRQ(ierr); 62a7e14dcfSSatish Balay ierr = VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 63a7e14dcfSSatish Balay ierr = VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = VecScatterDestroy(&scatter);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = ISDestroy(&ident);CHKERRQ(ierr); 66a7e14dcfSSatish Balay break; 67a7e14dcfSSatish Balay 68a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 69a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 70a7e14dcfSSatish Balay /* vr[i] = vf[i] if i in is 71a7e14dcfSSatish Balay vr[i] = 0 otherwise */ 726c23d075SBarry Smith if (*vreduced == NULL) { 73a7e14dcfSSatish Balay ierr = VecDuplicate(vfull,vreduced);CHKERRQ(ierr); 74a7e14dcfSSatish Balay } 75a7e14dcfSSatish Balay 76a7e14dcfSSatish Balay ierr = VecSet(*vreduced,maskvalue);CHKERRQ(ierr); 77a7e14dcfSSatish Balay ierr = ISGetLocalSize(is,&nlocal);CHKERRQ(ierr); 78a7e14dcfSSatish Balay ierr = VecGetOwnershipRange(vfull,&flow,&fhigh);CHKERRQ(ierr); 79a7e14dcfSSatish Balay ierr = VecGetArray(vfull,&fv);CHKERRQ(ierr); 80a7e14dcfSSatish Balay ierr = VecGetArray(*vreduced,&rv);CHKERRQ(ierr); 81a7e14dcfSSatish Balay ierr = ISGetIndices(is,&s);CHKERRQ(ierr); 8253506e15SBarry Smith if (nlocal > (fhigh-flow)) SETERRQ2(PETSC_COMM_WORLD,1,"IS local size %d > Vec local size %d",nlocal,fhigh-flow); 83a7e14dcfSSatish Balay for (i=0;i<nlocal;i++) { 84a7e14dcfSSatish Balay rv[s[i]-flow] = fv[s[i]-flow]; 85a7e14dcfSSatish Balay } 86a7e14dcfSSatish Balay ierr = ISRestoreIndices(is,&s);CHKERRQ(ierr); 87a7e14dcfSSatish Balay ierr = VecRestoreArray(vfull,&fv);CHKERRQ(ierr); 88a7e14dcfSSatish Balay ierr = VecRestoreArray(*vreduced,&rv);CHKERRQ(ierr); 89a7e14dcfSSatish Balay break; 90a7e14dcfSSatish Balay } 91a7e14dcfSSatish Balay } 92a7e14dcfSSatish Balay PetscFunctionReturn(0); 93a7e14dcfSSatish Balay } 94a7e14dcfSSatish Balay 95a7e14dcfSSatish Balay #undef __FUNCT__ 96a7e14dcfSSatish Balay #define __FUNCT__ "MatGetSubMat" 97a7e14dcfSSatish Balay /*@ 98a7e14dcfSSatish Balay MatGetSubMat - Gets a submatrix using the IS 99a7e14dcfSSatish Balay 100a7e14dcfSSatish Balay Input Parameters: 101a7e14dcfSSatish Balay + M - the full matrix (n x n) 102a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same) 103a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option 104a7e14dcfSSatish Balay - subset_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK, 105a7e14dcfSSatish Balay TAO_SUBSET_MATRIXFREE) 106a7e14dcfSSatish Balay 107a7e14dcfSSatish Balay Output Parameters: 108a7e14dcfSSatish Balay . Msub - the submatrix 109a7e14dcfSSatish Balay @*/ 110a7e14dcfSSatish Balay PetscErrorCode MatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub) 111a7e14dcfSSatish Balay { 112a7e14dcfSSatish Balay PetscErrorCode ierr; 113a7e14dcfSSatish Balay IS iscomp; 114a7e14dcfSSatish Balay PetscBool flg; 11553506e15SBarry Smith 116a7e14dcfSSatish Balay PetscFunctionBegin; 117a7e14dcfSSatish Balay PetscValidHeaderSpecific(M,MAT_CLASSID,1); 118a7e14dcfSSatish Balay PetscValidHeaderSpecific(is,IS_CLASSID,2); 119a7e14dcfSSatish Balay ierr = MatDestroy(Msub);CHKERRQ(ierr); 120a7e14dcfSSatish Balay switch (subset_type) { 121a7e14dcfSSatish Balay case TAO_SUBSET_SUBVEC: 122a7e14dcfSSatish Balay ierr = MatGetSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub);CHKERRQ(ierr); 123a7e14dcfSSatish Balay break; 124a7e14dcfSSatish Balay 125a7e14dcfSSatish Balay case TAO_SUBSET_MASK: 126a7e14dcfSSatish Balay /* Get Reduced Hessian 127a7e14dcfSSatish Balay Msub[i,j] = M[i,j] if i,j in Free_Local or i==j 128a7e14dcfSSatish Balay Msub[i,j] = 0 if i!=j and i or j not in Free_Local 129a7e14dcfSSatish Balay */ 1306c23d075SBarry Smith ierr = PetscOptionsBool("-different_submatrix","use separate hessian matrix when computing submatrices","TaoSubsetType",PETSC_FALSE,&flg,NULL); 131a7e14dcfSSatish Balay if (flg == PETSC_TRUE) { 132a7e14dcfSSatish Balay ierr = MatDuplicate(M, MAT_COPY_VALUES, Msub);CHKERRQ(ierr); 133a7e14dcfSSatish Balay } else { 134a7e14dcfSSatish Balay /* Act on hessian directly (default) */ 135a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)M);CHKERRQ(ierr); 136a7e14dcfSSatish Balay *Msub = M; 137a7e14dcfSSatish Balay } 138a7e14dcfSSatish Balay /* Save the diagonal to temporary vector */ 139a7e14dcfSSatish Balay ierr = MatGetDiagonal(*Msub,v1);CHKERRQ(ierr); 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay /* Zero out rows and columns */ 1424473680cSBarry Smith ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr); 143a7e14dcfSSatish Balay 144a7e14dcfSSatish Balay /* Use v1 instead of 0 here because of PETSc bug */ 145a7e14dcfSSatish Balay ierr = MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1);CHKERRQ(ierr); 146a7e14dcfSSatish Balay 147a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp);CHKERRQ(ierr); 148a7e14dcfSSatish Balay break; 149a7e14dcfSSatish Balay case TAO_SUBSET_MATRIXFREE: 1504473680cSBarry Smith ierr = ISComplementVec(is,v1,&iscomp);CHKERRQ(ierr); 151a7e14dcfSSatish Balay ierr = MatCreateSubMatrixFree(M,iscomp,iscomp,Msub);CHKERRQ(ierr); 152a7e14dcfSSatish Balay ierr = ISDestroy(&iscomp);CHKERRQ(ierr); 153a7e14dcfSSatish Balay break; 154a7e14dcfSSatish Balay } 155a7e14dcfSSatish Balay PetscFunctionReturn(0); 156a7e14dcfSSatish Balay } 157