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