1*dc2eb70dSMatthew G. Knepley static const char help[] = "Tests for injecting basis functions"; 2*dc2eb70dSMatthew G. Knepley 3*dc2eb70dSMatthew G. Knepley #include <petscdmplex.h> 4*dc2eb70dSMatthew G. Knepley #include <petscfe.h> 5*dc2eb70dSMatthew G. Knepley #include <petscds.h> 6*dc2eb70dSMatthew G. Knepley 7*dc2eb70dSMatthew G. Knepley typedef struct { 8*dc2eb70dSMatthew G. Knepley PetscInt its; /* Number of replications for timing */ 9*dc2eb70dSMatthew G. Knepley } AppCtx; 10*dc2eb70dSMatthew G. Knepley 11*dc2eb70dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 12*dc2eb70dSMatthew G. Knepley { 13*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 14*dc2eb70dSMatthew G. Knepley 15*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 16*dc2eb70dSMatthew G. Knepley options->its = 1; 17*dc2eb70dSMatthew G. Knepley 18*dc2eb70dSMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "FE Injection Options", "PETSCFE");CHKERRQ(ierr); 19*dc2eb70dSMatthew G. Knepley ierr = PetscOptionsInt("-its", "The number of replications for timing", "ex1.c", options->its, &options->its, NULL);CHKERRQ(ierr); 20*dc2eb70dSMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 21*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 22*dc2eb70dSMatthew G. Knepley } 23*dc2eb70dSMatthew G. Knepley 24*dc2eb70dSMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 25*dc2eb70dSMatthew G. Knepley { 26*dc2eb70dSMatthew G. Knepley PetscInt d; 27*dc2eb70dSMatthew G. Knepley *u = 0.0; 28*dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]); 29*dc2eb70dSMatthew G. Knepley return 0; 30*dc2eb70dSMatthew G. Knepley } 31*dc2eb70dSMatthew G. Knepley 32*dc2eb70dSMatthew G. Knepley static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 33*dc2eb70dSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 34*dc2eb70dSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 35*dc2eb70dSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 36*dc2eb70dSMatthew G. Knepley { 37*dc2eb70dSMatthew G. Knepley PetscInt d; 38*dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 39*dc2eb70dSMatthew G. Knepley } 40*dc2eb70dSMatthew G. Knepley 41*dc2eb70dSMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 42*dc2eb70dSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 43*dc2eb70dSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 44*dc2eb70dSMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 45*dc2eb70dSMatthew G. Knepley { 46*dc2eb70dSMatthew G. Knepley PetscInt d; 47*dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 48*dc2eb70dSMatthew G. Knepley } 49*dc2eb70dSMatthew G. Knepley 50*dc2eb70dSMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 51*dc2eb70dSMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 52*dc2eb70dSMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 53*dc2eb70dSMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 54*dc2eb70dSMatthew G. Knepley { 55*dc2eb70dSMatthew G. Knepley PetscInt d; 56*dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 57*dc2eb70dSMatthew G. Knepley } 58*dc2eb70dSMatthew G. Knepley 59*dc2eb70dSMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 60*dc2eb70dSMatthew G. Knepley { 61*dc2eb70dSMatthew G. Knepley PetscDS ds; 62*dc2eb70dSMatthew G. Knepley DMLabel label; 63*dc2eb70dSMatthew G. Knepley const PetscInt id = 1; 64*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 65*dc2eb70dSMatthew G. Knepley 66*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 67*dc2eb70dSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 68*dc2eb70dSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr); 69*dc2eb70dSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr); 70*dc2eb70dSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, trig_u, user);CHKERRQ(ierr); 71*dc2eb70dSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 72*dc2eb70dSMatthew G. Knepley if (label) {ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) trig_u, NULL, user, NULL);CHKERRQ(ierr);} 73*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 74*dc2eb70dSMatthew G. Knepley } 75*dc2eb70dSMatthew G. Knepley 76*dc2eb70dSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 77*dc2eb70dSMatthew G. Knepley { 78*dc2eb70dSMatthew G. Knepley DM cdm = dm; 79*dc2eb70dSMatthew G. Knepley PetscFE fe; 80*dc2eb70dSMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN]; 81*dc2eb70dSMatthew G. Knepley PetscInt dim; 82*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 83*dc2eb70dSMatthew G. Knepley 84*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 85*dc2eb70dSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 86*dc2eb70dSMatthew G. Knepley ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 87*dc2eb70dSMatthew G. Knepley ierr = DMCreateFEDefault(dm, dim, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 88*dc2eb70dSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 89*dc2eb70dSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 90*dc2eb70dSMatthew G. Knepley ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 91*dc2eb70dSMatthew G. Knepley ierr = DMCreateDS(dm);CHKERRQ(ierr); 92*dc2eb70dSMatthew G. Knepley ierr = (*setup)(dm, user);CHKERRQ(ierr); 93*dc2eb70dSMatthew G. Knepley while (cdm) { 94*dc2eb70dSMatthew G. Knepley ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 95*dc2eb70dSMatthew G. Knepley ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 96*dc2eb70dSMatthew G. Knepley } 97*dc2eb70dSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 98*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 99*dc2eb70dSMatthew G. Knepley } 100*dc2eb70dSMatthew G. Knepley 101*dc2eb70dSMatthew G. Knepley static PetscErrorCode PetscContainerUserDestroy_PetscFEGeom(void *ctx) 102*dc2eb70dSMatthew G. Knepley { 103*dc2eb70dSMatthew G. Knepley PetscFEGeom *geom = (PetscFEGeom *) ctx; 104*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 105*dc2eb70dSMatthew G. Knepley 106*dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 107*dc2eb70dSMatthew G. Knepley ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr); 108*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 109*dc2eb70dSMatthew G. Knepley } 110*dc2eb70dSMatthew G. Knepley 111*dc2eb70dSMatthew G. Knepley PetscErrorCode CellRangeGetFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 112*dc2eb70dSMatthew G. Knepley { 113*dc2eb70dSMatthew G. Knepley char composeStr[33] = {0}; 114*dc2eb70dSMatthew G. Knepley PetscObjectId id; 115*dc2eb70dSMatthew G. Knepley PetscContainer container; 116*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 117*dc2eb70dSMatthew G. Knepley 118*dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 119*dc2eb70dSMatthew G. Knepley ierr = PetscObjectGetId((PetscObject) quad, &id);CHKERRQ(ierr); 120*dc2eb70dSMatthew G. Knepley ierr = PetscSNPrintf(composeStr, 32, "CellRangeGetFEGeom_%x\n", id);CHKERRQ(ierr); 121*dc2eb70dSMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) cellIS, composeStr, (PetscObject *) &container);CHKERRQ(ierr); 122*dc2eb70dSMatthew G. Knepley if (container) { 123*dc2eb70dSMatthew G. Knepley ierr = PetscContainerGetPointer(container, (void **) geom);CHKERRQ(ierr); 124*dc2eb70dSMatthew G. Knepley } else { 125*dc2eb70dSMatthew G. Knepley ierr = DMFieldCreateFEGeom(coordField, cellIS, quad, faceData, geom);CHKERRQ(ierr); 126*dc2eb70dSMatthew G. Knepley ierr = PetscContainerCreate(PETSC_COMM_SELF, &container);CHKERRQ(ierr); 127*dc2eb70dSMatthew G. Knepley ierr = PetscContainerSetPointer(container, (void *) *geom);CHKERRQ(ierr); 128*dc2eb70dSMatthew G. Knepley ierr = PetscContainerSetUserDestroy(container, PetscContainerUserDestroy_PetscFEGeom);CHKERRQ(ierr); 129*dc2eb70dSMatthew G. Knepley ierr = PetscObjectCompose((PetscObject) cellIS, composeStr, (PetscObject) container);CHKERRQ(ierr); 130*dc2eb70dSMatthew G. Knepley ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 131*dc2eb70dSMatthew G. Knepley } 132*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 133*dc2eb70dSMatthew G. Knepley } 134*dc2eb70dSMatthew G. Knepley 135*dc2eb70dSMatthew G. Knepley PetscErrorCode CellRangeRestoreFEGeom(IS cellIS, DMField coordField, PetscQuadrature quad, PetscBool faceData, PetscFEGeom **geom) 136*dc2eb70dSMatthew G. Knepley { 137*dc2eb70dSMatthew G. Knepley PetscFunctionBegin; 138*dc2eb70dSMatthew G. Knepley *geom = NULL; 139*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 140*dc2eb70dSMatthew G. Knepley } 141*dc2eb70dSMatthew G. Knepley 142*dc2eb70dSMatthew G. Knepley static PetscErrorCode CreateFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 143*dc2eb70dSMatthew G. Knepley { 144*dc2eb70dSMatthew G. Knepley DMField coordField; 145*dc2eb70dSMatthew G. Knepley PetscInt Nf, f, maxDegree; 146*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 147*dc2eb70dSMatthew G. Knepley 148*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 149*dc2eb70dSMatthew G. Knepley *affineQuad = NULL; 150*dc2eb70dSMatthew G. Knepley *affineGeom = NULL; 151*dc2eb70dSMatthew G. Knepley *quads = NULL; 152*dc2eb70dSMatthew G. Knepley *geoms = NULL; 153*dc2eb70dSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 154*dc2eb70dSMatthew G. Knepley ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 155*dc2eb70dSMatthew G. Knepley ierr = DMFieldGetDegree(coordField, cellIS, NULL, &maxDegree);CHKERRQ(ierr); 156*dc2eb70dSMatthew G. Knepley if (maxDegree <= 1) { 157*dc2eb70dSMatthew G. Knepley ierr = DMFieldCreateDefaultQuadrature(coordField, cellIS, affineQuad);CHKERRQ(ierr); 158*dc2eb70dSMatthew G. Knepley if (*affineQuad) {ierr = CellRangeGetFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom);CHKERRQ(ierr);} 159*dc2eb70dSMatthew G. Knepley } else { 160*dc2eb70dSMatthew G. Knepley ierr = PetscCalloc2(Nf, quads, Nf, geoms);CHKERRQ(ierr); 161*dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 162*dc2eb70dSMatthew G. Knepley PetscFE fe; 163*dc2eb70dSMatthew G. Knepley 164*dc2eb70dSMatthew G. Knepley ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr); 165*dc2eb70dSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &(*quads)[f]);CHKERRQ(ierr); 166*dc2eb70dSMatthew G. Knepley ierr = PetscObjectReference((PetscObject) (*quads)[f]);CHKERRQ(ierr); 167*dc2eb70dSMatthew G. Knepley ierr = CellRangeGetFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]);CHKERRQ(ierr); 168*dc2eb70dSMatthew G. Knepley } 169*dc2eb70dSMatthew G. Knepley } 170*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 171*dc2eb70dSMatthew G. Knepley } 172*dc2eb70dSMatthew G. Knepley 173*dc2eb70dSMatthew G. Knepley static PetscErrorCode DestroyFEGeometry(DM dm, PetscDS ds, IS cellIS, PetscQuadrature *affineQuad, PetscFEGeom **affineGeom, PetscQuadrature **quads, PetscFEGeom ***geoms) 174*dc2eb70dSMatthew G. Knepley { 175*dc2eb70dSMatthew G. Knepley DMField coordField; 176*dc2eb70dSMatthew G. Knepley PetscInt Nf, f; 177*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 178*dc2eb70dSMatthew G. Knepley 179*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 180*dc2eb70dSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 181*dc2eb70dSMatthew G. Knepley ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr); 182*dc2eb70dSMatthew G. Knepley if (*affineQuad) { 183*dc2eb70dSMatthew G. Knepley ierr = CellRangeRestoreFEGeom(cellIS, coordField, *affineQuad, PETSC_FALSE, affineGeom);CHKERRQ(ierr); 184*dc2eb70dSMatthew G. Knepley ierr = PetscQuadratureDestroy(affineQuad);CHKERRQ(ierr); 185*dc2eb70dSMatthew G. Knepley } else { 186*dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 187*dc2eb70dSMatthew G. Knepley ierr = CellRangeRestoreFEGeom(cellIS, coordField, (*quads)[f], PETSC_FALSE, &(*geoms)[f]);CHKERRQ(ierr); 188*dc2eb70dSMatthew G. Knepley ierr = PetscQuadratureDestroy(&(*quads)[f]);CHKERRQ(ierr); 189*dc2eb70dSMatthew G. Knepley } 190*dc2eb70dSMatthew G. Knepley ierr = PetscFree2(*quads, *geoms);CHKERRQ(ierr); 191*dc2eb70dSMatthew G. Knepley } 192*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 193*dc2eb70dSMatthew G. Knepley } 194*dc2eb70dSMatthew G. Knepley 195*dc2eb70dSMatthew G. Knepley static PetscErrorCode TestEvaluation(DM dm) 196*dc2eb70dSMatthew G. Knepley { 197*dc2eb70dSMatthew G. Knepley PetscFE fe; 198*dc2eb70dSMatthew G. Knepley PetscSpace sp; 199*dc2eb70dSMatthew G. Knepley PetscReal *points; 200*dc2eb70dSMatthew G. Knepley PetscReal *B, *D, *H; 201*dc2eb70dSMatthew G. Knepley PetscInt dim, Nb, b, Nc, c, Np, p; 202*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 203*dc2eb70dSMatthew G. Knepley 204*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 205*dc2eb70dSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 206*dc2eb70dSMatthew G. Knepley ierr = DMGetField(dm, 0, NULL, (PetscObject *) &fe);CHKERRQ(ierr); 207*dc2eb70dSMatthew G. Knepley Np = 6; 208*dc2eb70dSMatthew G. Knepley ierr = PetscMalloc1(Np*dim, &points);CHKERRQ(ierr); 209*dc2eb70dSMatthew G. Knepley if (dim == 3) { 210*dc2eb70dSMatthew G. Knepley points[0] = -1.0; points[1] = -1.0; points[2] = -1.0; 211*dc2eb70dSMatthew G. Knepley points[3] = 1.0; points[4] = -1.0; points[5] = -1.0; 212*dc2eb70dSMatthew G. Knepley points[6] = -1.0; points[7] = 1.0; points[8] = -1.0; 213*dc2eb70dSMatthew G. Knepley points[9] = -1.0; points[10] = -1.0; points[11] = 1.0; 214*dc2eb70dSMatthew G. Knepley points[12] = 1.0; points[13] = -1.0; points[14] = 1.0; 215*dc2eb70dSMatthew G. Knepley points[15] = -1.0; points[16] = 1.0; points[17] = 1.0; 216*dc2eb70dSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only for 3D right now"); 217*dc2eb70dSMatthew G. Knepley ierr = PetscFEGetBasisSpace(fe, &sp);CHKERRQ(ierr); 218*dc2eb70dSMatthew G. Knepley ierr = PetscSpaceGetDimension(sp, &Nb);CHKERRQ(ierr); 219*dc2eb70dSMatthew G. Knepley ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 220*dc2eb70dSMatthew G. Knepley ierr = DMGetWorkArray(dm, Np*Nb*Nc, MPIU_REAL, &B);CHKERRQ(ierr); 221*dc2eb70dSMatthew G. Knepley ierr = DMGetWorkArray(dm, Np*Nb*Nc*dim, MPIU_REAL, &D);CHKERRQ(ierr); 222*dc2eb70dSMatthew G. Knepley ierr = DMGetWorkArray(dm, Np*Nb*Nc*dim*dim, MPIU_REAL, &H);CHKERRQ(ierr); 223*dc2eb70dSMatthew G. Knepley ierr = PetscSpaceEvaluate(sp, Np, points, B, NULL, NULL /*D, H*/);CHKERRQ(ierr); 224*dc2eb70dSMatthew G. Knepley for (p = 0; p < Np; ++p) { 225*dc2eb70dSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT "\n", p);CHKERRQ(ierr); 226*dc2eb70dSMatthew G. Knepley for (b = 0; b < Nb; ++b) { 227*dc2eb70dSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "B[%" PetscInt_FMT "]:", b);CHKERRQ(ierr); 228*dc2eb70dSMatthew G. Knepley for (c = 0; c < Nc; ++c) {ierr = PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[(p*Nb+b)*Nc+c]);CHKERRQ(ierr);} 229*dc2eb70dSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 230*dc2eb70dSMatthew G. Knepley #if 0 231*dc2eb70dSMatthew G. Knepley for (c = 0; c < Nc; ++c) { 232*dc2eb70dSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " D[%" PetscInt_FMT ",%" PetscInt_FMT "]:", b, c);CHKERRQ(ierr); 233*dc2eb70dSMatthew G. Knepley for (d = 0; d < dim; ++d) {ierr = PetscPrintf(PETSC_COMM_SELF, " %g", (double) B[((p*Nb+b)*Nc+c)*dim+d)]);CHKERRQ(ierr);} 234*dc2eb70dSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 235*dc2eb70dSMatthew G. Knepley } 236*dc2eb70dSMatthew G. Knepley #endif 237*dc2eb70dSMatthew G. Knepley } 238*dc2eb70dSMatthew G. Knepley } 239*dc2eb70dSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, Np*Nb, MPIU_REAL, &B);CHKERRQ(ierr); 240*dc2eb70dSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, Np*Nb*dim, MPIU_REAL, &D);CHKERRQ(ierr); 241*dc2eb70dSMatthew G. Knepley ierr = DMRestoreWorkArray(dm, Np*Nb*dim*dim, MPIU_REAL, &H);CHKERRQ(ierr); 242*dc2eb70dSMatthew G. Knepley ierr = PetscFree(points);CHKERRQ(ierr); 243*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 244*dc2eb70dSMatthew G. Knepley } 245*dc2eb70dSMatthew G. Knepley 246*dc2eb70dSMatthew G. Knepley static PetscErrorCode TestIntegration(DM dm, PetscInt cbs, PetscInt its) 247*dc2eb70dSMatthew G. Knepley { 248*dc2eb70dSMatthew G. Knepley PetscDS ds; 249*dc2eb70dSMatthew G. Knepley PetscFEGeom *chunkGeom = NULL; 250*dc2eb70dSMatthew G. Knepley PetscQuadrature affineQuad, *quads = NULL; 251*dc2eb70dSMatthew G. Knepley PetscFEGeom *affineGeom, **geoms = NULL; 252*dc2eb70dSMatthew G. Knepley PetscScalar *u, *elemVec; 253*dc2eb70dSMatthew G. Knepley IS cellIS; 254*dc2eb70dSMatthew G. Knepley PetscInt depth, cStart, cEnd, cell, chunkSize = cbs, Nch = 0, Nf, f, totDim, i, k; 255*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 256*dc2eb70dSMatthew G. Knepley 257*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 258*dc2eb70dSMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 259*dc2eb70dSMatthew G. Knepley ierr = DMGetStratumIS(dm, "depth", depth, &cellIS);CHKERRQ(ierr); 260*dc2eb70dSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 261*dc2eb70dSMatthew G. Knepley ierr = DMGetCellDS(dm, cStart, &ds);CHKERRQ(ierr); 262*dc2eb70dSMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 263*dc2eb70dSMatthew G. Knepley ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 264*dc2eb70dSMatthew G. Knepley ierr = CreateFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms);CHKERRQ(ierr); 265*dc2eb70dSMatthew G. Knepley ierr = PetscMalloc2(chunkSize*totDim, &u, chunkSize*totDim, &elemVec);CHKERRQ(ierr); 266*dc2eb70dSMatthew G. Knepley /* Assumptions: 267*dc2eb70dSMatthew G. Knepley - Single field 268*dc2eb70dSMatthew G. Knepley - No input data 269*dc2eb70dSMatthew G. Knepley - No auxiliary data 270*dc2eb70dSMatthew G. Knepley - No time-dependence 271*dc2eb70dSMatthew G. Knepley */ 272*dc2eb70dSMatthew G. Knepley for (i = 0; i < its; ++i) { 273*dc2eb70dSMatthew G. Knepley for (cell = cStart; cell < cEnd; cell += chunkSize, ++Nch) { 274*dc2eb70dSMatthew G. Knepley const PetscInt cS = cell, cE = PetscMin(cS + chunkSize, cEnd), Ne = cE - cS; 275*dc2eb70dSMatthew G. Knepley 276*dc2eb70dSMatthew G. Knepley ierr = PetscArrayzero(elemVec, chunkSize*totDim);CHKERRQ(ierr); 277*dc2eb70dSMatthew G. Knepley /* TODO Replace with DMPlexGetCellFields() */ 278*dc2eb70dSMatthew G. Knepley for (k = 0; k < chunkSize*totDim; ++k) u[k] = 1.0; 279*dc2eb70dSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 280*dc2eb70dSMatthew G. Knepley PetscFormKey key; 281*dc2eb70dSMatthew G. Knepley PetscFEGeom *geom = affineGeom ? affineGeom : geoms[f]; 282*dc2eb70dSMatthew G. Knepley /* PetscQuadrature quad = affineQuad ? affineQuad : quads[f]; */ 283*dc2eb70dSMatthew G. Knepley 284*dc2eb70dSMatthew G. Knepley key.label = NULL; key.value = 0; key.field = f; 285*dc2eb70dSMatthew G. Knepley ierr = PetscFEGeomGetChunk(geom, cS, cE, &chunkGeom);CHKERRQ(ierr); 286*dc2eb70dSMatthew G. Knepley ierr = PetscFEIntegrateResidual(ds, key, Ne, chunkGeom, u, NULL, NULL, NULL, 0.0, elemVec);CHKERRQ(ierr); 287*dc2eb70dSMatthew G. Knepley } 288*dc2eb70dSMatthew G. Knepley } 289*dc2eb70dSMatthew G. Knepley } 290*dc2eb70dSMatthew G. Knepley ierr = PetscFEGeomRestoreChunk(affineGeom, cStart, cEnd, &chunkGeom);CHKERRQ(ierr); 291*dc2eb70dSMatthew G. Knepley ierr = DestroyFEGeometry(dm, ds, cellIS, &affineQuad, &affineGeom, &quads, &geoms);CHKERRQ(ierr); 292*dc2eb70dSMatthew G. Knepley ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 293*dc2eb70dSMatthew G. Knepley ierr = PetscFree2(u, elemVec);CHKERRQ(ierr); 294*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 295*dc2eb70dSMatthew G. Knepley } 296*dc2eb70dSMatthew G. Knepley 297*dc2eb70dSMatthew G. Knepley static PetscErrorCode TestUnisolvence(DM dm) 298*dc2eb70dSMatthew G. Knepley { 299*dc2eb70dSMatthew G. Knepley Mat M; 300*dc2eb70dSMatthew G. Knepley Vec v; 301*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 302*dc2eb70dSMatthew G. Knepley 303*dc2eb70dSMatthew G. Knepley PetscFunctionBeginUser; 304*dc2eb70dSMatthew G. Knepley ierr = DMGetLocalVector(dm, &v);CHKERRQ(ierr); 305*dc2eb70dSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &v);CHKERRQ(ierr); 306*dc2eb70dSMatthew G. Knepley ierr = DMCreateMassMatrix(dm, dm, &M);CHKERRQ(ierr); 307*dc2eb70dSMatthew G. Knepley ierr = MatViewFromOptions(M, NULL, "-mass_view");CHKERRQ(ierr); 308*dc2eb70dSMatthew G. Knepley ierr = MatDestroy(&M);CHKERRQ(ierr); 309*dc2eb70dSMatthew G. Knepley PetscFunctionReturn(0); 310*dc2eb70dSMatthew G. Knepley } 311*dc2eb70dSMatthew G. Knepley 312*dc2eb70dSMatthew G. Knepley int main(int argc, char **argv) 313*dc2eb70dSMatthew G. Knepley { 314*dc2eb70dSMatthew G. Knepley DM dm; 315*dc2eb70dSMatthew G. Knepley AppCtx ctx; 316*dc2eb70dSMatthew G. Knepley PetscMPIInt size; 317*dc2eb70dSMatthew G. Knepley PetscErrorCode ierr; 318*dc2eb70dSMatthew G. Knepley 319*dc2eb70dSMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 320*dc2eb70dSMatthew G. Knepley ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 321*dc2eb70dSMatthew G. Knepley PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only."); 322*dc2eb70dSMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 323*dc2eb70dSMatthew G. Knepley ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr); 324*dc2eb70dSMatthew G. Knepley ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr); 325*dc2eb70dSMatthew G. Knepley ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 326*dc2eb70dSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) dm, "Mesh");CHKERRQ(ierr); 327*dc2eb70dSMatthew G. Knepley ierr = PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_view");CHKERRQ(ierr); 328*dc2eb70dSMatthew G. Knepley ierr = SetupDiscretization(dm, "field", SetupPrimalProblem, &ctx);CHKERRQ(ierr); 329*dc2eb70dSMatthew G. Knepley ierr = TestEvaluation(dm);CHKERRQ(ierr); 330*dc2eb70dSMatthew G. Knepley ierr = TestIntegration(dm, 1, ctx.its);CHKERRQ(ierr); 331*dc2eb70dSMatthew G. Knepley ierr = TestUnisolvence(dm);CHKERRQ(ierr); 332*dc2eb70dSMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 333*dc2eb70dSMatthew G. Knepley ierr = PetscFinalize(); 334*dc2eb70dSMatthew G. Knepley return ierr; 335*dc2eb70dSMatthew G. Knepley } 336*dc2eb70dSMatthew G. Knepley 337*dc2eb70dSMatthew G. Knepley /*TEST 338*dc2eb70dSMatthew G. Knepley test: 339*dc2eb70dSMatthew G. Knepley suffix: 0 340*dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -field_petscspace_degree 0 341*dc2eb70dSMatthew G. Knepley 342*dc2eb70dSMatthew G. Knepley test: 343*dc2eb70dSMatthew G. Knepley suffix: 2 344*dc2eb70dSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \ 345*dc2eb70dSMatthew G. Knepley -field_petscspace_type sum \ 346*dc2eb70dSMatthew G. Knepley -field_petscspace_variables 3 \ 347*dc2eb70dSMatthew G. Knepley -field_petscspace_components 3 \ 348*dc2eb70dSMatthew G. Knepley -field_petscspace_sum_spaces 2 \ 349*dc2eb70dSMatthew G. Knepley -field_petscspace_sum_concatenate false \ 350*dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_variables 3 \ 351*dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_components 3 \ 352*dc2eb70dSMatthew G. Knepley -field_sumcomp_0_petscspace_degree 1 \ 353*dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_variables 3 \ 354*dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_components 3 \ 355*dc2eb70dSMatthew G. Knepley -field_sumcomp_1_petscspace_type wxy \ 356*dc2eb70dSMatthew G. Knepley -field_petscdualspace_form_degree 0 \ 357*dc2eb70dSMatthew G. Knepley -field_petscdualspace_order 1 \ 358*dc2eb70dSMatthew G. Knepley -field_petscdualspace_components 3 359*dc2eb70dSMatthew G. Knepley 360*dc2eb70dSMatthew G. Knepley TEST*/ 361