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