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 25*9371c9d4SSatish Balay typedef enum { 26*9371c9d4SSatish Balay MOD_CONSTANT, 27*9371c9d4SSatish Balay MOD_OSCILLATORY, 28*9371c9d4SSatish Balay MOD_TANH, 29*9371c9d4SSatish Balay NUM_MOD_TYPES 30*9371c9d4SSatish Balay } ModType; 31f155c232SMatthew G. Knepley const char *modTypes[NUM_MOD_TYPES + 1] = {"constant", "oscillatory", "tanh", "unknown"}; 32f155c232SMatthew G. Knepley 33f155c232SMatthew G. Knepley /* Constants */ 34*9371c9d4SSatish Balay enum { 35*9371c9d4SSatish Balay EPSILON, 36*9371c9d4SSatish Balay NUM_CONSTANTS 37*9371c9d4SSatish 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*9371c9d4SSatish Balay static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 49f155c232SMatthew G. Knepley PetscInt d; 50f155c232SMatthew G. Knepley *u = 1.0; 51f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0 * PETSC_PI * x[d]); 52f155c232SMatthew G. Knepley return 0; 53f155c232SMatthew G. Knepley } 54f155c232SMatthew G. Knepley 55*9371c9d4SSatish Balay static PetscErrorCode oscillatory_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 56f155c232SMatthew G. Knepley Parameter *param = (Parameter *)ctx; 57f155c232SMatthew G. Knepley const PetscReal eps = param->epsilon; 58f155c232SMatthew G. Knepley 59f155c232SMatthew 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)); 60f155c232SMatthew G. Knepley return 0; 61f155c232SMatthew G. Knepley } 62f155c232SMatthew G. Knepley 63*9371c9d4SSatish Balay 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[]) { 64f155c232SMatthew G. Knepley PetscInt d; 65f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) { 66f155c232SMatthew G. Knepley PetscScalar v = 1.; 67f155c232SMatthew G. Knepley for (PetscInt e = 0; e < dim; e++) { 68f155c232SMatthew G. Knepley if (e == d) { 69f155c232SMatthew G. Knepley v *= -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 70f155c232SMatthew G. Knepley } else { 71f155c232SMatthew G. Knepley v *= PetscSinReal(2.0 * PETSC_PI * x[d]); 72f155c232SMatthew G. Knepley } 73f155c232SMatthew G. Knepley } 74f155c232SMatthew G. Knepley f0[0] += v; 75f155c232SMatthew G. Knepley } 76f155c232SMatthew G. Knepley } 77f155c232SMatthew G. Knepley 78*9371c9d4SSatish Balay 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[]) { 79f155c232SMatthew G. Knepley PetscInt d; 80f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 81f155c232SMatthew G. Knepley } 82f155c232SMatthew G. Knepley 83*9371c9d4SSatish Balay 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[]) { 84f155c232SMatthew G. Knepley PetscInt d; 85f155c232SMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 86f155c232SMatthew G. Knepley } 87f155c232SMatthew G. Knepley 88*9371c9d4SSatish Balay 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[]) { 89f155c232SMatthew G. Knepley f0[0] = -1.; 90f155c232SMatthew G. Knepley } 91f155c232SMatthew G. Knepley 92*9371c9d4SSatish Balay 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[]) { 93f155c232SMatthew G. Knepley const PetscReal eps = PetscRealPart(constants[EPSILON]); 94f155c232SMatthew G. Knepley 95f155c232SMatthew G. Knepley f1[0] = u_x[0] / (2. + PetscCosReal(2. * PETSC_PI * x[0] / eps)); 96f155c232SMatthew G. Knepley } 97f155c232SMatthew G. Knepley 98*9371c9d4SSatish Balay 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[]) { 99f155c232SMatthew G. Knepley const PetscReal eps = PetscRealPart(constants[EPSILON]); 100f155c232SMatthew G. Knepley 101f155c232SMatthew G. Knepley g3[0] = 1. / (2. + PetscCosReal(2. * PETSC_PI * x[0] / eps)); 102f155c232SMatthew G. Knepley } 103f155c232SMatthew G. Knepley 104*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 105f155c232SMatthew G. Knepley PetscInt mod; 106f155c232SMatthew G. Knepley 107f155c232SMatthew G. Knepley PetscFunctionBeginUser; 108f155c232SMatthew G. Knepley options->modType = MOD_CONSTANT; 109d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Homogenization Problem Options", "DMPLEX"); 110f155c232SMatthew G. Knepley mod = options->modType; 1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mod_type", "The model type", "ex36.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL)); 112f155c232SMatthew G. Knepley options->modType = (ModType)mod; 113d0609cedSBarry Smith PetscOptionsEnd(); 114f155c232SMatthew G. Knepley PetscFunctionReturn(0); 115f155c232SMatthew G. Knepley } 116f155c232SMatthew G. Knepley 117*9371c9d4SSatish Balay static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *user) { 118f155c232SMatthew G. Knepley PetscBag bag; 119f155c232SMatthew G. Knepley Parameter *p; 120f155c232SMatthew G. Knepley 121f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1229566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(comm, sizeof(Parameter), &user->bag)); 1239566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&p)); 1249566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(user->bag, "par", "Homogenization parameters")); 125f155c232SMatthew G. Knepley bag = user->bag; 1269566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->epsilon, 1.0, "epsilon", "Wavelength of fine scale oscillation")); 127f155c232SMatthew G. Knepley PetscFunctionReturn(0); 128f155c232SMatthew G. Knepley } 129f155c232SMatthew G. Knepley 130*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 131f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1359566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 1369566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 137f155c232SMatthew G. Knepley PetscFunctionReturn(0); 138f155c232SMatthew G. Knepley } 139f155c232SMatthew G. Knepley 140*9371c9d4SSatish Balay static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) { 141f155c232SMatthew G. Knepley PetscDS ds; 142f155c232SMatthew G. Knepley DMLabel label; 143f155c232SMatthew G. Knepley PetscSimplePointFunc ex; 144f155c232SMatthew G. Knepley const PetscInt id = 1; 145f155c232SMatthew G. Knepley void *ctx; 146f155c232SMatthew G. Knepley 147f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1489566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1499566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)&ctx)); 150f155c232SMatthew G. Knepley switch (user->modType) { 151f155c232SMatthew G. Knepley case MOD_CONSTANT: 1529566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_homogeneous_u, f1_u)); 1539566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1549566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 155f155c232SMatthew G. Knepley ex = trig_homogeneous_u; 1569566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))ex, NULL, ctx, NULL)); 157f155c232SMatthew G. Knepley break; 158f155c232SMatthew G. Knepley case MOD_OSCILLATORY: 1599566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_oscillatory_u, f1_oscillatory_u)); 1609566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_oscillatory_uu)); 1619566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 162f155c232SMatthew G. Knepley ex = oscillatory_u; 1639566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))ex, NULL, ctx, NULL)); 164f155c232SMatthew G. Knepley break; 16563a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Unsupported model type: %s (%d)", modTypes[PetscMin(user->modType, NUM_MOD_TYPES)], user->modType); 166f155c232SMatthew G. Knepley } 1679566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, ex, ctx)); 168f155c232SMatthew G. Knepley /* Setup constants */ 169f155c232SMatthew G. Knepley { 170f155c232SMatthew G. Knepley Parameter *param; 171f155c232SMatthew G. Knepley PetscScalar constants[NUM_CONSTANTS]; 172f155c232SMatthew G. Knepley 1739566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 174f155c232SMatthew G. Knepley 175f155c232SMatthew G. Knepley constants[EPSILON] = param->epsilon; 1769566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants)); 177f155c232SMatthew G. Knepley } 178f155c232SMatthew G. Knepley PetscFunctionReturn(0); 179f155c232SMatthew G. Knepley } 180f155c232SMatthew G. Knepley 181*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) { 182f155c232SMatthew G. Knepley DM cdm = dm; 183f155c232SMatthew G. Knepley PetscFE fe; 184f155c232SMatthew G. Knepley PetscBool simplex; 185f155c232SMatthew G. Knepley PetscInt dim; 186f155c232SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 187f155c232SMatthew G. Knepley 188f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1899566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1909566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 191f155c232SMatthew G. Knepley /* Create finite element */ 1929566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 1939566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 1949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name)); 195f155c232SMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 1969566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 1979566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1989566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 199f155c232SMatthew G. Knepley while (cdm) { 2009566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2019566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 202f155c232SMatthew G. Knepley } 2039566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 204f155c232SMatthew G. Knepley PetscFunctionReturn(0); 205f155c232SMatthew G. Knepley } 206f155c232SMatthew G. Knepley 207*9371c9d4SSatish Balay static PetscErrorCode CompareView(Vec u) { 208f155c232SMatthew G. Knepley DM dm; 209f155c232SMatthew G. Knepley Vec v[2], lv[2], exact; 210f155c232SMatthew G. Knepley PetscOptions options; 211f155c232SMatthew G. Knepley PetscViewer viewer; 212f155c232SMatthew G. Knepley PetscViewerFormat format; 213f155c232SMatthew G. Knepley PetscBool flg; 214f155c232SMatthew G. Knepley PetscInt i; 215f155c232SMatthew G. Knepley const char *name, *prefix; 216f155c232SMatthew G. Knepley 217f155c232SMatthew G. Knepley PetscFunctionBeginUser; 2189566063dSJacob Faibussowitsch PetscCall(VecGetDM(u, &dm)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptions((PetscObject)dm, &options)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)dm), options, prefix, "-compare_view", &viewer, &format, &flg)); 222f155c232SMatthew G. Knepley if (flg) { 2239566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &exact)); 2249566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, 0.0, exact, NULL)); 225f155c232SMatthew G. Knepley v[0] = u; 226f155c232SMatthew G. Knepley v[1] = exact; 227f155c232SMatthew G. Knepley for (i = 0; i < 2; ++i) { 2289566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lv[i])); 2299566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)v[i], &name)); 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)lv[i], name)); 2319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, v[i], INSERT_VALUES, lv[i])); 2329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, v[i], INSERT_VALUES, lv[i])); 2339566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, lv[i], 0., NULL, NULL, NULL)); 234f155c232SMatthew G. Knepley } 2359566063dSJacob Faibussowitsch PetscCall(DMPlexVecView1D(dm, 2, lv, viewer)); 2369566063dSJacob Faibussowitsch for (i = 0; i < 2; ++i) PetscCall(DMRestoreLocalVector(dm, &lv[i])); 2379566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &exact)); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 239f155c232SMatthew G. Knepley } 240f155c232SMatthew G. Knepley PetscFunctionReturn(0); 241f155c232SMatthew G. Knepley } 242f155c232SMatthew G. Knepley 243*9371c9d4SSatish Balay typedef struct { 244f155c232SMatthew G. Knepley Mat Mcoarse; /* Mass matrix on the coarse space */ 245f155c232SMatthew G. Knepley Mat Mfine; /* Mass matrix on the fine space */ 246f155c232SMatthew G. Knepley Mat Ifine; /* Interpolator from coarse to fine */ 247f155c232SMatthew G. Knepley Vec Iscale; /* Scaling vector for restriction */ 248f155c232SMatthew G. Knepley KSP kspCoarse; /* Solver for the coarse mass matrix */ 249f155c232SMatthew G. Knepley Vec tmpfine; /* Temporary vector in the fine space */ 250f155c232SMatthew G. Knepley Vec tmpcoarse; /* Temporary vector in the coarse space */ 251f155c232SMatthew G. Knepley } ProjStruct; 252f155c232SMatthew G. Knepley 253*9371c9d4SSatish Balay static PetscErrorCode DestroyCoarseProjection(Mat Pi) { 254f155c232SMatthew G. Knepley ProjStruct *ctx; 255f155c232SMatthew G. Knepley 256f155c232SMatthew G. Knepley PetscFunctionBegin; 2579566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(Pi, (void **)&ctx)); 2589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Mcoarse)); 2599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Mfine)); 2609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Ifine)); 2619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->Iscale)); 2629566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ctx->kspCoarse)); 2639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpcoarse)); 2649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpfine)); 2659566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 2669566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(Pi, NULL)); 267f155c232SMatthew G. Knepley PetscFunctionReturn(0); 268f155c232SMatthew G. Knepley } 269f155c232SMatthew G. Knepley 270*9371c9d4SSatish Balay static PetscErrorCode CoarseProjection(Mat Pi, Vec x, Vec y) { 271f155c232SMatthew G. Knepley ProjStruct *ctx; 272f155c232SMatthew G. Knepley 273f155c232SMatthew G. Knepley PetscFunctionBegin; 2749566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(Pi, (void **)&ctx)); 2759566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Mfine, x, ctx->tmpfine)); 2769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpfine, "Fine DG RHS")); 2779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpfine, "fine_dg_")); 2789566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpfine, NULL, "-rhs_view")); 2799566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->Ifine, ctx->tmpfine, ctx->tmpcoarse)); 2809566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpcoarse, "Coarse DG RHS")); 2819566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpcoarse, "coarse_dg_")); 2829566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 2839566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ctx->tmpcoarse, ctx->Iscale, ctx->tmpcoarse)); 2849566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 2859566063dSJacob Faibussowitsch PetscCall(KSPSolve(ctx->kspCoarse, ctx->tmpcoarse, y)); 286f155c232SMatthew G. Knepley PetscFunctionReturn(0); 287f155c232SMatthew G. Knepley } 288f155c232SMatthew G. Knepley 289*9371c9d4SSatish Balay static PetscErrorCode CreateCoarseProjection(DM dmc, DM dmf, Mat *Pi) { 290f155c232SMatthew G. Knepley ProjStruct *ctx; 291f155c232SMatthew G. Knepley PetscInt m, n, M, N; 292f155c232SMatthew G. Knepley 293f155c232SMatthew G. Knepley PetscFunctionBegin; 2949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &ctx)); 2959566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmc, &ctx->tmpcoarse)); 2969566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmf, &ctx->tmpfine)); 2979566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->tmpcoarse, &m)); 2989566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->tmpcoarse, &M)); 2999566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->tmpfine, &n)); 3009566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->tmpfine, &N)); 3019566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmc, dmc, &ctx->Mcoarse)); 3029566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmf, dmf, &ctx->Mfine)); 3039566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmc, dmf, &ctx->Ifine, &ctx->Iscale)); 3049566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)dmc), &ctx->kspCoarse)); 3059566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->kspCoarse, "coarse_")); 3069566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ctx->kspCoarse, ctx->Mcoarse, ctx->Mcoarse)); 3079566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ctx->kspCoarse)); 3089566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dmc), m, n, M, N, ctx, Pi)); 3099566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*Pi, MATOP_DESTROY, (void (*)(void))DestroyCoarseProjection)); 3109566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*Pi, MATOP_MULT, (void (*)(void))CoarseProjection)); 311f155c232SMatthew G. Knepley PetscFunctionReturn(0); 312f155c232SMatthew G. Knepley } 313f155c232SMatthew G. Knepley 314*9371c9d4SSatish Balay typedef struct { 315f155c232SMatthew G. Knepley Mat Ifdg; /* Embed the fine space into the DG version */ 316f155c232SMatthew G. Knepley Mat Pi; /* The L_2 stable projection to the DG coarse space */ 317f155c232SMatthew G. Knepley Vec tmpc; /* A temporary vector in the DG coarse space */ 318f155c232SMatthew G. Knepley Vec tmpf; /* A temporary vector in the DG fine space */ 319f155c232SMatthew G. Knepley } QuasiInterp; 320f155c232SMatthew G. Knepley 321*9371c9d4SSatish Balay static PetscErrorCode DestroyQuasiInterpolator(Mat P) { 322f155c232SMatthew G. Knepley QuasiInterp *ctx; 323f155c232SMatthew G. Knepley 324f155c232SMatthew G. Knepley PetscFunctionBegin; 3259566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(P, (void **)&ctx)); 3269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Ifdg)); 3279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Pi)); 3289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpc)); 3299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpf)); 3309566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 3319566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(P, NULL)); 332f155c232SMatthew G. Knepley PetscFunctionReturn(0); 333f155c232SMatthew G. Knepley } 334f155c232SMatthew G. Knepley 335*9371c9d4SSatish Balay static PetscErrorCode QuasiInterpolate(Mat P, Vec x, Vec y) { 336f155c232SMatthew G. Knepley QuasiInterp *ctx; 337f155c232SMatthew G. Knepley DM dmcdg, dmc; 338f155c232SMatthew G. Knepley Vec ly; 339f155c232SMatthew G. Knepley 340f155c232SMatthew G. Knepley PetscFunctionBegin; 3419566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(P, (void **)&ctx)); 3429566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Ifdg, x, ctx->tmpf)); 343f155c232SMatthew G. Knepley 3449566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpf, "Fine DG Potential")); 3459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpf, "fine_dg_")); 3469566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpf, NULL, "-vec_view")); 3479566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Pi, x, ctx->tmpc)); 348f155c232SMatthew G. Knepley 3499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ctx->tmpc, "Coarse DG Potential")); 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->tmpc, "coarse_dg_")); 3519566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpc, NULL, "-vec_view")); 3529566063dSJacob Faibussowitsch PetscCall(VecGetDM(ctx->tmpc, &dmcdg)); 353f155c232SMatthew G. Knepley 3549566063dSJacob Faibussowitsch PetscCall(VecGetDM(y, &dmc)); 3559566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmc, &ly)); 3569566063dSJacob Faibussowitsch PetscCall(DMPlexComputeClementInterpolant(dmcdg, ctx->tmpc, ly)); 3579566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dmc, ly, INSERT_VALUES, y)); 3589566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dmc, ly, INSERT_VALUES, y)); 3599566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmc, &ly)); 360f155c232SMatthew G. Knepley PetscFunctionReturn(0); 361f155c232SMatthew G. Knepley } 362f155c232SMatthew G. Knepley 363*9371c9d4SSatish Balay static PetscErrorCode CreateQuasiInterpolator(DM dmc, DM dmf, Mat *P) { 364f155c232SMatthew G. Knepley QuasiInterp *ctx; 365f155c232SMatthew G. Knepley DM dmcdg, dmfdg; 366f155c232SMatthew G. Knepley PetscFE fe; 367f155c232SMatthew G. Knepley Vec x, y; 368f155c232SMatthew G. Knepley DMPolytopeType ct; 369f155c232SMatthew G. Knepley PetscInt dim, cStart, m, n, M, N; 370f155c232SMatthew G. Knepley 371f155c232SMatthew G. Knepley PetscFunctionBegin; 3729566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(1, &ctx)); 3739566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmc, &x)); 3749566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmf, &y)); 3759566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &m)); 3769566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &M)); 3779566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 3789566063dSJacob Faibussowitsch PetscCall(VecGetSize(y, &N)); 3799566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmc, &x)); 3809566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmf, &y)); 381f155c232SMatthew G. Knepley 3829566063dSJacob Faibussowitsch PetscCall(DMClone(dmf, &dmfdg)); 3839566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmfdg, &dim)); 3849566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmfdg, 0, &cStart, NULL)); 3859566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dmfdg, cStart, &ct)); 3869566063dSJacob Faibussowitsch PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "fine_dg_", PETSC_DETERMINE, &fe)); 3879566063dSJacob Faibussowitsch PetscCall(DMSetField(dmfdg, 0, NULL, (PetscObject)fe)); 3889566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 3899566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmfdg)); 3909566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmf, dmfdg, &ctx->Ifdg, NULL)); 3919566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmfdg, &ctx->tmpf)); 3929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmfdg)); 393f155c232SMatthew G. Knepley 3949566063dSJacob Faibussowitsch PetscCall(DMClone(dmc, &dmcdg)); 3959566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmcdg, &dim)); 3969566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmcdg, 0, &cStart, NULL)); 3979566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dmcdg, cStart, &ct)); 3989566063dSJacob Faibussowitsch PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "coarse_dg_", PETSC_DETERMINE, &fe)); 3999566063dSJacob Faibussowitsch PetscCall(DMSetField(dmcdg, 0, NULL, (PetscObject)fe)); 4009566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 4019566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmcdg)); 402f155c232SMatthew G. Knepley 4039566063dSJacob Faibussowitsch PetscCall(CreateCoarseProjection(dmcdg, dmf, &ctx->Pi)); 4049566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmcdg, &ctx->tmpc)); 4059566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmcdg)); 406f155c232SMatthew G. Knepley 4079566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dmc), m, n, M, N, ctx, P)); 4089566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*P, MATOP_DESTROY, (void (*)(void))DestroyQuasiInterpolator)); 4099566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*P, MATOP_MULT, (void (*)(void))QuasiInterpolate)); 410f155c232SMatthew G. Knepley PetscFunctionReturn(0); 411f155c232SMatthew G. Knepley } 412f155c232SMatthew G. Knepley 413*9371c9d4SSatish Balay static PetscErrorCode CoarseTest(DM dm, Vec u, AppCtx *user) { 414f155c232SMatthew G. Knepley DM dmc; 415f155c232SMatthew G. Knepley Mat P; /* The quasi-interpolator to the coarse space */ 416f155c232SMatthew G. Knepley Vec uc; 417f155c232SMatthew G. Knepley 418f155c232SMatthew G. Knepley PetscFunctionBegin; 419f155c232SMatthew G. Knepley if (user->modType == MOD_CONSTANT) PetscFunctionReturn(0); 4209566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &dmc)); 4219566063dSJacob Faibussowitsch PetscCall(DMSetType(dmc, DMPLEX)); 4229566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dmc, "coarse_")); 4239566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dmc, user)); 4249566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dmc)); 4259566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmc, NULL, "-dm_view")); 426f155c232SMatthew G. Knepley 4279566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dmc, "potential", SetupPrimalProblem, user)); 4289566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmc, &uc)); 4299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)uc, "potential")); 4309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)uc, "coarse_")); 431f155c232SMatthew G. Knepley 4329566063dSJacob Faibussowitsch PetscCall(CreateQuasiInterpolator(dmc, dm, &P)); 433f155c232SMatthew G. Knepley #if 1 4349566063dSJacob Faibussowitsch PetscCall(MatMult(P, u, uc)); 435f155c232SMatthew G. Knepley #else 436f155c232SMatthew G. Knepley { 437f155c232SMatthew G. Knepley Mat In; 438f155c232SMatthew G. Knepley Vec sc; 439f155c232SMatthew G. Knepley 4409566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmc, dm, &In, &sc)); 4419566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(In, u, uc)); 4429566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(uc, sc, uc)); 4439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&In)); 4449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sc)); 445f155c232SMatthew G. Knepley } 446f155c232SMatthew G. Knepley #endif 4479566063dSJacob Faibussowitsch PetscCall(CompareView(uc)); 448f155c232SMatthew G. Knepley 4499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 4509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uc)); 4519566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmc)); 452f155c232SMatthew G. Knepley PetscFunctionReturn(0); 453f155c232SMatthew G. Knepley } 454f155c232SMatthew G. Knepley 455*9371c9d4SSatish Balay int main(int argc, char **argv) { 456f155c232SMatthew G. Knepley DM dm; /* Problem specification */ 457f155c232SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 458f155c232SMatthew G. Knepley Vec u; /* Solutions */ 459f155c232SMatthew G. Knepley AppCtx user; /* User-defined work context */ 460f155c232SMatthew G. Knepley 461327415f7SBarry Smith PetscFunctionBeginUser; 4629566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4639566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4649566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 465f155c232SMatthew G. Knepley /* Primal system */ 4669566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 4679566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4689566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 4699566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 4709566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 4719566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 4729566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "potential")); 4739566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 4749566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 4759566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 4769566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 4779566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 4789566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 4799566063dSJacob Faibussowitsch PetscCall(CompareView(u)); 480f155c232SMatthew G. Knepley /* Looking at a coarse problem */ 4819566063dSJacob Faibussowitsch PetscCall(CoarseTest(dm, u, &user)); 482f155c232SMatthew G. Knepley /* Cleanup */ 4839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 4849566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4859566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4869566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 4879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 488b122ec5aSJacob Faibussowitsch return 0; 489f155c232SMatthew G. Knepley } 490f155c232SMatthew G. Knepley 491f155c232SMatthew G. Knepley /*TEST 492f155c232SMatthew G. Knepley 493f155c232SMatthew G. Knepley test: 494f155c232SMatthew G. Knepley suffix: 1d_p1_constant 495f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dmsnes_check 496f155c232SMatthew G. Knepley 497f155c232SMatthew G. Knepley test: 498f155c232SMatthew G. Knepley suffix: 1d_p1_constant_conv 499f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 \ 500f155c232SMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 2 501f155c232SMatthew G. Knepley 502f155c232SMatthew G. Knepley test: 503f155c232SMatthew G. Knepley suffix: 1d_p1_oscillatory 504f155c232SMatthew G. Knepley args: -mod_type oscillatory -epsilon 0.03125 \ 505f155c232SMatthew G. Knepley -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dm_refine 2 -dmsnes_check \ 506f155c232SMatthew G. Knepley -coarse_dm_plex_dim 1 -coarse_dm_plex_box_faces 4 -coarse_dm_plex_hash_location \ 507f155c232SMatthew G. Knepley -fine_dg_petscspace_degree 1 -fine_dg_petscdualspace_lagrange_continuity 0 \ 508f155c232SMatthew G. Knepley -coarse_dg_petscspace_degree 1 -coarse_dg_petscdualspace_lagrange_continuity 0 509f155c232SMatthew G. Knepley 510f155c232SMatthew G. Knepley TEST*/ 511