xref: /petsc/src/snes/utils/convest.c (revision ebfe4b0dd6720ceca8680e82f5a07cba83b366b4)
1 #include <petscconvest.h>            /*I "petscconvest.h" I*/
2 #include <petscdmplex.h>
3 #include <petscds.h>
4 #include <petscblaslapack.h>
5 
6 #include <petsc/private/petscconvestimpl.h>
7 
8 static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
9 {
10   PetscInt c;
11   for (c = 0; c < Nc; ++c) u[c] = 0.0;
12   return 0;
13 }
14 
15 
16 /*@
17   PetscConvEstCreate - Create a PetscConvEst object
18 
19   Collective on MPI_Comm
20 
21   Input Parameter:
22 . comm - The communicator for the PetscConvEst object
23 
24   Output Parameter:
25 . ce   - The PetscConvEst object
26 
27   Level: beginner
28 
29 .keywords: PetscConvEst, convergence, create
30 .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate()
31 @*/
32 PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce)
33 {
34   PetscErrorCode ierr;
35 
36   PetscFunctionBegin;
37   PetscValidPointer(ce, 2);
38   ierr = PetscSysInitializePackage();CHKERRQ(ierr);
39   ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr);
40   (*ce)->monitor = PETSC_FALSE;
41   (*ce)->Nr      = 4;
42   PetscFunctionReturn(0);
43 }
44 
45 /*@
46   PetscConvEstDestroy - Destroys a PetscConvEst object
47 
48   Collective on PetscConvEst
49 
50   Input Parameter:
51 . ce - The PetscConvEst object
52 
53   Level: beginner
54 
55 .keywords: PetscConvEst, convergence, destroy
56 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
57 @*/
58 PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce)
59 {
60   PetscErrorCode ierr;
61 
62   PetscFunctionBegin;
63   if (!*ce) PetscFunctionReturn(0);
64   PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1);
65   if (--((PetscObject)(*ce))->refct > 0) {
66     *ce = NULL;
67     PetscFunctionReturn(0);
68   }
69   ierr = PetscFree2((*ce)->initGuess, (*ce)->exactSol);CHKERRQ(ierr);
70   ierr = PetscFree((*ce)->errors);CHKERRQ(ierr);
71   ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 /*@
76   PetscConvEstSetFromOptions - Sets a PetscConvEst object from options
77 
78   Collective on PetscConvEst
79 
80   Input Parameters:
81 . ce - The PetscConvEst object
82 
83   Level: beginner
84 
85 .keywords: PetscConvEst, convergence, options
86 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
87 @*/
88 PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce)
89 {
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr);
94   ierr = PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr);
95   ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr);
96   ierr = PetscOptionsEnd();
97   PetscFunctionReturn(0);
98 }
99 
100 /*@
101   PetscConvEstView - Views a PetscConvEst object
102 
103   Collective on PetscConvEst
104 
105   Input Parameters:
106 + ce     - The PetscConvEst object
107 - viewer - The PetscViewer object
108 
109   Level: beginner
110 
111 .keywords: PetscConvEst, convergence, view
112 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
113 @*/
114 PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer)
115 {
116   PetscErrorCode ierr;
117 
118   PetscFunctionBegin;
119   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr);
120   ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr);
121   PetscFunctionReturn(0);
122 }
123 
124 /*@
125   PetscConvEstGetSolver - Gets the solver used to produce discrete solutions
126 
127   Not collective
128 
129   Input Parameter:
130 . ce   - The PetscConvEst object
131 
132   Output Parameter:
133 . snes - The solver
134 
135   Level: intermediate
136 
137 .keywords: PetscConvEst, convergence
138 .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
139 @*/
140 PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, SNES *snes)
141 {
142   PetscFunctionBegin;
143   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
144   PetscValidPointer(snes, 2);
145   *snes = ce->snes;
146   PetscFunctionReturn(0);
147 }
148 
149 /*@
150   PetscConvEstSetSolver - Sets the solver used to produce discrete solutions
151 
152   Not collective
153 
154   Input Parameters:
155 + ce   - The PetscConvEst object
156 - snes - The solver
157 
158   Level: intermediate
159 
160   Note: The solver MUST have an attached DM/DS, so that we know the exact solution
161 
162 .keywords: PetscConvEst, convergence
163 .seealso: PetscConvEstGetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
164 @*/
165 PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, SNES snes)
166 {
167   PetscErrorCode ierr;
168 
169   PetscFunctionBegin;
170   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
171   PetscValidHeaderSpecific(snes, SNES_CLASSID, 2);
172   ce->snes = snes;
173   ierr = SNESGetDM(ce->snes, &ce->idm);CHKERRQ(ierr);
174   PetscFunctionReturn(0);
175 }
176 
177 /*@
178   PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence
179 
180   Collective on PetscConvEst
181 
182   Input Parameters:
183 . ce - The PetscConvEst object
184 
185   Level: beginner
186 
187 .keywords: PetscConvEst, convergence, setup
188 .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate()
189 @*/
190 PetscErrorCode PetscConvEstSetUp(PetscConvEst ce)
191 {
192   PetscDS        prob;
193   PetscInt       f;
194   PetscErrorCode ierr;
195 
196   PetscFunctionBegin;
197   ierr = DMGetDS(ce->idm, &prob);CHKERRQ(ierr);
198   ierr = PetscDSGetNumFields(prob, &ce->Nf);CHKERRQ(ierr);
199   ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr);
200   ierr = PetscMalloc2(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol);CHKERRQ(ierr);
201   for (f = 0; f < ce->Nf; ++f) ce->initGuess[f] = zero_private;
202   for (f = 0; f < ce->Nf; ++f) {
203     ierr = PetscDSGetExactSolution(prob, f, &ce->exactSol[f]);CHKERRQ(ierr);
204     if (!ce->exactSol[f]) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %D", f);
205   }
206   PetscFunctionReturn(0);
207 }
208 
209 static PetscErrorCode PetscConvEstLinearRegression_Private(PetscConvEst ce, PetscInt n, const PetscReal x[], const PetscReal y[], PetscReal *slope, PetscReal *intercept)
210 {
211   PetscScalar    H[4];
212   PetscReal     *X, *Y, beta[2];
213   PetscInt       i, j, k;
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   *slope = *intercept = 0.0;
218   ierr = PetscMalloc2(n*2, &X, n*2, &Y);CHKERRQ(ierr);
219   for (k = 0; k < n; ++k) {
220     /* X[n,2] = [1, x] */
221     X[k*2+0] = 1.0;
222     X[k*2+1] = x[k];
223   }
224   /* H = X^T X */
225   for (i = 0; i < 2; ++i) {
226     for (j = 0; j < 2; ++j) {
227       H[i*2+j] = 0.0;
228       for (k = 0; k < n; ++k) {
229         H[i*2+j] += X[k*2+i] * X[k*2+j];
230       }
231     }
232   }
233   /* H = (X^T X)^{-1} */
234   {
235     PetscBLASInt two = 2, ipiv[2], info;
236     PetscScalar  work[2];
237 
238     ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
239     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&two, &two, H, &two, ipiv, &info));
240     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&two, H, &two, ipiv, work, &two, &info));
241     ierr = PetscFPTrapPop();CHKERRQ(ierr);
242   }
243     /* Y = H X^T */
244   for (i = 0; i < 2; ++i) {
245     for (k = 0; k < n; ++k) {
246       Y[i*n+k] = 0.0;
247       for (j = 0; j < 2; ++j) {
248         Y[i*n+k] += PetscRealPart(H[i*2+j]) * X[k*2+j];
249       }
250     }
251   }
252   /* beta = Y error = [y-intercept, slope] */
253   for (i = 0; i < 2; ++i) {
254     beta[i] = 0.0;
255     for (k = 0; k < n; ++k) {
256       beta[i] += Y[i*n+k] * y[k];
257     }
258   }
259   ierr = PetscFree2(X, Y);CHKERRQ(ierr);
260   *intercept = beta[0];
261   *slope     = beta[1];
262   PetscFunctionReturn(0);
263 }
264 
265 /*@
266   PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization
267 
268   Not collective
269 
270   Input Parameter:
271 . ce   - The PetscConvEst object
272 
273   Output Parameter:
274 . alpha - The convergence rate for each field
275 
276   Note: The convergence rate alpha is defined by
277 $ || u_h - u_exact || < C h^alpha
278 where u_h is the discrete solution, and h is a measure of the discretization size.
279 
280 We solve a series of problems on refined meshes, calculate an error based upon the exact solution in the DS,
281 and then fit the result to our model above using linear regression.
282 
283   Options database keys:
284 . -snes_convergence_estimate : Execute convergence estimation and print out the rate
285 
286   Level: intermediate
287 
288 .keywords: PetscConvEst, convergence
289 .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate()
290 @*/
291 PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[])
292 {
293   DM            *dm;
294   PetscObject    disc;
295   MPI_Comm       comm;
296   const char    *uname, *dmname;
297   void          *ctx;
298   Vec            u;
299   PetscReal      t = 0.0, *x, *y, slope, intercept;
300   PetscInt      *dof, dim, Nr = ce->Nr, r, f;
301   PetscLogEvent  event;
302   PetscErrorCode ierr;
303 
304   PetscFunctionBegin;
305   ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr);
306   ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr);
307   ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr);
308   ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr);
309   ierr = PetscMalloc2((Nr+1), &dm, (Nr+1)*ce->Nf, &dof);CHKERRQ(ierr);
310   dm[0]  = ce->idm;
311   for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0;
312   /* Loop over meshes */
313   ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &event);CHKERRQ(ierr);
314   for (r = 0; r <= Nr; ++r) {
315     PetscLogStage stage;
316     char          stageName[PETSC_MAX_PATH_LEN];
317 
318     ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr);
319     ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);
320     ierr = PetscLogStagePush(stage);CHKERRQ(ierr);
321     if (r > 0) {
322       ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr);
323       ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr);
324       ierr = DMCopyDisc(ce->idm, dm[r]);CHKERRQ(ierr);
325       ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr);
326       ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr);
327       for (f = 0; f <= ce->Nf; ++f) {
328         PetscErrorCode (*nspconstr)(DM, PetscInt, MatNullSpace *);
329         ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr);
330         ierr = DMSetNullSpaceConstructor(dm[r],   f,  nspconstr);CHKERRQ(ierr);
331       }
332     }
333     ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr);
334     /* Create solution */
335     ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr);
336     ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr);
337     ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr);
338     ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr);
339     /* Setup solver */
340     ierr = SNESReset(ce->snes);CHKERRQ(ierr);
341     ierr = SNESSetDM(ce->snes, dm[r]);CHKERRQ(ierr);
342     ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr);
343     ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr);
344     /* Create initial guess */
345     ierr = DMProjectFunction(dm[r], t, ce->initGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
346     ierr = SNESSolve(ce->snes, NULL, u);CHKERRQ(ierr);
347     ierr = DMComputeL2FieldDiff(dm[r], t, ce->exactSol, NULL, u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr);
348     ierr = PetscLogEventEnd(event, ce, 0, 0, 0);CHKERRQ(ierr);
349     for (f = 0; f < ce->Nf; ++f) {
350       ierr = DMPlexGetHeightStratum(dm[r], 0, NULL, &dof[r*ce->Nf+f]);CHKERRQ(ierr);
351     }
352     /* Monitor */
353     if (ce->monitor) {
354       PetscReal *errors = &ce->errors[r*ce->Nf];
355 
356       ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr);
357       if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);}
358       for (f = 0; f < ce->Nf; ++f) {
359         if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);}
360         if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);}
361         else                     {ierr = PetscPrintf(comm, "%g", (double)errors[f]);CHKERRQ(ierr);}
362       }
363       if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);}
364       ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr);
365     }
366     /* Cleanup */
367     ierr = VecDestroy(&u);CHKERRQ(ierr);
368     ierr = PetscLogStagePop();CHKERRQ(ierr);
369   }
370   for (r = 1; r <= Nr; ++r) {
371     ierr = DMDestroy(&dm[r]);CHKERRQ(ierr);
372   }
373   /* Fit convergence rate */
374   ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr);
375   for (f = 0; f < ce->Nf; ++f) {
376     for (r = 0; r <= Nr; ++r) {
377       x[r] = PetscLog10Real(dof[r*ce->Nf+f]);
378       y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]);
379     }
380     ierr = PetscConvEstLinearRegression_Private(ce, Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr);
381     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
382     alpha[f] = -slope * dim;
383   }
384   ierr = PetscFree2(x, y);CHKERRQ(ierr);
385   ierr = PetscFree2(dm, dof);CHKERRQ(ierr);
386   /* Restore solver */
387   ierr = SNESReset(ce->snes);CHKERRQ(ierr);
388   ierr = SNESSetDM(ce->snes, ce->idm);CHKERRQ(ierr);
389   ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr);
390   ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr);
391   PetscFunctionReturn(0);
392 }
393 
394 /*@
395   PetscConvEstRateView - Displays the convergence rate to a viewer
396 
397    Collective on SNES
398 
399    Parameter:
400 +  snes - iterative context obtained from SNESCreate()
401 .  alpha - the convergence rate for each field
402 -  viewer - the viewer to display the reason
403 
404    Options Database Keys:
405 .  -snes_convergence_estimate - print the convergence rate
406 
407    Level: developer
408 
409 .seealso: PetscConvEstGetRate()
410 @*/
411 PetscErrorCode PetscConvEstRateView(PetscConvEst ce, PetscReal alpha[], PetscViewer viewer)
412 {
413   PetscBool      isAscii;
414   PetscInt       f;
415   PetscErrorCode ierr;
416 
417   PetscFunctionBegin;
418   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr);
419   if (isAscii) {
420     ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
421     ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr);
422     if (ce->Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);}
423     for (f = 0; f < ce->Nf; ++f) {
424       if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
425       ierr = PetscViewerASCIIPrintf(viewer, "%g", (double) alpha[f]);CHKERRQ(ierr);
426     }
427     if (ce->Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);}
428     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
429     ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr);
430   }
431   PetscFunctionReturn(0);
432 }
433