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 126f155c232SMatthew G. Knepley PetscFunctionBeginUser; 127f155c232SMatthew G. Knepley options->modType = MOD_CONSTANT; 128*d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Homogenization Problem Options", "DMPLEX"); 129f155c232SMatthew G. Knepley mod = options->modType; 1309566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mod_type", "The model type", "ex36.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL)); 131f155c232SMatthew G. Knepley options->modType = (ModType) mod; 132*d0609cedSBarry Smith PetscOptionsEnd(); 133f155c232SMatthew G. Knepley PetscFunctionReturn(0); 134f155c232SMatthew G. Knepley } 135f155c232SMatthew G. Knepley 136f155c232SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *user) 137f155c232SMatthew G. Knepley { 138f155c232SMatthew G. Knepley PetscBag bag; 139f155c232SMatthew G. Knepley Parameter *p; 140f155c232SMatthew G. Knepley 141f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1429566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(comm, sizeof(Parameter), &user->bag)); 1439566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **) &p)); 1449566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(user->bag, "par", "Homogenization parameters")); 145f155c232SMatthew G. Knepley bag = user->bag; 1469566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterReal(bag, &p->epsilon, 1.0, "epsilon", "Wavelength of fine scale oscillation")); 147f155c232SMatthew G. Knepley PetscFunctionReturn(0); 148f155c232SMatthew G. Knepley } 149f155c232SMatthew G. Knepley 150f155c232SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 151f155c232SMatthew G. Knepley { 152f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1539566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1549566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1559566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1569566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 1579566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 158f155c232SMatthew G. Knepley PetscFunctionReturn(0); 159f155c232SMatthew G. Knepley } 160f155c232SMatthew G. Knepley 161f155c232SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 162f155c232SMatthew G. Knepley { 163f155c232SMatthew G. Knepley PetscDS ds; 164f155c232SMatthew G. Knepley DMLabel label; 165f155c232SMatthew G. Knepley PetscSimplePointFunc ex; 166f155c232SMatthew G. Knepley const PetscInt id = 1; 167f155c232SMatthew G. Knepley void *ctx; 168f155c232SMatthew G. Knepley 169f155c232SMatthew G. Knepley PetscFunctionBeginUser; 1709566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1719566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **) &ctx)); 172f155c232SMatthew G. Knepley switch (user->modType) { 173f155c232SMatthew G. Knepley case MOD_CONSTANT: 1749566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_homogeneous_u, f1_u)); 1759566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1769566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 177f155c232SMatthew G. Knepley ex = trig_homogeneous_u; 1789566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, ctx, NULL)); 179f155c232SMatthew G. Knepley break; 180f155c232SMatthew G. Knepley case MOD_OSCILLATORY: 1819566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_oscillatory_u, f1_oscillatory_u)); 1829566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_oscillatory_uu)); 1839566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 184f155c232SMatthew G. Knepley ex = oscillatory_u; 1859566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, ctx, NULL)); 186f155c232SMatthew G. Knepley break; 187f155c232SMatthew 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); 188f155c232SMatthew G. Knepley } 1899566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, ex, ctx)); 190f155c232SMatthew G. Knepley /* Setup constants */ 191f155c232SMatthew G. Knepley { 192f155c232SMatthew G. Knepley Parameter *param; 193f155c232SMatthew G. Knepley PetscScalar constants[NUM_CONSTANTS]; 194f155c232SMatthew G. Knepley 1959566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 196f155c232SMatthew G. Knepley 197f155c232SMatthew G. Knepley constants[EPSILON] = param->epsilon; 1989566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants)); 199f155c232SMatthew G. Knepley } 200f155c232SMatthew G. Knepley PetscFunctionReturn(0); 201f155c232SMatthew G. Knepley } 202f155c232SMatthew G. Knepley 203f155c232SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 204f155c232SMatthew G. Knepley { 205f155c232SMatthew G. Knepley DM cdm = dm; 206f155c232SMatthew G. Knepley PetscFE fe; 207f155c232SMatthew G. Knepley PetscBool simplex; 208f155c232SMatthew G. Knepley PetscInt dim; 209f155c232SMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 210f155c232SMatthew G. Knepley 211f155c232SMatthew G. Knepley PetscFunctionBeginUser; 2129566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2139566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 214f155c232SMatthew G. Knepley /* Create finite element */ 2159566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 2169566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 218f155c232SMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 2199566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 2209566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2219566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 222f155c232SMatthew G. Knepley while (cdm) { 2239566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 2249566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 225f155c232SMatthew G. Knepley } 2269566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 227f155c232SMatthew G. Knepley PetscFunctionReturn(0); 228f155c232SMatthew G. Knepley } 229f155c232SMatthew G. Knepley 230f155c232SMatthew G. Knepley static PetscErrorCode CompareView(Vec u) 231f155c232SMatthew G. Knepley { 232f155c232SMatthew G. Knepley DM dm; 233f155c232SMatthew G. Knepley Vec v[2], lv[2], exact; 234f155c232SMatthew G. Knepley PetscOptions options; 235f155c232SMatthew G. Knepley PetscViewer viewer; 236f155c232SMatthew G. Knepley PetscViewerFormat format; 237f155c232SMatthew G. Knepley PetscBool flg; 238f155c232SMatthew G. Knepley PetscInt i; 239f155c232SMatthew G. Knepley const char *name, *prefix; 240f155c232SMatthew G. Knepley 241f155c232SMatthew G. Knepley PetscFunctionBeginUser; 2429566063dSJacob Faibussowitsch PetscCall(VecGetDM(u, &dm)); 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptions((PetscObject) dm, &options)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix)); 2459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject) dm), options, prefix, "-compare_view", &viewer, &format, &flg)); 246f155c232SMatthew G. Knepley if (flg) { 2479566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &exact)); 2489566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, 0.0, exact, NULL)); 249f155c232SMatthew G. Knepley v[0] = u; 250f155c232SMatthew G. Knepley v[1] = exact; 251f155c232SMatthew G. Knepley for (i = 0; i < 2; ++i) { 2529566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lv[i])); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) v[i], &name)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) lv[i], name)); 2559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, v[i], INSERT_VALUES, lv[i])); 2569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, v[i], INSERT_VALUES, lv[i])); 2579566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, lv[i], 0., NULL, NULL, NULL)); 258f155c232SMatthew G. Knepley } 2599566063dSJacob Faibussowitsch PetscCall(DMPlexVecView1D(dm, 2, lv, viewer)); 2609566063dSJacob Faibussowitsch for (i = 0; i < 2; ++i) PetscCall(DMRestoreLocalVector(dm, &lv[i])); 2619566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &exact)); 2629566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 263f155c232SMatthew G. Knepley } 264f155c232SMatthew G. Knepley PetscFunctionReturn(0); 265f155c232SMatthew G. Knepley } 266f155c232SMatthew G. Knepley 267f155c232SMatthew G. Knepley typedef struct 268f155c232SMatthew G. Knepley { 269f155c232SMatthew G. Knepley Mat Mcoarse; /* Mass matrix on the coarse space */ 270f155c232SMatthew G. Knepley Mat Mfine; /* Mass matrix on the fine space */ 271f155c232SMatthew G. Knepley Mat Ifine; /* Interpolator from coarse to fine */ 272f155c232SMatthew G. Knepley Vec Iscale; /* Scaling vector for restriction */ 273f155c232SMatthew G. Knepley KSP kspCoarse; /* Solver for the coarse mass matrix */ 274f155c232SMatthew G. Knepley Vec tmpfine; /* Temporary vector in the fine space */ 275f155c232SMatthew G. Knepley Vec tmpcoarse; /* Temporary vector in the coarse space */ 276f155c232SMatthew G. Knepley } ProjStruct; 277f155c232SMatthew G. Knepley 278f155c232SMatthew G. Knepley static PetscErrorCode DestroyCoarseProjection(Mat Pi) 279f155c232SMatthew G. Knepley { 280f155c232SMatthew G. Knepley ProjStruct *ctx; 281f155c232SMatthew G. Knepley 282f155c232SMatthew G. Knepley PetscFunctionBegin; 2839566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(Pi, (void **) &ctx)); 2849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Mcoarse)); 2859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Mfine)); 2869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Ifine)); 2879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->Iscale)); 2889566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ctx->kspCoarse)); 2899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpcoarse)); 2909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpfine)); 2919566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 2929566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(Pi, NULL)); 293f155c232SMatthew G. Knepley PetscFunctionReturn(0); 294f155c232SMatthew G. Knepley } 295f155c232SMatthew G. Knepley 296f155c232SMatthew G. Knepley static PetscErrorCode CoarseProjection(Mat Pi, Vec x, Vec y) 297f155c232SMatthew G. Knepley { 298f155c232SMatthew G. Knepley ProjStruct *ctx; 299f155c232SMatthew G. Knepley 300f155c232SMatthew G. Knepley PetscFunctionBegin; 3019566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(Pi, (void **) &ctx)); 3029566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Mfine, x, ctx->tmpfine)); 3039566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) ctx->tmpfine, "Fine DG RHS")); 3049566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpfine, "fine_dg_")); 3059566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpfine, NULL, "-rhs_view")); 3069566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->Ifine, ctx->tmpfine, ctx->tmpcoarse)); 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) ctx->tmpcoarse, "Coarse DG RHS")); 3089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpcoarse, "coarse_dg_")); 3099566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 3109566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(ctx->tmpcoarse, ctx->Iscale, ctx->tmpcoarse)); 3119566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpcoarse, NULL, "-rhs_view")); 3129566063dSJacob Faibussowitsch PetscCall(KSPSolve(ctx->kspCoarse, ctx->tmpcoarse, y)); 313f155c232SMatthew G. Knepley PetscFunctionReturn(0); 314f155c232SMatthew G. Knepley } 315f155c232SMatthew G. Knepley 316f155c232SMatthew G. Knepley static PetscErrorCode CreateCoarseProjection(DM dmc, DM dmf, Mat *Pi) 317f155c232SMatthew G. Knepley { 318f155c232SMatthew G. Knepley ProjStruct *ctx; 319f155c232SMatthew G. Knepley PetscInt m, n, M, N; 320f155c232SMatthew G. Knepley 321f155c232SMatthew G. Knepley PetscFunctionBegin; 3229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &ctx)); 3239566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmc, &ctx->tmpcoarse)); 3249566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmf, &ctx->tmpfine)); 3259566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->tmpcoarse, &m)); 3269566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->tmpcoarse, &M)); 3279566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ctx->tmpfine, &n)); 3289566063dSJacob Faibussowitsch PetscCall(VecGetSize(ctx->tmpfine, &N)); 3299566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmc, dmc, &ctx->Mcoarse)); 3309566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmf, dmf, &ctx->Mfine)); 3319566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmc, dmf, &ctx->Ifine, &ctx->Iscale)); 3329566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject) dmc), &ctx->kspCoarse)); 3339566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->kspCoarse, "coarse_")); 3349566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ctx->kspCoarse, ctx->Mcoarse, ctx->Mcoarse)); 3359566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ctx->kspCoarse)); 3369566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject) dmc), m, n, M, N, ctx, Pi)); 3379566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*Pi, MATOP_DESTROY, (void (*)(void)) DestroyCoarseProjection)); 3389566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*Pi, MATOP_MULT, (void (*)(void)) CoarseProjection)); 339f155c232SMatthew G. Knepley PetscFunctionReturn(0); 340f155c232SMatthew G. Knepley } 341f155c232SMatthew G. Knepley 342f155c232SMatthew G. Knepley typedef struct 343f155c232SMatthew G. Knepley { 344f155c232SMatthew G. Knepley Mat Ifdg; /* Embed the fine space into the DG version */ 345f155c232SMatthew G. Knepley Mat Pi; /* The L_2 stable projection to the DG coarse space */ 346f155c232SMatthew G. Knepley Vec tmpc; /* A temporary vector in the DG coarse space */ 347f155c232SMatthew G. Knepley Vec tmpf; /* A temporary vector in the DG fine space */ 348f155c232SMatthew G. Knepley } QuasiInterp; 349f155c232SMatthew G. Knepley 350f155c232SMatthew G. Knepley static PetscErrorCode DestroyQuasiInterpolator(Mat P) 351f155c232SMatthew G. Knepley { 352f155c232SMatthew G. Knepley QuasiInterp *ctx; 353f155c232SMatthew G. Knepley 354f155c232SMatthew G. Knepley PetscFunctionBegin; 3559566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(P, (void **) &ctx)); 3569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Ifdg)); 3579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Pi)); 3589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpc)); 3599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->tmpf)); 3609566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 3619566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(P, NULL)); 362f155c232SMatthew G. Knepley PetscFunctionReturn(0); 363f155c232SMatthew G. Knepley } 364f155c232SMatthew G. Knepley 365f155c232SMatthew G. Knepley static PetscErrorCode QuasiInterpolate(Mat P, Vec x, Vec y) 366f155c232SMatthew G. Knepley { 367f155c232SMatthew G. Knepley QuasiInterp *ctx; 368f155c232SMatthew G. Knepley DM dmcdg, dmc; 369f155c232SMatthew G. Knepley Vec ly; 370f155c232SMatthew G. Knepley 371f155c232SMatthew G. Knepley PetscFunctionBegin; 3729566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(P, (void **) &ctx)); 3739566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Ifdg, x, ctx->tmpf)); 374f155c232SMatthew G. Knepley 3759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) ctx->tmpf, "Fine DG Potential")); 3769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpf, "fine_dg_")); 3779566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpf, NULL, "-vec_view")); 3789566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->Pi, x, ctx->tmpc)); 379f155c232SMatthew G. Knepley 3809566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) ctx->tmpc, "Coarse DG Potential")); 3819566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) ctx->tmpc, "coarse_dg_")); 3829566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ctx->tmpc, NULL, "-vec_view")); 3839566063dSJacob Faibussowitsch PetscCall(VecGetDM(ctx->tmpc, &dmcdg)); 384f155c232SMatthew G. Knepley 3859566063dSJacob Faibussowitsch PetscCall(VecGetDM(y, &dmc)); 3869566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmc, &ly)); 3879566063dSJacob Faibussowitsch PetscCall(DMPlexComputeClementInterpolant(dmcdg, ctx->tmpc, ly)); 3889566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dmc, ly, INSERT_VALUES, y)); 3899566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dmc, ly, INSERT_VALUES, y)); 3909566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmc, &ly)); 391f155c232SMatthew G. Knepley PetscFunctionReturn(0); 392f155c232SMatthew G. Knepley } 393f155c232SMatthew G. Knepley 394f155c232SMatthew G. Knepley static PetscErrorCode CreateQuasiInterpolator(DM dmc, DM dmf, Mat *P) 395f155c232SMatthew G. Knepley { 396f155c232SMatthew G. Knepley QuasiInterp *ctx; 397f155c232SMatthew G. Knepley DM dmcdg, dmfdg; 398f155c232SMatthew G. Knepley PetscFE fe; 399f155c232SMatthew G. Knepley Vec x, y; 400f155c232SMatthew G. Knepley DMPolytopeType ct; 401f155c232SMatthew G. Knepley PetscInt dim, cStart, m, n, M, N; 402f155c232SMatthew G. Knepley 403f155c232SMatthew G. Knepley PetscFunctionBegin; 4049566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(1, &ctx)); 4059566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmc, &x)); 4069566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmf, &y)); 4079566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &m)); 4089566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &M)); 4099566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &n)); 4109566063dSJacob Faibussowitsch PetscCall(VecGetSize(y, &N)); 4119566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmc, &x)); 4129566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmf, &y)); 413f155c232SMatthew G. Knepley 4149566063dSJacob Faibussowitsch PetscCall(DMClone(dmf, &dmfdg)); 4159566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmfdg, &dim)); 4169566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmfdg, 0, &cStart, NULL)); 4179566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dmfdg, cStart, &ct)); 4189566063dSJacob Faibussowitsch PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "fine_dg_", PETSC_DETERMINE, &fe)); 4199566063dSJacob Faibussowitsch PetscCall(DMSetField(dmfdg, 0, NULL, (PetscObject) fe)); 4209566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 4219566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmfdg)); 4229566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmf, dmfdg, &ctx->Ifdg, NULL)); 4239566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmfdg, &ctx->tmpf)); 4249566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmfdg)); 425f155c232SMatthew G. Knepley 4269566063dSJacob Faibussowitsch PetscCall(DMClone(dmc, &dmcdg)); 4279566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmcdg, &dim)); 4289566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmcdg, 0, &cStart, NULL)); 4299566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dmcdg, cStart, &ct)); 4309566063dSJacob Faibussowitsch PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, "coarse_dg_", PETSC_DETERMINE, &fe)); 4319566063dSJacob Faibussowitsch PetscCall(DMSetField(dmcdg, 0, NULL, (PetscObject) fe)); 4329566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 4339566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmcdg)); 434f155c232SMatthew G. Knepley 4359566063dSJacob Faibussowitsch PetscCall(CreateCoarseProjection(dmcdg, dmf, &ctx->Pi)); 4369566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmcdg, &ctx->tmpc)); 4379566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmcdg)); 438f155c232SMatthew G. Knepley 4399566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject) dmc), m, n, M, N, ctx, P)); 4409566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*P, MATOP_DESTROY, (void (*)(void)) DestroyQuasiInterpolator)); 4419566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*P, MATOP_MULT, (void (*)(void)) QuasiInterpolate)); 442f155c232SMatthew G. Knepley PetscFunctionReturn(0); 443f155c232SMatthew G. Knepley } 444f155c232SMatthew G. Knepley 445f155c232SMatthew G. Knepley static PetscErrorCode CoarseTest(DM dm, Vec u, AppCtx *user) 446f155c232SMatthew G. Knepley { 447f155c232SMatthew G. Knepley DM dmc; 448f155c232SMatthew G. Knepley Mat P; /* The quasi-interpolator to the coarse space */ 449f155c232SMatthew G. Knepley Vec uc; 450f155c232SMatthew G. Knepley 451f155c232SMatthew G. Knepley PetscFunctionBegin; 452f155c232SMatthew G. Knepley if (user->modType == MOD_CONSTANT) PetscFunctionReturn(0); 4539566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject) dm), &dmc)); 4549566063dSJacob Faibussowitsch PetscCall(DMSetType(dmc, DMPLEX)); 4559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) dmc, "coarse_")); 4569566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dmc, user)); 4579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dmc)); 4589566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmc, NULL, "-dm_view")); 459f155c232SMatthew G. Knepley 4609566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dmc, "potential", SetupPrimalProblem, user)); 4619566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmc, &uc)); 4629566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) uc, "potential")); 4639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) uc, "coarse_")); 464f155c232SMatthew G. Knepley 4659566063dSJacob Faibussowitsch PetscCall(CreateQuasiInterpolator(dmc, dm, &P)); 466f155c232SMatthew G. Knepley #if 1 4679566063dSJacob Faibussowitsch PetscCall(MatMult(P, u, uc)); 468f155c232SMatthew G. Knepley #else 469f155c232SMatthew G. Knepley { 470f155c232SMatthew G. Knepley Mat In; 471f155c232SMatthew G. Knepley Vec sc; 472f155c232SMatthew G. Knepley 4739566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmc, dm, &In, &sc)); 4749566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(In, u, uc)); 4759566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(uc, sc, uc)); 4769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&In)); 4779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sc)); 478f155c232SMatthew G. Knepley } 479f155c232SMatthew G. Knepley #endif 4809566063dSJacob Faibussowitsch PetscCall(CompareView(uc)); 481f155c232SMatthew G. Knepley 4829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 4839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uc)); 4849566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmc)); 485f155c232SMatthew G. Knepley PetscFunctionReturn(0); 486f155c232SMatthew G. Knepley } 487f155c232SMatthew G. Knepley 488f155c232SMatthew G. Knepley int main(int argc, char **argv) 489f155c232SMatthew G. Knepley { 490f155c232SMatthew G. Knepley DM dm; /* Problem specification */ 491f155c232SMatthew G. Knepley SNES snes; /* Nonlinear solver */ 492f155c232SMatthew G. Knepley Vec u; /* Solutions */ 493f155c232SMatthew G. Knepley AppCtx user; /* User-defined work context */ 494f155c232SMatthew G. Knepley 4959566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4969566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4979566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 498f155c232SMatthew G. Knepley /* Primal system */ 4999566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 5009566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 5019566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 5029566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 5039566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 5049566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 5059566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 5069566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 5079566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 5089566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 5099566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 5109566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 5119566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 5129566063dSJacob Faibussowitsch PetscCall(CompareView(u)); 513f155c232SMatthew G. Knepley /* Looking at a coarse problem */ 5149566063dSJacob Faibussowitsch PetscCall(CoarseTest(dm, u, &user)); 515f155c232SMatthew G. Knepley /* Cleanup */ 5169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5179566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 5189566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5199566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 5209566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 521b122ec5aSJacob Faibussowitsch return 0; 522f155c232SMatthew G. Knepley } 523f155c232SMatthew G. Knepley 524f155c232SMatthew G. Knepley /*TEST 525f155c232SMatthew G. Knepley 526f155c232SMatthew G. Knepley test: 527f155c232SMatthew G. Knepley suffix: 1d_p1_constant 528f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dmsnes_check 529f155c232SMatthew G. Knepley 530f155c232SMatthew G. Knepley test: 531f155c232SMatthew G. Knepley suffix: 1d_p1_constant_conv 532f155c232SMatthew G. Knepley args: -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 \ 533f155c232SMatthew G. Knepley -snes_convergence_estimate -convest_num_refine 2 534f155c232SMatthew G. Knepley 535f155c232SMatthew G. Knepley test: 536f155c232SMatthew G. Knepley suffix: 1d_p1_oscillatory 537f155c232SMatthew G. Knepley args: -mod_type oscillatory -epsilon 0.03125 \ 538f155c232SMatthew G. Knepley -dm_plex_dim 1 -dm_plex_box_faces 4 -potential_petscspace_degree 1 -dm_refine 2 -dmsnes_check \ 539f155c232SMatthew G. Knepley -coarse_dm_plex_dim 1 -coarse_dm_plex_box_faces 4 -coarse_dm_plex_hash_location \ 540f155c232SMatthew G. Knepley -fine_dg_petscspace_degree 1 -fine_dg_petscdualspace_lagrange_continuity 0 \ 541f155c232SMatthew G. Knepley -coarse_dg_petscspace_degree 1 -coarse_dg_petscdualspace_lagrange_continuity 0 542f155c232SMatthew G. Knepley 543f155c232SMatthew G. Knepley TEST*/ 544