Actual source code: vecnest.c


  2: #include <../src/vec/vec/impls/nest/vecnestimpl.h>

  4: /* check all blocks are filled */
  5: static PetscErrorCode VecAssemblyBegin_Nest(Vec v)
  6: {
  7:   Vec_Nest *vs = (Vec_Nest *)v->data;
  8:   PetscInt  i;

 10:   PetscFunctionBegin;
 11:   for (i = 0; i < vs->nb; i++) {
 12:     PetscCheck(vs->v[i], PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Nest  vector cannot contain NULL blocks");
 13:     PetscCall(VecAssemblyBegin(vs->v[i]));
 14:   }
 15:   PetscFunctionReturn(PETSC_SUCCESS);
 16: }

 18: static PetscErrorCode VecAssemblyEnd_Nest(Vec v)
 19: {
 20:   Vec_Nest *vs = (Vec_Nest *)v->data;
 21:   PetscInt  i;

 23:   PetscFunctionBegin;
 24:   for (i = 0; i < vs->nb; i++) PetscCall(VecAssemblyEnd(vs->v[i]));
 25:   PetscFunctionReturn(PETSC_SUCCESS);
 26: }

 28: static PetscErrorCode VecDestroy_Nest(Vec v)
 29: {
 30:   Vec_Nest *vs = (Vec_Nest *)v->data;
 31:   PetscInt  i;

 33:   PetscFunctionBegin;
 34:   if (vs->v) {
 35:     for (i = 0; i < vs->nb; i++) PetscCall(VecDestroy(&vs->v[i]));
 36:     PetscCall(PetscFree(vs->v));
 37:   }
 38:   for (i = 0; i < vs->nb; i++) PetscCall(ISDestroy(&vs->is[i]));
 39:   PetscCall(PetscFree(vs->is));
 40:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "VecNestGetSubVec_C", NULL));
 41:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "VecNestGetSubVecs_C", NULL));
 42:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "VecNestSetSubVec_C", NULL));
 43:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "VecNestSetSubVecs_C", NULL));
 44:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "VecNestGetSize_C", NULL));

 46:   PetscCall(PetscFree(vs));
 47:   PetscFunctionReturn(PETSC_SUCCESS);
 48: }

 50: static PetscErrorCode VecCopy_Nest(Vec x, Vec y)
 51: {
 52:   Vec_Nest *bx = (Vec_Nest *)x->data;
 53:   Vec_Nest *by = (Vec_Nest *)y->data;
 54:   PetscInt  i;

 56:   PetscFunctionBegin;
 57:   PetscCheckTypeName(y, VECNEST);
 58:   VecNestCheckCompatible2(x, 1, y, 2);
 59:   for (i = 0; i < bx->nb; i++) PetscCall(VecCopy(bx->v[i], by->v[i]));
 60:   PetscFunctionReturn(PETSC_SUCCESS);
 61: }

 63: static PetscErrorCode VecDuplicate_Nest(Vec x, Vec *y)
 64: {
 65:   Vec_Nest *bx = (Vec_Nest *)x->data;
 66:   Vec       Y;
 67:   Vec      *sub;
 68:   PetscInt  i;

 70:   PetscFunctionBegin;
 71:   PetscCall(PetscMalloc1(bx->nb, &sub));
 72:   for (i = 0; i < bx->nb; i++) PetscCall(VecDuplicate(bx->v[i], &sub[i]));
 73:   PetscCall(VecCreateNest(PetscObjectComm((PetscObject)x), bx->nb, bx->is, sub, &Y));
 74:   for (i = 0; i < bx->nb; i++) PetscCall(VecDestroy(&sub[i]));
 75:   PetscCall(PetscFree(sub));
 76:   *y = Y;
 77:   PetscFunctionReturn(PETSC_SUCCESS);
 78: }

 80: static PetscErrorCode VecDot_Nest(Vec x, Vec y, PetscScalar *val)
 81: {
 82:   Vec_Nest   *bx = (Vec_Nest *)x->data;
 83:   Vec_Nest   *by = (Vec_Nest *)y->data;
 84:   PetscInt    i, nr;
 85:   PetscScalar x_dot_y, _val;

 87:   PetscFunctionBegin;
 88:   nr   = bx->nb;
 89:   _val = 0.0;
 90:   for (i = 0; i < nr; i++) {
 91:     PetscCall(VecDot(bx->v[i], by->v[i], &x_dot_y));
 92:     _val = _val + x_dot_y;
 93:   }
 94:   *val = _val;
 95:   PetscFunctionReturn(PETSC_SUCCESS);
 96: }

 98: static PetscErrorCode VecTDot_Nest(Vec x, Vec y, PetscScalar *val)
 99: {
100:   Vec_Nest   *bx = (Vec_Nest *)x->data;
101:   Vec_Nest   *by = (Vec_Nest *)y->data;
102:   PetscInt    i, nr;
103:   PetscScalar x_dot_y, _val;

105:   PetscFunctionBegin;
106:   nr   = bx->nb;
107:   _val = 0.0;
108:   for (i = 0; i < nr; i++) {
109:     PetscCall(VecTDot(bx->v[i], by->v[i], &x_dot_y));
110:     _val = _val + x_dot_y;
111:   }
112:   *val = _val;
113:   PetscFunctionReturn(PETSC_SUCCESS);
114: }

116: static PetscErrorCode VecDotNorm2_Nest(Vec x, Vec y, PetscScalar *dp, PetscScalar *nm)
117: {
118:   Vec_Nest   *bx = (Vec_Nest *)x->data;
119:   Vec_Nest   *by = (Vec_Nest *)y->data;
120:   PetscInt    i, nr;
121:   PetscScalar x_dot_y, _dp, _nm;
122:   PetscReal   norm2_y;

124:   PetscFunctionBegin;
125:   nr  = bx->nb;
126:   _dp = 0.0;
127:   _nm = 0.0;
128:   for (i = 0; i < nr; i++) {
129:     PetscCall(VecDotNorm2(bx->v[i], by->v[i], &x_dot_y, &norm2_y));
130:     _dp += x_dot_y;
131:     _nm += norm2_y;
132:   }
133:   *dp = _dp;
134:   *nm = _nm;
135:   PetscFunctionReturn(PETSC_SUCCESS);
136: }

138: static PetscErrorCode VecAXPY_Nest(Vec y, PetscScalar alpha, Vec x)
139: {
140:   Vec_Nest *bx = (Vec_Nest *)x->data;
141:   Vec_Nest *by = (Vec_Nest *)y->data;
142:   PetscInt  i, nr;

144:   PetscFunctionBegin;
145:   nr = bx->nb;
146:   for (i = 0; i < nr; i++) PetscCall(VecAXPY(by->v[i], alpha, bx->v[i]));
147:   PetscFunctionReturn(PETSC_SUCCESS);
148: }

