Actual source code: dmproject.c


  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:   DM  dm;
 10:   Vec mask;
 11: } projectConstraintsCtx;

 13: PetscErrorCode MatMult_GlobalToLocalNormal(Mat CtC, Vec x, Vec y)
 14: {
 15:   DM                     dm;
 16:   Vec                    local, mask;
 17:   projectConstraintsCtx *ctx;

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

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

 37:   for (f = 0; f < Nf; f++) u[f] = 1.;
 38:   return 0;
 39: }

 41: /*@
 42:   DMGlobalToLocalSolve - Solve for the global vector that is mapped to a given local vector by DMGlobalToLocalBegin()/DMGlobalToLocalEnd() with mode
 43:   = 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
 44:   a least-squares solution.  It is also assumed that DMLocalToGlobalBegin()/DMLocalToGlobalEnd() with mode = ADD_VALUES is the adjoint of the
 45:   global-to-local map, so that the least-squares solution may be found by the normal equations.

 47:   collective

 49:   Input Parameters:
 50: + dm - The DM object
 51: . x - The local vector
 52: - y - The global vector: the input value of globalVec is used as an initial guess

 54:   Output Parameters:
 55: . y - The least-squares solution

 57:   Level: advanced

 59:   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
 60:   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
 61:   closure of any local cell (see DMPlexGetAnchors()/DMPlexSetAnchors()).

 63: .seealso: `DMGlobalToLocalBegin()`, `DMGlobalToLocalEnd()`, `DMLocalToGlobalBegin()`, `DMLocalToGlobalEnd()`, `DMPlexGetAnchors()`, `DMPlexSetAnchors()`
 64: @*/
 65: PetscErrorCode DMGlobalToLocalSolve(DM dm, Vec x, Vec y)
 66: {
 67:   Mat                   CtC;
 68:   PetscInt              n, N, cStart, cEnd, c;
 69:   PetscBool             isPlex;
 70:   KSP                   ksp;
 71:   PC                    pc;
 72:   Vec                   global, mask = NULL;
 73:   projectConstraintsCtx ctx;

 75:   PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex);
 76:   if (isPlex) {
 77:     /* mark points in the closure */
 78:     DMCreateLocalVector(dm, &mask);
 79:     VecSet(mask, 0.0);
 80:     DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd);
 81:     if (cEnd > cStart) {
 82:       PetscScalar *ones;
 83:       PetscInt     numValues, i;

 85:       DMPlexVecGetClosure(dm, NULL, mask, cStart, &numValues, NULL);
 86:       PetscMalloc1(numValues, &ones);
 87:       for (i = 0; i < numValues; i++) ones[i] = 1.;
 88:       for (c = cStart; c < cEnd; c++) DMPlexVecSetClosure(dm, NULL, mask, c, ones, INSERT_VALUES);
 89:       PetscFree(ones);
 90:     }
 91:   } else {
 92:     PetscBool hasMask;

 94:     DMHasNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &hasMask);
 95:     if (!hasMask) {
 96:       PetscErrorCode (**func)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
 97:       void   **ctx;
 98:       PetscInt Nf, f;

100:       DMGetNumFields(dm, &Nf);
101:       PetscMalloc2(Nf, &func, Nf, &ctx);
102:       for (f = 0; f < Nf; ++f) {
103:         func[f] = DMGlobalToLocalSolve_project1;
104:         ctx[f]  = NULL;
105:       }
106:       DMGetNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
107:       DMProjectFunctionLocal(dm, 0.0, func, ctx, INSERT_ALL_VALUES, mask);
108:       DMRestoreNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
109:       PetscFree2(func, ctx);
110:     }
111:     DMGetNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
112:   }
113:   ctx.dm   = dm;
114:   ctx.mask = mask;
115:   VecGetSize(y, &N);
116:   VecGetLocalSize(y, &n);
117:   MatCreate(PetscObjectComm((PetscObject)dm), &CtC);
118:   MatSetSizes(CtC, n, n, N, N);
119:   MatSetType(CtC, MATSHELL);
120:   MatSetUp(CtC);
121:   MatShellSetContext(CtC, &ctx);
122:   MatShellSetOperation(CtC, MATOP_MULT, (void (*)(void))MatMult_GlobalToLocalNormal);
123:   KSPCreate(PetscObjectComm((PetscObject)dm), &ksp);
124:   KSPSetOperators(ksp, CtC, CtC);
125:   KSPSetType(ksp, KSPCG);
126:   KSPGetPC(ksp, &pc);
127:   PCSetType(pc, PCNONE);
128:   KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
129:   KSPSetUp(ksp);
130:   DMGetGlobalVector(dm, &global);
131:   VecSet(global, 0.);
132:   if (mask) VecPointwiseMult(x, mask, x);
133:   DMLocalToGlobalBegin(dm, x, ADD_VALUES, global);
134:   DMLocalToGlobalEnd(dm, x, ADD_VALUES, global);
135:   KSPSolve(ksp, global, y);
136:   DMRestoreGlobalVector(dm, &global);
137:   /* clean up */
138:   KSPDestroy(&ksp);
139:   MatDestroy(&CtC);
140:   if (isPlex) {
141:     VecDestroy(&mask);
142:   } else {
143:     DMRestoreNamedLocalVector(dm, "_DMGlobalToLocalSolve_mask", &mask);
144:   }

