Actual source code: taois_petsc.c
1: #include "tao_general.h" /*I "tao.h" I*/
2: #ifdef TAO_USE_PETSC
3: #include "taois_petsc.h"
4: #include "../vector/taovec_petsc.h"
6: int VecWhichBetween(Vec, Vec, Vec, IS *);
7: int VecWhichBetweenOrEqual(Vec, Vec, Vec, IS *);
8: int VecWhichGreaterThan(Vec, Vec, IS * );
9: int VecWhichLessThan(Vec, Vec, IS *);
10: int VecWhichEqual(Vec, Vec, IS *);
12: static int petscisevents=0;
13: static int tredistribute=0;
14: //static TaoPetscISType defaultistype=TaoRedistributeSubset;
15: //static TaoPetscISType defaultistype=TaoNoRedistributeSubset;
16: static TaoPetscISType defaultistype=TaoMaskFullSpace;
20: /* @C
21: TaoWrapPetscIS - Create a new TaoIndexSet object using PETSc.
23: Input Parameter:
24: + S - a PETSc IS
25: - imax - the maximum local length of the index set (Should be equal to or greater than the local length of the vector)
27: Output Parameter:
28: . SS - address of a new TaoIndexSet
30: Note:
31: A TaoIndexSetPetsc is an object with the methods of an abstract
32: TaoIndexSet object. A TaoIndexSetPetsc contains an implementation of the
33: TaoIndexSet methods. Routines using these vectors should declare a
34: pointer to a TaoIndexSet, assign this pointer to the address of a
35: TaoIndexSet object,use the pointer to invoke methods on the object, and use
36: this pointer as an argument when calling other routines. This usage is
37: different from the usage of a PETSc IS. In PETSc, applications will
38: typically declare a IS, and pass it as an argument into routines. That is,
39: applications will typically declare a pointer to a TaoIndexSet and use the
40: pointer, or declare a IS and use it directly.
42: .seealso TaoIndexSetGetPetscIS(), TaoIndexSetDestroy()
44: Level: developer
45: @ */
46: int TaoWrapPetscIS( IS S, int imax, TaoIndexSetPetsc ** SS){
49: // if (0==0) return 1;
50: // *SS = new TaoIndexSetPetsc(imax, S);
51: PetscFunctionReturn(1);
52: }
56: /* @C
57: TaoIndexSetPetsc - Create a new TaoIndexSet object using PETSc.
59: Input Parameter:
60: - V = A vector from the space that this indexset will describe.
61: + SS - an Index Set
63: Level: beginner
65: Note:
66: This index set will be deleted when the TaoIndexSet is deleted
68: Note:
69: The constuctor can also be called without the IndexSet
70:
71: Note:
72: The method TaoIndexSetPetsc::GetIS() returns a pointer to the current PETSc Index Set.
73: @ */
74: TaoIndexSetPetsc::TaoIndexSetPetsc(Vec V, IS SS):TaoIndexSet(){
75: Vec V2;
76: isp=0;isp2=0;ISGathered=0;ispComplement=0;scatter=0;
77: VecDuplicate(V,&V2);
78: this->SetUp(V2,SS);
79: return;
80: }
84: /* @C
85: TaoIndexSetPetsc - Create a new TaoIndexSet object using PETSc.
87: Input Parameter:
88: . V = A vector from the space that this indexset will describe.
90: Level: beginner
92: Note:
93: The method TaoIndexSetPetsc::GetIS() returns a pointer to the current PETSc Index Set.
94: @ */
95: TaoIndexSetPetsc::TaoIndexSetPetsc(Vec V):TaoIndexSet(){
96: int info,low,high;
97: IS is;
98: Vec V2;
99: isp=0;isp2=0;ISGathered=0;ispComplement=0;scatter=0;
100: VecDuplicate(V,&V2);
101: info=VecGetOwnershipRange(V,&low,&high);
102: info=ISCreateStride(PETSC_COMM_WORLD,1,low,1,&is);
103: this->SetUp(V2,is);
104: info = PetscObjectDereference((PetscObject)is);
105: return;
106: }
110: TaoIndexSetPetsc::~TaoIndexSetPetsc(){
111: clearit();
112: delete[] iptr;
113: if (this->VSpace){
114: VecDestroy(this->VSpace);
115: }
116: };
120: int TaoIndexSetPetsc::SetUp(Vec V2, IS SS){
122: int info,size,imax;
123: PetscTruth flg,ivalue;
124: MPI_Comm comm;
129: info=VecGetLocalSize(V2,&imax); CHKERRQ(info);
130: VSpace=V2;
132: iptr = new int[imax];
134: isp=0;isp2=0;ISGathered=0;ispComplement=0;scatter=0;
135: this->ispviewer=PETSC_VIEWER_STDOUT_WORLD;
136: reducedtype=defaultistype;
138: info=this->SetIS(SS); CHKERRQ(info);
139: nlocal=imax;
141: info = PetscObjectGetComm((PetscObject)V2,&comm);CHKERRQ(info);
142: MPI_Comm_size(comm,&size);
143: if (size==1 && defaultistype!=TaoMaskFullSpace &&defaultistype!=TaoMatrixFree ){
144: reducedtype=TaoSingleProcessor;
145: } else {
146: PetscOptionsGetTruth(0,"-redistribute",&ivalue,&flg);
147: if (flg && ivalue==PETSC_FALSE) reducedtype=TaoNoRedistributeSubset;
148: else if (flg && ivalue==PETSC_TRUE)reducedtype=TaoRedistributeSubset;
149: else reducedtype=defaultistype;
150: }
151: PetscOptionsHasName(0,"-taomask",&flg);
152: if (flg) { reducedtype= TaoMaskFullSpace; }
153: PetscOptionsHasName(0,"-taosubmatrixfree",&flg);
154: if (flg) { reducedtype= TaoMatrixFree; }
156: if (petscisevents==0){
157: PetscLogEventRegister(&petscisevents,"Identify Indices",IS_COOKIE);
158: }
159: if (tredistribute==0){
160: PetscLogEventRegister(&tredistribute,"IS Redistribute ",IS_COOKIE);
161: }
163: return(0);
164: }
168: int TaoIndexSetPetsc::SetIS(IS SS){
170: int info;
172: if (SS){
174: PetscObjectReference((PetscObject)SS);
175: }
176: info=this->clearit(); CHKERRQ(info);
177: this->isp=SS;
178: return(0);
179: }
183: /*@C
184: TaoSelectSubset - Select the manner in which subsets of variables in a matrix are selected
186: Input Parameter:
187: . type - the manner is which subsets are chosen.
189: Choice of types:
190: + TaoRedistributeSubset - When using multiple processors, redistribute the subset of variables
191: over the processors.
192: . TaoNoRedistributeSubset - When using multiple processors, keep variables on current processors
193: . TaoMaskFullSpace - Use the full set of variables, and mask the unneeded ones.
194: - TaoMatrixFree - Use the full set of variables, and mask the unneeded ones.
196: Note:
197: For use with active set methods in bound constrained optimization such as TRON and GPCG
199: Level: advanced
200: @*/
201: int TaoSelectSubset( TaoPetscISType type){
203: defaultistype=type;
204: return(0);
205: }
209: int TaoIndexSetPetsc::GetMask(Vec *VMask){
210: int info;
211: PetscScalar zero=0.0;
214: info=VecSet(this->VSpace, zero); CHKERRQ(info);
215: info=VecISSetToConstant(this->isp,1.0,this->VSpace); CHKERRQ(info);
216: *VMask=this->VSpace;
217: return(0);
218: }
222: int TaoIndexSetPetsc::UnionOf(TaoIndexSet* S1, TaoIndexSet* S2){
223: IS SS1=((TaoIndexSetPetsc*)S1)->GetIS();
224: IS SS2=((TaoIndexSetPetsc*)S2)->GetIS();
225: int info;
227: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
228: info=ISSum(SS1,SS2,&SS1); CHKERRQ(info);
229: info=PetscLogEventEnd(petscisevents,0,0,0,0); CHKERRQ(info);
230: return(0);
231: }
235: int TaoIndexSetPetsc::IntersectionOf(TaoIndexSet* S1, TaoIndexSet* S2){
236: IS SS1=((TaoIndexSetPetsc*)S1)->GetIS();
237: IS SS2=((TaoIndexSetPetsc*)S2)->GetIS();
238: IS SS3;
239: int info;
241: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
242: info=ISDifference(SS1,SS2,&SS3); CHKERRQ(info);
243: info=this->SetIS(SS3); CHKERRQ(info);
244: info=ISDestroy(SS3); CHKERRQ(info);
245: info=PetscLogEventEnd(petscisevents,0,0,0,0); CHKERRQ(info);
246: return(0);
247: }
251: int TaoIndexSetPetsc::ComplementOf(TaoIndexSet* SS){
252: IS SS1=((TaoIndexSetPetsc*)SS)->GetIS();
253: IS SS2;
254: int info;
256: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
257: info=ISCreateComplement(SS1,this->VSpace,&SS2); CHKERRQ(info);
258: info=this->SetIS(SS2); CHKERRQ(info);
259: info=ISDestroy(SS2); CHKERRQ(info);
260: info=PetscLogEventEnd(petscisevents,0,0,0,0); CHKERRQ(info);
261: return(0);
262: }
266: int TaoIndexSetPetsc::Duplicate(TaoIndexSet**SS){
267: int info;
268: TaoIndexSetPetsc *S;
269: IS is;
272: if (isp){
273: info = ISDuplicate(isp,&is); CHKERRQ(info);
274: S = new TaoIndexSetPetsc(this->VSpace, is);
275: info=ISDestroy(is); CHKERRQ(info);
276: } else {
277: S = new TaoIndexSetPetsc(this->VSpace);
278: }
279: *SS=S;
281: return(0);
282: }
286: int TaoIndexSetPetsc::IsSame(TaoIndexSet*SS, TaoTruth*flg){
287: int info;
288: PetscTruth pflg;
289: IS SS2=((TaoIndexSetPetsc*)SS)->GetIS();
291: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
292: info=ISEqual(isp,SS2,&pflg); CHKERRQ(info);
293: if (pflg==PETSC_TRUE) *flg=TAO_TRUE; else *flg=TAO_FALSE;
294: info=PetscLogEventEnd(petscisevents,0,0,0,0); CHKERRQ(info);
295: return(0);
296: }
300: int TaoIndexSetPetsc::View(){
301: int info;
303: info=ISView(this->isp,this->ispviewer);CHKERRQ(info); CHKERRQ(info);
304: return(0);
305: }
310: int TaoIndexSetPetsc::WhichEqual(TaoVec*v1,TaoVec*v2){
311: int info;
313: Vec V1=((TaoVecPetsc *)v1)->GetVec();
314: Vec V2=((TaoVecPetsc *)v2)->GetVec();
315: IS ISnew;
317: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
318: info=VecWhichEqual(V1,V2,&ISnew); CHKERRQ(info);
319: info=this->SetIS(ISnew); CHKERRQ(info);
320: info=ISDestroy(ISnew); CHKERRQ(info);
321: info=PetscLogEventEnd(petscisevents,0,0,0,0); CHKERRQ(info);
322: return(0);
323: }
327: int TaoIndexSetPetsc::WhichLessThan(TaoVec*v1,TaoVec*v2){
328: int info;
329: Vec V1=((TaoVecPetsc *)v1)->GetVec();
330: Vec V2=((TaoVecPetsc *)v2)->GetVec();
331: IS ISnew;
333: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
334: info=VecWhichLessThan(V1,V2,&ISnew); CHKERRQ(info);
335: info=this->SetIS(ISnew); CHKERRQ(info);
336: info=ISDestroy(ISnew); CHKERRQ(info);
337: info=PetscLogEventEnd(petscisevents,0,0,0,0); CHKERRQ(info);
338: return(0);
339: }
343: int TaoIndexSetPetsc::WhichGreaterThan(TaoVec*v1,TaoVec*v2){
344: int info;
345: Vec V1=((TaoVecPetsc *)v1)->GetVec();
346: Vec V2=((TaoVecPetsc *)v2)->GetVec();
347: IS ISnew;
349: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
350: info=VecWhichGreaterThan(V1,V2,&ISnew); CHKERRQ(info);
351: info=this->SetIS(ISnew); CHKERRQ(info);
352: info=ISDestroy(ISnew); CHKERRQ(info);
353: info=PetscLogEventEnd(petscisevents,0,0,0,0); CHKERRQ(info);
354: return(0);
355: }
359: int TaoIndexSetPetsc::WhichBetween(TaoVec*V1,TaoVec*V2,TaoVec*V3){
360: int info;
361: Vec VecLow=((TaoVecPetsc *)V1)->GetVec();
362: Vec VecHigh=((TaoVecPetsc *)V3)->GetVec();
363: Vec V=((TaoVecPetsc *)V2)->GetVec();
364: IS ISnew;
367: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
368: info=VecWhichBetween(VecLow,V,VecHigh,&ISnew); CHKERRQ(info);
369: info=this->SetIS(ISnew); CHKERRQ(info);
370: info=ISDestroy(ISnew); CHKERRQ(info);
371: info=PetscLogEventEnd(petscisevents,0,0,0,0);CHKERRQ(info);
373: return(0);
374: }
378: int TaoIndexSetPetsc::WhichBetweenOrEqual(TaoVec *V1, TaoVec *V2, TaoVec *V3)
379: {
380: int info;
381: Vec VecLow = ((TaoVecPetsc *)V1)->GetVec();
382: Vec VecHigh = ((TaoVecPetsc *)V3)->GetVec();
383: Vec V = ((TaoVecPetsc *)V2)->GetVec();
384: IS ISnew;
387: info=PetscLogEventBegin(petscisevents,0,0,0,0); CHKERRQ(info);
388: info=VecWhichBetweenOrEqual(VecLow,V,VecHigh,&ISnew); CHKERRQ(info);
389: info=this->SetIS(ISnew); CHKERRQ(info);
390: info=ISDestroy(ISnew); CHKERRQ(info);
391: info=PetscLogEventEnd(petscisevents,0,0,0,0);CHKERRQ(info);
392: return(0);
393: }
397: int TaoIndexSetPetsc::GetSize(int *nn){
399: int info=ISGetSize(isp,nn); CHKERRQ(info);
400: return(0);
401: }
405: int TaoIndexSetPetsc::clearit(){
406: int info;
408: if (scatter){
409: info=VecScatterDestroy(scatter); CHKERRQ(info);
410: scatter=0;
411: }
412: if (ISGathered){
413: info=ISDestroy(ISGathered); CHKERRQ(info);
414: ISGathered=0;
415: }
416: if (isp){
417: info=ISDestroy(isp); CHKERRQ(info);
418: isp=0;
419: }
420: if (isp2){
421: info=ISDestroy(isp2); CHKERRQ(info);
422: isp2=0;
423: }
424: if (ispComplement){
425: info=ISDestroy(ispComplement); CHKERRQ(info);
426: ispComplement=0;
427: }
428: isp=0;isp2=0;ISGathered=0;ispComplement=0;scatter=0;
430: return(0);
431: }
436: int TaoIndexSetPetsc::RedistributeIS(IS*ISNew){
437: int info,i,n,nn,low,high;
438: int *gl,*is;
439: Vec R;
440: IS ISAll;
441: MPI_Comm comm;
445: info=PetscLogEventBegin(petscisevents,0,0,0,0);CHKERRQ(info);
446: info = VecGetSize(VSpace,&n); CHKERRQ(info);
447: info = this->GetSize(&nn); CHKERRQ(info);
448: if (reducedtype==TaoMaskFullSpace){
449: *ISNew = isp;
450: } else if (reducedtype==TaoSingleProcessor){
451: *ISNew = isp;
452: } else if (reducedtype==TaoNoRedistributeSubset || n==nn){
453: *ISNew = isp;
454: } else if (reducedtype==TaoRedistributeSubset){
456: info=PetscLogEventBegin(tredistribute,0,0,0,0);CHKERRQ(info);
458: if (isp2==0){
459: info = this->GetWholeIS(&ISAll); CHKERRQ(info);
460:
461: info = PetscObjectGetComm((PetscObject)isp,&comm);CHKERRQ(info);
462: info = VecCreateMPI(comm,PETSC_DECIDE,nn,&R);CHKERRQ(info);
463: info = VecGetLocalSize(R,&n); CHKERRQ(info);
464: info = VecGetOwnershipRange(R, &low, &high); CHKERRQ(info);
465: info = ISGetIndices(ISAll, &is); CHKERRQ(info);
466: if (n>0){
467: info = PetscMalloc( n*sizeof(int),&gl ); CHKERRQ(info);
468: for (i=0; i<n; i++){
469: gl[i]= is[low+i];
470: }
471: } else gl=NULL;
472:
473: info = PetscObjectGetComm((PetscObject)isp,&comm);CHKERRQ(info);
474: info = ISCreateGeneral(comm,n,gl,ISNew); CHKERRQ(info);
475: info = ISRestoreIndices(ISGathered, &is); CHKERRQ(info);
476: if (n>0) {
477: info = PetscFree(gl); CHKERRQ(info);
478: }
480: isp2=*ISNew;
481: info = VecDestroy(R); CHKERRQ(info);
483: }
484: info=PetscLogEventEnd(tredistribute,0,0,0,0);CHKERRQ(info);
485: *ISNew = isp2;
486:
487: }
489: info=PetscLogEventEnd(petscisevents,0,0,0,0);CHKERRQ(info);
490: return(0);
491: }
496: int TaoIndexSetPetsc::GetWholeIS(IS*ISAll){
497: int info;
499: info=PetscLogEventBegin(petscisevents,0,0,0,0);CHKERRQ(info);
500: if (reducedtype==TaoSingleProcessor){
501: *ISAll = isp;
502: } else {
503: info=PetscLogEventBegin(tredistribute,0,0,0,0);CHKERRQ(info);
504: if (ISGathered==0){
505: info = ISAllGather(isp,&ISGathered); CHKERRQ(info);
506: info = ISSort(ISGathered); CHKERRQ(info);
507: } else if (ISGathered==0){
508: int *indices,*sizes,size,*offsets,n,*lindices,i,N,ierr;
509: MPI_Comm comm;
511: ierr=PetscObjectGetComm((PetscObject)isp,&comm);
512: MPI_Comm_size(comm,&size);
513: PetscMalloc(2*size*sizeof(int),&sizes);
514: offsets = sizes + size;
516: ierr=ISGetLocalSize(isp,&n);
517: MPI_Allgather(&n,1,MPI_INT,sizes,1,MPI_INT,comm);
518: offsets[0] = 0;
519: for (i=1;i<size; i++) offsets[i] = offsets[i-1] + sizes[i-1];
520: N = offsets[size-1] + sizes[size-1];
522: PetscMalloc((N+1)*sizeof(int),&indices);
523: ierr=ISGetIndices(isp,&lindices);
524: MPI_Allgatherv(lindices,n,MPI_INT,indices,sizes,offsets,MPI_INT,comm);
525: ierr=ISRestoreIndices(isp,&lindices);
527: ierr=ISCreateGeneral(comm,N,indices,&ISGathered);
528: PetscFree(indices);
529: PetscFree(sizes);
530: }
531: *ISAll = ISGathered;
532: info=PetscLogEventEnd(tredistribute,0,0,0,0);CHKERRQ(info);
533: }
535: info=PetscLogEventEnd(petscisevents,0,0,0,0);CHKERRQ(info);
536: return(0);
537: }
542: int TaoIndexSetPetsc::GetReducedVecScatter(Vec VFull, Vec VSub, VecScatter*scatterit){
543: int info,n,nlow,nhigh;
544: IS S1, Ident;
545: MPI_Comm comm;
549: if (scatter==0){
550: info = VecGetSize(VFull,&n); CHKERRQ(info);
551: info = VecGetSize(VSub,&nlow); CHKERRQ(info);
552: if (reducedtype==TaoMaskFullSpace || n==nlow){
554: info = VecScatterCreate(VFull,isp,VSub,isp,&scatter); CHKERRQ(info);
556: } else if (reducedtype==TaoSingleProcessor || reducedtype==TaoNoRedistributeSubset){
558: S1=isp;
559: info = VecGetOwnershipRange(VSub,&nlow,&nhigh); CHKERRQ(info);
560: n=nhigh-nlow;
562: info = PetscObjectGetComm((PetscObject)VFull,&comm);CHKERRQ(info);
563: info = ISCreateStride(comm,n,nlow,1,&Ident); CHKERRQ(info);
564: info = VecScatterCreate(VFull,S1,VSub,Ident,&scatter); CHKERRQ(info);
565: info = ISDestroy(Ident); CHKERRQ(info);
567: } else {
569: info = this->GetWholeIS(&S1); CHKERRQ(info);
570: nlow = 0;
571: info = VecGetSize(VSub,&n);CHKERRQ(info);
573: info = PetscObjectGetComm((PetscObject)VFull,&comm);CHKERRQ(info);
574: info = ISCreateStride(comm,n,nlow,1,&Ident); CHKERRQ(info);
575: info = VecScatterCreate(VFull,S1,VSub,Ident,&scatter); CHKERRQ(info);
576: info = ISDestroy(Ident); CHKERRQ(info);
578: }
579: }
580:
581: *scatterit=scatter;
583: return(0);
584: }
588: int TaoIndexSetPetsc::GetReducedType(TaoPetscISType* rtype){
590: *rtype=this->reducedtype;
591: return(0);
592: }
596: int TaoIndexSetPetsc::GetComplementIS(IS*isppp){
597: int info;
599: if (this->ispComplement==0){
600: info=PetscLogEventBegin(petscisevents,0,0,0,0);CHKERRQ(info);
601: info = ISCreateComplement(this->isp, this->VSpace, &this->ispComplement);CHKERRQ(info);
602: info=PetscLogEventEnd(petscisevents,0,0,0,0);CHKERRQ(info);
603: }
604: *isppp = this->ispComplement;
605: return(0);
606: }
610: int VecWhichEqual(Vec Vec1, Vec Vec2, IS * S)
611: {
612: /*
613: Create an index set containing the indices of
614: the vectors Vec1 and Vec2 with identical elements.
615: */
616: int i,n,low,high,low2,high2,n_same=0,info;
617: int *same;
618: PetscScalar *v1,*v2;
619: MPI_Comm comm;
627: info = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(info);
628: info = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(info);
630: if ( low != low2 || high != high2 )
631: SETERRQ(1,"Vectors must be identically loaded over processors");
633: info = VecGetLocalSize(Vec1,&n); CHKERRQ(info);
635: if (n>0){
636:
637: if (Vec1 == Vec2){
638: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
639: v2=v1;
640: } else {
641: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
642: info = VecGetArray(Vec2,&v2); CHKERRQ(info);
643: }
645: info = PetscMalloc( n*sizeof(int),&same ); CHKERRQ(info);
646:
647: for (i=0; i<n; i++){
648: if (v1[i] == v2[i]) {same[n_same]=low+i; n_same++;}
649: }
650:
651: if (Vec1 == Vec2){
652: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
653: } else {
654: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
655: info = VecRestoreArray(Vec2,&v2); CHKERRQ(info);
656: }
658: } else {
660: n_same = 0; same=NULL;
662: }
664: info = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(info);
665: info = ISCreateGeneral(comm,n_same,same,S);CHKERRQ(info);
667: if (same) {
668: info = PetscFree(same); CHKERRQ(info);
669: }
671: return(0);
672: }
676: int VecWhichLessThan(Vec Vec1, Vec Vec2, IS * S)
677: {
678: int i,n,low,high,low2,high2,n_lt=0,info;
679: int *lt;
680: PetscScalar *v1,*v2;
681: MPI_Comm comm;
688: info = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(info);
689: info = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(info);
691: if ( low != low2 || high != high2 )
692: SETERRQ(1,"Vectors must be identically loaded over processors");
694: info = VecGetLocalSize(Vec1,&n); CHKERRQ(info);
696: if (n>0){
698: if (Vec1 == Vec2){
699: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
700: v2=v1;
701: } else {
702: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
703: info = VecGetArray(Vec2,&v2); CHKERRQ(info);
704: }
705: info = PetscMalloc( n*sizeof(int),< ); CHKERRQ(info);
706:
707: for (i=0; i<n; i++){
708: if (v1[i] < v2[i]) {lt[n_lt]=high+i; n_lt++;}
709: }
711: if (Vec1 == Vec2){
712: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
713: } else {
714: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
715: info = VecRestoreArray(Vec2,&v2); CHKERRQ(info);
716: }
717:
718: } else {
719: n_lt=0; lt=NULL;
720: }
722: info = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(info);
723: info = ISCreateGeneral(comm,n_lt,lt,S);CHKERRQ(info);
725: if (lt) {
726: info = PetscFree(lt); CHKERRQ(info);
727: }
729: return(0);
730: }
734: int VecWhichGreaterThan(Vec Vec1, Vec Vec2, IS * S)
735: {
736: int i,n,low,high,low2,high2,n_gt=0,info;
737: int *gt=NULL;
738: PetscScalar *v1,*v2;
739: MPI_Comm comm;
746: info = VecGetOwnershipRange(Vec1, &low, &high); CHKERRQ(info);
747: info = VecGetOwnershipRange(Vec2, &low2, &high2); CHKERRQ(info);
749: if ( low != low2 || high != high2 )
750: SETERRQ(1,"Vectors must be identically loaded over processors");
752: info = VecGetLocalSize(Vec1,&n); CHKERRQ(info);
754: if (n>0){
756: if (Vec1 == Vec2){
757: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
758: v2=v1;
759: } else {
760: info = VecGetArray(Vec1,&v1); CHKERRQ(info);
761: info = VecGetArray(Vec2,&v2); CHKERRQ(info);
762: }
764: info = PetscMalloc( n*sizeof(int), > ); CHKERRQ(info);
765:
766: for (i=0; i<n; i++){
767: if (v1[i] > v2[i]) {gt[n_gt]=low+i; n_gt++;}
768: }
770: if (Vec1 == Vec2){
771: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
772: } else {
773: info = VecRestoreArray(Vec1,&v1); CHKERRQ(info);
774: info = VecRestoreArray(Vec2,&v2); CHKERRQ(info);
775: }
776:
777: } else{
778:
779: n_gt=0; gt=NULL;
781: }
783: info = PetscObjectGetComm((PetscObject)Vec1,&comm);CHKERRQ(info);
784: info = ISCreateGeneral(comm,n_gt,gt,S);CHKERRQ(info);
786: if (gt) {
787: info = PetscFree(gt); CHKERRQ(info);
788: }
789: return(0);
790: }
794: int VecWhichBetween(Vec VecLow, Vec V, Vec VecHigh, IS * S)
795: {
796: /*
797: Creates an index set with the indices of V whose
798: elements are stictly between the corresponding elements
799: of the vector VecLow and the Vector VecHigh
800: */
801: int i,n,low,high,low2,high2,low3,high3,n_vm=0,info;
802: int *vm;
803: PetscScalar *v1,*v2,*vmiddle;
804: MPI_Comm comm;
809: info = VecGetOwnershipRange(VecLow, &low, &high); CHKERRQ(info);
810: info = VecGetOwnershipRange(VecHigh, &low2, &high2); CHKERRQ(info);
811: info = VecGetOwnershipRange(V, &low3, &high3); CHKERRQ(info);
813: if ( low!=low2 || high!=high2 || low!=low3 || high!=high3 )
814: SETERRQ(1,"Vectors must be identically loaded over processors");
816: info = VecGetLocalSize(VecLow,&n); CHKERRQ(info);
818: if (n>0){
820: info = VecGetArray(VecLow,&v1); CHKERRQ(info);
821: if (VecLow != VecHigh){
822: info = VecGetArray(VecHigh,&v2); CHKERRQ(info);
823: } else {
824: v2=v1;
825: }
826: if ( V != VecLow && V != VecHigh){
827: info = VecGetArray(V,&vmiddle); CHKERRQ(info);
828: } else if ( V==VecLow ){
829: vmiddle=v1;
830: } else {
831: vmiddle =v2;
832: }
834: info = PetscMalloc( n*sizeof(int), &vm ); CHKERRQ(info);
835:
836: for (i=0; i<n; i++){
837: if (v1[i] < vmiddle[i] && vmiddle[i] < v2[i]) {vm[n_vm]=low+i; n_vm++;}
838: }
840: info = VecRestoreArray(VecLow,&v1); CHKERRQ(info);
841: if (VecLow != VecHigh){
842: info = VecRestoreArray(VecHigh,&v2); CHKERRQ(info);
843: }
844: if ( V != VecLow && V != VecHigh){
845: info = VecRestoreArray(V,&vmiddle); CHKERRQ(info);
846: }
848: } else {
850: n_vm=0; vm=NULL;
852: }
854: info = PetscObjectGetComm((PetscObject)V,&comm);CHKERRQ(info);
855: info = ISCreateGeneral(comm,n_vm,vm,S);CHKERRQ(info);
857: if (vm) {
858: info = PetscFree(vm); CHKERRQ(info);
859: }
861: return(0);
862: }
867: int VecWhichBetweenOrEqual(Vec VecLow, Vec V, Vec VecHigh, IS * S)
868: {
869: /*
870: Creates an index set with the indices of V whose
871: elements are stictly between the corresponding elements
872: of the vector VecLow and the Vector VecHigh
873: */
874: int i,n,low,high,low2,high2,low3,high3,n_vm=0,info;
875: int *vm;
876: PetscScalar *v1,*v2,*vmiddle;
877: MPI_Comm comm;
882: info = VecGetOwnershipRange(VecLow, &low, &high); CHKERRQ(info);
883: info = VecGetOwnershipRange(VecHigh, &low2, &high2); CHKERRQ(info);
884: info = VecGetOwnershipRange(V, &low3, &high3); CHKERRQ(info);
886: if ( low!=low2 || high!=high2 || low!=low3 || high!=high3 )
887: SETERRQ(1,"Vectors must be identically loaded over processors");
889: info = VecGetLocalSize(VecLow,&n); CHKERRQ(info);
891: if (n>0){
893: info = VecGetArray(VecLow,&v1); CHKERRQ(info);
894: if (VecLow != VecHigh){
895: info = VecGetArray(VecHigh,&v2); CHKERRQ(info);
896: } else {
897: v2=v1;
898: }
899: if ( V != VecLow && V != VecHigh){
900: info = VecGetArray(V,&vmiddle); CHKERRQ(info);
901: } else if ( V==VecLow ){
902: vmiddle=v1;
903: } else {
904: vmiddle =v2;
905: }
907: info = PetscMalloc( n*sizeof(int), &vm ); CHKERRQ(info);
908:
909: for (i=0; i<n; i++){
910: if (v1[i] <= vmiddle[i] && vmiddle[i] <= v2[i]) {vm[n_vm]=low+i; n_vm++;}
911: }
913: info = VecRestoreArray(VecLow,&v1); CHKERRQ(info);
914: if (VecLow != VecHigh){
915: info = VecRestoreArray(VecHigh,&v2); CHKERRQ(info);
916: }
917: if ( V != VecLow && V != VecHigh){
918: info = VecRestoreArray(V,&vmiddle); CHKERRQ(info);
919: }
921: } else {
923: n_vm=0; vm=NULL;
925: }
927: info = PetscObjectGetComm((PetscObject)V,&comm);CHKERRQ(info);
928: info = ISCreateGeneral(comm,n_vm,vm,S);CHKERRQ(info);
930: if (vm) {
931: info = PetscFree(vm); CHKERRQ(info);
932: }
934: return(0);
935: }
940: /*@C
941: ISCreateComplement - Creates the complement of the the index set
943: Input Parameter:
944: + S - a PETSc IS
945: - V - the reference vector space
947: Output Parameter:
948: . T - the complement of S
951: .seealso ISCreateGeneral()
953: Level: advanced
954: @*/
955: int ISCreateComplement(IS S, Vec V, IS *T){
956: int info,i,nis,nloc,high,low,n=0;
957: int *s,*tt,*ss;
958: MPI_Comm comm;
964: info = VecGetOwnershipRange(V,&low,&high); CHKERRQ(info);
965: info = VecGetLocalSize(V,&nloc); CHKERRQ(info);
966: info = ISGetLocalSize(S,&nis); CHKERRQ(info);
967: info = ISGetIndices(S, &s); CHKERRQ(info);
968: info = PetscMalloc( nloc*sizeof(int),&tt ); CHKERRQ(info);
969: info = PetscMalloc( nloc*sizeof(int),&ss ); CHKERRQ(info);
971: for (i=low; i<high; i++){ tt[i-low]=i; }
973: for (i=0; i<nis; i++){ tt[s[i]-low] = -2; }
974:
975: for (i=0; i<nloc; i++){
976: if (tt[i]>-1){ ss[n]=tt[i]; n++; }
977: }
979: info = ISRestoreIndices(S, &s); CHKERRQ(info);
980:
981: info = PetscObjectGetComm((PetscObject)S,&comm);CHKERRQ(info);
982: info = ISCreateGeneral(comm,n,ss,T);CHKERRQ(info);
983:
984: if (tt) {
985: info = PetscFree(tt); CHKERRQ(info);
986: }
987: if (ss) {
988: info = PetscFree(ss); CHKERRQ(info);
989: }
991: return(0);
992: }
993: #endif