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 */ 17c4762a1bSJed Brown } AppCtx; 18c4762a1bSJed Brown 19c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 20c4762a1bSJed Brown { 21c4762a1bSJed Brown *u = 0.0; 22c4762a1bSJed Brown return 0; 23c4762a1bSJed Brown } 24c4762a1bSJed Brown 25c4762a1bSJed Brown static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 26c4762a1bSJed Brown { 27c4762a1bSJed Brown PetscInt d; 28c4762a1bSJed Brown *u = 0.0; 29c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 30c4762a1bSJed Brown return 0; 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Compute integral of (residual of solution)*(adjoint solution - projection of adjoint solution) */ 34c4762a1bSJed Brown static void obj_error_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 35c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 36c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 37c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[]) 38c4762a1bSJed Brown { 39c4762a1bSJed Brown obj[0] = a[aOff[0]]*(u[0] - a[aOff[1]]); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 43c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 44c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 45c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 46c4762a1bSJed Brown { 47c4762a1bSJed Brown PetscInt d; 48c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown static void f0_unity_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 52c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 53c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 54c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 55c4762a1bSJed Brown { 56c4762a1bSJed Brown f0[0] = 1.0; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown static void f0_identityaux_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 60c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 61c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 62c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 63c4762a1bSJed Brown { 64c4762a1bSJed Brown f0[0] = a[0]; 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 68c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 69c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 70c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 71c4762a1bSJed Brown { 72c4762a1bSJed Brown PetscInt d; 73c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 77c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 78c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 79c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscInt d; 82c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown PetscErrorCode ierr; 88c4762a1bSJed Brown 89c4762a1bSJed Brown PetscFunctionBeginUser; 90c4762a1bSJed Brown options->shear = PETSC_FALSE; 91c4762a1bSJed Brown options->spectral = PETSC_FALSE; 92c4762a1bSJed Brown options->adjoint = PETSC_FALSE; 93c4762a1bSJed Brown 94c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = PetscOptionsEnd(); 99c4762a1bSJed Brown PetscFunctionReturn(0); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown static PetscErrorCode CreateSpectralPlanes(DM dm, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown PetscSection coordSection; 105c4762a1bSJed Brown Vec coordinates; 106c4762a1bSJed Brown const PetscScalar *coords; 107c4762a1bSJed Brown PetscInt dim, p, vStart, vEnd, v; 108c4762a1bSJed Brown PetscErrorCode ierr; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBeginUser; 111c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 116c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 117c4762a1bSJed Brown DMLabel label; 118c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 119c4762a1bSJed Brown 120c4762a1bSJed Brown ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = DMCreateLabel(dm, name);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = DMLabelAddStratum(label, 1);CHKERRQ(ierr); 124c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 125c4762a1bSJed Brown PetscInt off; 126c4762a1bSJed Brown 127c4762a1bSJed Brown ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 128c4762a1bSJed Brown if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off+planeDir[p]])) < PETSC_SMALL) { 129c4762a1bSJed Brown ierr = DMLabelSetValue(label, v, 1);CHKERRQ(ierr); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown } 132c4762a1bSJed Brown } 133c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 134c4762a1bSJed Brown PetscFunctionReturn(0); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 138c4762a1bSJed Brown { 139c4762a1bSJed Brown PetscErrorCode ierr; 140c4762a1bSJed Brown 141c4762a1bSJed Brown PetscFunctionBeginUser; 142c4762a1bSJed Brown /* Create box mesh */ 143*5be41b18SMatthew G. Knepley ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 144c4762a1bSJed Brown /* TODO: This should be pulled into the library */ 145c4762a1bSJed Brown { 146c4762a1bSJed Brown char convType[256]; 147c4762a1bSJed Brown PetscBool flg; 148c4762a1bSJed Brown 149c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = PetscOptionsEnd(); 152c4762a1bSJed Brown if (flg) { 153c4762a1bSJed Brown DM dmConv; 154c4762a1bSJed Brown 155c4762a1bSJed Brown ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr); 156c4762a1bSJed Brown if (dmConv) { 157c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 158c4762a1bSJed Brown *dm = dmConv; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown if (user->shear) {ierr = DMPlexShearGeometry(*dm, DM_X, NULL);CHKERRQ(ierr);} 163c4762a1bSJed Brown ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 164c4762a1bSJed Brown 165c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 169c4762a1bSJed Brown if (user->spectral) { 170c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 171c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 172c4762a1bSJed Brown 173c4762a1bSJed Brown ierr = CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user);CHKERRQ(ierr); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown PetscFunctionReturn(0); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown PetscDS prob; 181c4762a1bSJed Brown const PetscInt id = 1; 182c4762a1bSJed Brown PetscErrorCode ierr; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBeginUser; 185c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_trig_u, f1_u);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 0, trig_u, user);CHKERRQ(ierr); 18956cf3b9cSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) trig_u, NULL, 1, &id, user);CHKERRQ(ierr); 190c4762a1bSJed Brown PetscFunctionReturn(0); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user) 194c4762a1bSJed Brown { 195c4762a1bSJed Brown PetscDS prob; 196c4762a1bSJed Brown const PetscInt id = 1; 197c4762a1bSJed Brown PetscErrorCode ierr; 198c4762a1bSJed Brown 199c4762a1bSJed Brown PetscFunctionBeginUser; 200c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_unity_u, f1_u);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = PetscDSSetObjective(prob, 0, obj_error_u);CHKERRQ(ierr); 20456cf3b9cSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) zero, NULL, 1, &id, user);CHKERRQ(ierr); 205c4762a1bSJed Brown PetscFunctionReturn(0); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user) 209c4762a1bSJed Brown { 210c4762a1bSJed Brown PetscDS prob; 211c4762a1bSJed Brown PetscErrorCode ierr; 212c4762a1bSJed Brown 213c4762a1bSJed Brown PetscFunctionBeginUser; 214c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 215c4762a1bSJed Brown PetscFunctionReturn(0); 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 219c4762a1bSJed Brown { 220c4762a1bSJed Brown DM cdm = dm; 221c4762a1bSJed Brown PetscFE fe; 222*5be41b18SMatthew G. Knepley DMPolytopeType ct; 223*5be41b18SMatthew G. Knepley PetscBool simplex; 224*5be41b18SMatthew G. Knepley PetscInt dim, cStart; 225c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 226c4762a1bSJed Brown PetscErrorCode ierr; 227c4762a1bSJed Brown 228c4762a1bSJed Brown PetscFunctionBeginUser; 229*5be41b18SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 230*5be41b18SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 231*5be41b18SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 232*5be41b18SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 233c4762a1bSJed Brown /* Create finite element */ 234c4762a1bSJed Brown ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 235*5be41b18SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 237c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 238c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = (*setup)(dm, user);CHKERRQ(ierr); 241c4762a1bSJed Brown while (cdm) { 242c4762a1bSJed Brown ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 243c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 244c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 247c4762a1bSJed Brown PetscFunctionReturn(0); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 251c4762a1bSJed Brown { 252c4762a1bSJed Brown MPI_Comm comm; 253c4762a1bSJed Brown PetscSection coordSection, section; 254c4762a1bSJed Brown Vec coordinates, uloc; 255c4762a1bSJed Brown const PetscScalar *coords, *array; 256c4762a1bSJed Brown PetscInt p; 257c4762a1bSJed Brown PetscMPIInt size, rank; 258c4762a1bSJed Brown PetscErrorCode ierr; 259c4762a1bSJed Brown 260c4762a1bSJed Brown PetscFunctionBeginUser; 261c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 264c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &uloc);CHKERRQ(ierr); 265c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 268c4762a1bSJed Brown ierr = VecViewFromOptions(uloc, NULL, "-sol_view");CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = VecGetArrayRead(uloc, &array);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 274c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 275c4762a1bSJed Brown DMLabel label; 276c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 277c4762a1bSJed Brown Mat F; 278c4762a1bSJed Brown Vec x, y; 279c4762a1bSJed Brown IS stratum; 280c4762a1bSJed Brown PetscReal *ray, *gray; 281c4762a1bSJed Brown PetscScalar *rvals, *svals, *gsvals; 282c4762a1bSJed Brown PetscInt *perm, *nperm; 283c4762a1bSJed Brown PetscInt n, N, i, j, off, offu; 284c4762a1bSJed Brown const PetscInt *points; 285c4762a1bSJed Brown 286c4762a1bSJed Brown ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr); 287c4762a1bSJed Brown ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = DMLabelGetStratumIS(label, 1, &stratum);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = ISGetLocalSize(stratum, &n);CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = ISGetIndices(stratum, &points);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = PetscMalloc2(n, &ray, n, &svals);CHKERRQ(ierr); 292c4762a1bSJed Brown for (i = 0; i < n; ++i) { 293c4762a1bSJed Brown ierr = PetscSectionGetOffset(coordSection, points[i], &off);CHKERRQ(ierr); 294c4762a1bSJed Brown ierr = PetscSectionGetOffset(section, points[i], &offu);CHKERRQ(ierr); 295c4762a1bSJed Brown ray[i] = PetscRealPart(coords[off+((planeDir[p]+1)%2)]); 296c4762a1bSJed Brown svals[i] = array[offu]; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown /* Gather the ray data to proc 0 */ 299c4762a1bSJed Brown if (size > 1) { 300c4762a1bSJed Brown PetscMPIInt *cnt, *displs, p; 301c4762a1bSJed Brown 302c4762a1bSJed Brown ierr = PetscCalloc2(size, &cnt, size, &displs);CHKERRQ(ierr); 303c4762a1bSJed Brown ierr = MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 304c4762a1bSJed Brown for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1]; 305c4762a1bSJed Brown N = displs[size-1] + cnt[size-1]; 306c4762a1bSJed Brown ierr = PetscMalloc2(N, &gray, N, &gsvals);CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm);CHKERRQ(ierr); 308c4762a1bSJed Brown ierr = MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm);CHKERRQ(ierr); 309c4762a1bSJed Brown ierr = PetscFree2(cnt, displs);CHKERRQ(ierr); 310c4762a1bSJed Brown } else { 311c4762a1bSJed Brown N = n; 312c4762a1bSJed Brown gray = ray; 313c4762a1bSJed Brown gsvals = svals; 314c4762a1bSJed Brown } 315c4762a1bSJed Brown if (!rank) { 316c4762a1bSJed Brown /* Sort point along ray */ 317c4762a1bSJed Brown ierr = PetscMalloc2(N, &perm, N, &nperm);CHKERRQ(ierr); 318c4762a1bSJed Brown for (i = 0; i < N; ++i) {perm[i] = i;} 319c4762a1bSJed Brown ierr = PetscSortRealWithPermutation(N, gray, perm);CHKERRQ(ierr); 320c4762a1bSJed Brown /* Count duplicates and squish mapping */ 321c4762a1bSJed Brown nperm[0] = perm[0]; 322c4762a1bSJed Brown for (i = 1, j = 1; i < N; ++i) { 323c4762a1bSJed Brown if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i]; 324c4762a1bSJed Brown } 325c4762a1bSJed Brown /* Create FFT structs */ 326c4762a1bSJed Brown ierr = MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F);CHKERRQ(ierr); 327c4762a1bSJed Brown ierr = MatCreateVecs(F, &x, &y);CHKERRQ(ierr); 328c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y, name);CHKERRQ(ierr); 329c4762a1bSJed Brown ierr = VecGetArray(x, &rvals);CHKERRQ(ierr); 330c4762a1bSJed Brown for (i = 0, j = 0; i < N; ++i) { 331c4762a1bSJed Brown if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue; 332c4762a1bSJed Brown rvals[j] = gsvals[nperm[j]]; 333c4762a1bSJed Brown ++j; 334c4762a1bSJed Brown } 335c4762a1bSJed Brown ierr = PetscFree2(perm, nperm);CHKERRQ(ierr); 336c4762a1bSJed Brown if (size > 1) {ierr = PetscFree2(gray, gsvals);CHKERRQ(ierr);} 337c4762a1bSJed Brown ierr = VecRestoreArray(x, &rvals);CHKERRQ(ierr); 338c4762a1bSJed Brown /* Do FFT along the ray */ 339c4762a1bSJed Brown ierr = MatMult(F, x, y);CHKERRQ(ierr); 340c4762a1bSJed Brown /* Chop FFT */ 341c4762a1bSJed Brown ierr = VecChop(y, PETSC_SMALL);CHKERRQ(ierr); 342c4762a1bSJed Brown ierr = VecViewFromOptions(x, NULL, "-real_view");CHKERRQ(ierr); 343c4762a1bSJed Brown ierr = VecViewFromOptions(y, NULL, "-fft_view");CHKERRQ(ierr); 344c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 345c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 347c4762a1bSJed Brown } 348c4762a1bSJed Brown ierr = ISRestoreIndices(stratum, &points);CHKERRQ(ierr); 349c4762a1bSJed Brown ierr = ISDestroy(&stratum);CHKERRQ(ierr); 350c4762a1bSJed Brown ierr = PetscFree2(ray, svals);CHKERRQ(ierr); 351c4762a1bSJed Brown } 352c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 353c4762a1bSJed Brown ierr = VecRestoreArrayRead(uloc, &array);CHKERRQ(ierr); 354c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &uloc);CHKERRQ(ierr); 355c4762a1bSJed Brown PetscFunctionReturn(0); 356c4762a1bSJed Brown } 357c4762a1bSJed Brown 358c4762a1bSJed Brown int main(int argc, char **argv) 359c4762a1bSJed Brown { 360c4762a1bSJed Brown DM dm; /* Problem specification */ 361c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 362c4762a1bSJed Brown Vec u; /* Solutions */ 363c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 364c4762a1bSJed Brown PetscErrorCode ierr; 365c4762a1bSJed Brown 366c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 367c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 368c4762a1bSJed Brown /* Primal system */ 369c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = VecSet(u, 0.0);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 377c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr); 381c4762a1bSJed Brown if (user.spectral) { 382c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 383c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 384c4762a1bSJed Brown 385c4762a1bSJed Brown ierr = ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user);CHKERRQ(ierr); 386c4762a1bSJed Brown } 387c4762a1bSJed Brown /* Adjoint system */ 388c4762a1bSJed Brown if (user.adjoint) { 389c4762a1bSJed Brown DM dmAdj; 390c4762a1bSJed Brown SNES snesAdj; 391c4762a1bSJed Brown Vec uAdj; 392c4762a1bSJed Brown 393c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snesAdj);CHKERRQ(ierr); 394c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_");CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = DMClone(dm, &dmAdj);CHKERRQ(ierr); 396c4762a1bSJed Brown ierr = SNESSetDM(snesAdj, dmAdj);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user);CHKERRQ(ierr); 398c4762a1bSJed Brown ierr = DMCreateGlobalVector(dmAdj, &uAdj);CHKERRQ(ierr); 399c4762a1bSJed Brown ierr = VecSet(uAdj, 0.0);CHKERRQ(ierr); 400c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) uAdj, "adjoint");CHKERRQ(ierr); 401c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user);CHKERRQ(ierr); 402c4762a1bSJed Brown ierr = SNESSetFromOptions(snesAdj);CHKERRQ(ierr); 403c4762a1bSJed Brown ierr = SNESSolve(snesAdj, NULL, uAdj);CHKERRQ(ierr); 404c4762a1bSJed Brown ierr = SNESGetSolution(snesAdj, &uAdj);CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = VecViewFromOptions(uAdj, NULL, "-adjoint_view");CHKERRQ(ierr); 406c4762a1bSJed Brown /* Error representation */ 407c4762a1bSJed Brown { 408c4762a1bSJed Brown DM dmErr, dmErrAux, dms[2]; 409c4762a1bSJed Brown Vec errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj; 410c4762a1bSJed Brown IS *subis; 411c4762a1bSJed Brown PetscReal errorEstTot, errorL2Norm, errorL2Tot; 412c4762a1bSJed Brown PetscInt N, i; 413c4762a1bSJed Brown PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {trig_u}; 414c4762a1bSJed Brown void (*identity[1])(PetscInt, PetscInt, PetscInt, 415c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 416c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 417c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u}; 418c4762a1bSJed Brown void *ctxs[1] = {0}; 419c4762a1bSJed Brown 420c4762a1bSJed Brown ctxs[0] = &user; 421c4762a1bSJed Brown ierr = DMClone(dm, &dmErr);CHKERRQ(ierr); 422c4762a1bSJed Brown ierr = SetupDiscretization(dmErr, "error", SetupErrorProblem, &user);CHKERRQ(ierr); 423c4762a1bSJed Brown ierr = DMGetGlobalVector(dmErr, &errorEst);CHKERRQ(ierr); 424c4762a1bSJed Brown ierr = DMGetGlobalVector(dmErr, &errorL2);CHKERRQ(ierr); 425c4762a1bSJed Brown /* Compute auxiliary data (solution and projection of adjoint solution) */ 426c4762a1bSJed Brown ierr = DMGetLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr); 427c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr); 428c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr); 429c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &uAdjProj);CHKERRQ(ierr); 430c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAdj);CHKERRQ(ierr); 431c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) uAdjLoc);CHKERRQ(ierr); 432c4762a1bSJed Brown ierr = DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj);CHKERRQ(ierr); 433c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", NULL);CHKERRQ(ierr); 434c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", NULL);CHKERRQ(ierr); 435c4762a1bSJed Brown ierr = DMRestoreLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr); 436c4762a1bSJed Brown /* Attach auxiliary data */ 437c4762a1bSJed Brown dms[0] = dm; dms[1] = dm; 438c4762a1bSJed Brown ierr = DMCreateSuperDM(dms, 2, &subis, &dmErrAux);CHKERRQ(ierr); 439c4762a1bSJed Brown if (0) { 440c4762a1bSJed Brown PetscSection sec; 441c4762a1bSJed Brown 442c4762a1bSJed Brown ierr = DMGetLocalSection(dms[0], &sec);CHKERRQ(ierr); 443c4762a1bSJed Brown ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 444c4762a1bSJed Brown ierr = DMGetLocalSection(dms[1], &sec);CHKERRQ(ierr); 445c4762a1bSJed Brown ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = DMGetLocalSection(dmErrAux, &sec);CHKERRQ(ierr); 447c4762a1bSJed Brown ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 448c4762a1bSJed Brown } 449c4762a1bSJed Brown ierr = DMViewFromOptions(dmErrAux, NULL, "-dm_err_view");CHKERRQ(ierr); 450c4762a1bSJed Brown ierr = ISViewFromOptions(subis[0], NULL, "-super_is_view");CHKERRQ(ierr); 451c4762a1bSJed Brown ierr = ISViewFromOptions(subis[1], NULL, "-super_is_view");CHKERRQ(ierr); 452c4762a1bSJed Brown ierr = DMGetGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr); 453c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-map_vec_view");CHKERRQ(ierr); 454c4762a1bSJed Brown ierr = VecViewFromOptions(uAdjProj, NULL, "-map_vec_view");CHKERRQ(ierr); 455c4762a1bSJed Brown ierr = VecViewFromOptions(uErr, NULL, "-map_vec_view");CHKERRQ(ierr); 456c4762a1bSJed Brown ierr = VecISCopy(uErr, subis[0], SCATTER_FORWARD, u);CHKERRQ(ierr); 457c4762a1bSJed Brown ierr = VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &uAdjProj);CHKERRQ(ierr); 459c4762a1bSJed Brown for (i = 0; i < 2; ++i) {ierr = ISDestroy(&subis[i]);CHKERRQ(ierr);} 460c4762a1bSJed Brown ierr = PetscFree(subis);CHKERRQ(ierr); 461c4762a1bSJed Brown ierr = DMGetLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr); 462c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr); 463c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr); 464c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr); 465c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dmAdj, "dmAux", (PetscObject) dmErrAux);CHKERRQ(ierr); 466c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dmAdj, "A", (PetscObject) uErrLoc);CHKERRQ(ierr); 467c4762a1bSJed Brown /* Compute cellwise error estimate */ 468c4762a1bSJed Brown ierr = VecSet(errorEst, 0.0);CHKERRQ(ierr); 469c4762a1bSJed Brown ierr = DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user);CHKERRQ(ierr); 470c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dmAdj, "dmAux", NULL);CHKERRQ(ierr); 471c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dmAdj, "A", NULL);CHKERRQ(ierr); 472c4762a1bSJed Brown ierr = DMRestoreLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr); 473c4762a1bSJed Brown ierr = DMDestroy(&dmErrAux);CHKERRQ(ierr); 474c4762a1bSJed Brown /* Plot cellwise error vector */ 475c4762a1bSJed Brown ierr = VecViewFromOptions(errorEst, NULL, "-error_view");CHKERRQ(ierr); 476c4762a1bSJed Brown /* Compute ratio of estimate (sum over cells) with actual L_2 error */ 477c4762a1bSJed Brown ierr = DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm);CHKERRQ(ierr); 478c4762a1bSJed Brown ierr = DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2);CHKERRQ(ierr); 479c4762a1bSJed Brown ierr = VecViewFromOptions(errorL2, NULL, "-l2_error_view");CHKERRQ(ierr); 480c4762a1bSJed Brown ierr = VecNorm(errorL2, NORM_INFINITY, &errorL2Tot);CHKERRQ(ierr); 481c4762a1bSJed Brown ierr = VecNorm(errorEst, NORM_INFINITY, &errorEstTot);CHKERRQ(ierr); 482c4762a1bSJed Brown ierr = VecGetSize(errorEst, &N);CHKERRQ(ierr); 483c4762a1bSJed Brown ierr = VecPointwiseDivide(errorEst, errorEst, errorL2);CHKERRQ(ierr); 484c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) errorEst, "Error ratio");CHKERRQ(ierr); 485c4762a1bSJed Brown ierr = VecViewFromOptions(errorEst, NULL, "-error_ratio_view");CHKERRQ(ierr); 486c4762a1bSJed Brown ierr = 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));CHKERRQ(ierr); 487c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmErr, &errorEst);CHKERRQ(ierr); 488c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmErr, &errorL2);CHKERRQ(ierr); 489c4762a1bSJed Brown ierr = DMDestroy(&dmErr);CHKERRQ(ierr); 490c4762a1bSJed Brown } 491c4762a1bSJed Brown ierr = DMDestroy(&dmAdj);CHKERRQ(ierr); 492c4762a1bSJed Brown ierr = VecDestroy(&uAdj);CHKERRQ(ierr); 493c4762a1bSJed Brown ierr = SNESDestroy(&snesAdj);CHKERRQ(ierr); 494c4762a1bSJed Brown } 495c4762a1bSJed Brown /* Cleanup */ 496c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 498c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 499c4762a1bSJed Brown ierr = PetscFinalize(); 500c4762a1bSJed Brown return ierr; 501c4762a1bSJed Brown } 502c4762a1bSJed Brown 503c4762a1bSJed Brown /*TEST 504c4762a1bSJed Brown 505c4762a1bSJed Brown test: 506*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 507*5be41b18SMatthew G. Knepley suffix: 2d_p1_conv 508c4762a1bSJed Brown requires: triangle 509*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 510*5be41b18SMatthew G. Knepley test: 511*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 512*5be41b18SMatthew G. Knepley suffix: 2d_p2_conv 513*5be41b18SMatthew G. Knepley requires: triangle 514*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 515*5be41b18SMatthew G. Knepley test: 516*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 517*5be41b18SMatthew G. Knepley suffix: 2d_p3_conv 518*5be41b18SMatthew G. Knepley requires: triangle 519*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 520*5be41b18SMatthew G. Knepley test: 521*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 522*5be41b18SMatthew G. Knepley suffix: 2d_q1_conv 523*5be41b18SMatthew G. Knepley args: -dm_plex_box_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 524*5be41b18SMatthew G. Knepley test: 525*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 526*5be41b18SMatthew G. Knepley suffix: 2d_q2_conv 527*5be41b18SMatthew G. Knepley args: -dm_plex_box_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 528*5be41b18SMatthew G. Knepley test: 529*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 530*5be41b18SMatthew G. Knepley suffix: 2d_q3_conv 531*5be41b18SMatthew G. Knepley args: -dm_plex_box_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 532*5be41b18SMatthew G. Knepley test: 533*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 534*5be41b18SMatthew G. Knepley suffix: 2d_q1_shear_conv 535*5be41b18SMatthew G. Knepley args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 536*5be41b18SMatthew G. Knepley test: 537*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 538*5be41b18SMatthew G. Knepley suffix: 2d_q2_shear_conv 539*5be41b18SMatthew G. Knepley args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 540*5be41b18SMatthew G. Knepley test: 541*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 542*5be41b18SMatthew G. Knepley suffix: 2d_q3_shear_conv 543*5be41b18SMatthew G. Knepley args: -dm_plex_box_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 544*5be41b18SMatthew G. Knepley test: 545*5be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 1.7 546*5be41b18SMatthew G. Knepley suffix: 3d_p1_conv 547*5be41b18SMatthew G. Knepley requires: ctetgen 548*5be41b18SMatthew G. Knepley args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 549*5be41b18SMatthew G. Knepley test: 550*5be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8 551*5be41b18SMatthew G. Knepley suffix: 3d_p2_conv 552*5be41b18SMatthew G. Knepley requires: ctetgen 553*5be41b18SMatthew G. Knepley args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 554*5be41b18SMatthew G. Knepley test: 555*5be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0 556*5be41b18SMatthew G. Knepley suffix: 3d_p3_conv 557*5be41b18SMatthew G. Knepley requires: ctetgen 558*5be41b18SMatthew G. Knepley args: -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 559*5be41b18SMatthew G. Knepley test: 560*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8 561*5be41b18SMatthew G. Knepley suffix: 3d_q1_conv 562*5be41b18SMatthew G. Knepley args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 563*5be41b18SMatthew G. Knepley test: 564*5be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8 565*5be41b18SMatthew G. Knepley suffix: 3d_q2_conv 566*5be41b18SMatthew G. Knepley args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 567*5be41b18SMatthew G. Knepley test: 568*5be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8 569*5be41b18SMatthew G. Knepley suffix: 3d_q3_conv 570*5be41b18SMatthew G. Knepley args: -dm_plex_box_dim 3 -dm_plex_box_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 571*5be41b18SMatthew G. Knepley 572c4762a1bSJed Brown test: 573c4762a1bSJed Brown suffix: 2d_p1_scalable 574*5be41b18SMatthew G. Knepley requires: triangle 575*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 3 \ 576c4762a1bSJed Brown -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \ 577c4762a1bSJed Brown -pc_type gamg \ 578c4762a1bSJed Brown -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 579c4762a1bSJed Brown -pc_gamg_coarse_eq_limit 1000 \ 580c4762a1bSJed Brown -pc_gamg_square_graph 1 \ 581c4762a1bSJed Brown -pc_gamg_threshold 0.05 \ 582c4762a1bSJed Brown -pc_gamg_threshold_scale .0 \ 583c4762a1bSJed Brown -mg_levels_ksp_type chebyshev \ 584c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 585c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 586c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 587c4762a1bSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 588c4762a1bSJed Brown -mg_levels_pc_type jacobi \ 589c4762a1bSJed Brown -matptap_via scalable 590c4762a1bSJed Brown test: 591c4762a1bSJed Brown suffix: 2d_p1_gmg_vcycle 592c4762a1bSJed Brown requires: triangle 593*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 594c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg \ 595c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 596c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 597c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 598c4762a1bSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 599c4762a1bSJed Brown -mg_levels_pc_type jacobi 600c4762a1bSJed Brown test: 601c4762a1bSJed Brown suffix: 2d_p1_gmg_fcycle 602c4762a1bSJed Brown requires: triangle 603*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 604c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \ 605c4762a1bSJed Brown -mg_levels_ksp_max_it 2 \ 606c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 607c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 608c4762a1bSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 609c4762a1bSJed Brown -mg_levels_pc_type jacobi 610c4762a1bSJed Brown test: 611d6837840SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_adapt 612d6837840SMatthew G. Knepley requires: triangle bamg 613*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 614d6837840SMatthew 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 \ 615d6837840SMatthew G. Knepley -mg_levels_ksp_max_it 1 \ 616d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 617d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 618d6837840SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 619d6837840SMatthew G. Knepley -mg_levels_pc_type jacobi 620d6837840SMatthew G. Knepley test: 621c4762a1bSJed Brown suffix: 2d_p1_spectral_0 622c4762a1bSJed Brown requires: triangle fftw !complex 623*5be41b18SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view 624c4762a1bSJed Brown test: 625c4762a1bSJed Brown suffix: 2d_p1_spectral_1 626c4762a1bSJed Brown requires: triangle fftw !complex 627c4762a1bSJed Brown nsize: 2 628*5be41b18SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_distribute -potential_petscspace_degree 1 -spectral -fft_view 629c4762a1bSJed Brown test: 630c4762a1bSJed Brown suffix: 2d_p1_adj_0 631c4762a1bSJed Brown requires: triangle 632*5be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0 633c4762a1bSJed Brown 634c4762a1bSJed Brown TEST*/ 635