180aecfd0Sjosephpu static char help[] = "Tests L2 projection with DMSwarm using delta function particles and deposition of linear shape.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscfe.h> 5c4762a1bSJed Brown #include <petscdmswarm.h> 6c4762a1bSJed Brown #include <petscds.h> 7c4762a1bSJed Brown #include <petscksp.h> 8c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h> 9c4762a1bSJed Brown typedef struct { 10832e2b15SMatthew G. Knepley PetscReal L[3]; /* Dimensions of the mesh bounding box */ 11c4762a1bSJed Brown PetscInt particlesPerCell; /* The number of partices per cell */ 12c4762a1bSJed Brown PetscReal particleRelDx; /* Relative particle position perturbation compared to average cell diameter h */ 13c4762a1bSJed Brown PetscReal meshRelDx; /* Relative vertex position perturbation compared to average cell diameter h */ 14c4762a1bSJed Brown PetscInt k; /* Mode number for test function */ 15c4762a1bSJed Brown PetscReal momentTol; /* Tolerance for checking moment conservation */ 1680aecfd0Sjosephpu PetscBool shape; 17c4762a1bSJed Brown PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 18c4762a1bSJed Brown } AppCtx; 19c4762a1bSJed Brown 20c4762a1bSJed Brown /* const char *const ex2FunctionTypes[] = {"linear","x2_x4","sin","ex2FunctionTypes","EX2_FUNCTION_",0}; */ 21c4762a1bSJed Brown 22d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 23d71ae5a4SJacob Faibussowitsch { 24c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 25c4762a1bSJed Brown PetscInt d; 26c4762a1bSJed Brown 27c4762a1bSJed Brown u[0] = 0.0; 28832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d] / (ctx->L[d]); 293ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32d71ae5a4SJacob Faibussowitsch static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 33d71ae5a4SJacob Faibussowitsch { 34c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 35c4762a1bSJed Brown PetscInt d; 36c4762a1bSJed Brown 37c4762a1bSJed Brown u[0] = 1; 38832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d]) * PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4); 393ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42d71ae5a4SJacob Faibussowitsch static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 43d71ae5a4SJacob Faibussowitsch { 44c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 45c4762a1bSJed Brown 46832e2b15SMatthew G. Knepley u[0] = PetscSinScalar(2 * PETSC_PI * ctx->k * x[0] / (ctx->L[0])); 473ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 51d71ae5a4SJacob Faibussowitsch { 52c4762a1bSJed Brown char fstring[PETSC_MAX_PATH_LEN] = "linear"; 53c4762a1bSJed Brown PetscBool flag; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBeginUser; 56c4762a1bSJed Brown options->particlesPerCell = 1; 57c4762a1bSJed Brown options->k = 1; 58c4762a1bSJed Brown options->particleRelDx = 1.e-20; 59c4762a1bSJed Brown options->meshRelDx = 1.e-20; 60c4762a1bSJed Brown options->momentTol = 100. * PETSC_MACHINE_EPSILON; 6180aecfd0Sjosephpu options->shape = PETSC_FALSE; 62c4762a1bSJed Brown 63d0609cedSBarry Smith PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX"); 6480aecfd0Sjosephpu PetscCall(PetscOptionsBool("-shape", "Test linear function projection", "ex2.c", options->shape, &options->shape, NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL)); 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "linear", &flag)); 71c4762a1bSJed Brown if (flag) { 72c4762a1bSJed Brown options->func = linear; 73c4762a1bSJed Brown } else { 749566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "sin", &flag)); 75c4762a1bSJed Brown if (flag) { 76c4762a1bSJed Brown options->func = sinx; 77c4762a1bSJed Brown } else { 789566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "x2_x4", &flag)); 79c4762a1bSJed Brown options->func = x2_x4; 8028b400f6SJacob Faibussowitsch PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown function %s", fstring); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown } 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL)); 84d0609cedSBarry Smith PetscOptionsEnd(); 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode PerturbVertices(DM dm, AppCtx *user) 89d71ae5a4SJacob Faibussowitsch { 90c4762a1bSJed Brown PetscRandom rnd; 91c4762a1bSJed Brown PetscReal interval = user->meshRelDx; 92c4762a1bSJed Brown Vec coordinates; 93c4762a1bSJed Brown PetscScalar *coords; 94832e2b15SMatthew G. Knepley PetscReal *hh, low[3], high[3]; 95832e2b15SMatthew G. Knepley PetscInt d, cdim, cEnd, N, p, bs; 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscFunctionBeginUser; 989566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, low, high)); 999566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd)); 1009566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1019566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -interval, interval)); 1029566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 1039566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1049566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 1059566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(cdim, &hh)); 106832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) hh[d] = (user->L[d]) / PetscPowReal(cEnd, 1. / cdim); 1079566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &N)); 1089566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coordinates, &bs)); 10963a3b9bcSJacob Faibussowitsch PetscCheck(bs == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, cdim); 1109566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 111c4762a1bSJed Brown for (p = 0; p < N; p += cdim) { 112c4762a1bSJed Brown PetscScalar *coord = &coords[p], value; 113c4762a1bSJed Brown 114c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 1159566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rnd, &value)); 116832e2b15SMatthew G. Knepley coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value * hh[d]))); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown } 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1209566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 1219566063dSJacob Faibussowitsch PetscCall(PetscFree(hh)); 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 126d71ae5a4SJacob Faibussowitsch { 127832e2b15SMatthew G. Knepley PetscReal low[3], high[3]; 128832e2b15SMatthew G. Knepley PetscInt cdim, d; 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscFunctionBeginUser; 1319566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1329566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1339566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1349566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 13530602db0SMatthew G. Knepley 1369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*dm, &cdim)); 1379566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(*dm, low, high)); 138832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d]; 1399566063dSJacob Faibussowitsch PetscCall(PerturbVertices(*dm, user)); 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143d71ae5a4SJacob Faibussowitsch static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 144d71ae5a4SJacob Faibussowitsch { 145c4762a1bSJed Brown g0[0] = 1.0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEM(DM dm, AppCtx *user) 149d71ae5a4SJacob Faibussowitsch { 150c4762a1bSJed Brown PetscFE fe; 151c4762a1bSJed Brown PetscDS ds; 152832e2b15SMatthew G. Knepley DMPolytopeType ct; 153832e2b15SMatthew G. Knepley PetscInt dim, cStart; 154c4762a1bSJed Brown 155c4762a1bSJed Brown PetscFunctionBeginUser; 1569566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1579566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1589566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 159fc3eb83cSMatthew G. Knepley PetscCall(PetscFECreateByCell(PetscObjectComm((PetscObject)dm), dim, 1, ct, NULL, -1, &fe)); 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "fe")); 1619566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 1629566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 164c4762a1bSJed Brown /* Setup to form mass matrix */ 1659566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1669566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL)); 1673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 171d71ae5a4SJacob Faibussowitsch { 172c4762a1bSJed Brown PetscRandom rnd, rndp; 173c4762a1bSJed Brown PetscReal interval = user->particleRelDx; 174832e2b15SMatthew G. Knepley DMPolytopeType ct; 175832e2b15SMatthew G. Knepley PetscBool simplex; 176c4762a1bSJed Brown PetscScalar value, *vals; 177c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 178c4762a1bSJed Brown PetscInt *cellid; 179832e2b15SMatthew G. Knepley PetscInt Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d; 180c4762a1bSJed Brown 181c4762a1bSJed Brown PetscFunctionBeginUser; 1829566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1839566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1849566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 185832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1869566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 1879566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 1889566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1919566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 1929566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 1939566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rndp)); 1949566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rndp, -interval, interval)); 1959566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rndp)); 196c4762a1bSJed Brown 1979566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 1989566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 1999566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 2009566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 2019566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell)); 2029566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0)); 2039566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 2049566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2059566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2069566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals)); 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 209c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 210c4762a1bSJed Brown if (Np == 1) { 2119566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 212c4762a1bSJed Brown cellid[c] = c; 213c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d]; 214c4762a1bSJed Brown } else { 2159566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 216c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 217c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 218c4762a1bSJed Brown const PetscInt n = c * Np + p; 219c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 220c4762a1bSJed Brown 221c4762a1bSJed Brown cellid[n] = c; 2229371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 2239371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rnd, &value)); 2249371c9d4SSatish Balay refcoords[d] = PetscRealPart(value); 2259371c9d4SSatish Balay sum += refcoords[d]; 2269371c9d4SSatish Balay } 2279371c9d4SSatish Balay if (simplex && sum > 0.0) 2289371c9d4SSatish Balay for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 229c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown } 2339566063dSJacob Faibussowitsch PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 234c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 235c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 236c4762a1bSJed Brown const PetscInt n = c * Np + p; 237c4762a1bSJed Brown 2389371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 2399371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rndp, &value)); 2409371c9d4SSatish Balay coords[n * dim + d] += PetscRealPart(value); 2419371c9d4SSatish Balay } 2423ba16761SJacob Faibussowitsch PetscCall(user->func(dim, 0.0, &coords[n * dim], 1, &vals[c], user)); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } 2459566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2469566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2479566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals)); 2489566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 2499566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rndp)); 2509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 2519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 2523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 253c4762a1bSJed Brown } 254c4762a1bSJed Brown 25580aecfd0Sjosephpu static PetscErrorCode CreateParticles_Shape(DM dm, DM *sw, AppCtx *user) 25680aecfd0Sjosephpu { 25780aecfd0Sjosephpu PetscDS prob; 25880aecfd0Sjosephpu PetscFE fe; 25980aecfd0Sjosephpu PetscQuadrature quad; 26080aecfd0Sjosephpu PetscScalar *vals; 26180aecfd0Sjosephpu PetscReal *v0, *J, *invJ, detJ, *coords, *xi0; 26280aecfd0Sjosephpu PetscInt *cellid; 26380aecfd0Sjosephpu const PetscReal *qpoints, *qweights; 26480aecfd0Sjosephpu PetscInt cStart, cEnd, c, Nq, q, dim; 26580aecfd0Sjosephpu 26680aecfd0Sjosephpu PetscFunctionBeginUser; 26780aecfd0Sjosephpu PetscCall(DMGetDimension(dm, &dim)); 26880aecfd0Sjosephpu PetscCall(DMGetDS(dm, &prob)); 26980aecfd0Sjosephpu PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 27080aecfd0Sjosephpu PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe)); 27180aecfd0Sjosephpu PetscCall(PetscFEGetQuadrature(fe, &quad)); 27280aecfd0Sjosephpu PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, &qweights)); 27380aecfd0Sjosephpu PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 27480aecfd0Sjosephpu PetscCall(DMSetType(*sw, DMSWARM)); 27580aecfd0Sjosephpu PetscCall(DMSetDimension(*sw, dim)); 27680aecfd0Sjosephpu PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 27780aecfd0Sjosephpu PetscCall(DMSwarmSetCellDM(*sw, dm)); 27880aecfd0Sjosephpu PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_REAL)); 27980aecfd0Sjosephpu PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 28080aecfd0Sjosephpu PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Nq, 0)); 28180aecfd0Sjosephpu PetscCall(DMSetFromOptions(*sw)); 28280aecfd0Sjosephpu PetscCall(PetscMalloc4(dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 28380aecfd0Sjosephpu for (c = 0; c < dim; c++) xi0[c] = -1.; 28480aecfd0Sjosephpu PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 28580aecfd0Sjosephpu PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 28680aecfd0Sjosephpu PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals)); 28780aecfd0Sjosephpu for (c = cStart; c < cEnd; ++c) { 28880aecfd0Sjosephpu for (q = 0; q < Nq; ++q) { 28980aecfd0Sjosephpu PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 29080aecfd0Sjosephpu cellid[c * Nq + q] = c; 29180aecfd0Sjosephpu CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], &coords[(c * Nq + q) * dim]); 29280aecfd0Sjosephpu PetscCall(linear(dim, 0.0, &coords[(c * Nq + q) * dim], 1, &vals[c * Nq + q], (void *)user)); 29380aecfd0Sjosephpu vals[c * Nq + q] *= qweights[q] * detJ; 29480aecfd0Sjosephpu } 29580aecfd0Sjosephpu } 29680aecfd0Sjosephpu PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 29780aecfd0Sjosephpu PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 29880aecfd0Sjosephpu PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals)); 29980aecfd0Sjosephpu PetscCall(PetscFree4(xi0, v0, J, invJ)); 30080aecfd0Sjosephpu PetscCall(DMSwarmMigrate(*sw, PETSC_FALSE)); 30180aecfd0Sjosephpu PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 30280aecfd0Sjosephpu PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 30380aecfd0Sjosephpu PetscFunctionReturn(PETSC_SUCCESS); 30480aecfd0Sjosephpu } 30580aecfd0Sjosephpu 306d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user) 307d71ae5a4SJacob Faibussowitsch { 308c4762a1bSJed Brown DM dm; 309c4762a1bSJed Brown const PetscReal *coords; 310c4762a1bSJed Brown const PetscScalar *w; 311c4762a1bSJed Brown PetscReal mom[3] = {0.0, 0.0, 0.0}; 312c4762a1bSJed Brown PetscInt cell, cStart, cEnd, dim; 313c4762a1bSJed Brown 314c4762a1bSJed Brown PetscFunctionBeginUser; 3159566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 3169566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 3179566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3189566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 3199566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 3209566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 321c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 322c4762a1bSJed Brown PetscInt *pidx; 323c4762a1bSJed Brown PetscInt Np, p, d; 324c4762a1bSJed Brown 3259566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx)); 326c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 327c4762a1bSJed Brown const PetscInt idx = pidx[p]; 328c4762a1bSJed Brown const PetscReal *c = &coords[idx * dim]; 329c4762a1bSJed Brown 330c4762a1bSJed Brown mom[0] += PetscRealPart(w[idx]); 331c4762a1bSJed Brown mom[1] += PetscRealPart(w[idx]) * c[0]; 332c4762a1bSJed Brown for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d]; 333c4762a1bSJed Brown } 3349566063dSJacob Faibussowitsch PetscCall(PetscFree(pidx)); 335c4762a1bSJed Brown } 3369566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 3379566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 3389566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 339*462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw))); 3403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343d71ae5a4SJacob Faibussowitsch static void f0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 344d71ae5a4SJacob Faibussowitsch { 345c4762a1bSJed Brown f0[0] = u[0]; 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 348d71ae5a4SJacob Faibussowitsch static void f0_x(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 349d71ae5a4SJacob Faibussowitsch { 350c4762a1bSJed Brown f0[0] = x[0] * u[0]; 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 353d71ae5a4SJacob Faibussowitsch static void f0_r2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 354d71ae5a4SJacob Faibussowitsch { 355c4762a1bSJed Brown PetscInt d; 356c4762a1bSJed Brown 357c4762a1bSJed Brown f0[0] = 0.0; 358c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0]; 359c4762a1bSJed Brown } 360c4762a1bSJed Brown 361d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user) 362d71ae5a4SJacob Faibussowitsch { 363c4762a1bSJed Brown PetscDS prob; 364c4762a1bSJed Brown PetscScalar mom; 365c4762a1bSJed Brown 366c4762a1bSJed Brown PetscFunctionBeginUser; 3679566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 3689566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_1)); 3699566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 370c4762a1bSJed Brown moments[0] = PetscRealPart(mom); 3719566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_x)); 3729566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 373c4762a1bSJed Brown moments[1] = PetscRealPart(mom); 3749566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_r2)); 3759566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 376c4762a1bSJed Brown moments[2] = PetscRealPart(mom); 3773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 378c4762a1bSJed Brown } 379c4762a1bSJed Brown 380fc3eb83cSMatthew G. Knepley static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, Vec fhat, AppCtx *user) 381d71ae5a4SJacob Faibussowitsch { 382fc3eb83cSMatthew G. Knepley const char *fieldnames[1] = {"w_q"}; 383fc3eb83cSMatthew G. Knepley Vec fields[1] = {fhat}; 384fc3eb83cSMatthew G. Knepley PetscReal pmoments[3]; // \int f, \int x f, \int r^2 f 385fc3eb83cSMatthew G. Knepley PetscReal fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f 386c4762a1bSJed Brown 387c4762a1bSJed Brown PetscFunctionBeginUser; 388953494dbSMatthew G. Knepley PetscCall(DMSwarmProjectFields(sw, dm, 1, fieldnames, fields, SCATTER_FORWARD)); 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* Check moments of field */ 3919566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 3929566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 393fc3eb83cSMatthew G. Knepley PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 394fc3eb83cSMatthew G. Knepley for (PetscInt m = 0; m < 3; ++m) { 395fc3eb83cSMatthew G. Knepley PetscCheck(PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]) <= user->momentTol, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Moment %" PetscInt_FMT " error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]), 396fc3eb83cSMatthew G. Knepley user->momentTol); 397c4762a1bSJed Brown } 3983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown 401fc3eb83cSMatthew G. Knepley static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, Vec fhat, AppCtx *user) 402d71ae5a4SJacob Faibussowitsch { 403fc3eb83cSMatthew G. Knepley const char *fieldnames[1] = {"w_q"}; 404fc3eb83cSMatthew G. Knepley Vec fields[1] = {fhat}; 405fc3eb83cSMatthew G. Knepley PetscReal pmoments[3]; // \int f, \int x f, \int r^2 f 406fc3eb83cSMatthew G. Knepley PetscReal fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f 407c4762a1bSJed Brown 408c4762a1bSJed Brown PetscFunctionBeginUser; 409953494dbSMatthew G. Knepley PetscCall(DMSwarmProjectFields(sw, dm, 1, fieldnames, fields, SCATTER_REVERSE)); 410c4762a1bSJed Brown 411c4762a1bSJed Brown /* Check moments */ 4129566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 4139566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 414fc3eb83cSMatthew G. Knepley PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 415fc3eb83cSMatthew G. Knepley for (PetscInt m = 0; m < 3; ++m) { 416fc3eb83cSMatthew G. Knepley PetscCheck(PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]) <= user->momentTol, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Moment %" PetscInt_FMT " error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]), 417fc3eb83cSMatthew G. Knepley user->momentTol); 418c4762a1bSJed Brown } 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 420c4762a1bSJed Brown } 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */ 423d71ae5a4SJacob Faibussowitsch static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC) 424d71ae5a4SJacob Faibussowitsch { 425c4762a1bSJed Brown DM_Plex *mesh = (DM_Plex *)dm->data; 426c4762a1bSJed Brown PetscInt debug = mesh->printFEM; 427c4762a1bSJed Brown DM dmC; 428c4762a1bSJed Brown PetscSection section; 429c4762a1bSJed Brown PetscQuadrature quad = NULL; 430c4762a1bSJed Brown PetscScalar *interpolant, *gradsum; 431c4762a1bSJed Brown PetscFEGeom fegeom; 432c4762a1bSJed Brown PetscReal *coords; 433c4762a1bSJed Brown const PetscReal *quadPoints, *quadWeights; 434c4762a1bSJed Brown PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset; 435c4762a1bSJed Brown 436c4762a1bSJed Brown PetscFunctionBegin; 4379566063dSJacob Faibussowitsch PetscCall(VecGetDM(locC, &dmC)); 4389566063dSJacob Faibussowitsch PetscCall(VecSet(locC, 0.0)); 4399566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &coordDim)); 441c4762a1bSJed Brown fegeom.dimEmbed = coordDim; 4429566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 4439566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &numFields)); 444c4762a1bSJed Brown for (field = 0; field < numFields; ++field) { 445c4762a1bSJed Brown PetscObject obj; 446c4762a1bSJed Brown PetscClassId id; 447c4762a1bSJed Brown PetscInt Nc; 448c4762a1bSJed Brown 4499566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 451c4762a1bSJed Brown if (id == PETSCFE_CLASSID) { 452c4762a1bSJed Brown PetscFE fe = (PetscFE)obj; 453c4762a1bSJed Brown 4549566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 4559566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fe, &Nc)); 456c4762a1bSJed Brown } else if (id == PETSCFV_CLASSID) { 457c4762a1bSJed Brown PetscFV fv = (PetscFV)obj; 458c4762a1bSJed Brown 4599566063dSJacob Faibussowitsch PetscCall(PetscFVGetQuadrature(fv, &quad)); 4609566063dSJacob Faibussowitsch PetscCall(PetscFVGetNumComponents(fv, &Nc)); 46163a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 462c4762a1bSJed Brown numComponents += Nc; 463c4762a1bSJed Brown } 4649566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 46563a3b9bcSJacob Faibussowitsch PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents); 4669566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(coordDim * numComponents * 2, &gradsum, coordDim * numComponents, &interpolant, coordDim * Nq, &coords, Nq, &fegeom.detJ, coordDim * coordDim * Nq, &fegeom.J, coordDim * coordDim * Nq, &fegeom.invJ)); 4679566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4689566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 469c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 470c4762a1bSJed Brown PetscInt *star = NULL; 471c4762a1bSJed Brown PetscInt starSize, st, d, fc; 472c4762a1bSJed Brown 4739566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gradsum, coordDim * numComponents)); 4749566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 475c4762a1bSJed Brown for (st = 0; st < starSize * 2; st += 2) { 476c4762a1bSJed Brown const PetscInt cell = star[st]; 477c4762a1bSJed Brown PetscScalar *grad = &gradsum[coordDim * numComponents]; 478c4762a1bSJed Brown PetscScalar *x = NULL; 479c4762a1bSJed Brown 480c4762a1bSJed Brown if ((cell < cStart) || (cell >= cEnd)) continue; 4819566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ)); 4829566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x)); 483c4762a1bSJed Brown for (field = 0, fieldOffset = 0; field < numFields; ++field) { 484c4762a1bSJed Brown PetscObject obj; 485c4762a1bSJed Brown PetscClassId id; 486c4762a1bSJed Brown PetscInt Nb, Nc, q, qc = 0; 487c4762a1bSJed Brown 4889566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(grad, coordDim * numComponents)); 4899566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 4909566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 4919371c9d4SSatish Balay if (id == PETSCFE_CLASSID) { 4929371c9d4SSatish Balay PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 4939371c9d4SSatish Balay PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb)); 4949371c9d4SSatish Balay } else if (id == PETSCFV_CLASSID) { 4959371c9d4SSatish Balay PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 4969371c9d4SSatish Balay Nb = 1; 4979371c9d4SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 498c4762a1bSJed Brown for (q = 0; q < Nq; ++q) { 49963a3b9bcSJacob Faibussowitsch PetscCheck(fegeom.detJ[q] > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT ", quadrature points %" PetscInt_FMT, (double)fegeom.detJ[q], cell, q); 5009566063dSJacob Faibussowitsch if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &fegeom, q, interpolant)); 50163a3b9bcSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 502c4762a1bSJed Brown for (fc = 0; fc < Nc; ++fc) { 503c4762a1bSJed Brown const PetscReal wt = quadWeights[q * qNc + qc + fc]; 504c4762a1bSJed Brown 505c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q]; 506c4762a1bSJed Brown } 507c4762a1bSJed Brown } 508c4762a1bSJed Brown fieldOffset += Nb; 509c4762a1bSJed Brown } 5109566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x)); 511c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 512ad540459SPierre Jolivet for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d]; 513c4762a1bSJed Brown } 514c4762a1bSJed Brown if (debug) { 51563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " gradient: [", cell)); 516c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 517c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) { 5189566063dSJacob Faibussowitsch if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 5199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d]))); 520c4762a1bSJed Brown } 521c4762a1bSJed Brown } 5229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n")); 523c4762a1bSJed Brown } 524c4762a1bSJed Brown } 5259566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 5269566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES)); 527c4762a1bSJed Brown } 5289566063dSJacob Faibussowitsch PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ)); 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 530c4762a1bSJed Brown } 531c4762a1bSJed Brown 532d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user) 533d71ae5a4SJacob Faibussowitsch { 534c4762a1bSJed Brown MPI_Comm comm; 535c4762a1bSJed Brown KSP ksp; 536c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 537c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 538c4762a1bSJed Brown Vec f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */ 539c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 540c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 541c4762a1bSJed Brown PetscInt m; 542c4762a1bSJed Brown 543c4762a1bSJed Brown PetscFunctionBeginUser; 5449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 5459566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 5469566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ptof_")); 5479566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 5489566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 5499566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &grad)); 550c4762a1bSJed Brown 5519566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 5529566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 553c4762a1bSJed Brown 554c4762a1bSJed Brown /* make particle weight vector */ 5559566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 556c4762a1bSJed Brown 557c4762a1bSJed Brown /* create matrix RHS vector */ 5589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rhs)); 5599566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 5609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 5619566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 562c4762a1bSJed Brown 5639566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 5649566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 565c4762a1bSJed Brown 5669566063dSJacob Faibussowitsch PetscCall(InterpolateGradient(dm, fhat, grad)); 567c4762a1bSJed Brown 5689566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 5699566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M, M)); 5709566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 5719566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, grad)); 5729566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 5739566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 574c4762a1bSJed Brown 575c4762a1bSJed Brown /* Check moments of field */ 5769566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 5779566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, grad, fmoments, user)); 5789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 579c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 5801dca8a05SBarry Smith PetscCheck(PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]) <= user->momentTol, comm, PETSC_ERR_ARG_WRONG, "Moment %" PetscInt_FMT " error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]), user->momentTol); 581c4762a1bSJed Brown } 582c4762a1bSJed Brown 5839566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 5849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 5859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 5869566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 5879566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 5889566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &grad)); 5893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 590c4762a1bSJed Brown } 591c4762a1bSJed Brown 59280aecfd0Sjosephpu static PetscErrorCode TestL2Projection_Shape(DM dm, DM sw, AppCtx *user) 59380aecfd0Sjosephpu { 59480aecfd0Sjosephpu PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 59580aecfd0Sjosephpu KSP ksp; 59680aecfd0Sjosephpu Mat mass; 59780aecfd0Sjosephpu Vec u, rhs, uproj; 59880aecfd0Sjosephpu void *ctxs[1]; 59980aecfd0Sjosephpu PetscReal error; 60080aecfd0Sjosephpu 60180aecfd0Sjosephpu PetscFunctionBeginUser; 60280aecfd0Sjosephpu ctxs[0] = (void *)user; 60380aecfd0Sjosephpu funcs[0] = linear; 60480aecfd0Sjosephpu 60580aecfd0Sjosephpu PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &u)); 60680aecfd0Sjosephpu PetscCall(VecViewFromOptions(u, NULL, "-f_view")); 60780aecfd0Sjosephpu PetscCall(DMGetGlobalVector(dm, &rhs)); 60880aecfd0Sjosephpu PetscCall(DMCreateMassMatrix(sw, dm, &mass)); 60980aecfd0Sjosephpu PetscCall(MatMultTranspose(mass, u, rhs)); 61080aecfd0Sjosephpu PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 61180aecfd0Sjosephpu PetscCall(MatDestroy(&mass)); 61280aecfd0Sjosephpu PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &u)); 61380aecfd0Sjosephpu PetscCall(DMGetGlobalVector(dm, &uproj)); 61480aecfd0Sjosephpu PetscCall(DMCreateMatrix(dm, &mass)); 61580aecfd0Sjosephpu PetscCall(DMPlexSNESComputeJacobianFEM(dm, uproj, mass, mass, user)); 61680aecfd0Sjosephpu PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view")); 61780aecfd0Sjosephpu PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 61880aecfd0Sjosephpu PetscCall(KSPSetOperators(ksp, mass, mass)); 61980aecfd0Sjosephpu PetscCall(KSPSetFromOptions(ksp)); 62080aecfd0Sjosephpu PetscCall(KSPSolve(ksp, rhs, uproj)); 62180aecfd0Sjosephpu PetscCall(KSPDestroy(&ksp)); 62280aecfd0Sjosephpu PetscCall(PetscObjectSetName((PetscObject)uproj, "Full Projection")); 62380aecfd0Sjosephpu PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view")); 62480aecfd0Sjosephpu PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, uproj, &error)); 62580aecfd0Sjosephpu PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double)error)); 62680aecfd0Sjosephpu PetscCall(MatDestroy(&mass)); 62780aecfd0Sjosephpu PetscCall(DMRestoreGlobalVector(dm, &rhs)); 62880aecfd0Sjosephpu PetscCall(DMRestoreGlobalVector(dm, &uproj)); 62980aecfd0Sjosephpu PetscFunctionReturn(PETSC_SUCCESS); 63080aecfd0Sjosephpu } 63180aecfd0Sjosephpu 632d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 633d71ae5a4SJacob Faibussowitsch { 634c4762a1bSJed Brown MPI_Comm comm; 635c4762a1bSJed Brown DM dm, sw; 636c4762a1bSJed Brown AppCtx user; 637c4762a1bSJed Brown 638327415f7SBarry Smith PetscFunctionBeginUser; 6399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 640c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 6419566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 6429566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user)); 6439566063dSJacob Faibussowitsch PetscCall(CreateFEM(dm, &user)); 64480aecfd0Sjosephpu if (!user.shape) { 645fc3eb83cSMatthew G. Knepley Vec fhat; 646fc3eb83cSMatthew G. Knepley 6479566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user)); 648fc3eb83cSMatthew G. Knepley PetscCall(DMGetGlobalVector(dm, &fhat)); 649fc3eb83cSMatthew G. Knepley PetscCall(TestL2ProjectionParticlesToField(dm, sw, fhat, &user)); 650fc3eb83cSMatthew G. Knepley PetscCall(TestL2ProjectionFieldToParticles(dm, sw, fhat, &user)); 6519566063dSJacob Faibussowitsch PetscCall(TestFieldGradientProjection(dm, sw, &user)); 652fc3eb83cSMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dm, &fhat)); 65380aecfd0Sjosephpu } else { 65480aecfd0Sjosephpu PetscCall(CreateParticles_Shape(dm, &sw, &user)); 65580aecfd0Sjosephpu PetscCall(TestL2Projection_Shape(dm, sw, &user)); 65680aecfd0Sjosephpu } 6579566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6589566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 6599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 660b122ec5aSJacob Faibussowitsch return 0; 661c4762a1bSJed Brown } 662c4762a1bSJed Brown 663c4762a1bSJed Brown /*TEST 664c4762a1bSJed Brown 665832e2b15SMatthew G. Knepley # Swarm does not handle complex or quad 666832e2b15SMatthew G. Knepley build: 667832e2b15SMatthew G. Knepley requires: !complex double 668832e2b15SMatthew G. Knepley 669c4762a1bSJed Brown test: 670c4762a1bSJed Brown suffix: proj_tri_0 671832e2b15SMatthew G. Knepley requires: triangle 672832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 673c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 674c4762a1bSJed Brown 675c4762a1bSJed Brown test: 676c4762a1bSJed Brown suffix: proj_tri_2_faces 677832e2b15SMatthew G. Knepley requires: triangle 678832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 679c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 680c4762a1bSJed Brown 681c4762a1bSJed Brown test: 682c4762a1bSJed Brown suffix: proj_quad_0 683832e2b15SMatthew G. Knepley requires: triangle 68430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 685c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 686c4762a1bSJed Brown 687c4762a1bSJed Brown test: 688c4762a1bSJed Brown suffix: proj_quad_2_faces 689832e2b15SMatthew G. Knepley requires: triangle 69030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 691c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 692c4762a1bSJed Brown 693c4762a1bSJed Brown test: 694c4762a1bSJed Brown suffix: proj_tri_5P 695832e2b15SMatthew G. Knepley requires: triangle 696832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 697c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 698c4762a1bSJed Brown 699c4762a1bSJed Brown test: 700c4762a1bSJed Brown suffix: proj_quad_5P 701832e2b15SMatthew G. Knepley requires: triangle 70230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 703c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 704c4762a1bSJed Brown 705c4762a1bSJed Brown test: 706c4762a1bSJed Brown suffix: proj_tri_mdx 707832e2b15SMatthew G. Knepley requires: triangle 708832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 709c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 710c4762a1bSJed Brown 711c4762a1bSJed Brown test: 712c4762a1bSJed Brown suffix: proj_tri_mdx_5P 713832e2b15SMatthew G. Knepley requires: triangle 714832e2b15SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 715c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 716c4762a1bSJed Brown 717c4762a1bSJed Brown test: 718c4762a1bSJed Brown suffix: proj_tri_3d 719832e2b15SMatthew G. Knepley requires: ctetgen 72030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 721c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 722c4762a1bSJed Brown 723c4762a1bSJed Brown test: 724c4762a1bSJed Brown suffix: proj_tri_3d_2_faces 725832e2b15SMatthew G. Knepley requires: ctetgen 72630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 727c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 728c4762a1bSJed Brown 729c4762a1bSJed Brown test: 730c4762a1bSJed Brown suffix: proj_tri_3d_5P 731832e2b15SMatthew G. Knepley requires: ctetgen 73230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 733c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 734c4762a1bSJed Brown 735c4762a1bSJed Brown test: 736c4762a1bSJed Brown suffix: proj_tri_3d_mdx 737832e2b15SMatthew G. Knepley requires: ctetgen 73830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 739c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 740c4762a1bSJed Brown 741c4762a1bSJed Brown test: 742c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P 743832e2b15SMatthew G. Knepley requires: ctetgen 74430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 745c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 746c4762a1bSJed Brown 747c4762a1bSJed Brown test: 748c4762a1bSJed Brown suffix: proj_tri_3d_mdx_2_faces 749832e2b15SMatthew G. Knepley requires: ctetgen 75030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 751c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 752c4762a1bSJed Brown 753c4762a1bSJed Brown test: 754c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P_2_faces 755832e2b15SMatthew G. Knepley requires: ctetgen 75630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -particlesPerCell 5 -mesh_perturbation 1.0e-1 -dm_view -sw_view -petscspace_degree 2 -petscfe_default_quadrature_order {{2 3}} -ptof_pc_type lu -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none 757c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 758c4762a1bSJed Brown 759832e2b15SMatthew G. Knepley test: 760832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_scale 761832e2b15SMatthew G. Knepley requires: !complex 76230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 763832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 764832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 765832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 766832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none 767832e2b15SMatthew G. Knepley 768832e2b15SMatthew G. Knepley test: 769832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_prec_scale 770832e2b15SMatthew G. Knepley requires: !complex 77130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 772832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 773832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 774832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 775fc3eb83cSMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero 77680aecfd0Sjosephpu test: 77780aecfd0Sjosephpu suffix: proj_shape_linear_tri_2d 77880aecfd0Sjosephpu requires: ctetgen triangle 77980aecfd0Sjosephpu args: -shape -dm_plex_box_faces 2,2\ 78080aecfd0Sjosephpu -dm_view -sw_view\ 78180aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 78280aecfd0Sjosephpu -pc_type lu 78380aecfd0Sjosephpu test: 78480aecfd0Sjosephpu suffix: proj_shape_linear_quad_2d 78580aecfd0Sjosephpu requires: ctetgen triangle 78680aecfd0Sjosephpu args: -shape -dm_plex_simplex 0 -dm_plex_box_faces 2,2\ 78780aecfd0Sjosephpu -dm_view -sw_view\ 78880aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 78980aecfd0Sjosephpu -pc_type lu 79080aecfd0Sjosephpu test: 79180aecfd0Sjosephpu suffix: proj_shape_linear_tri_3d 79280aecfd0Sjosephpu requires: ctetgen triangle 79380aecfd0Sjosephpu args: -shape -dm_plex_dim 3 -dm_plex_box_faces 2,2,2\ 79480aecfd0Sjosephpu -dm_view -sw_view\ 79580aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 79680aecfd0Sjosephpu -pc_type lu 79780aecfd0Sjosephpu test: 79880aecfd0Sjosephpu suffix: proj_shape_linear_quad_3d 79980aecfd0Sjosephpu requires: ctetgen triangle 80080aecfd0Sjosephpu args: -shape -dm_plex_dim 3 -dm_plex_simplex 0\ 80180aecfd0Sjosephpu -dm_plex_box_faces 2,2,2\ 80280aecfd0Sjosephpu -dm_view -sw_view\ 80380aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 80480aecfd0Sjosephpu -pc_type lu 805c4762a1bSJed Brown TEST*/ 806