1c4762a1bSJed Brown static char help[] = "Poisson Problem in 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports automatic convergence estimation\n\ 5c4762a1bSJed Brown and eventually adaptivity.\n\n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petscdmplex.h> 8c4762a1bSJed Brown #include <petscsnes.h> 9c4762a1bSJed Brown #include <petscds.h> 10c4762a1bSJed Brown #include <petscconvest.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown typedef struct { 13c4762a1bSJed Brown /* Domain and mesh definition */ 14c4762a1bSJed Brown PetscBool spectral; /* Look at the spectrum along planes in the solution */ 15c4762a1bSJed Brown PetscBool shear; /* Shear the domain */ 16c4762a1bSJed Brown PetscBool adjoint; /* Solve the adjoint problem */ 17ad1a7c93SMatthew Knepley PetscBool homogeneous; /* Use homogeneous boudnary conditions */ 18ad1a7c93SMatthew Knepley PetscBool viewError; /* Output the solution error */ 19c4762a1bSJed Brown } AppCtx; 20c4762a1bSJed Brown 219371c9d4SSatish Balay static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 22c4762a1bSJed Brown *u = 0.0; 23c4762a1bSJed Brown return 0; 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 269371c9d4SSatish Balay static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 27c4762a1bSJed Brown PetscInt d; 28c4762a1bSJed Brown *u = 0.0; 29c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0 * PETSC_PI * x[d]); 30c4762a1bSJed Brown return 0; 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 339371c9d4SSatish Balay static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 34b724ec41SToby Isaac PetscInt d; 35b724ec41SToby Isaac *u = 1.0; 36b724ec41SToby Isaac for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0 * PETSC_PI * x[d]); 37b724ec41SToby Isaac return 0; 38b724ec41SToby Isaac } 39b724ec41SToby Isaac 40c4762a1bSJed Brown /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */ 419371c9d4SSatish Balay static void obj_error_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 obj[]) { 42c4762a1bSJed Brown obj[0] = a[aOff[0]] * (u[0] - a[aOff[1]]); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 459371c9d4SSatish Balay static void f0_trig_inhomogeneous_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[]) { 46c4762a1bSJed Brown PetscInt d; 47c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 509371c9d4SSatish Balay 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[]) { 51b724ec41SToby Isaac PetscInt d; 52b724ec41SToby Isaac for (d = 0; d < dim; ++d) { 53b724ec41SToby Isaac PetscScalar v = 1.; 54b724ec41SToby Isaac for (PetscInt e = 0; e < dim; e++) { 55b724ec41SToby Isaac if (e == d) { 56b724ec41SToby Isaac v *= -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 57b724ec41SToby Isaac } else { 58b724ec41SToby Isaac v *= PetscSinReal(2.0 * PETSC_PI * x[d]); 59b724ec41SToby Isaac } 60b724ec41SToby Isaac } 61b724ec41SToby Isaac f0[0] += v; 62b724ec41SToby Isaac } 63b724ec41SToby Isaac } 64b724ec41SToby Isaac 659371c9d4SSatish Balay static void f0_unity_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[]) { 66c4762a1bSJed Brown f0[0] = 1.0; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 699371c9d4SSatish Balay static void f0_identityaux_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[]) { 70c4762a1bSJed Brown f0[0] = a[0]; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 739371c9d4SSatish Balay 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[]) { 74c4762a1bSJed Brown PetscInt d; 75c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 789371c9d4SSatish Balay 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[]) { 79c4762a1bSJed Brown PetscInt d; 80c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 839371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 84c4762a1bSJed Brown PetscFunctionBeginUser; 85c4762a1bSJed Brown options->shear = PETSC_FALSE; 86c4762a1bSJed Brown options->spectral = PETSC_FALSE; 87c4762a1bSJed Brown options->adjoint = PETSC_FALSE; 88b724ec41SToby Isaac options->homogeneous = PETSC_FALSE; 89ad1a7c93SMatthew Knepley options->viewError = PETSC_FALSE; 90c4762a1bSJed Brown 91d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL)); 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL)); 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL)); 97d0609cedSBarry Smith PetscOptionsEnd(); 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 1019371c9d4SSatish Balay static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) { 102c4762a1bSJed Brown PetscSection coordSection; 103c4762a1bSJed Brown Vec coordinates; 104c4762a1bSJed Brown const PetscScalar *coords; 105c4762a1bSJed Brown PetscInt dim, p, vStart, vEnd, v; 106c4762a1bSJed Brown 107c4762a1bSJed Brown PetscFunctionBeginUser; 1089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dim)); 1099566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1109566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 113c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 114c4762a1bSJed Brown DMLabel label; 115c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 116c4762a1bSJed Brown 11763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 1189566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, name)); 1199566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 1209566063dSJacob Faibussowitsch PetscCall(DMLabelAddStratum(label, 1)); 121c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 122c4762a1bSJed Brown PetscInt off; 123c4762a1bSJed Brown 1249566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 1259371c9d4SSatish Balay if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off + planeDir[p]])) < PETSC_SMALL) { PetscCall(DMLabelSetValue(label, v, 1)); } 126c4762a1bSJed Brown } 127c4762a1bSJed Brown } 1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 129c4762a1bSJed Brown PetscFunctionReturn(0); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 1329371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 133c4762a1bSJed Brown PetscFunctionBeginUser; 1349566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1359566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1369566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1379566063dSJacob Faibussowitsch if (user->shear) PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL)); 1389566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 1399566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 140c4762a1bSJed Brown if (user->spectral) { 141c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 142c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 143c4762a1bSJed Brown 1449566063dSJacob Faibussowitsch PetscCall(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown PetscFunctionReturn(0); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 1499371c9d4SSatish Balay static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) { 15045480ffeSMatthew G. Knepley PetscDS ds; 15145480ffeSMatthew G. Knepley DMLabel label; 152c4762a1bSJed Brown const PetscInt id = 1; 153b724ec41SToby Isaac PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u; 154b724ec41SToby Isaac PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBeginUser; 1579566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1589566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0, f1_u)); 1599566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1609566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, ex, user)); 1619566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 1629566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))ex, NULL, user, NULL)); 163c4762a1bSJed Brown PetscFunctionReturn(0); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 1669371c9d4SSatish Balay static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user) { 16745480ffeSMatthew G. Knepley PetscDS ds; 16845480ffeSMatthew G. Knepley DMLabel label; 169c4762a1bSJed Brown const PetscInt id = 1; 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscFunctionBeginUser; 1729566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1739566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u)); 1749566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1759566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, obj_error_u)); 1769566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 1779566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, user, NULL)); 178c4762a1bSJed Brown PetscFunctionReturn(0); 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 1819371c9d4SSatish Balay static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user) { 182c4762a1bSJed Brown PetscDS prob; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBeginUser; 1859566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 186c4762a1bSJed Brown PetscFunctionReturn(0); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 1899371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) { 190c4762a1bSJed Brown DM cdm = dm; 191c4762a1bSJed Brown PetscFE fe; 1925be41b18SMatthew G. Knepley DMPolytopeType ct; 1935be41b18SMatthew G. Knepley PetscBool simplex; 1945be41b18SMatthew G. Knepley PetscInt dim, cStart; 195c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 196c4762a1bSJed Brown 197c4762a1bSJed Brown PetscFunctionBeginUser; 1989566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1999566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2009566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 2015be41b18SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 202c4762a1bSJed Brown /* Create finite element */ 2039566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 2049566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 2059566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name)); 206c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2079566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2089566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2099566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 210c4762a1bSJed Brown while (cdm) { 2119566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 212c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 2139566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 214c4762a1bSJed Brown } 2159566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 2199371c9d4SSatish Balay static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) { 220c4762a1bSJed Brown MPI_Comm comm; 221c4762a1bSJed Brown PetscSection coordSection, section; 222c4762a1bSJed Brown Vec coordinates, uloc; 223c4762a1bSJed Brown const PetscScalar *coords, *array; 224c4762a1bSJed Brown PetscInt p; 225c4762a1bSJed Brown PetscMPIInt size, rank; 226c4762a1bSJed Brown 227c4762a1bSJed Brown PetscFunctionBeginUser; 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2319566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &uloc)); 2329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc)); 2339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc)); 2349566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL)); 2359566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uloc, NULL, "-sol_view")); 2369566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 2379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(uloc, &array)); 2389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 241c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 242c4762a1bSJed Brown DMLabel label; 243c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 244c4762a1bSJed Brown Mat F; 245c4762a1bSJed Brown Vec x, y; 246c4762a1bSJed Brown IS stratum; 247c4762a1bSJed Brown PetscReal *ray, *gray; 248c4762a1bSJed Brown PetscScalar *rvals, *svals, *gsvals; 249c4762a1bSJed Brown PetscInt *perm, *nperm; 250c4762a1bSJed Brown PetscInt n, N, i, j, off, offu; 251c4762a1bSJed Brown const PetscInt *points; 252c4762a1bSJed Brown 25363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 2549566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 2559566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, 1, &stratum)); 2569566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(stratum, &n)); 2579566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratum, &points)); 2589566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ray, n, &svals)); 259c4762a1bSJed Brown for (i = 0; i < n; ++i) { 2609566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, points[i], &off)); 2619566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, points[i], &offu)); 262c4762a1bSJed Brown ray[i] = PetscRealPart(coords[off + ((planeDir[p] + 1) % 2)]); 263c4762a1bSJed Brown svals[i] = array[offu]; 264c4762a1bSJed Brown } 265c4762a1bSJed Brown /* Gather the ray data to proc 0 */ 266c4762a1bSJed Brown if (size > 1) { 267c4762a1bSJed Brown PetscMPIInt *cnt, *displs, p; 268c4762a1bSJed Brown 2699566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(size, &cnt, size, &displs)); 2709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm)); 271c4762a1bSJed Brown for (p = 1; p < size; ++p) displs[p] = displs[p - 1] + cnt[p - 1]; 272c4762a1bSJed Brown N = displs[size - 1] + cnt[size - 1]; 2739566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &gray, N, &gsvals)); 2749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm)); 2759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree2(cnt, displs)); 277c4762a1bSJed Brown } else { 278c4762a1bSJed Brown N = n; 279c4762a1bSJed Brown gray = ray; 280c4762a1bSJed Brown gsvals = svals; 281c4762a1bSJed Brown } 282dd400576SPatrick Sanan if (rank == 0) { 283c4762a1bSJed Brown /* Sort point along ray */ 2849566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &perm, N, &nperm)); 285c4762a1bSJed Brown for (i = 0; i < N; ++i) { perm[i] = i; } 2869566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithPermutation(N, gray, perm)); 287c4762a1bSJed Brown /* Count duplicates and squish mapping */ 288c4762a1bSJed Brown nperm[0] = perm[0]; 289c4762a1bSJed Brown for (i = 1, j = 1; i < N; ++i) { 290c4762a1bSJed Brown if (PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) > PETSC_SMALL) nperm[j++] = perm[i]; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown /* Create FFT structs */ 2939566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F)); 2949566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(F, &x, &y)); 2959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, name)); 2969566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &rvals)); 297c4762a1bSJed Brown for (i = 0, j = 0; i < N; ++i) { 298c4762a1bSJed Brown if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) < PETSC_SMALL) continue; 299c4762a1bSJed Brown rvals[j] = gsvals[nperm[j]]; 300c4762a1bSJed Brown ++j; 301c4762a1bSJed Brown } 3029566063dSJacob Faibussowitsch PetscCall(PetscFree2(perm, nperm)); 3039566063dSJacob Faibussowitsch if (size > 1) PetscCall(PetscFree2(gray, gsvals)); 3049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &rvals)); 305c4762a1bSJed Brown /* Do FFT along the ray */ 3069566063dSJacob Faibussowitsch PetscCall(MatMult(F, x, y)); 307c4762a1bSJed Brown /* Chop FFT */ 3089566063dSJacob Faibussowitsch PetscCall(VecChop(y, PETSC_SMALL)); 3099566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(x, NULL, "-real_view")); 3109566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(y, NULL, "-fft_view")); 3119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 3139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 314c4762a1bSJed Brown } 3159566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratum, &points)); 3169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratum)); 3179566063dSJacob Faibussowitsch PetscCall(PetscFree2(ray, svals)); 318c4762a1bSJed Brown } 3199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 3209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(uloc, &array)); 3219566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &uloc)); 322c4762a1bSJed Brown PetscFunctionReturn(0); 323c4762a1bSJed Brown } 324c4762a1bSJed Brown 3259371c9d4SSatish Balay int main(int argc, char **argv) { 326c4762a1bSJed Brown DM dm; /* Problem specification */ 327c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 328c4762a1bSJed Brown Vec u; /* Solutions */ 329c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 330c4762a1bSJed Brown 331327415f7SBarry Smith PetscFunctionBeginUser; 3329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3339566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 334c4762a1bSJed Brown /* Primal system */ 3359566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 3369566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 3379566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 3389566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 3399566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3409566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "potential")); 3429566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 3439566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 3449566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3459566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 3469566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 347ad1a7c93SMatthew Knepley if (user.viewError) { 348ad1a7c93SMatthew Knepley PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 349ad1a7c93SMatthew Knepley void *ctx; 350ad1a7c93SMatthew Knepley PetscDS ds; 351ad1a7c93SMatthew Knepley PetscReal error; 352ad1a7c93SMatthew Knepley PetscInt N; 353ad1a7c93SMatthew Knepley 3549566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3559566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &sol, &ctx)); 3569566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &N)); 3579566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error)); 35863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g\n", N, (double)error)); 359ad1a7c93SMatthew Knepley } 360c4762a1bSJed Brown if (user.spectral) { 361c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 362c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 363c4762a1bSJed Brown 3649566063dSJacob Faibussowitsch PetscCall(ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user)); 365c4762a1bSJed Brown } 366c4762a1bSJed Brown /* Adjoint system */ 367c4762a1bSJed Brown if (user.adjoint) { 368c4762a1bSJed Brown DM dmAdj; 369c4762a1bSJed Brown SNES snesAdj; 370c4762a1bSJed Brown Vec uAdj; 371c4762a1bSJed Brown 3729566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snesAdj)); 3739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)snesAdj, "adjoint_")); 3749566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAdj)); 3759566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snesAdj, dmAdj)); 3769566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user)); 3779566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmAdj, &uAdj)); 3789566063dSJacob Faibussowitsch PetscCall(VecSet(uAdj, 0.0)); 3799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)uAdj, "adjoint")); 3809566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user)); 3819566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snesAdj)); 3829566063dSJacob Faibussowitsch PetscCall(SNESSolve(snesAdj, NULL, uAdj)); 3839566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snesAdj, &uAdj)); 3849566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uAdj, NULL, "-adjoint_view")); 385c4762a1bSJed Brown /* Error representation */ 386c4762a1bSJed Brown { 387c4762a1bSJed Brown DM dmErr, dmErrAux, dms[2]; 388c4762a1bSJed Brown Vec errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj; 389c4762a1bSJed Brown IS *subis; 390c4762a1bSJed Brown PetscReal errorEstTot, errorL2Norm, errorL2Tot; 391c4762a1bSJed Brown PetscInt N, i; 392b724ec41SToby Isaac PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u}; 3939371c9d4SSatish Balay void (*identity[1])(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u}; 394c4762a1bSJed Brown void *ctxs[1] = {0}; 395c4762a1bSJed Brown 396c4762a1bSJed Brown ctxs[0] = &user; 3979566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmErr)); 3989566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dmErr, "error", SetupErrorProblem, &user)); 3999566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErr, &errorEst)); 4009566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErr, &errorL2)); 401c4762a1bSJed Brown /* Compute auxiliary data (solution and projection of adjoint solution) */ 4029566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmAdj, &uAdjLoc)); 4039566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 4049566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 4059566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &uAdjProj)); 4069566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc)); 4079566063dSJacob Faibussowitsch PetscCall(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj)); 4089566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 4099566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmAdj, &uAdjLoc)); 410c4762a1bSJed Brown /* Attach auxiliary data */ 4119371c9d4SSatish Balay dms[0] = dm; 4129371c9d4SSatish Balay dms[1] = dm; 4139566063dSJacob Faibussowitsch PetscCall(DMCreateSuperDM(dms, 2, &subis, &dmErrAux)); 414c4762a1bSJed Brown if (0) { 415c4762a1bSJed Brown PetscSection sec; 416c4762a1bSJed Brown 4179566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[0], &sec)); 4189566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 4199566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[1], &sec)); 4209566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 4219566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmErrAux, &sec)); 4229566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 423c4762a1bSJed Brown } 4249566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view")); 4259566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(subis[0], NULL, "-super_is_view")); 4269566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(subis[1], NULL, "-super_is_view")); 4279566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErrAux, &uErr)); 4289566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-map_vec_view")); 4299566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view")); 4309566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uErr, NULL, "-map_vec_view")); 4319566063dSJacob Faibussowitsch PetscCall(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u)); 4329566063dSJacob Faibussowitsch PetscCall(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj)); 4339566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &uAdjProj)); 4349566063dSJacob Faibussowitsch for (i = 0; i < 2; ++i) PetscCall(ISDestroy(&subis[i])); 4359566063dSJacob Faibussowitsch PetscCall(PetscFree(subis)); 4369566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmErrAux, &uErrLoc)); 4379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc)); 4389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc)); 4399566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErrAux, &uErr)); 4409566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc)); 441c4762a1bSJed Brown /* Compute cellwise error estimate */ 4429566063dSJacob Faibussowitsch PetscCall(VecSet(errorEst, 0.0)); 4439566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user)); 4449566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL)); 4459566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmErrAux, &uErrLoc)); 4469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmErrAux)); 447c4762a1bSJed Brown /* Plot cellwise error vector */ 4489566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorEst, NULL, "-error_view")); 449c4762a1bSJed Brown /* Compute ratio of estimate (sum over cells) with actual L_2 error */ 4509566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm)); 4519566063dSJacob Faibussowitsch PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2)); 4529566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorL2, NULL, "-l2_error_view")); 4539566063dSJacob Faibussowitsch PetscCall(VecNorm(errorL2, NORM_INFINITY, &errorL2Tot)); 4549566063dSJacob Faibussowitsch PetscCall(VecNorm(errorEst, NORM_INFINITY, &errorEstTot)); 4559566063dSJacob Faibussowitsch PetscCall(VecGetSize(errorEst, &N)); 4569566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(errorEst, errorEst, errorL2)); 4579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)errorEst, "Error ratio")); 4589566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorEst, NULL, "-error_ratio_view")); 45963a3b9bcSJacob Faibussowitsch 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)))); 4609566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErr, &errorEst)); 4619566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErr, &errorL2)); 4629566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmErr)); 463c4762a1bSJed Brown } 4649566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAdj)); 4659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uAdj)); 4669566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snesAdj)); 467c4762a1bSJed Brown } 468c4762a1bSJed Brown /* Cleanup */ 4699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 4709566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4719566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4729566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 473b122ec5aSJacob Faibussowitsch return 0; 474c4762a1bSJed Brown } 475c4762a1bSJed Brown 476c4762a1bSJed Brown /*TEST 477c4762a1bSJed Brown 478c4762a1bSJed Brown test: 4795be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 4805be41b18SMatthew G. Knepley suffix: 2d_p1_conv 481c4762a1bSJed Brown requires: triangle 4825be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 4835be41b18SMatthew G. Knepley test: 4845be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 4855be41b18SMatthew G. Knepley suffix: 2d_p2_conv 4865be41b18SMatthew G. Knepley requires: triangle 4875be41b18SMatthew G. Knepley args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 4885be41b18SMatthew G. Knepley test: 4895be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 4905be41b18SMatthew G. Knepley suffix: 2d_p3_conv 4915be41b18SMatthew G. Knepley requires: triangle 4925be41b18SMatthew G. Knepley args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 4935be41b18SMatthew G. Knepley test: 4945be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 4955be41b18SMatthew G. Knepley suffix: 2d_q1_conv 49630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 4975be41b18SMatthew G. Knepley test: 4985be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 4995be41b18SMatthew G. Knepley suffix: 2d_q2_conv 50030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5015be41b18SMatthew G. Knepley test: 5025be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5035be41b18SMatthew G. Knepley suffix: 2d_q3_conv 50430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5055be41b18SMatthew G. Knepley test: 5065be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5075be41b18SMatthew G. Knepley suffix: 2d_q1_shear_conv 50830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5095be41b18SMatthew G. Knepley test: 5105be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5115be41b18SMatthew G. Knepley suffix: 2d_q2_shear_conv 51230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5135be41b18SMatthew G. Knepley test: 5145be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5155be41b18SMatthew G. Knepley suffix: 2d_q3_shear_conv 51630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5175be41b18SMatthew G. Knepley test: 5180fdc7489SMatthew Knepley # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7 5195be41b18SMatthew G. Knepley suffix: 3d_p1_conv 5205be41b18SMatthew G. Knepley requires: ctetgen 52130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 5225be41b18SMatthew G. Knepley test: 5235be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5245be41b18SMatthew G. Knepley suffix: 3d_p2_conv 5255be41b18SMatthew G. Knepley requires: ctetgen 52630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 5275be41b18SMatthew G. Knepley test: 5285be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0 5295be41b18SMatthew G. Knepley suffix: 3d_p3_conv 5305be41b18SMatthew G. Knepley requires: ctetgen 53130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 5325be41b18SMatthew G. Knepley test: 5335be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8 5345be41b18SMatthew G. Knepley suffix: 3d_q1_conv 53530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 5365be41b18SMatthew G. Knepley test: 5375be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5385be41b18SMatthew G. Knepley suffix: 3d_q2_conv 53930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 5405be41b18SMatthew G. Knepley test: 5415be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8 5425be41b18SMatthew G. Knepley suffix: 3d_q3_conv 54330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 5449139b27eSToby Isaac test: 5459139b27eSToby Isaac suffix: 2d_p1_fas_full 5469139b27eSToby Isaac requires: triangle 547e600fa54SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 5489139b27eSToby Isaac -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \ 5499139b27eSToby Isaac -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 5509139b27eSToby Isaac -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 5519139b27eSToby Isaac -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 55273f7197eSJed Brown -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10 5539139b27eSToby Isaac test: 5549139b27eSToby Isaac suffix: 2d_p1_fas_full_homogeneous 5559139b27eSToby Isaac requires: triangle 556e600fa54SMatthew G. Knepley args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 5579139b27eSToby Isaac -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \ 5589139b27eSToby Isaac -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 5599139b27eSToby Isaac -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 5609139b27eSToby Isaac -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 56173f7197eSJed Brown -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.1 -fas_levels_esteig_ksp_max_it 10 5625be41b18SMatthew G. Knepley 563c4762a1bSJed Brown test: 564c4762a1bSJed Brown suffix: 2d_p1_scalable 5655be41b18SMatthew G. Knepley requires: triangle 5665be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 3 \ 567c4762a1bSJed Brown -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \ 56873f7197eSJed Brown -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \ 569c4762a1bSJed Brown -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 570c4762a1bSJed Brown -pc_gamg_coarse_eq_limit 1000 \ 571c4762a1bSJed Brown -pc_gamg_threshold 0.05 \ 572c4762a1bSJed Brown -pc_gamg_threshold_scale .0 \ 573c4762a1bSJed Brown -mg_levels_ksp_type chebyshev \ 574c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 575c4762a1bSJed Brown -mg_levels_pc_type jacobi \ 576c4762a1bSJed Brown -matptap_via scalable 577c4762a1bSJed Brown test: 578c4762a1bSJed Brown suffix: 2d_p1_gmg_vcycle 579c4762a1bSJed Brown requires: triangle 5805be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 581c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg \ 582c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 583c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 584c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 58573f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 586c4762a1bSJed Brown -mg_levels_pc_type jacobi 587c4762a1bSJed Brown test: 588c4762a1bSJed Brown suffix: 2d_p1_gmg_fcycle 589c4762a1bSJed Brown requires: triangle 5905be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 591c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \ 592c4762a1bSJed Brown -mg_levels_ksp_max_it 2 \ 593c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 594c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 59573f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 596c4762a1bSJed Brown -mg_levels_pc_type jacobi 597c4762a1bSJed Brown test: 598d6837840SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_adapt 5992b3cbbdaSStefano Zampini requires: triangle 600908b9b43SStefano Zampini args: -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 6012b3cbbdaSStefano Zampini -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \ 602d6837840SMatthew G. Knepley -mg_levels_ksp_max_it 1 \ 603d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 604d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 60573f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 606d6837840SMatthew G. Knepley -mg_levels_pc_type jacobi 607d6837840SMatthew G. Knepley test: 608c4762a1bSJed Brown suffix: 2d_p1_spectral_0 609c4762a1bSJed Brown requires: triangle fftw !complex 6105be41b18SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view 611c4762a1bSJed Brown test: 612c4762a1bSJed Brown suffix: 2d_p1_spectral_1 613c4762a1bSJed Brown requires: triangle fftw !complex 614c4762a1bSJed Brown nsize: 2 615e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view 616c4762a1bSJed Brown test: 617c4762a1bSJed Brown suffix: 2d_p1_adj_0 618c4762a1bSJed Brown requires: triangle 6195be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0 620b1ad28e3SJunchao Zhang test: 621b1ad28e3SJunchao Zhang nsize: 2 622*dcfd994dSJunchao Zhang requires: kokkos_kernels 623b1ad28e3SJunchao Zhang suffix: kokkos 624e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,3,6 -petscpartitioner_type simple -dm_plex_simplex 0 -potential_petscspace_degree 1 \ 625b1ad28e3SJunchao Zhang -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 \ 62673f7197eSJed Brown -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 \ 62773f7197eSJed Brown -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 628c4762a1bSJed Brown 629c4762a1bSJed Brown TEST*/ 630