146:   return 0;
147: }

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

152:   Collective on DM

154:   Input Parameters:
155: + dm      - The DM
156: . time    - The time
157: . U       - The input field vector
158: . funcs   - The functions to evaluate, one per field
159: - mode    - The insertion mode for values

161:   Output Parameter:
162: . X       - The output vector

164:    Calling sequence of func:
165: $    func(PetscInt dim, PetscInt Nf, PetscInt NfAux,
166: $         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
167: $         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
168: $         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]);

170: +  dim          - The spatial dimension
171: .  Nf           - The number of input fields
172: .  NfAux        - The number of input auxiliary fields
173: .  uOff         - The offset of each field in u[]
174: .  uOff_x       - The offset of each field in u_x[]
175: .  u            - The field values at this point in space
176: .  u_t          - The field time derivative at this point in space (or NULL)
177: .  u_x          - The field derivatives at this point in space
178: .  aOff         - The offset of each auxiliary field in u[]
179: .  aOff_x       - The offset of each auxiliary field in u_x[]
180: .  a            - The auxiliary field values at this point in space
181: .  a_t          - The auxiliary field time derivative at this point in space (or NULL)
182: .  a_x          - The auxiliary field derivatives at this point in space
183: .  t            - The current time
184: .  x            - The coordinates of this point
185: .  numConstants - The number of constants
186: .  constants    - The value of each constant
187: -  f            - The value of the function at this point in space

189:   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.
190:   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
191:   a subdomain. You can also output a different number of fields than the input, with different discretizations. Last the auxiliary DM, attached to the
192:   auxiliary field vector, which is attached to dm, can also be different. It can have a different topology, number of fields, and discretizations.

194:   Level: intermediate

196: .seealso: `DMProjectFieldLocal()`, `DMProjectFieldLabelLocal()`, `DMProjectFunction()`, `DMComputeL2Diff()`
197: @*/
198: PetscErrorCode DMProjectField(DM dm, PetscReal time, Vec U, void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode mode, Vec X)
199: {
200:   Vec localX, localU;
201:   DM  dmIn;

204:   DMGetLocalVector(dm, &localX);
205:   /* We currently check whether locU == locX to see if we need to apply BC */
206:   if (U != X) {
207:     VecGetDM(U, &dmIn);
208:     DMGetLocalVector(dmIn, &localU);
209:   } else {
210:     dmIn   = dm;
211:     localU = localX;
212:   }
213:   DMGlobalToLocalBegin(dmIn, U, INSERT_VALUES, localU);
214:   DMGlobalToLocalEnd(dmIn, U, INSERT_VALUES, localU);
215:   DMProjectFieldLocal(dm, time, localU, funcs, mode, localX);
216:   DMLocalToGlobalBegin(dm, localX, mode, X);
217:   DMLocalToGlobalEnd(dm, localX, mode, X);
218:   if (mode == INSERT_VALUES || mode == INSERT_ALL_VALUES || mode == INSERT_BC_VALUES) {
219:     Mat cMat;

221:     DMGetDefaultConstraints(dm, NULL, &cMat, NULL);
222:     if (cMat) DMGlobalToLocalSolve(dm, localX, X);
223:   }
224:   DMRestoreLocalVector(dm, &localX);
225:   if (U != X) DMRestoreLocalVector(dmIn, &localU);
226:   return 0;
227: }

