1*80aecfd0Sjosephpu 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 */ 16832e2b15SMatthew G. Knepley PetscBool useBlockDiagPrec; /* Use the block diagonal of the normal equations as a preconditioner */ 17*80aecfd0Sjosephpu PetscBool shape; 18c4762a1bSJed Brown PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 19c4762a1bSJed Brown } AppCtx; 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* const char *const ex2FunctionTypes[] = {"linear","x2_x4","sin","ex2FunctionTypes","EX2_FUNCTION_",0}; */ 22c4762a1bSJed Brown 23d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 24d71ae5a4SJacob Faibussowitsch { 25c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 26c4762a1bSJed Brown PetscInt d; 27c4762a1bSJed Brown 28c4762a1bSJed Brown u[0] = 0.0; 29832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d] / (ctx->L[d]); 303ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 34d71ae5a4SJacob Faibussowitsch { 35c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 36c4762a1bSJed Brown PetscInt d; 37c4762a1bSJed Brown 38c4762a1bSJed Brown u[0] = 1; 39832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d]) * PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4); 403ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 44d71ae5a4SJacob Faibussowitsch { 45c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 46c4762a1bSJed Brown 47832e2b15SMatthew G. Knepley u[0] = PetscSinScalar(2 * PETSC_PI * ctx->k * x[0] / (ctx->L[0])); 483ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 52d71ae5a4SJacob Faibussowitsch { 53c4762a1bSJed Brown char fstring[PETSC_MAX_PATH_LEN] = "linear"; 54c4762a1bSJed Brown PetscBool flag; 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscFunctionBeginUser; 57c4762a1bSJed Brown options->particlesPerCell = 1; 58c4762a1bSJed Brown options->k = 1; 59c4762a1bSJed Brown options->particleRelDx = 1.e-20; 60c4762a1bSJed Brown options->meshRelDx = 1.e-20; 61c4762a1bSJed Brown options->momentTol = 100. * PETSC_MACHINE_EPSILON; 62832e2b15SMatthew G. Knepley options->useBlockDiagPrec = PETSC_FALSE; 63*80aecfd0Sjosephpu options->shape = PETSC_FALSE; 64c4762a1bSJed Brown 65d0609cedSBarry Smith PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX"); 66*80aecfd0Sjosephpu PetscCall(PetscOptionsBool("-shape", "Test linear function projection", "ex2.c", options->shape, &options->shape, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL)); 719566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL)); 729566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "linear", &flag)); 73c4762a1bSJed Brown if (flag) { 74c4762a1bSJed Brown options->func = linear; 75c4762a1bSJed Brown } else { 769566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "sin", &flag)); 77c4762a1bSJed Brown if (flag) { 78c4762a1bSJed Brown options->func = sinx; 79c4762a1bSJed Brown } else { 809566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "x2_x4", &flag)); 81c4762a1bSJed Brown options->func = x2_x4; 8228b400f6SJacob Faibussowitsch PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown function %s", fstring); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown } 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL)); 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-block_diag_prec", "Use the block diagonal of the normal equations to precondition the particle projection", "ex2.c", options->useBlockDiagPrec, &options->useBlockDiagPrec, NULL)); 87d0609cedSBarry Smith PetscOptionsEnd(); 88c4762a1bSJed Brown 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92d71ae5a4SJacob Faibussowitsch static PetscErrorCode PerturbVertices(DM dm, AppCtx *user) 93d71ae5a4SJacob Faibussowitsch { 94c4762a1bSJed Brown PetscRandom rnd; 95c4762a1bSJed Brown PetscReal interval = user->meshRelDx; 96c4762a1bSJed Brown Vec coordinates; 97c4762a1bSJed Brown PetscScalar *coords; 98832e2b15SMatthew G. Knepley PetscReal *hh, low[3], high[3]; 99832e2b15SMatthew G. Knepley PetscInt d, cdim, cEnd, N, p, bs; 100c4762a1bSJed Brown 101c4762a1bSJed Brown PetscFunctionBeginUser; 1029566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, low, high)); 1039566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd)); 1049566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1059566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -interval, interval)); 1069566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 1079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 1099566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(cdim, &hh)); 110832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) hh[d] = (user->L[d]) / PetscPowReal(cEnd, 1. / cdim); 1119566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &N)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coordinates, &bs)); 11363a3b9bcSJacob Faibussowitsch PetscCheck(bs == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, cdim); 1149566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 115c4762a1bSJed Brown for (p = 0; p < N; p += cdim) { 116c4762a1bSJed Brown PetscScalar *coord = &coords[p], value; 117c4762a1bSJed Brown 118c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 1199566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rnd, &value)); 120832e2b15SMatthew G. Knepley coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value * hh[d]))); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown } 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1249566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 1259566063dSJacob Faibussowitsch PetscCall(PetscFree(hh)); 1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 130d71ae5a4SJacob Faibussowitsch { 131832e2b15SMatthew G. Knepley PetscReal low[3], high[3]; 132832e2b15SMatthew G. Knepley PetscInt cdim, d; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBeginUser; 1359566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1369566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1379566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1389566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 13930602db0SMatthew G. Knepley 1409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*dm, &cdim)); 1419566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(*dm, low, high)); 142832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d]; 1439566063dSJacob Faibussowitsch PetscCall(PerturbVertices(*dm, user)); 1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147d71ae5a4SJacob 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[]) 148d71ae5a4SJacob Faibussowitsch { 149c4762a1bSJed Brown g0[0] = 1.0; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEM(DM dm, AppCtx *user) 153d71ae5a4SJacob Faibussowitsch { 154c4762a1bSJed Brown PetscFE fe; 155c4762a1bSJed Brown PetscDS ds; 156832e2b15SMatthew G. Knepley DMPolytopeType ct; 157832e2b15SMatthew G. Knepley PetscBool simplex; 158832e2b15SMatthew G. Knepley PetscInt dim, cStart; 159c4762a1bSJed Brown 160c4762a1bSJed Brown PetscFunctionBeginUser; 1619566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1629566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1639566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 164832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1659566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe)); 1669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "fe")); 1679566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 1689566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1699566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 170c4762a1bSJed Brown /* Setup to form mass matrix */ 1719566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1729566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL)); 1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 177d71ae5a4SJacob Faibussowitsch { 178c4762a1bSJed Brown PetscRandom rnd, rndp; 179c4762a1bSJed Brown PetscReal interval = user->particleRelDx; 180832e2b15SMatthew G. Knepley DMPolytopeType ct; 181832e2b15SMatthew G. Knepley PetscBool simplex; 182c4762a1bSJed Brown PetscScalar value, *vals; 183c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 184c4762a1bSJed Brown PetscInt *cellid; 185832e2b15SMatthew G. Knepley PetscInt Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d; 186c4762a1bSJed Brown 187c4762a1bSJed Brown PetscFunctionBeginUser; 1889566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1899566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1909566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 191832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1929566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 1939566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 1949566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 195c4762a1bSJed Brown 1969566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1979566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 1989566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 1999566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rndp)); 2009566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rndp, -interval, interval)); 2019566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rndp)); 202c4762a1bSJed Brown 2039566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 2049566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 2059566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 2069566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 2079566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell)); 2089566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0)); 2099566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 2109566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2119566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2129566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals)); 213c4762a1bSJed Brown 2149566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 215c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 216c4762a1bSJed Brown if (Np == 1) { 2179566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 218c4762a1bSJed Brown cellid[c] = c; 219c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d]; 220c4762a1bSJed Brown } else { 2219566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 222c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 223c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 224c4762a1bSJed Brown const PetscInt n = c * Np + p; 225c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 226c4762a1bSJed Brown 227c4762a1bSJed Brown cellid[n] = c; 2289371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 2299371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rnd, &value)); 2309371c9d4SSatish Balay refcoords[d] = PetscRealPart(value); 2319371c9d4SSatish Balay sum += refcoords[d]; 2329371c9d4SSatish Balay } 2339371c9d4SSatish Balay if (simplex && sum > 0.0) 2349371c9d4SSatish Balay for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 235c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown } 238c4762a1bSJed Brown } 2399566063dSJacob Faibussowitsch PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 240c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 241c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 242c4762a1bSJed Brown const PetscInt n = c * Np + p; 243c4762a1bSJed Brown 2449371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 2459371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rndp, &value)); 2469371c9d4SSatish Balay coords[n * dim + d] += PetscRealPart(value); 2479371c9d4SSatish Balay } 2483ba16761SJacob Faibussowitsch PetscCall(user->func(dim, 0.0, &coords[n * dim], 1, &vals[c], user)); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown } 2519566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2529566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2539566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals)); 2549566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 2559566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rndp)); 2569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 2579566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261*80aecfd0Sjosephpu static PetscErrorCode CreateParticles_Shape(DM dm, DM *sw, AppCtx *user) 262*80aecfd0Sjosephpu { 263*80aecfd0Sjosephpu PetscDS prob; 264*80aecfd0Sjosephpu PetscFE fe; 265*80aecfd0Sjosephpu PetscQuadrature quad; 266*80aecfd0Sjosephpu PetscScalar *vals; 267*80aecfd0Sjosephpu PetscReal *v0, *J, *invJ, detJ, *coords, *xi0; 268*80aecfd0Sjosephpu PetscInt *cellid; 269*80aecfd0Sjosephpu const PetscReal *qpoints, *qweights; 270*80aecfd0Sjosephpu PetscInt cStart, cEnd, c, Nq, q, dim; 271*80aecfd0Sjosephpu 272*80aecfd0Sjosephpu PetscFunctionBeginUser; 273*80aecfd0Sjosephpu PetscCall(DMGetDimension(dm, &dim)); 274*80aecfd0Sjosephpu PetscCall(DMGetDS(dm, &prob)); 275*80aecfd0Sjosephpu PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 276*80aecfd0Sjosephpu PetscCall(PetscDSGetDiscretization(prob, 0, (PetscObject *)&fe)); 277*80aecfd0Sjosephpu PetscCall(PetscFEGetQuadrature(fe, &quad)); 278*80aecfd0Sjosephpu PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, &qpoints, &qweights)); 279*80aecfd0Sjosephpu PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 280*80aecfd0Sjosephpu PetscCall(DMSetType(*sw, DMSWARM)); 281*80aecfd0Sjosephpu PetscCall(DMSetDimension(*sw, dim)); 282*80aecfd0Sjosephpu PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 283*80aecfd0Sjosephpu PetscCall(DMSwarmSetCellDM(*sw, dm)); 284*80aecfd0Sjosephpu PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_REAL)); 285*80aecfd0Sjosephpu PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 286*80aecfd0Sjosephpu PetscCall(DMSwarmSetLocalSizes(*sw, (cEnd - cStart) * Nq, 0)); 287*80aecfd0Sjosephpu PetscCall(DMSetFromOptions(*sw)); 288*80aecfd0Sjosephpu PetscCall(PetscMalloc4(dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 289*80aecfd0Sjosephpu for (c = 0; c < dim; c++) xi0[c] = -1.; 290*80aecfd0Sjosephpu PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 291*80aecfd0Sjosephpu PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 292*80aecfd0Sjosephpu PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals)); 293*80aecfd0Sjosephpu for (c = cStart; c < cEnd; ++c) { 294*80aecfd0Sjosephpu for (q = 0; q < Nq; ++q) { 295*80aecfd0Sjosephpu PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 296*80aecfd0Sjosephpu cellid[c * Nq + q] = c; 297*80aecfd0Sjosephpu CoordinatesRefToReal(dim, dim, xi0, v0, J, &qpoints[q * dim], &coords[(c * Nq + q) * dim]); 298*80aecfd0Sjosephpu PetscCall(linear(dim, 0.0, &coords[(c * Nq + q) * dim], 1, &vals[c * Nq + q], (void *)user)); 299*80aecfd0Sjosephpu vals[c * Nq + q] *= qweights[q] * detJ; 300*80aecfd0Sjosephpu } 301*80aecfd0Sjosephpu } 302*80aecfd0Sjosephpu PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 303*80aecfd0Sjosephpu PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 304*80aecfd0Sjosephpu PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals)); 305*80aecfd0Sjosephpu PetscCall(PetscFree4(xi0, v0, J, invJ)); 306*80aecfd0Sjosephpu PetscCall(DMSwarmMigrate(*sw, PETSC_FALSE)); 307*80aecfd0Sjosephpu PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 308*80aecfd0Sjosephpu PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 309*80aecfd0Sjosephpu PetscFunctionReturn(PETSC_SUCCESS); 310*80aecfd0Sjosephpu } 311*80aecfd0Sjosephpu 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user) 313d71ae5a4SJacob Faibussowitsch { 314c4762a1bSJed Brown DM dm; 315c4762a1bSJed Brown const PetscReal *coords; 316c4762a1bSJed Brown const PetscScalar *w; 317c4762a1bSJed Brown PetscReal mom[3] = {0.0, 0.0, 0.0}; 318c4762a1bSJed Brown PetscInt cell, cStart, cEnd, dim; 319c4762a1bSJed Brown 320c4762a1bSJed Brown PetscFunctionBeginUser; 3219566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 3229566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 3239566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3249566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 3259566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 3269566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 327c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 328c4762a1bSJed Brown PetscInt *pidx; 329c4762a1bSJed Brown PetscInt Np, p, d; 330c4762a1bSJed Brown 3319566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx)); 332c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 333c4762a1bSJed Brown const PetscInt idx = pidx[p]; 334c4762a1bSJed Brown const PetscReal *c = &coords[idx * dim]; 335c4762a1bSJed Brown 336c4762a1bSJed Brown mom[0] += PetscRealPart(w[idx]); 337c4762a1bSJed Brown mom[1] += PetscRealPart(w[idx]) * c[0]; 338c4762a1bSJed Brown for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d]; 339c4762a1bSJed Brown } 3409566063dSJacob Faibussowitsch PetscCall(PetscFree(pidx)); 341c4762a1bSJed Brown } 3429566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 3439566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 3449566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 345712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw))); 3463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 347c4762a1bSJed Brown } 348c4762a1bSJed Brown 349d71ae5a4SJacob 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[]) 350d71ae5a4SJacob Faibussowitsch { 351c4762a1bSJed Brown f0[0] = u[0]; 352c4762a1bSJed Brown } 353c4762a1bSJed Brown 354d71ae5a4SJacob 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[]) 355d71ae5a4SJacob Faibussowitsch { 356c4762a1bSJed Brown f0[0] = x[0] * u[0]; 357c4762a1bSJed Brown } 358c4762a1bSJed Brown 359d71ae5a4SJacob 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[]) 360d71ae5a4SJacob Faibussowitsch { 361c4762a1bSJed Brown PetscInt d; 362c4762a1bSJed Brown 363c4762a1bSJed Brown f0[0] = 0.0; 364c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0]; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user) 368d71ae5a4SJacob Faibussowitsch { 369c4762a1bSJed Brown PetscDS prob; 370c4762a1bSJed Brown PetscScalar mom; 371c4762a1bSJed Brown 372c4762a1bSJed Brown PetscFunctionBeginUser; 3739566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 3749566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_1)); 3759566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 376c4762a1bSJed Brown moments[0] = PetscRealPart(mom); 3779566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_x)); 3789566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 379c4762a1bSJed Brown moments[1] = PetscRealPart(mom); 3809566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_r2)); 3819566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 382c4762a1bSJed Brown moments[2] = PetscRealPart(mom); 3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 384c4762a1bSJed Brown } 385c4762a1bSJed Brown 386d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user) 387d71ae5a4SJacob Faibussowitsch { 388c4762a1bSJed Brown MPI_Comm comm; 389c4762a1bSJed Brown KSP ksp; 390c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 391c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 392c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */ 393c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 394c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 395c4762a1bSJed Brown PetscInt m; 396c4762a1bSJed Brown 397c4762a1bSJed Brown PetscFunctionBeginUser; 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3999566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 4009566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ptof_")); 4019566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 4029566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 403c4762a1bSJed Brown 4049566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 4059566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 406c4762a1bSJed Brown 407c4762a1bSJed Brown /* make particle weight vector */ 4089566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 409c4762a1bSJed Brown 410c4762a1bSJed Brown /* create matrix RHS vector */ 4119566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rhs)); 4129566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 4139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 4149566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 415c4762a1bSJed Brown 4169566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 4179566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 4189566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 4199566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M, M)); 4209566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 4219566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, fhat)); 4229566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 4239566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 424c4762a1bSJed Brown 425c4762a1bSJed Brown /* Check moments of field */ 4269566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 4279566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 4289566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 429c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 4301dca8a05SBarry 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); 431c4762a1bSJed Brown } 432c4762a1bSJed Brown 4339566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 4349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 4359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 4369566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 4379566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 438c4762a1bSJed Brown 4393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 440c4762a1bSJed Brown } 441c4762a1bSJed Brown 442d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user) 443d71ae5a4SJacob Faibussowitsch { 444c4762a1bSJed Brown MPI_Comm comm; 445c4762a1bSJed Brown KSP ksp; 446c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 447832e2b15SMatthew G. Knepley Mat M_p, PM_p; /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */ 448c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */ 449c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 450c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 451c4762a1bSJed Brown PetscInt m; 452c4762a1bSJed Brown 453c4762a1bSJed Brown PetscFunctionBeginUser; 4549566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 455c4762a1bSJed Brown 4569566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 4579566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ftop_")); 4589566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 4599566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 4609566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 4619566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 4629566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 463c4762a1bSJed Brown /* make particle weight vector */ 4649566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 465c4762a1bSJed Brown /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */ 4669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 4679566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 4689566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 4699566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 4709566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 4719566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M, fhat, rhs)); 4729566063dSJacob Faibussowitsch if (user->useBlockDiagPrec) PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 4739371c9d4SSatish Balay else { 4749371c9d4SSatish Balay PetscCall(PetscObjectReference((PetscObject)M_p)); 4759371c9d4SSatish Balay PM_p = M_p; 4769371c9d4SSatish Balay } 4779566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M_p, PM_p)); 4789566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, rhs, f)); 4799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 4809566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 4819566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 482c4762a1bSJed Brown /* Check moments */ 4839566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 4849566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 4859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 486c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 4871dca8a05SBarry 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); 488c4762a1bSJed Brown } 4899566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 4909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 4919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 4929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PM_p)); 4939566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 4949566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 4953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 496c4762a1bSJed Brown } 497c4762a1bSJed Brown 498c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */ 499d71ae5a4SJacob Faibussowitsch static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC) 500d71ae5a4SJacob Faibussowitsch { 501c4762a1bSJed Brown DM_Plex *mesh = (DM_Plex *)dm->data; 502c4762a1bSJed Brown PetscInt debug = mesh->printFEM; 503c4762a1bSJed Brown DM dmC; 504c4762a1bSJed Brown PetscSection section; 505c4762a1bSJed Brown PetscQuadrature quad = NULL; 506c4762a1bSJed Brown PetscScalar *interpolant, *gradsum; 507c4762a1bSJed Brown PetscFEGeom fegeom; 508c4762a1bSJed Brown PetscReal *coords; 509c4762a1bSJed Brown const PetscReal *quadPoints, *quadWeights; 510c4762a1bSJed Brown PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset; 511c4762a1bSJed Brown 512c4762a1bSJed Brown PetscFunctionBegin; 5139566063dSJacob Faibussowitsch PetscCall(VecGetDM(locC, &dmC)); 5149566063dSJacob Faibussowitsch PetscCall(VecSet(locC, 0.0)); 5159566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &coordDim)); 517c4762a1bSJed Brown fegeom.dimEmbed = coordDim; 5189566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 5199566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &numFields)); 520c4762a1bSJed Brown for (field = 0; field < numFields; ++field) { 521c4762a1bSJed Brown PetscObject obj; 522c4762a1bSJed Brown PetscClassId id; 523c4762a1bSJed Brown PetscInt Nc; 524c4762a1bSJed Brown 5259566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 5269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 527c4762a1bSJed Brown if (id == PETSCFE_CLASSID) { 528c4762a1bSJed Brown PetscFE fe = (PetscFE)obj; 529c4762a1bSJed Brown 5309566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 5319566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fe, &Nc)); 532c4762a1bSJed Brown } else if (id == PETSCFV_CLASSID) { 533c4762a1bSJed Brown PetscFV fv = (PetscFV)obj; 534c4762a1bSJed Brown 5359566063dSJacob Faibussowitsch PetscCall(PetscFVGetQuadrature(fv, &quad)); 5369566063dSJacob Faibussowitsch PetscCall(PetscFVGetNumComponents(fv, &Nc)); 53763a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 538c4762a1bSJed Brown numComponents += Nc; 539c4762a1bSJed Brown } 5409566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 54163a3b9bcSJacob Faibussowitsch PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents); 5429566063dSJacob 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)); 5439566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 5449566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 545c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 546c4762a1bSJed Brown PetscInt *star = NULL; 547c4762a1bSJed Brown PetscInt starSize, st, d, fc; 548c4762a1bSJed Brown 5499566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gradsum, coordDim * numComponents)); 5509566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 551c4762a1bSJed Brown for (st = 0; st < starSize * 2; st += 2) { 552c4762a1bSJed Brown const PetscInt cell = star[st]; 553c4762a1bSJed Brown PetscScalar *grad = &gradsum[coordDim * numComponents]; 554c4762a1bSJed Brown PetscScalar *x = NULL; 555c4762a1bSJed Brown 556c4762a1bSJed Brown if ((cell < cStart) || (cell >= cEnd)) continue; 5579566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ)); 5589566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x)); 559c4762a1bSJed Brown for (field = 0, fieldOffset = 0; field < numFields; ++field) { 560c4762a1bSJed Brown PetscObject obj; 561c4762a1bSJed Brown PetscClassId id; 562c4762a1bSJed Brown PetscInt Nb, Nc, q, qc = 0; 563c4762a1bSJed Brown 5649566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(grad, coordDim * numComponents)); 5659566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 5669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 5679371c9d4SSatish Balay if (id == PETSCFE_CLASSID) { 5689371c9d4SSatish Balay PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 5699371c9d4SSatish Balay PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb)); 5709371c9d4SSatish Balay } else if (id == PETSCFV_CLASSID) { 5719371c9d4SSatish Balay PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 5729371c9d4SSatish Balay Nb = 1; 5739371c9d4SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 574c4762a1bSJed Brown for (q = 0; q < Nq; ++q) { 57563a3b9bcSJacob 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); 5769566063dSJacob Faibussowitsch if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &fegeom, q, interpolant)); 57763a3b9bcSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 578c4762a1bSJed Brown for (fc = 0; fc < Nc; ++fc) { 579c4762a1bSJed Brown const PetscReal wt = quadWeights[q * qNc + qc + fc]; 580c4762a1bSJed Brown 581c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q]; 582c4762a1bSJed Brown } 583c4762a1bSJed Brown } 584c4762a1bSJed Brown fieldOffset += Nb; 585c4762a1bSJed Brown } 5869566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x)); 587c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 588ad540459SPierre Jolivet for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d]; 589c4762a1bSJed Brown } 590c4762a1bSJed Brown if (debug) { 59163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " gradient: [", cell)); 592c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 593c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) { 5949566063dSJacob Faibussowitsch if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 5959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d]))); 596c4762a1bSJed Brown } 597c4762a1bSJed Brown } 5989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n")); 599c4762a1bSJed Brown } 600c4762a1bSJed Brown } 6019566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 6029566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES)); 603c4762a1bSJed Brown } 6049566063dSJacob Faibussowitsch PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ)); 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 606c4762a1bSJed Brown } 607c4762a1bSJed Brown 608d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user) 609d71ae5a4SJacob Faibussowitsch { 610c4762a1bSJed Brown MPI_Comm comm; 611c4762a1bSJed Brown KSP ksp; 612c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 613c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 614c4762a1bSJed Brown Vec f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */ 615c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 616c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 617c4762a1bSJed Brown PetscInt m; 618c4762a1bSJed Brown 619c4762a1bSJed Brown PetscFunctionBeginUser; 6209566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 6219566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 6229566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ptof_")); 6239566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 6249566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 6259566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &grad)); 626c4762a1bSJed Brown 6279566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 6289566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 629c4762a1bSJed Brown 630c4762a1bSJed Brown /* make particle weight vector */ 6319566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 632c4762a1bSJed Brown 633c4762a1bSJed Brown /* create matrix RHS vector */ 6349566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rhs)); 6359566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 6369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 6379566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 638c4762a1bSJed Brown 6399566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 6409566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 641c4762a1bSJed Brown 6429566063dSJacob Faibussowitsch PetscCall(InterpolateGradient(dm, fhat, grad)); 643c4762a1bSJed Brown 6449566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 6459566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M, M)); 6469566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 6479566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, grad)); 6489566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 6499566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 650c4762a1bSJed Brown 651c4762a1bSJed Brown /* Check moments of field */ 6529566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 6539566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, grad, fmoments, user)); 6549566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 655c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 6561dca8a05SBarry 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); 657c4762a1bSJed Brown } 658c4762a1bSJed Brown 6599566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 6609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 6619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 6629566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 6639566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 6649566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &grad)); 665c4762a1bSJed Brown 6663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 667c4762a1bSJed Brown } 668c4762a1bSJed Brown 669*80aecfd0Sjosephpu static PetscErrorCode TestL2Projection_Shape(DM dm, DM sw, AppCtx *user) 670*80aecfd0Sjosephpu { 671*80aecfd0Sjosephpu PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 672*80aecfd0Sjosephpu KSP ksp; 673*80aecfd0Sjosephpu Mat mass; 674*80aecfd0Sjosephpu Vec u, rhs, uproj; 675*80aecfd0Sjosephpu void *ctxs[1]; 676*80aecfd0Sjosephpu PetscReal error; 677*80aecfd0Sjosephpu 678*80aecfd0Sjosephpu PetscFunctionBeginUser; 679*80aecfd0Sjosephpu ctxs[0] = (void *)user; 680*80aecfd0Sjosephpu funcs[0] = linear; 681*80aecfd0Sjosephpu 682*80aecfd0Sjosephpu PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &u)); 683*80aecfd0Sjosephpu PetscCall(VecViewFromOptions(u, NULL, "-f_view")); 684*80aecfd0Sjosephpu PetscCall(DMGetGlobalVector(dm, &rhs)); 685*80aecfd0Sjosephpu PetscCall(DMCreateMassMatrix(sw, dm, &mass)); 686*80aecfd0Sjosephpu PetscCall(MatMultTranspose(mass, u, rhs)); 687*80aecfd0Sjosephpu PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 688*80aecfd0Sjosephpu PetscCall(MatDestroy(&mass)); 689*80aecfd0Sjosephpu PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &u)); 690*80aecfd0Sjosephpu PetscCall(DMGetGlobalVector(dm, &uproj)); 691*80aecfd0Sjosephpu PetscCall(DMCreateMatrix(dm, &mass)); 692*80aecfd0Sjosephpu PetscCall(DMPlexSNESComputeJacobianFEM(dm, uproj, mass, mass, user)); 693*80aecfd0Sjosephpu PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view")); 694*80aecfd0Sjosephpu PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 695*80aecfd0Sjosephpu PetscCall(KSPSetOperators(ksp, mass, mass)); 696*80aecfd0Sjosephpu PetscCall(KSPSetFromOptions(ksp)); 697*80aecfd0Sjosephpu PetscCall(KSPSolve(ksp, rhs, uproj)); 698*80aecfd0Sjosephpu PetscCall(KSPDestroy(&ksp)); 699*80aecfd0Sjosephpu PetscCall(PetscObjectSetName((PetscObject)uproj, "Full Projection")); 700*80aecfd0Sjosephpu PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view")); 701*80aecfd0Sjosephpu PetscCall(DMComputeL2Diff(dm, 0.0, funcs, ctxs, uproj, &error)); 702*80aecfd0Sjosephpu PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double)error)); 703*80aecfd0Sjosephpu PetscCall(MatDestroy(&mass)); 704*80aecfd0Sjosephpu PetscCall(DMRestoreGlobalVector(dm, &rhs)); 705*80aecfd0Sjosephpu PetscCall(DMRestoreGlobalVector(dm, &uproj)); 706*80aecfd0Sjosephpu PetscFunctionReturn(PETSC_SUCCESS); 707*80aecfd0Sjosephpu } 708*80aecfd0Sjosephpu 709d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 710d71ae5a4SJacob Faibussowitsch { 711c4762a1bSJed Brown MPI_Comm comm; 712c4762a1bSJed Brown DM dm, sw; 713c4762a1bSJed Brown AppCtx user; 714c4762a1bSJed Brown 715327415f7SBarry Smith PetscFunctionBeginUser; 7169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 717c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 7189566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 7199566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user)); 7209566063dSJacob Faibussowitsch PetscCall(CreateFEM(dm, &user)); 721*80aecfd0Sjosephpu if (!user.shape) { 7229566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user)); 7239566063dSJacob Faibussowitsch PetscCall(TestL2ProjectionParticlesToField(dm, sw, &user)); 7249566063dSJacob Faibussowitsch PetscCall(TestL2ProjectionFieldToParticles(dm, sw, &user)); 7259566063dSJacob Faibussowitsch PetscCall(TestFieldGradientProjection(dm, sw, &user)); 726*80aecfd0Sjosephpu } else { 727*80aecfd0Sjosephpu PetscCall(CreateParticles_Shape(dm, &sw, &user)); 728*80aecfd0Sjosephpu PetscCall(TestL2Projection_Shape(dm, sw, &user)); 729*80aecfd0Sjosephpu } 7309566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 7319566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 7329566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 733b122ec5aSJacob Faibussowitsch return 0; 734c4762a1bSJed Brown } 735c4762a1bSJed Brown 736c4762a1bSJed Brown /*TEST 737c4762a1bSJed Brown 738832e2b15SMatthew G. Knepley # Swarm does not handle complex or quad 739832e2b15SMatthew G. Knepley build: 740832e2b15SMatthew G. Knepley requires: !complex double 741832e2b15SMatthew G. Knepley 742c4762a1bSJed Brown test: 743c4762a1bSJed Brown suffix: proj_tri_0 744832e2b15SMatthew G. Knepley requires: triangle 745832e2b15SMatthew 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 746c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 747c4762a1bSJed Brown 748c4762a1bSJed Brown test: 749c4762a1bSJed Brown suffix: proj_tri_2_faces 750832e2b15SMatthew G. Knepley requires: triangle 751832e2b15SMatthew 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 752c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 753c4762a1bSJed Brown 754c4762a1bSJed Brown test: 755c4762a1bSJed Brown suffix: proj_quad_0 756832e2b15SMatthew G. Knepley requires: triangle 75730602db0SMatthew 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 758c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 759c4762a1bSJed Brown 760c4762a1bSJed Brown test: 761c4762a1bSJed Brown suffix: proj_quad_2_faces 762832e2b15SMatthew G. Knepley requires: triangle 76330602db0SMatthew 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 764c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 765c4762a1bSJed Brown 766c4762a1bSJed Brown test: 767c4762a1bSJed Brown suffix: proj_tri_5P 768832e2b15SMatthew G. Knepley requires: triangle 769832e2b15SMatthew 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 770c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 771c4762a1bSJed Brown 772c4762a1bSJed Brown test: 773c4762a1bSJed Brown suffix: proj_quad_5P 774832e2b15SMatthew G. Knepley requires: triangle 77530602db0SMatthew 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 776c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 777c4762a1bSJed Brown 778c4762a1bSJed Brown test: 779c4762a1bSJed Brown suffix: proj_tri_mdx 780832e2b15SMatthew G. Knepley requires: triangle 781832e2b15SMatthew 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 782c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 783c4762a1bSJed Brown 784c4762a1bSJed Brown test: 785c4762a1bSJed Brown suffix: proj_tri_mdx_5P 786832e2b15SMatthew G. Knepley requires: triangle 787832e2b15SMatthew 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 788c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 789c4762a1bSJed Brown 790c4762a1bSJed Brown test: 791c4762a1bSJed Brown suffix: proj_tri_3d 792832e2b15SMatthew G. Knepley requires: ctetgen 79330602db0SMatthew 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 794c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 795c4762a1bSJed Brown 796c4762a1bSJed Brown test: 797c4762a1bSJed Brown suffix: proj_tri_3d_2_faces 798832e2b15SMatthew G. Knepley requires: ctetgen 79930602db0SMatthew 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 800c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 801c4762a1bSJed Brown 802c4762a1bSJed Brown test: 803c4762a1bSJed Brown suffix: proj_tri_3d_5P 804832e2b15SMatthew G. Knepley requires: ctetgen 80530602db0SMatthew 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 806c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 807c4762a1bSJed Brown 808c4762a1bSJed Brown test: 809c4762a1bSJed Brown suffix: proj_tri_3d_mdx 810832e2b15SMatthew G. Knepley requires: ctetgen 81130602db0SMatthew 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 812c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 813c4762a1bSJed Brown 814c4762a1bSJed Brown test: 815c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P 816832e2b15SMatthew G. Knepley requires: ctetgen 81730602db0SMatthew 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 818c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 819c4762a1bSJed Brown 820c4762a1bSJed Brown test: 821c4762a1bSJed Brown suffix: proj_tri_3d_mdx_2_faces 822832e2b15SMatthew G. Knepley requires: ctetgen 82330602db0SMatthew 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 824c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 825c4762a1bSJed Brown 826c4762a1bSJed Brown test: 827c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P_2_faces 828832e2b15SMatthew G. Knepley requires: ctetgen 82930602db0SMatthew 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 830c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 831c4762a1bSJed Brown 832832e2b15SMatthew G. Knepley test: 833832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_scale 834832e2b15SMatthew G. Knepley requires: !complex 83530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 836832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 837832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 838832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 839832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none 840832e2b15SMatthew G. Knepley 841832e2b15SMatthew G. Knepley test: 842832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_prec_scale 843832e2b15SMatthew G. Knepley requires: !complex 84430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 845832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 846832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 847832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 848832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec 849*80aecfd0Sjosephpu test: 850*80aecfd0Sjosephpu suffix: proj_shape_linear_tri_2d 851*80aecfd0Sjosephpu requires: ctetgen triangle 852*80aecfd0Sjosephpu args: -shape -dm_plex_box_faces 2,2\ 853*80aecfd0Sjosephpu -dm_view -sw_view\ 854*80aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 855*80aecfd0Sjosephpu -pc_type lu 856*80aecfd0Sjosephpu test: 857*80aecfd0Sjosephpu suffix: proj_shape_linear_quad_2d 858*80aecfd0Sjosephpu requires: ctetgen triangle 859*80aecfd0Sjosephpu args: -shape -dm_plex_simplex 0 -dm_plex_box_faces 2,2\ 860*80aecfd0Sjosephpu -dm_view -sw_view\ 861*80aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 862*80aecfd0Sjosephpu -pc_type lu 863*80aecfd0Sjosephpu test: 864*80aecfd0Sjosephpu suffix: proj_shape_linear_tri_3d 865*80aecfd0Sjosephpu requires: ctetgen triangle 866*80aecfd0Sjosephpu args: -shape -dm_plex_dim 3 -dm_plex_box_faces 2,2,2\ 867*80aecfd0Sjosephpu -dm_view -sw_view\ 868*80aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 869*80aecfd0Sjosephpu -pc_type lu 870*80aecfd0Sjosephpu test: 871*80aecfd0Sjosephpu suffix: proj_shape_linear_quad_3d 872*80aecfd0Sjosephpu requires: ctetgen triangle 873*80aecfd0Sjosephpu args: -shape -dm_plex_dim 3 -dm_plex_simplex 0\ 874*80aecfd0Sjosephpu -dm_plex_box_faces 2,2,2\ 875*80aecfd0Sjosephpu -dm_view -sw_view\ 876*80aecfd0Sjosephpu -petscspace_degree 1 -petscfe_default_quadrature_order 1\ 877*80aecfd0Sjosephpu -pc_type lu 878c4762a1bSJed Brown TEST*/ 879