1f155c232SMatthew G. Knepley static char help[] = "First example in homogenization book\n\n"; 2f155c232SMatthew G. Knepley 3f155c232SMatthew G. Knepley #include <petscsnes.h> 4f155c232SMatthew G. Knepley #include <petscdmplex.h> 5f155c232SMatthew G. Knepley #include <petscds.h> 6f155c232SMatthew G. Knepley #include <petscconvest.h> 7f155c232SMatthew G. Knepley #include <petscbag.h> 8f155c232SMatthew G. Knepley 9f155c232SMatthew G. Knepley /* 10f155c232SMatthew G. Knepley To control the refinement use -dm_plex_box_faces <n> or -dm_refine <k>, or both 11f155c232SMatthew G. Knepley 12f155c232SMatthew G. Knepley To see the exact and computed solutions 13f155c232SMatthew G. Knepley 14f155c232SMatthew G. Knepley -compare_view draw -draw_size 500,500 -draw_pause -1 15f155c232SMatthew G. Knepley 16f155c232SMatthew G. Knepley To see the delay in convergence of the discretization use 17f155c232SMatthew G. Knepley 18f155c232SMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 7 -convest_monitor 19f155c232SMatthew G. Knepley 20f155c232SMatthew G. Knepley and to see the proper rate use 21f155c232SMatthew G. Knepley 22f155c232SMatthew G. Knepley -dm_refine 5 -snes_convergence_estimate -convest_num_refine 2 -convest_monitor 23f155c232SMatthew G. Knepley */ 24f155c232SMatthew G. Knepley 259371c9d4SSatish Balay typedef enum { 269371c9d4SSatish Balay MOD_CONSTANT, 279371c9d4SSatish Balay MOD_OSCILLATORY, 289371c9d4SSatish Balay MOD_TANH, 299371c9d4SSatish Balay NUM_MOD_TYPES 309371c9d4SSatish Balay } ModType; 31f155c232SMatthew G. Knepley const char *modTypes[NUM_MOD_TYPES + 1] = {"constant", "oscillatory", "tanh", "unknown"}; 32f155c232SMatthew G. Knepley 33f155c232SMatthew G. Knepley /* Constants */ 349371c9d4SSatish Balay enum { 359371c9d4SSatish Balay EPSILON, 369371c9d4SSatish Balay NUM_CONSTANTS 379371c9d4SSatish Balay }; 38f155c232SMatthew G. Knepley 39f155c232SMatthew G. Knepley typedef struct { 40f155c232SMatthew G. Knepley PetscReal epsilon; /* Wavelength of fine scale oscillation */ 41f155c232SMatthew G. Knepley } Parameter; 42f155c232SMatthew G. Knepley 43f155c232SMatthew G. Knepley typedef struct { 44f155c232SMatthew G. Knepley PetscBag bag; /* Holds problem parameters */ 45f155c232SMatthew G. Knepley ModType modType; /* Model type */ 46f155c232SMatthew G. Knepley } AppCtx; 47f155c232SMatthew G. Knepley 48*2a8381b2SBarry Smith static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx) 49d71ae5a4SJacob Faibussowitsch { 50f155c232SMatthew G. Knepley PetscInt d; 51f155c232SMatthew G. Knepley *u = 1.0; 52f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0 * PETSC_PI * x[d]); 533ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 54f155c232SMatthew G. Knepley } 55f155c232SMatthew G. Knepley 56*2a8381b2SBarry Smith static PetscErrorCode oscillatory_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx) 57d71ae5a4SJacob Faibussowitsch { 58f155c232SMatthew G. Knepley Parameter *param = (Parameter *)ctx; 59f155c232SMatthew G. Knepley const PetscReal eps = param->epsilon; 60f155c232SMatthew G. Knepley 61f155c232SMatthew G. Knepley u[0] = x[0] - x[0] * x[0] + (eps / (2. * PETSC_PI)) * (0.5 - x[0]) * PetscSinReal(2. * PETSC_PI * x[0] / eps) + PetscSqr(eps / (2. * PETSC_PI)) * (1. - PetscCosReal(2. * PETSC_PI * x[0] / eps)); 623ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 63f155c232SMatthew G. Knepley } 64f155c232SMatthew G. Knepley 65d71ae5a4SJacob Faibussowitsch static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 66d71ae5a4SJacob Faibussowitsch { 67f155c232SMatthew G. Knepley PetscInt d; 68f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) { 69f155c232SMatthew G. Knepley PetscScalar v = 1.; 70f155c232SMatthew G. Knepley for (PetscInt e = 0; e < dim; e++) { 71f155c232SMatthew G. Knepley if (e == d) { 72f155c232SMatthew G. Knepley v *= -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 73f155c232SMatthew G. Knepley } else { 74f155c232SMatthew G. Knepley v *= PetscSinReal(2.0 * PETSC_PI * x[d]); 75f155c232SMatthew G. Knepley } 76f155c232SMatthew G. Knepley } 77f155c232SMatthew G. Knepley f0[0] += v; 78f155c232SMatthew G. Knepley } 79f155c232SMatthew G. Knepley } 80f155c232SMatthew G. Knepley 81d71ae5a4SJacob Faibussowitsch static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 82d71ae5a4SJacob Faibussowitsch { 83f155c232SMatthew G. Knepley PetscInt d; 84f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 85f155c232SMatthew G. Knepley } 86f155c232SMatthew G. Knepley 87d71ae5a4SJacob Faibussowitsch static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 88d71ae5a4SJacob Faibussowitsch { 89f155c232SMatthew G. Knepley PetscInt d; 90f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 91f155c232SMatthew G. Knepley } 92f155c232SMatthew G. Knepley 93d71ae5a4SJacob Faibussowitsch static void f0_oscillatory_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 94d71ae5a4SJacob Faibussowitsch { 95f155c232SMatthew G. Knepley f0[0] = -1.; 96f155c232SMatthew G. Knepley } 97f155c232SMatthew G. Knepley 98d71ae5a4SJacob Faibussowitsch static void f1_oscillatory_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 99d71ae5a4SJacob Faibussowitsch { 100f155c232SMatthew G. Knepley const PetscReal eps = PetscRealPart(constants[EPSILON]); 101f155c232SMatthew G. Knepley 102f155c232SMatthew G. Knepley f1[0] = u_x[0] / (2. + PetscCosReal(2. * PETSC_PI * x[0] / eps)); 103f155c232SMatthew G. Knepley } 104f155c232SMatthew G. Knepley 105d71ae5a4SJacob Faibussowitsch static void g3_oscillatory_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 106d71ae5a4SJacob Faibussowitsch { 107f155c232SMatthew G. Knepley const PetscReal eps = PetscRealPart(constants[EPSILON]); 108f155c232SMatthew G. Knepley 109f155c232SMatthew G. Knepley g3[0] = 1. / (2. + PetscCosReal(2. * PETSC_PI * x[0] / eps)); 110f155c232SMatthew G. Knepley } 111f155c232SMatthew G. Knepley 112d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 113d71ae5a4SJacob Faibussowitsch { 114f155c232SMatthew G. Knepley PetscInt mod; 115f155c232SMatthew G. Knepley 116f155c232SMatthew G. Knepley PetscFunctionBeginUser; 117f155c232SMatthew G. Knepley options->modType = MOD_CONSTANT; 118d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Homogenization Problem Options", "DMPLEX"); 119f155c232SMatthew G. Knepley mod = options->modType; 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mod_type", "The model type", "ex36.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL)); 121f155c232SMatthew G. Knepley options->modType = (ModType)mod; 122d0609cedSBarry Smith PetscOptionsEnd(); 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124f155c232SMatthew G. Knepley } 125f155c232SMatthew G. Knepley 126d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *user) 127d71ae5a4SJacob Faibussowitsch { 128f155c232SMatthew G. Knepley PetscBag bag; 129f155c232SMatthew G. Knepley Parameter *p; 130f155c232SMatthew G. Knepley 131f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1329566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(comm, sizeof(Parameter), &user->bag)); 133*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, &p)); 1349566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(user->bag, "par", "Homogenization parameters")); 135f155c232SMatthew G. Knepley bag = user->bag; 1369566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->epsilon, 1.0, "epsilon", "Wavelength of fine scale oscillation")); 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 138f155c232SMatthew G. Knepley } 139f155c232SMatthew G. Knepley 140d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 141d71ae5a4SJacob Faibussowitsch { 142f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1439566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1449566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1459566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1469566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 1479566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 149f155c232SMatthew G. Knepley } 150f155c232SMatthew G. Knepley 151d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 152d71ae5a4SJacob Faibussowitsch { 153f155c232SMatthew G. Knepley PetscDS ds; 154f155c232SMatthew G. Knepley DMLabel label; 1558434afd1SBarry Smith PetscSimplePointFn *ex; 156f155c232SMatthew G. Knepley const PetscInt id = 1; 157f155c232SMatthew G. Knepley void *ctx; 158f155c232SMatthew G. Knepley 159f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1609566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 161300f1712SStefano Zampini PetscCall(PetscBagGetData(user->bag, &ctx)); 162f155c232SMatthew G. Knepley switch (user->modType) { 163f155c232SMatthew G. Knepley case MOD_CONSTANT: 1649566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_homogeneous_u, f1_u)); 1659566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1669566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 167f155c232SMatthew G. Knepley ex = trig_homogeneous_u; 16857d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)ex, NULL, ctx, NULL)); 169f155c232SMatthew G. Knepley break; 170f155c232SMatthew G. Knepley case MOD_OSCILLATORY: 1719566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_oscillatory_u, f1_oscillatory_u)); 1729566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_oscillatory_uu)); 1739566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 174f155c232SMatthew G. Knepley ex = oscillatory_u; 17557d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)ex, NULL, ctx, NULL)); 176f155c232SMatthew G. Knepley break; 177d71ae5a4SJacob Faibussowitsch default: 178d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Unsupported model type: %s (%d)", modTypes[PetscMin(user->modType, NUM_MOD_TYPES)], user->modType); 179f155c232SMatthew G. Knepley } 1809566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, ex, ctx)); 181f155c232SMatthew G. Knepley /* Setup constants */ 182f155c232SMatthew G. Knepley { 183f155c232SMatthew G. Knepley Parameter *param; 184f155c232SMatthew G. Knepley PetscScalar constants[NUM_CONSTANTS]; 185f155c232SMatthew G. Knepley 186*2a8381b2SBarry Smith PetscCall(PetscBagGetData(user->bag, ¶m)); 187f155c232SMatthew G. Knepley 188f155c232SMatthew G. Knepley constants[EPSILON] = param->epsilon; 1899566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants)); 190f155c232SMatthew G. Knepley } 1913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 192f155c232SMatthew G. Knepley } 193f155c232SMatthew G. Knepley 194d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 195d71ae5a4SJacob Faibussowitsch { 196f155c232SMatthew G. Knepley DM cdm = dm; 197f155c232SMatthew G. Knepley PetscFE fe; 198f155c232SMatthew G. Knepley PetscBool simplex; 199f155c232SMatthew G. Knepley PetscInt dim; 200f155c232SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 201f155c232SMatthew G. Knepley 202f155c232SMatthew G. Knepley PetscFunctionBeginUser; 2039566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 205f155c232SMatthew G. Knepley /* Create finite element */ 2069566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 2079566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 2089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name)); 209f155c232SMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 2109566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2119566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2129566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 213f155c232SMatthew G. Knepley while (cdm) { 2149566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2159566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 216f155c232SMatthew G. Knepley } 2179566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 219f155c232SMatthew G. Knepley } 220f155c232SMatthew G. Knepley 221d71ae5a4SJacob Faibussowitsch static PetscErrorCode CompareView(Vec u) 222d71ae5a4SJacob Faibussowitsch { 223f155c232SMatthew G. Knepley DM dm; 224f155c232SMatthew G. Knepley Vec v[2], lv[2], exact; 225f155c232SMatthew G. Knepley PetscOptions options; 226f155c232SMatthew G. Knepley PetscViewer viewer; 227f155c232SMatthew G. Knepley PetscViewerFormat format; 228f155c232SMatthew G. Knepley PetscBool flg; 229f155c232SMatthew G. Knepley PetscInt i; 230f155c232SMatthew G. Knepley const char *name, *prefix; 231f155c232SMatthew G. Knepley 232f155c232SMatthew G. Knepley PetscFunctionBeginUser; 2339566063dSJacob Faibussowitsch PetscCall(VecGetDM(u, &dm)); 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptions((PetscObject)dm, &options)); 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 236648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)dm), options, prefix, "-compare_view", &viewer, &format, &flg)); 237f155c232SMatthew G. Knepley if (flg) { 2389566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &exact)); 2399566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, 0.0, exact, NULL)); 240f155c232SMatthew G. Knepley v[0] = u; 241f155c232SMatthew G. Knepley v[1] = exact; 242f155c232SMatthew G. Knepley for (i = 0; i < 2; ++i) { 2439566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lv[i])); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)v[i], &name)); 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)lv[i], name)); 2469566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, v[i], INSERT_VALUES, lv[i])); 2479566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, v[i], INSERT_VALUES, lv[i])); 2489566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, lv[i], 0., NULL, NULL, NULL)); 249f155c232SMatthew G. Knepley } 2509566063dSJacob Faibussowitsch PetscCall(DMPlexVecView1D(dm, 2, lv, viewer)); 2519566063dSJacob Faibussowitsch for (i = 0; i < 2; ++i) PetscCall(DMRestoreLocalVector(dm, &lv[i])); 2529566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &exact)); 253648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 254f155c232SMatthew G. Knepley } 2553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 256f155c232SMatthew G. Knepley } 257f155c232SMatthew G. Knepley 2589371c9d4SSatish Balay typedef struct { 259f155c232SMatthew G. Knepley Mat Mcoarse; /* Mass matrix on the coarse space */ 260f155c232SMatthew G. Knepley Mat Mfine; /* Mass matrix on the fine space */ 261f155c232SMatthew G. Knepley Mat Ifine; /* Interpolator from coarse to fine */ 262f155c232SMatthew G. Knepley Vec Iscale; /* Scaling vector for restriction */ 263f155c232SMatthew G. Knepley KSP kspCoarse; /* Solver for the coarse mass matrix */ 264f155c232SMatthew G. Knepley Vec tmpfine; /* Temporary vector in the fine space */ 265f155c232SMatthew G. Knepley Vec tmpcoarse; /* Temporary vector in the coarse space */ 266f155c232SMatthew G. Knepley } ProjStruct; 267f155c232SMatthew G. Knepley 268d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyCoarseProjection(Mat Pi) 269d71ae5a4SJacob Faibussowitsch { 270f155c232SMatthew G. Knepley ProjStruct *ctx; 271f155c232SMatthew G. Knepley 272f155c232SMatthew G. Knepley PetscFunctionBegin; 273*2a8381b2SBarry Smith PetscCall(MatShellGetContext(Pi, &ctx)); 2749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Mcoarse)); 2759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Mfine)); 2769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Ifine)); 2779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->Iscale)); 2789566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ctx->kspCoarse)); 2799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpcoarse)); 2809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpfine)); 2819566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 2829566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(Pi, NULL)); 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284f155c232SMatthew G. Knepley } 285f155c232SMatthew G. Knepley 286d71ae5a4SJacob Faibussowitsch static PetscErrorCode CoarseProjection(Mat Pi, Vec x, Vec y) 287d71ae5a4SJacob Faibussowitsch { 288f155c232SMatthew G. Knepley ProjStruct *ctx; 289f155c232SMatthew G. Knepley 290f155c232SMatthew G. Knepley PetscFunctionBegin; 291*2a8381b2SBarry Smith PetscCall(MatShellGetContext(Pi, &ctx)); 2929566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Mfine, x, ctx->tmpfine)); 2939566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpfine, "Fine DG RHS")); 2949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpfine, "fine_dg_")); 2959566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpfine, NULL, "-rhs_view")); 2969566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->Ifine, ctx->tmpfine, ctx->tmpcoarse)); 2979566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpcoarse, "Coarse DG RHS")); 2989566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpcoarse, "coarse_dg_")); 2999566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 3009566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ctx->tmpcoarse, ctx->Iscale, ctx->tmpcoarse)); 3019566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 3029566063dSJacob Faibussowitsch PetscCall(KSPSolve(ctx->kspCoarse, ctx->tmpcoarse, y)); 3033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 304f155c232SMatthew G. Knepley } 305f155c232SMatthew G. Knepley 306d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateCoarseProjection(DM dmc, DM dmf, Mat *Pi) 307d71ae5a4SJacob Faibussowitsch { 308f155c232SMatthew G. Knepley ProjStruct *ctx; 309f155c232SMatthew G. Knepley PetscInt m, n, M, N; 310f155c232SMatthew G. Knepley 311f155c232SMatthew G. Knepley PetscFunctionBegin; 3129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &ctx)); 3139566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmc, &ctx->tmpcoarse)); 3149566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmf, &ctx->tmpfine)); 3159566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->tmpcoarse, &m)); 3169566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->tmpcoarse, &M)); 3179566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->tmpfine, &n)); 3189566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->tmpfine, &N)); 3199566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmc, dmc, &ctx->Mcoarse)); 3209566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmf, dmf, &ctx->Mfine)); 3219566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmc, dmf, &ctx->Ifine, &ctx->Iscale)); 3229566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)dmc), &ctx->kspCoarse)); 3239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->kspCoarse, "coarse_")); 3249566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ctx->kspCoarse, ctx->Mcoarse, ctx->Mcoarse)); 3259566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ctx->kspCoarse)); 3269566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dmc), m, n, M, N, ctx, Pi)); 32757d50842SBarry Smith PetscCall(MatShellSetOperation(*Pi, MATOP_DESTROY, (PetscErrorCodeFn *)DestroyCoarseProjection)); 32857d50842SBarry Smith PetscCall(MatShellSetOperation(*Pi, MATOP_MULT, (PetscErrorCodeFn *)CoarseProjection)); 3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 330f155c232SMatthew G. Knepley } 331f155c232SMatthew G. Knepley 3329371c9d4SSatish Balay typedef struct { 333f155c232SMatthew G. Knepley Mat Ifdg; /* Embed the fine space into the DG version */ 334f155c232SMatthew G. Knepley Mat Pi; /* The L_2 stable projection to the DG coarse space */ 335f155c232SMatthew G. Knepley Vec tmpc; /* A temporary vector in the DG coarse space */ 336f155c232SMatthew G. Knepley Vec tmpf; /* A temporary vector in the DG fine space */ 337f155c232SMatthew G. Knepley } QuasiInterp; 338f155c232SMatthew G. Knepley 339d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyQuasiInterpolator(Mat P) 340d71ae5a4SJacob Faibussowitsch { 341f155c232SMatthew G. Knepley QuasiInterp *ctx; 342f155c232SMatthew G. Knepley 343f155c232SMatthew G. Knepley PetscFunctionBegin; 344*2a8381b2SBarry Smith PetscCall(MatShellGetContext(P, &ctx)); 3459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Ifdg)); 3469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Pi)); 3479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpc)); 3489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpf)); 3499566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 3509566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(P, NULL)); 3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 352f155c232SMatthew G. Knepley } 353f155c232SMatthew G. Knepley 354d71ae5a4SJacob Faibussowitsch static PetscErrorCode QuasiInterpolate(Mat P, Vec x, Vec y) 355d71ae5a4SJacob Faibussowitsch { 356f155c232SMatthew G. Knepley QuasiInterp *ctx; 357f155c232SMatthew G. Knepley DM dmcdg, dmc; 358f155c232SMatthew G. Knepley Vec ly; 359f155c232SMatthew G. Knepley 360f155c232SMatthew G. Knepley PetscFunctionBegin; 361*2a8381b2SBarry Smith PetscCall(MatShellGetContext(P, &ctx)); 3629566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Ifdg, x, ctx->tmpf)); 363f155c232SMatthew G. Knepley 3649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpf, "Fine DG Potential")); 3659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpf, "fine_dg_")); 3669566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpf, NULL, "-vec_view")); 3679566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Pi, x, ctx->tmpc)); 368f155c232SMatthew G. Knepley 3699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpc, "Coarse DG Potential")); 3709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpc, "coarse_dg_")); 3719566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpc, NULL, "-vec_view")); 3729566063dSJacob Faibussowitsch PetscCall(VecGetDM(ctx->tmpc, &dmcdg)); 373f155c232SMatthew G. Knepley 3749566063dSJacob Faibussowitsch PetscCall(VecGetDM(y, &dmc)); 3759566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmc, &ly)); 3769566063dSJacob Faibussowitsch PetscCall(DMPlexComputeClementInterpolant(dmcdg, ctx->tmpc, ly)); 3779566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dmc, ly, INSERT_VALUES, y)); 3789566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dmc, ly, INSERT_VALUES, y)); 3799566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmc, &ly)); 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381f155c232SMatthew G. Knepley } 382f155c232SMatthew G. Knepley 383d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateQuasiInterpolator(DM dmc, DM dmf, Mat *P) 384d71ae5a4SJacob Faibussowitsch { 385f155c232SMatthew G. Knepley QuasiInterp *ctx; 386f155c232SMatthew G. Knepley DM dmcdg, dmfdg; 387f155c232SMatthew G. Knepley PetscFE fe; 388f155c232SMatthew G. Knepley Vec x, y; 389f155c232SMatthew G. Knepley DMPolytopeType ct; 390f155c232SMatthew G. Knepley PetscInt dim, cStart, m, n, M, N; 391f155c232SMatthew G. Knepley 392f155c232SMatthew G. Knepley PetscFunctionBegin; 3939566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(1, &ctx)); 3949566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmc, &x)); 3959566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmf, &y)); 3969566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &m)); 3979566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &M)); 3989566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 3999566063dSJacob Faibussowitsch PetscCall(VecGetSize(y, &N)); 4009566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmc, &x)); 4019566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmf, &y)); 402f155c232SMatthew G. Knepley 4039566063dSJacob Faibussowitsch PetscCall(DMClone(dmf, &dmfdg)); 4049566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmfdg, &dim)); 4059566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmfdg, 0, &cStart, NULL)); 4069566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dmfdg, cStart, &ct)); 4079566063dSJacob Faibussowitsch PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "fine_dg_", PETSC_DETERMINE, &fe)); 4089566063dSJacob Faibussowitsch PetscCall(DMSetField(dmfdg, 0, NULL, (PetscObject)fe)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 4109566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmfdg)); 4119566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmf, dmfdg, &ctx->Ifdg, NULL)); 4129566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmfdg, &ctx->tmpf)); 4139566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmfdg)); 414f155c232SMatthew G. Knepley 4159566063dSJacob Faibussowitsch PetscCall(DMClone(dmc, &dmcdg)); 4169566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmcdg, &dim)); 4179566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmcdg, 0, &cStart, NULL)); 4189566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dmcdg, cStart, &ct)); 4199566063dSJacob Faibussowitsch PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "coarse_dg_", PETSC_DETERMINE, &fe)); 4209566063dSJacob Faibussowitsch PetscCall(DMSetField(dmcdg, 0, NULL, (PetscObject)fe)); 4219566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 4229566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmcdg)); 423f155c232SMatthew G. Knepley 4249566063dSJacob Faibussowitsch PetscCall(CreateCoarseProjection(dmcdg, dmf, &ctx->Pi)); 4259566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmcdg, &ctx->tmpc)); 4269566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmcdg)); 427f155c232SMatthew G. Knepley 4289566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dmc), m, n, M, N, ctx, P)); 42957d50842SBarry Smith PetscCall(MatShellSetOperation(*P, MATOP_DESTROY, (PetscErrorCodeFn *)DestroyQuasiInterpolator)); 43057d50842SBarry Smith PetscCall(MatShellSetOperation(*P, MATOP_MULT, (PetscErrorCodeFn *)QuasiInterpolate)); 4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 432f155c232SMatthew G. Knepley } 433f155c232SMatthew G. Knepley 434d71ae5a4SJacob Faibussowitsch static PetscErrorCode CoarseTest(DM dm, Vec u, AppCtx *user) 435d71ae5a4SJacob Faibussowitsch { 436f155c232SMatthew G. Knepley DM dmc; 437f155c232SMatthew G. Knepley Mat P; /* The quasi-interpolator to the coarse space */ 438f155c232SMatthew G. Knepley Vec uc; 439f155c232SMatthew G. Knepley 440f155c232SMatthew G. Knepley PetscFunctionBegin; 4413ba16761SJacob Faibussowitsch if (user->modType == MOD_CONSTANT) PetscFunctionReturn(PETSC_SUCCESS); 4429566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &dmc)); 4439566063dSJacob Faibussowitsch PetscCall(DMSetType(dmc, DMPLEX)); 4449566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dmc, "coarse_")); 4459566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dmc, user)); 4469566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dmc)); 4479566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmc, NULL, "-dm_view")); 448f155c232SMatthew G. Knepley 4499566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dmc, "potential", SetupPrimalProblem, user)); 4509566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmc, &uc)); 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)uc, "potential")); 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)uc, "coarse_")); 453f155c232SMatthew G. Knepley 4549566063dSJacob Faibussowitsch PetscCall(CreateQuasiInterpolator(dmc, dm, &P)); 455f155c232SMatthew G. Knepley #if 1 4569566063dSJacob Faibussowitsch PetscCall(MatMult(P, u, uc)); 457f155c232SMatthew G. Knepley #else 458f155c232SMatthew G. Knepley { 459f155c232SMatthew G. Knepley Mat In; 460f155c232SMatthew G. Knepley Vec sc; 461f155c232SMatthew G. Knepley 4629566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmc, dm, &In, &sc)); 4639566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(In, u, uc)); 4649566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(uc, sc, uc)); 4659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&In)); 4669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sc)); 467f155c232SMatthew G. Knepley } 468f155c232SMatthew G. Knepley #endif 4699566063dSJacob Faibussowitsch PetscCall(CompareView(uc)); 470f155c232SMatthew G. Knepley 4719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 4729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uc)); 4739566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmc)); 4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 475f155c232SMatthew G. Knepley } 476f155c232SMatthew G. Knepley 477d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 478d71ae5a4SJacob Faibussowitsch { 479f155c232SMatthew G. Knepley DM dm; /* Problem specification */ 480f155c232SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 481f155c232SMatthew G. Knepley Vec u; /* Solutions */ 482f155c232SMatthew G. Knepley AppCtx user; /* User-defined work context */ 483f155c232SMatthew G. Knepley 484327415f7SBarry Smith PetscFunctionBeginUser; 4859566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4869566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4879566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 488f155c232SMatthew G. Knepley /* Primal system */ 4899566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 4909566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4919566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 4929566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 4939566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 4949566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 4959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "potential")); 4966493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 4979566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 4989566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 4999566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 5009566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 5019566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 5029566063dSJacob Faibussowitsch PetscCall(CompareView(u)); 503f155c232SMatthew G. Knepley /* Looking at a coarse problem */ 5049566063dSJacob Faibussowitsch PetscCall(CoarseTest(dm, u, &user)); 505f155c232SMatthew G. Knepley /* Cleanup */ 5069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5079566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 5089566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5099566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 5109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 511b122ec5aSJacob Faibussowitsch return 0; 512f155c232SMatthew G. Knepley } 513f155c232SMatthew G. Knepley 514f155c232SMatthew G. Knepley /*TEST 515f155c232SMatthew G. Knepley 516f155c232SMatthew G. Knepley test: 517f155c232SMatthew G. Knepley suffix: 1d_p1_constant 518f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dmsnes_check 519f155c232SMatthew G. Knepley 520f155c232SMatthew G. Knepley test: 521f155c232SMatthew G. Knepley suffix: 1d_p1_constant_conv 522f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 \ 523f155c232SMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 2 524f155c232SMatthew G. Knepley 525f155c232SMatthew G. Knepley test: 526f155c232SMatthew G. Knepley suffix: 1d_p1_oscillatory 527f155c232SMatthew G. Knepley args: -mod_type oscillatory -epsilon 0.03125 \ 528f155c232SMatthew G. Knepley -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dm_refine 2 -dmsnes_check \ 529f155c232SMatthew G. Knepley -coarse_dm_plex_dim 1 -coarse_dm_plex_box_faces 4 -coarse_dm_plex_hash_location \ 530f155c232SMatthew G. Knepley -fine_dg_petscspace_degree 1 -fine_dg_petscdualspace_lagrange_continuity 0 \ 531f155c232SMatthew G. Knepley -coarse_dg_petscspace_degree 1 -coarse_dg_petscdualspace_lagrange_continuity 0 532f155c232SMatthew G. Knepley 533f155c232SMatthew G. Knepley TEST*/ 534