229: /********************* Adaptive Interpolation **************************/

231: /* See the discussion of Adaptive Interpolation in manual/high_level_mg.rst */
232: PetscErrorCode DMAdaptInterpolator(DM dmc, DM dmf, Mat In, KSP smoother, Mat MF, Mat MC, Mat *InAdapt, void *user)
233: {
234:   Mat                globalA, AF;
235:   Vec                tmp;
236:   const PetscScalar *af, *ac;
237:   PetscScalar       *A, *b, *x, *workscalar;
238:   PetscReal         *w, *sing, *workreal, rcond = PETSC_SMALL;
239:   PetscBLASInt       M, N, one = 1, irank, lwrk, info;
240:   PetscInt           debug = 0, rStart, rEnd, r, maxcols = 0, k, Nc, ldac, ldaf;
241:   PetscBool          allocVc = PETSC_FALSE;

243:   PetscLogEventBegin(DM_AdaptInterpolator, dmc, dmf, 0, 0);
244:   PetscOptionsGetInt(NULL, NULL, "-dm_interpolator_adapt_debug", &debug, NULL);
245:   MatGetSize(MF, NULL, &Nc);
246:   MatDuplicate(In, MAT_SHARE_NONZERO_PATTERN, InAdapt);
247:   MatGetOwnershipRange(In, &rStart, &rEnd);
248: #if 0
249:   MatGetMaxRowLen(In, &maxcols);
250: #else
251:   for (r = rStart; r < rEnd; ++r) {
252:     PetscInt ncols;

254:     MatGetRow(In, r, &ncols, NULL, NULL);
255:     maxcols = PetscMax(maxcols, ncols);
256:     MatRestoreRow(In, r, &ncols, NULL, NULL);
257:   }
258: #endif
259:   if (Nc < maxcols) PetscPrintf(PETSC_COMM_SELF, "The number of input vectors %" PetscInt_FMT " < %" PetscInt_FMT " the maximum number of column entries\n", Nc, maxcols);
260:   for (k = 0; k < Nc && debug; ++k) {
261:     char        name[PETSC_MAX_PATH_LEN];
262:     const char *prefix;
263:     Vec         vc, vf;

265:     PetscObjectGetOptionsPrefix((PetscObject)smoother, &prefix);

267:     if (MC) {
268:       PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "%sCoarse Vector %" PetscInt_FMT, prefix ? prefix : NULL, k);
269:       MatDenseGetColumnVecRead(MC, k, &vc);
270:       PetscObjectSetName((PetscObject)vc, name);
271:       VecViewFromOptions(vc, NULL, "-dm_adapt_interp_view_coarse");
272:       MatDenseRestoreColumnVecRead(MC, k, &vc);
273:     }
274:     PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "%sFine Vector %" PetscInt_FMT, prefix ? prefix : NULL, k);
275:     MatDenseGetColumnVecRead(MF, k, &vf);
276:     PetscObjectSetName((PetscObject)vf, name);
277:     VecViewFromOptions(vf, NULL, "-dm_adapt_interp_view_fine");
278:     MatDenseRestoreColumnVecRead(MF, k, &vf);
279:   }
280:   PetscBLASIntCast(3 * PetscMin(Nc, maxcols) + PetscMax(2 * PetscMin(Nc, maxcols), PetscMax(Nc, maxcols)), &lwrk);
281:   PetscMalloc7(Nc * maxcols, &A, PetscMax(Nc, maxcols), &b, Nc, &w, maxcols, &x, maxcols, &sing, lwrk, &workscalar, 5 * PetscMin(Nc, maxcols), &workreal);
282:   /* 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} */
283:   KSPGetOperators(smoother, &globalA, NULL);

