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 25f155c232SMatthew G. Knepley typedef enum {MOD_CONSTANT, MOD_OSCILLATORY, MOD_TANH, NUM_MOD_TYPES} ModType; 26f155c232SMatthew G. Knepley const char *modTypes[NUM_MOD_TYPES+1] = {"constant", "oscillatory", "tanh", "unknown"}; 27f155c232SMatthew G. Knepley 28f155c232SMatthew G. Knepley /* Constants */ 29f155c232SMatthew G. Knepley enum {EPSILON, NUM_CONSTANTS}; 30f155c232SMatthew G. Knepley 31f155c232SMatthew G. Knepley typedef struct { 32f155c232SMatthew G. Knepley PetscReal epsilon; /* Wavelength of fine scale oscillation */ 33f155c232SMatthew G. Knepley } Parameter; 34f155c232SMatthew G. Knepley 35f155c232SMatthew G. Knepley typedef struct { 36f155c232SMatthew G. Knepley PetscBag bag; /* Holds problem parameters */ 37f155c232SMatthew G. Knepley ModType modType; /* Model type */ 38f155c232SMatthew G. Knepley } AppCtx; 39f155c232SMatthew G. Knepley 40f155c232SMatthew G. Knepley static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 41f155c232SMatthew G. Knepley { 42f155c232SMatthew G. Knepley PetscInt d; 43f155c232SMatthew G. Knepley *u = 1.0; 44f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0*PETSC_PI*x[d]); 45f155c232SMatthew G. Knepley return 0; 46f155c232SMatthew G. Knepley } 47f155c232SMatthew G. Knepley 48f155c232SMatthew G. Knepley static PetscErrorCode oscillatory_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 49f155c232SMatthew G. Knepley { 50f155c232SMatthew G. Knepley Parameter *param = (Parameter *) ctx; 51f155c232SMatthew G. Knepley const PetscReal eps = param->epsilon; 52f155c232SMatthew G. Knepley 53f155c232SMatthew 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)); 54f155c232SMatthew G. Knepley return 0; 55f155c232SMatthew G. Knepley } 56f155c232SMatthew G. Knepley 57f155c232SMatthew G. Knepley static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 58f155c232SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 59f155c232SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 60f155c232SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 61f155c232SMatthew G. Knepley { 62f155c232SMatthew G. Knepley PetscInt d; 63f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) { 64f155c232SMatthew G. Knepley PetscScalar v = 1.; 65f155c232SMatthew G. Knepley for (PetscInt e = 0; e < dim; e++) { 66f155c232SMatthew G. Knepley if (e == d) { 67f155c232SMatthew G. Knepley v *= -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 68f155c232SMatthew G. Knepley } else { 69f155c232SMatthew G. Knepley v *= PetscSinReal(2.0*PETSC_PI*x[d]); 70f155c232SMatthew G. Knepley } 71f155c232SMatthew G. Knepley } 72f155c232SMatthew G. Knepley f0[0] += v; 73f155c232SMatthew G. Knepley } 74f155c232SMatthew G. Knepley } 75f155c232SMatthew G. Knepley 76f155c232SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 77f155c232SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 78f155c232SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 79f155c232SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 80f155c232SMatthew G. Knepley { 81f155c232SMatthew G. Knepley PetscInt d; 82f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 83f155c232SMatthew G. Knepley } 84f155c232SMatthew G. Knepley 85f155c232SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 86f155c232SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 87f155c232SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 88f155c232SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 89f155c232SMatthew G. Knepley { 90f155c232SMatthew G. Knepley PetscInt d; 91f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 92f155c232SMatthew G. Knepley } 93f155c232SMatthew G. Knepley 94f155c232SMatthew G. Knepley static void f0_oscillatory_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 95f155c232SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 96f155c232SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 97f155c232SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 98f155c232SMatthew G. Knepley { 99f155c232SMatthew G. Knepley f0[0] = -1.; 100f155c232SMatthew G. Knepley } 101f155c232SMatthew G. Knepley 102f155c232SMatthew G. Knepley static void f1_oscillatory_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 103f155c232SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 104f155c232SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 105f155c232SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 106f155c232SMatthew G. Knepley { 107f155c232SMatthew G. Knepley const PetscReal eps = PetscRealPart(constants[EPSILON]); 108f155c232SMatthew G. Knepley 109f155c232SMatthew G. Knepley f1[0] = u_x[0] / (2. + PetscCosReal(2.*PETSC_PI*x[0]/eps)); 110f155c232SMatthew G. Knepley } 111f155c232SMatthew G. Knepley 112f155c232SMatthew G. Knepley static void g3_oscillatory_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 113f155c232SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 114f155c232SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 115f155c232SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 116f155c232SMatthew G. Knepley { 117f155c232SMatthew G. Knepley const PetscReal eps = PetscRealPart(constants[EPSILON]); 118f155c232SMatthew G. Knepley 119f155c232SMatthew G. Knepley g3[0] = 1. / (2. + PetscCosReal(2.*PETSC_PI*x[0]/eps)); 120f155c232SMatthew G. Knepley } 121f155c232SMatthew G. Knepley 122f155c232SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 123f155c232SMatthew G. Knepley { 124f155c232SMatthew G. Knepley PetscInt mod; 125f155c232SMatthew G. Knepley PetscErrorCode ierr; 126f155c232SMatthew G. Knepley 127f155c232SMatthew G. Knepley PetscFunctionBeginUser; 128f155c232SMatthew G. Knepley options->modType = MOD_CONSTANT; 129f155c232SMatthew G. Knepley 130f155c232SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Homogenization Problem Options", "DMPLEX");CHKERRQ(ierr); 131f155c232SMatthew G. Knepley mod = options->modType; 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-mod_type", "The model type", "ex36.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL)); 133f155c232SMatthew G. Knepley options->modType = (ModType) mod; 134f155c232SMatthew G. Knepley ierr = PetscOptionsEnd(); 135f155c232SMatthew G. Knepley PetscFunctionReturn(0); 136f155c232SMatthew G. Knepley } 137f155c232SMatthew G. Knepley 138f155c232SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *user) 139f155c232SMatthew G. Knepley { 140f155c232SMatthew G. Knepley PetscBag bag; 141f155c232SMatthew G. Knepley Parameter *p; 142f155c232SMatthew G. Knepley 143f155c232SMatthew G. Knepley PetscFunctionBeginUser; 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagCreate(comm, sizeof(Parameter), &user->bag)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) &p)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagSetName(user->bag, "par", "Homogenization parameters")); 147f155c232SMatthew G. Knepley bag = user->bag; 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterReal(bag, &p->epsilon, 1.0, "epsilon", "Wavelength of fine scale oscillation")); 149f155c232SMatthew G. Knepley PetscFunctionReturn(0); 150f155c232SMatthew G. Knepley } 151f155c232SMatthew G. Knepley 152f155c232SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 153f155c232SMatthew G. Knepley { 154f155c232SMatthew G. Knepley PetscFunctionBeginUser; 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(*dm, user)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 160f155c232SMatthew G. Knepley PetscFunctionReturn(0); 161f155c232SMatthew G. Knepley } 162f155c232SMatthew G. Knepley 163f155c232SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 164f155c232SMatthew G. Knepley { 165f155c232SMatthew G. Knepley PetscDS ds; 166f155c232SMatthew G. Knepley DMLabel label; 167f155c232SMatthew G. Knepley PetscSimplePointFunc ex; 168f155c232SMatthew G. Knepley const PetscInt id = 1; 169f155c232SMatthew G. Knepley void *ctx; 170f155c232SMatthew G. Knepley 171f155c232SMatthew G. Knepley PetscFunctionBeginUser; 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) &ctx)); 174f155c232SMatthew G. Knepley switch (user->modType) { 175f155c232SMatthew G. Knepley case MOD_CONSTANT: 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_homogeneous_u, f1_u)); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 179f155c232SMatthew G. Knepley ex = trig_homogeneous_u; 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, ctx, NULL)); 181f155c232SMatthew G. Knepley break; 182f155c232SMatthew G. Knepley case MOD_OSCILLATORY: 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_oscillatory_u, f1_oscillatory_u)); 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_oscillatory_uu)); 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 186f155c232SMatthew G. Knepley ex = oscillatory_u; 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, ctx, NULL)); 188f155c232SMatthew G. Knepley break; 189f155c232SMatthew G. Knepley default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported model type: %s (%D)", modTypes[PetscMin(user->modType, NUM_MOD_TYPES)], user->modType); 190f155c232SMatthew G. Knepley } 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, ex, ctx)); 192f155c232SMatthew G. Knepley /* Setup constants */ 193f155c232SMatthew G. Knepley { 194f155c232SMatthew G. Knepley Parameter *param; 195f155c232SMatthew G. Knepley PetscScalar constants[NUM_CONSTANTS]; 196f155c232SMatthew G. Knepley 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 198f155c232SMatthew G. Knepley 199f155c232SMatthew G. Knepley constants[EPSILON] = param->epsilon; 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetConstants(ds, NUM_CONSTANTS, constants)); 201f155c232SMatthew G. Knepley } 202f155c232SMatthew G. Knepley PetscFunctionReturn(0); 203f155c232SMatthew G. Knepley } 204f155c232SMatthew G. Knepley 205f155c232SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 206f155c232SMatthew G. Knepley { 207f155c232SMatthew G. Knepley DM cdm = dm; 208f155c232SMatthew G. Knepley PetscFE fe; 209f155c232SMatthew G. Knepley PetscBool simplex; 210f155c232SMatthew G. Knepley PetscInt dim; 211f155c232SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 212f155c232SMatthew G. Knepley 213f155c232SMatthew G. Knepley PetscFunctionBeginUser; 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 216f155c232SMatthew G. Knepley /* Create finite element */ 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 219*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, name)); 220f155c232SMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 222*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 223*5f80ce2aSJacob Faibussowitsch CHKERRQ((*setup)(dm, user)); 224f155c232SMatthew G. Knepley while (cdm) { 225*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm,cdm)); 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 227f155c232SMatthew G. Knepley } 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 229f155c232SMatthew G. Knepley PetscFunctionReturn(0); 230f155c232SMatthew G. Knepley } 231f155c232SMatthew G. Knepley 232f155c232SMatthew G. Knepley static PetscErrorCode CompareView(Vec u) 233f155c232SMatthew G. Knepley { 234f155c232SMatthew G. Knepley DM dm; 235f155c232SMatthew G. Knepley Vec v[2], lv[2], exact; 236f155c232SMatthew G. Knepley PetscOptions options; 237f155c232SMatthew G. Knepley PetscViewer viewer; 238f155c232SMatthew G. Knepley PetscViewerFormat format; 239f155c232SMatthew G. Knepley PetscBool flg; 240f155c232SMatthew G. Knepley PetscInt i; 241f155c232SMatthew G. Knepley const char *name, *prefix; 242f155c232SMatthew G. Knepley 243f155c232SMatthew G. Knepley PetscFunctionBeginUser; 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(u, &dm)); 245*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptions((PetscObject) dm, &options)); 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetViewer(PetscObjectComm((PetscObject) dm), options, prefix, "-compare_view", &viewer, &format, &flg)); 248f155c232SMatthew G. Knepley if (flg) { 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &exact)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeExactSolution(dm, 0.0, exact, NULL)); 251f155c232SMatthew G. Knepley v[0] = u; 252f155c232SMatthew G. Knepley v[1] = exact; 253f155c232SMatthew G. Knepley for (i = 0; i < 2; ++i) { 254*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lv[i])); 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject) v[i], &name)); 256*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lv[i], name)); 257*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(dm, v[i], INSERT_VALUES, lv[i])); 258*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(dm, v[i], INSERT_VALUES, lv[i])); 259*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, lv[i], 0., NULL, NULL, NULL)); 260f155c232SMatthew G. Knepley } 261*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecView1D(dm, 2, lv, viewer)); 262*5f80ce2aSJacob Faibussowitsch for (i = 0; i < 2; ++i) CHKERRQ(DMRestoreLocalVector(dm, &lv[i])); 263*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &exact)); 264*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 265f155c232SMatthew G. Knepley } 266f155c232SMatthew G. Knepley PetscFunctionReturn(0); 267f155c232SMatthew G. Knepley } 268f155c232SMatthew G. Knepley 269f155c232SMatthew G. Knepley typedef struct 270f155c232SMatthew G. Knepley { 271f155c232SMatthew G. Knepley Mat Mcoarse; /* Mass matrix on the coarse space */ 272f155c232SMatthew G. Knepley Mat Mfine; /* Mass matrix on the fine space */ 273f155c232SMatthew G. Knepley Mat Ifine; /* Interpolator from coarse to fine */ 274f155c232SMatthew G. Knepley Vec Iscale; /* Scaling vector for restriction */ 275f155c232SMatthew G. Knepley KSP kspCoarse; /* Solver for the coarse mass matrix */ 276f155c232SMatthew G. Knepley Vec tmpfine; /* Temporary vector in the fine space */ 277f155c232SMatthew G. Knepley Vec tmpcoarse; /* Temporary vector in the coarse space */ 278f155c232SMatthew G. Knepley } ProjStruct; 279f155c232SMatthew G. Knepley 280f155c232SMatthew G. Knepley static PetscErrorCode DestroyCoarseProjection(Mat Pi) 281f155c232SMatthew G. Knepley { 282f155c232SMatthew G. Knepley ProjStruct *ctx; 283f155c232SMatthew G. Knepley 284f155c232SMatthew G. Knepley PetscFunctionBegin; 285*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(Pi, (void **) &ctx)); 286*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ctx->Mcoarse)); 287*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ctx->Mfine)); 288*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ctx->Ifine)); 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx->Iscale)); 290*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&ctx->kspCoarse)); 291*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx->tmpcoarse)); 292*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx->tmpfine)); 293*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ctx)); 294*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(Pi, NULL)); 295f155c232SMatthew G. Knepley PetscFunctionReturn(0); 296f155c232SMatthew G. Knepley } 297f155c232SMatthew G. Knepley 298f155c232SMatthew G. Knepley static PetscErrorCode CoarseProjection(Mat Pi, Vec x, Vec y) 299f155c232SMatthew G. Knepley { 300f155c232SMatthew G. Knepley ProjStruct *ctx; 301f155c232SMatthew G. Knepley 302f155c232SMatthew G. Knepley PetscFunctionBegin; 303*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(Pi, (void **) &ctx)); 304*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(ctx->Mfine, x, ctx->tmpfine)); 305*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) ctx->tmpfine, "Fine DG RHS")); 306*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpfine, "fine_dg_")); 307*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(ctx->tmpfine, NULL, "-rhs_view")); 308*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(ctx->Ifine, ctx->tmpfine, ctx->tmpcoarse)); 309*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) ctx->tmpcoarse, "Coarse DG RHS")); 310*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpcoarse, "coarse_dg_")); 311*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 312*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(ctx->tmpcoarse, ctx->Iscale, ctx->tmpcoarse)); 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(ctx->kspCoarse, ctx->tmpcoarse, y)); 315f155c232SMatthew G. Knepley PetscFunctionReturn(0); 316f155c232SMatthew G. Knepley } 317f155c232SMatthew G. Knepley 318f155c232SMatthew G. Knepley static PetscErrorCode CreateCoarseProjection(DM dmc, DM dmf, Mat *Pi) 319f155c232SMatthew G. Knepley { 320f155c232SMatthew G. Knepley ProjStruct *ctx; 321f155c232SMatthew G. Knepley PetscInt m, n, M, N; 322f155c232SMatthew G. Knepley 323f155c232SMatthew G. Knepley PetscFunctionBegin; 324*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(1, &ctx)); 325*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dmc, &ctx->tmpcoarse)); 326*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dmf, &ctx->tmpfine)); 327*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(ctx->tmpcoarse, &m)); 328*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(ctx->tmpcoarse, &M)); 329*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(ctx->tmpfine, &n)); 330*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(ctx->tmpfine, &N)); 331*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMassMatrix(dmc, dmc, &ctx->Mcoarse)); 332*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMassMatrix(dmf, dmf, &ctx->Mfine)); 333*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(dmc, dmf, &ctx->Ifine, &ctx->Iscale)); 334*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(PetscObjectComm((PetscObject) dmc), &ctx->kspCoarse)); 335*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) ctx->kspCoarse, "coarse_")); 336*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(ctx->kspCoarse, ctx->Mcoarse, ctx->Mcoarse)); 337*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(ctx->kspCoarse)); 338*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PetscObjectComm((PetscObject) dmc), m, n, M, N, ctx, Pi)); 339*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*Pi, MATOP_DESTROY, (void (*)(void)) DestroyCoarseProjection)); 340*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*Pi, MATOP_MULT, (void (*)(void)) CoarseProjection)); 341f155c232SMatthew G. Knepley PetscFunctionReturn(0); 342f155c232SMatthew G. Knepley } 343f155c232SMatthew G. Knepley 344f155c232SMatthew G. Knepley typedef struct 345f155c232SMatthew G. Knepley { 346f155c232SMatthew G. Knepley Mat Ifdg; /* Embed the fine space into the DG version */ 347f155c232SMatthew G. Knepley Mat Pi; /* The L_2 stable projection to the DG coarse space */ 348f155c232SMatthew G. Knepley Vec tmpc; /* A temporary vector in the DG coarse space */ 349f155c232SMatthew G. Knepley Vec tmpf; /* A temporary vector in the DG fine space */ 350f155c232SMatthew G. Knepley } QuasiInterp; 351f155c232SMatthew G. Knepley 352f155c232SMatthew G. Knepley static PetscErrorCode DestroyQuasiInterpolator(Mat P) 353f155c232SMatthew G. Knepley { 354f155c232SMatthew G. Knepley QuasiInterp *ctx; 355f155c232SMatthew G. Knepley 356f155c232SMatthew G. Knepley PetscFunctionBegin; 357*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(P, (void **) &ctx)); 358*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ctx->Ifdg)); 359*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&ctx->Pi)); 360*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx->tmpc)); 361*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctx->tmpf)); 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ctx)); 363*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetContext(P, NULL)); 364f155c232SMatthew G. Knepley PetscFunctionReturn(0); 365f155c232SMatthew G. Knepley } 366f155c232SMatthew G. Knepley 367f155c232SMatthew G. Knepley static PetscErrorCode QuasiInterpolate(Mat P, Vec x, Vec y) 368f155c232SMatthew G. Knepley { 369f155c232SMatthew G. Knepley QuasiInterp *ctx; 370f155c232SMatthew G. Knepley DM dmcdg, dmc; 371f155c232SMatthew G. Knepley Vec ly; 372f155c232SMatthew G. Knepley 373f155c232SMatthew G. Knepley PetscFunctionBegin; 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(P, (void **) &ctx)); 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(ctx->Ifdg, x, ctx->tmpf)); 376f155c232SMatthew G. Knepley 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) ctx->tmpf, "Fine DG Potential")); 378*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpf, "fine_dg_")); 379*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(ctx->tmpf, NULL, "-vec_view")); 380*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(ctx->Pi, x, ctx->tmpc)); 381f155c232SMatthew G. Knepley 382*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) ctx->tmpc, "Coarse DG Potential")); 383*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpc, "coarse_dg_")); 384*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(ctx->tmpc, NULL, "-vec_view")); 385*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(ctx->tmpc, &dmcdg)); 386f155c232SMatthew G. Knepley 387*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(y, &dmc)); 388*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dmc, &ly)); 389*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeClementInterpolant(dmcdg, ctx->tmpc, ly)); 390*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(dmc, ly, INSERT_VALUES, y)); 391*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(dmc, ly, INSERT_VALUES, y)); 392*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dmc, &ly)); 393f155c232SMatthew G. Knepley PetscFunctionReturn(0); 394f155c232SMatthew G. Knepley } 395f155c232SMatthew G. Knepley 396f155c232SMatthew G. Knepley static PetscErrorCode CreateQuasiInterpolator(DM dmc, DM dmf, Mat *P) 397f155c232SMatthew G. Knepley { 398f155c232SMatthew G. Knepley QuasiInterp *ctx; 399f155c232SMatthew G. Knepley DM dmcdg, dmfdg; 400f155c232SMatthew G. Knepley PetscFE fe; 401f155c232SMatthew G. Knepley Vec x, y; 402f155c232SMatthew G. Knepley DMPolytopeType ct; 403f155c232SMatthew G. Knepley PetscInt dim, cStart, m, n, M, N; 404f155c232SMatthew G. Knepley 405f155c232SMatthew G. Knepley PetscFunctionBegin; 406*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(1, &ctx)); 407*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmc, &x)); 408*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmf, &y)); 409*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(x, &m)); 410*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x, &M)); 411*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(y, &n)); 412*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(y, &N)); 413*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmc, &x)); 414*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmf, &y)); 415f155c232SMatthew G. Knepley 416*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dmf, &dmfdg)); 417*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dmfdg, &dim)); 418*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dmfdg, 0, &cStart, NULL)); 419*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dmfdg, cStart, &ct)); 420*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "fine_dg_", PETSC_DETERMINE, &fe)); 421*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dmfdg, 0, NULL, (PetscObject) fe)); 422*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 423*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dmfdg)); 424*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(dmf, dmfdg, &ctx->Ifdg, NULL)); 425*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dmfdg, &ctx->tmpf)); 426*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmfdg)); 427f155c232SMatthew G. Knepley 428*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dmc, &dmcdg)); 429*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dmcdg, &dim)); 430*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dmcdg, 0, &cStart, NULL)); 431*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dmcdg, cStart, &ct)); 432*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "coarse_dg_", PETSC_DETERMINE, &fe)); 433*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dmcdg, 0, NULL, (PetscObject) fe)); 434*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 435*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dmcdg)); 436f155c232SMatthew G. Knepley 437*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateCoarseProjection(dmcdg, dmf, &ctx->Pi)); 438*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dmcdg, &ctx->tmpc)); 439*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmcdg)); 440f155c232SMatthew G. Knepley 441*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PetscObjectComm((PetscObject) dmc), m, n, M, N, ctx, P)); 442*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*P, MATOP_DESTROY, (void (*)(void)) DestroyQuasiInterpolator)); 443*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(*P, MATOP_MULT, (void (*)(void)) QuasiInterpolate)); 444f155c232SMatthew G. Knepley PetscFunctionReturn(0); 445f155c232SMatthew G. Knepley } 446f155c232SMatthew G. Knepley 447f155c232SMatthew G. Knepley static PetscErrorCode CoarseTest(DM dm, Vec u, AppCtx *user) 448f155c232SMatthew G. Knepley { 449f155c232SMatthew G. Knepley DM dmc; 450f155c232SMatthew G. Knepley Mat P; /* The quasi-interpolator to the coarse space */ 451f155c232SMatthew G. Knepley Vec uc; 452f155c232SMatthew G. Knepley 453f155c232SMatthew G. Knepley PetscFunctionBegin; 454f155c232SMatthew G. Knepley if (user->modType == MOD_CONSTANT) PetscFunctionReturn(0); 455*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PetscObjectComm((PetscObject) dm), &dmc)); 456*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dmc, DMPLEX)); 457*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) dmc, "coarse_")); 458*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(dmc, user)); 459*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dmc)); 460*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dmc, NULL, "-dm_view")); 461f155c232SMatthew G. Knepley 462*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dmc, "potential", SetupPrimalProblem, user)); 463*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dmc, &uc)); 464*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) uc, "potential")); 465*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) uc, "coarse_")); 466f155c232SMatthew G. Knepley 467*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateQuasiInterpolator(dmc, dm, &P)); 468f155c232SMatthew G. Knepley #if 1 469*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(P, u, uc)); 470f155c232SMatthew G. Knepley #else 471f155c232SMatthew G. Knepley { 472f155c232SMatthew G. Knepley Mat In; 473f155c232SMatthew G. Knepley Vec sc; 474f155c232SMatthew G. Knepley 475*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(dmc, dm, &In, &sc)); 476*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(In, u, uc)); 477*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(uc, sc, uc)); 478*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&In)); 479*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&sc)); 480f155c232SMatthew G. Knepley } 481f155c232SMatthew G. Knepley #endif 482*5f80ce2aSJacob Faibussowitsch CHKERRQ(CompareView(uc)); 483f155c232SMatthew G. Knepley 484*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 485*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&uc)); 486*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmc)); 487f155c232SMatthew G. Knepley PetscFunctionReturn(0); 488f155c232SMatthew G. Knepley } 489f155c232SMatthew G. Knepley 490f155c232SMatthew G. Knepley int main(int argc, char **argv) 491f155c232SMatthew G. Knepley { 492f155c232SMatthew G. Knepley DM dm; /* Problem specification */ 493f155c232SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 494f155c232SMatthew G. Knepley Vec u; /* Solutions */ 495f155c232SMatthew G. Knepley AppCtx user; /* User-defined work context */ 496f155c232SMatthew G. Knepley PetscErrorCode ierr; 497f155c232SMatthew G. Knepley 498f155c232SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 499*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 500*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupParameters(PETSC_COMM_WORLD, &user)); 501f155c232SMatthew G. Knepley /* Primal system */ 502*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 503*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 504*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes, dm)); 505*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 506*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 507*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(u, 0.0)); 508*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "potential")); 509*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 510*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 511*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSNESCheckFromOptions(snes, u)); 512*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes, NULL, u)); 513*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetSolution(snes, &u)); 514*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-potential_view")); 515*5f80ce2aSJacob Faibussowitsch CHKERRQ(CompareView(u)); 516f155c232SMatthew G. Knepley /* Looking at a coarse problem */ 517*5f80ce2aSJacob Faibussowitsch CHKERRQ(CoarseTest(dm, u, &user)); 518f155c232SMatthew G. Knepley /* Cleanup */ 519*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 520*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 521*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 522*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagDestroy(&user.bag)); 523f155c232SMatthew G. Knepley ierr = PetscFinalize(); 524f155c232SMatthew G. Knepley return ierr; 525f155c232SMatthew G. Knepley } 526f155c232SMatthew G. Knepley 527f155c232SMatthew G. Knepley /*TEST 528f155c232SMatthew G. Knepley 529f155c232SMatthew G. Knepley test: 530f155c232SMatthew G. Knepley suffix: 1d_p1_constant 531f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dmsnes_check 532f155c232SMatthew G. Knepley 533f155c232SMatthew G. Knepley test: 534f155c232SMatthew G. Knepley suffix: 1d_p1_constant_conv 535f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 \ 536f155c232SMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 2 537f155c232SMatthew G. Knepley 538f155c232SMatthew G. Knepley test: 539f155c232SMatthew G. Knepley suffix: 1d_p1_oscillatory 540f155c232SMatthew G. Knepley args: -mod_type oscillatory -epsilon 0.03125 \ 541f155c232SMatthew G. Knepley -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dm_refine 2 -dmsnes_check \ 542f155c232SMatthew G. Knepley -coarse_dm_plex_dim 1 -coarse_dm_plex_box_faces 4 -coarse_dm_plex_hash_location \ 543f155c232SMatthew G. Knepley -fine_dg_petscspace_degree 1 -fine_dg_petscdualspace_lagrange_continuity 0 \ 544f155c232SMatthew G. Knepley -coarse_dg_petscspace_degree 1 -coarse_dg_petscdualspace_lagrange_continuity 0 545f155c232SMatthew G. Knepley 546f155c232SMatthew G. Knepley TEST*/ 547