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