285:   MatMatMult(globalA, MF, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &AF);
286:   for (k = 0; k < Nc; ++k) {
287:     PetscScalar vnorm, vAnorm;
288:     Vec         vf;

290:     w[k] = 1.0;
291:     MatDenseGetColumnVecRead(MF, k, &vf);
292:     MatDenseGetColumnVecRead(AF, k, &tmp);
293:     VecDot(vf, vf, &vnorm);
294: #if 0
295:     DMGetGlobalVector(dmf, &tmp2);
296:     KSPSolve(smoother, tmp, tmp2);
297:     VecDot(vf, tmp2, &vAnorm);
298:     DMRestoreGlobalVector(dmf, &tmp2);
299: #else
300:     VecDot(vf, tmp, &vAnorm);
301: #endif
302:     w[k] = PetscRealPart(vnorm) / PetscRealPart(vAnorm);
303:     MatDenseRestoreColumnVecRead(MF, k, &vf);
304:     MatDenseRestoreColumnVecRead(AF, k, &tmp);
305:   }
306:   MatDestroy(&AF);
307:   if (!MC) {
308:     allocVc = PETSC_TRUE;
309:     MatTransposeMatMult(In, MF, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &MC);
310:   }
311:   /* Solve a LS system for each fine row
312:      MATT: Can we generalize to the case where Nc for the fine space
313:      is different for Nc for the coarse? */
314:   MatDenseGetArrayRead(MF, &af);
315:   MatDenseGetLDA(MF, &ldaf);
316:   MatDenseGetArrayRead(MC, &ac);
317:   MatDenseGetLDA(MC, &ldac);
318:   for (r = rStart; r < rEnd; ++r) {
319:     PetscInt           ncols, c;
320:     const PetscInt    *cols;
321:     const PetscScalar *vals;

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

328:       /* b_k = \sqrt{w_k} f^{F,k}_r */
329:       b[k] = wk * af[r - rStart + k * ldaf];
330:       /* A_{kc} = \sqrt{w_k} f^{C,k}_c */
331:       /* TODO Must pull out VecScatter from In, scatter in vc[k] values up front, and access them indirectly just as in MatMult() */
332:       for (c = 0; c < ncols; ++c) {
333:         /* This is element (k, c) of A */
334:         A[c * Nc + k] = wk * ac[cols[c] - rStart + k * ldac];
335:       }
336:     }
337:     PetscBLASIntCast(Nc, &M);
338:     PetscBLASIntCast(ncols, &N);
339:     if (debug) {
340: #if defined(PETSC_USE_COMPLEX)
341:       PetscScalar *tmp;
342:       PetscInt     j;

344:       DMGetWorkArray(dmc, Nc, MPIU_SCALAR, (void *)&tmp);
345:       for (j = 0; j < Nc; ++j) tmp[j] = w[j];
346:       DMPrintCellMatrix(r, "Interpolator Row LS weights", Nc, 1, tmp);
347:       DMPrintCellMatrix(r, "Interpolator Row LS matrix", Nc, ncols, A);
348:       for (j = 0; j < Nc; ++j) tmp[j] = b[j];
349:       DMPrintCellMatrix(r, "Interpolator Row LS rhs", Nc, 1, tmp);
350:       DMRestoreWorkArray(dmc, Nc, MPIU_SCALAR, (void *)&tmp);
351: #else
352:       DMPrintCellMatrix(r, "Interpolator Row LS weights", Nc, 1, w);
353:       DMPrintCellMatrix(r, "Interpolator Row LS matrix", Nc, ncols, A);
354:       DMPrintCellMatrix(r, "Interpolator Row LS rhs", Nc, 1, b);
355: #endif
356:     }
357: #if defined(PETSC_USE_COMPLEX)
358:     /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */
359:     PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&M, &N, &one, A, &M, b, M > N ? &M : &N, sing, &rcond, &irank, workscalar, &lwrk, workreal, &info));
360: #else
361:     /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */
362:     PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&M, &N, &one, A, &M, b, M > N ? &M : &N, sing, &rcond, &irank, workscalar, &lwrk, &info));
363: #endif
366:     if (debug) {
367:       PetscPrintf(PETSC_COMM_SELF, "rank %" PetscBLASInt_FMT " rcond %g\n", irank, (double)rcond);
368: #if defined(PETSC_USE_COMPLEX)
369:       {
370:         PetscScalar *tmp;
371:         PetscInt     j;

373:         DMGetWorkArray(dmc, Nc, MPIU_SCALAR, (void *)&tmp);
374:         for (j = 0; j < PetscMin(Nc, ncols); ++j) tmp[j] = sing[j];
375:         DMPrintCellMatrix(r, "Interpolator Row LS singular values", PetscMin(Nc, ncols), 1, tmp);
376:         DMRestoreWorkArray(dmc, Nc, MPIU_SCALAR, (void *)&tmp);
377:       }
378: #else
379:       DMPrintCellMatrix(r, "Interpolator Row LS singular values", PetscMin(Nc, ncols), 1, sing);
380: #endif
381:       DMPrintCellMatrix(r, "Interpolator Row LS old P", ncols, 1, vals);
382:       DMPrintCellMatrix(r, "Interpolator Row LS sol", ncols, 1, b);
383:     }
384:     MatSetValues(*InAdapt, 1, &r, ncols, cols, b, INSERT_VALUES);
385:     MatRestoreRow(In, r, &ncols, &cols, &vals);
386:   }
387:   MatDenseRestoreArrayRead(MF, &af);
388:   MatDenseRestoreArrayRead(MC, &ac);
389:   PetscFree7(A, b, w, x, sing, workscalar, workreal);
390:   if (allocVc) MatDestroy(&MC);
391:   MatAssemblyBegin(*InAdapt, MAT_FINAL_ASSEMBLY);
392:   MatAssemblyEnd(*InAdapt, MAT_FINAL_ASSEMBLY);
393:   PetscLogEventEnd(DM_AdaptInterpolator, dmc, dmf, 0, 0);
394:   return 0;
395: }

