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; 1721b3e102SJoe Wallwork DMLabel bdLabel = NULL, rgLabel = NULL; 1800d473e3SJoe Wallwork MPI_Comm comm; 19e5697243SJoe Wallwork PetscBool uniform = PETSC_FALSE, isotropic = PETSC_FALSE; 2000d473e3SJoe Wallwork PetscErrorCode ierr; 21f49e87b0SJoe Wallwork PetscInt *faces, dim = 3, 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 = PetscOptionsBool("-uniform", "Should the metric be assumed uniform?", "ex60.c", uniform, &uniform, NULL);CHKERRQ(ierr); 31e5697243SJoe Wallwork ierr = PetscOptionsBool("-isotropic", "Should the metric be assumed isotropic, or computed as a recovered Hessian?", "ex60.c", isotropic, &isotropic, NULL);CHKERRQ(ierr); 3200d473e3SJoe Wallwork ierr = PetscOptionsEnd(); 3300d473e3SJoe Wallwork 3400d473e3SJoe Wallwork /* Create box mesh */ 3500d473e3SJoe Wallwork ierr = PetscMalloc1(dim, &faces);CHKERRQ(ierr); 36f49e87b0SJoe Wallwork for (d = 0; d < dim; ++d) faces[d] = 4; 3700d473e3SJoe Wallwork ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr); 3800d473e3SJoe Wallwork ierr = PetscFree(faces);CHKERRQ(ierr); 39e5697243SJoe Wallwork ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 4000d473e3SJoe Wallwork 4100d473e3SJoe Wallwork /* Distribute mesh over processes */ 4200d473e3SJoe Wallwork ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr); 4300d473e3SJoe Wallwork if (dmDist) { 4400d473e3SJoe Wallwork ierr = DMDestroy(&dm);CHKERRQ(ierr); 4500d473e3SJoe Wallwork dm = dmDist; 4600d473e3SJoe Wallwork } 4700d473e3SJoe Wallwork ierr = PetscObjectSetName((PetscObject) dm, "DM_init");CHKERRQ(ierr); 4800d473e3SJoe Wallwork ierr = DMViewFromOptions(dm, NULL, "-initial_mesh_view");CHKERRQ(ierr); 4900d473e3SJoe Wallwork 5000d473e3SJoe Wallwork /* Construct metric */ 5100d473e3SJoe Wallwork if (uniform) { 52e5697243SJoe Wallwork if (isotropic) { ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric);CHKERRQ(ierr); } 5300d473e3SJoe Wallwork else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported."); 5400d473e3SJoe Wallwork } 5500d473e3SJoe Wallwork else { 5600d473e3SJoe Wallwork DM dmIndi; 5700d473e3SJoe Wallwork PetscFE fe; 5800d473e3SJoe Wallwork Vec indicator; 5900d473e3SJoe Wallwork 6000d473e3SJoe Wallwork /* Construct "error indicator" */ 6100d473e3SJoe Wallwork ierr = DMClone(dm, &dmIndi);CHKERRQ(ierr); 6200d473e3SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 635767498bSJoe Wallwork ierr = DMSetField(dmIndi, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 6400d473e3SJoe Wallwork ierr = DMCreateDS(dmIndi);CHKERRQ(ierr); 6500d473e3SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 6600d473e3SJoe Wallwork ierr = DMCreateLocalVector(dmIndi, &indicator);CHKERRQ(ierr); 67e5697243SJoe Wallwork if (isotropic) { 6800d473e3SJoe Wallwork 6900d473e3SJoe Wallwork /* Isotropic case: just specify unity */ 7000d473e3SJoe Wallwork ierr = VecSet(indicator, scaling);CHKERRQ(ierr); 7100d473e3SJoe Wallwork ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric);CHKERRQ(ierr); 7200d473e3SJoe Wallwork 7300d473e3SJoe Wallwork } else { 7400d473e3SJoe Wallwork 7500d473e3SJoe Wallwork /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */ 7600d473e3SJoe Wallwork DM dmGrad; 7700d473e3SJoe Wallwork PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl}; 7800d473e3SJoe Wallwork Vec gradient; 7900d473e3SJoe Wallwork 8000d473e3SJoe Wallwork /* Project the parabola into P1 space */ 8100d473e3SJoe Wallwork ierr = DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator);CHKERRQ(ierr); 8200d473e3SJoe Wallwork 8300d473e3SJoe Wallwork /* Approximate the gradient */ 8400d473e3SJoe Wallwork ierr = DMClone(dmIndi, &dmGrad);CHKERRQ(ierr); 8500d473e3SJoe Wallwork ierr = PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 865767498bSJoe Wallwork ierr = DMSetField(dmGrad, 0, NULL, (PetscObject)fe);CHKERRQ(ierr); 8700d473e3SJoe Wallwork ierr = DMCreateDS(dmGrad);CHKERRQ(ierr); 8800d473e3SJoe Wallwork ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 8900d473e3SJoe Wallwork ierr = DMCreateLocalVector(dmGrad, &gradient);CHKERRQ(ierr); 9000d473e3SJoe Wallwork ierr = DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient);CHKERRQ(ierr); 9100d473e3SJoe Wallwork ierr = VecViewFromOptions(gradient, NULL, "-adapt_gradient_view");CHKERRQ(ierr); 9200d473e3SJoe Wallwork 9300d473e3SJoe Wallwork /* Approximate the Hessian */ 945767498bSJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, &metric);CHKERRQ(ierr); 9500d473e3SJoe Wallwork ierr = DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric);CHKERRQ(ierr); 9600d473e3SJoe Wallwork ierr = VecViewFromOptions(metric, NULL, "-adapt_hessian_view");CHKERRQ(ierr); 9700d473e3SJoe Wallwork ierr = VecDestroy(&gradient);CHKERRQ(ierr); 9800d473e3SJoe Wallwork ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 9900d473e3SJoe Wallwork } 10000d473e3SJoe Wallwork ierr = VecDestroy(&indicator);CHKERRQ(ierr); 10100d473e3SJoe Wallwork ierr = DMDestroy(&dmIndi);CHKERRQ(ierr); 10200d473e3SJoe Wallwork } 10300d473e3SJoe Wallwork 10400d473e3SJoe Wallwork /* Test metric routines */ 10500d473e3SJoe Wallwork { 106f49e87b0SJoe Wallwork DM dmDet; 10700d473e3SJoe Wallwork PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2}; 108f49e87b0SJoe Wallwork Vec metric1, metric2, metricComb, determinant; 10900d473e3SJoe Wallwork Vec metrics[2]; 11000d473e3SJoe Wallwork 11100d473e3SJoe Wallwork ierr = VecDuplicate(metric, &metric1);CHKERRQ(ierr); 11200d473e3SJoe Wallwork ierr = VecSet(metric1, 0);CHKERRQ(ierr); 11300d473e3SJoe Wallwork ierr = VecAXPY(metric1, 0.625, metric);CHKERRQ(ierr); 11400d473e3SJoe Wallwork ierr = VecDuplicate(metric, &metric2);CHKERRQ(ierr); 11500d473e3SJoe Wallwork ierr = VecSet(metric2, 0);CHKERRQ(ierr); 11600d473e3SJoe Wallwork ierr = VecAXPY(metric2, 2.5, metric);CHKERRQ(ierr); 11700d473e3SJoe Wallwork metrics[0] = metric1; 11800d473e3SJoe Wallwork metrics[1] = metric2; 11900d473e3SJoe Wallwork 12000d473e3SJoe Wallwork /* Test metric average */ 12100d473e3SJoe Wallwork ierr = DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb);CHKERRQ(ierr); 12200d473e3SJoe Wallwork ierr = VecAXPY(metricComb, -1, metric);CHKERRQ(ierr); 12300d473e3SJoe Wallwork ierr = VecNorm(metric, NORM_2, &norm);CHKERRQ(ierr); 12400d473e3SJoe Wallwork ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr); 12500d473e3SJoe Wallwork errornorm /= norm; 126*7cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 127*7cdd5544SJoe Wallwork if (errornorm > tol) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed"); 12800d473e3SJoe Wallwork ierr = VecDestroy(&metricComb);CHKERRQ(ierr); 12900d473e3SJoe Wallwork 13000d473e3SJoe Wallwork /* Test metric intersection */ 131e5697243SJoe Wallwork if (isotropic) { 13200d473e3SJoe Wallwork ierr = DMPlexMetricIntersection(dm, 2, metrics, &metricComb);CHKERRQ(ierr); 13300d473e3SJoe Wallwork ierr = VecAXPY(metricComb, -1, metric1);CHKERRQ(ierr); 13400d473e3SJoe Wallwork ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr); 13500d473e3SJoe Wallwork errornorm /= norm; 136*7cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 137*7cdd5544SJoe Wallwork if (errornorm > tol) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed"); 13800d473e3SJoe Wallwork } 139f49e87b0SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 14000d473e3SJoe Wallwork ierr = VecDestroy(&metric2);CHKERRQ(ierr); 14100d473e3SJoe Wallwork ierr = VecDestroy(&metricComb);CHKERRQ(ierr); 14200d473e3SJoe Wallwork 14300d473e3SJoe Wallwork /* Test metric SPD enforcement */ 144f49e87b0SJoe Wallwork ierr = DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1, &determinant);CHKERRQ(ierr); 145e5697243SJoe Wallwork if (isotropic) { 146f49e87b0SJoe Wallwork Vec err; 147f49e87b0SJoe Wallwork 148f49e87b0SJoe Wallwork ierr = VecDuplicate(determinant, &err);CHKERRQ(ierr); 149f49e87b0SJoe Wallwork ierr = VecSet(err, 1.0);CHKERRQ(ierr); 150f49e87b0SJoe Wallwork ierr = VecNorm(err, NORM_2, &norm);CHKERRQ(ierr); 151f49e87b0SJoe Wallwork ierr = VecAXPY(err, -1, determinant);CHKERRQ(ierr); 152f49e87b0SJoe Wallwork ierr = VecNorm(err, NORM_2, &errornorm);CHKERRQ(ierr); 153f49e87b0SJoe Wallwork ierr = VecDestroy(&err);CHKERRQ(ierr); 154f49e87b0SJoe Wallwork errornorm /= norm; 155*7cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 156*7cdd5544SJoe Wallwork if (errornorm > tol) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit"); 15700d473e3SJoe Wallwork ierr = VecAXPY(metric1, -1, metric);CHKERRQ(ierr); 15800d473e3SJoe Wallwork ierr = VecNorm(metric1, NORM_2, &errornorm);CHKERRQ(ierr); 15900d473e3SJoe Wallwork errornorm /= norm; 160*7cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 161*7cdd5544SJoe Wallwork if (errornorm > tol) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed"); 16200d473e3SJoe Wallwork } 16300d473e3SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 164f49e87b0SJoe Wallwork ierr = VecGetDM(determinant, &dmDet);CHKERRQ(ierr); 165f49e87b0SJoe Wallwork ierr = VecDestroy(&determinant);CHKERRQ(ierr); 166f49e87b0SJoe Wallwork ierr = DMDestroy(&dmDet);CHKERRQ(ierr); 16700d473e3SJoe Wallwork 16800d473e3SJoe Wallwork /* Test metric normalization */ 169e5697243SJoe Wallwork ierr = DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1);CHKERRQ(ierr); 170e5697243SJoe Wallwork if (isotropic) { 171e5697243SJoe Wallwork PetscReal target; 172e5697243SJoe Wallwork 173e5697243SJoe Wallwork ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr); 174e5697243SJoe Wallwork scaling = PetscPowReal(target, 2.0/dim); 17500d473e3SJoe Wallwork ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric2);CHKERRQ(ierr); 17600d473e3SJoe Wallwork ierr = VecAXPY(metric2, -1, metric1);CHKERRQ(ierr); 17700d473e3SJoe Wallwork ierr = VecNorm(metric2, NORM_2, &errornorm);CHKERRQ(ierr); 17800d473e3SJoe Wallwork errornorm /= norm; 179*7cdd5544SJoe Wallwork ierr = PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", 100*errornorm);CHKERRQ(ierr); 180*7cdd5544SJoe Wallwork if (errornorm > tol) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed"); 18100d473e3SJoe Wallwork } 18200d473e3SJoe Wallwork ierr = VecCopy(metric1, metric);CHKERRQ(ierr); 18300d473e3SJoe Wallwork ierr = VecDestroy(&metric2);CHKERRQ(ierr); 18400d473e3SJoe Wallwork ierr = VecDestroy(&metric1);CHKERRQ(ierr); 18500d473e3SJoe Wallwork } 18600d473e3SJoe Wallwork 18700d473e3SJoe Wallwork /* Adapt the mesh */ 18821b3e102SJoe Wallwork ierr = DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt);CHKERRQ(ierr); 18921b3e102SJoe Wallwork ierr = DMDestroy(&dm);CHKERRQ(ierr); 19000d473e3SJoe Wallwork ierr = PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted");CHKERRQ(ierr); 19100d473e3SJoe Wallwork ierr = VecDestroy(&metric);CHKERRQ(ierr); 19200d473e3SJoe Wallwork ierr = DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view");CHKERRQ(ierr); 19300d473e3SJoe Wallwork 19400d473e3SJoe Wallwork /* Clean up */ 19500d473e3SJoe Wallwork ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr); 19600d473e3SJoe Wallwork ierr = PetscFinalize(); 19700d473e3SJoe Wallwork return 0; 19800d473e3SJoe Wallwork } 19900d473e3SJoe Wallwork 20000d473e3SJoe Wallwork /*TEST 20100d473e3SJoe Wallwork 202dd0671eeSJoe Wallwork test: 203dd0671eeSJoe Wallwork suffix: uniform_2d_pragmatic 20400d473e3SJoe Wallwork requires: pragmatic 205dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor pragmatic -uniform -isotropic 20600d473e3SJoe Wallwork test: 207dd0671eeSJoe Wallwork suffix: uniform_3d_pragmatic 208dd0671eeSJoe Wallwork requires: pragmatic tetgen 209dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor pragmatic -uniform -isotropic 21000d473e3SJoe Wallwork test: 211dd0671eeSJoe Wallwork suffix: iso_2d_pragmatic 212dd0671eeSJoe Wallwork requires: pragmatic 213dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor pragmatic -isotropic 21400d473e3SJoe Wallwork test: 215dd0671eeSJoe Wallwork suffix: iso_3d_pragmatic 216dd0671eeSJoe Wallwork requires: pragmatic tetgen 217dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor pragmatic -isotropic 21800d473e3SJoe Wallwork test: 219dd0671eeSJoe Wallwork suffix: hessian_2d_pragmatic 220dd0671eeSJoe Wallwork requires: pragmatic 221dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor pragmatic 22200d473e3SJoe Wallwork test: 223dd0671eeSJoe Wallwork suffix: hessian_3d_pragmatic 224dd0671eeSJoe Wallwork requires: pragmatic tetgen 225dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor pragmatic 22600d473e3SJoe Wallwork test: 227dd0671eeSJoe Wallwork suffix: uniform_2d_mmg 228dd0671eeSJoe Wallwork requires: mmg 229dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor mmg -uniform -isotropic 230dd0671eeSJoe Wallwork test: 231dd0671eeSJoe Wallwork suffix: uniform_3d_mmg 232dd0671eeSJoe Wallwork requires: mmg tetgen 233dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor mmg -uniform -isotropic 234dd0671eeSJoe Wallwork test: 235dd0671eeSJoe Wallwork suffix: iso_2d_mmg 236dd0671eeSJoe Wallwork requires: mmg 237dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor mmg -isotropic 238dd0671eeSJoe Wallwork test: 239dd0671eeSJoe Wallwork suffix: iso_3d_mmg 240dd0671eeSJoe Wallwork requires: mmg tetgen 241dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor mmg -isotropic 242dd0671eeSJoe Wallwork test: 243dd0671eeSJoe Wallwork suffix: hessian_2d_mmg 244dd0671eeSJoe Wallwork requires: mmg 245dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor mmg 246dd0671eeSJoe Wallwork test: 247dd0671eeSJoe Wallwork suffix: hessian_3d_mmg 248dd0671eeSJoe Wallwork requires: mmg tetgen 249dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor mmg 250dd0671eeSJoe Wallwork test: 251dd0671eeSJoe Wallwork suffix: uniform_3d_parmmg 252dd0671eeSJoe Wallwork requires: parmmg tetgen 253dd0671eeSJoe Wallwork nsize: 2 254dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor parmmg -uniform -isotropic 255dd0671eeSJoe Wallwork test: 256dd0671eeSJoe Wallwork suffix: iso_3d_parmmg 257dd0671eeSJoe Wallwork requires: parmmg tetgen 258dd0671eeSJoe Wallwork nsize: 2 259dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor parmmg -isotropic 260dd0671eeSJoe Wallwork test: 261dd0671eeSJoe Wallwork suffix: hessian_3d_parmmg 262dd0671eeSJoe Wallwork requires: parmmg tetgen 263dd0671eeSJoe Wallwork nsize: 2 264dd0671eeSJoe Wallwork args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor parmmg 26500d473e3SJoe Wallwork 26600d473e3SJoe Wallwork TEST*/ 267