Actual source code: dmproject.c

petsc-3.14.3 2021-01-09
Report Typos and Errors

  2: #include <petsc/private/dmimpl.h>
  3: #include <petscdm.h>
  4: #include <petscdmplex.h>
  5: #include <petscksp.h>
  6: #include <petscblaslapack.h>

  8: typedef struct _projectConstraintsCtx
  9: {
 10:   DM  dm;
 11:   Vec mask;
 12: }
 13: projectConstraintsCtx;

 15: PetscErrorCode MatMult_GlobalToLocalNormal(Mat CtC, Vec x, Vec y)
 16: {
 17:   DM                    dm;
 18:   Vec                   local, mask;
 19:   projectConstraintsCtx *ctx;
 20:   PetscErrorCode        ierr;

 23:   MatShellGetContext(CtC,&ctx);
 24:   dm   = ctx->dm;
 25:   mask = ctx->mask;
 26:   DMGetLocalVector(dm,&local);
 27:   DMGlobalToLocalBegin(dm,x,INSERT_VALUES,local);
 28:   DMGlobalToLocalEnd(dm,x,INSERT_VALUES,local);
 29:   if (mask) {VecPointwiseMult(local,mask,local);}
 30:   VecSet(y,0.);
 31:   DMLocalToGlobalBegin(dm,local,ADD_VALUES,y);
 32:   DMLocalToGlobalEnd(dm,local,ADD_VALUES,y);
 33:   DMRestoreLocalVector(dm,&local);
 34:   return(0);
 35: }

 37: static PetscErrorCode DMGlobalToLocalSolve_project1 (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
 38: {
 39:   PetscInt f;

 42:   for (f = 0; f < Nf; f++) {
 43:     u[f] = 1.;
 44:   }
 45:   return(0);
 46: }

 48: /*@
 49:   DMGlobalToLocalSolve - Solve for the global vector that is mapped to a given local vector by DMGlobalToLocalBegin()/DMGlobalToLocalEnd() with mode
 50:   = INSERT_VALUES.  It is assumed that the sum of all the local vector sizes is greater than or equal to the global vector size, so the solution is
 51:   a least-squares solution.  It is also assumed that DMLocalToGlobalBegin()/DMLocalToGlobalEnd() with mode = ADD_VALUES is the adjoint of the
 52:   global-to-local map, so that the least-squares solution may be found by the normal equations.

 54:   collective

 56:   Input Parameters:
 57: + dm - The DM object
 58: . x - The local vector
 59: - y - The global vector: the input value of globalVec is used as an initial guess

 61:   Output Parameters:
 62: . y - The least-squares solution

 64:   Level: advanced

 66:   Note: If the DM is of type DMPLEX, then y is the solution of L' * D * L * y = L' * D * x, where D is a diagonal mask that is 1 for every point in
 67:   the union of the closures of the local cells and 0 otherwise.  This difference is only relevant if there are anchor points that are not in the
 68:   closure of any local cell (see DMPlexGetAnchors()/DMPlexSetAnchors()).

 70: .seealso: DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToGlobalBegin(), DMLocalToGlobalEnd(), DMPlexGetAnchors(), DMPlexSetAnchors()
 71: @*/
 72: PetscErrorCode DMGlobalToLocalSolve(DM dm, Vec x, Vec y)
 73: {
 74:   Mat                   CtC;
 75:   PetscInt              n, N, cStart, cEnd, c;
 76:   PetscBool             isPlex;
 77:   KSP                   ksp;
 78:   PC                    pc;
 79:   Vec                   global, mask=NULL;
 80:   projectConstraintsCtx ctx;

 84:   PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&isPlex);
 85:   if (isPlex) {
 86:     /* mark points in the closure */
 87:     DMCreateLocalVector(dm,&mask);
 88:     VecSet(mask,0.0);
 89:     DMPlexGetSimplexOrBoxCells(dm,0,&cStart,&cEnd);
 90:     if (cEnd > cStart) {
 91:       PetscScalar *ones;
 92:       PetscInt numValues, i;

 94:       DMPlexVecGetClosure(dm,NULL,mask,cStart,&numValues,NULL);
 95:       PetscMalloc1(numValues,&ones);
 96:       for (i = 0; i < numValues; i++) {
 97:         ones[i] = 1.;
 98:       }
 99:       for (c = cStart; c < cEnd; c++) {
100:         DMPlexVecSetClosure(dm,NULL,mask,c,ones,INSERT_VALUES);
101:       }
102:       PetscFree(ones);
103:     }
104:   }
105:   else {
106:     PetscBool hasMask;

108:     DMHasNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &hasMask);
109:     if (!hasMask) {
110:       PetscErrorCode (**func) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
111:       void            **ctx;
112:       PetscInt          Nf, f;

114:       DMGetNumFields(dm, &Nf);
115:       PetscMalloc2(Nf, &func, Nf, &ctx);
116:       for (f = 0; f < Nf; ++f) {
117:         func[f] = DMGlobalToLocalSolve_project1;
118:         ctx[f]  = NULL;
119:       }
120:       DMGetNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
121:       DMProjectFunctionLocal(dm,0.0,func,ctx,INSERT_ALL_VALUES,mask);
122:       DMRestoreNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
123:       PetscFree2(func, ctx);
124:     }
125:     DMGetNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
126:   }
127:   ctx.dm   = dm;
128:   ctx.mask = mask;
129:   VecGetSize(y,&N);
130:   VecGetLocalSize(y,&n);
131:   MatCreate(PetscObjectComm((PetscObject)dm),&CtC);
132:   MatSetSizes(CtC,n,n,N,N);
133:   MatSetType(CtC,MATSHELL);
134:   MatSetUp(CtC);
135:   MatShellSetContext(CtC,&ctx);
136:   MatShellSetOperation(CtC,MATOP_MULT,(void(*)(void))MatMult_GlobalToLocalNormal);
137:   KSPCreate(PetscObjectComm((PetscObject)dm),&ksp);
138:   KSPSetOperators(ksp,CtC,CtC);
139:   KSPSetType(ksp,KSPCG);
140:   KSPGetPC(ksp,&pc);
141:   PCSetType(pc,PCNONE);
142:   KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
143:   KSPSetUp(ksp);
144:   DMGetGlobalVector(dm,&global);
145:   VecSet(global,0.);
146:   if (mask) {VecPointwiseMult(x,mask,x);}
147:   DMLocalToGlobalBegin(dm,x,ADD_VALUES,global);
148:   DMLocalToGlobalEnd(dm,x,ADD_VALUES,global);
149:   KSPSolve(ksp,global,y);
150:   DMRestoreGlobalVector(dm,&global);
151:   /* clean up */
152:   KSPDestroy(&ksp);
153:   MatDestroy(&CtC);
154:   if (isPlex) {
155:     VecDestroy(&mask);
156:   }
157:   else {
158:     DMRestoreNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
159:   }

