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