150: static PetscErrorCode VecAYPX_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;

156:   PetscFunctionBegin;
157:   nr = bx->nb;
158:   for (i = 0; i < nr; i++) PetscCall(VecAYPX(by->v[i], alpha, bx->v[i]));
159:   PetscFunctionReturn(PETSC_SUCCESS);
160: }

162: static PetscErrorCode VecAXPBY_Nest(Vec y, PetscScalar alpha, PetscScalar beta, Vec x)
163: {
164:   Vec_Nest *bx = (Vec_Nest *)x->data;
165:   Vec_Nest *by = (Vec_Nest *)y->data;
166:   PetscInt  i, nr;

168:   PetscFunctionBegin;
169:   nr = bx->nb;
170:   for (i = 0; i < nr; i++) PetscCall(VecAXPBY(by->v[i], alpha, beta, bx->v[i]));
171:   PetscFunctionReturn(PETSC_SUCCESS);
172: }

174: static PetscErrorCode VecAXPBYPCZ_Nest(Vec z, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec x, Vec y)
175: {
176:   Vec_Nest *bx = (Vec_Nest *)x->data;
177:   Vec_Nest *by = (Vec_Nest *)y->data;
178:   Vec_Nest *bz = (Vec_Nest *)z->data;
179:   PetscInt  i, nr;

181:   PetscFunctionBegin;
182:   nr = bx->nb;
183:   for (i = 0; i < nr; i++) PetscCall(VecAXPBYPCZ(bz->v[i], alpha, beta, gamma, bx->v[i], by->v[i]));
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: static PetscErrorCode VecScale_Nest(Vec x, PetscScalar alpha)
188: {
189:   Vec_Nest *bx = (Vec_Nest *)x->data;
190:   PetscInt  i, nr;

192:   PetscFunctionBegin;
193:   nr = bx->nb;
194:   for (i = 0; i < nr; i++) PetscCall(VecScale(bx->v[i], alpha));
195:   PetscFunctionReturn(PETSC_SUCCESS);
196: }

198: static PetscErrorCode VecPointwiseMult_Nest(Vec w, Vec x, Vec y)
199: {
200:   Vec_Nest *bx = (Vec_Nest *)x->data;
201:   Vec_Nest *by = (Vec_Nest *)y->data;
202:   Vec_Nest *bw = (Vec_Nest *)w->data;
203:   PetscInt  i, nr;

205:   PetscFunctionBegin;
206:   VecNestCheckCompatible3(w, 1, x, 2, y, 3);
207:   nr = bx->nb;
208:   for (i = 0; i < nr; i++) PetscCall(VecPointwiseMult(bw->v[i], bx->v[i], by->v[i]));
209:   PetscFunctionReturn(PETSC_SUCCESS);
210: }

212: static PetscErrorCode VecPointwiseDivide_Nest(Vec w, Vec x, Vec y)
213: {
214:   Vec_Nest *bx = (Vec_Nest *)x->data;
215:   Vec_Nest *by = (Vec_Nest *)y->data;
216:   Vec_Nest *bw = (Vec_Nest *)w->data;
217:   PetscInt  i, nr;

219:   PetscFunctionBegin;
220:   VecNestCheckCompatible3(w, 1, x, 2, y, 3);

222:   nr = bx->nb;
223:   for (i = 0; i < nr; i++) PetscCall(VecPointwiseDivide(bw->v[i], bx->v[i], by->v[i]));
224:   PetscFunctionReturn(PETSC_SUCCESS);
225: }

227: static PetscErrorCode VecReciprocal_Nest(Vec x)
228: {
229:   Vec_Nest *bx = (Vec_Nest *)x->data;
230:   PetscInt  i, nr;

232:   PetscFunctionBegin;
233:   nr = bx->nb;
234:   for (i = 0; i < nr; i++) PetscCall(VecReciprocal(bx->v[i]));
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: static PetscErrorCode VecNorm_Nest(Vec xin, NormType type, PetscReal *z)
239: {
240:   Vec_Nest *bx = (Vec_Nest *)xin->data;
241:   PetscInt  i, nr;
242:   PetscReal z_i;
243:   PetscReal _z;

245:   PetscFunctionBegin;
246:   nr = bx->nb;
247:   _z = 0.0;

249:   if (type == NORM_2) {
250:     PetscScalar dot;
251:     PetscCall(VecDot(xin, xin, &dot));
252:     _z = PetscAbsScalar(PetscSqrtScalar(dot));
253:   } else if (type == NORM_1) {
254:     for (i = 0; i < nr; i++) {
255:       PetscCall(VecNorm(bx->v[i], type, &z_i));
256:       _z = _z + z_i;
257:     }
258:   } else if (type == NORM_INFINITY) {
259:     for (i = 0; i < nr; i++) {
260:       PetscCall(VecNorm(bx->v[i], type, &z_i));
261:       if (z_i > _z) _z = z_i;
262:     }
263:   }

265:   *z = _z;
266:   PetscFunctionReturn(PETSC_SUCCESS);
267: }

269: static PetscErrorCode VecMAXPY_Nest(Vec y, PetscInt nv, const PetscScalar alpha[], Vec *x)
270: {
271:   PetscInt v;

273:   PetscFunctionBegin;
274:   for (v = 0; v < nv; v++) {
275:     /* Do axpy on each vector,v */
276:     PetscCall(VecAXPY(y, alpha[v], x[v]));
277:   }
278:   PetscFunctionReturn(PETSC_SUCCESS);
279: }

281: static PetscErrorCode VecMDot_Nest(Vec x, PetscInt nv, const Vec y[], PetscScalar *val)
282: {
283:   PetscInt j;

285:   PetscFunctionBegin;
286:   for (j = 0; j < nv; j++) PetscCall(VecDot(x, y[j], &val[j]));
287:   PetscFunctionReturn(PETSC_SUCCESS);
288: }

290: static PetscErrorCode VecMTDot_Nest(Vec x, PetscInt nv, const Vec y[], PetscScalar *val)
291: {
292:   PetscInt j;

294:   PetscFunctionBegin;
295:   for (j = 0; j < nv; j++) PetscCall(VecTDot(x, y[j], &val[j]));
296:   PetscFunctionReturn(PETSC_SUCCESS);
297: }

299: static PetscErrorCode VecSet_Nest(Vec x, PetscScalar alpha)
300: {
301:   Vec_Nest *bx = (Vec_Nest *)x->data;
302:   PetscInt  i, nr;

304:   PetscFunctionBegin;
305:   nr = bx->nb;
306:   for (i = 0; i < nr; i++) PetscCall(VecSet(bx->v[i], alpha));
307:   PetscFunctionReturn(PETSC_SUCCESS);
308: }

310: static PetscErrorCode VecConjugate_Nest(Vec x)
311: {
312:   Vec_Nest *bx = (Vec_Nest *)x->data;
313:   PetscInt  j, nr;

315:   PetscFunctionBegin;
316:   nr = bx->nb;
317:   for (j = 0; j < nr; j++) PetscCall(VecConjugate(bx->v[j]));
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }

321: static PetscErrorCode VecSwap_Nest(Vec x, Vec y)
322: {
323:   Vec_Nest *bx = (Vec_Nest *)x->data;
324:   Vec_Nest *by = (Vec_Nest *)y->data;
325:   PetscInt  i, nr;

327:   PetscFunctionBegin;
328:   VecNestCheckCompatible2(x, 1, y, 2);
329:   nr = bx->nb;
330:   for (i = 0; i < nr; i++) PetscCall(VecSwap(bx->v[i], by->v[i]));
331:   PetscFunctionReturn(PETSC_SUCCESS);
332: }

334: static PetscErrorCode VecWAXPY_Nest(Vec w, PetscScalar alpha, Vec x, Vec y)
335: {
336:   Vec_Nest *bx = (Vec_Nest *)x->data;
337:   Vec_Nest *by = (Vec_Nest *)y->data;
338:   Vec_Nest *bw = (Vec_Nest *)w->data;
339:   PetscInt  i, nr;

341:   PetscFunctionBegin;
342:   VecNestCheckCompatible3(w, 1, x, 3, y, 4);

344:   nr = bx->nb;
345:   for (i = 0; i < nr; i++) PetscCall(VecWAXPY(bw->v[i], alpha, bx->v[i], by->v[i]));
346:   PetscFunctionReturn(PETSC_SUCCESS);
347: }

349: static PetscErrorCode VecMin_Nest(Vec x, PetscInt *p, PetscReal *v)
350: {
351:   PetscInt  i, lp = -1, lw = -1;
352:   PetscReal lv;
353:   Vec_Nest *bx = (Vec_Nest *)x->data;

355:   PetscFunctionBegin;
356:   *v = PETSC_MAX_REAL;
357:   for (i = 0; i < bx->nb; i++) {
358:     PetscInt tp;
359:     PetscCall(VecMin(bx->v[i], &tp, &lv));
360:     if (lv < *v) {
361:       *v = lv;
362:       lw = i;
363:       lp = tp;
364:     }
365:   }
366:   if (p && lw > -1) {
367:     PetscInt        st, en;
368:     const PetscInt *idxs;

370:     *p = -1;
371:     PetscCall(VecGetOwnershipRange(bx->v[lw], &st, &en));
372:     if (lp >= st && lp < en) {
373:       PetscCall(ISGetIndices(bx->is[lw], &idxs));
374:       *p = idxs[lp - st];
375:       PetscCall(ISRestoreIndices(bx->is[lw], &idxs));
376:     }
377:     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, p, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)x)));
378:   }
379:   PetscFunctionReturn(PETSC_SUCCESS);
380: }

