1 static char help[] = "Poisson Problem in 2d and 3d with finite elements.\n\ 2 We solve the Poisson problem in a rectangular\n\ 3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4 This example supports automatic convergence estimation\n\ 5 and eventually adaptivity.\n\n\n"; 6 7 #include <petscdmplex.h> 8 #include <petscsnes.h> 9 #include <petscds.h> 10 #include <petscconvest.h> 11 12 typedef struct { 13 /* Domain and mesh definition */ 14 PetscBool spectral; /* Look at the spectrum along planes in the solution */ 15 PetscBool shear; /* Shear the domain */ 16 PetscBool adjoint; /* Solve the adjoint problem */ 17 PetscBool homogeneous; /* Use homogeneous boudnary conditions */ 18 PetscBool viewError; /* Output the solution error */ 19 } AppCtx; 20 21 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 22 { 23 *u = 0.0; 24 return 0; 25 } 26 27 static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 28 { 29 PetscInt d; 30 *u = 0.0; 31 for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 32 return 0; 33 } 34 35 static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 36 { 37 PetscInt d; 38 *u = 1.0; 39 for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0*PETSC_PI*x[d]); 40 return 0; 41 } 42 43 /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */ 44 static void obj_error_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 45 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 46 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 47 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[]) 48 { 49 obj[0] = a[aOff[0]]*(u[0] - a[aOff[1]]); 50 } 51 52 static void f0_trig_inhomogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 53 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 54 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 55 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 56 { 57 PetscInt d; 58 for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 59 } 60 61 static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 62 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 63 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 64 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 65 { 66 PetscInt d; 67 for (d = 0; d < dim; ++d) { 68 PetscScalar v = 1.; 69 for (PetscInt e = 0; e < dim; e++) { 70 if (e == d) { 71 v *= -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 72 } else { 73 v *= PetscSinReal(2.0*PETSC_PI*x[d]); 74 } 75 } 76 f0[0] += v; 77 } 78 } 79 80 static void f0_unity_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 81 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 82 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 83 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 84 { 85 f0[0] = 1.0; 86 } 87 88 static void f0_identityaux_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 89 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 90 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 91 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 92 { 93 f0[0] = a[0]; 94 } 95 96 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 97 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 98 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 99 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 100 { 101 PetscInt d; 102 for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 103 } 104 105 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 106 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 107 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 108 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 109 { 110 PetscInt d; 111 for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 112 } 113 114 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 115 { 116 PetscFunctionBeginUser; 117 options->shear = PETSC_FALSE; 118 options->spectral = PETSC_FALSE; 119 options->adjoint = PETSC_FALSE; 120 options->homogeneous = PETSC_FALSE; 121 options->viewError = PETSC_FALSE; 122 123 PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 124 PetscCall(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL)); 125 PetscCall(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL)); 126 PetscCall(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL)); 127 PetscCall(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL)); 128 PetscCall(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL)); 129 PetscOptionsEnd(); 130 PetscFunctionReturn(0); 131 } 132 133 static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 134 { 135 PetscSection coordSection; 136 Vec coordinates; 137 const PetscScalar *coords; 138 PetscInt dim, p, vStart, vEnd, v; 139 140 PetscFunctionBeginUser; 141 PetscCall(DMGetCoordinateDim(dm, &dim)); 142 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 143 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 144 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 145 PetscCall(VecGetArrayRead(coordinates, &coords)); 146 for (p = 0; p < numPlanes; ++p) { 147 DMLabel label; 148 char name[PETSC_MAX_PATH_LEN]; 149 150 PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 151 PetscCall(DMCreateLabel(dm, name)); 152 PetscCall(DMGetLabel(dm, name, &label)); 153 PetscCall(DMLabelAddStratum(label, 1)); 154 for (v = vStart; v < vEnd; ++v) { 155 PetscInt off; 156 157 PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 158 if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off+planeDir[p]])) < PETSC_SMALL) { 159 PetscCall(DMLabelSetValue(label, v, 1)); 160 } 161 } 162 } 163 PetscCall(VecRestoreArrayRead(coordinates, &coords)); 164 PetscFunctionReturn(0); 165 } 166 167 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 168 { 169 PetscFunctionBeginUser; 170 PetscCall(DMCreate(comm, dm)); 171 PetscCall(DMSetType(*dm, DMPLEX)); 172 PetscCall(DMSetFromOptions(*dm)); 173 if (user->shear) PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL)); 174 PetscCall(DMSetApplicationContext(*dm, user)); 175 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 176 if (user->spectral) { 177 PetscInt planeDir[2] = {0, 1}; 178 PetscReal planeCoord[2] = {0., 1.}; 179 180 PetscCall(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user)); 181 } 182 PetscFunctionReturn(0); 183 } 184 185 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 186 { 187 PetscDS ds; 188 DMLabel label; 189 const PetscInt id = 1; 190 PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u; 191 PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u; 192 193 PetscFunctionBeginUser; 194 PetscCall(DMGetDS(dm, &ds)); 195 PetscCall(PetscDSSetResidual(ds, 0, f0, f1_u)); 196 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 197 PetscCall(PetscDSSetExactSolution(ds, 0, ex, user)); 198 PetscCall(DMGetLabel(dm, "marker", &label)); 199 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, user, NULL)); 200 PetscFunctionReturn(0); 201 } 202 203 static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user) 204 { 205 PetscDS ds; 206 DMLabel label; 207 const PetscInt id = 1; 208 209 PetscFunctionBeginUser; 210 PetscCall(DMGetDS(dm, &ds)); 211 PetscCall(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u)); 212 PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 213 PetscCall(PetscDSSetObjective(ds, 0, obj_error_u)); 214 PetscCall(DMGetLabel(dm, "marker", &label)); 215 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, user, NULL)); 216 PetscFunctionReturn(0); 217 } 218 219 static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user) 220 { 221 PetscDS prob; 222 223 PetscFunctionBeginUser; 224 PetscCall(DMGetDS(dm, &prob)); 225 PetscFunctionReturn(0); 226 } 227 228 static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 229 { 230 DM cdm = dm; 231 PetscFE fe; 232 DMPolytopeType ct; 233 PetscBool simplex; 234 PetscInt dim, cStart; 235 char prefix[PETSC_MAX_PATH_LEN]; 236 237 PetscFunctionBeginUser; 238 PetscCall(DMGetDimension(dm, &dim)); 239 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 240 PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 241 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 242 /* Create finite element */ 243 PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 244 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 245 PetscCall(PetscObjectSetName((PetscObject) fe, name)); 246 /* Set discretization and boundary conditions for each mesh */ 247 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 248 PetscCall(DMCreateDS(dm)); 249 PetscCall((*setup)(dm, user)); 250 while (cdm) { 251 PetscCall(DMCopyDisc(dm,cdm)); 252 /* TODO: Check whether the boundary of coarse meshes is marked */ 253 PetscCall(DMGetCoarseDM(cdm, &cdm)); 254 } 255 PetscCall(PetscFEDestroy(&fe)); 256 PetscFunctionReturn(0); 257 } 258 259 static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 260 { 261 MPI_Comm comm; 262 PetscSection coordSection, section; 263 Vec coordinates, uloc; 264 const PetscScalar *coords, *array; 265 PetscInt p; 266 PetscMPIInt size, rank; 267 268 PetscFunctionBeginUser; 269 PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 270 PetscCallMPI(MPI_Comm_size(comm, &size)); 271 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 272 PetscCall(DMGetLocalVector(dm, &uloc)); 273 PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc)); 274 PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc)); 275 PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL)); 276 PetscCall(VecViewFromOptions(uloc, NULL, "-sol_view")); 277 PetscCall(DMGetLocalSection(dm, §ion)); 278 PetscCall(VecGetArrayRead(uloc, &array)); 279 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 280 PetscCall(DMGetCoordinateSection(dm, &coordSection)); 281 PetscCall(VecGetArrayRead(coordinates, &coords)); 282 for (p = 0; p < numPlanes; ++p) { 283 DMLabel label; 284 char name[PETSC_MAX_PATH_LEN]; 285 Mat F; 286 Vec x, y; 287 IS stratum; 288 PetscReal *ray, *gray; 289 PetscScalar *rvals, *svals, *gsvals; 290 PetscInt *perm, *nperm; 291 PetscInt n, N, i, j, off, offu; 292 const PetscInt *points; 293 294 PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 295 PetscCall(DMGetLabel(dm, name, &label)); 296 PetscCall(DMLabelGetStratumIS(label, 1, &stratum)); 297 PetscCall(ISGetLocalSize(stratum, &n)); 298 PetscCall(ISGetIndices(stratum, &points)); 299 PetscCall(PetscMalloc2(n, &ray, n, &svals)); 300 for (i = 0; i < n; ++i) { 301 PetscCall(PetscSectionGetOffset(coordSection, points[i], &off)); 302 PetscCall(PetscSectionGetOffset(section, points[i], &offu)); 303 ray[i] = PetscRealPart(coords[off+((planeDir[p]+1)%2)]); 304 svals[i] = array[offu]; 305 } 306 /* Gather the ray data to proc 0 */ 307 if (size > 1) { 308 PetscMPIInt *cnt, *displs, p; 309 310 PetscCall(PetscCalloc2(size, &cnt, size, &displs)); 311 PetscCallMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm)); 312 for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1]; 313 N = displs[size-1] + cnt[size-1]; 314 PetscCall(PetscMalloc2(N, &gray, N, &gsvals)); 315 PetscCallMPI(MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm)); 316 PetscCallMPI(MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm)); 317 PetscCall(PetscFree2(cnt, displs)); 318 } else { 319 N = n; 320 gray = ray; 321 gsvals = svals; 322 } 323 if (rank == 0) { 324 /* Sort point along ray */ 325 PetscCall(PetscMalloc2(N, &perm, N, &nperm)); 326 for (i = 0; i < N; ++i) {perm[i] = i;} 327 PetscCall(PetscSortRealWithPermutation(N, gray, perm)); 328 /* Count duplicates and squish mapping */ 329 nperm[0] = perm[0]; 330 for (i = 1, j = 1; i < N; ++i) { 331 if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i]; 332 } 333 /* Create FFT structs */ 334 PetscCall(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F)); 335 PetscCall(MatCreateVecs(F, &x, &y)); 336 PetscCall(PetscObjectSetName((PetscObject) y, name)); 337 PetscCall(VecGetArray(x, &rvals)); 338 for (i = 0, j = 0; i < N; ++i) { 339 if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue; 340 rvals[j] = gsvals[nperm[j]]; 341 ++j; 342 } 343 PetscCall(PetscFree2(perm, nperm)); 344 if (size > 1) PetscCall(PetscFree2(gray, gsvals)); 345 PetscCall(VecRestoreArray(x, &rvals)); 346 /* Do FFT along the ray */ 347 PetscCall(MatMult(F, x, y)); 348 /* Chop FFT */ 349 PetscCall(VecChop(y, PETSC_SMALL)); 350 PetscCall(VecViewFromOptions(x, NULL, "-real_view")); 351 PetscCall(VecViewFromOptions(y, NULL, "-fft_view")); 352 PetscCall(VecDestroy(&x)); 353 PetscCall(VecDestroy(&y)); 354 PetscCall(MatDestroy(&F)); 355 } 356 PetscCall(ISRestoreIndices(stratum, &points)); 357 PetscCall(ISDestroy(&stratum)); 358 PetscCall(PetscFree2(ray, svals)); 359 } 360 PetscCall(VecRestoreArrayRead(coordinates, &coords)); 361 PetscCall(VecRestoreArrayRead(uloc, &array)); 362 PetscCall(DMRestoreLocalVector(dm, &uloc)); 363 PetscFunctionReturn(0); 364 } 365 366 int main(int argc, char **argv) 367 { 368 DM dm; /* Problem specification */ 369 SNES snes; /* Nonlinear solver */ 370 Vec u; /* Solutions */ 371 AppCtx user; /* User-defined work context */ 372 373 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 374 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 375 /* Primal system */ 376 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 377 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 378 PetscCall(SNESSetDM(snes, dm)); 379 PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 380 PetscCall(DMCreateGlobalVector(dm, &u)); 381 PetscCall(VecSet(u, 0.0)); 382 PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 383 PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 384 PetscCall(SNESSetFromOptions(snes)); 385 PetscCall(SNESSolve(snes, NULL, u)); 386 PetscCall(SNESGetSolution(snes, &u)); 387 PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 388 if (user.viewError) { 389 PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 390 void *ctx; 391 PetscDS ds; 392 PetscReal error; 393 PetscInt N; 394 395 PetscCall(DMGetDS(dm, &ds)); 396 PetscCall(PetscDSGetExactSolution(ds, 0, &sol, &ctx)); 397 PetscCall(VecGetSize(u, &N)); 398 PetscCall(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error)); 399 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g\n", N, (double)error)); 400 } 401 if (user.spectral) { 402 PetscInt planeDir[2] = {0, 1}; 403 PetscReal planeCoord[2] = {0., 1.}; 404 405 PetscCall(ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user)); 406 } 407 /* Adjoint system */ 408 if (user.adjoint) { 409 DM dmAdj; 410 SNES snesAdj; 411 Vec uAdj; 412 413 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snesAdj)); 414 PetscCall(PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_")); 415 PetscCall(DMClone(dm, &dmAdj)); 416 PetscCall(SNESSetDM(snesAdj, dmAdj)); 417 PetscCall(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user)); 418 PetscCall(DMCreateGlobalVector(dmAdj, &uAdj)); 419 PetscCall(VecSet(uAdj, 0.0)); 420 PetscCall(PetscObjectSetName((PetscObject) uAdj, "adjoint")); 421 PetscCall(DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user)); 422 PetscCall(SNESSetFromOptions(snesAdj)); 423 PetscCall(SNESSolve(snesAdj, NULL, uAdj)); 424 PetscCall(SNESGetSolution(snesAdj, &uAdj)); 425 PetscCall(VecViewFromOptions(uAdj, NULL, "-adjoint_view")); 426 /* Error representation */ 427 { 428 DM dmErr, dmErrAux, dms[2]; 429 Vec errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj; 430 IS *subis; 431 PetscReal errorEstTot, errorL2Norm, errorL2Tot; 432 PetscInt N, i; 433 PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u}; 434 void (*identity[1])(PetscInt, PetscInt, PetscInt, 435 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 436 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 437 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u}; 438 void *ctxs[1] = {0}; 439 440 ctxs[0] = &user; 441 PetscCall(DMClone(dm, &dmErr)); 442 PetscCall(SetupDiscretization(dmErr, "error", SetupErrorProblem, &user)); 443 PetscCall(DMGetGlobalVector(dmErr, &errorEst)); 444 PetscCall(DMGetGlobalVector(dmErr, &errorL2)); 445 /* Compute auxiliary data (solution and projection of adjoint solution) */ 446 PetscCall(DMGetLocalVector(dmAdj, &uAdjLoc)); 447 PetscCall(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 448 PetscCall(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 449 PetscCall(DMGetGlobalVector(dm, &uAdjProj)); 450 PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc)); 451 PetscCall(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj)); 452 PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 453 PetscCall(DMRestoreLocalVector(dmAdj, &uAdjLoc)); 454 /* Attach auxiliary data */ 455 dms[0] = dm; dms[1] = dm; 456 PetscCall(DMCreateSuperDM(dms, 2, &subis, &dmErrAux)); 457 if (0) { 458 PetscSection sec; 459 460 PetscCall(DMGetLocalSection(dms[0], &sec)); 461 PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 462 PetscCall(DMGetLocalSection(dms[1], &sec)); 463 PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 464 PetscCall(DMGetLocalSection(dmErrAux, &sec)); 465 PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 466 } 467 PetscCall(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view")); 468 PetscCall(ISViewFromOptions(subis[0], NULL, "-super_is_view")); 469 PetscCall(ISViewFromOptions(subis[1], NULL, "-super_is_view")); 470 PetscCall(DMGetGlobalVector(dmErrAux, &uErr)); 471 PetscCall(VecViewFromOptions(u, NULL, "-map_vec_view")); 472 PetscCall(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view")); 473 PetscCall(VecViewFromOptions(uErr, NULL, "-map_vec_view")); 474 PetscCall(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u)); 475 PetscCall(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj)); 476 PetscCall(DMRestoreGlobalVector(dm, &uAdjProj)); 477 for (i = 0; i < 2; ++i) PetscCall(ISDestroy(&subis[i])); 478 PetscCall(PetscFree(subis)); 479 PetscCall(DMGetLocalVector(dmErrAux, &uErrLoc)); 480 PetscCall(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc)); 481 PetscCall(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc)); 482 PetscCall(DMRestoreGlobalVector(dmErrAux, &uErr)); 483 PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc)); 484 /* Compute cellwise error estimate */ 485 PetscCall(VecSet(errorEst, 0.0)); 486 PetscCall(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user)); 487 PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL)); 488 PetscCall(DMRestoreLocalVector(dmErrAux, &uErrLoc)); 489 PetscCall(DMDestroy(&dmErrAux)); 490 /* Plot cellwise error vector */ 491 PetscCall(VecViewFromOptions(errorEst, NULL, "-error_view")); 492 /* Compute ratio of estimate (sum over cells) with actual L_2 error */ 493 PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm)); 494 PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2)); 495 PetscCall(VecViewFromOptions(errorL2, NULL, "-l2_error_view")); 496 PetscCall(VecNorm(errorL2, NORM_INFINITY, &errorL2Tot)); 497 PetscCall(VecNorm(errorEst, NORM_INFINITY, &errorEstTot)); 498 PetscCall(VecGetSize(errorEst, &N)); 499 PetscCall(VecPointwiseDivide(errorEst, errorEst, errorL2)); 500 PetscCall(PetscObjectSetName((PetscObject) errorEst, "Error ratio")); 501 PetscCall(VecViewFromOptions(errorEst, NULL, "-error_ratio_view")); 502 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g Error Ratio: %g/%g = %g\n", N, (double) errorL2Norm, (double) errorEstTot, (double) PetscSqrtReal(errorL2Tot), (double)(errorEstTot/PetscSqrtReal(errorL2Tot)))); 503 PetscCall(DMRestoreGlobalVector(dmErr, &errorEst)); 504 PetscCall(DMRestoreGlobalVector(dmErr, &errorL2)); 505 PetscCall(DMDestroy(&dmErr)); 506 } 507 PetscCall(DMDestroy(&dmAdj)); 508 PetscCall(VecDestroy(&uAdj)); 509 PetscCall(SNESDestroy(&snesAdj)); 510 } 511 /* Cleanup */ 512 PetscCall(VecDestroy(&u)); 513 PetscCall(SNESDestroy(&snes)); 514 PetscCall(DMDestroy(&dm)); 515 PetscCall(PetscFinalize()); 516 return 0; 517 } 518 519 /*TEST 520 521 test: 522 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 523 suffix: 2d_p1_conv 524 requires: triangle 525 args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 526 test: 527 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 528 suffix: 2d_p2_conv 529 requires: triangle 530 args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 531 test: 532 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 533 suffix: 2d_p3_conv 534 requires: triangle 535 args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 536 test: 537 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 538 suffix: 2d_q1_conv 539 args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 540 test: 541 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 542 suffix: 2d_q2_conv 543 args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 544 test: 545 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 546 suffix: 2d_q3_conv 547 args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 548 test: 549 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 550 suffix: 2d_q1_shear_conv 551 args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 552 test: 553 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 554 suffix: 2d_q2_shear_conv 555 args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 556 test: 557 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 558 suffix: 2d_q3_shear_conv 559 args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 560 test: 561 # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7 562 suffix: 3d_p1_conv 563 requires: ctetgen 564 args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 565 test: 566 # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8 567 suffix: 3d_p2_conv 568 requires: ctetgen 569 args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 570 test: 571 # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0 572 suffix: 3d_p3_conv 573 requires: ctetgen 574 args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 575 test: 576 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8 577 suffix: 3d_q1_conv 578 args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 579 test: 580 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8 581 suffix: 3d_q2_conv 582 args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 583 test: 584 # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8 585 suffix: 3d_q3_conv 586 args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 587 test: 588 suffix: 2d_p1_fas_full 589 requires: triangle 590 args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 591 -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \ 592 -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 593 -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 594 -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 595 -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10 596 test: 597 suffix: 2d_p1_fas_full_homogeneous 598 requires: triangle 599 args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 600 -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \ 601 -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 602 -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 603 -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 604 -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10 605 606 test: 607 suffix: 2d_p1_scalable 608 requires: triangle 609 args: -potential_petscspace_degree 1 -dm_refine 3 \ 610 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \ 611 -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \ 612 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 613 -pc_gamg_coarse_eq_limit 1000 \ 614 -pc_gamg_square_graph 1 \ 615 -pc_gamg_threshold 0.05 \ 616 -pc_gamg_threshold_scale .0 \ 617 -mg_levels_ksp_type chebyshev \ 618 -mg_levels_ksp_max_it 1 \ 619 -mg_levels_pc_type jacobi \ 620 -matptap_via scalable 621 test: 622 suffix: 2d_p1_gmg_vcycle 623 requires: triangle 624 args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 625 -ksp_rtol 5e-10 -pc_type mg \ 626 -mg_levels_ksp_max_it 1 \ 627 -mg_levels_esteig_ksp_type cg \ 628 -mg_levels_esteig_ksp_max_it 10 \ 629 -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 630 -mg_levels_pc_type jacobi 631 test: 632 suffix: 2d_p1_gmg_fcycle 633 requires: triangle 634 args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 635 -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \ 636 -mg_levels_ksp_max_it 2 \ 637 -mg_levels_esteig_ksp_type cg \ 638 -mg_levels_esteig_ksp_max_it 10 \ 639 -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 640 -mg_levels_pc_type jacobi 641 test: 642 suffix: 2d_p1_gmg_vcycle_adapt 643 requires: triangle bamg 644 args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 645 -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \ 646 -mg_levels_ksp_max_it 1 \ 647 -mg_levels_esteig_ksp_type cg \ 648 -mg_levels_esteig_ksp_max_it 10 \ 649 -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 650 -mg_levels_pc_type jacobi 651 test: 652 suffix: 2d_p1_spectral_0 653 requires: triangle fftw !complex 654 args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view 655 test: 656 suffix: 2d_p1_spectral_1 657 requires: triangle fftw !complex 658 nsize: 2 659 args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view 660 test: 661 suffix: 2d_p1_adj_0 662 requires: triangle 663 args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0 664 test: 665 nsize: 2 666 requires: !sycl kokkos_kernels 667 suffix: kokkos 668 args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,6 -petscpartitioner_type simple -dm_plex_simplex 0 -potential_petscspace_degree 1 \ 669 -dm_refine 0 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_coarse_eq_limit 1000 -pc_gamg_threshold 0.0 \ 670 -pc_gamg_threshold_scale .5 -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 2 -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \ 671 -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 672 673 TEST*/ 674