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