161:   return(0);
162: }

164: /*@C
165:   DMProjectField - This projects the given function of the input fields into the function space provided, putting the coefficients in a global vector.

167:   Collective on DM

169:   Input Parameters:
170: + dm      - The DM
171: . time    - The time
172: . U       - The input field vector
173: . funcs   - The functions to evaluate, one per field
174: - mode    - The insertion mode for values

176:   Output Parameter:
177: . X       - The output vector

179:    Calling sequence of func:
180: $    func(PetscInt dim, PetscInt Nf, PetscInt NfAux,
181: $         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
182: $         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
183: $         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]);

185: +  dim          - The spatial dimension
186: .  Nf           - The number of input fields
187: .  NfAux        - The number of input auxiliary fields
188: .  uOff         - The offset of each field in u[]
189: .  uOff_x       - The offset of each field in u_x[]
190: .  u            - The field values at this point in space
191: .  u_t          - The field time derivative at this point in space (or NULL)
192: .  u_x          - The field derivatives at this point in space
193: .  aOff         - The offset of each auxiliary field in u[]
194: .  aOff_x       - The offset of each auxiliary field in u_x[]
195: .  a            - The auxiliary field values at this point in space
196: .  a_t          - The auxiliary field time derivative at this point in space (or NULL)
197: .  a_x          - The auxiliary field derivatives at this point in space
198: .  t            - The current time
199: .  x            - The coordinates of this point
200: .  numConstants - The number of constants
201: .  constants    - The value of each constant
202: -  f            - The value of the function at this point in space

204:   Note: There are three different DMs that potentially interact in this function. The output DM, dm, specifies the layout of the values calculates by funcs.
205:   The input DM, attached to U, may be different. For example, you can input the solution over the full domain, but output over a piece of the boundary, or
206:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary DM, attached to the
207:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

209:   Level: intermediate

211: .seealso: DMProjectFieldLocal(), DMProjectFieldLabelLocal(), DMProjectFunction(), DMComputeL2Diff()
212: @*/
213: PetscErrorCode DMProjectField(DM dm, PetscReal time, Vec U,
214:                               void (**funcs)(PetscInt, PetscInt, PetscInt,
215:                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
216:                                              const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
217:                                              PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
218:                               InsertMode mode, Vec X)
219: {
220:   Vec            localX, localU;

225:   DMGetLocalVector(dm, &localX);
226:   /* We currently check whether locU == locX to see if we need to apply BC */
227:   if (U != X) {DMGetLocalVector(dm, &localU);}
228:   else        {localU = localX;}
229:   DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);
230:   DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);
231:   DMProjectFieldLocal(dm, time, localU, funcs, mode, localX);
232:   DMLocalToGlobalBegin(dm, localX, mode, X);
233:   DMLocalToGlobalEnd(dm, localX, mode, X);
234:   if (mode == INSERT_VALUES || mode == INSERT_ALL_VALUES || mode == INSERT_BC_VALUES) {
235:     Mat cMat;

237:     DMGetDefaultConstraints(dm, NULL, &cMat);
238:     if (cMat) {
239:       DMGlobalToLocalSolve(dm, localX, X);
240:     }
241:   }
242:   DMRestoreLocalVector(dm, &localX);
243:   if (U != X) {DMRestoreLocalVector(dm, &localU);}
244:   return(0);
245: }