397: PetscErrorCode DMCheckInterpolator(DM dmf, Mat In, Mat MC, Mat MF, PetscReal tol)
398: {
399:   Vec       tmp;
400:   PetscReal norminf, norm2, maxnorminf = 0.0, maxnorm2 = 0.0;
401:   PetscInt  k, Nc;

403:   DMGetGlobalVector(dmf, &tmp);
404:   MatViewFromOptions(In, NULL, "-dm_interpolator_adapt_error");
405:   MatGetSize(MF, NULL, &Nc);
406:   for (k = 0; k < Nc; ++k) {
407:     Vec vc, vf;

409:     MatDenseGetColumnVecRead(MC, k, &vc);
410:     MatDenseGetColumnVecRead(MF, k, &vf);
411:     MatMult(In, vc, tmp);
412:     VecAXPY(tmp, -1.0, vf);
413:     VecViewFromOptions(vc, NULL, "-dm_interpolator_adapt_error");
414:     VecViewFromOptions(vf, NULL, "-dm_interpolator_adapt_error");
415:     VecViewFromOptions(tmp, NULL, "-dm_interpolator_adapt_error");
416:     VecNorm(tmp, NORM_INFINITY, &norminf);
417:     VecNorm(tmp, NORM_2, &norm2);
418:     maxnorminf = PetscMax(maxnorminf, norminf);
419:     maxnorm2   = PetscMax(maxnorm2, norm2);
420:     PetscPrintf(PetscObjectComm((PetscObject)dmf), "Coarse vec %" PetscInt_FMT " ||vf - P vc||_\\infty %g, ||vf - P vc||_2 %g\n", k, (double)norminf, (double)norm2);
421:     MatDenseRestoreColumnVecRead(MC, k, &vc);
422:     MatDenseRestoreColumnVecRead(MF, k, &vf);
423:   }
424:   DMRestoreGlobalVector(dmf, &tmp);
426:   return 0;
427: }