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