247: /********************* Adaptive Interpolation **************************/

249: /* See the discussion of Adaptive Interpolation in manual/high_level_mg.rst */
250: PetscErrorCode DMAdaptInterpolator(DM dmc, DM dmf, Mat In, KSP smoother, PetscInt Nc, Vec vf[], Vec vc[], Mat *InAdapt, void *user)
251: {
252:   Mat            globalA;
253:   Vec            tmp, tmp2;
254:   PetscScalar   *A, *b, *x, *workscalar;
255:   PetscReal     *w, *sing, *workreal, rcond = PETSC_SMALL;
256:   PetscBLASInt   M, N, one = 1, irank, lwrk, info;
257:   PetscInt       debug = 0, rStart, rEnd, r, maxcols = 0, k;
258:   PetscBool      allocVc = PETSC_FALSE;

262:   PetscLogEventBegin(DM_AdaptInterpolator,dmc,dmf,0,0);
263:   PetscOptionsGetInt(NULL, NULL, "-dm_interpolator_adapt_debug", &debug, NULL);
264:   MatDuplicate(In, MAT_SHARE_NONZERO_PATTERN, InAdapt);
265:   MatGetOwnershipRange(In, &rStart, &rEnd);
266:   #if 0
267:   MatGetMaxRowLen(In, &maxcols);
268:   #else
269:   for (r = rStart; r < rEnd; ++r) {
270:     PetscInt           ncols;
271:     const PetscInt    *cols;
272:     const PetscScalar *vals;

274:     MatGetRow(In, r, &ncols, &cols, &vals);
275:     maxcols = PetscMax(maxcols, ncols);
276:     MatRestoreRow(In, r, &ncols, &cols, &vals);
277:   }
278:   #endif
279:   if (Nc < maxcols) PetscPrintf(PETSC_COMM_SELF, "The number of input vectors %D < %D the maximum number of column entries\n", Nc, maxcols);
280:   for (k = 0; k < Nc; ++k) {
281:     char        name[PETSC_MAX_PATH_LEN];
282:     const char *prefix;

284:     PetscObjectGetOptionsPrefix((PetscObject) smoother, &prefix);
285:     PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "%sCoarse Vector %D", prefix ? prefix : NULL, k);
286:     PetscObjectSetName((PetscObject) vc[k], name);
287:     VecViewFromOptions(vc[k], NULL, "-dm_adapt_interp_view_coarse");
288:     PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "%sFine Vector %D", prefix ? prefix : NULL, k);
289:     PetscObjectSetName((PetscObject) vf[k], name);
290:     VecViewFromOptions(vf[k], NULL, "-dm_adapt_interp_view_fine");
291:   }
292:   PetscBLASIntCast(3*PetscMin(Nc, maxcols) + PetscMax(2*PetscMin(Nc, maxcols), PetscMax(Nc, maxcols)), &lwrk);
293:   PetscMalloc7(Nc*maxcols, &A, PetscMax(Nc, maxcols), &b, Nc, &w, maxcols, &x, maxcols, &sing, lwrk, &workscalar, 5*PetscMin(Nc, maxcols), &workreal);
294:   /* w_k = \frac{\HC{v_k} B_l v_k}{\HC{v_k} A_l v_k} or the inverse Rayleigh quotient, which we calculate using \frac{\HC{v_k} v_k}{\HC{v_k} B^{-1}_l A_l v_k} */
295:   KSPGetOperators(smoother, &globalA, NULL);
296:   DMGetGlobalVector(dmf, &tmp);
297:   DMGetGlobalVector(dmf, &tmp2);
298:   for (k = 0; k < Nc; ++k) {
299:     PetscScalar vnorm, vAnorm;
300:     PetscBool   canMult = PETSC_FALSE;
301:     const char *type;

303:     w[k] = 1.0;
304:     PetscObjectGetType((PetscObject) globalA, &type);
305:     if (type) {MatAssembled(globalA, &canMult);}
306:     if (type && canMult) {
307:       VecDot(vf[k], vf[k], &vnorm);
308:       MatMult(globalA, vf[k], tmp);
309: #if 0
310:       KSPSolve(smoother, tmp, tmp2);
311:       VecDot(vf[k], tmp2, &vAnorm);
312: #else
313:       VecDot(vf[k], tmp, &vAnorm);
314: #endif
315:       w[k] = PetscRealPart(vnorm) / PetscRealPart(vAnorm);
316:     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "System matrix is not assembled.");
317:   }
318:   DMRestoreGlobalVector(dmf, &tmp);
319:   DMRestoreGlobalVector(dmf, &tmp2);
320:   if (!vc) {
321:     allocVc = PETSC_TRUE;
322:     PetscMalloc1(Nc, &vc);
323:     for (k = 0; k < Nc; ++k) {
324:       DMGetGlobalVector(dmc, &vc[k]);
325:       MatMultTranspose(In, vf[k], vc[k]);
326:     }
327:   }
328:   /* Solve a LS system for each fine row */
329:   for (r = rStart; r < rEnd; ++r) {
330:     PetscInt           ncols, c;
331:     const PetscInt    *cols;
332:     const PetscScalar *vals, *a;

334:     MatGetRow(In, r, &ncols, &cols, &vals);
335:     for (k = 0; k < Nc; ++k) {
336:       /* Need to fit lowest mode exactly */
337:       const PetscReal wk = ((ncols == 1) && (k > 0)) ? 0.0 : PetscSqrtReal(w[k]);

339:       /* b_k = \sqrt{w_k} f^{F,k}_r */
340:       VecGetArrayRead(vf[k], &a);
341:       b[k] = wk * a[r-rStart];
342:       VecRestoreArrayRead(vf[k], &a);
343:       /* A_{kc} = \sqrt{w_k} f^{C,k}_c */
344:       /* TODO Must pull out VecScatter from In, scatter in vc[k] values up front, and access them indirectly just as in MatMult() */
345:       VecGetArrayRead(vc[k], &a);
346:       for (c = 0; c < ncols; ++c) {
347:         /* This is element (k, c) of A */
348:         A[c*Nc+k] = wk * a[cols[c]-rStart];
349:       }
350:       VecRestoreArrayRead(vc[k], &a);
351:     }
352:     PetscBLASIntCast(Nc,    &M);
353:     PetscBLASIntCast(ncols, &N);
354:     if (debug) {
355: #if defined(PETSC_USE_COMPLEX)
356:       PetscScalar *tmp;
357:       PetscInt     j;

359:       DMGetWorkArray(dmc, Nc, MPIU_SCALAR, (void *) &tmp);
360:       for (j = 0; j < Nc; ++j) tmp[j] = w[j];
361:       DMPrintCellMatrix(r, "Interpolator Row LS weights", Nc, 1, tmp);
362:       DMPrintCellMatrix(r, "Interpolator Row LS matrix", Nc, ncols, A);
363:       for (j = 0; j < Nc; ++j) tmp[j] = b[j];
364:       DMPrintCellMatrix(r, "Interpolator Row LS rhs", Nc, 1, tmp);
365:       DMRestoreWorkArray(dmc, Nc, MPIU_SCALAR, (void *) &tmp);
366: #else
367:       DMPrintCellMatrix(r, "Interpolator Row LS weights", Nc, 1, w);
368:       DMPrintCellMatrix(r, "Interpolator Row LS matrix", Nc, ncols, A);
369:       DMPrintCellMatrix(r, "Interpolator Row LS rhs", Nc, 1, b);
370: #endif
371:     }
372: #if defined(PETSC_USE_COMPLEX)
373:     /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */
374:     PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&M, &N, &one, A, &M, b, M > N ? &M : &N, sing, &rcond, &irank, workscalar, &lwrk, workreal, &info));
375: #else
376:     /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */
377:     PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&M, &N, &one, A, &M, b, M > N ? &M : &N, sing, &rcond, &irank, workscalar, &lwrk, &info));
378: #endif
379:     if (info < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GELSS");
380:     if (info > 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "SVD failed to converge");
381:     if (debug) {
382:       PetscPrintf(PETSC_COMM_SELF, "rank %d rcond %g\n", irank, (double) rcond);
383: #if defined(PETSC_USE_COMPLEX)
384:       {
385:         PetscScalar *tmp;
386:         PetscInt     j;

388:         DMGetWorkArray(dmc, Nc, MPIU_SCALAR, (void *) &tmp);
389:         for (j = 0; j < PetscMin(Nc, ncols); ++j) tmp[j] = sing[j];
390:         DMPrintCellMatrix(r, "Interpolator Row LS singular values", PetscMin(Nc, ncols), 1, tmp);
391:         DMRestoreWorkArray(dmc, Nc, MPIU_SCALAR, (void *) &tmp);
392:       }
393: #else
394:       DMPrintCellMatrix(r, "Interpolator Row LS singular values", PetscMin(Nc, ncols), 1, sing);
395: #endif
396:       DMPrintCellMatrix(r, "Interpolator Row LS old P", ncols, 1, vals);
397:       DMPrintCellMatrix(r, "Interpolator Row LS sol", ncols, 1, b);
398:     }
399:     MatSetValues(*InAdapt, 1, &r, ncols, cols, b, INSERT_VALUES);
400:     MatRestoreRow(In, r, &ncols, &cols, &vals);
401:   }
402:   PetscFree7(A, b, w, x, sing, workscalar, workreal);
403:   if (allocVc) {
404:     for (k = 0; k < Nc; ++k) {DMRestoreGlobalVector(dmc, &vc[k]);}
405:     PetscFree(vc);
406:   }
407:   MatAssemblyBegin(*InAdapt, MAT_FINAL_ASSEMBLY);
408:   MatAssemblyEnd(*InAdapt, MAT_FINAL_ASSEMBLY);
409:   PetscLogEventEnd(DM_AdaptInterpolator,dmc,dmf,0,0);
410:   return(0);
411: }

