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 21c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 22c4762a1bSJed Brown { 23c4762a1bSJed Brown *u = 0.0; 24c4762a1bSJed Brown return 0; 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27b724ec41SToby Isaac static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown PetscInt d; 30c4762a1bSJed Brown *u = 0.0; 31c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 32c4762a1bSJed Brown return 0; 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35b724ec41SToby Isaac static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 36b724ec41SToby Isaac { 37b724ec41SToby Isaac PetscInt d; 38b724ec41SToby Isaac *u = 1.0; 39b724ec41SToby Isaac for (d = 0; d < dim; ++d) *u *= PetscSinReal(2.0*PETSC_PI*x[d]); 40b724ec41SToby Isaac return 0; 41b724ec41SToby Isaac } 42b724ec41SToby Isaac 43c4762a1bSJed Brown /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */ 44c4762a1bSJed Brown static void obj_error_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 45c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 46c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 47c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[]) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown obj[0] = a[aOff[0]]*(u[0] - a[aOff[1]]); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52b724ec41SToby Isaac static void f0_trig_inhomogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 53c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 54c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 55c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 56c4762a1bSJed Brown { 57c4762a1bSJed Brown PetscInt d; 58c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61b724ec41SToby Isaac static void f0_trig_homogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 62b724ec41SToby Isaac const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 63b724ec41SToby Isaac const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 64b724ec41SToby Isaac PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 65b724ec41SToby Isaac { 66b724ec41SToby Isaac PetscInt d; 67b724ec41SToby Isaac for (d = 0; d < dim; ++d) { 68b724ec41SToby Isaac PetscScalar v = 1.; 69b724ec41SToby Isaac for (PetscInt e = 0; e < dim; e++) { 70b724ec41SToby Isaac if (e == d) { 71b724ec41SToby Isaac v *= -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 72b724ec41SToby Isaac } else { 73b724ec41SToby Isaac v *= PetscSinReal(2.0*PETSC_PI*x[d]); 74b724ec41SToby Isaac } 75b724ec41SToby Isaac } 76b724ec41SToby Isaac f0[0] += v; 77b724ec41SToby Isaac } 78b724ec41SToby Isaac } 79b724ec41SToby Isaac 80c4762a1bSJed Brown static void f0_unity_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 81c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 82c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 83c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 84c4762a1bSJed Brown { 85c4762a1bSJed Brown f0[0] = 1.0; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown static void f0_identityaux_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 89c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 90c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 91c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 92c4762a1bSJed Brown { 93c4762a1bSJed Brown f0[0] = a[0]; 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 97c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 98c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 99c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 100c4762a1bSJed Brown { 101c4762a1bSJed Brown PetscInt d; 102c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 106c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 107c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 108c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 109c4762a1bSJed Brown { 110c4762a1bSJed Brown PetscInt d; 111c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 115c4762a1bSJed Brown { 116c4762a1bSJed Brown PetscFunctionBeginUser; 117c4762a1bSJed Brown options->shear = PETSC_FALSE; 118c4762a1bSJed Brown options->spectral = PETSC_FALSE; 119c4762a1bSJed Brown options->adjoint = PETSC_FALSE; 120b724ec41SToby Isaac options->homogeneous = PETSC_FALSE; 121ad1a7c93SMatthew Knepley options->viewError = PETSC_FALSE; 122c4762a1bSJed Brown 123d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 1249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL)); 1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL)); 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL)); 1279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL)); 1289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL)); 129d0609cedSBarry Smith PetscOptionsEnd(); 130c4762a1bSJed Brown PetscFunctionReturn(0); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 134c4762a1bSJed Brown { 135c4762a1bSJed Brown PetscSection coordSection; 136c4762a1bSJed Brown Vec coordinates; 137c4762a1bSJed Brown const PetscScalar *coords; 138c4762a1bSJed Brown PetscInt dim, p, vStart, vEnd, v; 139c4762a1bSJed Brown 140c4762a1bSJed Brown PetscFunctionBeginUser; 1419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dim)); 1429566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 146c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 147c4762a1bSJed Brown DMLabel label; 148c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 149c4762a1bSJed Brown 15063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 1519566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, name)); 1529566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 1539566063dSJacob Faibussowitsch PetscCall(DMLabelAddStratum(label, 1)); 154c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 155c4762a1bSJed Brown PetscInt off; 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 158c4762a1bSJed Brown if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off+planeDir[p]])) < PETSC_SMALL) { 1599566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(label, v, 1)); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown } 1639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown PetscFunctionBeginUser; 1709566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1719566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1729566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1739566063dSJacob Faibussowitsch if (user->shear) PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL)); 1749566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 1759566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 176c4762a1bSJed Brown if (user->spectral) { 177c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 178c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 179c4762a1bSJed Brown 1809566063dSJacob Faibussowitsch PetscCall(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user)); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown PetscFunctionReturn(0); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 186c4762a1bSJed Brown { 18745480ffeSMatthew G. Knepley PetscDS ds; 18845480ffeSMatthew G. Knepley DMLabel label; 189c4762a1bSJed Brown const PetscInt id = 1; 190b724ec41SToby Isaac PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u; 191b724ec41SToby Isaac PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u; 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscFunctionBeginUser; 1949566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1959566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0, f1_u)); 1969566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, ex, user)); 1989566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 1999566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, user, NULL)); 200c4762a1bSJed Brown PetscFunctionReturn(0); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203c4762a1bSJed Brown static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user) 204c4762a1bSJed Brown { 20545480ffeSMatthew G. Knepley PetscDS ds; 20645480ffeSMatthew G. Knepley DMLabel label; 207c4762a1bSJed Brown const PetscInt id = 1; 208c4762a1bSJed Brown 209c4762a1bSJed Brown PetscFunctionBeginUser; 2109566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2119566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u)); 2129566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 2139566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, obj_error_u)); 2149566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2159566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, user, NULL)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown PetscDS prob; 222c4762a1bSJed Brown 223c4762a1bSJed Brown PetscFunctionBeginUser; 2249566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 225c4762a1bSJed Brown PetscFunctionReturn(0); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 229c4762a1bSJed Brown { 230c4762a1bSJed Brown DM cdm = dm; 231c4762a1bSJed Brown PetscFE fe; 2325be41b18SMatthew G. Knepley DMPolytopeType ct; 2335be41b18SMatthew G. Knepley PetscBool simplex; 2345be41b18SMatthew G. Knepley PetscInt dim, cStart; 235c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 236c4762a1bSJed Brown 237c4762a1bSJed Brown PetscFunctionBeginUser; 2389566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2399566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2409566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 2415be41b18SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 242c4762a1bSJed Brown /* Create finite element */ 2439566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 2449566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 246c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2479566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 2489566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2499566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 250c4762a1bSJed Brown while (cdm) { 2519566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm,cdm)); 252c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 2539566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 254c4762a1bSJed Brown } 2559566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 256c4762a1bSJed Brown PetscFunctionReturn(0); 257c4762a1bSJed Brown } 258c4762a1bSJed Brown 259c4762a1bSJed Brown static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 260c4762a1bSJed Brown { 261c4762a1bSJed Brown MPI_Comm comm; 262c4762a1bSJed Brown PetscSection coordSection, section; 263c4762a1bSJed Brown Vec coordinates, uloc; 264c4762a1bSJed Brown const PetscScalar *coords, *array; 265c4762a1bSJed Brown PetscInt p; 266c4762a1bSJed Brown PetscMPIInt size, rank; 267c4762a1bSJed Brown 268c4762a1bSJed Brown PetscFunctionBeginUser; 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 2709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2729566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &uloc)); 2739566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc)); 2749566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc)); 2759566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL)); 2769566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uloc, NULL, "-sol_view")); 2779566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 2789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(uloc, &array)); 2799566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2809566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 282c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 283c4762a1bSJed Brown DMLabel label; 284c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 285c4762a1bSJed Brown Mat F; 286c4762a1bSJed Brown Vec x, y; 287c4762a1bSJed Brown IS stratum; 288c4762a1bSJed Brown PetscReal *ray, *gray; 289c4762a1bSJed Brown PetscScalar *rvals, *svals, *gsvals; 290c4762a1bSJed Brown PetscInt *perm, *nperm; 291c4762a1bSJed Brown PetscInt n, N, i, j, off, offu; 292c4762a1bSJed Brown const PetscInt *points; 293c4762a1bSJed Brown 29463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 2959566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 2969566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, 1, &stratum)); 2979566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(stratum, &n)); 2989566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratum, &points)); 2999566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ray, n, &svals)); 300c4762a1bSJed Brown for (i = 0; i < n; ++i) { 3019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, points[i], &off)); 3029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, points[i], &offu)); 303c4762a1bSJed Brown ray[i] = PetscRealPart(coords[off+((planeDir[p]+1)%2)]); 304c4762a1bSJed Brown svals[i] = array[offu]; 305c4762a1bSJed Brown } 306c4762a1bSJed Brown /* Gather the ray data to proc 0 */ 307c4762a1bSJed Brown if (size > 1) { 308c4762a1bSJed Brown PetscMPIInt *cnt, *displs, p; 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(size, &cnt, size, &displs)); 3119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm)); 312c4762a1bSJed Brown for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1]; 313c4762a1bSJed Brown N = displs[size-1] + cnt[size-1]; 3149566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &gray, N, &gsvals)); 3159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm)); 3169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm)); 3179566063dSJacob Faibussowitsch PetscCall(PetscFree2(cnt, displs)); 318c4762a1bSJed Brown } else { 319c4762a1bSJed Brown N = n; 320c4762a1bSJed Brown gray = ray; 321c4762a1bSJed Brown gsvals = svals; 322c4762a1bSJed Brown } 323dd400576SPatrick Sanan if (rank == 0) { 324c4762a1bSJed Brown /* Sort point along ray */ 3259566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &perm, N, &nperm)); 326c4762a1bSJed Brown for (i = 0; i < N; ++i) {perm[i] = i;} 3279566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithPermutation(N, gray, perm)); 328c4762a1bSJed Brown /* Count duplicates and squish mapping */ 329c4762a1bSJed Brown nperm[0] = perm[0]; 330c4762a1bSJed Brown for (i = 1, j = 1; i < N; ++i) { 331c4762a1bSJed Brown if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i]; 332c4762a1bSJed Brown } 333c4762a1bSJed Brown /* Create FFT structs */ 3349566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F)); 3359566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(F, &x, &y)); 3369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) y, name)); 3379566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &rvals)); 338c4762a1bSJed Brown for (i = 0, j = 0; i < N; ++i) { 339c4762a1bSJed Brown if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue; 340c4762a1bSJed Brown rvals[j] = gsvals[nperm[j]]; 341c4762a1bSJed Brown ++j; 342c4762a1bSJed Brown } 3439566063dSJacob Faibussowitsch PetscCall(PetscFree2(perm, nperm)); 3449566063dSJacob Faibussowitsch if (size > 1) PetscCall(PetscFree2(gray, gsvals)); 3459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &rvals)); 346c4762a1bSJed Brown /* Do FFT along the ray */ 3479566063dSJacob Faibussowitsch PetscCall(MatMult(F, x, y)); 348c4762a1bSJed Brown /* Chop FFT */ 3499566063dSJacob Faibussowitsch PetscCall(VecChop(y, PETSC_SMALL)); 3509566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(x, NULL, "-real_view")); 3519566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(y, NULL, "-fft_view")); 3529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 3549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 355c4762a1bSJed Brown } 3569566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratum, &points)); 3579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratum)); 3589566063dSJacob Faibussowitsch PetscCall(PetscFree2(ray, svals)); 359c4762a1bSJed Brown } 3609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 3619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(uloc, &array)); 3629566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &uloc)); 363c4762a1bSJed Brown PetscFunctionReturn(0); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown 366c4762a1bSJed Brown int main(int argc, char **argv) 367c4762a1bSJed Brown { 368c4762a1bSJed Brown DM dm; /* Problem specification */ 369c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 370c4762a1bSJed Brown Vec u; /* Solutions */ 371c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 372c4762a1bSJed Brown 3739566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3749566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 375c4762a1bSJed Brown /* Primal system */ 3769566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 3779566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 3789566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 3799566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 3809566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3819566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 3839566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 3849566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 3859566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3869566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 3879566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 388ad1a7c93SMatthew Knepley if (user.viewError) { 389ad1a7c93SMatthew Knepley PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 390ad1a7c93SMatthew Knepley void *ctx; 391ad1a7c93SMatthew Knepley PetscDS ds; 392ad1a7c93SMatthew Knepley PetscReal error; 393ad1a7c93SMatthew Knepley PetscInt N; 394ad1a7c93SMatthew Knepley 3959566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3969566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &sol, &ctx)); 3979566063dSJacob Faibussowitsch PetscCall(VecGetSize(u, &N)); 3989566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error)); 39963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g\n", N, (double)error)); 400ad1a7c93SMatthew Knepley } 401c4762a1bSJed Brown if (user.spectral) { 402c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 403c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 404c4762a1bSJed Brown 4059566063dSJacob Faibussowitsch PetscCall(ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user)); 406c4762a1bSJed Brown } 407c4762a1bSJed Brown /* Adjoint system */ 408c4762a1bSJed Brown if (user.adjoint) { 409c4762a1bSJed Brown DM dmAdj; 410c4762a1bSJed Brown SNES snesAdj; 411c4762a1bSJed Brown Vec uAdj; 412c4762a1bSJed Brown 4139566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snesAdj)); 4149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_")); 4159566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAdj)); 4169566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snesAdj, dmAdj)); 4179566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user)); 4189566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmAdj, &uAdj)); 4199566063dSJacob Faibussowitsch PetscCall(VecSet(uAdj, 0.0)); 4209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) uAdj, "adjoint")); 4219566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user)); 4229566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snesAdj)); 4239566063dSJacob Faibussowitsch PetscCall(SNESSolve(snesAdj, NULL, uAdj)); 4249566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snesAdj, &uAdj)); 4259566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uAdj, NULL, "-adjoint_view")); 426c4762a1bSJed Brown /* Error representation */ 427c4762a1bSJed Brown { 428c4762a1bSJed Brown DM dmErr, dmErrAux, dms[2]; 429c4762a1bSJed Brown Vec errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj; 430c4762a1bSJed Brown IS *subis; 431c4762a1bSJed Brown PetscReal errorEstTot, errorL2Norm, errorL2Tot; 432c4762a1bSJed Brown PetscInt N, i; 433b724ec41SToby Isaac PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u}; 434c4762a1bSJed Brown void (*identity[1])(PetscInt, PetscInt, PetscInt, 435c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 436c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 437c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u}; 438c4762a1bSJed Brown void *ctxs[1] = {0}; 439c4762a1bSJed Brown 440c4762a1bSJed Brown ctxs[0] = &user; 4419566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmErr)); 4429566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dmErr, "error", SetupErrorProblem, &user)); 4439566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErr, &errorEst)); 4449566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErr, &errorL2)); 445c4762a1bSJed Brown /* Compute auxiliary data (solution and projection of adjoint solution) */ 4469566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmAdj, &uAdjLoc)); 4479566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 4489566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 4499566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &uAdjProj)); 4509566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc)); 4519566063dSJacob Faibussowitsch PetscCall(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj)); 4529566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 4539566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmAdj, &uAdjLoc)); 454c4762a1bSJed Brown /* Attach auxiliary data */ 455c4762a1bSJed Brown dms[0] = dm; dms[1] = dm; 4569566063dSJacob Faibussowitsch PetscCall(DMCreateSuperDM(dms, 2, &subis, &dmErrAux)); 457c4762a1bSJed Brown if (0) { 458c4762a1bSJed Brown PetscSection sec; 459c4762a1bSJed Brown 4609566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[0], &sec)); 4619566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 4629566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[1], &sec)); 4639566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 4649566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmErrAux, &sec)); 4659566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 466c4762a1bSJed Brown } 4679566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view")); 4689566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(subis[0], NULL, "-super_is_view")); 4699566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(subis[1], NULL, "-super_is_view")); 4709566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErrAux, &uErr)); 4719566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-map_vec_view")); 4729566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view")); 4739566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uErr, NULL, "-map_vec_view")); 4749566063dSJacob Faibussowitsch PetscCall(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u)); 4759566063dSJacob Faibussowitsch PetscCall(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj)); 4769566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &uAdjProj)); 4779566063dSJacob Faibussowitsch for (i = 0; i < 2; ++i) PetscCall(ISDestroy(&subis[i])); 4789566063dSJacob Faibussowitsch PetscCall(PetscFree(subis)); 4799566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmErrAux, &uErrLoc)); 4809566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc)); 4819566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc)); 4829566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErrAux, &uErr)); 4839566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc)); 484c4762a1bSJed Brown /* Compute cellwise error estimate */ 4859566063dSJacob Faibussowitsch PetscCall(VecSet(errorEst, 0.0)); 4869566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user)); 4879566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL)); 4889566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmErrAux, &uErrLoc)); 4899566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmErrAux)); 490c4762a1bSJed Brown /* Plot cellwise error vector */ 4919566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorEst, NULL, "-error_view")); 492c4762a1bSJed Brown /* Compute ratio of estimate (sum over cells) with actual L_2 error */ 4939566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm)); 4949566063dSJacob Faibussowitsch PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2)); 4959566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorL2, NULL, "-l2_error_view")); 4969566063dSJacob Faibussowitsch PetscCall(VecNorm(errorL2, NORM_INFINITY, &errorL2Tot)); 4979566063dSJacob Faibussowitsch PetscCall(VecNorm(errorEst, NORM_INFINITY, &errorEstTot)); 4989566063dSJacob Faibussowitsch PetscCall(VecGetSize(errorEst, &N)); 4999566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(errorEst, errorEst, errorL2)); 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) errorEst, "Error ratio")); 5019566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorEst, NULL, "-error_ratio_view")); 50263a3b9bcSJacob 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)))); 5039566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErr, &errorEst)); 5049566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErr, &errorL2)); 5059566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmErr)); 506c4762a1bSJed Brown } 5079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAdj)); 5089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uAdj)); 5099566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snesAdj)); 510c4762a1bSJed Brown } 511c4762a1bSJed Brown /* Cleanup */ 5129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5139566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 5149566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5159566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 516b122ec5aSJacob Faibussowitsch return 0; 517c4762a1bSJed Brown } 518c4762a1bSJed Brown 519c4762a1bSJed Brown /*TEST 520c4762a1bSJed Brown 521c4762a1bSJed Brown test: 5225be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5235be41b18SMatthew G. Knepley suffix: 2d_p1_conv 524c4762a1bSJed Brown requires: triangle 5255be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5265be41b18SMatthew G. Knepley test: 5275be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5285be41b18SMatthew G. Knepley suffix: 2d_p2_conv 5295be41b18SMatthew G. Knepley requires: triangle 5305be41b18SMatthew G. Knepley args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5315be41b18SMatthew G. Knepley test: 5325be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5335be41b18SMatthew G. Knepley suffix: 2d_p3_conv 5345be41b18SMatthew G. Knepley requires: triangle 5355be41b18SMatthew G. Knepley args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5365be41b18SMatthew G. Knepley test: 5375be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5385be41b18SMatthew G. Knepley suffix: 2d_q1_conv 53930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5405be41b18SMatthew G. Knepley test: 5415be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5425be41b18SMatthew G. Knepley suffix: 2d_q2_conv 54330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5445be41b18SMatthew G. Knepley test: 5455be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5465be41b18SMatthew G. Knepley suffix: 2d_q3_conv 54730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5485be41b18SMatthew G. Knepley test: 5495be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5505be41b18SMatthew G. Knepley suffix: 2d_q1_shear_conv 55130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5525be41b18SMatthew G. Knepley test: 5535be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5545be41b18SMatthew G. Knepley suffix: 2d_q2_shear_conv 55530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5565be41b18SMatthew G. Knepley test: 5575be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5585be41b18SMatthew G. Knepley suffix: 2d_q3_shear_conv 55930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5605be41b18SMatthew G. Knepley test: 5610fdc7489SMatthew Knepley # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7 5625be41b18SMatthew G. Knepley suffix: 3d_p1_conv 5635be41b18SMatthew G. Knepley requires: ctetgen 56430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 5655be41b18SMatthew G. Knepley test: 5665be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5675be41b18SMatthew G. Knepley suffix: 3d_p2_conv 5685be41b18SMatthew G. Knepley requires: ctetgen 56930602db0SMatthew 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 5705be41b18SMatthew G. Knepley test: 5715be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0 5725be41b18SMatthew G. Knepley suffix: 3d_p3_conv 5735be41b18SMatthew G. Knepley requires: ctetgen 57430602db0SMatthew 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 5755be41b18SMatthew G. Knepley test: 5765be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8 5775be41b18SMatthew G. Knepley suffix: 3d_q1_conv 57830602db0SMatthew 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 5795be41b18SMatthew G. Knepley test: 5805be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5815be41b18SMatthew G. Knepley suffix: 3d_q2_conv 58230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 5835be41b18SMatthew G. Knepley test: 5845be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8 5855be41b18SMatthew G. Knepley suffix: 3d_q3_conv 58630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 5879139b27eSToby Isaac test: 5889139b27eSToby Isaac suffix: 2d_p1_fas_full 5899139b27eSToby Isaac requires: triangle 590e600fa54SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 5919139b27eSToby Isaac -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \ 5929139b27eSToby Isaac -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 5939139b27eSToby Isaac -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 5949139b27eSToby Isaac -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 59573f7197eSJed 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 5969139b27eSToby Isaac test: 5979139b27eSToby Isaac suffix: 2d_p1_fas_full_homogeneous 5989139b27eSToby Isaac requires: triangle 599e600fa54SMatthew G. Knepley args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 6009139b27eSToby Isaac -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \ 6019139b27eSToby Isaac -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 6029139b27eSToby Isaac -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 6039139b27eSToby Isaac -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 60473f7197eSJed 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 6055be41b18SMatthew G. Knepley 606c4762a1bSJed Brown test: 607c4762a1bSJed Brown suffix: 2d_p1_scalable 6085be41b18SMatthew G. Knepley requires: triangle 6095be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 3 \ 610c4762a1bSJed Brown -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \ 61173f7197eSJed Brown -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \ 612c4762a1bSJed Brown -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 613c4762a1bSJed Brown -pc_gamg_coarse_eq_limit 1000 \ 614c4762a1bSJed Brown -pc_gamg_square_graph 1 \ 615c4762a1bSJed Brown -pc_gamg_threshold 0.05 \ 616c4762a1bSJed Brown -pc_gamg_threshold_scale .0 \ 617c4762a1bSJed Brown -mg_levels_ksp_type chebyshev \ 618c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 619c4762a1bSJed Brown -mg_levels_pc_type jacobi \ 620c4762a1bSJed Brown -matptap_via scalable 621c4762a1bSJed Brown test: 622c4762a1bSJed Brown suffix: 2d_p1_gmg_vcycle 623c4762a1bSJed Brown requires: triangle 6245be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 625c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg \ 626c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 627c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 628c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 62973f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 630c4762a1bSJed Brown -mg_levels_pc_type jacobi 631c4762a1bSJed Brown test: 632c4762a1bSJed Brown suffix: 2d_p1_gmg_fcycle 633c4762a1bSJed Brown requires: triangle 6345be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 635c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \ 636c4762a1bSJed Brown -mg_levels_ksp_max_it 2 \ 637c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 638c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 63973f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 640c4762a1bSJed Brown -mg_levels_pc_type jacobi 641c4762a1bSJed Brown test: 642d6837840SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_adapt 643*2b3cbbdaSStefano Zampini requires: triangle 6445be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 645*2b3cbbdaSStefano Zampini -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \ 646d6837840SMatthew G. Knepley -mg_levels_ksp_max_it 1 \ 647d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 648d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 64973f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 650d6837840SMatthew G. Knepley -mg_levels_pc_type jacobi 651d6837840SMatthew G. Knepley test: 652c4762a1bSJed Brown suffix: 2d_p1_spectral_0 653c4762a1bSJed Brown requires: triangle fftw !complex 6545be41b18SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view 655c4762a1bSJed Brown test: 656c4762a1bSJed Brown suffix: 2d_p1_spectral_1 657c4762a1bSJed Brown requires: triangle fftw !complex 658c4762a1bSJed Brown nsize: 2 659e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view 660c4762a1bSJed Brown test: 661c4762a1bSJed Brown suffix: 2d_p1_adj_0 662c4762a1bSJed Brown requires: triangle 6635be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0 664b1ad28e3SJunchao Zhang test: 665b1ad28e3SJunchao Zhang nsize: 2 6663078479eSJunchao Zhang requires: !sycl kokkos_kernels 667b1ad28e3SJunchao Zhang suffix: kokkos 668e600fa54SMatthew 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 \ 669b1ad28e3SJunchao 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 \ 67073f7197eSJed 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 \ 67173f7197eSJed Brown -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 672c4762a1bSJed Brown 673c4762a1bSJed Brown TEST*/ 674