100d473e3SJoe Wallwork static char help[] = "Test metric utils in the uniform, isotropic case.\n\n"; 200d473e3SJoe Wallwork 300d473e3SJoe Wallwork #include <petscdmplex.h> 400d473e3SJoe Wallwork 500d473e3SJoe Wallwork static PetscErrorCode bowl(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 600d473e3SJoe Wallwork { 700d473e3SJoe Wallwork PetscInt d; 800d473e3SJoe Wallwork 900d473e3SJoe Wallwork *u = 0.0; 1000d473e3SJoe Wallwork for (d = 0; d < dim; d++) *u += 0.5*(x[d] - 0.5)*(x[d] - 0.5); 1100d473e3SJoe Wallwork 1200d473e3SJoe Wallwork return 0; 1300d473e3SJoe Wallwork } 1400d473e3SJoe Wallwork 1500d473e3SJoe Wallwork int main(int argc, char **argv) { 1600d473e3SJoe Wallwork DM dm, dmDist, dmAdapt; 1700d473e3SJoe Wallwork DMLabel bdLabel = NULL; 1800d473e3SJoe Wallwork MPI_Comm comm; 19*e5697243SJoe Wallwork PetscBool uniform = PETSC_FALSE, isotropic = PETSC_FALSE; 2000d473e3SJoe Wallwork PetscErrorCode ierr; 2100d473e3SJoe Wallwork PetscInt *faces, dim = 3, numEdges = 4, d; 2200d473e3SJoe Wallwork PetscReal scaling = 1.0; 2300d473e3SJoe Wallwork Vec metric; 2400d473e3SJoe Wallwork 2500d473e3SJoe Wallwork /* Set up */ 2600d473e3SJoe Wallwork ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 2700d473e3SJoe Wallwork comm = PETSC_COMM_WORLD; 2800d473e3SJoe Wallwork ierr = PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");CHKERRQ(ierr); 2900d473e3SJoe Wallwork ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex60.c", dim, &dim, NULL, 2, 3);CHKERRQ(ierr); 3000d473e3SJoe Wallwork ierr = PetscOptionsInt("-num_edges", "Number of edges on each boundary of the initial mesh", "ex60.c", numEdges, &numEdges, NULL);CHKERRQ(ierr); 3100d473e3SJoe Wallwork ierr = PetscOptionsBool("-uniform", "Should the metric be assumed uniform?", "ex60.c", uniform, &uniform, NULL);CHKERRQ(ierr); 32*e5697243SJoe Wallwork ierr = PetscOptionsBool("-isotropic", "Should the metric be assumed isotropic, or computed as a recovered Hessian?", "ex60.c", isotropic, &isotropic, NULL);CHKERRQ(ierr); 3300d473e3SJoe Wallwork ierr = PetscOptionsEnd(); 3400d473e3SJoe Wallwork 3500d473e3SJoe Wallwork /* Create box mesh */ 3600d473e3SJoe Wallwork ierr = PetscMalloc1(dim, &faces);CHKERRQ(ierr); 3700d473e3SJoe Wallwork for (d = 0; d < dim; ++d) faces[d] = numEdges; 3800d473e3SJoe Wallwork ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr); 3900d473e3SJoe Wallwork ierr = PetscFree(faces);CHKERRQ(ierr); 40*e5697243SJoe Wallwork ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 4100d473e3SJoe Wallwork 4200d473e3SJoe Wallwork /* Distribute mesh over processes */ 4300d473e3SJoe Wallwork ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr); 4400d473e3SJoe Wallwork if (dmDist) { 4500d473e3SJoe Wallwork ierr = DMDestroy(&dm);CHKERRQ(ierr); 4600d473e3SJoe Wallwork dm = dmDist; 4700d473e3SJoe Wallwork } 4800d473e3SJoe Wallwork ierr = PetscObjectSetName((PetscObject) dm, "DM_init");CHKERRQ(ierr); 4900d473e3SJoe Wallwork ierr = DMViewFromOptions(dm, NULL, "-initial_mesh_view");CHKERRQ(ierr); 5000d473e3SJoe Wallwork 5100d473e3SJoe Wallwork /* Construct metric */ 5200d473e3SJoe Wallwork if (uniform) { 53*e5697243SJoe Wallwork if (isotropic) { ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric);CHKERRQ(ierr); } 5400d473e3SJoe Wallwork else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported."); 5500d473e3SJoe Wallwork } 5600d473e3SJoe Wallwork else { 5700d473e3SJoe Wallwork DM dmIndi; 5800d473e3SJoe Wallwork PetscFE fe; 5900d473e3SJoe Wallwork Vec indicator; 6000d473e3SJoe Wallwork 6100d473e3SJoe Wallwork /* Construct "error indicator" */ 6200d473e3SJoe Wallwork ierr = DMClone(dm, &dmIndi);CHKERRQ(ierr); 6300d473e3SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 6400d473e3SJoe Wallwork ierr = DMAddField(dmIndi, NULL, (PetscObject)fe);CHKERRQ(ierr); 6500d473e3SJoe Wallwork ierr = DMCreateDS(dmIndi);CHKERRQ(ierr); 6600d473e3SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 6700d473e3SJoe Wallwork ierr = DMCreateLocalVector(dmIndi, &indicator);CHKERRQ(ierr); 68*e5697243SJoe Wallwork if (isotropic) { 6900d473e3SJoe Wallwork 7000d473e3SJoe Wallwork /* Isotropic case: just specify unity */ 7100d473e3SJoe Wallwork ierr = VecSet(indicator, scaling);CHKERRQ(ierr); 7200d473e3SJoe Wallwork ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric);CHKERRQ(ierr); 7300d473e3SJoe Wallwork 7400d473e3SJoe Wallwork } else { 7500d473e3SJoe Wallwork 7600d473e3SJoe Wallwork /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */ 7700d473e3SJoe Wallwork DM dmGrad; 7800d473e3SJoe Wallwork PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl}; 7900d473e3SJoe Wallwork PetscInt pStart, pEnd, p; 8000d473e3SJoe Wallwork PetscSection sec, secCoord; 8100d473e3SJoe Wallwork Vec gradient; 8200d473e3SJoe Wallwork 8300d473e3SJoe Wallwork /* Project the parabola into P1 space */ 8400d473e3SJoe Wallwork ierr = DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator);CHKERRQ(ierr); 8500d473e3SJoe Wallwork 8600d473e3SJoe Wallwork /* Approximate the gradient */ 8700d473e3SJoe Wallwork ierr = DMClone(dmIndi, &dmGrad);CHKERRQ(ierr); 8800d473e3SJoe Wallwork ierr = DMGetCoordinateSection(dmGrad, &secCoord);CHKERRQ(ierr); 8900d473e3SJoe Wallwork ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dmGrad), &sec);CHKERRQ(ierr); 9000d473e3SJoe Wallwork ierr = PetscSectionSetNumFields(sec, 1);CHKERRQ(ierr); 9100d473e3SJoe Wallwork ierr = PetscSectionSetFieldComponents(sec, 0, dim);CHKERRQ(ierr); 9200d473e3SJoe Wallwork ierr = PetscSectionGetChart(secCoord, &pStart, &pEnd);CHKERRQ(ierr); 9300d473e3SJoe Wallwork ierr = PetscSectionSetChart(sec, pStart, pEnd);CHKERRQ(ierr); 9400d473e3SJoe Wallwork for (p = pStart; p < pEnd; ++p) { 9500d473e3SJoe Wallwork ierr = PetscSectionSetDof(sec, p, dim);CHKERRQ(ierr); 9600d473e3SJoe Wallwork ierr = PetscSectionSetFieldDof(sec, p, 0, dim);CHKERRQ(ierr); 9700d473e3SJoe Wallwork } 9800d473e3SJoe Wallwork ierr = PetscSectionSetUp(sec);CHKERRQ(ierr); 9900d473e3SJoe Wallwork ierr = DMSetLocalSection(dmGrad, sec);CHKERRQ(ierr); 10000d473e3SJoe Wallwork ierr = PetscSectionDestroy(&sec);CHKERRQ(ierr); 10100d473e3SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 10200d473e3SJoe Wallwork ierr = DMClearFields(dmGrad);CHKERRQ(ierr); 10300d473e3SJoe Wallwork ierr = DMAddField(dmGrad, NULL, (PetscObject)fe);CHKERRQ(ierr); 10400d473e3SJoe Wallwork ierr = DMCreateDS(dmGrad);CHKERRQ(ierr); 10500d473e3SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 10600d473e3SJoe Wallwork ierr = DMCreateLocalVector(dmGrad, &gradient);CHKERRQ(ierr); 10700d473e3SJoe Wallwork ierr = DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient);CHKERRQ(ierr); 10800d473e3SJoe Wallwork ierr = VecViewFromOptions(gradient, NULL, "-adapt_gradient_view");CHKERRQ(ierr); 10900d473e3SJoe Wallwork 11000d473e3SJoe Wallwork /* Approximate the Hessian */ 11100d473e3SJoe Wallwork ierr = DMGetCoordinateSection(dm, &secCoord);CHKERRQ(ierr); 11200d473e3SJoe Wallwork ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sec);CHKERRQ(ierr); 11300d473e3SJoe Wallwork ierr = PetscSectionSetNumFields(sec, 1);CHKERRQ(ierr); 11400d473e3SJoe Wallwork ierr = PetscSectionSetFieldComponents(sec, 0, dim*dim);CHKERRQ(ierr); 11500d473e3SJoe Wallwork ierr = PetscSectionGetChart(secCoord, &pStart, &pEnd);CHKERRQ(ierr); 11600d473e3SJoe Wallwork ierr = PetscSectionSetChart(sec, pStart, pEnd);CHKERRQ(ierr); 11700d473e3SJoe Wallwork for (p = pStart; p < pEnd; ++p) { 11800d473e3SJoe Wallwork ierr = PetscSectionSetDof(sec, p, dim*dim);CHKERRQ(ierr); 11900d473e3SJoe Wallwork ierr = PetscSectionSetFieldDof(sec, p, 0, dim*dim);CHKERRQ(ierr); 12000d473e3SJoe Wallwork } 12100d473e3SJoe Wallwork ierr = PetscSectionSetUp(sec);CHKERRQ(ierr); 12200d473e3SJoe Wallwork ierr = DMSetLocalSection(dm, sec);CHKERRQ(ierr); 12300d473e3SJoe Wallwork ierr = PetscSectionDestroy(&sec);CHKERRQ(ierr); 12400d473e3SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, dim*dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 12500d473e3SJoe Wallwork ierr = DMClearFields(dm);CHKERRQ(ierr); 12600d473e3SJoe Wallwork ierr = DMAddField(dm, NULL, (PetscObject)fe);CHKERRQ(ierr); 12700d473e3SJoe Wallwork ierr = DMCreateDS(dm);CHKERRQ(ierr); 12800d473e3SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 12900d473e3SJoe Wallwork ierr = DMCreateLocalVector(dm, &metric);CHKERRQ(ierr); 13000d473e3SJoe Wallwork ierr = DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric);CHKERRQ(ierr); 13100d473e3SJoe Wallwork ierr = VecViewFromOptions(metric, NULL, "-adapt_hessian_view");CHKERRQ(ierr); 13200d473e3SJoe Wallwork ierr = VecDestroy(&gradient);CHKERRQ(ierr); 13300d473e3SJoe Wallwork ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 13400d473e3SJoe Wallwork } 13500d473e3SJoe Wallwork ierr = VecDestroy(&indicator);CHKERRQ(ierr); 13600d473e3SJoe Wallwork ierr = DMDestroy(&dmIndi);CHKERRQ(ierr); 13700d473e3SJoe Wallwork } 13800d473e3SJoe Wallwork 13900d473e3SJoe Wallwork /* Test metric routines */ 14000d473e3SJoe Wallwork { 14100d473e3SJoe Wallwork PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2}; 14200d473e3SJoe Wallwork Vec metric1, metric2, metricComb; 14300d473e3SJoe Wallwork Vec metrics[2]; 14400d473e3SJoe Wallwork 14500d473e3SJoe Wallwork ierr = VecDuplicate(metric, &metric1);CHKERRQ(ierr); 14600d473e3SJoe Wallwork ierr = VecSet(metric1, 0);CHKERRQ(ierr); 14700d473e3SJoe Wallwork ierr = VecAXPY(metric1, 0.625, metric);CHKERRQ(ierr); 14800d473e3SJoe Wallwork ierr = VecDuplicate(metric, &metric2);CHKERRQ(ierr); 14900d473e3SJoe Wallwork ierr = VecSet(metric2, 0);CHKERRQ(ierr); 15000d473e3SJoe Wallwork ierr = VecAXPY(metric2, 2.5, metric);CHKERRQ(ierr); 15100d473e3SJoe Wallwork metrics[0] = metric1; 15200d473e3SJoe Wallwork metrics[1] = metric2; 15300d473e3SJoe Wallwork 15400d473e3SJoe Wallwork /* Test metric average */ 15500d473e3SJoe Wallwork ierr = DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb);CHKERRQ(ierr); 15600d473e3SJoe Wallwork ierr = VecAXPY(metricComb, -1, metric);CHKERRQ(ierr); 15700d473e3SJoe Wallwork ierr = VecNorm(metric, NORM_2, &norm);CHKERRQ(ierr); 15800d473e3SJoe Wallwork ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr); 15900d473e3SJoe Wallwork errornorm /= norm; 16000d473e3SJoe Wallwork if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed (L2 error %f)", errornorm); 16100d473e3SJoe Wallwork ierr = VecDestroy(&metricComb);CHKERRQ(ierr); 16200d473e3SJoe Wallwork 16300d473e3SJoe Wallwork /* Test metric intersection */ 164*e5697243SJoe Wallwork if (isotropic) { 16500d473e3SJoe Wallwork ierr = DMPlexMetricIntersection(dm, 2, metrics, &metricComb);CHKERRQ(ierr); 16600d473e3SJoe Wallwork ierr = VecAXPY(metricComb, -1, metric1);CHKERRQ(ierr); 16700d473e3SJoe Wallwork ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr); 16800d473e3SJoe Wallwork errornorm /= norm; 16900d473e3SJoe Wallwork if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed (L2 error %f)", errornorm); 17000d473e3SJoe Wallwork } 17100d473e3SJoe Wallwork ierr = VecDestroy(&metric2);CHKERRQ(ierr); 17200d473e3SJoe Wallwork ierr = VecDestroy(&metricComb);CHKERRQ(ierr); 17300d473e3SJoe Wallwork ierr = VecCopy(metric, metric1);CHKERRQ(ierr); 17400d473e3SJoe Wallwork 17500d473e3SJoe Wallwork /* Test metric SPD enforcement */ 176*e5697243SJoe Wallwork ierr = DMPlexMetricEnforceSPD(dm, PETSC_TRUE, PETSC_TRUE, metric);CHKERRQ(ierr); 177*e5697243SJoe Wallwork if (isotropic) { 17800d473e3SJoe Wallwork ierr = VecAXPY(metric1, -1, metric);CHKERRQ(ierr); 17900d473e3SJoe Wallwork ierr = VecNorm(metric1, NORM_2, &errornorm);CHKERRQ(ierr); 18000d473e3SJoe Wallwork errornorm /= norm; 18100d473e3SJoe Wallwork if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed (L2 error %f)", errornorm); 18200d473e3SJoe Wallwork } 18300d473e3SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 18400d473e3SJoe Wallwork 18500d473e3SJoe Wallwork /* Test metric normalization */ 186*e5697243SJoe Wallwork ierr = DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1);CHKERRQ(ierr); 187*e5697243SJoe Wallwork if (isotropic) { 188*e5697243SJoe Wallwork PetscReal target; 189*e5697243SJoe Wallwork 190*e5697243SJoe Wallwork ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr); 191*e5697243SJoe Wallwork scaling = PetscPowReal(target, 2.0/dim); 19200d473e3SJoe Wallwork ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric2);CHKERRQ(ierr); 19300d473e3SJoe Wallwork ierr = VecAXPY(metric2, -1, metric1);CHKERRQ(ierr); 19400d473e3SJoe Wallwork ierr = VecNorm(metric2, NORM_2, &errornorm);CHKERRQ(ierr); 19500d473e3SJoe Wallwork errornorm /= norm; 19600d473e3SJoe Wallwork if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed (L2 error %f)", errornorm); 19700d473e3SJoe Wallwork } 19800d473e3SJoe Wallwork ierr = VecCopy(metric1, metric);CHKERRQ(ierr); 19900d473e3SJoe Wallwork ierr = VecDestroy(&metric2);CHKERRQ(ierr); 20000d473e3SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 20100d473e3SJoe Wallwork } 20200d473e3SJoe Wallwork 20300d473e3SJoe Wallwork /* Adapt the mesh */ 20400d473e3SJoe Wallwork ierr = DMAdaptMetric(dm, metric, bdLabel, &dmAdapt);CHKERRQ(ierr); 20500d473e3SJoe Wallwork ierr = PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted");CHKERRQ(ierr); 20600d473e3SJoe Wallwork ierr = VecDestroy(&metric);CHKERRQ(ierr); 20700d473e3SJoe Wallwork ierr = DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view");CHKERRQ(ierr); 20800d473e3SJoe Wallwork 20900d473e3SJoe Wallwork /* Compare DMs */ 21000d473e3SJoe Wallwork ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 21100d473e3SJoe Wallwork ierr = DMView(dmAdapt, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 21200d473e3SJoe Wallwork 21300d473e3SJoe Wallwork /* Clean up */ 21400d473e3SJoe Wallwork ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr); 21500d473e3SJoe Wallwork ierr = DMDestroy(&dm);CHKERRQ(ierr); 21600d473e3SJoe Wallwork ierr = PetscFinalize(); 21700d473e3SJoe Wallwork return 0; 21800d473e3SJoe Wallwork } 21900d473e3SJoe Wallwork 22000d473e3SJoe Wallwork /*TEST 22100d473e3SJoe Wallwork 22200d473e3SJoe Wallwork build: 22300d473e3SJoe Wallwork requires: pragmatic 22400d473e3SJoe Wallwork 22500d473e3SJoe Wallwork test: 22600d473e3SJoe Wallwork suffix: uniform_2d 227*e5697243SJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -uniform -isotropic 22800d473e3SJoe Wallwork test: 22900d473e3SJoe Wallwork suffix: uniform_3d 230*e5697243SJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -uniform -isotropic 23100d473e3SJoe Wallwork test: 23200d473e3SJoe Wallwork suffix: iso_2d 233*e5697243SJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -isotropic 23400d473e3SJoe Wallwork test: 23500d473e3SJoe Wallwork suffix: iso_3d 236*e5697243SJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -isotropic 23700d473e3SJoe Wallwork test: 23800d473e3SJoe Wallwork suffix: hessian_2d 239*e5697243SJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 24000d473e3SJoe Wallwork test: 24100d473e3SJoe Wallwork suffix: hessian_3d 242*e5697243SJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 24300d473e3SJoe Wallwork 24400d473e3SJoe Wallwork TEST*/ 245