413: PetscErrorCode DMCheckInterpolator(DM dmf, Mat In, PetscInt Nc, Vec vc[], Vec vf[], PetscReal tol)
414: {
415:   Vec            tmp;
416:   PetscReal      norminf, norm2, maxnorminf = 0.0, maxnorm2 = 0.0;
417:   PetscInt       k;

421:   DMGetGlobalVector(dmf, &tmp);
422:   MatViewFromOptions(In, NULL, "-dm_interpolator_adapt_error");
423:   for (k = 0; k < Nc; ++k) {
424:     MatMult(In, vc[k], tmp);
425:     VecAXPY(tmp, -1.0, vf[k]);
426:     VecViewFromOptions(vc[k], NULL, "-dm_interpolator_adapt_error");
427:     VecViewFromOptions(vf[k], NULL, "-dm_interpolator_adapt_error");
428:     VecViewFromOptions(tmp, NULL, "-dm_interpolator_adapt_error");
429:     VecNorm(tmp, NORM_INFINITY, &norminf);
430:     VecNorm(tmp, NORM_2, &norm2);
431:     maxnorminf = PetscMax(maxnorminf, norminf);
432:     maxnorm2   = PetscMax(maxnorm2,   norm2);
433:     PetscPrintf(PetscObjectComm((PetscObject) dmf), "Coarse vec %D ||vf - P vc||_\\infty %g, ||vf - P vc||_2 %g\n", k, norminf, norm2);
434:   }
435:   DMRestoreGlobalVector(dmf, &tmp);
436:   if (maxnorm2 > tol) SETERRQ2(PetscObjectComm((PetscObject) dmf), PETSC_ERR_ARG_WRONG, "max_k ||vf_k - P vc_k||_2 %g > tol %g", maxnorm2, tol);
437:   return(0);
438: }