382: static PetscErrorCode VecMax_Nest(Vec x, PetscInt *p, PetscReal *v)
383: {
384:   PetscInt  i, lp = -1, lw = -1;
385:   PetscReal lv;
386:   Vec_Nest *bx = (Vec_Nest *)x->data;

388:   PetscFunctionBegin;
389:   *v = PETSC_MIN_REAL;
390:   for (i = 0; i < bx->nb; i++) {
391:     PetscInt tp;
392:     PetscCall(VecMax(bx->v[i], &tp, &lv));
393:     if (lv > *v) {
394:       *v = lv;
395:       lw = i;
396:       lp = tp;
397:     }
398:   }
399:   if (p && lw > -1) {
400:     PetscInt        st, en;
401:     const PetscInt *idxs;

403:     *p = -1;
404:     PetscCall(VecGetOwnershipRange(bx->v[lw], &st, &en));
405:     if (lp >= st && lp < en) {
406:       PetscCall(ISGetIndices(bx->is[lw], &idxs));
407:       *p = idxs[lp - st];
408:       PetscCall(ISRestoreIndices(bx->is[lw], &idxs));
409:     }
410:     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, p, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)x)));
411:   }
412:   PetscFunctionReturn(PETSC_SUCCESS);
413: }

415: static PetscErrorCode VecView_Nest(Vec x, PetscViewer viewer)
416: {
417:   Vec_Nest *bx = (Vec_Nest *)x->data;
418:   PetscBool isascii;
419:   PetscInt  i;

421:   PetscFunctionBegin;
422:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
423:   if (isascii) {
424:     PetscCall(PetscViewerASCIIPushTab(viewer));
425:     PetscCall(PetscViewerASCIIPrintf(viewer, "VecNest, rows=%" PetscInt_FMT ",  structure: \n", bx->nb));
426:     for (i = 0; i < bx->nb; i++) {
427:       VecType  type;
428:       char     name[256] = "", prefix[256] = "";
429:       PetscInt NR;

431:       PetscCall(VecGetSize(bx->v[i], &NR));
432:       PetscCall(VecGetType(bx->v[i], &type));
433:       if (((PetscObject)bx->v[i])->name) PetscCall(PetscSNPrintf(name, sizeof(name), "name=\"%s\", ", ((PetscObject)bx->v[i])->name));
434:       if (((PetscObject)bx->v[i])->prefix) PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "prefix=\"%s\", ", ((PetscObject)bx->v[i])->prefix));

436:       PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") : %s%stype=%s, rows=%" PetscInt_FMT " \n", i, name, prefix, type, NR));

438:       PetscCall(PetscViewerASCIIPushTab(viewer)); /* push1 */
439:       PetscCall(VecView(bx->v[i], viewer));
440:       PetscCall(PetscViewerASCIIPopTab(viewer)); /* pop1 */
441:     }
442:     PetscCall(PetscViewerASCIIPopTab(viewer)); /* pop0 */
443:   }
444:   PetscFunctionReturn(PETSC_SUCCESS);
445: }

447: static PetscErrorCode VecSize_Nest_Recursive(Vec x, PetscBool globalsize, PetscInt *L)
448: {
449:   Vec_Nest *bx;
450:   PetscInt  size, i, nr;
451:   PetscBool isnest;

453:   PetscFunctionBegin;
454:   PetscCall(PetscObjectTypeCompare((PetscObject)x, VECNEST, &isnest));
455:   if (!isnest) {
456:     /* Not nest */
457:     if (globalsize) PetscCall(VecGetSize(x, &size));
458:     else PetscCall(VecGetLocalSize(x, &size));
459:     *L = *L + size;
460:     PetscFunctionReturn(PETSC_SUCCESS);
461:   }

463:   /* Otherwise we have a nest */
464:   bx = (Vec_Nest *)x->data;
465:   nr = bx->nb;

467:   /* now descend recursively */
468:   for (i = 0; i < nr; i++) PetscCall(VecSize_Nest_Recursive(bx->v[i], globalsize, L));
469:   PetscFunctionReturn(PETSC_SUCCESS);
470: }

