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 125c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = PetscOptionsBool("-shear", "Shear the domain", "ex13.c", options->shear, &options->shear, NULL);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = PetscOptionsBool("-spectral", "Look at the spectrum along planes of the solution", "ex13.c", options->spectral, &options->spectral, NULL);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = PetscOptionsBool("-adjoint", "Solve the adjoint problem", "ex13.c", options->adjoint, &options->adjoint, NULL);CHKERRQ(ierr); 129b724ec41SToby Isaac ierr = PetscOptionsBool("-homogeneous", "Use homogeneous boundary conditions", "ex13.c", options->homogeneous, &options->homogeneous, NULL);CHKERRQ(ierr); 130ad1a7c93SMatthew Knepley ierr = PetscOptionsBool("-error_view", "Output the solution error", "ex13.c", options->viewError, &options->viewError, NULL);CHKERRQ(ierr); 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 PetscErrorCode ierr; 142c4762a1bSJed Brown 143c4762a1bSJed Brown PetscFunctionBeginUser; 144c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm, &dim);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 149c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 150c4762a1bSJed Brown DMLabel label; 151c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 152c4762a1bSJed Brown 153c4762a1bSJed Brown ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = DMCreateLabel(dm, name);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = DMLabelAddStratum(label, 1);CHKERRQ(ierr); 157c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 158c4762a1bSJed Brown PetscInt off; 159c4762a1bSJed Brown 160c4762a1bSJed Brown ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 161c4762a1bSJed Brown if (PetscAbsReal(planeCoord[p] - PetscRealPart(coords[off+planeDir[p]])) < PETSC_SMALL) { 162c4762a1bSJed Brown ierr = DMLabelSetValue(label, v, 1);CHKERRQ(ierr); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 167c4762a1bSJed Brown PetscFunctionReturn(0); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 171c4762a1bSJed Brown { 172c4762a1bSJed Brown PetscErrorCode ierr; 173c4762a1bSJed Brown 174c4762a1bSJed Brown PetscFunctionBeginUser; 17530602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 17630602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 17830602db0SMatthew G. Knepley if (user->shear) {ierr = DMPlexShearGeometry(*dm, DM_X, NULL);CHKERRQ(ierr);} 17930602db0SMatthew G. Knepley ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 181c4762a1bSJed Brown if (user->spectral) { 182c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 183c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 184c4762a1bSJed Brown 185c4762a1bSJed Brown ierr = CreateSpectralPlanes(*dm, 2, planeDir, planeCoord, user);CHKERRQ(ierr); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown PetscFunctionReturn(0); 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 191c4762a1bSJed Brown { 19245480ffeSMatthew G. Knepley PetscDS ds; 19345480ffeSMatthew G. Knepley DMLabel label; 194c4762a1bSJed Brown const PetscInt id = 1; 195b724ec41SToby Isaac PetscPointFunc f0 = user->homogeneous ? f0_trig_homogeneous_u : f0_trig_inhomogeneous_u; 196b724ec41SToby Isaac PetscErrorCode (*ex)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = user->homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u; 197c4762a1bSJed Brown PetscErrorCode ierr; 198c4762a1bSJed Brown 199c4762a1bSJed Brown PetscFunctionBeginUser; 20045480ffeSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 20145480ffeSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0, f1_u);CHKERRQ(ierr); 20245480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 20345480ffeSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, ex, user);CHKERRQ(ierr); 20445480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 20545480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ex, NULL, user, NULL);CHKERRQ(ierr); 206c4762a1bSJed Brown PetscFunctionReturn(0); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209c4762a1bSJed Brown static PetscErrorCode SetupAdjointProblem(DM dm, AppCtx *user) 210c4762a1bSJed Brown { 21145480ffeSMatthew G. Knepley PetscDS ds; 21245480ffeSMatthew G. Knepley DMLabel label; 213c4762a1bSJed Brown const PetscInt id = 1; 214c4762a1bSJed Brown PetscErrorCode ierr; 215c4762a1bSJed Brown 216c4762a1bSJed Brown PetscFunctionBeginUser; 21745480ffeSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 21845480ffeSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_unity_u, f1_u);CHKERRQ(ierr); 21945480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 22045480ffeSMatthew G. Knepley ierr = PetscDSSetObjective(ds, 0, obj_error_u);CHKERRQ(ierr); 22145480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 22245480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr); 223c4762a1bSJed Brown PetscFunctionReturn(0); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown static PetscErrorCode SetupErrorProblem(DM dm, AppCtx *user) 227c4762a1bSJed Brown { 228c4762a1bSJed Brown PetscDS prob; 229c4762a1bSJed Brown PetscErrorCode ierr; 230c4762a1bSJed Brown 231c4762a1bSJed Brown PetscFunctionBeginUser; 232c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 233c4762a1bSJed Brown PetscFunctionReturn(0); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 237c4762a1bSJed Brown { 238c4762a1bSJed Brown DM cdm = dm; 239c4762a1bSJed Brown PetscFE fe; 2405be41b18SMatthew G. Knepley DMPolytopeType ct; 2415be41b18SMatthew G. Knepley PetscBool simplex; 2425be41b18SMatthew G. Knepley PetscInt dim, cStart; 243c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 244c4762a1bSJed Brown PetscErrorCode ierr; 245c4762a1bSJed Brown 246c4762a1bSJed Brown PetscFunctionBeginUser; 2475be41b18SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 2485be41b18SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 2495be41b18SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 2505be41b18SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 251c4762a1bSJed Brown /* Create finite element */ 252c4762a1bSJed Brown ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 253260ffc31SMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 255c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 256c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = (*setup)(dm, user);CHKERRQ(ierr); 259c4762a1bSJed Brown while (cdm) { 260c4762a1bSJed Brown ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 261c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 262c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 265c4762a1bSJed Brown PetscFunctionReturn(0); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown static PetscErrorCode ComputeSpectral(DM dm, Vec u, PetscInt numPlanes, const PetscInt planeDir[], const PetscReal planeCoord[], AppCtx *user) 269c4762a1bSJed Brown { 270c4762a1bSJed Brown MPI_Comm comm; 271c4762a1bSJed Brown PetscSection coordSection, section; 272c4762a1bSJed Brown Vec coordinates, uloc; 273c4762a1bSJed Brown const PetscScalar *coords, *array; 274c4762a1bSJed Brown PetscInt p; 275c4762a1bSJed Brown PetscMPIInt size, rank; 276c4762a1bSJed Brown PetscErrorCode ierr; 277c4762a1bSJed Brown 278c4762a1bSJed Brown PetscFunctionBeginUser; 279c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 280ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 281ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 282c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &uloc);CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm, u, INSERT_VALUES, uloc);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, uloc, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = VecViewFromOptions(uloc, NULL, "-sol_view");CHKERRQ(ierr); 287c4762a1bSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = VecGetArrayRead(uloc, &array);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 292c4762a1bSJed Brown for (p = 0; p < numPlanes; ++p) { 293c4762a1bSJed Brown DMLabel label; 294c4762a1bSJed Brown char name[PETSC_MAX_PATH_LEN]; 295c4762a1bSJed Brown Mat F; 296c4762a1bSJed Brown Vec x, y; 297c4762a1bSJed Brown IS stratum; 298c4762a1bSJed Brown PetscReal *ray, *gray; 299c4762a1bSJed Brown PetscScalar *rvals, *svals, *gsvals; 300c4762a1bSJed Brown PetscInt *perm, *nperm; 301c4762a1bSJed Brown PetscInt n, N, i, j, off, offu; 302c4762a1bSJed Brown const PetscInt *points; 303c4762a1bSJed Brown 304c4762a1bSJed Brown ierr = PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "spectral_plane_%D", p);CHKERRQ(ierr); 305c4762a1bSJed Brown ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = DMLabelGetStratumIS(label, 1, &stratum);CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = ISGetLocalSize(stratum, &n);CHKERRQ(ierr); 308c4762a1bSJed Brown ierr = ISGetIndices(stratum, &points);CHKERRQ(ierr); 309c4762a1bSJed Brown ierr = PetscMalloc2(n, &ray, n, &svals);CHKERRQ(ierr); 310c4762a1bSJed Brown for (i = 0; i < n; ++i) { 311c4762a1bSJed Brown ierr = PetscSectionGetOffset(coordSection, points[i], &off);CHKERRQ(ierr); 312c4762a1bSJed Brown ierr = PetscSectionGetOffset(section, points[i], &offu);CHKERRQ(ierr); 313c4762a1bSJed Brown ray[i] = PetscRealPart(coords[off+((planeDir[p]+1)%2)]); 314c4762a1bSJed Brown svals[i] = array[offu]; 315c4762a1bSJed Brown } 316c4762a1bSJed Brown /* Gather the ray data to proc 0 */ 317c4762a1bSJed Brown if (size > 1) { 318c4762a1bSJed Brown PetscMPIInt *cnt, *displs, p; 319c4762a1bSJed Brown 320c4762a1bSJed Brown ierr = PetscCalloc2(size, &cnt, size, &displs);CHKERRQ(ierr); 321ffc4695bSBarry Smith ierr = MPI_Gather(&n, 1, MPIU_INT, cnt, 1, MPIU_INT, 0, comm);CHKERRMPI(ierr); 322c4762a1bSJed Brown for (p = 1; p < size; ++p) displs[p] = displs[p-1] + cnt[p-1]; 323c4762a1bSJed Brown N = displs[size-1] + cnt[size-1]; 324c4762a1bSJed Brown ierr = PetscMalloc2(N, &gray, N, &gsvals);CHKERRQ(ierr); 325ffc4695bSBarry Smith ierr = MPI_Gatherv(ray, n, MPIU_REAL, gray, cnt, displs, MPIU_REAL, 0, comm);CHKERRMPI(ierr); 326ffc4695bSBarry Smith ierr = MPI_Gatherv(svals, n, MPIU_SCALAR, gsvals, cnt, displs, MPIU_SCALAR, 0, comm);CHKERRMPI(ierr); 327c4762a1bSJed Brown ierr = PetscFree2(cnt, displs);CHKERRQ(ierr); 328c4762a1bSJed Brown } else { 329c4762a1bSJed Brown N = n; 330c4762a1bSJed Brown gray = ray; 331c4762a1bSJed Brown gsvals = svals; 332c4762a1bSJed Brown } 333dd400576SPatrick Sanan if (rank == 0) { 334c4762a1bSJed Brown /* Sort point along ray */ 335c4762a1bSJed Brown ierr = PetscMalloc2(N, &perm, N, &nperm);CHKERRQ(ierr); 336c4762a1bSJed Brown for (i = 0; i < N; ++i) {perm[i] = i;} 337c4762a1bSJed Brown ierr = PetscSortRealWithPermutation(N, gray, perm);CHKERRQ(ierr); 338c4762a1bSJed Brown /* Count duplicates and squish mapping */ 339c4762a1bSJed Brown nperm[0] = perm[0]; 340c4762a1bSJed Brown for (i = 1, j = 1; i < N; ++i) { 341c4762a1bSJed Brown if (PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) > PETSC_SMALL) nperm[j++] = perm[i]; 342c4762a1bSJed Brown } 343c4762a1bSJed Brown /* Create FFT structs */ 344c4762a1bSJed Brown ierr = MatCreateFFT(PETSC_COMM_SELF, 1, &j, MATFFTW, &F);CHKERRQ(ierr); 345c4762a1bSJed Brown ierr = MatCreateVecs(F, &x, &y);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y, name);CHKERRQ(ierr); 347c4762a1bSJed Brown ierr = VecGetArray(x, &rvals);CHKERRQ(ierr); 348c4762a1bSJed Brown for (i = 0, j = 0; i < N; ++i) { 349c4762a1bSJed Brown if (i > 0 && PetscAbsReal(gray[perm[i]] - gray[perm[i-1]]) < PETSC_SMALL) continue; 350c4762a1bSJed Brown rvals[j] = gsvals[nperm[j]]; 351c4762a1bSJed Brown ++j; 352c4762a1bSJed Brown } 353c4762a1bSJed Brown ierr = PetscFree2(perm, nperm);CHKERRQ(ierr); 354c4762a1bSJed Brown if (size > 1) {ierr = PetscFree2(gray, gsvals);CHKERRQ(ierr);} 355c4762a1bSJed Brown ierr = VecRestoreArray(x, &rvals);CHKERRQ(ierr); 356c4762a1bSJed Brown /* Do FFT along the ray */ 357c4762a1bSJed Brown ierr = MatMult(F, x, y);CHKERRQ(ierr); 358c4762a1bSJed Brown /* Chop FFT */ 359c4762a1bSJed Brown ierr = VecChop(y, PETSC_SMALL);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = VecViewFromOptions(x, NULL, "-real_view");CHKERRQ(ierr); 361c4762a1bSJed Brown ierr = VecViewFromOptions(y, NULL, "-fft_view");CHKERRQ(ierr); 362c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 364c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 365c4762a1bSJed Brown } 366c4762a1bSJed Brown ierr = ISRestoreIndices(stratum, &points);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = ISDestroy(&stratum);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = PetscFree2(ray, svals);CHKERRQ(ierr); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = VecRestoreArrayRead(uloc, &array);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &uloc);CHKERRQ(ierr); 373c4762a1bSJed Brown PetscFunctionReturn(0); 374c4762a1bSJed Brown } 375c4762a1bSJed Brown 376c4762a1bSJed Brown int main(int argc, char **argv) 377c4762a1bSJed Brown { 378c4762a1bSJed Brown DM dm; /* Problem specification */ 379c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 380c4762a1bSJed Brown Vec u; /* Solutions */ 381c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 382c4762a1bSJed Brown PetscErrorCode ierr; 383c4762a1bSJed Brown 384c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 385c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 386c4762a1bSJed Brown /* Primal system */ 387c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 388c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 389c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 390c4762a1bSJed Brown ierr = SetupDiscretization(dm, "potential", SetupPrimalProblem, &user);CHKERRQ(ierr); 391c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 392c4762a1bSJed Brown ierr = VecSet(u, 0.0);CHKERRQ(ierr); 393c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr); 394c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 396c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 398c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr); 399ad1a7c93SMatthew Knepley if (user.viewError) { 400ad1a7c93SMatthew Knepley PetscErrorCode (*sol)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 401ad1a7c93SMatthew Knepley void *ctx; 402ad1a7c93SMatthew Knepley PetscDS ds; 403ad1a7c93SMatthew Knepley PetscReal error; 404ad1a7c93SMatthew Knepley PetscInt N; 405ad1a7c93SMatthew Knepley 406ad1a7c93SMatthew Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 407ad1a7c93SMatthew Knepley ierr = PetscDSGetExactSolution(ds, 0, &sol, &ctx);CHKERRQ(ierr); 408ad1a7c93SMatthew Knepley ierr = VecGetSize(u, &N);CHKERRQ(ierr); 409ad1a7c93SMatthew Knepley ierr = DMComputeL2Diff(dm, 0.0, &sol, &ctx, u, &error);CHKERRQ(ierr); 410ad1a7c93SMatthew Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "N: %D L2 error: %g\n", N, (double)error);CHKERRQ(ierr); 411ad1a7c93SMatthew Knepley } 412c4762a1bSJed Brown if (user.spectral) { 413c4762a1bSJed Brown PetscInt planeDir[2] = {0, 1}; 414c4762a1bSJed Brown PetscReal planeCoord[2] = {0., 1.}; 415c4762a1bSJed Brown 416c4762a1bSJed Brown ierr = ComputeSpectral(dm, u, 2, planeDir, planeCoord, &user);CHKERRQ(ierr); 417c4762a1bSJed Brown } 418c4762a1bSJed Brown /* Adjoint system */ 419c4762a1bSJed Brown if (user.adjoint) { 420c4762a1bSJed Brown DM dmAdj; 421c4762a1bSJed Brown SNES snesAdj; 422c4762a1bSJed Brown Vec uAdj; 423c4762a1bSJed Brown 424c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snesAdj);CHKERRQ(ierr); 425c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) snesAdj, "adjoint_");CHKERRQ(ierr); 426c4762a1bSJed Brown ierr = DMClone(dm, &dmAdj);CHKERRQ(ierr); 427c4762a1bSJed Brown ierr = SNESSetDM(snesAdj, dmAdj);CHKERRQ(ierr); 428c4762a1bSJed Brown ierr = SetupDiscretization(dmAdj, "adjoint", SetupAdjointProblem, &user);CHKERRQ(ierr); 429c4762a1bSJed Brown ierr = DMCreateGlobalVector(dmAdj, &uAdj);CHKERRQ(ierr); 430c4762a1bSJed Brown ierr = VecSet(uAdj, 0.0);CHKERRQ(ierr); 431c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) uAdj, "adjoint");CHKERRQ(ierr); 432c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dmAdj, &user, &user, &user);CHKERRQ(ierr); 433c4762a1bSJed Brown ierr = SNESSetFromOptions(snesAdj);CHKERRQ(ierr); 434c4762a1bSJed Brown ierr = SNESSolve(snesAdj, NULL, uAdj);CHKERRQ(ierr); 435c4762a1bSJed Brown ierr = SNESGetSolution(snesAdj, &uAdj);CHKERRQ(ierr); 436c4762a1bSJed Brown ierr = VecViewFromOptions(uAdj, NULL, "-adjoint_view");CHKERRQ(ierr); 437c4762a1bSJed Brown /* Error representation */ 438c4762a1bSJed Brown { 439c4762a1bSJed Brown DM dmErr, dmErrAux, dms[2]; 440c4762a1bSJed Brown Vec errorEst, errorL2, uErr, uErrLoc, uAdjLoc, uAdjProj; 441c4762a1bSJed Brown IS *subis; 442c4762a1bSJed Brown PetscReal errorEstTot, errorL2Norm, errorL2Tot; 443c4762a1bSJed Brown PetscInt N, i; 444b724ec41SToby Isaac PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {user.homogeneous ? trig_homogeneous_u : trig_inhomogeneous_u}; 445c4762a1bSJed Brown void (*identity[1])(PetscInt, PetscInt, PetscInt, 446c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 447c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 448c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]) = {f0_identityaux_u}; 449c4762a1bSJed Brown void *ctxs[1] = {0}; 450c4762a1bSJed Brown 451c4762a1bSJed Brown ctxs[0] = &user; 452c4762a1bSJed Brown ierr = DMClone(dm, &dmErr);CHKERRQ(ierr); 453c4762a1bSJed Brown ierr = SetupDiscretization(dmErr, "error", SetupErrorProblem, &user);CHKERRQ(ierr); 454c4762a1bSJed Brown ierr = DMGetGlobalVector(dmErr, &errorEst);CHKERRQ(ierr); 455c4762a1bSJed Brown ierr = DMGetGlobalVector(dmErr, &errorL2);CHKERRQ(ierr); 456c4762a1bSJed Brown /* Compute auxiliary data (solution and projection of adjoint solution) */ 457c4762a1bSJed Brown ierr = DMGetLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr); 459c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dmAdj, uAdj, INSERT_VALUES, uAdjLoc);CHKERRQ(ierr); 460c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &uAdjProj);CHKERRQ(ierr); 4619a2a23afSMatthew G. Knepley ierr = DMSetAuxiliaryVec(dm, NULL, 0, uAdjLoc);CHKERRQ(ierr); 462c4762a1bSJed Brown ierr = DMProjectField(dm, 0.0, u, identity, INSERT_VALUES, uAdjProj);CHKERRQ(ierr); 4639a2a23afSMatthew G. Knepley ierr = DMSetAuxiliaryVec(dm, NULL, 0, NULL);CHKERRQ(ierr); 464c4762a1bSJed Brown ierr = DMRestoreLocalVector(dmAdj, &uAdjLoc);CHKERRQ(ierr); 465c4762a1bSJed Brown /* Attach auxiliary data */ 466c4762a1bSJed Brown dms[0] = dm; dms[1] = dm; 467c4762a1bSJed Brown ierr = DMCreateSuperDM(dms, 2, &subis, &dmErrAux);CHKERRQ(ierr); 468c4762a1bSJed Brown if (0) { 469c4762a1bSJed Brown PetscSection sec; 470c4762a1bSJed Brown 471c4762a1bSJed Brown ierr = DMGetLocalSection(dms[0], &sec);CHKERRQ(ierr); 472c4762a1bSJed Brown ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 473c4762a1bSJed Brown ierr = DMGetLocalSection(dms[1], &sec);CHKERRQ(ierr); 474c4762a1bSJed Brown ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 475c4762a1bSJed Brown ierr = DMGetLocalSection(dmErrAux, &sec);CHKERRQ(ierr); 476c4762a1bSJed Brown ierr = PetscSectionView(sec, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 477c4762a1bSJed Brown } 478c4762a1bSJed Brown ierr = DMViewFromOptions(dmErrAux, NULL, "-dm_err_view");CHKERRQ(ierr); 479c4762a1bSJed Brown ierr = ISViewFromOptions(subis[0], NULL, "-super_is_view");CHKERRQ(ierr); 480c4762a1bSJed Brown ierr = ISViewFromOptions(subis[1], NULL, "-super_is_view");CHKERRQ(ierr); 481c4762a1bSJed Brown ierr = DMGetGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr); 482c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-map_vec_view");CHKERRQ(ierr); 483c4762a1bSJed Brown ierr = VecViewFromOptions(uAdjProj, NULL, "-map_vec_view");CHKERRQ(ierr); 484c4762a1bSJed Brown ierr = VecViewFromOptions(uErr, NULL, "-map_vec_view");CHKERRQ(ierr); 485c4762a1bSJed Brown ierr = VecISCopy(uErr, subis[0], SCATTER_FORWARD, u);CHKERRQ(ierr); 486c4762a1bSJed Brown ierr = VecISCopy(uErr, subis[1], SCATTER_FORWARD, uAdjProj);CHKERRQ(ierr); 487c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &uAdjProj);CHKERRQ(ierr); 488c4762a1bSJed Brown for (i = 0; i < 2; ++i) {ierr = ISDestroy(&subis[i]);CHKERRQ(ierr);} 489c4762a1bSJed Brown ierr = PetscFree(subis);CHKERRQ(ierr); 490c4762a1bSJed Brown ierr = DMGetLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr); 491c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr); 492c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm, uErr, INSERT_VALUES, uErrLoc);CHKERRQ(ierr); 493c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmErrAux, &uErr);CHKERRQ(ierr); 4949a2a23afSMatthew G. Knepley ierr = DMSetAuxiliaryVec(dmAdj, NULL, 0, uErrLoc);CHKERRQ(ierr); 495c4762a1bSJed Brown /* Compute cellwise error estimate */ 496c4762a1bSJed Brown ierr = VecSet(errorEst, 0.0);CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = DMPlexComputeCellwiseIntegralFEM(dmAdj, uAdj, errorEst, &user);CHKERRQ(ierr); 4989a2a23afSMatthew G. Knepley ierr = DMSetAuxiliaryVec(dmAdj, NULL, 0, NULL);CHKERRQ(ierr); 499c4762a1bSJed Brown ierr = DMRestoreLocalVector(dmErrAux, &uErrLoc);CHKERRQ(ierr); 500c4762a1bSJed Brown ierr = DMDestroy(&dmErrAux);CHKERRQ(ierr); 501c4762a1bSJed Brown /* Plot cellwise error vector */ 502c4762a1bSJed Brown ierr = VecViewFromOptions(errorEst, NULL, "-error_view");CHKERRQ(ierr); 503c4762a1bSJed Brown /* Compute ratio of estimate (sum over cells) with actual L_2 error */ 504c4762a1bSJed Brown ierr = DMComputeL2Diff(dm, 0.0, funcs, ctxs, u, &errorL2Norm);CHKERRQ(ierr); 505c4762a1bSJed Brown ierr = DMPlexComputeL2DiffVec(dm, 0.0, funcs, ctxs, u, errorL2);CHKERRQ(ierr); 506c4762a1bSJed Brown ierr = VecViewFromOptions(errorL2, NULL, "-l2_error_view");CHKERRQ(ierr); 507c4762a1bSJed Brown ierr = VecNorm(errorL2, NORM_INFINITY, &errorL2Tot);CHKERRQ(ierr); 508c4762a1bSJed Brown ierr = VecNorm(errorEst, NORM_INFINITY, &errorEstTot);CHKERRQ(ierr); 509c4762a1bSJed Brown ierr = VecGetSize(errorEst, &N);CHKERRQ(ierr); 510c4762a1bSJed Brown ierr = VecPointwiseDivide(errorEst, errorEst, errorL2);CHKERRQ(ierr); 511c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) errorEst, "Error ratio");CHKERRQ(ierr); 512c4762a1bSJed Brown ierr = VecViewFromOptions(errorEst, NULL, "-error_ratio_view");CHKERRQ(ierr); 513c4762a1bSJed 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); 514c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmErr, &errorEst);CHKERRQ(ierr); 515c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmErr, &errorL2);CHKERRQ(ierr); 516c4762a1bSJed Brown ierr = DMDestroy(&dmErr);CHKERRQ(ierr); 517c4762a1bSJed Brown } 518c4762a1bSJed Brown ierr = DMDestroy(&dmAdj);CHKERRQ(ierr); 519c4762a1bSJed Brown ierr = VecDestroy(&uAdj);CHKERRQ(ierr); 520c4762a1bSJed Brown ierr = SNESDestroy(&snesAdj);CHKERRQ(ierr); 521c4762a1bSJed Brown } 522c4762a1bSJed Brown /* Cleanup */ 523c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 524c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 525c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 526c4762a1bSJed Brown ierr = PetscFinalize(); 527c4762a1bSJed Brown return ierr; 528c4762a1bSJed Brown } 529c4762a1bSJed Brown 530c4762a1bSJed Brown /*TEST 531c4762a1bSJed Brown 532c4762a1bSJed Brown test: 5335be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5345be41b18SMatthew G. Knepley suffix: 2d_p1_conv 535c4762a1bSJed Brown requires: triangle 5365be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5375be41b18SMatthew G. Knepley test: 5385be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5395be41b18SMatthew G. Knepley suffix: 2d_p2_conv 5405be41b18SMatthew G. Knepley requires: triangle 5415be41b18SMatthew G. Knepley args: -potential_petscspace_degree 2 -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: 3.9 5445be41b18SMatthew G. Knepley suffix: 2d_p3_conv 5455be41b18SMatthew G. Knepley requires: triangle 5465be41b18SMatthew G. Knepley args: -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5475be41b18SMatthew G. Knepley test: 5485be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5495be41b18SMatthew G. Knepley suffix: 2d_q1_conv 55030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5515be41b18SMatthew G. Knepley test: 5525be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5535be41b18SMatthew G. Knepley suffix: 2d_q2_conv 55430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5555be41b18SMatthew G. Knepley test: 5565be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5575be41b18SMatthew G. Knepley suffix: 2d_q3_conv 55830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5595be41b18SMatthew G. Knepley test: 5605be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.9 5615be41b18SMatthew G. Knepley suffix: 2d_q1_shear_conv 56230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 5635be41b18SMatthew G. Knepley test: 5645be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.9 5655be41b18SMatthew G. Knepley suffix: 2d_q2_shear_conv 56630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 5675be41b18SMatthew G. Knepley test: 5685be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 3.9 5695be41b18SMatthew G. Knepley suffix: 2d_q3_shear_conv 57030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -shear -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 2 5715be41b18SMatthew G. Knepley test: 5720fdc7489SMatthew Knepley # Using -convest_num_refine 3 we get L_2 convergence rate: 1.7 5735be41b18SMatthew G. Knepley suffix: 3d_p1_conv 5745be41b18SMatthew G. Knepley requires: ctetgen 57530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -potential_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 5765be41b18SMatthew G. Knepley test: 5775be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5785be41b18SMatthew G. Knepley suffix: 3d_p2_conv 5795be41b18SMatthew G. Knepley requires: ctetgen 58030602db0SMatthew 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 5815be41b18SMatthew G. Knepley test: 5825be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 4.0 5835be41b18SMatthew G. Knepley suffix: 3d_p3_conv 5845be41b18SMatthew G. Knepley requires: ctetgen 58530602db0SMatthew 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 5865be41b18SMatthew G. Knepley test: 5875be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 1.8 5885be41b18SMatthew G. Knepley suffix: 3d_q1_conv 58930602db0SMatthew 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 5905be41b18SMatthew G. Knepley test: 5915be41b18SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: 2.8 5925be41b18SMatthew G. Knepley suffix: 3d_q2_conv 59330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 1 5945be41b18SMatthew G. Knepley test: 5955be41b18SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: 3.8 5965be41b18SMatthew G. Knepley suffix: 3d_q3_conv 59730602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -potential_petscspace_degree 3 -snes_convergence_estimate -convest_num_refine 1 5989139b27eSToby Isaac test: 5999139b27eSToby Isaac suffix: 2d_p1_fas_full 6009139b27eSToby Isaac requires: triangle 601*e600fa54SMatthew G. Knepley args: -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 -snes_fas_full_total \ 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 \ 6069139b27eSToby Isaac -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.05 -fas_levels_esteig_ksp_max_it 10 6079139b27eSToby Isaac test: 6089139b27eSToby Isaac suffix: 2d_p1_fas_full_homogeneous 6099139b27eSToby Isaac requires: triangle 610*e600fa54SMatthew G. Knepley args: -homogeneous -potential_petscspace_degree 1 -dm_refine_hierarchy 5 \ 6119139b27eSToby Isaac -snes_max_it 1 -snes_type fas -snes_fas_levels 5 -snes_fas_type full \ 6129139b27eSToby Isaac -fas_coarse_snes_monitor -fas_coarse_snes_max_it 1 -fas_coarse_ksp_atol 1.e-13 \ 6139139b27eSToby Isaac -fas_levels_snes_monitor -fas_levels_snes_max_it 1 -fas_levels_snes_type newtonls \ 6149139b27eSToby Isaac -fas_levels_pc_type none -fas_levels_ksp_max_it 2 -fas_levels_ksp_converged_maxits -fas_levels_ksp_type chebyshev \ 6159139b27eSToby Isaac -fas_levels_esteig_ksp_type cg -fas_levels_ksp_chebyshev_esteig 0,0.25,0,1.05 -fas_levels_esteig_ksp_max_it 10 6165be41b18SMatthew G. Knepley 617c4762a1bSJed Brown test: 618c4762a1bSJed Brown suffix: 2d_p1_scalable 6195be41b18SMatthew G. Knepley requires: triangle 6205be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 3 \ 621c4762a1bSJed Brown -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned \ 622c4762a1bSJed Brown -pc_type gamg \ 623c4762a1bSJed Brown -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 624c4762a1bSJed Brown -pc_gamg_coarse_eq_limit 1000 \ 625c4762a1bSJed Brown -pc_gamg_square_graph 1 \ 626c4762a1bSJed Brown -pc_gamg_threshold 0.05 \ 627c4762a1bSJed Brown -pc_gamg_threshold_scale .0 \ 628c4762a1bSJed Brown -mg_levels_ksp_type chebyshev \ 629c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 630c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 631c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 632c4762a1bSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 633c4762a1bSJed Brown -mg_levels_pc_type jacobi \ 634c4762a1bSJed Brown -matptap_via scalable 635c4762a1bSJed Brown test: 636c4762a1bSJed Brown suffix: 2d_p1_gmg_vcycle 637c4762a1bSJed Brown requires: triangle 6385be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 639c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg \ 640c4762a1bSJed Brown -mg_levels_ksp_max_it 1 \ 641c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 642c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 643c4762a1bSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 644c4762a1bSJed Brown -mg_levels_pc_type jacobi 645c4762a1bSJed Brown test: 646c4762a1bSJed Brown suffix: 2d_p1_gmg_fcycle 647c4762a1bSJed Brown requires: triangle 6485be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 649c4762a1bSJed Brown -ksp_rtol 5e-10 -pc_type mg -pc_mg_type full \ 650c4762a1bSJed Brown -mg_levels_ksp_max_it 2 \ 651c4762a1bSJed Brown -mg_levels_esteig_ksp_type cg \ 652c4762a1bSJed Brown -mg_levels_esteig_ksp_max_it 10 \ 653c4762a1bSJed Brown -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 654c4762a1bSJed Brown -mg_levels_pc_type jacobi 655c4762a1bSJed Brown test: 656d6837840SMatthew G. Knepley suffix: 2d_p1_gmg_vcycle_adapt 657d6837840SMatthew G. Knepley requires: triangle bamg 6585be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 659d6837840SMatthew 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 \ 660d6837840SMatthew G. Knepley -mg_levels_ksp_max_it 1 \ 661d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_type cg \ 662d6837840SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 \ 663d6837840SMatthew G. Knepley -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 \ 664d6837840SMatthew G. Knepley -mg_levels_pc_type jacobi 665d6837840SMatthew G. Knepley test: 666c4762a1bSJed Brown suffix: 2d_p1_spectral_0 667c4762a1bSJed Brown requires: triangle fftw !complex 6685be41b18SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -potential_petscspace_degree 1 -dm_refine 6 -spectral -fft_view 669c4762a1bSJed Brown test: 670c4762a1bSJed Brown suffix: 2d_p1_spectral_1 671c4762a1bSJed Brown requires: triangle fftw !complex 672c4762a1bSJed Brown nsize: 2 673*e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 4,4 -potential_petscspace_degree 1 -spectral -fft_view 674c4762a1bSJed Brown test: 675c4762a1bSJed Brown suffix: 2d_p1_adj_0 676c4762a1bSJed Brown requires: triangle 6775be41b18SMatthew G. Knepley args: -potential_petscspace_degree 1 -dm_refine 1 -adjoint -adjoint_petscspace_degree 1 -error_petscspace_degree 0 678b1ad28e3SJunchao Zhang test: 679b1ad28e3SJunchao Zhang nsize: 2 6803078479eSJunchao Zhang requires: !sycl kokkos_kernels 681b1ad28e3SJunchao Zhang suffix: kokkos 682*e600fa54SMatthew 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 \ 683b1ad28e3SJunchao 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 \ 684b1ad28e3SJunchao Zhang -pc_gamg_threshold_scale .5 -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 2 -mg_levels_esteig_ksp_type cg -mg_levels_esteig_ksp_max_it 10 \ 685b1ad28e3SJunchao Zhang -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -ksp_monitor -snes_monitor -dm_view -dm_mat_type aijkokkos -dm_vec_type kokkos 686c4762a1bSJed Brown 687c4762a1bSJed Brown TEST*/ 688