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