472: /* Returns the sum of the global size of all the constituent vectors in the nest */
473: static PetscErrorCode VecGetSize_Nest(Vec x, PetscInt *N)
474: {
475:   PetscFunctionBegin;
476:   *N = x->map->N;
477:   PetscFunctionReturn(PETSC_SUCCESS);
478: }

480: /* Returns the sum of the local size of all the constituent vectors in the nest */
481: static PetscErrorCode VecGetLocalSize_Nest(Vec x, PetscInt *n)
482: {
483:   PetscFunctionBegin;
484:   *n = x->map->n;
485:   PetscFunctionReturn(PETSC_SUCCESS);
486: }

488: static PetscErrorCode VecMaxPointwiseDivide_Nest(Vec x, Vec y, PetscReal *max)
489: {
490:   Vec_Nest *bx = (Vec_Nest *)x->data;
491:   Vec_Nest *by = (Vec_Nest *)y->data;
492:   PetscInt  i, nr;
493:   PetscReal local_max, m;

495:   PetscFunctionBegin;
496:   VecNestCheckCompatible2(x, 1, y, 2);
497:   nr = bx->nb;
498:   m  = 0.0;
499:   for (i = 0; i < nr; i++) {
500:     PetscCall(VecMaxPointwiseDivide(bx->v[i], by->v[i], &local_max));
501:     if (local_max > m) m = local_max;
502:   }
503:   *max = m;
504:   PetscFunctionReturn(PETSC_SUCCESS);
505: }

507: static PetscErrorCode VecGetSubVector_Nest(Vec X, IS is, Vec *x)
508: {
509:   Vec_Nest *bx = (Vec_Nest *)X->data;
510:   PetscInt  i;

512:   PetscFunctionBegin;
513:   *x = NULL;
514:   for (i = 0; i < bx->nb; i++) {
515:     PetscBool issame = PETSC_FALSE;
516:     PetscCall(ISEqual(is, bx->is[i], &issame));
517:     if (issame) {
518:       *x = bx->v[i];
519:       PetscCall(PetscObjectReference((PetscObject)(*x)));
520:       break;
521:     }
522:   }
523:   PetscCheck(*x, PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_OUTOFRANGE, "Index set not found in nested Vec");
524:   PetscFunctionReturn(PETSC_SUCCESS);
525: }

527: static PetscErrorCode VecRestoreSubVector_Nest(Vec X, IS is, Vec *x)
528: {
529:   PetscFunctionBegin;
530:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
531:   PetscCall(VecDestroy(x));
532:   PetscFunctionReturn(PETSC_SUCCESS);
533: }

