Actual source code: vecnest.c
2: #include vecnestimpl.h
4: /* check all blocks are filled */
7: static PetscErrorCode VecAssemblyBegin_Nest(Vec v)
8: {
9: Vec_Nest *vs = (Vec_Nest*)v->data;
10: PetscInt i;
14: for (i=0;i<vs->nb;i++) {
15: if (!vs->v[i]) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Nest vector cannot contain NULL blocks");
16: VecAssemblyBegin(vs->v[i]);
17: }
18: return(0);
19: }
23: static PetscErrorCode VecAssemblyEnd_Nest(Vec v)
24: {
25: Vec_Nest *vs = (Vec_Nest*)v->data;
26: PetscInt i;
30: for (i=0;i<vs->nb;i++) {
31: VecAssemblyEnd(vs->v[i]);
32: }
33: return(0);
34: }
38: static PetscErrorCode VecDestroy_Nest(Vec v)
39: {
40: Vec_Nest *vs = (Vec_Nest*)v->data;
41: PetscInt i;
45: if (vs->v) {
46: for (i=0; i<vs->nb; i++) {
47: VecDestroy(&vs->v[i]);
48: }
49: PetscFree(vs->v);
50: }
51: for (i=0; i<vs->nb; i++) {
52: ISDestroy(&vs->is[i]);
53: }
54: PetscFree(vs->is);
56: PetscFree(vs);
57: return(0);
58: }
60: /* supports nested blocks */
63: static PetscErrorCode VecCopy_Nest(Vec x,Vec y)
64: {
65: Vec_Nest *bx = (Vec_Nest*)x->data;
66: Vec_Nest *by = (Vec_Nest*)y->data;
67: PetscInt i;
71: VecNestCheckCompatible2(x,1,y,2);
72: for (i=0; i<bx->nb; i++) {
73: VecCopy(bx->v[i],by->v[i]);
74: }
75: return(0);
76: }
78: /* supports nested blocks */
81: static PetscErrorCode VecDuplicate_Nest(Vec x,Vec *y)
82: {
83: Vec_Nest *bx = (Vec_Nest*)x->data;
84: Vec Y;
85: Vec *sub;
86: PetscInt i;
91: PetscMalloc(sizeof(Vec)*bx->nb,&sub);
92: for (i=0; i<bx->nb; i++ ) {
93: VecDuplicate(bx->v[i],&sub[i]);
94: }
95: VecCreateNest(((PetscObject)x)->comm,bx->nb,bx->is,sub,&Y);
96: for (i=0; i<bx->nb; i++ ) {
97: VecDestroy(&sub[i]);
98: }
99: PetscFree(sub);
100: *y = Y;
101: return(0);
102: }
104: /* supports nested blocks */
107: static PetscErrorCode VecDot_Nest(Vec x,Vec y,PetscScalar *val)
108: {
109: Vec_Nest *bx = (Vec_Nest*)x->data;
110: Vec_Nest *by = (Vec_Nest*)y->data;
111: PetscInt i,nr;
112: PetscScalar x_dot_y,_val;
116: nr = bx->nb;
117: _val = 0.0;
118: for (i=0; i<nr; i++) {
119: VecDot(bx->v[i],by->v[i],&x_dot_y);
120: _val = _val + x_dot_y;
121: }
122: *val = _val;
123: return(0);
124: }
126: /* supports nested blocks */
129: static PetscErrorCode VecTDot_Nest(Vec x,Vec y,PetscScalar *val)
130: {
131: Vec_Nest *bx = (Vec_Nest*)x->data;
132: Vec_Nest *by = (Vec_Nest*)y->data;
133: PetscInt i,nr;
134: PetscScalar x_dot_y,_val;
138: nr = bx->nb;
139: _val = 0.0;
140: for (i=0; i<nr; i++) {
141: VecTDot(bx->v[i],by->v[i],&x_dot_y);
142: _val = _val + x_dot_y;
143: }
144: *val = _val;
145: return(0);
146: }
150: static PetscErrorCode VecAXPY_Nest(Vec y,PetscScalar alpha,Vec x)
151: {
152: Vec_Nest *bx = (Vec_Nest*)x->data;
153: Vec_Nest *by = (Vec_Nest*)y->data;
154: PetscInt i,nr;
158: nr = bx->nb;
159: for (i=0; i<nr; i++) {
160: VecAXPY(by->v[i],alpha,bx->v[i]);
161: }
162: return(0);
163: }
167: static PetscErrorCode VecAYPX_Nest(Vec y,PetscScalar alpha,Vec x)
168: {
169: Vec_Nest *bx = (Vec_Nest*)x->data;
170: Vec_Nest *by = (Vec_Nest*)y->data;
171: PetscInt i,nr;
175: nr = bx->nb;
176: for (i=0; i<nr; i++) {
177: VecAYPX(by->v[i],alpha,bx->v[i]);
178: }
179: return(0);
180: }
184: static PetscErrorCode VecAXPBY_Nest(Vec y,PetscScalar alpha,PetscScalar beta,Vec x)
185: {
186: Vec_Nest *bx = (Vec_Nest*)x->data;
187: Vec_Nest *by = (Vec_Nest*)y->data;
188: PetscInt i,nr;
192: nr = bx->nb;
193: for (i=0; i<nr; i++) {
194: VecAXPBY(by->v[i],alpha,beta,bx->v[i]);
195: }
196: return(0);
197: }
201: static PetscErrorCode VecScale_Nest(Vec x,PetscScalar alpha)
202: {
203: Vec_Nest *bx = (Vec_Nest*)x->data;
204: PetscInt i,nr;
208: nr = bx->nb;
209: for (i=0; i<nr; i++) {
210: VecScale(bx->v[i],alpha);
211: }
212: return(0);
213: }
217: static PetscErrorCode VecPointwiseMult_Nest(Vec w,Vec x,Vec y)
218: {
219: Vec_Nest *bx = (Vec_Nest*)x->data;
220: Vec_Nest *by = (Vec_Nest*)y->data;
221: Vec_Nest *bw = (Vec_Nest*)w->data;
222: PetscInt i,nr;
226: VecNestCheckCompatible3(w,1,x,3,y,4);
227: nr = bx->nb;
228: for (i=0; i<nr; i++) {
229: VecPointwiseMult(bw->v[i],bx->v[i],by->v[i]);
230: }
231: return(0);
232: }
236: static PetscErrorCode VecPointwiseDivide_Nest(Vec w,Vec x,Vec y)
237: {
238: Vec_Nest *bx = (Vec_Nest*)x->data;
239: Vec_Nest *by = (Vec_Nest*)y->data;
240: Vec_Nest *bw = (Vec_Nest*)w->data;
241: PetscInt i,nr;
245: VecNestCheckCompatible3(w,1,x,2,y,3);
247: nr = bx->nb;
248: for (i=0; i<nr; i++) {
249: VecPointwiseDivide(bw->v[i],bx->v[i],by->v[i]);
250: }
251: return(0);
252: }
256: static PetscErrorCode VecReciprocal_Nest(Vec x)
257: {
258: Vec_Nest *bx = (Vec_Nest*)x->data;
259: PetscInt i,nr;
263: nr = bx->nb;
264: for (i=0; i<nr; i++) {
265: VecReciprocal(bx->v[i]);
266: }
268: return(0);
269: }
273: static PetscErrorCode VecNorm_Nest(Vec xin,NormType type,PetscReal* z)
274: {
275: Vec_Nest *bx = (Vec_Nest*)xin->data;
276: PetscInt i,nr;
277: PetscReal z_i;
278: PetscReal _z;
282: nr = bx->nb;
283: _z = 0.0;
285: if (type == NORM_2) {
286: PetscScalar dot;
287: VecDot(xin,xin,&dot);
288: _z = PetscAbsScalar(PetscSqrtScalar(dot));
289: } else if (type == NORM_1) {
290: for (i=0; i<nr; i++) {
291: VecNorm(bx->v[i],type,&z_i);
292: _z = _z + z_i;
293: }
294: } else if (type == NORM_INFINITY) {
295: for (i=0; i<nr; i++) {
296: VecNorm(bx->v[i],type,&z_i);
297: if (z_i > _z) {
298: _z = z_i;
299: }
300: }
301: }
303: *z = _z;
304: return(0);
305: }
309: static PetscErrorCode VecMAXPY_Nest(Vec y,PetscInt nv,const PetscScalar alpha[],Vec *x)
310: {
311: PetscInt v;
315: for (v=0; v<nv; v++) {
316: /* Do axpy on each vector,v */
317: VecAXPY(y,alpha[v],x[v]);
318: }
319: return(0);
320: }
324: static PetscErrorCode VecMDot_Nest(Vec x,PetscInt nv,const Vec y[],PetscScalar *val)
325: {
326: PetscInt j;
330: for (j=0; j<nv; j++) {
331: VecDot(x,y[j],&val[j]);
332: }
333: return(0);
334: }
338: static PetscErrorCode VecMTDot_Nest(Vec x,PetscInt nv,const Vec y[],PetscScalar *val)
339: {
340: PetscInt j;
344: for (j=0; j<nv; j++) {
345: VecTDot(x,y[j],&val[j]);
346: }
347: return(0);
348: }
352: static PetscErrorCode VecSet_Nest(Vec x,PetscScalar alpha)
353: {
354: Vec_Nest *bx = (Vec_Nest*)x->data;
355: PetscInt i,nr;
359: nr = bx->nb;
360: for (i=0; i<nr; i++) {
361: VecSet(bx->v[i],alpha);
362: }
363: return(0);
364: }
368: static PetscErrorCode VecConjugate_Nest(Vec x)
369: {
370: Vec_Nest *bx = (Vec_Nest*)x->data;
371: PetscInt j,nr;
375: nr = bx->nb;
376: for (j=0; j<nr; j++) {
377: VecConjugate(bx->v[j]);
378: }
379: return(0);
380: }
384: static PetscErrorCode VecSwap_Nest(Vec x,Vec y)
385: {
386: Vec_Nest *bx = (Vec_Nest*)x->data;
387: Vec_Nest *by = (Vec_Nest*)y->data;
388: PetscInt i,nr;
392: VecNestCheckCompatible2(x,1,y,2);
393: nr = bx->nb;
394: for (i=0; i<nr; i++) {
395: VecSwap(bx->v[i],by->v[i]);
396: }
397: return(0);
398: }
402: static PetscErrorCode VecWAXPY_Nest(Vec w,PetscScalar alpha,Vec x,Vec y)
403: {
404: Vec_Nest *bx = (Vec_Nest*)x->data;
405: Vec_Nest *by = (Vec_Nest*)y->data;
406: Vec_Nest *bw = (Vec_Nest*)w->data;
407: PetscInt i,nr;
411: VecNestCheckCompatible3(w,1,x,3,y,4);
413: nr = bx->nb;
414: for (i=0; i<nr; i++) {
415: VecWAXPY(bw->v[i],alpha,bx->v[i],by->v[i]);
416: }
417: return(0);
418: }
422: static PetscErrorCode VecMax_Nest_Recursive(Vec x,PetscInt *cnt,PetscInt *p,PetscReal *max)
423: {
424: Vec_Nest *bx = (Vec_Nest*)x->data;
426: PetscInt i,nr;
427: PetscBool isnest;
428: PetscInt L;
429: PetscInt _entry_loc;
430: PetscReal _entry_val;
434: PetscTypeCompare((PetscObject)x,VECNEST,&isnest);
435: if (!isnest) {
436: /* Not nest */
437: VecMax(x,&_entry_loc,&_entry_val);
438: if (_entry_val > *max) {
439: *max = _entry_val;
440: *p = _entry_loc + *cnt;
441: }
442: VecGetSize(x,&L);
443: *cnt = *cnt + L;
444: return(0);
445: }
447: /* Otherwise we have a nest */
448: bx = (Vec_Nest*)x->data;
449: nr = bx->nb;
451: /* now descend recursively */
452: for (i=0; i<nr; i++) {
453: VecMax_Nest_Recursive(bx->v[i],cnt,p,max);
454: }
455: return(0);
456: }
458: /* supports nested blocks */
461: static PetscErrorCode VecMax_Nest(Vec x,PetscInt *p,PetscReal *max)
462: {
463: PetscInt cnt;
467: cnt = 0;
468: *p = 0;
469: *max = 0.0;
470: VecMax_Nest_Recursive(x,&cnt,p,max);
471: return(0);
472: }
476: static PetscErrorCode VecMin_Nest_Recursive(Vec x,PetscInt *cnt,PetscInt *p,PetscReal *min)
477: {
478: Vec_Nest *bx = (Vec_Nest*)x->data;
479: PetscInt i,nr,L,_entry_loc;
480: PetscBool isnest;
481: PetscReal _entry_val;
485: PetscTypeCompare((PetscObject)x,VECNEST,&isnest);
486: if (!isnest) {
487: /* Not nest */
488: VecMin(x,&_entry_loc,&_entry_val);
489: if (_entry_val < *min) {
490: *min = _entry_val;
491: *p = _entry_loc + *cnt;
492: }
493: VecGetSize(x,&L);
494: *cnt = *cnt + L;
495: return(0);
496: }
498: /* Otherwise we have a nest */
499: bx = (Vec_Nest*)x->data;
500: nr = bx->nb;
502: /* now descend recursively */
503: for (i=0; i<nr; i++) {
504: VecMin_Nest_Recursive(bx->v[i],cnt,p,min);
505: }
506: return(0);
507: }
511: static PetscErrorCode VecMin_Nest(Vec x,PetscInt *p,PetscReal *min)
512: {
513: PetscInt cnt;
517: cnt = 0;
518: *p = 0;
519: *min = 1.0e308;
520: VecMin_Nest_Recursive(x,&cnt,p,min);
521: return(0);
522: }
524: /* supports nested blocks */
527: static PetscErrorCode VecView_Nest(Vec x,PetscViewer viewer)
528: {
529: Vec_Nest *bx = (Vec_Nest*)x->data;
530: PetscBool isascii;
531: PetscInt i;
535: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
536: if (isascii) {
537: PetscViewerASCIIPrintf(viewer,"Vector Object:\n");
538: PetscViewerASCIIPushTab(viewer); /* push0 */
539: PetscViewerASCIIPrintf(viewer,"type=nest, rows=%d \n",bx->nb);
541: PetscViewerASCIIPrintf(viewer,"VecNest structure: \n");
542: for (i=0; i<bx->nb; i++) {
543: const VecType type;
544: char name[256] = "",prefix[256] = "";
545: PetscInt NR;
547: VecGetSize(bx->v[i],&NR);
548: VecGetType(bx->v[i],&type);
549: if (((PetscObject)bx->v[i])->name) {PetscSNPrintf(name,sizeof name,"name=\"%s\", ",((PetscObject)bx->v[i])->name);}
550: if (((PetscObject)bx->v[i])->prefix) {PetscSNPrintf(prefix,sizeof prefix,"prefix=\"%s\", ",((PetscObject)bx->v[i])->prefix);}
552: PetscViewerASCIIPrintf(viewer,"(%D) : %s%stype=%s, rows=%D \n",i,name,prefix,type,NR);
554: PetscViewerASCIIPushTab(viewer); /* push1 */
555: VecView(bx->v[i],viewer);
556: PetscViewerASCIIPopTab(viewer); /* pop1 */
557: }
558: PetscViewerASCIIPopTab(viewer); /* pop0 */
559: }
560: return(0);
561: }
565: static PetscErrorCode VecSize_Nest_Recursive(Vec x,PetscBool globalsize,PetscInt *L)
566: {
567: Vec_Nest *bx = (Vec_Nest*)x->data;
568: PetscInt size,i,nr;
569: PetscBool isnest;
573: PetscTypeCompare((PetscObject)x,VECNEST,&isnest);
574: if (!isnest) {
575: /* Not nest */
576: if (globalsize) { VecGetSize(x,&size); }
577: else { VecGetLocalSize(x,&size); }
578: *L = *L + size;
579: return(0);
580: }
582: /* Otherwise we have a nest */
583: bx = (Vec_Nest*)x->data;
584: nr = bx->nb;
586: /* now descend recursively */
587: for (i=0; i<nr; i++) {
588: VecSize_Nest_Recursive(bx->v[i],globalsize,L);
589: }
590: return(0);
591: }
593: /* Returns the sum of the global size of all the consituent vectors in the nest */
596: static PetscErrorCode VecGetSize_Nest(Vec x,PetscInt *N)
597: {
599: *N = x->map->N;
600: return(0);
601: }
603: /* Returns the sum of the local size of all the consituent vectors in the nest */
606: static PetscErrorCode VecGetLocalSize_Nest(Vec x,PetscInt *n)
607: {
609: *n = x->map->n;
610: return(0);
611: }
615: static PetscErrorCode VecMaxPointwiseDivide_Nest(Vec x,Vec y,PetscReal *max)
616: {
617: Vec_Nest *bx = (Vec_Nest*)x->data;
618: Vec_Nest *by = (Vec_Nest*)y->data;
619: PetscInt i,nr;
620: PetscReal local_max,m;
624: VecNestCheckCompatible2(x,1,y,2);
625: nr = bx->nb;
626: m = 0.0;
627: for (i=0; i<nr; i++) {
628: VecMaxPointwiseDivide(bx->v[i],by->v[i],&local_max);
629: if (local_max > m) {
630: m = local_max;
631: }
632: }
633: *max = m;
634: return(0);
635: }
639: static PetscErrorCode VecGetSubVector_Nest(Vec X,IS is,Vec *x)
640: {
641: Vec_Nest *bx = (Vec_Nest*)X->data;
642: PetscInt i;
646: *x = PETSC_NULL;
647: for (i=0;i<bx->nb; i++) {
648: PetscBool issame = PETSC_FALSE;
649: ISEqual(is,bx->is[i],&issame);
650: if (issame) {
651: *x = bx->v[i];
652: PetscObjectReference((PetscObject)(*x));
653: break;
654: }
655: }
656: if (!*x) SETERRQ(((PetscObject)is)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Index set not found in nested Vec");
657: return(0);
658: }
662: static PetscErrorCode VecRestoreSubVector_Nest(Vec X,IS is,Vec *x)
663: {
667: VecDestroy(x);
668: return(0);
669: }
673: static PetscErrorCode VecNestSetOps_Private(struct _VecOps *ops)
674: {
677: /* 0 */
678: ops->duplicate = VecDuplicate_Nest;
679: ops->duplicatevecs = VecDuplicateVecs_Default;
680: ops->destroyvecs = VecDestroyVecs_Default;
681: ops->dot = VecDot_Nest;
682: ops->mdot = VecMDot_Nest;
684: /* 5 */
685: ops->norm = VecNorm_Nest;
686: ops->tdot = VecTDot_Nest;
687: ops->mtdot = VecMTDot_Nest;
688: ops->scale = VecScale_Nest;
689: ops->copy = VecCopy_Nest;
691: /* 10 */
692: ops->set = VecSet_Nest;
693: ops->swap = VecSwap_Nest;
694: ops->axpy = VecAXPY_Nest;
695: ops->axpby = VecAXPBY_Nest;
696: ops->maxpy = VecMAXPY_Nest;
698: /* 15 */
699: ops->aypx = VecAYPX_Nest;
700: ops->waxpy = VecWAXPY_Nest;
701: ops->axpbypcz = 0;
702: ops->pointwisemult = VecPointwiseMult_Nest;
703: ops->pointwisedivide = VecPointwiseDivide_Nest;
704: /* 20 */
705: ops->setvalues = 0;
706: ops->assemblybegin = VecAssemblyBegin_Nest;
707: ops->assemblyend = VecAssemblyEnd_Nest;
708: ops->getarray = 0;
709: ops->getsize = VecGetSize_Nest;
711: /* 25 */
712: ops->getlocalsize = VecGetLocalSize_Nest;
713: ops->restorearray = 0;
714: ops->max = VecMax_Nest;
715: ops->min = VecMin_Nest;
716: ops->setrandom = 0;
718: /* 30 */
719: ops->setoption = 0;
720: ops->setvaluesblocked = 0;
721: ops->destroy = VecDestroy_Nest;
722: ops->view = VecView_Nest;
723: ops->placearray = 0;
725: /* 35 */
726: ops->replacearray = 0;
727: ops->dot_local = VecDot_Nest;
728: ops->tdot_local = VecTDot_Nest;
729: ops->norm_local = VecNorm_Nest;
730: ops->mdot_local = VecMDot_Nest;
732: /* 40 */
733: ops->mtdot_local = VecMTDot_Nest;
734: ops->load = 0;
735: ops->reciprocal = VecReciprocal_Nest;
736: ops->conjugate = VecConjugate_Nest;
737: ops->setlocaltoglobalmapping = 0;
739: /* 45 */
740: ops->setvalueslocal = 0;
741: ops->resetarray = 0;
742: ops->setfromoptions = 0;
743: ops->maxpointwisedivide = VecMaxPointwiseDivide_Nest;
744: ops->load = 0;
746: /* 50 */
747: ops->pointwisemax = 0;
748: ops->pointwisemaxabs = 0;
749: ops->pointwisemin = 0;
750: ops->getvalues = 0;
751: ops->sqrt = 0;
753: /* 55 */
754: ops->abs = 0;
755: ops->exp = 0;
756: ops->shift = 0;
757: ops->create = 0;
758: ops->stridegather = 0;
760: /* 60 */
761: ops->stridescatter = 0;
762: ops->dotnorm2 = 0;
763: ops->getsubvector = VecGetSubVector_Nest;
764: ops->restoresubvector = VecRestoreSubVector_Nest;
766: return(0);
767: }
771: static PetscErrorCode VecNestGetSubVecs_Private(Vec x,PetscInt m,const PetscInt idxm[],Vec vec[])
772: {
773: Vec_Nest *b = (Vec_Nest*)x->data;
774: PetscInt i;
775: PetscInt row;
778: if (!m) return(0);
779: for (i=0; i<m; i++) {
780: row = idxm[i];
781: if (row >= b->nb) SETERRQ2(((PetscObject)x)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,b->nb-1);
782: vec[i] = b->v[row];
783: }
784: return(0);
785: }
790: PetscErrorCode VecNestGetSubVec_Nest(Vec X,PetscInt idxm,Vec *sx)
791: {
795: VecNestGetSubVecs_Private(X,1,&idxm,sx);
796: return(0);
797: }
802: /*@C
803: VecNestGetSubVec - Returns a single, sub-vector from a nest vector.
805: Not collective
807: Input Parameters:
808: . X - nest vector
809: . idxm - index of the vector within the nest
811: Output Parameter:
812: . sx - vector at index idxm within the nest
814: Notes:
816: Level: developer
818: .seealso: VecNestGetSize(), VecNestGetSubVecs()
819: @*/
820: PetscErrorCode VecNestGetSubVec(Vec X,PetscInt idxm,Vec *sx)
821: {
825: PetscUseMethod(X,"VecNestGetSubVec_C",(Vec,PetscInt,Vec*),(X,idxm,sx));
826: return(0);
827: }
832: PetscErrorCode VecNestGetSubVecs_Nest(Vec X,PetscInt *N,Vec **sx)
833: {
834: Vec_Nest *b = (Vec_Nest*)X->data;
836: if (N) { *N = b->nb; }
837: if (sx) { *sx = b->v; }
838: return(0);
839: }
844: /*@C
845: VecNestGetSubVecs - Returns the entire array of vectors defining a nest vector.
847: Not collective
849: Input Parameters:
850: . X - nest vector
852: Output Parameter:
853: . N - number of nested vecs
854: . sx - array of vectors
856: Notes:
858: The user should not free the array sx.
860: Level: developer
862: .seealso: VecNestGetSize(), VecNestGetSubVec()
863: @*/
864: PetscErrorCode VecNestGetSubVecs(Vec X,PetscInt *N,Vec **sx)
865: {
869: PetscUseMethod(X,"VecNestGetSubVecs_C",(Vec,PetscInt*,Vec**),(X,N,sx));
870: return(0);
871: }
876: PetscErrorCode VecNestGetSize_Nest(Vec X,PetscInt *N)
877: {
878: Vec_Nest *b = (Vec_Nest*)X->data;
881: *N = b->nb;
882: return(0);
883: }
888: /*@C
889: VecNestGetSize - Returns the size of the nest vector.
891: Not collective
893: Input Parameters:
894: . X - nest vector
896: Output Parameter:
897: . N - number of nested vecs
899: Notes:
901: Level: developer
903: .seealso: VecNestGetSubVec(), VecNestGetSubVecs()
904: @*/
905: PetscErrorCode VecNestGetSize(Vec X,PetscInt *N)
906: {
912: PetscUseMethod(X,"VecNestGetSize_C",(Vec,PetscInt*),(X,N));
913: return(0);
914: }
918: static PetscErrorCode VecSetUp_Nest_Private(Vec V,PetscInt nb,Vec x[])
919: {
920: Vec_Nest *ctx = (Vec_Nest*)V->data;
921: PetscInt i;
925: if (ctx->setup_called) return(0);
927: ctx->nb = nb;
928: if (ctx->nb < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Cannot create VECNEST with < 0 blocks.");
930: /* Create space */
931: PetscMalloc(ctx->nb*sizeof(Vec),&ctx->v);
932: for (i=0; i<ctx->nb; i++) {
933: ctx->v[i] = x[i];
934: PetscObjectReference((PetscObject)x[i]);
935: /* Do not allocate memory for internal sub blocks */
936: }
938: PetscMalloc(ctx->nb*sizeof(IS),&ctx->is);
940: ctx->setup_called = PETSC_TRUE;
941: return(0);
942: }
946: static PetscErrorCode VecSetUp_NestIS_Private(Vec V,PetscInt nb,IS is[])
947: {
948: Vec_Nest *ctx = (Vec_Nest*)V->data;
949: PetscInt i,offset,m,n,M,N;
953: if (is) { /* Do some consistency checks and reference the is */
954: offset = V->map->rstart;
955: for (i=0; i<ctx->nb; i++) {
956: ISGetSize(is[i],&M);
957: VecGetSize(ctx->v[i],&N);
958: if (M != N) SETERRQ3(((PetscObject)V)->comm,PETSC_ERR_ARG_INCOMP,"In slot %D, IS of size %D is not compatible with Vec of size %D",i,M,N);
959: ISGetLocalSize(is[i],&m);
960: VecGetLocalSize(ctx->v[i],&n);
961: if (m != n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"In slot %D, IS of local size %D is not compatible with Vec of local size %D",i,m,n);
962: #if defined(PETSC_USE_DEBUG)
963: { /* This test can be expensive */
964: PetscInt start;
965: PetscBool contiguous;
966: ISContiguousLocal(is[i],offset,offset+n,&start,&contiguous);
967: if (!contiguous) SETERRQ1(((PetscObject)V)->comm,PETSC_ERR_SUP,"Index set %D is not contiguous with layout of matching vector",i);
968: if (start != 0) SETERRQ1(((PetscObject)V)->comm,PETSC_ERR_SUP,"Index set %D introduces overlap or a hole",i);
969: }
970: #endif
971: PetscObjectReference((PetscObject)is[i]);
972: ctx->is[i] = is[i];
973: offset += n;
974: }
975: } else { /* Create a contiguous ISStride for each entry */
976: offset = V->map->rstart;
977: for (i=0; i<ctx->nb; i++) {
978: PetscInt bs;
979: VecGetLocalSize(ctx->v[i],&n);
980: VecGetBlockSize(ctx->v[i],&bs);
981: ISCreateStride(((PetscObject)ctx->v[i])->comm,n,offset,1,&ctx->is[i]);
982: ISSetBlockSize(ctx->is[i],bs);
983: offset += n;
984: }
985: }
986: return(0);
987: }
991: /*@
992: VecCreateNest - Creates a new vector containing several nested subvectors, each stored separately
994: Collective on Vec
996: Input Parameter:
997: + comm - Communicator for the new Vec
998: . nb - number of nested blocks
999: . is - array of nb index sets describing each nested block, or PETSC_NULL to pack subvectors contiguously
1000: - x - array of nb sub-vectors
1002: Output Parameter:
1003: . Y - new vector
1005: Level: advanced
1007: .seealso: VecCreate(), MatCreateNest(), DMSetVecType(), VECNEST
1008: @*/
1009: PetscErrorCode VecCreateNest(MPI_Comm comm,PetscInt nb,IS is[],Vec x[],Vec *Y)
1010: {
1011: Vec V;
1012: Vec_Nest *s;
1013: PetscInt n,N;
1017: VecCreate(comm,&V);
1018: PetscObjectChangeTypeName((PetscObject)V,VECNEST);
1020: /* allocate and set pointer for implememtation data */
1021: PetscMalloc(sizeof(Vec_Nest),&s);
1022: PetscMemzero(s,sizeof(Vec_Nest));
1023: V->data = (void*)s;
1024: s->setup_called = PETSC_FALSE;
1025: s->nb = -1;
1026: s->v = PETSC_NULL;
1028: VecSetUp_Nest_Private(V,nb,x);
1030: n = N = 0;
1031: VecSize_Nest_Recursive(V,PETSC_TRUE,&N);
1032: VecSize_Nest_Recursive(V,PETSC_FALSE,&n);
1033: PetscLayoutSetSize(V->map,N);
1034: PetscLayoutSetLocalSize(V->map,n);
1035: PetscLayoutSetBlockSize(V->map,1);
1036: PetscLayoutSetUp(V->map);
1038: VecSetUp_NestIS_Private(V,nb,is);
1040: VecNestSetOps_Private(V->ops);
1041: V->petscnative = PETSC_TRUE;
1044: /* expose block api's */
1045: PetscObjectComposeFunctionDynamic((PetscObject)V,"VecNestGetSubVec_C","VecNestGetSubVec_Nest",VecNestGetSubVec_Nest);
1046: PetscObjectComposeFunctionDynamic((PetscObject)V,"VecNestGetSubVecs_C","VecNestGetSubVecs_Nest",VecNestGetSubVecs_Nest);
1047: PetscObjectComposeFunctionDynamic((PetscObject)V,"VecNestGetSize_C","VecNestGetSize_Nest",VecNestGetSize_Nest);
1049: *Y = V;
1050: return(0);
1051: }
1053: /*MC
1054: VECNEST - VECNEST = "nest" - Vector type consisting of nested subvectors, each stored separately.
1056: Level: intermediate
1058: Notes:
1059: This vector type reduces the number of copies for certain solvers applied to multi-physics problems.
1060: It is usually used with MATNEST and DMComposite via DMSetVecType().
1062: .seealso: VecCreate(), VecType, VecCreateNest(), MatCreateNest()
1063: M*/