1c4762a1bSJed Brown static char help[] = "Tests L2 projection with DMSwarm using delta function particles.\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 */ 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 22*d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 23*d71ae5a4SJacob 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]); 29c4762a1bSJed Brown return 0; 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32*d71ae5a4SJacob Faibussowitsch static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 33*d71ae5a4SJacob 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); 39c4762a1bSJed Brown return 0; 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42*d71ae5a4SJacob Faibussowitsch static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) 43*d71ae5a4SJacob 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])); 47c4762a1bSJed Brown return 0; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 51*d71ae5a4SJacob 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; 61832e2b15SMatthew G. Knepley options->useBlockDiagPrec = PETSC_FALSE; 62c4762a1bSJed Brown 63d0609cedSBarry Smith PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX"); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL)); 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL)); 679566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "linear", &flag)); 70c4762a1bSJed Brown if (flag) { 71c4762a1bSJed Brown options->func = linear; 72c4762a1bSJed Brown } else { 739566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "sin", &flag)); 74c4762a1bSJed Brown if (flag) { 75c4762a1bSJed Brown options->func = sinx; 76c4762a1bSJed Brown } else { 779566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "x2_x4", &flag)); 78c4762a1bSJed Brown options->func = x2_x4; 7928b400f6SJacob Faibussowitsch PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown function %s", fstring); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL)); 839566063dSJacob 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)); 84d0609cedSBarry Smith PetscOptionsEnd(); 85c4762a1bSJed Brown 86c4762a1bSJed Brown PetscFunctionReturn(0); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PerturbVertices(DM dm, AppCtx *user) 90*d71ae5a4SJacob Faibussowitsch { 91c4762a1bSJed Brown PetscRandom rnd; 92c4762a1bSJed Brown PetscReal interval = user->meshRelDx; 93c4762a1bSJed Brown Vec coordinates; 94c4762a1bSJed Brown PetscScalar *coords; 95832e2b15SMatthew G. Knepley PetscReal *hh, low[3], high[3]; 96832e2b15SMatthew G. Knepley PetscInt d, cdim, cEnd, N, p, bs; 97c4762a1bSJed Brown 98c4762a1bSJed Brown PetscFunctionBeginUser; 999566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, low, high)); 1009566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd)); 1019566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1029566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -interval, interval)); 1039566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 1049566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1059566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 1069566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(cdim, &hh)); 107832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) hh[d] = (user->L[d]) / PetscPowReal(cEnd, 1. / cdim); 1089566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &N)); 1099566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coordinates, &bs)); 11063a3b9bcSJacob Faibussowitsch PetscCheck(bs == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, cdim); 1119566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 112c4762a1bSJed Brown for (p = 0; p < N; p += cdim) { 113c4762a1bSJed Brown PetscScalar *coord = &coords[p], value; 114c4762a1bSJed Brown 115c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 1169566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rnd, &value)); 117832e2b15SMatthew G. Knepley coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value * hh[d]))); 118c4762a1bSJed Brown } 119c4762a1bSJed Brown } 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1219566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 1229566063dSJacob Faibussowitsch PetscCall(PetscFree(hh)); 123c4762a1bSJed Brown PetscFunctionReturn(0); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) 127*d71ae5a4SJacob Faibussowitsch { 128832e2b15SMatthew G. Knepley PetscReal low[3], high[3]; 129832e2b15SMatthew G. Knepley PetscInt cdim, d; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBeginUser; 1329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1359566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 13630602db0SMatthew G. Knepley 1379566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*dm, &cdim)); 1389566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(*dm, low, high)); 139832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d]; 1409566063dSJacob Faibussowitsch PetscCall(PerturbVertices(*dm, user)); 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144*d71ae5a4SJacob 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[]) 145*d71ae5a4SJacob Faibussowitsch { 146c4762a1bSJed Brown g0[0] = 1.0; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateFEM(DM dm, AppCtx *user) 150*d71ae5a4SJacob Faibussowitsch { 151c4762a1bSJed Brown PetscFE fe; 152c4762a1bSJed Brown PetscDS ds; 153832e2b15SMatthew G. Knepley DMPolytopeType ct; 154832e2b15SMatthew G. Knepley PetscBool simplex; 155832e2b15SMatthew G. Knepley PetscInt dim, cStart; 156c4762a1bSJed Brown 157c4762a1bSJed Brown PetscFunctionBeginUser; 1589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1599566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1609566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 161832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1629566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe)); 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "fe")); 1649566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 1659566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 167c4762a1bSJed Brown /* Setup to form mass matrix */ 1689566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1699566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL)); 170c4762a1bSJed Brown PetscFunctionReturn(0); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) 174*d71ae5a4SJacob Faibussowitsch { 175c4762a1bSJed Brown PetscRandom rnd, rndp; 176c4762a1bSJed Brown PetscReal interval = user->particleRelDx; 177832e2b15SMatthew G. Knepley DMPolytopeType ct; 178832e2b15SMatthew G. Knepley PetscBool simplex; 179c4762a1bSJed Brown PetscScalar value, *vals; 180c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 181c4762a1bSJed Brown PetscInt *cellid; 182832e2b15SMatthew G. Knepley PetscInt Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBeginUser; 1859566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1869566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1879566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 188832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1899566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 1909566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 1919566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 192c4762a1bSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1949566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 1959566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 1969566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rndp)); 1979566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rndp, -interval, interval)); 1989566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rndp)); 199c4762a1bSJed Brown 2009566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 2019566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 2029566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 2039566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell)); 2059566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0)); 2069566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 2079566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2089566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2099566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals)); 210c4762a1bSJed Brown 2119566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 212c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 213c4762a1bSJed Brown if (Np == 1) { 2149566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 215c4762a1bSJed Brown cellid[c] = c; 216c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d]; 217c4762a1bSJed Brown } else { 2189566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 219c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 220c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 221c4762a1bSJed Brown const PetscInt n = c * Np + p; 222c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 223c4762a1bSJed Brown 224c4762a1bSJed Brown cellid[n] = c; 2259371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 2269371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rnd, &value)); 2279371c9d4SSatish Balay refcoords[d] = PetscRealPart(value); 2289371c9d4SSatish Balay sum += refcoords[d]; 2299371c9d4SSatish Balay } 2309371c9d4SSatish Balay if (simplex && sum > 0.0) 2319371c9d4SSatish Balay for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 232c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown } 235c4762a1bSJed Brown } 2369566063dSJacob Faibussowitsch PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 237c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 238c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 239c4762a1bSJed Brown const PetscInt n = c * Np + p; 240c4762a1bSJed Brown 2419371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 2429371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rndp, &value)); 2439371c9d4SSatish Balay coords[n * dim + d] += PetscRealPart(value); 2449371c9d4SSatish Balay } 245c4762a1bSJed Brown user->func(dim, 0.0, &coords[n * dim], 1, &vals[c], user); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown } 2489566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2499566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2509566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals)); 2519566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 2529566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rndp)); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 2549566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 255c4762a1bSJed Brown PetscFunctionReturn(0); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258*d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user) 259*d71ae5a4SJacob Faibussowitsch { 260c4762a1bSJed Brown DM dm; 261c4762a1bSJed Brown const PetscReal *coords; 262c4762a1bSJed Brown const PetscScalar *w; 263c4762a1bSJed Brown PetscReal mom[3] = {0.0, 0.0, 0.0}; 264c4762a1bSJed Brown PetscInt cell, cStart, cEnd, dim; 265c4762a1bSJed Brown 266c4762a1bSJed Brown PetscFunctionBeginUser; 2679566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 2689566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 2699566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2709566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 2719566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2729566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 273c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 274c4762a1bSJed Brown PetscInt *pidx; 275c4762a1bSJed Brown PetscInt Np, p, d; 276c4762a1bSJed Brown 2779566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx)); 278c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 279c4762a1bSJed Brown const PetscInt idx = pidx[p]; 280c4762a1bSJed Brown const PetscReal *c = &coords[idx * dim]; 281c4762a1bSJed Brown 282c4762a1bSJed Brown mom[0] += PetscRealPart(w[idx]); 283c4762a1bSJed Brown mom[1] += PetscRealPart(w[idx]) * c[0]; 284c4762a1bSJed Brown for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d]; 285c4762a1bSJed Brown } 2869566063dSJacob Faibussowitsch PetscCall(PetscFree(pidx)); 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2899566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 2909566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 2919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw))); 292c4762a1bSJed Brown PetscFunctionReturn(0); 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295*d71ae5a4SJacob 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[]) 296*d71ae5a4SJacob Faibussowitsch { 297c4762a1bSJed Brown f0[0] = u[0]; 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 300*d71ae5a4SJacob 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[]) 301*d71ae5a4SJacob Faibussowitsch { 302c4762a1bSJed Brown f0[0] = x[0] * u[0]; 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305*d71ae5a4SJacob 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[]) 306*d71ae5a4SJacob Faibussowitsch { 307c4762a1bSJed Brown PetscInt d; 308c4762a1bSJed Brown 309c4762a1bSJed Brown f0[0] = 0.0; 310c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0]; 311c4762a1bSJed Brown } 312c4762a1bSJed Brown 313*d71ae5a4SJacob Faibussowitsch static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user) 314*d71ae5a4SJacob Faibussowitsch { 315c4762a1bSJed Brown PetscDS prob; 316c4762a1bSJed Brown PetscScalar mom; 317c4762a1bSJed Brown 318c4762a1bSJed Brown PetscFunctionBeginUser; 3199566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 3209566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_1)); 3219566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 322c4762a1bSJed Brown moments[0] = PetscRealPart(mom); 3239566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_x)); 3249566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 325c4762a1bSJed Brown moments[1] = PetscRealPart(mom); 3269566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_r2)); 3279566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 328c4762a1bSJed Brown moments[2] = PetscRealPart(mom); 329c4762a1bSJed Brown PetscFunctionReturn(0); 330c4762a1bSJed Brown } 331c4762a1bSJed Brown 332*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user) 333*d71ae5a4SJacob Faibussowitsch { 334c4762a1bSJed Brown MPI_Comm comm; 335c4762a1bSJed Brown KSP ksp; 336c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 337c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 338c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */ 339c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 340c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 341c4762a1bSJed Brown PetscInt m; 342c4762a1bSJed Brown 343c4762a1bSJed Brown PetscFunctionBeginUser; 3449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3459566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 3469566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ptof_")); 3479566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 3489566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 349c4762a1bSJed Brown 3509566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 3519566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 352c4762a1bSJed Brown 353c4762a1bSJed Brown /* make particle weight vector */ 3549566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* create matrix RHS vector */ 3579566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rhs)); 3589566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 3609566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 361c4762a1bSJed Brown 3629566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 3639566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 3649566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 3659566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M, M)); 3669566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 3679566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, fhat)); 3689566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 3699566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 370c4762a1bSJed Brown 371c4762a1bSJed Brown /* Check moments of field */ 3729566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 3739566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 3749566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 375c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 3761dca8a05SBarry 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); 377c4762a1bSJed Brown } 378c4762a1bSJed Brown 3799566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 3809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 3819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 3829566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 3839566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 384c4762a1bSJed Brown 385c4762a1bSJed Brown PetscFunctionReturn(0); 386c4762a1bSJed Brown } 387c4762a1bSJed Brown 388*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user) 389*d71ae5a4SJacob Faibussowitsch { 390c4762a1bSJed Brown MPI_Comm comm; 391c4762a1bSJed Brown KSP ksp; 392c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 393832e2b15SMatthew G. Knepley Mat M_p, PM_p; /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */ 394c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */ 395c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 396c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 397c4762a1bSJed Brown PetscInt m; 398c4762a1bSJed Brown 399c4762a1bSJed Brown PetscFunctionBeginUser; 4009566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 401c4762a1bSJed Brown 4029566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 4039566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ftop_")); 4049566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 405c4762a1bSJed Brown 4069566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 4079566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 408c4762a1bSJed Brown 4099566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 4109566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 411c4762a1bSJed Brown 412c4762a1bSJed Brown /* make particle weight vector */ 4139566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */ 4169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 4179566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 4189566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 4199566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 4209566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 4219566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M, fhat, rhs)); 4229566063dSJacob Faibussowitsch if (user->useBlockDiagPrec) PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 4239371c9d4SSatish Balay else { 4249371c9d4SSatish Balay PetscCall(PetscObjectReference((PetscObject)M_p)); 4259371c9d4SSatish Balay PM_p = M_p; 4269371c9d4SSatish Balay } 427c4762a1bSJed Brown 4289566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M_p, PM_p)); 4299566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, rhs, f)); 4309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 4319566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 432c4762a1bSJed Brown 4339566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* Check moments */ 4369566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 4379566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 4389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 439c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 4401dca8a05SBarry 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); 441c4762a1bSJed Brown } 4429566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 4439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 4449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 4459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PM_p)); 4469566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 4479566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 448c4762a1bSJed Brown PetscFunctionReturn(0); 449c4762a1bSJed Brown } 450c4762a1bSJed Brown 451c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */ 452*d71ae5a4SJacob Faibussowitsch static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC) 453*d71ae5a4SJacob Faibussowitsch { 454c4762a1bSJed Brown DM_Plex *mesh = (DM_Plex *)dm->data; 455c4762a1bSJed Brown PetscInt debug = mesh->printFEM; 456c4762a1bSJed Brown DM dmC; 457c4762a1bSJed Brown PetscSection section; 458c4762a1bSJed Brown PetscQuadrature quad = NULL; 459c4762a1bSJed Brown PetscScalar *interpolant, *gradsum; 460c4762a1bSJed Brown PetscFEGeom fegeom; 461c4762a1bSJed Brown PetscReal *coords; 462c4762a1bSJed Brown const PetscReal *quadPoints, *quadWeights; 463c4762a1bSJed Brown PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset; 464c4762a1bSJed Brown 465c4762a1bSJed Brown PetscFunctionBegin; 4669566063dSJacob Faibussowitsch PetscCall(VecGetDM(locC, &dmC)); 4679566063dSJacob Faibussowitsch PetscCall(VecSet(locC, 0.0)); 4689566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4699566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &coordDim)); 470c4762a1bSJed Brown fegeom.dimEmbed = coordDim; 4719566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 4729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &numFields)); 473c4762a1bSJed Brown for (field = 0; field < numFields; ++field) { 474c4762a1bSJed Brown PetscObject obj; 475c4762a1bSJed Brown PetscClassId id; 476c4762a1bSJed Brown PetscInt Nc; 477c4762a1bSJed Brown 4789566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 4799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 480c4762a1bSJed Brown if (id == PETSCFE_CLASSID) { 481c4762a1bSJed Brown PetscFE fe = (PetscFE)obj; 482c4762a1bSJed Brown 4839566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 4849566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fe, &Nc)); 485c4762a1bSJed Brown } else if (id == PETSCFV_CLASSID) { 486c4762a1bSJed Brown PetscFV fv = (PetscFV)obj; 487c4762a1bSJed Brown 4889566063dSJacob Faibussowitsch PetscCall(PetscFVGetQuadrature(fv, &quad)); 4899566063dSJacob Faibussowitsch PetscCall(PetscFVGetNumComponents(fv, &Nc)); 49063a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 491c4762a1bSJed Brown numComponents += Nc; 492c4762a1bSJed Brown } 4939566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 49463a3b9bcSJacob Faibussowitsch PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents); 4959566063dSJacob 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)); 4969566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4979566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 498c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 499c4762a1bSJed Brown PetscInt *star = NULL; 500c4762a1bSJed Brown PetscInt starSize, st, d, fc; 501c4762a1bSJed Brown 5029566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gradsum, coordDim * numComponents)); 5039566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 504c4762a1bSJed Brown for (st = 0; st < starSize * 2; st += 2) { 505c4762a1bSJed Brown const PetscInt cell = star[st]; 506c4762a1bSJed Brown PetscScalar *grad = &gradsum[coordDim * numComponents]; 507c4762a1bSJed Brown PetscScalar *x = NULL; 508c4762a1bSJed Brown 509c4762a1bSJed Brown if ((cell < cStart) || (cell >= cEnd)) continue; 5109566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ)); 5119566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x)); 512c4762a1bSJed Brown for (field = 0, fieldOffset = 0; field < numFields; ++field) { 513c4762a1bSJed Brown PetscObject obj; 514c4762a1bSJed Brown PetscClassId id; 515c4762a1bSJed Brown PetscInt Nb, Nc, q, qc = 0; 516c4762a1bSJed Brown 5179566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(grad, coordDim * numComponents)); 5189566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 5199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 5209371c9d4SSatish Balay if (id == PETSCFE_CLASSID) { 5219371c9d4SSatish Balay PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 5229371c9d4SSatish Balay PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb)); 5239371c9d4SSatish Balay } else if (id == PETSCFV_CLASSID) { 5249371c9d4SSatish Balay PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 5259371c9d4SSatish Balay Nb = 1; 5269371c9d4SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 527c4762a1bSJed Brown for (q = 0; q < Nq; ++q) { 52863a3b9bcSJacob 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); 5299566063dSJacob Faibussowitsch if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &fegeom, q, interpolant)); 53063a3b9bcSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 531c4762a1bSJed Brown for (fc = 0; fc < Nc; ++fc) { 532c4762a1bSJed Brown const PetscReal wt = quadWeights[q * qNc + qc + fc]; 533c4762a1bSJed Brown 534c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q]; 535c4762a1bSJed Brown } 536c4762a1bSJed Brown } 537c4762a1bSJed Brown fieldOffset += Nb; 538c4762a1bSJed Brown } 5399566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x)); 540c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 541ad540459SPierre Jolivet for (d = 0; d < coordDim; ++d) gradsum[fc * coordDim + d] += grad[fc * coordDim + d]; 542c4762a1bSJed Brown } 543c4762a1bSJed Brown if (debug) { 54463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " gradient: [", cell)); 545c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 546c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) { 5479566063dSJacob Faibussowitsch if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 5489566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d]))); 549c4762a1bSJed Brown } 550c4762a1bSJed Brown } 5519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n")); 552c4762a1bSJed Brown } 553c4762a1bSJed Brown } 5549566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 5559566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES)); 556c4762a1bSJed Brown } 5579566063dSJacob Faibussowitsch PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ)); 558c4762a1bSJed Brown PetscFunctionReturn(0); 559c4762a1bSJed Brown } 560c4762a1bSJed Brown 561*d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user) 562*d71ae5a4SJacob Faibussowitsch { 563c4762a1bSJed Brown MPI_Comm comm; 564c4762a1bSJed Brown KSP ksp; 565c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 566c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 567c4762a1bSJed Brown Vec f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */ 568c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 569c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 570c4762a1bSJed Brown PetscInt m; 571c4762a1bSJed Brown 572c4762a1bSJed Brown PetscFunctionBeginUser; 5739566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 5749566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 5759566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ptof_")); 5769566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 5779566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 5789566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &grad)); 579c4762a1bSJed Brown 5809566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 5819566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 582c4762a1bSJed Brown 583c4762a1bSJed Brown /* make particle weight vector */ 5849566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 585c4762a1bSJed Brown 586c4762a1bSJed Brown /* create matrix RHS vector */ 5879566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rhs)); 5889566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 5899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 5909566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 591c4762a1bSJed Brown 5929566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 5939566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 594c4762a1bSJed Brown 5959566063dSJacob Faibussowitsch PetscCall(InterpolateGradient(dm, fhat, grad)); 596c4762a1bSJed Brown 5979566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 5989566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M, M)); 5999566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 6009566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, grad)); 6019566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 6029566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 603c4762a1bSJed Brown 604c4762a1bSJed Brown /* Check moments of field */ 6059566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 6069566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, grad, fmoments, user)); 6079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 608c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 6091dca8a05SBarry 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); 610c4762a1bSJed Brown } 611c4762a1bSJed Brown 6129566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 6139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 6149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 6159566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 6169566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 6179566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &grad)); 618c4762a1bSJed Brown 619c4762a1bSJed Brown PetscFunctionReturn(0); 620c4762a1bSJed Brown } 621c4762a1bSJed Brown 622*d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 623*d71ae5a4SJacob Faibussowitsch { 624c4762a1bSJed Brown MPI_Comm comm; 625c4762a1bSJed Brown DM dm, sw; 626c4762a1bSJed Brown AppCtx user; 627c4762a1bSJed Brown 628327415f7SBarry Smith PetscFunctionBeginUser; 6299566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 630c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 6319566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 6329566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user)); 6339566063dSJacob Faibussowitsch PetscCall(CreateFEM(dm, &user)); 6349566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user)); 6359566063dSJacob Faibussowitsch PetscCall(TestL2ProjectionParticlesToField(dm, sw, &user)); 6369566063dSJacob Faibussowitsch PetscCall(TestL2ProjectionFieldToParticles(dm, sw, &user)); 6379566063dSJacob Faibussowitsch PetscCall(TestFieldGradientProjection(dm, sw, &user)); 6389566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6399566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 6409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 641b122ec5aSJacob Faibussowitsch return 0; 642c4762a1bSJed Brown } 643c4762a1bSJed Brown 644c4762a1bSJed Brown /*TEST 645c4762a1bSJed Brown 646832e2b15SMatthew G. Knepley # Swarm does not handle complex or quad 647832e2b15SMatthew G. Knepley build: 648832e2b15SMatthew G. Knepley requires: !complex double 649832e2b15SMatthew G. Knepley 650c4762a1bSJed Brown test: 651c4762a1bSJed Brown suffix: proj_tri_0 652832e2b15SMatthew G. Knepley requires: triangle 653832e2b15SMatthew 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 654c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 655c4762a1bSJed Brown 656c4762a1bSJed Brown test: 657c4762a1bSJed Brown suffix: proj_tri_2_faces 658832e2b15SMatthew G. Knepley requires: triangle 659832e2b15SMatthew 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 660c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 661c4762a1bSJed Brown 662c4762a1bSJed Brown test: 663c4762a1bSJed Brown suffix: proj_quad_0 664832e2b15SMatthew G. Knepley requires: triangle 66530602db0SMatthew 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 666c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 667c4762a1bSJed Brown 668c4762a1bSJed Brown test: 669c4762a1bSJed Brown suffix: proj_quad_2_faces 670832e2b15SMatthew G. Knepley requires: triangle 67130602db0SMatthew 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 672c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 673c4762a1bSJed Brown 674c4762a1bSJed Brown test: 675c4762a1bSJed Brown suffix: proj_tri_5P 676832e2b15SMatthew G. Knepley requires: triangle 677832e2b15SMatthew 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 678c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 679c4762a1bSJed Brown 680c4762a1bSJed Brown test: 681c4762a1bSJed Brown suffix: proj_quad_5P 682832e2b15SMatthew G. Knepley requires: triangle 68330602db0SMatthew 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 684c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 685c4762a1bSJed Brown 686c4762a1bSJed Brown test: 687c4762a1bSJed Brown suffix: proj_tri_mdx 688832e2b15SMatthew G. Knepley requires: triangle 689832e2b15SMatthew 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 690c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 691c4762a1bSJed Brown 692c4762a1bSJed Brown test: 693c4762a1bSJed Brown suffix: proj_tri_mdx_5P 694832e2b15SMatthew G. Knepley requires: triangle 695832e2b15SMatthew 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 696c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 697c4762a1bSJed Brown 698c4762a1bSJed Brown test: 699c4762a1bSJed Brown suffix: proj_tri_3d 700832e2b15SMatthew G. Knepley requires: ctetgen 70130602db0SMatthew 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 702c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 703c4762a1bSJed Brown 704c4762a1bSJed Brown test: 705c4762a1bSJed Brown suffix: proj_tri_3d_2_faces 706832e2b15SMatthew G. Knepley requires: ctetgen 70730602db0SMatthew 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 708c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 709c4762a1bSJed Brown 710c4762a1bSJed Brown test: 711c4762a1bSJed Brown suffix: proj_tri_3d_5P 712832e2b15SMatthew G. Knepley requires: ctetgen 71330602db0SMatthew 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 714c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 715c4762a1bSJed Brown 716c4762a1bSJed Brown test: 717c4762a1bSJed Brown suffix: proj_tri_3d_mdx 718832e2b15SMatthew G. Knepley requires: ctetgen 71930602db0SMatthew 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 720c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 721c4762a1bSJed Brown 722c4762a1bSJed Brown test: 723c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P 724832e2b15SMatthew G. Knepley requires: ctetgen 72530602db0SMatthew 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 726c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 727c4762a1bSJed Brown 728c4762a1bSJed Brown test: 729c4762a1bSJed Brown suffix: proj_tri_3d_mdx_2_faces 730832e2b15SMatthew G. Knepley requires: ctetgen 73130602db0SMatthew 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 732c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 733c4762a1bSJed Brown 734c4762a1bSJed Brown test: 735c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P_2_faces 736832e2b15SMatthew G. Knepley requires: ctetgen 73730602db0SMatthew 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 738c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 739c4762a1bSJed Brown 740832e2b15SMatthew G. Knepley test: 741832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_scale 742832e2b15SMatthew G. Knepley requires: !complex 74330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 744832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 745832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 746832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 747832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none 748832e2b15SMatthew G. Knepley 749832e2b15SMatthew G. Knepley test: 750832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_prec_scale 751832e2b15SMatthew G. Knepley requires: !complex 75230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 753832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 754832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 755832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 756832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec 757832e2b15SMatthew G. Knepley 758c4762a1bSJed Brown TEST*/ 759