535: static PetscErrorCode VecGetArray_Nest(Vec X, PetscScalar **x)
536: {
537:   Vec_Nest *bx = (Vec_Nest *)X->data;
538:   PetscInt  i, m, rstart, rend;

540:   PetscFunctionBegin;
541:   PetscCall(VecGetOwnershipRange(X, &rstart, &rend));
542:   PetscCall(VecGetLocalSize(X, &m));
543:   PetscCall(PetscMalloc1(m, x));
544:   for (i = 0; i < bx->nb; i++) {
545:     Vec                subvec = bx->v[i];
546:     IS                 isy    = bx->is[i];
547:     PetscInt           j, sm;
548:     const PetscInt    *ixy;
549:     const PetscScalar *y;
550:     PetscCall(VecGetLocalSize(subvec, &sm));
551:     PetscCall(VecGetArrayRead(subvec, &y));
552:     PetscCall(ISGetIndices(isy, &ixy));
553:     for (j = 0; j < sm; j++) {
554:       PetscInt ix = ixy[j];
555:       PetscCheck(ix >= rstart && rend > ix, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for getting array from nonlocal subvec");
556:       (*x)[ix - rstart] = y[j];
557:     }
558:     PetscCall(ISRestoreIndices(isy, &ixy));
559:     PetscCall(VecRestoreArrayRead(subvec, &y));
560:   }
561:   PetscFunctionReturn(PETSC_SUCCESS);
562: }

564: static PetscErrorCode VecRestoreArray_Nest(Vec X, PetscScalar **x)
565: {
566:   Vec_Nest *bx = (Vec_Nest *)X->data;
567:   PetscInt  i, m, rstart, rend;

569:   PetscFunctionBegin;
570:   PetscCall(VecGetOwnershipRange(X, &rstart, &rend));
571:   PetscCall(VecGetLocalSize(X, &m));
572:   for (i = 0; i < bx->nb; i++) {
573:     Vec             subvec = bx->v[i];
574:     IS              isy    = bx->is[i];
575:     PetscInt        j, sm;
576:     const PetscInt *ixy;
577:     PetscScalar    *y;
578:     PetscCall(VecGetLocalSize(subvec, &sm));
579:     PetscCall(VecGetArray(subvec, &y));
580:     PetscCall(ISGetIndices(isy, &ixy));
581:     for (j = 0; j < sm; j++) {
582:       PetscInt ix = ixy[j];
583:       PetscCheck(ix >= rstart && rend > ix, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for getting array from nonlocal subvec");
584:       y[j] = (*x)[ix - rstart];
585:     }
586:     PetscCall(ISRestoreIndices(isy, &ixy));
587:     PetscCall(VecRestoreArray(subvec, &y));
588:   }
589:   PetscCall(PetscFree(*x));
590:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
591:   PetscFunctionReturn(PETSC_SUCCESS);
592: }

594: static PetscErrorCode VecRestoreArrayRead_Nest(Vec X, const PetscScalar **x)
595: {
596:   PetscFunctionBegin;
597:   PetscCall(PetscFree(*x));
598:   PetscFunctionReturn(PETSC_SUCCESS);
599: }

601: static PetscErrorCode VecConcatenate_Nest(PetscInt nx, const Vec X[], Vec *Y, IS *x_is[])
602: {
603:   PetscFunctionBegin;
604:   PetscCheck(nx <= 0, PetscObjectComm((PetscObject)(*X)), PETSC_ERR_SUP, "VecConcatenate() is not supported for VecNest");
605:   PetscFunctionReturn(PETSC_SUCCESS);
606: }

608: static PetscErrorCode VecCreateLocalVector_Nest(Vec v, Vec *w)
609: {
610:   Vec      *ww;
611:   IS       *wis;
612:   Vec_Nest *bv = (Vec_Nest *)v->data;
613:   PetscInt  i;

615:   PetscFunctionBegin;
616:   PetscCall(PetscMalloc2(bv->nb, &ww, bv->nb, &wis));
617:   for (i = 0; i < bv->nb; i++) PetscCall(VecCreateLocalVector(bv->v[i], &ww[i]));
618:   for (i = 0; i < bv->nb; i++) {
619:     PetscInt off;

621:     PetscCall(VecGetOwnershipRange(v, &off, NULL));
622:     PetscCall(ISOnComm(bv->is[i], PetscObjectComm((PetscObject)ww[i]), PETSC_COPY_VALUES, &wis[i]));
623:     PetscCall(ISShift(wis[i], -off, wis[i]));
624:   }
625:   PetscCall(VecCreateNest(PETSC_COMM_SELF, bv->nb, wis, ww, w));
626:   for (i = 0; i < bv->nb; i++) {
627:     PetscCall(VecDestroy(&ww[i]));
628:     PetscCall(ISDestroy(&wis[i]));
629:   }
630:   PetscCall(PetscFree2(ww, wis));
631:   PetscFunctionReturn(PETSC_SUCCESS);
632: }

634: static PetscErrorCode VecGetLocalVector_Nest(Vec v, Vec w)
635: {
636:   Vec_Nest *bv = (Vec_Nest *)v->data;
637:   Vec_Nest *bw = (Vec_Nest *)w->data;
638:   PetscInt  i;

640:   PetscFunctionBegin;
641:   PetscCheckSameType(v, 1, w, 2);
642:   PetscCheck(bv->nb = bw->nb, PetscObjectComm((PetscObject)w), PETSC_ERR_ARG_WRONG, "Invalid local vector");
643:   for (i = 0; i < bv->nb; i++) PetscCall(VecGetLocalVector(bv->v[i], bw->v[i]));
644:   PetscFunctionReturn(PETSC_SUCCESS);
645: }

647: static PetscErrorCode VecRestoreLocalVector_Nest(Vec v, Vec w)
648: {
649:   Vec_Nest *bv = (Vec_Nest *)v->data;
650:   Vec_Nest *bw = (Vec_Nest *)w->data;
651:   PetscInt  i;

653:   PetscFunctionBegin;
654:   PetscCheckSameType(v, 1, w, 2);
655:   PetscCheck(bv->nb = bw->nb, PetscObjectComm((PetscObject)w), PETSC_ERR_ARG_WRONG, "Invalid local vector");
656:   for (i = 0; i < bv->nb; i++) PetscCall(VecRestoreLocalVector(bv->v[i], bw->v[i]));
657:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
658:   PetscFunctionReturn(PETSC_SUCCESS);
659: }

661: static PetscErrorCode VecGetLocalVectorRead_Nest(Vec v, Vec w)
662: {
663:   Vec_Nest *bv = (Vec_Nest *)v->data;
664:   Vec_Nest *bw = (Vec_Nest *)w->data;
665:   PetscInt  i;

667:   PetscFunctionBegin;
668:   PetscCheckSameType(v, 1, w, 2);
669:   PetscCheck(bv->nb = bw->nb, PetscObjectComm((PetscObject)w), PETSC_ERR_ARG_WRONG, "Invalid local vector");
670:   for (i = 0; i < bv->nb; i++) PetscCall(VecGetLocalVectorRead(bv->v[i], bw->v[i]));
671:   PetscFunctionReturn(PETSC_SUCCESS);
672: }

674: static PetscErrorCode VecRestoreLocalVectorRead_Nest(Vec v, Vec w)
675: {
676:   Vec_Nest *bv = (Vec_Nest *)v->data;
677:   Vec_Nest *bw = (Vec_Nest *)w->data;
678:   PetscInt  i;

680:   PetscFunctionBegin;
681:   PetscCheckSameType(v, 1, w, 2);
682:   PetscCheck(bv->nb = bw->nb, PetscObjectComm((PetscObject)w), PETSC_ERR_ARG_WRONG, "Invalid local vector");
683:   for (i = 0; i < bv->nb; i++) PetscCall(VecRestoreLocalVectorRead(bv->v[i], bw->v[i]));
684:   PetscFunctionReturn(PETSC_SUCCESS);
685: }

687: static PetscErrorCode VecSetRandom_Nest(Vec v, PetscRandom r)
688: {
689:   Vec_Nest *bv = (Vec_Nest *)v->data;
690:   PetscInt  i;

692:   PetscFunctionBegin;
693:   for (i = 0; i < bv->nb; i++) PetscCall(VecSetRandom(bv->v[i], r));
694:   PetscFunctionReturn(PETSC_SUCCESS);
695: }

697: static PetscErrorCode VecNestSetOps_Private(struct _VecOps *ops)
698: {
699:   PetscFunctionBegin;
700:   ops->duplicate               = VecDuplicate_Nest;
701:   ops->duplicatevecs           = VecDuplicateVecs_Default;
702:   ops->destroyvecs             = VecDestroyVecs_Default;
703:   ops->dot                     = VecDot_Nest;
704:   ops->mdot                    = VecMDot_Nest;
705:   ops->norm                    = VecNorm_Nest;
706:   ops->tdot                    = VecTDot_Nest;
707:   ops->mtdot                   = VecMTDot_Nest;
708:   ops->scale                   = VecScale_Nest;
709:   ops->copy                    = VecCopy_Nest;
710:   ops->set                     = VecSet_Nest;
711:   ops->swap                    = VecSwap_Nest;
712:   ops->axpy                    = VecAXPY_Nest;
713:   ops->axpby                   = VecAXPBY_Nest;
714:   ops->maxpy                   = VecMAXPY_Nest;
715:   ops->aypx                    = VecAYPX_Nest;
716:   ops->waxpy                   = VecWAXPY_Nest;
717:   ops->axpbypcz                = NULL;
718:   ops->pointwisemult           = VecPointwiseMult_Nest;
719:   ops->pointwisedivide         = VecPointwiseDivide_Nest;
720:   ops->setvalues               = NULL;
721:   ops->assemblybegin           = VecAssemblyBegin_Nest;
722:   ops->assemblyend             = VecAssemblyEnd_Nest;
723:   ops->getarray                = VecGetArray_Nest;
724:   ops->getsize                 = VecGetSize_Nest;
725:   ops->getlocalsize            = VecGetLocalSize_Nest;
726:   ops->restorearray            = VecRestoreArray_Nest;
727:   ops->restorearrayread        = VecRestoreArrayRead_Nest;
728:   ops->max                     = VecMax_Nest;
729:   ops->min                     = VecMin_Nest;
730:   ops->setrandom               = NULL;
731:   ops->setoption               = NULL;
732:   ops->setvaluesblocked        = NULL;
733:   ops->destroy                 = VecDestroy_Nest;
734:   ops->view                    = VecView_Nest;
735:   ops->placearray              = NULL;
736:   ops->replacearray            = NULL;
737:   ops->dot_local               = VecDot_Nest;
738:   ops->tdot_local              = VecTDot_Nest;
739:   ops->norm_local              = VecNorm_Nest;
740:   ops->mdot_local              = VecMDot_Nest;
741:   ops->mtdot_local             = VecMTDot_Nest;
742:   ops->load                    = NULL;
743:   ops->reciprocal              = VecReciprocal_Nest;
744:   ops->conjugate               = VecConjugate_Nest;
745:   ops->setlocaltoglobalmapping = NULL;
746:   ops->setvalueslocal          = NULL;
747:   ops->resetarray              = NULL;
748:   ops->setfromoptions          = NULL;
749:   ops->maxpointwisedivide      = VecMaxPointwiseDivide_Nest;
750:   ops->load                    = NULL;
751:   ops->pointwisemax            = NULL;
752:   ops->pointwisemaxabs         = NULL;
753:   ops->pointwisemin            = NULL;
754:   ops->getvalues               = NULL;
755:   ops->sqrt                    = NULL;
756:   ops->abs                     = NULL;
757:   ops->exp                     = NULL;
758:   ops->shift                   = NULL;
759:   ops->create                  = NULL;
760:   ops->stridegather            = NULL;
761:   ops->stridescatter           = NULL;
762:   ops->dotnorm2                = VecDotNorm2_Nest;
763:   ops->getsubvector            = VecGetSubVector_Nest;
764:   ops->restoresubvector        = VecRestoreSubVector_Nest;
765:   ops->axpbypcz                = VecAXPBYPCZ_Nest;
766:   ops->concatenate             = VecConcatenate_Nest;
767:   ops->createlocalvector       = VecCreateLocalVector_Nest;
768:   ops->getlocalvector          = VecGetLocalVector_Nest;
769:   ops->getlocalvectorread      = VecGetLocalVectorRead_Nest;
770:   ops->restorelocalvector      = VecRestoreLocalVector_Nest;
771:   ops->restorelocalvectorread  = VecRestoreLocalVectorRead_Nest;
772:   ops->setrandom               = VecSetRandom_Nest;
773:   PetscFunctionReturn(PETSC_SUCCESS);
774: }

776: static PetscErrorCode VecNestGetSubVecs_Private(Vec x, PetscInt m, const PetscInt idxm[], Vec vec[])
777: {
778:   Vec_Nest *b = (Vec_Nest *)x->data;
779:   PetscInt  i;
780:   PetscInt  row;

782:   PetscFunctionBegin;
783:   if (!m) PetscFunctionReturn(PETSC_SUCCESS);
784:   for (i = 0; i < m; i++) {
785:     row = idxm[i];
786:     PetscCheck(row < b->nb, PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, row, b->nb - 1);
787:     vec[i] = b->v[row];
788:   }
789:   PetscFunctionReturn(PETSC_SUCCESS);
790: }

792: PetscErrorCode VecNestGetSubVec_Nest(Vec X, PetscInt idxm, Vec *sx)
793: {
794:   PetscFunctionBegin;
795:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
796:   PetscCall(VecNestGetSubVecs_Private(X, 1, &idxm, sx));
797:   PetscFunctionReturn(PETSC_SUCCESS);
798: }

800: /*@
801:  VecNestGetSubVec - Returns a single, sub-vector from a nest vector.

803:  Not Collective

805:  Input Parameters:
806: +  X  - nest vector
807: -  idxm - index of the vector within the nest

809:  Output Parameter:
810: .  sx - vector at index `idxm` within the nest

812:  Level: developer

814: .seealso: `VECNEST`,  [](chapter_vectors), `Vec`, `VecType`, `VecNestGetSize()`, `VecNestGetSubVecs()`
815: @*/
816: PetscErrorCode VecNestGetSubVec(Vec X, PetscInt idxm, Vec *sx)
817: {
818:   PetscFunctionBegin;
819:   PetscUseMethod(X, "VecNestGetSubVec_C", (Vec, PetscInt, Vec *), (X, idxm, sx));
820:   PetscFunctionReturn(PETSC_SUCCESS);
821: }

823: PetscErrorCode VecNestGetSubVecs_Nest(Vec X, PetscInt *N, Vec **sx)
824: {
825:   Vec_Nest *b = (Vec_Nest *)X->data;

827:   PetscFunctionBegin;
828:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
829:   if (N) *N = b->nb;
830:   if (sx) *sx = b->v;
831:   PetscFunctionReturn(PETSC_SUCCESS);
832: }

834: /*@C
835:  VecNestGetSubVecs - Returns the entire array of vectors defining a nest vector.

837:  Not collective

839:  Input Parameter:
840: .  X  - nest vector

842:  Output Parameters:
843: +  N - number of nested vecs
844: -  sx - array of vectors

846:  Level: developer

848:  Note:
849:  The user should not free the array `sx`.

851:  Fortran Note:
852:  The caller must allocate the array to hold the subvectors.

854: .seealso: `VECNEST`,  [](chapter_vectors), `Vec`, `VecType`, `VecNestGetSize()`, `VecNestGetSubVec()`
855: @*/
856: PetscErrorCode VecNestGetSubVecs(Vec X, PetscInt *N, Vec **sx)
857: {
858:   PetscFunctionBegin;
859:   PetscUseMethod(X, "VecNestGetSubVecs_C", (Vec, PetscInt *, Vec **), (X, N, sx));
860:   PetscFunctionReturn(PETSC_SUCCESS);
861: }

863: static PetscErrorCode VecNestSetSubVec_Private(Vec X, PetscInt idxm, Vec x)
864: {
865:   Vec_Nest *bx = (Vec_Nest *)X->data;
866:   PetscInt  i, offset = 0, n = 0, bs;
867:   IS        is;
868:   PetscBool issame = PETSC_FALSE;
869:   PetscInt  N      = 0;

871:   /* check if idxm < bx->nb */
872:   PetscCheck(idxm < bx->nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " maximum %" PetscInt_FMT, idxm, bx->nb);

874:   PetscFunctionBegin;
875:   PetscCall(VecDestroy(&bx->v[idxm]));      /* destroy the existing vector */
876:   PetscCall(VecDuplicate(x, &bx->v[idxm])); /* duplicate the layout of given vector */
877:   PetscCall(VecCopy(x, bx->v[idxm]));       /* copy the contents of the given vector */

879:   /* check if we need to update the IS for the block */
880:   offset = X->map->rstart;
881:   for (i = 0; i < idxm; i++) {
882:     n = 0;
883:     PetscCall(VecGetLocalSize(bx->v[i], &n));
884:     offset += n;
885:   }

887:   /* get the local size and block size */
888:   PetscCall(VecGetLocalSize(x, &n));
889:   PetscCall(VecGetBlockSize(x, &bs));

891:   /* create the new IS */
892:   PetscCall(ISCreateStride(PetscObjectComm((PetscObject)x), n, offset, 1, &is));
893:   PetscCall(ISSetBlockSize(is, bs));

895:   /* check if they are equal */
896:   PetscCall(ISEqual(is, bx->is[idxm], &issame));

898:   if (!issame) {
899:     /* The IS of given vector has a different layout compared to the existing block vector.
900:      Destroy the existing reference and update the IS. */
901:     PetscCall(ISDestroy(&bx->is[idxm]));
902:     PetscCall(ISDuplicate(is, &bx->is[idxm]));
903:     PetscCall(ISCopy(is, bx->is[idxm]));

905:     offset += n;
906:     /* Since the current IS[idxm] changed, we need to update all the subsequent IS */
907:     for (i = idxm + 1; i < bx->nb; i++) {
908:       /* get the local size and block size */
909:       PetscCall(VecGetLocalSize(bx->v[i], &n));
910:       PetscCall(VecGetBlockSize(bx->v[i], &bs));

912:       /* destroy the old and create the new IS */
913:       PetscCall(ISDestroy(&bx->is[i]));
914:       PetscCall(ISCreateStride(((PetscObject)bx->v[i])->comm, n, offset, 1, &bx->is[i]));
915:       PetscCall(ISSetBlockSize(bx->is[i], bs));

917:       offset += n;
918:     }

920:     n = 0;
921:     PetscCall(VecSize_Nest_Recursive(X, PETSC_TRUE, &N));
922:     PetscCall(VecSize_Nest_Recursive(X, PETSC_FALSE, &n));
923:     PetscCall(PetscLayoutSetSize(X->map, N));
924:     PetscCall(PetscLayoutSetLocalSize(X->map, n));
925:   }

927:   PetscCall(ISDestroy(&is));
928:   PetscFunctionReturn(PETSC_SUCCESS);
929: }

931: PetscErrorCode VecNestSetSubVec_Nest(Vec X, PetscInt idxm, Vec sx)
932: {
933:   PetscFunctionBegin;
934:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
935:   PetscCall(VecNestSetSubVec_Private(X, idxm, sx));
936:   PetscFunctionReturn(PETSC_SUCCESS);
937: }

939: /*@
940:    VecNestSetSubVec - Set a single component vector in a nest vector at specified index.

942:    Not collective

944:    Input Parameters:
945: +  X  - nest vector
946: .  idxm - index of the vector within the nest vector
947: -  sx - vector at index `idxm` within the nest vector

949:    Level: developer

951:    Note:
952:    The new vector `sx` does not have to be of same size as X[idxm]. Arbitrary vector layouts are allowed.

954: .seealso: `VECNEST`,  [](chapter_vectors), `Vec`, `VecType`, `VecNestSetSubVecs()`, `VecNestGetSubVec()`
955: @*/
956: PetscErrorCode VecNestSetSubVec(Vec X, PetscInt idxm, Vec sx)
957: {
958:   PetscFunctionBegin;
959:   PetscUseMethod(X, "VecNestSetSubVec_C", (Vec, PetscInt, Vec), (X, idxm, sx));
960:   PetscFunctionReturn(PETSC_SUCCESS);
961: }

963: PetscErrorCode VecNestSetSubVecs_Nest(Vec X, PetscInt N, PetscInt *idxm, Vec *sx)
964: {
965:   PetscInt i;

967:   PetscFunctionBegin;
968:   PetscCall(PetscObjectStateIncrease((PetscObject)X));
969:   for (i = 0; i < N; i++) PetscCall(VecNestSetSubVec_Private(X, idxm[i], sx[i]));
970:   PetscFunctionReturn(PETSC_SUCCESS);
971: }

973: /*@C
974:    VecNestSetSubVecs - Sets the component vectors at the specified indices in a nest vector.

976:    Not collective

978:    Input Parameters:
979: +  X  - nest vector
980: .  N - number of component vecs in `sx`
981: .  idxm - indices of component vectors that are to be replaced
982: -  sx - array of vectors

984:    Level: developer

986:    Note:
987:    The components in the vector array `sx` do not have to be of the same size as corresponding
988:    components in `X`. The user can also free the array `sx` after the call.

990: .seealso: `VECNEST`,  [](chapter_vectors), `Vec`, `VecType`, `VecNestGetSize()`, `VecNestGetSubVec()`
991: @*/
992: PetscErrorCode VecNestSetSubVecs(Vec X, PetscInt N, PetscInt *idxm, Vec *sx)
993: {
994:   PetscFunctionBegin;
995:   PetscUseMethod(X, "VecNestSetSubVecs_C", (Vec, PetscInt, PetscInt *, Vec *), (X, N, idxm, sx));
996:   PetscFunctionReturn(PETSC_SUCCESS);
997: }

999: PetscErrorCode VecNestGetSize_Nest(Vec X, PetscInt *N)
1000: {
1001:   Vec_Nest *b = (Vec_Nest *)X->data;

1003:   PetscFunctionBegin;
1004:   *N = b->nb;
1005:   PetscFunctionReturn(PETSC_SUCCESS);
1006: }

1008: /*@
1009:  VecNestGetSize - Returns the size of the nest vector.

1011:  Not collective

1013:  Input Parameter:
1014: .  X  - nest vector

1016:  Output Parameter:
1017: .  N - number of nested vecs

1019:  Level: developer

1021: .seealso: `VECNEST`,  [](chapter_vectors), `Vec`, `VecType`, `VecNestGetSubVec()`, `VecNestGetSubVecs()`
1022: @*/
1023: PetscErrorCode VecNestGetSize(Vec X, PetscInt *N)
1024: {
1025:   PetscFunctionBegin;
1028:   PetscUseMethod(X, "VecNestGetSize_C", (Vec, PetscInt *), (X, N));
1029:   PetscFunctionReturn(PETSC_SUCCESS);
1030: }

1032: static PetscErrorCode VecSetUp_Nest_Private(Vec V, PetscInt nb, Vec x[])
1033: {
1034:   Vec_Nest *ctx = (Vec_Nest *)V->data;
1035:   PetscInt  i;

1037:   PetscFunctionBegin;
1038:   if (ctx->setup_called) PetscFunctionReturn(PETSC_SUCCESS);

1040:   ctx->nb = nb;
1041:   PetscCheck(ctx->nb >= 0, PetscObjectComm((PetscObject)V), PETSC_ERR_ARG_WRONG, "Cannot create VECNEST with < 0 blocks.");

1043:   /* Create space */
1044:   PetscCall(PetscMalloc1(ctx->nb, &ctx->v));
1045:   for (i = 0; i < ctx->nb; i++) {
1046:     ctx->v[i] = x[i];
1047:     PetscCall(PetscObjectReference((PetscObject)x[i]));
1048:     /* Do not allocate memory for internal sub blocks */
1049:   }

1051:   PetscCall(PetscMalloc1(ctx->nb, &ctx->is));

1053:   ctx->setup_called = PETSC_TRUE;
1054:   PetscFunctionReturn(PETSC_SUCCESS);
1055: }

1057: static PetscErrorCode VecSetUp_NestIS_Private(Vec V, PetscInt nb, IS is[])
1058: {
1059:   Vec_Nest *ctx = (Vec_Nest *)V->data;
1060:   PetscInt  i, offset, m, n, M, N;

1062:   PetscFunctionBegin;
1063:   if (is) { /* Do some consistency checks and reference the is */
1064:     offset = V->map->rstart;
1065:     for (i = 0; i < ctx->nb; i++) {
1066:       PetscCall(ISGetSize(is[i], &M));
1067:       PetscCall(VecGetSize(ctx->v[i], &N));
1068:       PetscCheck(M == N, PetscObjectComm((PetscObject)V), PETSC_ERR_ARG_INCOMP, "In slot %" PetscInt_FMT ", IS of size %" PetscInt_FMT " is not compatible with Vec of size %" PetscInt_FMT, i, M, N);
1069:       PetscCall(ISGetLocalSize(is[i], &m));
1070:       PetscCall(VecGetLocalSize(ctx->v[i], &n));
1071:       PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "In slot %" PetscInt_FMT ", IS of local size %" PetscInt_FMT " is not compatible with Vec of local size %" PetscInt_FMT, i, m, n);
1072:       PetscCall(PetscObjectReference((PetscObject)is[i]));
1073:       ctx->is[i] = is[i];
1074:       offset += n;
1075:     }
1076:   } else { /* Create a contiguous ISStride for each entry */
1077:     offset = V->map->rstart;
1078:     for (i = 0; i < ctx->nb; i++) {
1079:       PetscInt bs;
1080:       PetscCall(VecGetLocalSize(ctx->v[i], &n));
1081:       PetscCall(VecGetBlockSize(ctx->v[i], &bs));
1082:       PetscCall(ISCreateStride(((PetscObject)ctx->v[i])->comm, n, offset, 1, &ctx->is[i]));
1083:       PetscCall(ISSetBlockSize(ctx->is[i], bs));
1084:       offset += n;
1085:     }
1086:   }
1087:   PetscFunctionReturn(PETSC_SUCCESS);
1088: }

1090: /*@C
1091:    VecCreateNest - Creates a new vector containing several nested subvectors, each stored separately

1093:    Collective

1095:    Input Parameters:
1096: +  comm - Communicator for the new `Vec`
1097: .  nb - number of nested blocks
1098: .  is - array of `nb` index sets describing each nested block, or `NULL` to pack subvectors contiguously
1099: -  x - array of `nb` sub-vectors

1101:    Output Parameter:
1102: .  Y - new vector

1104:    Level: advanced

1106: .seealso: `VECNEST`,  [](chapter_vectors), `Vec`, `VecType`, `VecCreate()`, `MatCreateNest()`, `DMSetVecType()`, `VECNEST`
1107: @*/
1108: PetscErrorCode VecCreateNest(MPI_Comm comm, PetscInt nb, IS is[], Vec x[], Vec *Y)
1109: {
1110:   Vec       V;
1111:   Vec_Nest *s;
1112:   PetscInt  n, N;

1114:   PetscFunctionBegin;
1115:   PetscCall(VecCreate(comm, &V));
1116:   PetscCall(PetscObjectChangeTypeName((PetscObject)V, VECNEST));

1118:   /* allocate and set pointer for implementation data */
1119:   PetscCall(PetscNew(&s));
1120:   V->data         = (void *)s;
1121:   s->setup_called = PETSC_FALSE;
1122:   s->nb           = -1;
1123:   s->v            = NULL;

1125:   PetscCall(VecSetUp_Nest_Private(V, nb, x));

1127:   n = N = 0;
1128:   PetscCall(VecSize_Nest_Recursive(V, PETSC_TRUE, &N));
1129:   PetscCall(VecSize_Nest_Recursive(V, PETSC_FALSE, &n));
1130:   PetscCall(PetscLayoutSetSize(V->map, N));
1131:   PetscCall(PetscLayoutSetLocalSize(V->map, n));
1132:   PetscCall(PetscLayoutSetBlockSize(V->map, 1));
1133:   PetscCall(PetscLayoutSetUp(V->map));

1135:   PetscCall(VecSetUp_NestIS_Private(V, nb, is));

1137:   PetscCall(VecNestSetOps_Private(V->ops));
1138:   V->petscnative = PETSC_FALSE;

1140:   /* expose block api's */
1141:   PetscCall(PetscObjectComposeFunction((PetscObject)V, "VecNestGetSubVec_C", VecNestGetSubVec_Nest));
1142:   PetscCall(PetscObjectComposeFunction((PetscObject)V, "VecNestGetSubVecs_C", VecNestGetSubVecs_Nest));
1143:   PetscCall(PetscObjectComposeFunction((PetscObject)V, "VecNestSetSubVec_C", VecNestSetSubVec_Nest));
1144:   PetscCall(PetscObjectComposeFunction((PetscObject)V, "VecNestSetSubVecs_C", VecNestSetSubVecs_Nest));
1145:   PetscCall(PetscObjectComposeFunction((PetscObject)V, "VecNestGetSize_C", VecNestGetSize_Nest));

1147:   *Y = V;
1148:   PetscFunctionReturn(PETSC_SUCCESS);
1149: }

1151: /*MC
1152:   VECNEST - VECNEST = "nest" - Vector type consisting of nested subvectors, each stored separately.

1154:   Level: intermediate

1156:   Notes:
1157:   This vector type reduces the number of copies for certain solvers applied to multi-physics problems.
1158:   It is usually used with `MATNEST` and `DMCOMPOSITE` via `DMSetVecType()`.

1160: .seealso: [](chapter_vectors), `Vec`, `VecType`, `VecCreate()`, `VecType`, `VecCreateNest()`, `MatCreateNest()`
1161: M*/