Actual source code: vecnest.c
petsc-3.7.5 2017-01-01
2: #include <../src/vec/vec/impls/nest/vecnestimpl.h> /*I "petscvec.h" I*/
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(PetscObjectComm((PetscObject)v),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);
55: PetscObjectComposeFunction((PetscObject)v,"",NULL);
56: PetscObjectComposeFunction((PetscObject)v,"",NULL);
57: PetscObjectComposeFunction((PetscObject)v,"",NULL);
58: PetscObjectComposeFunction((PetscObject)v,"",NULL);
59: PetscObjectComposeFunction((PetscObject)v,"",NULL);
61: PetscFree(vs);
62: return(0);
63: }
65: /* supports nested blocks */
68: static PetscErrorCode VecCopy_Nest(Vec x,Vec y)
69: {
70: Vec_Nest *bx = (Vec_Nest*)x->data;
71: Vec_Nest *by = (Vec_Nest*)y->data;
72: PetscInt i;
76: VecNestCheckCompatible2(x,1,y,2);
77: for (i=0; i<bx->nb; i++) {
78: VecCopy(bx->v[i],by->v[i]);
79: }
80: return(0);
81: }
83: /* supports nested blocks */
86: static PetscErrorCode VecDuplicate_Nest(Vec x,Vec *y)
87: {
88: Vec_Nest *bx = (Vec_Nest*)x->data;
89: Vec Y;
90: Vec *sub;
91: PetscInt i;
95: PetscMalloc(sizeof(Vec)*bx->nb,&sub);
96: for (i=0; i<bx->nb; i++) {
97: VecDuplicate(bx->v[i],&sub[i]);
98: }
99: VecCreateNest(PetscObjectComm((PetscObject)x),bx->nb,bx->is,sub,&Y);
100: for (i=0; i<bx->nb; i++) {
101: VecDestroy(&sub[i]);
102: }
103: PetscFree(sub);
104: *y = Y;
105: return(0);
106: }
108: /* supports nested blocks */
111: static PetscErrorCode VecDot_Nest(Vec x,Vec y,PetscScalar *val)
112: {
113: Vec_Nest *bx = (Vec_Nest*)x->data;
114: Vec_Nest *by = (Vec_Nest*)y->data;
115: PetscInt i,nr;
116: PetscScalar x_dot_y,_val;
120: nr = bx->nb;
121: _val = 0.0;
122: for (i=0; i<nr; i++) {
123: VecDot(bx->v[i],by->v[i],&x_dot_y);
124: _val = _val + x_dot_y;
125: }
126: *val = _val;
127: return(0);
128: }
130: /* supports nested blocks */
133: static PetscErrorCode VecTDot_Nest(Vec x,Vec y,PetscScalar *val)
134: {
135: Vec_Nest *bx = (Vec_Nest*)x->data;
136: Vec_Nest *by = (Vec_Nest*)y->data;
137: PetscInt i,nr;
138: PetscScalar x_dot_y,_val;
142: nr = bx->nb;
143: _val = 0.0;
144: for (i=0; i<nr; i++) {
145: VecTDot(bx->v[i],by->v[i],&x_dot_y);
146: _val = _val + x_dot_y;
147: }
148: *val = _val;
149: return(0);
150: }
154: static PetscErrorCode VecDotNorm2_Nest(Vec x,Vec y,PetscScalar *dp, PetscScalar *nm)
155: {
156: Vec_Nest *bx = (Vec_Nest*)x->data;
157: Vec_Nest *by = (Vec_Nest*)y->data;
158: PetscInt i,nr;
159: PetscScalar x_dot_y,_dp,_nm;
160: PetscReal norm2_y;
164: nr = bx->nb;
165: _dp = 0.0;
166: _nm = 0.0;
167: for (i=0; i<nr; i++) {
168: VecDotNorm2(bx->v[i],by->v[i],&x_dot_y,&norm2_y);
169: _dp += x_dot_y;
170: _nm += norm2_y;
171: }
172: *dp = _dp;
173: *nm = _nm;
174: return(0);
175: }
179: static PetscErrorCode VecAXPY_Nest(Vec y,PetscScalar alpha,Vec x)
180: {
181: Vec_Nest *bx = (Vec_Nest*)x->data;
182: Vec_Nest *by = (Vec_Nest*)y->data;
183: PetscInt i,nr;
187: nr = bx->nb;
188: for (i=0; i<nr; i++) {
189: VecAXPY(by->v[i],alpha,bx->v[i]);
190: }
191: return(0);
192: }
196: static PetscErrorCode VecAYPX_Nest(Vec y,PetscScalar alpha,Vec x)
197: {
198: Vec_Nest *bx = (Vec_Nest*)x->data;
199: Vec_Nest *by = (Vec_Nest*)y->data;
200: PetscInt i,nr;
204: nr = bx->nb;
205: for (i=0; i<nr; i++) {
206: VecAYPX(by->v[i],alpha,bx->v[i]);
207: }
208: return(0);
209: }
213: static PetscErrorCode VecAXPBY_Nest(Vec y,PetscScalar alpha,PetscScalar beta,Vec x)
214: {
215: Vec_Nest *bx = (Vec_Nest*)x->data;
216: Vec_Nest *by = (Vec_Nest*)y->data;
217: PetscInt i,nr;
221: nr = bx->nb;
222: for (i=0; i<nr; i++) {
223: VecAXPBY(by->v[i],alpha,beta,bx->v[i]);
224: }
225: return(0);
226: }
230: static PetscErrorCode VecScale_Nest(Vec x,PetscScalar alpha)
231: {
232: Vec_Nest *bx = (Vec_Nest*)x->data;
233: PetscInt i,nr;
237: nr = bx->nb;
238: for (i=0; i<nr; i++) {
239: VecScale(bx->v[i],alpha);
240: }
241: return(0);
242: }
246: static PetscErrorCode VecPointwiseMult_Nest(Vec w,Vec x,Vec y)
247: {
248: Vec_Nest *bx = (Vec_Nest*)x->data;
249: Vec_Nest *by = (Vec_Nest*)y->data;
250: Vec_Nest *bw = (Vec_Nest*)w->data;
251: PetscInt i,nr;
255: VecNestCheckCompatible3(w,1,x,3,y,4);
256: nr = bx->nb;
257: for (i=0; i<nr; i++) {
258: VecPointwiseMult(bw->v[i],bx->v[i],by->v[i]);
259: }
260: return(0);
261: }
265: static PetscErrorCode VecPointwiseDivide_Nest(Vec w,Vec x,Vec y)
266: {
267: Vec_Nest *bx = (Vec_Nest*)x->data;
268: Vec_Nest *by = (Vec_Nest*)y->data;
269: Vec_Nest *bw = (Vec_Nest*)w->data;
270: PetscInt i,nr;
274: VecNestCheckCompatible3(w,1,x,2,y,3);
276: nr = bx->nb;
277: for (i=0; i<nr; i++) {
278: VecPointwiseDivide(bw->v[i],bx->v[i],by->v[i]);
279: }
280: return(0);
281: }
285: static PetscErrorCode VecReciprocal_Nest(Vec x)
286: {
287: Vec_Nest *bx = (Vec_Nest*)x->data;
288: PetscInt i,nr;
292: nr = bx->nb;
293: for (i=0; i<nr; i++) {
294: VecReciprocal(bx->v[i]);
295: }
296: return(0);
297: }
301: static PetscErrorCode VecNorm_Nest(Vec xin,NormType type,PetscReal *z)
302: {
303: Vec_Nest *bx = (Vec_Nest*)xin->data;
304: PetscInt i,nr;
305: PetscReal z_i;
306: PetscReal _z;
310: nr = bx->nb;
311: _z = 0.0;
313: if (type == NORM_2) {
314: PetscScalar dot;
315: VecDot(xin,xin,&dot);
316: _z = PetscAbsScalar(PetscSqrtScalar(dot));
317: } else if (type == NORM_1) {
318: for (i=0; i<nr; i++) {
319: VecNorm(bx->v[i],type,&z_i);
320: _z = _z + z_i;
321: }
322: } else if (type == NORM_INFINITY) {
323: for (i=0; i<nr; i++) {
324: VecNorm(bx->v[i],type,&z_i);
325: if (z_i > _z) _z = z_i;
326: }
327: }
329: *z = _z;
330: return(0);
331: }
335: static PetscErrorCode VecMAXPY_Nest(Vec y,PetscInt nv,const PetscScalar alpha[],Vec *x)
336: {
337: PetscInt v;
341: for (v=0; v<nv; v++) {
342: /* Do axpy on each vector,v */
343: VecAXPY(y,alpha[v],x[v]);
344: }
345: return(0);
346: }
350: static PetscErrorCode VecMDot_Nest(Vec x,PetscInt nv,const Vec y[],PetscScalar *val)
351: {
352: PetscInt j;
356: for (j=0; j<nv; j++) {
357: VecDot(x,y[j],&val[j]);
358: }
359: return(0);
360: }
364: static PetscErrorCode VecMTDot_Nest(Vec x,PetscInt nv,const Vec y[],PetscScalar *val)
365: {
366: PetscInt j;
370: for (j=0; j<nv; j++) {
371: VecTDot(x,y[j],&val[j]);
372: }
373: return(0);
374: }
378: static PetscErrorCode VecSet_Nest(Vec x,PetscScalar alpha)
379: {
380: Vec_Nest *bx = (Vec_Nest*)x->data;
381: PetscInt i,nr;
385: nr = bx->nb;
386: for (i=0; i<nr; i++) {
387: VecSet(bx->v[i],alpha);
388: }
389: return(0);
390: }
394: static PetscErrorCode VecConjugate_Nest(Vec x)
395: {
396: Vec_Nest *bx = (Vec_Nest*)x->data;
397: PetscInt j,nr;
401: nr = bx->nb;
402: for (j=0; j<nr; j++) {
403: VecConjugate(bx->v[j]);
404: }
405: return(0);
406: }
410: static PetscErrorCode VecSwap_Nest(Vec x,Vec y)
411: {
412: Vec_Nest *bx = (Vec_Nest*)x->data;
413: Vec_Nest *by = (Vec_Nest*)y->data;
414: PetscInt i,nr;
418: VecNestCheckCompatible2(x,1,y,2);
419: nr = bx->nb;
420: for (i=0; i<nr; i++) {
421: VecSwap(bx->v[i],by->v[i]);
422: }
423: return(0);
424: }
428: static PetscErrorCode VecWAXPY_Nest(Vec w,PetscScalar alpha,Vec x,Vec y)
429: {
430: Vec_Nest *bx = (Vec_Nest*)x->data;
431: Vec_Nest *by = (Vec_Nest*)y->data;
432: Vec_Nest *bw = (Vec_Nest*)w->data;
433: PetscInt i,nr;
437: VecNestCheckCompatible3(w,1,x,3,y,4);
439: nr = bx->nb;
440: for (i=0; i<nr; i++) {
441: VecWAXPY(bw->v[i],alpha,bx->v[i],by->v[i]);
442: }
443: return(0);
444: }
448: static PetscErrorCode VecMax_Nest_Recursive(Vec x,PetscInt *cnt,PetscInt *p,PetscReal *max)
449: {
450: Vec_Nest *bx;
451: PetscInt i,nr;
452: PetscBool isnest;
453: PetscInt L;
454: PetscInt _entry_loc;
455: PetscReal _entry_val;
459: PetscObjectTypeCompare((PetscObject)x,VECNEST,&isnest);
460: if (!isnest) {
461: /* Not nest */
462: VecMax(x,&_entry_loc,&_entry_val);
463: if (_entry_val > *max) {
464: *max = _entry_val;
465: *p = _entry_loc + *cnt;
466: }
467: VecGetSize(x,&L);
468: *cnt = *cnt + L;
469: return(0);
470: }
472: /* Otherwise we have a nest */
473: bx = (Vec_Nest*)x->data;
474: nr = bx->nb;
476: /* now descend recursively */
477: for (i=0; i<nr; i++) {
478: VecMax_Nest_Recursive(bx->v[i],cnt,p,max);
479: }
480: return(0);
481: }
483: /* supports nested blocks */
486: static PetscErrorCode VecMax_Nest(Vec x,PetscInt *p,PetscReal *max)
487: {
488: PetscInt cnt;
492: cnt = 0;
493: *p = 0;
494: *max = PETSC_MIN_REAL;
495: VecMax_Nest_Recursive(x,&cnt,p,max);
496: return(0);
497: }
501: static PetscErrorCode VecMin_Nest_Recursive(Vec x,PetscInt *cnt,PetscInt *p,PetscReal *min)
502: {
503: Vec_Nest *bx = (Vec_Nest*)x->data;
504: PetscInt i,nr,L,_entry_loc;
505: PetscBool isnest;
506: PetscReal _entry_val;
510: PetscObjectTypeCompare((PetscObject)x,VECNEST,&isnest);
511: if (!isnest) {
512: /* Not nest */
513: VecMin(x,&_entry_loc,&_entry_val);
514: if (_entry_val < *min) {
515: *min = _entry_val;
516: *p = _entry_loc + *cnt;
517: }
518: VecGetSize(x,&L);
519: *cnt = *cnt + L;
520: return(0);
521: }
523: /* Otherwise we have a nest */
524: nr = bx->nb;
526: /* now descend recursively */
527: for (i=0; i<nr; i++) {
528: VecMin_Nest_Recursive(bx->v[i],cnt,p,min);
529: }
530: return(0);
531: }
535: static PetscErrorCode VecMin_Nest(Vec x,PetscInt *p,PetscReal *min)
536: {
537: PetscInt cnt;
541: cnt = 0;
542: *p = 0;
543: *min = PETSC_MAX_REAL;
544: VecMin_Nest_Recursive(x,&cnt,p,min);
545: return(0);
546: }
548: /* supports nested blocks */
551: static PetscErrorCode VecView_Nest(Vec x,PetscViewer viewer)
552: {
553: Vec_Nest *bx = (Vec_Nest*)x->data;
554: PetscBool isascii;
555: PetscInt i;
559: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
560: if (isascii) {
561: PetscViewerASCIIPushTab(viewer);
562: PetscViewerASCIIPrintf(viewer,"VecNest, rows=%D, structure: \n",bx->nb);
563: for (i=0; i<bx->nb; i++) {
564: VecType type;
565: char name[256] = "",prefix[256] = "";
566: PetscInt NR;
568: VecGetSize(bx->v[i],&NR);
569: VecGetType(bx->v[i],&type);
570: if (((PetscObject)bx->v[i])->name) {PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bx->v[i])->name);}
571: if (((PetscObject)bx->v[i])->prefix) {PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bx->v[i])->prefix);}
573: PetscViewerASCIIPrintf(viewer,"(%D) : %s%stype=%s, rows=%D \n",i,name,prefix,type,NR);
575: PetscViewerASCIIPushTab(viewer); /* push1 */
576: VecView(bx->v[i],viewer);
577: PetscViewerASCIIPopTab(viewer); /* pop1 */
578: }
579: PetscViewerASCIIPopTab(viewer); /* pop0 */
580: }
581: return(0);
582: }
586: static PetscErrorCode VecSize_Nest_Recursive(Vec x,PetscBool globalsize,PetscInt *L)
587: {
588: Vec_Nest *bx;
589: PetscInt size,i,nr;
590: PetscBool isnest;
594: PetscObjectTypeCompare((PetscObject)x,VECNEST,&isnest);
595: if (!isnest) {
596: /* Not nest */
597: if (globalsize) { VecGetSize(x,&size); }
598: else { VecGetLocalSize(x,&size); }
599: *L = *L + size;
600: return(0);
601: }
603: /* Otherwise we have a nest */
604: bx = (Vec_Nest*)x->data;
605: nr = bx->nb;
607: /* now descend recursively */
608: for (i=0; i<nr; i++) {
609: VecSize_Nest_Recursive(bx->v[i],globalsize,L);
610: }
611: return(0);
612: }
614: /* Returns the sum of the global size of all the consituent vectors in the nest */
617: static PetscErrorCode VecGetSize_Nest(Vec x,PetscInt *N)
618: {
620: *N = x->map->N;
621: return(0);
622: }
624: /* Returns the sum of the local size of all the consituent vectors in the nest */
627: static PetscErrorCode VecGetLocalSize_Nest(Vec x,PetscInt *n)
628: {
630: *n = x->map->n;
631: return(0);
632: }
636: static PetscErrorCode VecMaxPointwiseDivide_Nest(Vec x,Vec y,PetscReal *max)
637: {
638: Vec_Nest *bx = (Vec_Nest*)x->data;
639: Vec_Nest *by = (Vec_Nest*)y->data;
640: PetscInt i,nr;
641: PetscReal local_max,m;
645: VecNestCheckCompatible2(x,1,y,2);
646: nr = bx->nb;
647: m = 0.0;
648: for (i=0; i<nr; i++) {
649: VecMaxPointwiseDivide(bx->v[i],by->v[i],&local_max);
650: if (local_max > m) m = local_max;
651: }
652: *max = m;
653: return(0);
654: }
658: static PetscErrorCode VecGetSubVector_Nest(Vec X,IS is,Vec *x)
659: {
660: Vec_Nest *bx = (Vec_Nest*)X->data;
661: PetscInt i;
665: *x = NULL;
666: for (i=0; i<bx->nb; i++) {
667: PetscBool issame = PETSC_FALSE;
668: ISEqual(is,bx->is[i],&issame);
669: if (issame) {
670: *x = bx->v[i];
671: PetscObjectReference((PetscObject)(*x));
672: break;
673: }
674: }
675: if (!*x) SETERRQ(PetscObjectComm((PetscObject)is),PETSC_ERR_ARG_OUTOFRANGE,"Index set not found in nested Vec");
676: return(0);
677: }
681: static PetscErrorCode VecRestoreSubVector_Nest(Vec X,IS is,Vec *x)
682: {
686: VecDestroy(x);
687: return(0);
688: }
692: static PetscErrorCode VecGetArray_Nest(Vec X,PetscScalar **x)
693: {
694: Vec_Nest *bx = (Vec_Nest*)X->data;
695: PetscInt i,m,rstart,rend;
699: VecGetOwnershipRange(X,&rstart,&rend);
700: VecGetLocalSize(X,&m);
701: PetscMalloc1(m,x);
702: for (i=0; i<bx->nb; i++) {
703: Vec subvec = bx->v[i];
704: IS isy = bx->is[i];
705: PetscInt j,sm;
706: const PetscInt *ixy;
707: const PetscScalar *y;
708: VecGetLocalSize(subvec,&sm);
709: VecGetArrayRead(subvec,&y);
710: ISGetIndices(isy,&ixy);
711: for (j=0; j<sm; j++) {
712: PetscInt ix = ixy[j];
713: if (ix < rstart || rend <= ix) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for getting array from nonlocal subvec");
714: (*x)[ix-rstart] = y[j];
715: }
716: ISRestoreIndices(isy,&ixy);
717: VecRestoreArrayRead(subvec,&y);
718: }
719: return(0);
720: }
724: static PetscErrorCode VecRestoreArray_Nest(Vec X,PetscScalar **x)
725: {
726: Vec_Nest *bx = (Vec_Nest*)X->data;
727: PetscInt i,m,rstart,rend;
731: VecGetOwnershipRange(X,&rstart,&rend);
732: VecGetLocalSize(X,&m);
733: for (i=0; i<bx->nb; i++) {
734: Vec subvec = bx->v[i];
735: IS isy = bx->is[i];
736: PetscInt j,sm;
737: const PetscInt *ixy;
738: PetscScalar *y;
739: VecGetLocalSize(subvec,&sm);
740: VecGetArray(subvec,&y);
741: ISGetIndices(isy,&ixy);
742: for (j=0; j<sm; j++) {
743: PetscInt ix = ixy[j];
744: if (ix < rstart || rend <= ix) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for getting array from nonlocal subvec");
745: y[j] = (*x)[ix-rstart];
746: }
747: ISRestoreIndices(isy,&ixy);
748: VecRestoreArray(subvec,&y);
749: }
750: PetscFree(*x);
751: return(0);
752: }
757: static PetscErrorCode VecNestSetOps_Private(struct _VecOps *ops)
758: {
760: /* 0 */
761: ops->duplicate = VecDuplicate_Nest;
762: ops->duplicatevecs = VecDuplicateVecs_Default;
763: ops->destroyvecs = VecDestroyVecs_Default;
764: ops->dot = VecDot_Nest;
765: ops->mdot = VecMDot_Nest;
767: /* 5 */
768: ops->norm = VecNorm_Nest;
769: ops->tdot = VecTDot_Nest;
770: ops->mtdot = VecMTDot_Nest;
771: ops->scale = VecScale_Nest;
772: ops->copy = VecCopy_Nest;
774: /* 10 */
775: ops->set = VecSet_Nest;
776: ops->swap = VecSwap_Nest;
777: ops->axpy = VecAXPY_Nest;
778: ops->axpby = VecAXPBY_Nest;
779: ops->maxpy = VecMAXPY_Nest;
781: /* 15 */
782: ops->aypx = VecAYPX_Nest;
783: ops->waxpy = VecWAXPY_Nest;
784: ops->axpbypcz = 0;
785: ops->pointwisemult = VecPointwiseMult_Nest;
786: ops->pointwisedivide = VecPointwiseDivide_Nest;
787: /* 20 */
788: ops->setvalues = 0;
789: ops->assemblybegin = VecAssemblyBegin_Nest;
790: ops->assemblyend = VecAssemblyEnd_Nest;
791: ops->getarray = VecGetArray_Nest;
792: ops->getsize = VecGetSize_Nest;
794: /* 25 */
795: ops->getlocalsize = VecGetLocalSize_Nest;
796: ops->restorearray = VecRestoreArray_Nest;
797: ops->max = VecMax_Nest;
798: ops->min = VecMin_Nest;
799: ops->setrandom = 0;
801: /* 30 */
802: ops->setoption = 0;
803: ops->setvaluesblocked = 0;
804: ops->destroy = VecDestroy_Nest;
805: ops->view = VecView_Nest;
806: ops->placearray = 0;
808: /* 35 */
809: ops->replacearray = 0;
810: ops->dot_local = VecDot_Nest;
811: ops->tdot_local = VecTDot_Nest;
812: ops->norm_local = VecNorm_Nest;
813: ops->mdot_local = VecMDot_Nest;
815: /* 40 */
816: ops->mtdot_local = VecMTDot_Nest;
817: ops->load = 0;
818: ops->reciprocal = VecReciprocal_Nest;
819: ops->conjugate = VecConjugate_Nest;
820: ops->setlocaltoglobalmapping = 0;
822: /* 45 */
823: ops->setvalueslocal = 0;
824: ops->resetarray = 0;
825: ops->setfromoptions = 0;
826: ops->maxpointwisedivide = VecMaxPointwiseDivide_Nest;
827: ops->load = 0;
829: /* 50 */
830: ops->pointwisemax = 0;
831: ops->pointwisemaxabs = 0;
832: ops->pointwisemin = 0;
833: ops->getvalues = 0;
834: ops->sqrt = 0;
836: /* 55 */
837: ops->abs = 0;
838: ops->exp = 0;
839: ops->shift = 0;
840: ops->create = 0;
841: ops->stridegather = 0;
843: /* 60 */
844: ops->stridescatter = 0;
845: ops->dotnorm2 = VecDotNorm2_Nest;
846: ops->getsubvector = VecGetSubVector_Nest;
847: ops->restoresubvector = VecRestoreSubVector_Nest;
848: return(0);
849: }
853: static PetscErrorCode VecNestGetSubVecs_Private(Vec x,PetscInt m,const PetscInt idxm[],Vec vec[])
854: {
855: Vec_Nest *b = (Vec_Nest*)x->data;
856: PetscInt i;
857: PetscInt row;
860: if (!m) return(0);
861: for (i=0; i<m; i++) {
862: row = idxm[i];
863: if (row >= b->nb) SETERRQ2(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,b->nb-1);
864: vec[i] = b->v[row];
865: }
866: return(0);
867: }
871: PetscErrorCode VecNestGetSubVec_Nest(Vec X,PetscInt idxm,Vec *sx)
872: {
876: VecNestGetSubVecs_Private(X,1,&idxm,sx);
877: return(0);
878: }
882: /*@
883: VecNestGetSubVec - Returns a single, sub-vector from a nest vector.
885: Not collective
887: Input Parameters:
888: . X - nest vector
889: . idxm - index of the vector within the nest
891: Output Parameter:
892: . sx - vector at index idxm within the nest
894: Notes:
896: Level: developer
898: .seealso: VecNestGetSize(), VecNestGetSubVecs()
899: @*/
900: PetscErrorCode VecNestGetSubVec(Vec X,PetscInt idxm,Vec *sx)
901: {
905: PetscUseMethod(X,"VecNestGetSubVec_C",(Vec,PetscInt,Vec*),(X,idxm,sx));
906: return(0);
907: }
911: PetscErrorCode VecNestGetSubVecs_Nest(Vec X,PetscInt *N,Vec **sx)
912: {
913: Vec_Nest *b = (Vec_Nest*)X->data;
916: if (N) *N = b->nb;
917: if (sx) *sx = b->v;
918: return(0);
919: }
923: /*@C
924: VecNestGetSubVecs - Returns the entire array of vectors defining a nest vector.
926: Not collective
928: Input Parameters:
929: . X - nest vector
931: Output Parameter:
932: + N - number of nested vecs
933: - sx - array of vectors
935: Notes:
936: The user should not free the array sx.
938: Fortran Notes:
939: The caller must allocate the array to hold the subvectors.
941: Level: developer
943: .seealso: VecNestGetSize(), VecNestGetSubVec()
944: @*/
945: PetscErrorCode VecNestGetSubVecs(Vec X,PetscInt *N,Vec **sx)
946: {
950: PetscUseMethod(X,"VecNestGetSubVecs_C",(Vec,PetscInt*,Vec**),(X,N,sx));
951: return(0);
952: }
956: static PetscErrorCode VecNestSetSubVec_Private(Vec X,PetscInt idxm,Vec x)
957: {
958: Vec_Nest *bx = (Vec_Nest*)X->data;
959: PetscInt i,offset=0,n=0,bs;
960: IS is;
962: PetscBool issame = PETSC_FALSE;
963: PetscInt N=0;
965: /* check if idxm < bx->nb */
966: if (idxm >= bx->nb) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Out of range index value %D maximum %D",idxm,bx->nb);
969: VecDestroy(&bx->v[idxm]); /* destroy the existing vector */
970: VecDuplicate(x,&bx->v[idxm]); /* duplicate the layout of given vector */
971: VecCopy(x,bx->v[idxm]); /* copy the contents of the given vector */
973: /* check if we need to update the IS for the block */
974: offset = X->map->rstart;
975: for (i=0; i<idxm; i++) {
976: n=0;
977: VecGetLocalSize(bx->v[i],&n);
978: offset += n;
979: }
981: /* get the local size and block size */
982: VecGetLocalSize(x,&n);
983: VecGetBlockSize(x,&bs);
985: /* create the new IS */
986: ISCreateStride(PetscObjectComm((PetscObject)x),n,offset,1,&is);
987: ISSetBlockSize(is,bs);
989: /* check if they are equal */
990: ISEqual(is,bx->is[idxm],&issame);
992: if (!issame) {
993: /* The IS of given vector has a different layout compared to the existing block vector.
994: Destroy the existing reference and update the IS. */
995: ISDestroy(&bx->is[idxm]);
996: ISDuplicate(is,&bx->is[idxm]);
997: ISCopy(is,bx->is[idxm]);
999: offset += n;
1000: /* Since the current IS[idxm] changed, we need to update all the subsequent IS */
1001: for (i=idxm+1; i<bx->nb; i++) {
1002: /* get the local size and block size */
1003: VecGetLocalSize(bx->v[i],&n);
1004: VecGetBlockSize(bx->v[i],&bs);
1006: /* destroy the old and create the new IS */
1007: ISDestroy(&bx->is[i]);
1008: ISCreateStride(((PetscObject)bx->v[i])->comm,n,offset,1,&bx->is[i]);
1009: ISSetBlockSize(bx->is[i],bs);
1011: offset += n;
1012: }
1014: n=0;
1015: VecSize_Nest_Recursive(X,PETSC_TRUE,&N);
1016: VecSize_Nest_Recursive(X,PETSC_FALSE,&n);
1017: PetscLayoutSetSize(X->map,N);
1018: PetscLayoutSetLocalSize(X->map,n);
1019: }
1021: ISDestroy(&is);
1022: return(0);
1023: }
1027: PetscErrorCode VecNestSetSubVec_Nest(Vec X,PetscInt idxm,Vec sx)
1028: {
1032: VecNestSetSubVec_Private(X,idxm,sx);
1033: return(0);
1034: }
1038: /*@
1039: VecNestSetSubVec - Set a single component vector in a nest vector at specified index.
1041: Not collective
1043: Input Parameters:
1044: + X - nest vector
1045: . idxm - index of the vector within the nest vector
1046: - sx - vector at index idxm within the nest vector
1048: Notes:
1049: The new vector sx does not have to be of same size as X[idxm]. Arbitrary vector layouts are allowed.
1051: Level: developer
1053: .seealso: VecNestSetSubVecs(), VecNestGetSubVec()
1054: @*/
1055: PetscErrorCode VecNestSetSubVec(Vec X,PetscInt idxm,Vec sx)
1056: {
1060: PetscUseMethod(X,"VecNestSetSubVec_C",(Vec,PetscInt,Vec),(X,idxm,sx));
1061: return(0);
1062: }
1066: PetscErrorCode VecNestSetSubVecs_Nest(Vec X,PetscInt N,PetscInt *idxm,Vec *sx)
1067: {
1068: PetscInt i;
1072: for (i=0; i<N; i++) {
1073: VecNestSetSubVec_Private(X,idxm[i],sx[i]);
1074: }
1075: return(0);
1076: }
1080: /*@C
1081: VecNestSetSubVecs - Sets the component vectors at the specified indices in a nest vector.
1083: Not collective
1085: Input Parameters:
1086: + X - nest vector
1087: . N - number of component vecs in sx
1088: . idxm - indices of component vecs that are to be replaced
1089: - sx - array of vectors
1091: Notes:
1092: The components in the vector array sx do not have to be of the same size as corresponding
1093: components in X. The user can also free the array "sx" after the call.
1095: Level: developer
1097: .seealso: VecNestGetSize(), VecNestGetSubVec()
1098: @*/
1099: PetscErrorCode VecNestSetSubVecs(Vec X,PetscInt N,PetscInt *idxm,Vec *sx)
1100: {
1104: PetscUseMethod(X,"VecNestSetSubVecs_C",(Vec,PetscInt,PetscInt*,Vec*),(X,N,idxm,sx));
1105: return(0);
1106: }
1110: PetscErrorCode VecNestGetSize_Nest(Vec X,PetscInt *N)
1111: {
1112: Vec_Nest *b = (Vec_Nest*)X->data;
1115: *N = b->nb;
1116: return(0);
1117: }
1121: /*@
1122: VecNestGetSize - Returns the size of the nest vector.
1124: Not collective
1126: Input Parameters:
1127: . X - nest vector
1129: Output Parameter:
1130: . N - number of nested vecs
1132: Notes:
1134: Level: developer
1136: .seealso: VecNestGetSubVec(), VecNestGetSubVecs()
1137: @*/
1138: PetscErrorCode VecNestGetSize(Vec X,PetscInt *N)
1139: {
1145: PetscUseMethod(X,"VecNestGetSize_C",(Vec,PetscInt*),(X,N));
1146: return(0);
1147: }
1151: static PetscErrorCode VecSetUp_Nest_Private(Vec V,PetscInt nb,Vec x[])
1152: {
1153: Vec_Nest *ctx = (Vec_Nest*)V->data;
1154: PetscInt i;
1158: if (ctx->setup_called) return(0);
1160: ctx->nb = nb;
1161: if (ctx->nb < 0) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_WRONG,"Cannot create VECNEST with < 0 blocks.");
1163: /* Create space */
1164: PetscMalloc1(ctx->nb,&ctx->v);
1165: for (i=0; i<ctx->nb; i++) {
1166: ctx->v[i] = x[i];
1167: PetscObjectReference((PetscObject)x[i]);
1168: /* Do not allocate memory for internal sub blocks */
1169: }
1171: PetscMalloc1(ctx->nb,&ctx->is);
1173: ctx->setup_called = PETSC_TRUE;
1174: return(0);
1175: }
1179: static PetscErrorCode VecSetUp_NestIS_Private(Vec V,PetscInt nb,IS is[])
1180: {
1181: Vec_Nest *ctx = (Vec_Nest*)V->data;
1182: PetscInt i,offset,m,n,M,N;
1186: if (is) { /* Do some consistency checks and reference the is */
1187: offset = V->map->rstart;
1188: for (i=0; i<ctx->nb; i++) {
1189: ISGetSize(is[i],&M);
1190: VecGetSize(ctx->v[i],&N);
1191: if (M != N) SETERRQ3(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"In slot %D, IS of size %D is not compatible with Vec of size %D",i,M,N);
1192: ISGetLocalSize(is[i],&m);
1193: VecGetLocalSize(ctx->v[i],&n);
1194: 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);
1195: #if defined(PETSC_USE_DEBUG)
1196: { /* This test can be expensive */
1197: PetscInt start;
1198: PetscBool contiguous;
1199: ISContiguousLocal(is[i],offset,offset+n,&start,&contiguous);
1200: if (!contiguous) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Index set %D is not contiguous with layout of matching vector",i);
1201: if (start != 0) SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Index set %D introduces overlap or a hole",i);
1202: }
1203: #endif
1204: PetscObjectReference((PetscObject)is[i]);
1205: ctx->is[i] = is[i];
1206: offset += n;
1207: }
1208: } else { /* Create a contiguous ISStride for each entry */
1209: offset = V->map->rstart;
1210: for (i=0; i<ctx->nb; i++) {
1211: PetscInt bs;
1212: VecGetLocalSize(ctx->v[i],&n);
1213: VecGetBlockSize(ctx->v[i],&bs);
1214: ISCreateStride(((PetscObject)ctx->v[i])->comm,n,offset,1,&ctx->is[i]);
1215: ISSetBlockSize(ctx->is[i],bs);
1216: offset += n;
1217: }
1218: }
1219: return(0);
1220: }
1224: /*@C
1225: VecCreateNest - Creates a new vector containing several nested subvectors, each stored separately
1227: Collective on Vec
1229: Input Parameter:
1230: + comm - Communicator for the new Vec
1231: . nb - number of nested blocks
1232: . is - array of nb index sets describing each nested block, or NULL to pack subvectors contiguously
1233: - x - array of nb sub-vectors
1235: Output Parameter:
1236: . Y - new vector
1238: Level: advanced
1240: .seealso: VecCreate(), MatCreateNest(), DMSetVecType(), VECNEST
1241: @*/
1242: PetscErrorCode VecCreateNest(MPI_Comm comm,PetscInt nb,IS is[],Vec x[],Vec *Y)
1243: {
1244: Vec V;
1245: Vec_Nest *s;
1246: PetscInt n,N;
1250: VecCreate(comm,&V);
1251: PetscObjectChangeTypeName((PetscObject)V,VECNEST);
1253: /* allocate and set pointer for implememtation data */
1254: PetscMalloc(sizeof(Vec_Nest),&s);
1255: PetscMemzero(s,sizeof(Vec_Nest));
1256: V->data = (void*)s;
1257: s->setup_called = PETSC_FALSE;
1258: s->nb = -1;
1259: s->v = NULL;
1261: VecSetUp_Nest_Private(V,nb,x);
1263: n = N = 0;
1264: VecSize_Nest_Recursive(V,PETSC_TRUE,&N);
1265: VecSize_Nest_Recursive(V,PETSC_FALSE,&n);
1266: PetscLayoutSetSize(V->map,N);
1267: PetscLayoutSetLocalSize(V->map,n);
1268: PetscLayoutSetBlockSize(V->map,1);
1269: PetscLayoutSetUp(V->map);
1271: VecSetUp_NestIS_Private(V,nb,is);
1273: VecNestSetOps_Private(V->ops);
1274: V->petscnative = PETSC_FALSE;
1277: /* expose block api's */
1278: PetscObjectComposeFunction((PetscObject)V,"VecNestGetSubVec_C",VecNestGetSubVec_Nest);
1279: PetscObjectComposeFunction((PetscObject)V,"VecNestGetSubVecs_C",VecNestGetSubVecs_Nest);
1280: PetscObjectComposeFunction((PetscObject)V,"VecNestSetSubVec_C",VecNestSetSubVec_Nest);
1281: PetscObjectComposeFunction((PetscObject)V,"VecNestSetSubVecs_C",VecNestSetSubVecs_Nest);
1282: PetscObjectComposeFunction((PetscObject)V,"VecNestGetSize_C",VecNestGetSize_Nest);
1284: *Y = V;
1285: return(0);
1286: }
1288: /*MC
1289: VECNEST - VECNEST = "nest" - Vector type consisting of nested subvectors, each stored separately.
1291: Level: intermediate
1293: Notes:
1294: This vector type reduces the number of copies for certain solvers applied to multi-physics problems.
1295: It is usually used with MATNEST and DMComposite via DMSetVecType().
1297: .seealso: VecCreate(), VecType, VecCreateNest(), MatCreateNest()
1298: M*/