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