1a12d352dSMatthew G. Knepley static char help[] = "Tests 1D nested mesh refinement.\n\n"; 2a12d352dSMatthew G. Knepley 3a12d352dSMatthew G. Knepley #include <petscdmplex.h> 4a12d352dSMatthew G. Knepley #include <petscds.h> 5a12d352dSMatthew G. Knepley 6a12d352dSMatthew G. Knepley typedef struct { 7a12d352dSMatthew G. Knepley PetscInt Nr; /* Number of refinements */ 8a12d352dSMatthew G. Knepley PetscSimplePointFunc funcs[1]; /* Functions to test */ 9a12d352dSMatthew G. Knepley } AppCtx; 10a12d352dSMatthew G. Knepley 11a12d352dSMatthew G. Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 12a12d352dSMatthew G. Knepley { 13a12d352dSMatthew G. Knepley u[0] = 1.; 14a12d352dSMatthew G. Knepley return 0; 15a12d352dSMatthew G. Knepley } 16a12d352dSMatthew G. Knepley 17a12d352dSMatthew G. Knepley static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 18a12d352dSMatthew G. Knepley { 19a12d352dSMatthew G. Knepley u[0] = 2.*x[0] + 1.; 20a12d352dSMatthew G. Knepley return 0; 21a12d352dSMatthew G. Knepley } 22a12d352dSMatthew G. Knepley 23a12d352dSMatthew G. Knepley static PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 24a12d352dSMatthew G. Knepley { 25a12d352dSMatthew G. Knepley u[0] = 3.*x[0]*x[0] + 2.*x[0] + 1.; 26a12d352dSMatthew G. Knepley return 0; 27a12d352dSMatthew G. Knepley } 28a12d352dSMatthew G. Knepley 29a12d352dSMatthew G. Knepley static PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 30a12d352dSMatthew G. Knepley { 31a12d352dSMatthew G. Knepley u[0] = 4.*x[0]*x[0]*x[0] + 3.*x[0]*x[0] + 2.*x[0] + 1.; 32a12d352dSMatthew G. Knepley return 0; 33a12d352dSMatthew G. Knepley } 34a12d352dSMatthew G. Knepley 35a12d352dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 36a12d352dSMatthew G. Knepley { 37a12d352dSMatthew G. Knepley PetscErrorCode ierr; 38a12d352dSMatthew G. Knepley 39a12d352dSMatthew G. Knepley PetscFunctionBeginUser; 40a12d352dSMatthew G. Knepley options->Nr = 1; 41*9566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "1D Refinement Options", "DMPLEX");PetscCall(ierr); 42*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-num_refine", "Refine cycles", "ex46.c", options->Nr, &options->Nr, NULL)); 43*9566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 44a12d352dSMatthew G. Knepley PetscFunctionReturn(0); 45a12d352dSMatthew G. Knepley } 46a12d352dSMatthew G. Knepley 47a12d352dSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 48a12d352dSMatthew G. Knepley { 49a12d352dSMatthew G. Knepley PetscFunctionBeginUser; 50*9566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 51*9566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 52*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 53*9566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 54*9566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 55a12d352dSMatthew G. Knepley PetscFunctionReturn(0); 56a12d352dSMatthew G. Knepley } 57a12d352dSMatthew G. Knepley 58a12d352dSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 59a12d352dSMatthew G. Knepley { 60a12d352dSMatthew G. Knepley DM cdm = dm; 61a12d352dSMatthew G. Knepley PetscFE fe; 62a12d352dSMatthew G. Knepley PetscSpace sp; 63a12d352dSMatthew G. Knepley PetscInt dim, deg; 64a12d352dSMatthew G. Knepley 65a12d352dSMatthew G. Knepley PetscFunctionBeginUser; 66*9566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 67*9566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, -1, &fe)); 68*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, "scalar")); 69*9566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 70*9566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 71a12d352dSMatthew G. Knepley while (cdm) { 72*9566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 73*9566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 74a12d352dSMatthew G. Knepley } 75*9566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &sp)); 76*9566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(sp, °, NULL)); 77a12d352dSMatthew G. Knepley switch (deg) { 78a12d352dSMatthew G. Knepley case 0: user->funcs[0] = constant;break; 79a12d352dSMatthew G. Knepley case 1: user->funcs[0] = linear;break; 80a12d352dSMatthew G. Knepley case 2: user->funcs[0] = quadratic;break; 81a12d352dSMatthew G. Knepley case 3: user->funcs[0] = cubic;break; 8298921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Could not determine function to test for degree %D", deg); 83a12d352dSMatthew G. Knepley } 84*9566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 85a12d352dSMatthew G. Knepley PetscFunctionReturn(0); 86a12d352dSMatthew G. Knepley } 87a12d352dSMatthew G. Knepley 88a12d352dSMatthew G. Knepley static PetscErrorCode CheckError(DM dm, Vec u, PetscSimplePointFunc funcs[]) 89a12d352dSMatthew G. Knepley { 90a12d352dSMatthew G. Knepley PetscReal error, tol = PETSC_SMALL; 91a12d352dSMatthew G. Knepley MPI_Comm comm; 92a12d352dSMatthew G. Knepley 93a12d352dSMatthew G. Knepley PetscFunctionBeginUser; 94*9566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, funcs, NULL, u, &error)); 95*9566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 96*9566063dSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL at tolerance %g error %g\n", (double)tol,(double) error)); 97*9566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass at tolerance %g\n", (double)tol)); 98a12d352dSMatthew G. Knepley PetscFunctionReturn(0); 99a12d352dSMatthew G. Knepley } 100a12d352dSMatthew G. Knepley 101a12d352dSMatthew G. Knepley int main(int argc, char **argv) 102a12d352dSMatthew G. Knepley { 103a12d352dSMatthew G. Knepley DM dm; 104a12d352dSMatthew G. Knepley Vec u; 105a12d352dSMatthew G. Knepley AppCtx user; 106a12d352dSMatthew G. Knepley PetscInt cStart, cEnd, c, r; 107a12d352dSMatthew G. Knepley 108*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 109*9566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 110*9566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 111*9566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 112*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u)); 113*9566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, user.funcs, NULL, INSERT_ALL_VALUES, u)); 114*9566063dSJacob Faibussowitsch PetscCall(CheckError(dm, u, user.funcs)); 115a12d352dSMatthew G. Knepley for (r = 0; r < user.Nr; ++r) { 116a12d352dSMatthew G. Knepley DM adm; 117a12d352dSMatthew G. Knepley DMLabel adapt; 118a12d352dSMatthew G. Knepley Vec au; 119a12d352dSMatthew G. Knepley Mat Interp; 120a12d352dSMatthew G. Knepley 121*9566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adapt)); 122*9566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(adapt, DM_ADAPT_COARSEN)); 123*9566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 124a12d352dSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 125*9566063dSJacob Faibussowitsch if (c % 2) PetscCall(DMLabelSetValue(adapt, c, DM_ADAPT_REFINE)); 126a12d352dSMatthew G. Knepley } 127*9566063dSJacob Faibussowitsch PetscCall(DMAdaptLabel(dm, adapt, &adm)); 128*9566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adapt)); 129*9566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) adm, "Adapted Mesh")); 130*9566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(adm, NULL, "-dm_view")); 131a12d352dSMatthew G. Knepley 132*9566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dm, adm, &Interp, NULL)); 133*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(adm, &au)); 134*9566063dSJacob Faibussowitsch PetscCall(MatInterpolate(Interp, u, au)); 135*9566063dSJacob Faibussowitsch PetscCall(CheckError(adm, au, user.funcs)); 136*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp)); 137*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 138*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 139a12d352dSMatthew G. Knepley dm = adm; 140a12d352dSMatthew G. Knepley u = au; 141a12d352dSMatthew G. Knepley } 142*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 143*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 144*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 145b122ec5aSJacob Faibussowitsch return 0; 146a12d352dSMatthew G. Knepley } 147a12d352dSMatthew G. Knepley 148a12d352dSMatthew G. Knepley /*TEST 149a12d352dSMatthew G. Knepley 150a12d352dSMatthew G. Knepley test: 151a12d352dSMatthew G. Knepley suffix: 0 152a12d352dSMatthew G. Knepley args: -num_refine 4 -petscspace_degree 3 \ 153a12d352dSMatthew G. Knepley -dm_plex_dim 1 -dm_plex_box_faces 5 -dm_plex_transform_type refine_1d -dm_plex_hash_location -dm_view 154a12d352dSMatthew G. Knepley 155a12d352dSMatthew G. Knepley TEST*/ 156