xref: /petsc/src/dm/impls/plex/plexmetric.c (revision ddb78f6b74cf1ff5a4ff9836096f0e738d54f39e)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petscblaslapack.h>
3 
4 PetscErrorCode DMPlexP1FieldCreate_Private(DM dm, PetscInt f, PetscInt size, Vec *metric)
5 {
6   MPI_Comm       comm;
7   PetscErrorCode ierr;
8   PetscFE        fe;
9   PetscInt       dim;
10 
11   PetscFunctionBegin;
12 
13   /* Extract metadata from dm */
14   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
15   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
16 
17   /* Create a P1 field of the requested size */
18   ierr = PetscFECreateLagrange(comm, dim, size, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
19   ierr = DMSetField(dm, f, NULL, (PetscObject)fe);CHKERRQ(ierr);
20   ierr = DMCreateDS(dm);CHKERRQ(ierr);
21   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
22   ierr = DMCreateLocalVector(dm, metric);CHKERRQ(ierr);
23 
24   PetscFunctionReturn(0);
25 }
26 
27 /*
28   DMPlexMetricCreate - Create a Riemannian metric field
29 
30   Input parameters:
31 + dm     - The DM
32 - f      - The field number to use
33 
34   Output parameter:
35 . metric - The metric
36 
37   Level: beginner
38 
39   Note: It is assumed that the DM is comprised of simplices.
40 
41 .seealso: DMPlexMetricCreateUniform(), DMPlexMetricCreateIsotropic()
42 */
43 PetscErrorCode DMPlexMetricCreate(DM dm, PetscInt f, Vec *metric)
44 {
45   PetscErrorCode ierr;
46   PetscInt       coordDim, Nd;
47 
48   PetscFunctionBegin;
49   ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr);
50   Nd = coordDim*coordDim;
51   ierr = DMPlexP1FieldCreate_Private(dm, f, Nd, metric);CHKERRQ(ierr);
52   PetscFunctionReturn(0);
53 }
54 
55 typedef struct {
56   PetscReal scaling;  /* Scaling for uniform metric diagonal */
57 } DMPlexMetricUniformCtx;
58 
59 static PetscErrorCode diagonal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
60 {
61   DMPlexMetricUniformCtx *user = (DMPlexMetricUniformCtx*)ctx;
62   PetscInt                i, j;
63 
64   for (i = 0; i < dim; ++i) {
65     for (j = 0; j < dim; ++j) {
66       if (i == j) u[i+dim*j] = user->scaling;
67       else u[i+dim*j] = 0.0;
68     }
69   }
70   return 0;
71 }
72 
73 /*
74   DMPlexMetricCreateUniform - Construct a uniform isotropic metric
75 
76   Input parameters:
77 + dm     - The DM
78 . f      - The field number to use
79 - alpha  - Scaling parameter for the diagonal
80 
81   Output parameter:
82 . metric - The uniform metric
83 
84   Level: beginner
85 
86   Note: It is assumed that the DM is comprised of simplices.
87 
88 .seealso: DMPlexMetricCreate(), DMPlexMetricCreateIsotropic()
89 */
90 PetscErrorCode DMPlexMetricCreateUniform(DM dm, PetscInt f, PetscReal alpha, Vec *metric)
91 {
92   DMPlexMetricUniformCtx user;
93   PetscErrorCode       (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
94   PetscErrorCode         ierr;
95   void                  *ctxs[1];
96 
97   PetscFunctionBegin;
98   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
99   if (!alpha) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling is undefined");
100   if (alpha < 1.0e-30) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Uniform metric scaling %e should be positive", alpha);
101   else user.scaling = alpha;
102   funcs[0] = diagonal;
103   ctxs[0] = &user;
104   ierr = DMProjectFunctionLocal(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, *metric);CHKERRQ(ierr);
105   PetscFunctionReturn(0);
106 }
107 
108 /*
109   DMPlexMetricCreateIsotropic - Construct an isotropic metric from an error indicator
110 
111   Input parameters:
112 + dm        - The DM
113 . f         - The field number to use
114 - indicator - The error indicator
115 
116   Output parameter:
117 . metric    - The isotropic metric
118 
119   Level: beginner
120 
121   Notes:
122 
123   It is assumed that the DM is comprised of simplices.
124 
125   The indicator needs to be a scalar field defined at *vertices*.
126 
127 .seealso: DMPlexMetricCreate(), DMPlexMetricCreateUniform()
128 */
129 PetscErrorCode DMPlexMetricCreateIsotropic(DM dm, PetscInt f, Vec indicator, Vec *metric)
130 {
131   DM                 dmIndi;
132   PetscErrorCode     ierr;
133   PetscInt           dim, d, vStart, vEnd, v;
134   const PetscScalar *indi;
135   PetscScalar       *met;
136 
137   PetscFunctionBegin;
138   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
139   ierr = DMPlexMetricCreate(dm, f, metric);CHKERRQ(ierr);
140   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
141   ierr = VecGetArrayRead(indicator, &indi);CHKERRQ(ierr);
142   ierr = VecGetArrayWrite(*metric, &met);CHKERRQ(ierr);
143   ierr = VecGetDM(indicator, &dmIndi);CHKERRQ(ierr);
144   for (v = vStart; v < vEnd; ++v) {
145     PetscScalar *vindi, *vmet;
146     ierr = DMPlexPointLocalRead(dmIndi, v, indi, &vindi);CHKERRQ(ierr);
147     ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
148     for (d = 0; d < dim; ++d) vmet[d*(dim+1)] = vindi[0];
149   }
150   ierr = VecRestoreArrayWrite(*metric, &met);CHKERRQ(ierr);
151   ierr = VecRestoreArrayRead(indicator, &indi);CHKERRQ(ierr);
152   PetscFunctionReturn(0);
153 }
154 
155 static PetscErrorCode DMPlexMetricModify_Private(PetscInt dim, PetscReal h_min, PetscReal h_max, PetscReal a_max, PetscScalar Mp[])
156 {
157   PetscErrorCode ierr;
158   PetscInt       i, j, k;
159   PetscReal     *eigs, max_eig, l_min = 1.0/(h_max*h_max), l_max = 1.0/(h_min*h_min), la_min = 1.0/(a_max*a_max);
160   PetscScalar   *Mpos;
161 
162   PetscFunctionBegin;
163   ierr = PetscMalloc2(dim*dim, &Mpos, dim, &eigs);CHKERRQ(ierr);
164 
165   /* Symmetrize */
166   for (i = 0; i < dim; ++i) {
167     Mpos[i*dim+i] = Mp[i*dim+i];
168     for (j = i+1; j < dim; ++j) {
169       Mpos[i*dim+j] = 0.5*(Mp[i*dim+j] + Mp[j*dim+i]);
170       Mpos[j*dim+i] = Mpos[i*dim+j];
171     }
172   }
173 
174   /* Compute eigendecomposition */
175   {
176     PetscScalar  *work;
177     PetscBLASInt lwork;
178 
179     lwork = 5*dim;
180     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
181     {
182       PetscBLASInt lierr;
183       PetscBLASInt nb;
184 
185       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
186       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
187 #if defined(PETSC_USE_COMPLEX)
188       {
189         PetscReal *rwork;
190         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
191         PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,rwork,&lierr));
192         ierr = PetscFree(rwork);CHKERRQ(ierr);
193       }
194 #else
195       PetscStackCallBLAS("LAPACKsyev",LAPACKsyev_("V","U",&nb,Mpos,&nb,eigs,work,&lwork,&lierr));
196 #endif
197       if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
198       ierr = PetscFPTrapPop();CHKERRQ(ierr);
199     }
200     ierr = PetscFree(work);CHKERRQ(ierr);
201   }
202 
203   /* Reflect to positive orthant and enforce maximum and minimum size */
204   max_eig = 0.0;
205   for (i = 0; i < dim; ++i) {
206     eigs[i] = PetscMin(l_max, PetscMax(l_min, PetscAbsReal(eigs[i])));
207     max_eig = PetscMax(eigs[i], max_eig);
208   }
209 
210   /* Enforce maximum anisotropy */
211   for (i = 0; i < dim; ++i) {
212     if (a_max > 1.0) eigs[i] = PetscMax(eigs[i], max_eig*la_min);
213   }
214 
215   /* Reconstruct Hessian */
216   for (i = 0; i < dim; ++i) {
217     for (j = 0; j < dim; ++j) {
218       Mp[i*dim+j] = 0.0;
219       for (k = 0; k < dim; ++k) {
220         Mp[i*dim+j] += Mpos[k*dim+i] * eigs[k] * Mpos[k*dim+j];
221       }
222     }
223   }
224   ierr = PetscFree2(Mpos, eigs);CHKERRQ(ierr);
225 
226   PetscFunctionReturn(0);
227 }
228 
229 /*
230   DMPlexMetricEnforceSPD - Enforce symmetric positive-definiteness of a metric
231 
232   Input parameters:
233 + dm            - The DM
234 . restrictSizes - Should maximum/minimum metric magnitudes and anisotropy be enforced?
235 - metric        - The metric
236 
237   Output parameter:
238 . metric        - The metric
239 
240   Level: beginner
241 
242 .seealso: DMPlexMetricNormalize(), DMPlexMetricIntersection()
243 */
244 PetscErrorCode DMPlexMetricEnforceSPD(DM dm, PetscBool restrictSizes, Vec metric)
245 {
246   DMPlexMetricCtx *user;
247   PetscErrorCode   ierr;
248   PetscInt         dim, vStart, vEnd, v;
249   PetscScalar     *met;
250   PetscReal        h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0;
251 
252   PetscFunctionBegin;
253 
254   /* Extract metadata from dm */
255   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
256   ierr = DMGetApplicationContext(dm, (void**)&user);CHKERRQ(ierr);
257   if (restrictSizes) {
258     if (user->h_max > h_min) h_max = PetscMin(h_max, user->h_max);
259     if (user->h_min > 0.0) h_min = PetscMax(h_min, user->h_min);
260     if (user->a_max > 1.0) a_max = user->a_max;
261   }
262 
263   /* Enforce SPD */
264   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
265   ierr = VecGetArray(metric, &met);CHKERRQ(ierr);
266   for (v = vStart; v < vEnd; ++v) {
267     PetscScalar *vmet;
268     ierr = DMPlexPointLocalRef(dm, v, met, &vmet);CHKERRQ(ierr);
269     ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, vmet);CHKERRQ(ierr);
270   }
271   ierr = VecRestoreArray(metric, &met);CHKERRQ(ierr);
272 
273   PetscFunctionReturn(0);
274 }
275 
276 static void detMFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
277                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
278                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
279                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
280 {
281   const PetscScalar p = constants[0];
282   PetscReal         detH = 0.0;
283 
284   if      (dim == 2) DMPlex_Det2D_Scalar_Internal(&detH, u);
285   else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detH, u);
286   f0[0] = PetscPowReal(detH, p/(2.0*p + dim));
287 }
288 
289 /*
290   DMPlexMetricNormalize - Apply L-p normalization to a metric
291 
292   Input parameters:
293 + dm            - The DM
294 . metricIn      - The unnormalized metric
295 - restrictSizes - Should maximum/minimum metric magnitudes and anisotropy be enforced?
296 
297   Output parameter:
298 . metricOut     - The normalized metric
299 
300   Level: beginner
301 
302 .seealso: DMPlexMetricEnforceSPD(), DMPlexMetricIntersection()
303 */
304 PetscErrorCode DMPlexMetricNormalize(DM dm, Vec metricIn, PetscBool restrictSizes, Vec *metricOut)
305 {
306   DMPlexMetricCtx *user;
307   MPI_Comm         comm;
308   PetscDS          ds;
309   PetscErrorCode   ierr;
310   PetscInt         dim, Nd, vStart, vEnd, v, i;
311   PetscScalar     *met, integral, constants[1];
312   PetscReal        p, h_min = 1.0e-30, h_max = 1.0e+30, a_max = 0.0, factGlob, target;
313 
314   PetscFunctionBegin;
315 
316   /* Extract metadata from dm */
317   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
318   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
319   Nd = dim*dim;
320   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
321   ierr = DMGetApplicationContext(dm, (void**)&user);CHKERRQ(ierr);
322   if (restrictSizes && user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max;
323   if (PetscAbsReal(user->p) >= 1.0) p = user->p;
324   else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Metric normalization order %f should be greater than one.", user->p);
325   constants[0] = p;
326   if (user->targetComplexity > 0.0) target = user->targetComplexity;
327   else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Target metric complexity %f should be positive.", user->targetComplexity);
328 
329   /* Set up metric and ensure it is SPD */
330   ierr = DMPlexMetricCreate(dm, 0, metricOut);CHKERRQ(ierr);
331   ierr = VecCopy(metricIn, *metricOut);CHKERRQ(ierr);
332   ierr = DMPlexMetricEnforceSPD(dm, PETSC_FALSE, *metricOut);CHKERRQ(ierr);
333 
334   /* Compute global normalization factor */
335   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
336   ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
337   ierr = PetscDSSetObjective(ds, 0, detMFunc);CHKERRQ(ierr);
338   ierr = DMPlexComputeIntegralFEM(dm, *metricOut, &integral, NULL);CHKERRQ(ierr);
339   factGlob = PetscPowReal(target/PetscRealPart(integral), 2.0/dim);
340 
341   /* Apply local scaling */
342   a_max = 0.0;
343   if (restrictSizes) {
344     if (user->h_max > h_min) h_max = PetscMin(h_max, user->h_max);
345     if (user->h_min > 0.0) h_min = PetscMax(h_min, user->h_min);
346     if (!user->restrictAnisotropyFirst && user->a_max > 1.0) a_max = user->a_max;
347   }
348   ierr = VecGetArray(*metricOut, &met);CHKERRQ(ierr);
349   for (v = vStart; v < vEnd; ++v) {
350     PetscScalar       *Mp;
351     PetscReal          detM, fact;
352 
353     ierr = DMPlexPointLocalRef(dm, v, met, &Mp);CHKERRQ(ierr);
354     if      (dim == 2) DMPlex_Det2D_Scalar_Internal(&detM, Mp);
355     else if (dim == 3) DMPlex_Det3D_Scalar_Internal(&detM, Mp);
356     else SETERRQ1(comm, PETSC_ERR_SUP, "Dimension %d not supported", dim);
357     fact = factGlob * PetscPowReal(PetscAbsReal(detM), -1.0/(2*p+dim));
358     for (i = 0; i < Nd; ++i) Mp[i] *= fact;
359     if (restrictSizes) { ierr = DMPlexMetricModify_Private(dim, h_min, h_max, a_max, Mp);CHKERRQ(ierr); }
360   }
361   ierr = VecRestoreArray(*metricOut, &met);CHKERRQ(ierr);
362 
363   PetscFunctionReturn(0);
364 }
365 
366 /*
367   DMPlexMetricAverage - Compute the average of a list of metrics
368 
369   Input Parameter:
370 + dm         - The DM
371 . numMetrics - The number of metrics to be averaged
372 . weights    - Weights for the average
373 - metrics    - The metrics to be averaged
374 
375   Output Parameter:
376 . metricAvg  - The averaged metric
377 
378   Level: beginner
379 
380   Notes:
381   The weights should sum to unity.
382 
383   If weights are not provided then an unweighted average is used.
384 
385 .seealso: DMPlexMetricAverage2(), DMPlexMetricAverage3(), DMPlexMetricIntersection()
386 */
387 PetscErrorCode DMPlexMetricAverage(DM dm, PetscInt numMetrics, PetscReal weights[], Vec metrics[], Vec *metricAvg)
388 {
389   PetscBool      haveWeights = PETSC_TRUE;
390   PetscErrorCode ierr;
391   PetscInt       i;
392   PetscReal      sum = 0.0, tol = 1.0e-10;
393 
394   PetscFunctionBegin;
395   if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot average %d < 1 metrics", numMetrics); }
396   ierr = DMPlexMetricCreate(dm, 0, metricAvg);CHKERRQ(ierr);
397   ierr = VecSet(*metricAvg, 0.0);CHKERRQ(ierr);
398 
399   /* Default to the unweighted case */
400   if (!weights) {
401     ierr = PetscMalloc1(numMetrics, &weights);CHKERRQ(ierr);
402     haveWeights = PETSC_FALSE;
403     for (i = 0; i < numMetrics; ++i) {weights[i] = 1.0/numMetrics; }
404   }
405 
406   /* Check weights sum to unity */
407   for (i = 0; i < numMetrics; ++i) { sum += weights[i]; }
408   if (PetscAbsReal(sum - 1) > tol) { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Weights do not sum to unity"); }
409 
410   /* Compute metric average */
411   for (i = 0; i < numMetrics; ++i) { ierr = VecAXPY(*metricAvg, weights[i], metrics[i]);CHKERRQ(ierr); }
412   if (!haveWeights) {ierr = PetscFree(weights); }
413   PetscFunctionReturn(0);
414 }
415 
416 /*
417   DMPlexMetricAverage2 - Compute the unweighted average of two metrics
418 
419   Input Parameter:
420 + dm         - The DM
421 . metric1    - The first metric to be averaged
422 - metric2    - The second metric to be averaged
423 
424   Output Parameter:
425 . metricAvg  - The averaged metric
426 
427   Level: beginner
428 
429 .seealso: DMPlexMetricAverage(), DMPlexMetricAverage3()
430 */
431 PetscErrorCode DMPlexMetricAverage2(DM dm, Vec metric1, Vec metric2, Vec *metricAvg)
432 {
433   PetscErrorCode ierr;
434   PetscReal      weights[2] = {0.5, 0.5};
435   Vec            metrics[2] = {metric1, metric2};
436 
437   PetscFunctionBegin;
438   ierr = DMPlexMetricAverage(dm, 2, weights, metrics, metricAvg);CHKERRQ(ierr);
439   PetscFunctionReturn(0);
440 }
441 
442 /*
443   DMPlexMetricAverage3 - Compute the unweighted average of three metrics
444 
445   Input Parameter:
446 + dm         - The DM
447 . metric1    - The first metric to be averaged
448 . metric2    - The second metric to be averaged
449 - metric3    - The third metric to be averaged
450 
451   Output Parameter:
452 . metricAvg  - The averaged metric
453 
454   Level: beginner
455 
456 .seealso: DMPlexMetricAverage(), DMPlexMetricAverage2()
457 */
458 PetscErrorCode DMPlexMetricAverage3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricAvg)
459 {
460   PetscErrorCode ierr;
461   PetscReal      weights[3] = {1.0/3.0, 1.0/3.0, 1.0/3.0};
462   Vec            metrics[3] = {metric1, metric2, metric3};
463 
464   PetscFunctionBegin;
465   ierr = DMPlexMetricAverage(dm, 3, weights, metrics, metricAvg);CHKERRQ(ierr);
466   PetscFunctionReturn(0);
467 }
468 
469 static PetscErrorCode DMPlexMetricIntersection_Private(PetscInt dim, PetscScalar M1[], PetscScalar M2[])
470 {
471   PetscErrorCode ierr;
472   PetscInt       i, j, k, l, m;
473   PetscReal     *evals, *evals1;
474   PetscScalar   *evecs, *sqrtM1, *isqrtM1;
475 
476   PetscFunctionBegin;
477   ierr = PetscMalloc5(dim*dim, &evecs, dim*dim, &sqrtM1, dim*dim, &isqrtM1, dim, &evals, dim, &evals1);CHKERRQ(ierr);
478   for (i = 0; i < dim; ++i) {
479     for (j = 0; j < dim; ++j) {
480       evecs[i*dim+j] = M1[i*dim+j];
481     }
482   }
483   {
484     PetscScalar *work;
485     PetscBLASInt lwork;
486 
487     lwork = 5*dim;
488     ierr = PetscMalloc1(5*dim, &work);CHKERRQ(ierr);
489     {
490       PetscBLASInt lierr, nb;
491       PetscReal    sqrtk;
492 
493       /* Compute eigendecomposition of M1 */
494       ierr = PetscBLASIntCast(dim, &nb);CHKERRQ(ierr);
495       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
496 #if defined(PETSC_USE_COMPLEX)
497       {
498         PetscReal *rwork;
499         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
500         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, rwork, &lierr));
501         ierr = PetscFree(rwork);CHKERRQ(ierr);
502       }
503 #else
504       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals1, work, &lwork, &lierr));
505 #endif
506       if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
507       ierr = PetscFPTrapPop();
508 
509       /* Compute square root and reciprocal */
510       for (i = 0; i < dim; ++i) {
511         for (j = 0; j < dim; ++j) {
512           sqrtM1[i*dim+j] = 0.0;
513           isqrtM1[i*dim+j] = 0.0;
514           for (k = 0; k < dim; ++k) {
515             sqrtk = PetscSqrtReal(evals1[k]);
516             sqrtM1[i*dim+j] += evecs[k*dim+i] * sqrtk * evecs[k*dim+j];
517             isqrtM1[i*dim+j] += evecs[k*dim+i] * (1.0/sqrtk) * evecs[k*dim+j];
518           }
519         }
520       }
521 
522       /* Map into the space spanned by the eigenvectors of M1 */
523       for (i = 0; i < dim; ++i) {
524         for (j = 0; j < dim; ++j) {
525           evecs[i*dim+j] = 0.0;
526           for (k = 0; k < dim; ++k) {
527             for (l = 0; l < dim; ++l) {
528               evecs[i*dim+j] += isqrtM1[i*dim+k] * M2[l*dim+k] * isqrtM1[j*dim+l];
529             }
530           }
531         }
532       }
533 
534       /* Compute eigendecomposition */
535       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
536 #if defined(PETSC_USE_COMPLEX)
537       {
538         PetscReal *rwork;
539         ierr = PetscMalloc1(3*dim, &rwork);CHKERRQ(ierr);
540         PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, rwork, &lierr));
541         ierr = PetscFree(rwork);CHKERRQ(ierr);
542       }
543 #else
544       PetscStackCallBLAS("LAPACKsyev", LAPACKsyev_("V", "U", &nb, evecs, &nb, evals, work, &lwork, &lierr));
545 #endif
546       if (lierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in LAPACK routine %d", (int) lierr);
547       ierr = PetscFPTrapPop();
548 
549       /* Modify eigenvalues */
550       for (i = 0; i < dim; ++i) evals[i] = PetscMin(evals[i], evals1[i]);
551 
552       /* Map back to get the intersection */
553       for (i = 0; i < dim; ++i) {
554         for (j = 0; j < dim; ++j) {
555           M2[i*dim+j] = 0.0;
556           for (k = 0; k < dim; ++k) {
557             for (l = 0; l < dim; ++l) {
558               for (m = 0; m < dim; ++m) {
559                 M2[i*dim+j] += sqrtM1[i*dim+k] * evecs[l*dim+k] * evals[l] * evecs[l*dim+m] * sqrtM1[j*dim+m];
560               }
561             }
562           }
563         }
564       }
565     }
566     ierr = PetscFree(work);CHKERRQ(ierr);
567   }
568   ierr = PetscFree5(evecs, sqrtM1, isqrtM1, evals, evals1);CHKERRQ(ierr);
569   PetscFunctionReturn(0);
570 }
571 
572 /*
573   DMPlexMetricIntersection - Compute the intersection of a list of metrics
574 
575   Input Parameter:
576 + dm         - The DM
577 . numMetrics - The number of metrics to be intersected
578 - metrics    - The metrics to be intersected
579 
580   Output Parameter:
581 . metricInt  - The intersected metric
582 
583   Level: beginner
584 
585   Notes:
586 
587   The intersection of a list of metrics has the maximal ellipsoid which fits within the ellipsoids of the component metrics.
588 
589   The implementation used here is only consistent with the maximal ellipsoid definition in the case numMetrics = 2.
590 
591 .seealso: DMPlexMetricIntersection2(), DMPlexMetricIntersection3(), DMPlexMetricAverage()
592 */
593 PetscErrorCode DMPlexMetricIntersection(DM dm, PetscInt numMetrics, Vec metrics[], Vec *metricInt)
594 {
595   PetscErrorCode ierr;
596   PetscInt       dim, vStart, vEnd, v, i;
597   PetscScalar   *met, *meti, *M, *Mi;
598 
599   PetscFunctionBegin;
600   if (numMetrics < 1) { SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot intersect %d < 1 metrics", numMetrics); }
601 
602   /* Extract metadata from dm */
603   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
604   ierr = DMPlexMetricCreate(dm, 0, metricInt);CHKERRQ(ierr);
605   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
606 
607   /* Copy over the first metric */
608   ierr = VecCopy(metrics[0], *metricInt);CHKERRQ(ierr);
609 
610   /* Intersect subsequent metrics in turn */
611   if (numMetrics > 1) {
612     ierr = VecGetArray(*metricInt, &met);CHKERRQ(ierr);
613     for (i = 1; i < numMetrics; ++i) {
614       ierr = VecGetArray(metrics[i], &meti);CHKERRQ(ierr);
615       for (v = vStart; v < vEnd; ++v) {
616         ierr = DMPlexPointLocalRef(dm, v, met, &M);CHKERRQ(ierr);
617         ierr = DMPlexPointLocalRef(dm, v, meti, &Mi);CHKERRQ(ierr);
618         ierr = DMPlexMetricIntersection_Private(dim, Mi, M);CHKERRQ(ierr);
619       }
620       ierr = VecRestoreArray(metrics[i], &meti);CHKERRQ(ierr);
621     }
622     ierr = VecRestoreArray(*metricInt, &met);CHKERRQ(ierr);
623   }
624 
625   PetscFunctionReturn(0);
626 }
627 
628 /*
629   DMPlexMetricIntersection2 - Compute the intersection of two metrics
630 
631   Input Parameter:
632 + dm        - The DM
633 . metric1   - The first metric to be intersected
634 - metric2   - The second metric to be intersected
635 
636   Output Parameter:
637 . metricInt - The intersected metric
638 
639   Level: beginner
640 
641 .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection3()
642 */
643 PetscErrorCode DMPlexMetricIntersection2(DM dm, Vec metric1, Vec metric2, Vec *metricInt)
644 {
645   PetscErrorCode ierr;
646   Vec            metrics[2] = {metric1, metric2};
647 
648   PetscFunctionBegin;
649   ierr = DMPlexMetricIntersection(dm, 2, metrics, metricInt);CHKERRQ(ierr);
650   PetscFunctionReturn(0);
651 }
652 
653 /*
654   DMPlexMetricIntersection3 - Compute the intersection of three metrics
655 
656   Input Parameter:
657 + dm        - The DM
658 . metric1   - The first metric to be intersected
659 . metric2   - The second metric to be intersected
660 - metric3   - The third metric to be intersected
661 
662   Output Parameter:
663 . metricInt - The intersected metric
664 
665   Level: beginner
666 
667 .seealso: DMPlexMetricIntersection(), DMPlexMetricIntersection2()
668 */
669 PetscErrorCode DMPlexMetricIntersection3(DM dm, Vec metric1, Vec metric2, Vec metric3, Vec *metricInt)
670 {
671   PetscErrorCode ierr;
672   Vec            metrics[3] = {metric1, metric2, metric3};
673 
674   PetscFunctionBegin;
675   ierr = DMPlexMetricIntersection(dm, 3, metrics, metricInt);CHKERRQ(ierr);
676   PetscFunctionReturn(0);
677 }
678