xref: /petsc/src/tao/bound/utils/isutil.c (revision f89ca46fb01025fa5f21ef09d10cb4723982ea5b)
1a7e14dcfSSatish Balay #include "tao.h" /*I "tao.h" I*/
2a7e14dcfSSatish Balay #include "petsc-private/matimpl.h"
3*f89ca46fSSatish Balay #include "../src/tao/matrix/submatfree.h"
4a7e14dcfSSatish Balay #include "tao_util.h" /*I "tao_util.h" I*/
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay 
7a7e14dcfSSatish Balay #undef __FUNCT__
8a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichEqual"
9a7e14dcfSSatish Balay /*@
10a7e14dcfSSatish Balay   VecWhichEqual - Creates an index set containing the indices
11a7e14dcfSSatish Balay   where the vectors Vec1 and Vec2 have identical elements.
12a7e14dcfSSatish Balay 
13a7e14dcfSSatish Balay   Collective on S
14a7e14dcfSSatish Balay 
15a7e14dcfSSatish Balay   Input Parameters:
16a7e14dcfSSatish Balay . Vec1, Vec2 - the two vectors to compare
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay   OutputParameter:
19a7e14dcfSSatish Balay . S - The index set containing the indices i where vec1[i] == vec2[i]
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay   Level: advanced
22a7e14dcfSSatish Balay @*/
23a7e14dcfSSatish Balay PetscErrorCode VecWhichEqual(Vec Vec1, Vec Vec2, IS * S)
24a7e14dcfSSatish Balay {
25a7e14dcfSSatish Balay   PetscErrorCode    ierr;
26a7e14dcfSSatish Balay   PetscInt i,n_same=0;
27a7e14dcfSSatish Balay   PetscInt n,low,high,low2,high2;
28a7e14dcfSSatish Balay   PetscInt    *same;
29a7e14dcfSSatish Balay   PetscReal *v1,*v2;
30a7e14dcfSSatish Balay   MPI_Comm comm;
31a7e14dcfSSatish Balay 
32a7e14dcfSSatish Balay   PetscFunctionBegin;
33a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec1,VEC_CLASSID,1);
34a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec2,VEC_CLASSID,2);
35a7e14dcfSSatish Balay   PetscCheckSameComm(Vec1,1,Vec2,2);
36a7e14dcfSSatish Balay 
37a7e14dcfSSatish Balay 
38a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(ierr);
39a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(ierr);
40a7e14dcfSSatish Balay 
41a7e14dcfSSatish Balay   if ( low != low2 || high != high2 )
42a7e14dcfSSatish Balay     SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay   ierr = VecGetLocalSize(Vec1,&n); CHKERRQ(ierr);
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay   if (n>0){
47a7e14dcfSSatish Balay 
48a7e14dcfSSatish Balay     if (Vec1 == Vec2){
49a7e14dcfSSatish Balay       ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr);
50a7e14dcfSSatish Balay       v2=v1;
51a7e14dcfSSatish Balay     } else {
52a7e14dcfSSatish Balay       ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr);
53a7e14dcfSSatish Balay       ierr = VecGetArray(Vec2,&v2); CHKERRQ(ierr);
54a7e14dcfSSatish Balay     }
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay     ierr = PetscMalloc( n*sizeof(PetscInt),&same ); CHKERRQ(ierr);
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay     for (i=0; i<n; i++){
59a7e14dcfSSatish Balay       if (v1[i] == v2[i]) {same[n_same]=low+i; n_same++;}
60a7e14dcfSSatish Balay     }
61a7e14dcfSSatish Balay 
62a7e14dcfSSatish Balay     if (Vec1 == Vec2){
63a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr);
64a7e14dcfSSatish Balay     } else {
65a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr);
66a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec2,&v2); CHKERRQ(ierr);
67a7e14dcfSSatish Balay     }
68a7e14dcfSSatish Balay 
69a7e14dcfSSatish Balay   } else {
70a7e14dcfSSatish Balay 
71a7e14dcfSSatish Balay     n_same = 0; same=NULL;
72a7e14dcfSSatish Balay 
73a7e14dcfSSatish Balay   }
74a7e14dcfSSatish Balay 
75a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(ierr);
76a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,n_same,same,PETSC_COPY_VALUES,S);CHKERRQ(ierr);
77a7e14dcfSSatish Balay 
78a7e14dcfSSatish Balay   if (same) {
79a7e14dcfSSatish Balay     ierr = PetscFree(same); CHKERRQ(ierr);
80a7e14dcfSSatish Balay   }
81a7e14dcfSSatish Balay 
82a7e14dcfSSatish Balay   PetscFunctionReturn(0);
83a7e14dcfSSatish Balay }
84a7e14dcfSSatish Balay 
85a7e14dcfSSatish Balay #undef __FUNCT__
86a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichLessThan"
87a7e14dcfSSatish Balay /*@
88a7e14dcfSSatish Balay   VecWhichLessThan - Creates an index set containing the indices
89a7e14dcfSSatish Balay   where the vectors Vec1 < Vec2
90a7e14dcfSSatish Balay 
91a7e14dcfSSatish Balay   Collective on S
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay   Input Parameters:
94a7e14dcfSSatish Balay . Vec1, Vec2 - the two vectors to compare
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay   OutputParameter:
97a7e14dcfSSatish Balay . S - The index set containing the indices i where vec1[i] < vec2[i]
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay   Level: advanced
100a7e14dcfSSatish Balay @*/
101a7e14dcfSSatish Balay PetscErrorCode VecWhichLessThan(Vec Vec1, Vec Vec2, IS * S)
102a7e14dcfSSatish Balay {
103a7e14dcfSSatish Balay   int ierr;
104a7e14dcfSSatish Balay   PetscInt i;
105a7e14dcfSSatish Balay   PetscInt n,low,high,low2,high2,n_lt=0;
106a7e14dcfSSatish Balay   PetscInt *lt;
107a7e14dcfSSatish Balay   PetscReal *v1,*v2;
108a7e14dcfSSatish Balay   MPI_Comm comm;
109a7e14dcfSSatish Balay 
110a7e14dcfSSatish Balay   PetscFunctionBegin;
111a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec1,VEC_CLASSID,1);
112a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec2,VEC_CLASSID,2);
113a7e14dcfSSatish Balay   PetscCheckSameComm(Vec1,1,Vec2,2);
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(ierr);
116a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(ierr);
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay   if ( low != low2 || high != high2 )
119a7e14dcfSSatish Balay     SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
120a7e14dcfSSatish Balay 
121a7e14dcfSSatish Balay   ierr = VecGetLocalSize(Vec1,&n); CHKERRQ(ierr);
122a7e14dcfSSatish Balay 
123a7e14dcfSSatish Balay   if (n>0){
124a7e14dcfSSatish Balay 
125a7e14dcfSSatish Balay     if (Vec1 == Vec2){
126a7e14dcfSSatish Balay       ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr);
127a7e14dcfSSatish Balay       v2=v1;
128a7e14dcfSSatish Balay     } else {
129a7e14dcfSSatish Balay       ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr);
130a7e14dcfSSatish Balay       ierr = VecGetArray(Vec2,&v2); CHKERRQ(ierr);
131a7e14dcfSSatish Balay     }
132a7e14dcfSSatish Balay     ierr = PetscMalloc( n*sizeof(PetscInt),&lt ); CHKERRQ(ierr);
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay     for (i=0; i<n; i++){
135a7e14dcfSSatish Balay       if (v1[i] < v2[i]) {lt[n_lt]=low+i; n_lt++;}
136a7e14dcfSSatish Balay     }
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay     if (Vec1 == Vec2){
139a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr);
140a7e14dcfSSatish Balay     } else {
141a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr);
142a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec2,&v2); CHKERRQ(ierr);
143a7e14dcfSSatish Balay     }
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   } else {
146a7e14dcfSSatish Balay     n_lt=0; lt=NULL;
147a7e14dcfSSatish Balay   }
148a7e14dcfSSatish Balay 
149a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(ierr);
150a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,n_lt,lt,PETSC_COPY_VALUES,S);CHKERRQ(ierr);
151a7e14dcfSSatish Balay 
152a7e14dcfSSatish Balay   if (lt) {
153a7e14dcfSSatish Balay     ierr = PetscFree(lt); CHKERRQ(ierr);
154a7e14dcfSSatish Balay   }
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay   PetscFunctionReturn(0);
157a7e14dcfSSatish Balay }
158a7e14dcfSSatish Balay 
159a7e14dcfSSatish Balay #undef __FUNCT__
160a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichGreaterThan"
161a7e14dcfSSatish Balay /*@
162a7e14dcfSSatish Balay   VecWhichGreaterThan - Creates an index set containing the indices
163a7e14dcfSSatish Balay   where the vectors Vec1 > Vec2
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay   Collective on S
166a7e14dcfSSatish Balay 
167a7e14dcfSSatish Balay   Input Parameters:
168a7e14dcfSSatish Balay . Vec1, Vec2 - the two vectors to compare
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay   OutputParameter:
171a7e14dcfSSatish Balay . S - The index set containing the indices i where vec1[i] > vec2[i]
172a7e14dcfSSatish Balay 
173a7e14dcfSSatish Balay   Level: advanced
174a7e14dcfSSatish Balay @*/
175a7e14dcfSSatish Balay PetscErrorCode VecWhichGreaterThan(Vec Vec1, Vec Vec2, IS * S)
176a7e14dcfSSatish Balay {
177a7e14dcfSSatish Balay   int    ierr;
178a7e14dcfSSatish Balay   PetscInt n,low,high,low2,high2,n_gt=0,i;
179a7e14dcfSSatish Balay   PetscInt    *gt=NULL;
180a7e14dcfSSatish Balay   PetscReal *v1,*v2;
181a7e14dcfSSatish Balay   MPI_Comm comm;
182a7e14dcfSSatish Balay 
183a7e14dcfSSatish Balay   PetscFunctionBegin;
184a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec1,VEC_CLASSID,1);
185a7e14dcfSSatish Balay   PetscValidHeaderSpecific(Vec2,VEC_CLASSID,2);
186a7e14dcfSSatish Balay   PetscCheckSameComm(Vec1,1,Vec2,2);
187a7e14dcfSSatish Balay 
188a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(ierr);
189a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(ierr);
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay   if ( low != low2 || high != high2 )
192a7e14dcfSSatish Balay     SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay   ierr = VecGetLocalSize(Vec1,&n); CHKERRQ(ierr);
195a7e14dcfSSatish Balay 
196a7e14dcfSSatish Balay   if (n>0){
197a7e14dcfSSatish Balay 
198a7e14dcfSSatish Balay     if (Vec1 == Vec2){
199a7e14dcfSSatish Balay       ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr);
200a7e14dcfSSatish Balay       v2=v1;
201a7e14dcfSSatish Balay     } else {
202a7e14dcfSSatish Balay       ierr = VecGetArray(Vec1,&v1); CHKERRQ(ierr);
203a7e14dcfSSatish Balay       ierr = VecGetArray(Vec2,&v2); CHKERRQ(ierr);
204a7e14dcfSSatish Balay     }
205a7e14dcfSSatish Balay 
206a7e14dcfSSatish Balay     ierr = PetscMalloc( n*sizeof(PetscInt), &gt ); CHKERRQ(ierr);
207a7e14dcfSSatish Balay 
208a7e14dcfSSatish Balay     for (i=0; i<n; i++){
209a7e14dcfSSatish Balay       if (v1[i] > v2[i]) {gt[n_gt]=low+i; n_gt++;}
210a7e14dcfSSatish Balay     }
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay     if (Vec1 == Vec2){
213a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr);
214a7e14dcfSSatish Balay     } else {
215a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec1,&v1); CHKERRQ(ierr);
216a7e14dcfSSatish Balay       ierr = VecRestoreArray(Vec2,&v2); CHKERRQ(ierr);
217a7e14dcfSSatish Balay     }
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay   } else{
220a7e14dcfSSatish Balay 
221a7e14dcfSSatish Balay     n_gt=0; gt=NULL;
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay   }
224a7e14dcfSSatish Balay 
225a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(ierr);
226a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,n_gt,gt,PETSC_COPY_VALUES,S);CHKERRQ(ierr);
227a7e14dcfSSatish Balay 
228a7e14dcfSSatish Balay   if (gt) {
229a7e14dcfSSatish Balay     ierr = PetscFree(gt); CHKERRQ(ierr);
230a7e14dcfSSatish Balay   }
231a7e14dcfSSatish Balay   PetscFunctionReturn(0);
232a7e14dcfSSatish Balay }
233a7e14dcfSSatish Balay 
234a7e14dcfSSatish Balay #undef __FUNCT__
235a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichBetween"
236a7e14dcfSSatish Balay /*@
237a7e14dcfSSatish Balay   VecWhichBetween - Creates an index set containing the indices
238a7e14dcfSSatish Balay   where  VecLow < V < VecHigh
239a7e14dcfSSatish Balay 
240a7e14dcfSSatish Balay   Collective on S
241a7e14dcfSSatish Balay 
242a7e14dcfSSatish Balay   Input Parameters:
243a7e14dcfSSatish Balay + VecLow - lower bound
244a7e14dcfSSatish Balay . V - Vector to compare
245a7e14dcfSSatish Balay - VecHigh - higher bound
246a7e14dcfSSatish Balay 
247a7e14dcfSSatish Balay   OutputParameter:
248a7e14dcfSSatish Balay . S - The index set containing the indices i where veclow[i] < v[i] < vechigh[i]
249a7e14dcfSSatish Balay 
250a7e14dcfSSatish Balay   Level: advanced
251a7e14dcfSSatish Balay @*/
252a7e14dcfSSatish Balay PetscErrorCode VecWhichBetween(Vec VecLow, Vec V, Vec VecHigh, IS *S)
253a7e14dcfSSatish Balay {
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay   PetscErrorCode ierr;
256a7e14dcfSSatish Balay   PetscInt n,low,high,low2,high2,low3,high3,n_vm=0;
257a7e14dcfSSatish Balay   PetscInt *vm,i;
258a7e14dcfSSatish Balay   PetscReal *v1,*v2,*vmiddle;
259a7e14dcfSSatish Balay   MPI_Comm comm;
260a7e14dcfSSatish Balay 
261a7e14dcfSSatish Balay   PetscValidHeaderSpecific(V,VEC_CLASSID,2);
262a7e14dcfSSatish Balay   PetscCheckSameComm(V,2,VecLow,1); PetscCheckSameComm(V,2,VecHigh,3);
263a7e14dcfSSatish Balay 
264a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(VecLow, &low, &high); CHKERRQ(ierr);
265a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(VecHigh, &low2, &high2); CHKERRQ(ierr);
266a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(V, &low3, &high3); CHKERRQ(ierr);
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay   if ( low!=low2 || high!=high2 || low!=low3 || high!=high3 )
269a7e14dcfSSatish Balay     SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
270a7e14dcfSSatish Balay 
271a7e14dcfSSatish Balay   ierr = VecGetLocalSize(VecLow,&n); CHKERRQ(ierr);
272a7e14dcfSSatish Balay 
273a7e14dcfSSatish Balay   if (n>0){
274a7e14dcfSSatish Balay 
275a7e14dcfSSatish Balay     ierr = VecGetArray(VecLow,&v1); CHKERRQ(ierr);
276a7e14dcfSSatish Balay     if (VecLow != VecHigh){
277a7e14dcfSSatish Balay       ierr = VecGetArray(VecHigh,&v2); CHKERRQ(ierr);
278a7e14dcfSSatish Balay     } else {
279a7e14dcfSSatish Balay       v2=v1;
280a7e14dcfSSatish Balay     }
281a7e14dcfSSatish Balay     if ( V != VecLow && V != VecHigh){
282a7e14dcfSSatish Balay       ierr = VecGetArray(V,&vmiddle); CHKERRQ(ierr);
283a7e14dcfSSatish Balay     } else if ( V==VecLow ){
284a7e14dcfSSatish Balay       vmiddle=v1;
285a7e14dcfSSatish Balay     } else {
286a7e14dcfSSatish Balay       vmiddle =v2;
287a7e14dcfSSatish Balay     }
288a7e14dcfSSatish Balay 
289a7e14dcfSSatish Balay     ierr = PetscMalloc( n*sizeof(PetscInt), &vm ); CHKERRQ(ierr);
290a7e14dcfSSatish Balay 
291a7e14dcfSSatish Balay     for (i=0; i<n; i++){
292a7e14dcfSSatish Balay       if (v1[i] < vmiddle[i] && vmiddle[i] < v2[i]) {vm[n_vm]=low+i; n_vm++;}
293a7e14dcfSSatish Balay     }
294a7e14dcfSSatish Balay 
295a7e14dcfSSatish Balay     ierr = VecRestoreArray(VecLow,&v1); CHKERRQ(ierr);
296a7e14dcfSSatish Balay     if (VecLow != VecHigh){
297a7e14dcfSSatish Balay       ierr = VecRestoreArray(VecHigh,&v2); CHKERRQ(ierr);
298a7e14dcfSSatish Balay     }
299a7e14dcfSSatish Balay     if ( V != VecLow && V != VecHigh){
300a7e14dcfSSatish Balay       ierr = VecRestoreArray(V,&vmiddle); CHKERRQ(ierr);
301a7e14dcfSSatish Balay     }
302a7e14dcfSSatish Balay 
303a7e14dcfSSatish Balay   } else {
304a7e14dcfSSatish Balay 
305a7e14dcfSSatish Balay     n_vm=0; vm=NULL;
306a7e14dcfSSatish Balay 
307a7e14dcfSSatish Balay   }
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)V,&comm);CHKERRQ(ierr);
310a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,n_vm,vm,PETSC_COPY_VALUES,S);CHKERRQ(ierr);
311a7e14dcfSSatish Balay 
312a7e14dcfSSatish Balay   if (vm) {
313a7e14dcfSSatish Balay     ierr = PetscFree(vm); CHKERRQ(ierr);
314a7e14dcfSSatish Balay   }
315a7e14dcfSSatish Balay 
316a7e14dcfSSatish Balay   PetscFunctionReturn(0);
317a7e14dcfSSatish Balay }
318a7e14dcfSSatish Balay 
319a7e14dcfSSatish Balay 
320a7e14dcfSSatish Balay #undef __FUNCT__
321a7e14dcfSSatish Balay #define __FUNCT__ "VecWhichBetweenOrEqual"
322a7e14dcfSSatish Balay /*@
323a7e14dcfSSatish Balay   VecWhichBetweenOrEqual - Creates an index set containing the indices
324a7e14dcfSSatish Balay   where  VecLow <= V <= VecHigh
325a7e14dcfSSatish Balay 
326a7e14dcfSSatish Balay   Collective on S
327a7e14dcfSSatish Balay 
328a7e14dcfSSatish Balay   Input Parameters:
329a7e14dcfSSatish Balay + VecLow - lower bound
330a7e14dcfSSatish Balay . V - Vector to compare
331a7e14dcfSSatish Balay - VecHigh - higher bound
332a7e14dcfSSatish Balay 
333a7e14dcfSSatish Balay   OutputParameter:
334a7e14dcfSSatish Balay . S - The index set containing the indices i where veclow[i] <= v[i] <= vechigh[i]
335a7e14dcfSSatish Balay 
336a7e14dcfSSatish Balay   Level: advanced
337a7e14dcfSSatish Balay @*/
338a7e14dcfSSatish Balay 
339a7e14dcfSSatish Balay PetscErrorCode VecWhichBetweenOrEqual(Vec VecLow, Vec V, Vec VecHigh, IS * S)
340a7e14dcfSSatish Balay {
341a7e14dcfSSatish Balay   PetscErrorCode ierr;
342a7e14dcfSSatish Balay   PetscInt n,low,high,low2,high2,low3,high3,n_vm=0,i;
343a7e14dcfSSatish Balay   PetscInt *vm;
344a7e14dcfSSatish Balay   PetscReal *v1,*v2,*vmiddle;
345a7e14dcfSSatish Balay   MPI_Comm comm;
346a7e14dcfSSatish Balay 
347a7e14dcfSSatish Balay   PetscValidHeaderSpecific(V,VEC_CLASSID,2);
348a7e14dcfSSatish Balay   PetscCheckSameComm(V,2,VecLow,1); PetscCheckSameComm(V,2,VecHigh,3);
349a7e14dcfSSatish Balay 
350a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(VecLow, &low, &high); CHKERRQ(ierr);
351a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(VecHigh, &low2, &high2); CHKERRQ(ierr);
352a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(V, &low3, &high3); CHKERRQ(ierr);
353a7e14dcfSSatish Balay 
354a7e14dcfSSatish Balay   if ( low!=low2 || high!=high2 || low!=low3 || high!=high3 )
355a7e14dcfSSatish Balay     SETERRQ(PETSC_COMM_SELF,1,"Vectors must be identically loaded over processors");
356a7e14dcfSSatish Balay 
357a7e14dcfSSatish Balay   ierr = VecGetLocalSize(VecLow,&n); CHKERRQ(ierr);
358a7e14dcfSSatish Balay 
359a7e14dcfSSatish Balay   if (n>0){
360a7e14dcfSSatish Balay 
361a7e14dcfSSatish Balay     ierr = VecGetArray(VecLow,&v1); CHKERRQ(ierr);
362a7e14dcfSSatish Balay     if (VecLow != VecHigh){
363a7e14dcfSSatish Balay       ierr = VecGetArray(VecHigh,&v2); CHKERRQ(ierr);
364a7e14dcfSSatish Balay     } else {
365a7e14dcfSSatish Balay       v2=v1;
366a7e14dcfSSatish Balay     }
367a7e14dcfSSatish Balay     if ( V != VecLow && V != VecHigh){
368a7e14dcfSSatish Balay       ierr = VecGetArray(V,&vmiddle); CHKERRQ(ierr);
369a7e14dcfSSatish Balay     } else if ( V==VecLow ){
370a7e14dcfSSatish Balay       vmiddle=v1;
371a7e14dcfSSatish Balay     } else {
372a7e14dcfSSatish Balay       vmiddle =v2;
373a7e14dcfSSatish Balay     }
374a7e14dcfSSatish Balay 
375a7e14dcfSSatish Balay     ierr = PetscMalloc( n*sizeof(PetscInt), &vm ); CHKERRQ(ierr);
376a7e14dcfSSatish Balay 
377a7e14dcfSSatish Balay     for (i=0; i<n; i++){
378a7e14dcfSSatish Balay       if (v1[i] <= vmiddle[i] && vmiddle[i] <= v2[i]) {vm[n_vm]=low+i; n_vm++;}
379a7e14dcfSSatish Balay     }
380a7e14dcfSSatish Balay 
381a7e14dcfSSatish Balay     ierr = VecRestoreArray(VecLow,&v1); CHKERRQ(ierr);
382a7e14dcfSSatish Balay     if (VecLow != VecHigh){
383a7e14dcfSSatish Balay       ierr = VecRestoreArray(VecHigh,&v2); CHKERRQ(ierr);
384a7e14dcfSSatish Balay     }
385a7e14dcfSSatish Balay     if ( V != VecLow && V != VecHigh){
386a7e14dcfSSatish Balay       ierr = VecRestoreArray(V,&vmiddle); CHKERRQ(ierr);
387a7e14dcfSSatish Balay     }
388a7e14dcfSSatish Balay 
389a7e14dcfSSatish Balay   } else {
390a7e14dcfSSatish Balay 
391a7e14dcfSSatish Balay     n_vm=0; vm=NULL;
392a7e14dcfSSatish Balay 
393a7e14dcfSSatish Balay   }
394a7e14dcfSSatish Balay 
395a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)V,&comm);CHKERRQ(ierr);
396a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,n_vm,vm,PETSC_COPY_VALUES,S);CHKERRQ(ierr);
397a7e14dcfSSatish Balay 
398a7e14dcfSSatish Balay   if (vm) {
399a7e14dcfSSatish Balay     ierr = PetscFree(vm); CHKERRQ(ierr);
400a7e14dcfSSatish Balay   }
401a7e14dcfSSatish Balay 
402a7e14dcfSSatish Balay   PetscFunctionReturn(0);
403a7e14dcfSSatish Balay }
404a7e14dcfSSatish Balay 
405a7e14dcfSSatish Balay 
406a7e14dcfSSatish Balay #undef __FUNCT__
407a7e14dcfSSatish Balay #define __FUNCT__ "VecGetSubVec"
408a7e14dcfSSatish Balay /*@
409a7e14dcfSSatish Balay   VecGetSubVec - Gets a subvector using the IS
410a7e14dcfSSatish Balay 
411a7e14dcfSSatish Balay   Input Parameters:
412a7e14dcfSSatish Balay + vfull - the full matrix
413a7e14dcfSSatish Balay . is - the index set for the subvector
414a7e14dcfSSatish Balay . reduced_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK,  TAO_SUBSET_MATRIXFREE)
415a7e14dcfSSatish Balay - maskvalue - the value to set the unused vector elements to (for TAO_SUBSET_MASK or TAO_SUBSET_MATRIXFREE)
416a7e14dcfSSatish Balay 
417a7e14dcfSSatish Balay 
418a7e14dcfSSatish Balay   Output Parameters:
419a7e14dcfSSatish Balay . vreduced - the subvector
420a7e14dcfSSatish Balay 
421a7e14dcfSSatish Balay   Note:
422a7e14dcfSSatish Balay   maskvalue should usually be 0.0, unless a pointwise divide will be used.
423a7e14dcfSSatish Balay @*/
424a7e14dcfSSatish Balay PetscErrorCode VecGetSubVec(Vec vfull, IS is, PetscInt reduced_type, PetscReal maskvalue, Vec *vreduced)
425a7e14dcfSSatish Balay {
426a7e14dcfSSatish Balay     PetscErrorCode ierr;
427a7e14dcfSSatish Balay     PetscInt nfull,nreduced,nreduced_local,rlow,rhigh,flow,fhigh;
428a7e14dcfSSatish Balay     PetscInt i,nlocal;
429a7e14dcfSSatish Balay     PetscReal *fv,*rv;
430a7e14dcfSSatish Balay     const PetscInt *s;
431a7e14dcfSSatish Balay     IS ident;
432a7e14dcfSSatish Balay     VecType vtype;
433a7e14dcfSSatish Balay     VecScatter scatter;
434a7e14dcfSSatish Balay     MPI_Comm comm;
435a7e14dcfSSatish Balay 
436a7e14dcfSSatish Balay     PetscFunctionBegin;
437a7e14dcfSSatish Balay     PetscValidHeaderSpecific(vfull,VEC_CLASSID,1);
438a7e14dcfSSatish Balay     PetscValidHeaderSpecific(is,IS_CLASSID,2);
439a7e14dcfSSatish Balay 
440a7e14dcfSSatish Balay     ierr = VecGetSize(vfull, &nfull); CHKERRQ(ierr);
441a7e14dcfSSatish Balay     ierr = ISGetSize(is, &nreduced); CHKERRQ(ierr);
442a7e14dcfSSatish Balay 
443a7e14dcfSSatish Balay     if (nreduced == nfull) {
444a7e14dcfSSatish Balay 
445a7e14dcfSSatish Balay       ierr = VecDestroy(vreduced); CHKERRQ(ierr);
446a7e14dcfSSatish Balay       ierr = VecDuplicate(vfull,vreduced); CHKERRQ(ierr);
447a7e14dcfSSatish Balay       ierr = VecCopy(vfull,*vreduced); CHKERRQ(ierr);
448a7e14dcfSSatish Balay 
449a7e14dcfSSatish Balay     } else {
450a7e14dcfSSatish Balay 
451a7e14dcfSSatish Balay 	switch (reduced_type) {
452a7e14dcfSSatish Balay 	case TAO_SUBSET_SUBVEC:
453a7e14dcfSSatish Balay 	  ierr = VecGetType(vfull,&vtype); CHKERRQ(ierr);
454a7e14dcfSSatish Balay 	  ierr = VecGetOwnershipRange(vfull,&flow,&fhigh); CHKERRQ(ierr);
455a7e14dcfSSatish Balay 	  ierr = ISGetLocalSize(is,&nreduced_local); CHKERRQ(ierr);
456a7e14dcfSSatish Balay 	  ierr = PetscObjectGetComm((PetscObject)vfull,&comm); CHKERRQ(ierr);
457a7e14dcfSSatish Balay 	  if (*vreduced) {
458a7e14dcfSSatish Balay 	    ierr = VecDestroy(vreduced); CHKERRQ(ierr);
459a7e14dcfSSatish Balay 	  }
460a7e14dcfSSatish Balay 	  ierr = VecCreate(comm,vreduced); CHKERRQ(ierr);
461a7e14dcfSSatish Balay 	  ierr = VecSetType(*vreduced,vtype); CHKERRQ(ierr);
462a7e14dcfSSatish Balay 
463a7e14dcfSSatish Balay 
464a7e14dcfSSatish Balay 	  ierr = VecSetSizes(*vreduced,nreduced_local,nreduced); CHKERRQ(ierr);
465a7e14dcfSSatish Balay 
466a7e14dcfSSatish Balay 	  ierr = VecGetOwnershipRange(*vreduced,&rlow,&rhigh); CHKERRQ(ierr);
467a7e14dcfSSatish Balay 
468a7e14dcfSSatish Balay 	  ierr = ISCreateStride(comm,nreduced_local,rlow,1,&ident); CHKERRQ(ierr);
469a7e14dcfSSatish Balay 	  ierr = VecScatterCreate(vfull,is,*vreduced,ident,&scatter); CHKERRQ(ierr);
470a7e14dcfSSatish Balay 	  ierr = VecScatterBegin(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
471a7e14dcfSSatish Balay 	  ierr = VecScatterEnd(scatter,vfull,*vreduced,INSERT_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
472a7e14dcfSSatish Balay 	  ierr = VecScatterDestroy(&scatter); CHKERRQ(ierr);
473a7e14dcfSSatish Balay 	  ierr = ISDestroy(&ident); CHKERRQ(ierr);
474a7e14dcfSSatish Balay 	  break;
475a7e14dcfSSatish Balay 
476a7e14dcfSSatish Balay 	case TAO_SUBSET_MASK:
477a7e14dcfSSatish Balay 	case TAO_SUBSET_MATRIXFREE:
478a7e14dcfSSatish Balay 	  /* vr[i] = vf[i]   if i in is
479a7e14dcfSSatish Balay              vr[i] = 0       otherwise */
480a7e14dcfSSatish Balay 	  if (*vreduced == PETSC_NULL) {
481a7e14dcfSSatish Balay 	    ierr = VecDuplicate(vfull,vreduced); CHKERRQ(ierr);
482a7e14dcfSSatish Balay 	  }
483a7e14dcfSSatish Balay 
484a7e14dcfSSatish Balay 	  ierr = VecSet(*vreduced,maskvalue); CHKERRQ(ierr);
485a7e14dcfSSatish Balay 	  ierr = ISGetLocalSize(is,&nlocal); CHKERRQ(ierr);
486a7e14dcfSSatish Balay 	  ierr = VecGetOwnershipRange(vfull,&flow,&fhigh); CHKERRQ(ierr);
487a7e14dcfSSatish Balay 	  ierr = VecGetArray(vfull,&fv); CHKERRQ(ierr);
488a7e14dcfSSatish Balay 	  ierr = VecGetArray(*vreduced,&rv); CHKERRQ(ierr);
489a7e14dcfSSatish Balay 	  ierr = ISGetIndices(is,&s); CHKERRQ(ierr);
490a7e14dcfSSatish Balay 	  if (nlocal > (fhigh-flow)) {
491a7e14dcfSSatish Balay 	    SETERRQ2(PETSC_COMM_WORLD,1,"IS local size %d > Vec local size %d",nlocal,fhigh-flow);
492a7e14dcfSSatish Balay 	  }
493a7e14dcfSSatish Balay 	  for (i=0;i<nlocal;i++) {
494a7e14dcfSSatish Balay 	    rv[s[i]-flow] = fv[s[i]-flow];
495a7e14dcfSSatish Balay 	  }
496a7e14dcfSSatish Balay 	  ierr = ISRestoreIndices(is,&s); CHKERRQ(ierr);
497a7e14dcfSSatish Balay 	  ierr = VecRestoreArray(vfull,&fv); CHKERRQ(ierr);
498a7e14dcfSSatish Balay 	  ierr = VecRestoreArray(*vreduced,&rv); CHKERRQ(ierr);
499a7e14dcfSSatish Balay 	  break;
500a7e14dcfSSatish Balay 	}
501a7e14dcfSSatish Balay     }
502a7e14dcfSSatish Balay     PetscFunctionReturn(0);
503a7e14dcfSSatish Balay 
504a7e14dcfSSatish Balay }
505a7e14dcfSSatish Balay 
506a7e14dcfSSatish Balay #undef __FUNCT__
507a7e14dcfSSatish Balay #define __FUNCT__ "VecReducedXPY"
508a7e14dcfSSatish Balay /*@
509a7e14dcfSSatish Balay   VecReducedXPY - Adds a reduced vector to the appropriate elements of a full-space vector.
510a7e14dcfSSatish Balay 
511a7e14dcfSSatish Balay   Input Parameters:
512a7e14dcfSSatish Balay + vfull - the full-space vector
513a7e14dcfSSatish Balay . vreduced - the reduced-space vector
514a7e14dcfSSatish Balay - is - the index set for the reduced space
515a7e14dcfSSatish Balay 
516a7e14dcfSSatish Balay   Output Parameters:
517a7e14dcfSSatish Balay . vfull - the sum of the full-space vector and reduced-space vector
518a7e14dcfSSatish Balay @*/
519a7e14dcfSSatish Balay PetscErrorCode VecReducedXPY(Vec vfull, Vec vreduced, IS is)
520a7e14dcfSSatish Balay {
521a7e14dcfSSatish Balay     VecScatter scatter;
522a7e14dcfSSatish Balay     IS ident;
523a7e14dcfSSatish Balay     PetscInt nfull,nreduced,rlow,rhigh;
524a7e14dcfSSatish Balay     MPI_Comm comm;
525a7e14dcfSSatish Balay     PetscErrorCode ierr;
526a7e14dcfSSatish Balay 
527a7e14dcfSSatish Balay     PetscFunctionBegin;
528a7e14dcfSSatish Balay     PetscValidHeaderSpecific(vfull,VEC_CLASSID,1);
529a7e14dcfSSatish Balay     PetscValidHeaderSpecific(vreduced,VEC_CLASSID,2);
530a7e14dcfSSatish Balay     PetscValidHeaderSpecific(is,IS_CLASSID,3);
531a7e14dcfSSatish Balay     ierr = VecGetSize(vfull,&nfull); CHKERRQ(ierr);
532a7e14dcfSSatish Balay     ierr = VecGetSize(vreduced,&nreduced); CHKERRQ(ierr);
533a7e14dcfSSatish Balay 
534a7e14dcfSSatish Balay     if (nfull == nreduced) { /* Also takes care of masked vectors */
535a7e14dcfSSatish Balay       ierr = VecAXPY(vfull,1.0,vreduced); CHKERRQ(ierr);
536a7e14dcfSSatish Balay     } else {
537a7e14dcfSSatish Balay       ierr = PetscObjectGetComm((PetscObject)vfull,&comm); CHKERRQ(ierr);
538a7e14dcfSSatish Balay       ierr = VecGetOwnershipRange(vreduced,&rlow,&rhigh); CHKERRQ(ierr);
539a7e14dcfSSatish Balay       ierr = ISCreateStride(comm,rhigh-rlow,rlow,1,&ident); CHKERRQ(ierr);
540a7e14dcfSSatish Balay       ierr = VecScatterCreate(vreduced,ident,vfull,is,&scatter); CHKERRQ(ierr);
541a7e14dcfSSatish Balay       ierr = VecScatterBegin(scatter,vreduced,vfull,ADD_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
542a7e14dcfSSatish Balay       ierr = VecScatterEnd(scatter,vreduced,vfull,ADD_VALUES,SCATTER_FORWARD); CHKERRQ(ierr);
543a7e14dcfSSatish Balay       ierr = VecScatterDestroy(&scatter); CHKERRQ(ierr);
544a7e14dcfSSatish Balay       ierr = ISDestroy(&ident); CHKERRQ(ierr);
545a7e14dcfSSatish Balay     }
546a7e14dcfSSatish Balay 
547a7e14dcfSSatish Balay     PetscFunctionReturn(0);
548a7e14dcfSSatish Balay }
549a7e14dcfSSatish Balay 
550a7e14dcfSSatish Balay 
551a7e14dcfSSatish Balay #undef __FUNCT__
552a7e14dcfSSatish Balay #define __FUNCT__ "ISCreateComplement"
553a7e14dcfSSatish Balay /*@
554a7e14dcfSSatish Balay    ISCreateComplement - Creates the complement of the the index set
555a7e14dcfSSatish Balay 
556a7e14dcfSSatish Balay    Collective on IS
557a7e14dcfSSatish Balay 
558a7e14dcfSSatish Balay    Input Parameter:
559a7e14dcfSSatish Balay +  S -  a PETSc IS
560a7e14dcfSSatish Balay -  V - the reference vector space
561a7e14dcfSSatish Balay 
562a7e14dcfSSatish Balay    Output Parameter:
563a7e14dcfSSatish Balay .  T -  the complement of S
564a7e14dcfSSatish Balay 
565a7e14dcfSSatish Balay 
566a7e14dcfSSatish Balay .seealso ISCreateGeneral()
567a7e14dcfSSatish Balay 
568a7e14dcfSSatish Balay    Level: advanced
569a7e14dcfSSatish Balay @*/
570a7e14dcfSSatish Balay PetscErrorCode ISCreateComplement(IS S, Vec V, IS *T){
571a7e14dcfSSatish Balay   PetscErrorCode ierr;
572a7e14dcfSSatish Balay   PetscInt i,nis,nloc,high,low,n=0;
573a7e14dcfSSatish Balay   const PetscInt *s;
574a7e14dcfSSatish Balay   PetscInt *tt,*ss;
575a7e14dcfSSatish Balay   MPI_Comm comm;
576a7e14dcfSSatish Balay 
577a7e14dcfSSatish Balay   PetscFunctionBegin;
578a7e14dcfSSatish Balay   PetscValidHeaderSpecific(S,IS_CLASSID,1);
579a7e14dcfSSatish Balay   PetscValidHeaderSpecific(V,VEC_CLASSID,2);
580a7e14dcfSSatish Balay 
581a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(V,&low,&high); CHKERRQ(ierr);
582a7e14dcfSSatish Balay   ierr = VecGetLocalSize(V,&nloc); CHKERRQ(ierr);
583a7e14dcfSSatish Balay   ierr = ISGetLocalSize(S,&nis); CHKERRQ(ierr);
584a7e14dcfSSatish Balay   ierr = ISGetIndices(S, &s); CHKERRQ(ierr);
585a7e14dcfSSatish Balay   ierr = PetscMalloc( nloc*sizeof(PetscInt),&tt ); CHKERRQ(ierr);
586a7e14dcfSSatish Balay   ierr = PetscMalloc( nloc*sizeof(PetscInt),&ss ); CHKERRQ(ierr);
587a7e14dcfSSatish Balay 
588a7e14dcfSSatish Balay   for (i=low; i<high; i++){ tt[i-low]=i; }
589a7e14dcfSSatish Balay 
590a7e14dcfSSatish Balay   for (i=0; i<nis; i++){ tt[s[i]-low] = -2; }
591a7e14dcfSSatish Balay 
592a7e14dcfSSatish Balay   for (i=0; i<nloc; i++){
593a7e14dcfSSatish Balay     if (tt[i]>-1){ ss[n]=tt[i]; n++; }
594a7e14dcfSSatish Balay   }
595a7e14dcfSSatish Balay 
596a7e14dcfSSatish Balay   ierr = ISRestoreIndices(S, &s); CHKERRQ(ierr);
597a7e14dcfSSatish Balay 
598a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)S,&comm);CHKERRQ(ierr);
599a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,n,ss,PETSC_COPY_VALUES,T);CHKERRQ(ierr);
600a7e14dcfSSatish Balay 
601a7e14dcfSSatish Balay   if (tt) {
602a7e14dcfSSatish Balay     ierr = PetscFree(tt); CHKERRQ(ierr);
603a7e14dcfSSatish Balay   }
604a7e14dcfSSatish Balay   if (ss) {
605a7e14dcfSSatish Balay     ierr = PetscFree(ss); CHKERRQ(ierr);
606a7e14dcfSSatish Balay   }
607a7e14dcfSSatish Balay 
608a7e14dcfSSatish Balay   PetscFunctionReturn(0);
609a7e14dcfSSatish Balay }
610a7e14dcfSSatish Balay 
611a7e14dcfSSatish Balay #undef __FUNCT__
612a7e14dcfSSatish Balay #define __FUNCT__ "VecISSetToConstant"
613a7e14dcfSSatish Balay /*@
614a7e14dcfSSatish Balay    VecISSetToConstant - Sets the elements of a vector, specified by an index set, to a constant
615a7e14dcfSSatish Balay 
616a7e14dcfSSatish Balay    Input Parameter:
617a7e14dcfSSatish Balay +  S -  a PETSc IS
618a7e14dcfSSatish Balay .  c - the constant
619a7e14dcfSSatish Balay -  V - a Vec
620a7e14dcfSSatish Balay 
621a7e14dcfSSatish Balay .seealso VecSet()
622a7e14dcfSSatish Balay 
623a7e14dcfSSatish Balay    Level: advanced
624a7e14dcfSSatish Balay @*/
625a7e14dcfSSatish Balay PetscErrorCode VecISSetToConstant(IS S, PetscReal c, Vec V){
626a7e14dcfSSatish Balay   PetscErrorCode ierr;
627a7e14dcfSSatish Balay   PetscInt nloc,low,high,i;
628a7e14dcfSSatish Balay   const PetscInt *s;
629a7e14dcfSSatish Balay   PetscReal *v;
630a7e14dcfSSatish Balay   PetscFunctionBegin;
631a7e14dcfSSatish Balay   PetscValidHeaderSpecific(V,VEC_CLASSID,3);
632a7e14dcfSSatish Balay   PetscValidHeaderSpecific(S,IS_CLASSID,1);
633a7e14dcfSSatish Balay   PetscValidType(V,3);
634a7e14dcfSSatish Balay   PetscCheckSameComm(V,3,S,1);
635a7e14dcfSSatish Balay 
636a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(V, &low, &high); CHKERRQ(ierr);
637a7e14dcfSSatish Balay   ierr = ISGetLocalSize(S,&nloc);CHKERRQ(ierr);
638a7e14dcfSSatish Balay 
639a7e14dcfSSatish Balay   ierr = ISGetIndices(S, &s); CHKERRQ(ierr);
640a7e14dcfSSatish Balay   ierr = VecGetArray(V,&v); CHKERRQ(ierr);
641a7e14dcfSSatish Balay   for (i=0; i<nloc; i++){
642a7e14dcfSSatish Balay     v[s[i]-low] = c;
643a7e14dcfSSatish Balay   }
644a7e14dcfSSatish Balay 
645a7e14dcfSSatish Balay   ierr = ISRestoreIndices(S, &s); CHKERRQ(ierr);
646a7e14dcfSSatish Balay   ierr = VecRestoreArray(V,&v); CHKERRQ(ierr);
647a7e14dcfSSatish Balay 
648a7e14dcfSSatish Balay   PetscFunctionReturn(0);
649a7e14dcfSSatish Balay 
650a7e14dcfSSatish Balay }
651a7e14dcfSSatish Balay 
652a7e14dcfSSatish Balay #undef __FUNCT__
653a7e14dcfSSatish Balay #define __FUNCT__ "MatGetSubMat"
654a7e14dcfSSatish Balay /*@
655a7e14dcfSSatish Balay   MatGetSubMat - Gets a submatrix using the IS
656a7e14dcfSSatish Balay 
657a7e14dcfSSatish Balay   Input Parameters:
658a7e14dcfSSatish Balay + M - the full matrix (n x n)
659a7e14dcfSSatish Balay . is - the index set for the submatrix (both row and column index sets need to be the same)
660a7e14dcfSSatish Balay . v1 - work vector of dimension n, needed for TAO_SUBSET_MASK option
661a7e14dcfSSatish Balay - subset_type - the method TAO is using for subsetting (TAO_SUBSET_SUBVEC, TAO_SUBSET_MASK,
662a7e14dcfSSatish Balay   TAO_SUBSET_MATRIXFREE)
663a7e14dcfSSatish Balay 
664a7e14dcfSSatish Balay   Output Parameters:
665a7e14dcfSSatish Balay . Msub - the submatrix
666a7e14dcfSSatish Balay @*/
667a7e14dcfSSatish Balay PetscErrorCode MatGetSubMat(Mat M, IS is, Vec v1, TaoSubsetType subset_type, Mat *Msub)
668a7e14dcfSSatish Balay {
669a7e14dcfSSatish Balay   PetscErrorCode ierr;
670a7e14dcfSSatish Balay   IS iscomp;
671a7e14dcfSSatish Balay   PetscBool flg;
672a7e14dcfSSatish Balay   PetscFunctionBegin;
673a7e14dcfSSatish Balay   PetscValidHeaderSpecific(M,MAT_CLASSID,1);
674a7e14dcfSSatish Balay   PetscValidHeaderSpecific(is,IS_CLASSID,2);
675a7e14dcfSSatish Balay   if (*Msub) {
676a7e14dcfSSatish Balay     ierr = MatDestroy(Msub); CHKERRQ(ierr);
677a7e14dcfSSatish Balay   }
678a7e14dcfSSatish Balay   switch (subset_type) {
679a7e14dcfSSatish Balay     case TAO_SUBSET_SUBVEC:
680a7e14dcfSSatish Balay       ierr = MatGetSubMatrix(M, is, is, MAT_INITIAL_MATRIX, Msub); CHKERRQ(ierr);
681a7e14dcfSSatish Balay       break;
682a7e14dcfSSatish Balay 
683a7e14dcfSSatish Balay     case TAO_SUBSET_MASK:
684a7e14dcfSSatish Balay       /* Get Reduced Hessian
685a7e14dcfSSatish Balay 	 Msub[i,j] = M[i,j] if i,j in Free_Local or i==j
686a7e14dcfSSatish Balay 	 Msub[i,j] = 0      if i!=j and i or j not in Free_Local
687a7e14dcfSSatish Balay       */
688a7e14dcfSSatish Balay       ierr = PetscOptionsBool("-different_submatrix","use separate hessian matrix when computing submatrices","TaoSubsetType",PETSC_FALSE,&flg,PETSC_NULL);
689a7e14dcfSSatish Balay       if (flg == PETSC_TRUE) {
690a7e14dcfSSatish Balay 	ierr = MatDuplicate(M, MAT_COPY_VALUES, Msub); CHKERRQ(ierr);
691a7e14dcfSSatish Balay       } else {
692a7e14dcfSSatish Balay 	/* Act on hessian directly (default) */
693a7e14dcfSSatish Balay 	ierr = PetscObjectReference((PetscObject)M); CHKERRQ(ierr);
694a7e14dcfSSatish Balay 	*Msub = M;
695a7e14dcfSSatish Balay       }
696a7e14dcfSSatish Balay       /* Save the diagonal to temporary vector */
697a7e14dcfSSatish Balay       ierr = MatGetDiagonal(*Msub,v1); CHKERRQ(ierr);
698a7e14dcfSSatish Balay 
699a7e14dcfSSatish Balay       /* Zero out rows and columns */
700a7e14dcfSSatish Balay       ierr = ISCreateComplement(is,v1,&iscomp); CHKERRQ(ierr);
701a7e14dcfSSatish Balay 
702a7e14dcfSSatish Balay       /* Use v1 instead of 0 here because of PETSc bug */
703a7e14dcfSSatish Balay       ierr = MatZeroRowsColumnsIS(*Msub,iscomp,1.0,v1,v1); CHKERRQ(ierr);
704a7e14dcfSSatish Balay 
705a7e14dcfSSatish Balay       ierr = ISDestroy(&iscomp); CHKERRQ(ierr);
706a7e14dcfSSatish Balay       break;
707a7e14dcfSSatish Balay     case TAO_SUBSET_MATRIXFREE:
708a7e14dcfSSatish Balay       ierr = ISCreateComplement(is,v1,&iscomp); CHKERRQ(ierr);
709a7e14dcfSSatish Balay       ierr = MatCreateSubMatrixFree(M,iscomp,iscomp,Msub); CHKERRQ(ierr);
710a7e14dcfSSatish Balay       ierr = ISDestroy(&iscomp); CHKERRQ(ierr);
711a7e14dcfSSatish Balay       break;
712a7e14dcfSSatish Balay   }
713a7e14dcfSSatish Balay   PetscFunctionReturn(0);
714a7e14dcfSSatish Balay }
715