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