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 */ 17da81f932SPierre Jolivet PetscBool homogeneous; /* Use homogeneous boundary conditions */ 18ad1a7c93SMatthew Knepley PetscBool viewError; /* Output the solution error */ 19c4762a1bSJed Brown } AppCtx; 20c4762a1bSJed Brown 21d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 22d71ae5a4SJacob Faibussowitsch { 23c4762a1bSJed Brown *u = 0.0; 243ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 28d71ae5a4SJacob Faibussowitsch { 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]); 323ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_homogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 36d71ae5a4SJacob Faibussowitsch { 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]); 403ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 41b724ec41SToby Isaac } 42b724ec41SToby Isaac 43c4762a1bSJed Brown /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */ 44d71ae5a4SJacob Faibussowitsch 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[]) 45d71ae5a4SJacob Faibussowitsch { 46c4762a1bSJed Brown obj[0] = a[aOff[0]] * (u[0] - a[aOff[1]]); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49d71ae5a4SJacob Faibussowitsch 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[]) 50d71ae5a4SJacob Faibussowitsch { 51c4762a1bSJed Brown PetscInt d; 52c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55d71ae5a4SJacob Faibussowitsch 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[]) 56d71ae5a4SJacob Faibussowitsch { 57b724ec41SToby Isaac PetscInt d; 58b724ec41SToby Isaac for (d = 0; d < dim; ++d) { 59b724ec41SToby Isaac PetscScalar v = 1.; 60b724ec41SToby Isaac for (PetscInt e = 0; e < dim; e++) { 61b724ec41SToby Isaac if (e == d) { 62b724ec41SToby Isaac v *= -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 63b724ec41SToby Isaac } else { 64b724ec41SToby Isaac v *= PetscSinReal(2.0 * PETSC_PI * x[d]); 65b724ec41SToby Isaac } 66b724ec41SToby Isaac } 67b724ec41SToby Isaac f0[0] += v; 68b724ec41SToby Isaac } 69b724ec41SToby Isaac } 70b724ec41SToby Isaac 71d71ae5a4SJacob Faibussowitsch 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[]) 72d71ae5a4SJacob Faibussowitsch { 73c4762a1bSJed Brown f0[0] = 1.0; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76d71ae5a4SJacob Faibussowitsch 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[]) 77d71ae5a4SJacob Faibussowitsch { 78c4762a1bSJed Brown f0[0] = a[0]; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81d71ae5a4SJacob Faibussowitsch 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[]) 82d71ae5a4SJacob Faibussowitsch { 83c4762a1bSJed Brown PetscInt d; 84c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87d71ae5a4SJacob Faibussowitsch 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[]) 88d71ae5a4SJacob Faibussowitsch { 89c4762a1bSJed Brown PetscInt d; 90c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 94d71ae5a4SJacob Faibussowitsch { 95c4762a1bSJed Brown PetscFunctionBeginUser; 96c4762a1bSJed Brown options->shear = PETSC_FALSE; 97c4762a1bSJed Brown options->spectral = PETSC_FALSE; 98c4762a1bSJed Brown options->adjoint = PETSC_FALSE; 99b724ec41SToby Isaac options->homogeneous = PETSC_FALSE; 100ad1a7c93SMatthew Knepley options->viewError = PETSC_FALSE; 101c4762a1bSJed Brown 102d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL)); 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL)); 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL)); 108d0609cedSBarry Smith PetscOptionsEnd(); 1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 113d71ae5a4SJacob Faibussowitsch { 114c4762a1bSJed Brown PetscSection coordSection; 115c4762a1bSJed Brown Vec coordinates; 116c4762a1bSJed Brown const PetscScalar *coords; 117c4762a1bSJed Brown PetscInt dim, p, vStart, vEnd, v; 118c4762a1bSJed Brown 119c4762a1bSJed Brown PetscFunctionBeginUser; 1209566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dim)); 1219566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1229566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1239566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 125c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 126c4762a1bSJed Brown DMLabel label; 127c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 128c4762a1bSJed Brown 12963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 1309566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, name)); 1319566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 1329566063dSJacob Faibussowitsch PetscCall(DMLabelAddStratum(label, 1)); 133c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 134c4762a1bSJed Brown PetscInt off; 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 13748a46eb9SPierre Jolivet if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off + planeDir[p]])) < PETSC_SMALL) PetscCall(DMLabelSetValue(label, v, 1)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown } 1409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 145d71ae5a4SJacob Faibussowitsch { 146c4762a1bSJed Brown PetscFunctionBeginUser; 1479566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1489566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1499566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1509566063dSJacob Faibussowitsch if (user->shear) PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL)); 1519566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 1529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 153c4762a1bSJed Brown if (user->spectral) { 154c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 155c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user)); 158c4762a1bSJed Brown } 1593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 163d71ae5a4SJacob Faibussowitsch { 16445480ffeSMatthew G. Knepley PetscDS ds; 16545480ffeSMatthew G. Knepley DMLabel label; 166c4762a1bSJed Brown const PetscInt id = 1; 167b724ec41SToby Isaac PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u; 168b724ec41SToby Isaac PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBeginUser; 1719566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1729566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0, f1_u)); 1739566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1749566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, ex, user)); 1759566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 176c78f73d1SMatthew G. Knepley if (label) PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))ex, NULL, user, NULL)); 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178c4762a1bSJed Brown } 179c4762a1bSJed Brown 180d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user) 181d71ae5a4SJacob Faibussowitsch { 18245480ffeSMatthew G. Knepley PetscDS ds; 18345480ffeSMatthew G. Knepley DMLabel label; 184c4762a1bSJed Brown const PetscInt id = 1; 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscFunctionBeginUser; 1879566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1889566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_unity_u, f1_u)); 1899566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 1909566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, obj_error_u)); 1919566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 1929566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, user, NULL)); 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user) 197d71ae5a4SJacob Faibussowitsch { 198c4762a1bSJed Brown PetscDS prob; 199c4762a1bSJed Brown 200c4762a1bSJed Brown PetscFunctionBeginUser; 2019566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 206d71ae5a4SJacob Faibussowitsch { 207c4762a1bSJed Brown DM cdm = dm; 208c4762a1bSJed Brown PetscFE fe; 2095be41b18SMatthew G. Knepley DMPolytopeType ct; 2105be41b18SMatthew G. Knepley PetscBool simplex; 2115be41b18SMatthew G. Knepley PetscInt dim, cStart; 212c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 213c4762a1bSJed Brown 214c4762a1bSJed Brown PetscFunctionBeginUser; 2159566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2169566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2179566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 2185be41b18SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 219c4762a1bSJed Brown /* Create finite element */ 2209566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 2219566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name)); 223c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2249566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2259566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2269566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 227c4762a1bSJed Brown while (cdm) { 2289566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 229c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 2309566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 231c4762a1bSJed Brown } 2329566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236*e8e188d2SZach Atkins static PetscErrorCode ComputeSpectral(Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 237d71ae5a4SJacob Faibussowitsch { 238c4762a1bSJed Brown MPI_Comm comm; 239*e8e188d2SZach Atkins DM dm; 240c4762a1bSJed Brown PetscSection coordSection, section; 241c4762a1bSJed Brown Vec coordinates, uloc; 242c4762a1bSJed Brown const PetscScalar *coords, *array; 243c4762a1bSJed Brown PetscInt p; 244c4762a1bSJed Brown PetscMPIInt size, rank; 245c4762a1bSJed Brown 246c4762a1bSJed Brown PetscFunctionBeginUser; 247*e8e188d2SZach Atkins if (!user->spectral) PetscFunctionReturn(PETSC_SUCCESS); 248*e8e188d2SZach Atkins PetscCall(VecGetDM(u, &dm)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2529566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &uloc)); 2539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc)); 2549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc)); 2559566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL)); 2569566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uloc, NULL, "-sol_view")); 2579566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(uloc, &array)); 2599566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 2619566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 262c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 263c4762a1bSJed Brown DMLabel label; 264c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 265c4762a1bSJed Brown Mat F; 266c4762a1bSJed Brown Vec x, y; 267c4762a1bSJed Brown IS stratum; 268c4762a1bSJed Brown PetscReal *ray, *gray; 269c4762a1bSJed Brown PetscScalar *rvals, *svals, *gsvals; 270c4762a1bSJed Brown PetscInt *perm, *nperm; 271c4762a1bSJed Brown PetscInt n, N, i, j, off, offu; 272c4762a1bSJed Brown const PetscInt *points; 273c4762a1bSJed Brown 27463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%" PetscInt_FMT, p)); 2759566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 2769566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, 1, &stratum)); 2779566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(stratum, &n)); 2789566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratum, &points)); 2799566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &ray, n, &svals)); 280c4762a1bSJed Brown for (i = 0; i < n; ++i) { 2819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, points[i], &off)); 2829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, points[i], &offu)); 283c4762a1bSJed Brown ray[i] = PetscRealPart(coords[off + ((planeDir[p] + 1) % 2)]); 284c4762a1bSJed Brown svals[i] = array[offu]; 285c4762a1bSJed Brown } 286c4762a1bSJed Brown /* Gather the ray data to proc 0 */ 287c4762a1bSJed Brown if (size > 1) { 288c4762a1bSJed Brown PetscMPIInt *cnt, *displs, p; 289c4762a1bSJed Brown 2909566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(size, &cnt, size, &displs)); 2919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm)); 292c4762a1bSJed Brown for (p = 1; p < size; ++p) displs[p] = displs[p - 1] + cnt[p - 1]; 293c4762a1bSJed Brown N = displs[size - 1] + cnt[size - 1]; 2949566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &gray, N, &gsvals)); 2959566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm)); 2969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm)); 2979566063dSJacob Faibussowitsch PetscCall(PetscFree2(cnt, displs)); 298c4762a1bSJed Brown } else { 299c4762a1bSJed Brown N = n; 300c4762a1bSJed Brown gray = ray; 301c4762a1bSJed Brown gsvals = svals; 302c4762a1bSJed Brown } 303dd400576SPatrick Sanan if (rank == 0) { 304c4762a1bSJed Brown /* Sort point along ray */ 3059566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &perm, N, &nperm)); 306ad540459SPierre Jolivet for (i = 0; i < N; ++i) perm[i] = i; 3079566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithPermutation(N, gray, perm)); 308c4762a1bSJed Brown /* Count duplicates and squish mapping */ 309c4762a1bSJed Brown nperm[0] = perm[0]; 310c4762a1bSJed Brown for (i = 1, j = 1; i < N; ++i) { 311c4762a1bSJed Brown if (PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) > PETSC_SMALL) nperm[j++] = perm[i]; 312c4762a1bSJed Brown } 313c4762a1bSJed Brown /* Create FFT structs */ 3149566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F)); 3159566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(F, &x, &y)); 3169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, name)); 3179566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &rvals)); 318c4762a1bSJed Brown for (i = 0, j = 0; i < N; ++i) { 319c4762a1bSJed Brown if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i - 1]]) < PETSC_SMALL) continue; 320c4762a1bSJed Brown rvals[j] = gsvals[nperm[j]]; 321c4762a1bSJed Brown ++j; 322c4762a1bSJed Brown } 3239566063dSJacob Faibussowitsch PetscCall(PetscFree2(perm, nperm)); 3249566063dSJacob Faibussowitsch if (size > 1) PetscCall(PetscFree2(gray, gsvals)); 3259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &rvals)); 326c4762a1bSJed Brown /* Do FFT along the ray */ 3279566063dSJacob Faibussowitsch PetscCall(MatMult(F, x, y)); 328c4762a1bSJed Brown /* Chop FFT */ 32993ec0da9SPierre Jolivet PetscCall(VecFilter(y, PETSC_SMALL)); 3309566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(x, NULL, "-real_view")); 3319566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(y, NULL, "-fft_view")); 3329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 3349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 335c4762a1bSJed Brown } 3369566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratum, &points)); 3379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratum)); 3389566063dSJacob Faibussowitsch PetscCall(PetscFree2(ray, svals)); 339c4762a1bSJed Brown } 3409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 3419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(uloc, &array)); 3429566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &uloc)); 3433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 346*e8e188d2SZach Atkins static PetscErrorCode ComputeAdjoint(Vec u, AppCtx *user) 347d71ae5a4SJacob Faibussowitsch { 348*e8e188d2SZach Atkins PetscFunctionBegin; 349*e8e188d2SZach Atkins if (!user->adjoint) PetscFunctionReturn(PETSC_SUCCESS); 350*e8e188d2SZach Atkins DM dm, dmAdj; 351c4762a1bSJed Brown SNES snesAdj; 352c4762a1bSJed Brown Vec uAdj; 353c4762a1bSJed Brown 354*e8e188d2SZach Atkins PetscCall(VecGetDM(u, &dm)); 3559566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snesAdj)); 3569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)snesAdj, "adjoint_")); 3579566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAdj)); 3589566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snesAdj, dmAdj)); 359*e8e188d2SZach Atkins PetscCall(SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, user)); 3609566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dmAdj, &uAdj)); 3619566063dSJacob Faibussowitsch PetscCall(VecSet(uAdj, 0.0)); 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)uAdj, "adjoint")); 363*e8e188d2SZach Atkins PetscCall(DMPlexSetSNESLocalFEM(dmAdj, user, user, user)); 3649566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snesAdj)); 3659566063dSJacob Faibussowitsch PetscCall(SNESSolve(snesAdj, NULL, uAdj)); 3669566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snesAdj, &uAdj)); 3679566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uAdj, NULL, "-adjoint_view")); 368c4762a1bSJed Brown /* Error representation */ 369c4762a1bSJed Brown { 370c4762a1bSJed Brown DM dmErr, dmErrAux, dms[2]; 371c4762a1bSJed Brown Vec errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj; 372c4762a1bSJed Brown IS *subis; 373c4762a1bSJed Brown PetscReal errorEstTot, errorL2Norm, errorL2Tot; 374c4762a1bSJed Brown PetscInt N, i; 375*e8e188d2SZach Atkins PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u}; 3769371c9d4SSatish 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}; 377c4762a1bSJed Brown void *ctxs[1] = {0}; 378c4762a1bSJed Brown 379*e8e188d2SZach Atkins ctxs[0] = user; 3809566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmErr)); 381*e8e188d2SZach Atkins PetscCall(SetupDiscretization(dmErr, "error", SetupErrorProblem, user)); 3829566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErr, &errorEst)); 3839566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErr, &errorL2)); 384c4762a1bSJed Brown /* Compute auxiliary data (solution and projection of adjoint solution) */ 3859566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmAdj, &uAdjLoc)); 3869566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 3879566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc)); 3889566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &uAdjProj)); 3899566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, uAdjLoc)); 3909566063dSJacob Faibussowitsch PetscCall(DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj)); 3919566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 3929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmAdj, &uAdjLoc)); 393c4762a1bSJed Brown /* Attach auxiliary data */ 3949371c9d4SSatish Balay dms[0] = dm; 3959371c9d4SSatish Balay dms[1] = dm; 3969566063dSJacob Faibussowitsch PetscCall(DMCreateSuperDM(dms, 2, &subis, &dmErrAux)); 397c4762a1bSJed Brown if (0) { 398c4762a1bSJed Brown PetscSection sec; 399c4762a1bSJed Brown 4009566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[0], &sec)); 4019566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 4029566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dms[1], &sec)); 4039566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 4049566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmErrAux, &sec)); 4059566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD)); 406c4762a1bSJed Brown } 4079566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmErrAux, NULL, "-dm_err_view")); 4089566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(subis[0], NULL, "-super_is_view")); 4099566063dSJacob Faibussowitsch PetscCall(ISViewFromOptions(subis[1], NULL, "-super_is_view")); 4109566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmErrAux, &uErr)); 4119566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-map_vec_view")); 4129566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uAdjProj, NULL, "-map_vec_view")); 4139566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uErr, NULL, "-map_vec_view")); 4149566063dSJacob Faibussowitsch PetscCall(VecISCopy(uErr, subis[0], SCATTER_FORWARD, u)); 4159566063dSJacob Faibussowitsch PetscCall(VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj)); 4169566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &uAdjProj)); 417*e8e188d2SZach Atkins for (i = 0; i < 2; ++i) { PetscCall(ISDestroy(&subis[i])); } 4189566063dSJacob Faibussowitsch PetscCall(PetscFree(subis)); 4199566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmErrAux, &uErrLoc)); 4209566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc)); 4219566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc)); 4229566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErrAux, &uErr)); 4239566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, uErrLoc)); 424c4762a1bSJed Brown /* Compute cellwise error estimate */ 4259566063dSJacob Faibussowitsch PetscCall(VecSet(errorEst, 0.0)); 426*e8e188d2SZach Atkins PetscCall(DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, user)); 4279566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dmAdj, NULL, 0, 0, NULL)); 4289566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmErrAux, &uErrLoc)); 4299566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmErrAux)); 430c4762a1bSJed Brown /* Plot cellwise error vector */ 4319566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorEst, NULL, "-error_view")); 432c4762a1bSJed Brown /* Compute ratio of estimate (sum over cells) with actual L_2 error */ 4339566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm)); 4349566063dSJacob Faibussowitsch PetscCall(DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2)); 4359566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorL2, NULL, "-l2_error_view")); 4369566063dSJacob Faibussowitsch PetscCall(VecNorm(errorL2, NORM_INFINITY, &errorL2Tot)); 4379566063dSJacob Faibussowitsch PetscCall(VecNorm(errorEst, NORM_INFINITY, &errorEstTot)); 4389566063dSJacob Faibussowitsch PetscCall(VecGetSize(errorEst, &N)); 4399566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(errorEst, errorEst, errorL2)); 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)errorEst, "Error ratio")); 4419566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(errorEst, NULL, "-error_ratio_view")); 44263a3b9bcSJacob 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)))); 4439566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErr, &errorEst)); 4449566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmErr, &errorL2)); 4459566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmErr)); 446c4762a1bSJed Brown } 4479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAdj)); 4489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uAdj)); 4499566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snesAdj)); 450*e8e188d2SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 451c4762a1bSJed Brown } 452*e8e188d2SZach Atkins 453*e8e188d2SZach Atkins static PetscErrorCode ErrorView(Vec u, AppCtx *user) 454*e8e188d2SZach Atkins { 455*e8e188d2SZach Atkins PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 456*e8e188d2SZach Atkins void *ctx; 457*e8e188d2SZach Atkins DM dm; 458*e8e188d2SZach Atkins PetscDS ds; 459*e8e188d2SZach Atkins PetscReal error; 460*e8e188d2SZach Atkins PetscInt N; 461*e8e188d2SZach Atkins 462*e8e188d2SZach Atkins PetscFunctionBegin; 463*e8e188d2SZach Atkins if (!user->viewError) PetscFunctionReturn(PETSC_SUCCESS); 464*e8e188d2SZach Atkins PetscCall(VecGetDM(u, &dm)); 465*e8e188d2SZach Atkins PetscCall(DMGetDS(dm, &ds)); 466*e8e188d2SZach Atkins PetscCall(PetscDSGetExactSolution(ds, 0, &sol, &ctx)); 467*e8e188d2SZach Atkins PetscCall(VecGetSize(u, &N)); 468*e8e188d2SZach Atkins PetscCall(DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error)); 469*e8e188d2SZach Atkins PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " L2 error: %g\n", N, (double)error)); 470*e8e188d2SZach Atkins PetscFunctionReturn(PETSC_SUCCESS); 471*e8e188d2SZach Atkins } 472*e8e188d2SZach Atkins 473*e8e188d2SZach Atkins int main(int argc, char **argv) 474*e8e188d2SZach Atkins { 475*e8e188d2SZach Atkins DM dm; /* Problem specification */ 476*e8e188d2SZach Atkins SNES snes; /* Nonlinear solver */ 477*e8e188d2SZach Atkins Vec u; /* Solutions */ 478*e8e188d2SZach Atkins AppCtx user; /* User-defined work context */ 479*e8e188d2SZach Atkins PetscInt planeDir[2] = {0, 1}; 480*e8e188d2SZach Atkins PetscReal planeCoord[2] = {0., 1.}; 481*e8e188d2SZach Atkins 482*e8e188d2SZach Atkins PetscFunctionBeginUser; 483*e8e188d2SZach Atkins PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 484*e8e188d2SZach Atkins PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 485*e8e188d2SZach Atkins /* Primal system */ 486*e8e188d2SZach Atkins PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 487*e8e188d2SZach Atkins PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 488*e8e188d2SZach Atkins PetscCall(SNESSetDM(snes, dm)); 489*e8e188d2SZach Atkins PetscCall(SetupDiscretization(dm, "potential", SetupPrimalProblem, &user)); 490*e8e188d2SZach Atkins PetscCall(DMCreateGlobalVector(dm, &u)); 491*e8e188d2SZach Atkins PetscCall(VecSet(u, 0.0)); 492*e8e188d2SZach Atkins PetscCall(PetscObjectSetName((PetscObject)u, "potential")); 493*e8e188d2SZach Atkins PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 494*e8e188d2SZach Atkins PetscCall(SNESSetFromOptions(snes)); 495*e8e188d2SZach Atkins PetscCall(SNESSolve(snes, NULL, u)); 496*e8e188d2SZach Atkins PetscCall(SNESGetSolution(snes, &u)); 497*e8e188d2SZach Atkins PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 498*e8e188d2SZach Atkins PetscCall(ErrorView(u, &user)); 499*e8e188d2SZach Atkins PetscCall(ComputeSpectral(u, 2, planeDir, planeCoord, &user)); 500*e8e188d2SZach Atkins PetscCall(ComputeAdjoint(u, &user)); 501c4762a1bSJed Brown /* Cleanup */ 5029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5039566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 5049566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5059566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 506b122ec5aSJacob Faibussowitsch return 0; 507c4762a1bSJed Brown } 508c4762a1bSJed Brown 509c4762a1bSJed Brown /*TEST 510c4762a1bSJed Brown 511c4762a1bSJed Brown test: 5125be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5135be41b18SMatthew G. Knepley suffix: 2d_p1_conv 514c4762a1bSJed Brown requires: triangle 5155be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5165be41b18SMatthew G. Knepley test: 5175be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5185be41b18SMatthew G. Knepley suffix: 2d_p2_conv 5195be41b18SMatthew G. Knepley requires: triangle 5205be41b18SMatthew G. Knepley args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5215be41b18SMatthew G. Knepley test: 5225be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5235be41b18SMatthew G. Knepley suffix: 2d_p3_conv 5245be41b18SMatthew G. Knepley requires: triangle 5255be41b18SMatthew G. Knepley args: -potential_petscspace_degree 3 -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: 1.9 5285be41b18SMatthew G. Knepley suffix: 2d_q1_conv 52930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5305be41b18SMatthew G. Knepley test: 5315be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5325be41b18SMatthew G. Knepley suffix: 2d_q2_conv 53330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5345be41b18SMatthew G. Knepley test: 5355be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5365be41b18SMatthew G. Knepley suffix: 2d_q3_conv 53730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5385be41b18SMatthew G. Knepley test: 5395be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5405be41b18SMatthew G. Knepley suffix: 2d_q1_shear_conv 54130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5425be41b18SMatthew G. Knepley test: 5435be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5445be41b18SMatthew G. Knepley suffix: 2d_q2_shear_conv 54530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5465be41b18SMatthew G. Knepley test: 5475be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5485be41b18SMatthew G. Knepley suffix: 2d_q3_shear_conv 54930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5505be41b18SMatthew G. Knepley test: 5510fdc7489SMatthew Knepley # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7 5525be41b18SMatthew G. Knepley suffix: 3d_p1_conv 5535be41b18SMatthew G. Knepley requires: ctetgen 55430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 5555be41b18SMatthew G. Knepley test: 5565be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5575be41b18SMatthew G. Knepley suffix: 3d_p2_conv 5585be41b18SMatthew G. Knepley requires: ctetgen 55930602db0SMatthew 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 5605be41b18SMatthew G. Knepley test: 5615be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0 5625be41b18SMatthew G. Knepley suffix: 3d_p3_conv 5635be41b18SMatthew G. Knepley requires: ctetgen 56430602db0SMatthew 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 5655be41b18SMatthew G. Knepley test: 5665be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8 5675be41b18SMatthew G. Knepley suffix: 3d_q1_conv 56830602db0SMatthew 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 5695be41b18SMatthew G. Knepley test: 5705be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5715be41b18SMatthew G. Knepley suffix: 3d_q2_conv 57230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 5735be41b18SMatthew G. Knepley test: 5745be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8 5755be41b18SMatthew G. Knepley suffix: 3d_q3_conv 57630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 5779139b27eSToby Isaac test: 5789139b27eSToby Isaac suffix: 2d_p1_fas_full 5799139b27eSToby Isaac requires: triangle 580e600fa54SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 5819139b27eSToby Isaac -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full -snes_fas_full_total \ 5829139b27eSToby Isaac -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 5839139b27eSToby Isaac -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 5849139b27eSToby Isaac -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 58573f7197eSJed 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 5869139b27eSToby Isaac test: 5879139b27eSToby Isaac suffix: 2d_p1_fas_full_homogeneous 5889139b27eSToby Isaac requires: triangle 589e600fa54SMatthew G. Knepley args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 5909139b27eSToby Isaac -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \ 5919139b27eSToby Isaac -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 5929139b27eSToby Isaac -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 5939139b27eSToby Isaac -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 59473f7197eSJed 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 5955be41b18SMatthew G. Knepley 596c4762a1bSJed Brown test: 597c4762a1bSJed Brown suffix: 2d_p1_scalable 5985be41b18SMatthew G. Knepley requires: triangle 5995be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 3 \ 600c4762a1bSJed Brown -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \ 60173f7197eSJed Brown -pc_type gamg -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 \ 602c4762a1bSJed Brown -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 603c4762a1bSJed Brown -pc_gamg_coarse_eq_limit 1000 \ 604c4762a1bSJed Brown -pc_gamg_threshold 0.05 \ 605c4762a1bSJed Brown -pc_gamg_threshold_scale .0 \ 606c4762a1bSJed Brown -mg_levels_ksp_type chebyshev \ 607c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 608c4762a1bSJed Brown -mg_levels_pc_type jacobi \ 609c4762a1bSJed Brown -matptap_via scalable 610c4762a1bSJed Brown test: 611c4762a1bSJed Brown suffix: 2d_p1_gmg_vcycle 612c4762a1bSJed Brown requires: triangle 6135be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 614c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg \ 615c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 616c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 617c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 61873f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 619c4762a1bSJed Brown -mg_levels_pc_type jacobi 620c4762a1bSJed Brown test: 621c4762a1bSJed Brown suffix: 2d_p1_gmg_fcycle 622c4762a1bSJed Brown requires: triangle 6235be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 624c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \ 625c4762a1bSJed Brown -mg_levels_ksp_max_it 2 \ 626c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 627c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 62873f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 629c4762a1bSJed Brown -mg_levels_pc_type jacobi 630c4762a1bSJed Brown test: 631d6837840SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_adapt 6322b3cbbdaSStefano Zampini requires: triangle 633908b9b43SStefano Zampini args: -petscpartitioner_type simple -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 6342b3cbbdaSStefano Zampini -ksp_rtol 5e-10 -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space harmonic -pc_mg_adapt_interp_n 8 \ 635d6837840SMatthew G. Knepley -mg_levels_ksp_max_it 1 \ 636d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 637d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 63873f7197eSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 639d6837840SMatthew G. Knepley -mg_levels_pc_type jacobi 640d6837840SMatthew G. Knepley test: 641c4762a1bSJed Brown suffix: 2d_p1_spectral_0 642c4762a1bSJed Brown requires: triangle fftw !complex 6435be41b18SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view 644c4762a1bSJed Brown test: 645c4762a1bSJed Brown suffix: 2d_p1_spectral_1 646c4762a1bSJed Brown requires: triangle fftw !complex 647c4762a1bSJed Brown nsize: 2 648e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view 649c4762a1bSJed Brown test: 650c4762a1bSJed Brown suffix: 2d_p1_adj_0 651c4762a1bSJed Brown requires: triangle 6525be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0 653b1ad28e3SJunchao Zhang test: 654b1ad28e3SJunchao Zhang nsize: 2 655dcfd994dSJunchao Zhang requires: kokkos_kernels 656b1ad28e3SJunchao Zhang suffix: kokkos 657e600fa54SMatthew 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 \ 658b1ad28e3SJunchao 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 \ 65973f7197eSJed 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 \ 66073f7197eSJed Brown -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 661c4762a1bSJed Brown 662c4762a1bSJed Brown TEST*/ 663