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