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*9371c9d4SSatish Balay static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) { 23c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 24c4762a1bSJed Brown PetscInt d; 25c4762a1bSJed Brown 26c4762a1bSJed Brown u[0] = 0.0; 27832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d] / (ctx->L[d]); 28c4762a1bSJed Brown return 0; 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 31*9371c9d4SSatish Balay static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) { 32c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 33c4762a1bSJed Brown PetscInt d; 34c4762a1bSJed Brown 35c4762a1bSJed Brown u[0] = 1; 36832e2b15SMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d]) * PetscSqr(ctx->L[d]) - PetscPowRealInt(x[d], 4); 37c4762a1bSJed Brown return 0; 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40*9371c9d4SSatish Balay static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *a_ctx) { 41c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)a_ctx; 42c4762a1bSJed Brown 43832e2b15SMatthew G. Knepley u[0] = PetscSinScalar(2 * PETSC_PI * ctx->k * x[0] / (ctx->L[0])); 44c4762a1bSJed Brown return 0; 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 48c4762a1bSJed Brown char fstring[PETSC_MAX_PATH_LEN] = "linear"; 49c4762a1bSJed Brown PetscBool flag; 50c4762a1bSJed Brown 51c4762a1bSJed Brown PetscFunctionBeginUser; 52c4762a1bSJed Brown options->particlesPerCell = 1; 53c4762a1bSJed Brown options->k = 1; 54c4762a1bSJed Brown options->particleRelDx = 1.e-20; 55c4762a1bSJed Brown options->meshRelDx = 1.e-20; 56c4762a1bSJed Brown options->momentTol = 100. * PETSC_MACHINE_EPSILON; 57832e2b15SMatthew G. Knepley options->useBlockDiagPrec = PETSC_FALSE; 58c4762a1bSJed Brown 59d0609cedSBarry Smith PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX"); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-k", "Mode number of test", "ex2.c", options->k, &options->k, NULL)); 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-particlesPerCell", "Number of particles per cell", "ex2.c", options->particlesPerCell, &options->particlesPerCell, NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-particle_perturbation", "Relative perturbation of particles (0,1)", "ex2.c", options->particleRelDx, &options->particleRelDx, NULL)); 639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mesh_perturbation", "Relative perturbation of mesh points (0,1)", "ex2.c", options->meshRelDx, &options->meshRelDx, NULL)); 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-function", "Name of test function", "ex2.c", fstring, fstring, sizeof(fstring), NULL)); 659566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "linear", &flag)); 66c4762a1bSJed Brown if (flag) { 67c4762a1bSJed Brown options->func = linear; 68c4762a1bSJed Brown } else { 699566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "sin", &flag)); 70c4762a1bSJed Brown if (flag) { 71c4762a1bSJed Brown options->func = sinx; 72c4762a1bSJed Brown } else { 739566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fstring, "x2_x4", &flag)); 74c4762a1bSJed Brown options->func = x2_x4; 7528b400f6SJacob Faibussowitsch PetscCheck(flag, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown function %s", fstring); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown } 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-moment_tol", "Tolerance for moment checks", "ex2.c", options->momentTol, &options->momentTol, NULL)); 799566063dSJacob 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)); 80d0609cedSBarry Smith PetscOptionsEnd(); 81c4762a1bSJed Brown 82c4762a1bSJed Brown PetscFunctionReturn(0); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85*9371c9d4SSatish Balay static PetscErrorCode PerturbVertices(DM dm, AppCtx *user) { 86c4762a1bSJed Brown PetscRandom rnd; 87c4762a1bSJed Brown PetscReal interval = user->meshRelDx; 88c4762a1bSJed Brown Vec coordinates; 89c4762a1bSJed Brown PetscScalar *coords; 90832e2b15SMatthew G. Knepley PetscReal *hh, low[3], high[3]; 91832e2b15SMatthew G. Knepley PetscInt d, cdim, cEnd, N, p, bs; 92c4762a1bSJed Brown 93c4762a1bSJed Brown PetscFunctionBeginUser; 949566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, low, high)); 959566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &cEnd)); 969566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 979566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -interval, interval)); 989566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 999566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1009566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 1019566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(cdim, &hh)); 102832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) hh[d] = (user->L[d]) / PetscPowReal(cEnd, 1. / cdim); 1039566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &N)); 1049566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(coordinates, &bs)); 10563a3b9bcSJacob Faibussowitsch PetscCheck(bs == cdim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Coordinate vector has wrong block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, cdim); 1069566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 107c4762a1bSJed Brown for (p = 0; p < N; p += cdim) { 108c4762a1bSJed Brown PetscScalar *coord = &coords[p], value; 109c4762a1bSJed Brown 110c4762a1bSJed Brown for (d = 0; d < cdim; ++d) { 1119566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rnd, &value)); 112832e2b15SMatthew G. Knepley coord[d] = PetscMax(low[d], PetscMin(high[d], PetscRealPart(coord[d] + value * hh[d]))); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 1159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1169566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 1179566063dSJacob Faibussowitsch PetscCall(PetscFree(hh)); 118c4762a1bSJed Brown PetscFunctionReturn(0); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user) { 122832e2b15SMatthew G. Knepley PetscReal low[3], high[3]; 123832e2b15SMatthew G. Knepley PetscInt cdim, d; 124c4762a1bSJed Brown 125c4762a1bSJed Brown PetscFunctionBeginUser; 1269566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1279566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1289566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1299566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 13030602db0SMatthew G. Knepley 1319566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*dm, &cdim)); 1329566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(*dm, low, high)); 133832e2b15SMatthew G. Knepley for (d = 0; d < cdim; ++d) user->L[d] = high[d] - low[d]; 1349566063dSJacob Faibussowitsch PetscCall(PerturbVertices(*dm, user)); 135c4762a1bSJed Brown PetscFunctionReturn(0); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138*9371c9d4SSatish Balay 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[]) { 139c4762a1bSJed Brown g0[0] = 1.0; 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142*9371c9d4SSatish Balay static PetscErrorCode CreateFEM(DM dm, AppCtx *user) { 143c4762a1bSJed Brown PetscFE fe; 144c4762a1bSJed Brown PetscDS ds; 145832e2b15SMatthew G. Knepley DMPolytopeType ct; 146832e2b15SMatthew G. Knepley PetscBool simplex; 147832e2b15SMatthew G. Knepley PetscInt dim, cStart; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBeginUser; 1509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1519566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1529566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 153832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1549566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, NULL, -1, &fe)); 1559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "fe")); 1569566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 1579566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1589566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 159c4762a1bSJed Brown /* Setup to form mass matrix */ 1609566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1619566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, identity, NULL, NULL, NULL)); 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165*9371c9d4SSatish Balay static PetscErrorCode CreateParticles(DM dm, DM *sw, AppCtx *user) { 166c4762a1bSJed Brown PetscRandom rnd, rndp; 167c4762a1bSJed Brown PetscReal interval = user->particleRelDx; 168832e2b15SMatthew G. Knepley DMPolytopeType ct; 169832e2b15SMatthew G. Knepley PetscBool simplex; 170c4762a1bSJed Brown PetscScalar value, *vals; 171c4762a1bSJed Brown PetscReal *centroid, *coords, *xi0, *v0, *J, *invJ, detJ; 172c4762a1bSJed Brown PetscInt *cellid; 173832e2b15SMatthew G. Knepley PetscInt Ncell, Np = user->particlesPerCell, p, cStart, c, dim, d; 174c4762a1bSJed Brown 175c4762a1bSJed Brown PetscFunctionBeginUser; 1769566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1779566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 1789566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 179832e2b15SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 1809566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw)); 1819566063dSJacob Faibussowitsch PetscCall(DMSetType(*sw, DMSWARM)); 1829566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*sw, dim)); 183c4762a1bSJed Brown 1849566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rnd)); 1859566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, -1.0, 1.0)); 1869566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd)); 1879566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm), &rndp)); 1889566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rndp, -interval, interval)); 1899566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rndp)); 190c4762a1bSJed Brown 1919566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); 1929566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(*sw, dm)); 1939566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR)); 1949566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(*sw)); 1959566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, NULL, &Ncell)); 1969566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(*sw, Ncell * Np, 0)); 1979566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*sw)); 1989566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 1999566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2009566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&vals)); 201c4762a1bSJed Brown 2029566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dim, ¢roid, dim, &xi0, dim, &v0, dim * dim, &J, dim * dim, &invJ)); 203c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 204c4762a1bSJed Brown if (Np == 1) { 2059566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 206c4762a1bSJed Brown cellid[c] = c; 207c4762a1bSJed Brown for (d = 0; d < dim; ++d) coords[c * dim + d] = centroid[d]; 208c4762a1bSJed Brown } else { 2099566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); /* affine */ 210c4762a1bSJed Brown for (d = 0; d < dim; ++d) xi0[d] = -1.0; 211c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 212c4762a1bSJed Brown const PetscInt n = c * Np + p; 213c4762a1bSJed Brown PetscReal sum = 0.0, refcoords[3]; 214c4762a1bSJed Brown 215c4762a1bSJed Brown cellid[n] = c; 216*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 217*9371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rnd, &value)); 218*9371c9d4SSatish Balay refcoords[d] = PetscRealPart(value); 219*9371c9d4SSatish Balay sum += refcoords[d]; 220*9371c9d4SSatish Balay } 221*9371c9d4SSatish Balay if (simplex && sum > 0.0) 222*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) refcoords[d] -= PetscSqrtReal(dim) * sum; 223c4762a1bSJed Brown CoordinatesRefToReal(dim, dim, xi0, v0, J, refcoords, &coords[n * dim]); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown } 226c4762a1bSJed Brown } 2279566063dSJacob Faibussowitsch PetscCall(PetscFree5(centroid, xi0, v0, J, invJ)); 228c4762a1bSJed Brown for (c = 0; c < Ncell; ++c) { 229c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 230c4762a1bSJed Brown const PetscInt n = c * Np + p; 231c4762a1bSJed Brown 232*9371c9d4SSatish Balay for (d = 0; d < dim; ++d) { 233*9371c9d4SSatish Balay PetscCall(PetscRandomGetValue(rndp, &value)); 234*9371c9d4SSatish Balay coords[n * dim + d] += PetscRealPart(value); 235*9371c9d4SSatish Balay } 236c4762a1bSJed Brown user->func(dim, 0.0, &coords[n * dim], 1, &vals[c], user); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown } 2399566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2409566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_cellid, NULL, NULL, (void **)&cellid)); 2419566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&vals)); 2429566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 2439566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rndp)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles")); 2459566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view")); 246c4762a1bSJed Brown PetscFunctionReturn(0); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249*9371c9d4SSatish Balay static PetscErrorCode computeParticleMoments(DM sw, PetscReal moments[3], AppCtx *user) { 250c4762a1bSJed Brown DM dm; 251c4762a1bSJed Brown const PetscReal *coords; 252c4762a1bSJed Brown const PetscScalar *w; 253c4762a1bSJed Brown PetscReal mom[3] = {0.0, 0.0, 0.0}; 254c4762a1bSJed Brown PetscInt cell, cStart, cEnd, dim; 255c4762a1bSJed Brown 256c4762a1bSJed Brown PetscFunctionBeginUser; 2579566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sw, &dim)); 2589566063dSJacob Faibussowitsch PetscCall(DMSwarmGetCellDM(sw, &dm)); 2599566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2609566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetAccess(sw)); 2619566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2629566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&w)); 263c4762a1bSJed Brown for (cell = cStart; cell < cEnd; ++cell) { 264c4762a1bSJed Brown PetscInt *pidx; 265c4762a1bSJed Brown PetscInt Np, p, d; 266c4762a1bSJed Brown 2679566063dSJacob Faibussowitsch PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx)); 268c4762a1bSJed Brown for (p = 0; p < Np; ++p) { 269c4762a1bSJed Brown const PetscInt idx = pidx[p]; 270c4762a1bSJed Brown const PetscReal *c = &coords[idx * dim]; 271c4762a1bSJed Brown 272c4762a1bSJed Brown mom[0] += PetscRealPart(w[idx]); 273c4762a1bSJed Brown mom[1] += PetscRealPart(w[idx]) * c[0]; 274c4762a1bSJed Brown for (d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx]) * c[d] * c[d]; 275c4762a1bSJed Brown } 2769566063dSJacob Faibussowitsch PetscCall(PetscFree(pidx)); 277c4762a1bSJed Brown } 2789566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords)); 2799566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&w)); 2809566063dSJacob Faibussowitsch PetscCall(DMSwarmSortRestoreAccess(sw)); 2819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw))); 282c4762a1bSJed Brown PetscFunctionReturn(0); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285*9371c9d4SSatish Balay 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[]) { 286c4762a1bSJed Brown f0[0] = u[0]; 287c4762a1bSJed Brown } 288c4762a1bSJed Brown 289*9371c9d4SSatish Balay 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[]) { 290c4762a1bSJed Brown f0[0] = x[0] * u[0]; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 293*9371c9d4SSatish Balay 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[]) { 294c4762a1bSJed Brown PetscInt d; 295c4762a1bSJed Brown 296c4762a1bSJed Brown f0[0] = 0.0; 297c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[0]; 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 300*9371c9d4SSatish Balay static PetscErrorCode computeFEMMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user) { 301c4762a1bSJed Brown PetscDS prob; 302c4762a1bSJed Brown PetscScalar mom; 303c4762a1bSJed Brown 304c4762a1bSJed Brown PetscFunctionBeginUser; 3059566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 3069566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_1)); 3079566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 308c4762a1bSJed Brown moments[0] = PetscRealPart(mom); 3099566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_x)); 3109566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 311c4762a1bSJed Brown moments[1] = PetscRealPart(mom); 3129566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_r2)); 3139566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user)); 314c4762a1bSJed Brown moments[2] = PetscRealPart(mom); 315c4762a1bSJed Brown PetscFunctionReturn(0); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318*9371c9d4SSatish Balay static PetscErrorCode TestL2ProjectionParticlesToField(DM dm, DM sw, AppCtx *user) { 319c4762a1bSJed Brown MPI_Comm comm; 320c4762a1bSJed Brown KSP ksp; 321c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 322c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 323c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i f, FEM field */ 324c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 325c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 326c4762a1bSJed Brown PetscInt m; 327c4762a1bSJed Brown 328c4762a1bSJed Brown PetscFunctionBeginUser; 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3309566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 3319566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ptof_")); 3329566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 3339566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 334c4762a1bSJed Brown 3359566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 3369566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* make particle weight vector */ 3399566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* create matrix RHS vector */ 3429566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rhs)); 3439566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 3449566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 3459566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 346c4762a1bSJed Brown 3479566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 3489566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 3499566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 3509566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M, M)); 3519566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 3529566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, fhat)); 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 3549566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* Check moments of field */ 3579566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 3589566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 3599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 360c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 3611dca8a05SBarry 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); 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 3649566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 3659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 3669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 3679566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 3689566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 369c4762a1bSJed Brown 370c4762a1bSJed Brown PetscFunctionReturn(0); 371c4762a1bSJed Brown } 372c4762a1bSJed Brown 373*9371c9d4SSatish Balay static PetscErrorCode TestL2ProjectionFieldToParticles(DM dm, DM sw, AppCtx *user) { 374c4762a1bSJed Brown MPI_Comm comm; 375c4762a1bSJed Brown KSP ksp; 376c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 377832e2b15SMatthew G. Knepley Mat M_p, PM_p; /* Particle mass matrix M_p, and the preconditioning matrix, e.g. M_p M^T_p */ 378c4762a1bSJed Brown Vec f, rhs, fhat; /* Particle field f, \int phi_i fhat, FEM field */ 379c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 380c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 381c4762a1bSJed Brown PetscInt m; 382c4762a1bSJed Brown 383c4762a1bSJed Brown PetscFunctionBeginUser; 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 385c4762a1bSJed Brown 3869566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 3879566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ftop_")); 3889566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 389c4762a1bSJed Brown 3909566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 3919566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 392c4762a1bSJed Brown 3939566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 3949566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 395c4762a1bSJed Brown 396c4762a1bSJed Brown /* make particle weight vector */ 3979566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 398c4762a1bSJed Brown 399c4762a1bSJed Brown /* create matrix RHS vector, in this case the FEM field fhat with the coefficients vector #alpha */ 4009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 4019566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 4029566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 4039566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 4049566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 4059566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M, fhat, rhs)); 4069566063dSJacob Faibussowitsch if (user->useBlockDiagPrec) PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 407*9371c9d4SSatish Balay else { 408*9371c9d4SSatish Balay PetscCall(PetscObjectReference((PetscObject)M_p)); 409*9371c9d4SSatish Balay PM_p = M_p; 410*9371c9d4SSatish Balay } 411c4762a1bSJed Brown 4129566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M_p, PM_p)); 4139566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, rhs, f)); 4149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 4159566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 416c4762a1bSJed Brown 4179566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 418c4762a1bSJed Brown 419c4762a1bSJed Brown /* Check moments */ 4209566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 4219566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, fhat, fmoments, user)); 4229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 423c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 4241dca8a05SBarry 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); 425c4762a1bSJed Brown } 4269566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 4279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 4289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 4299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PM_p)); 4309566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 4319566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 432c4762a1bSJed Brown PetscFunctionReturn(0); 433c4762a1bSJed Brown } 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* Interpolate the gradient of an FEM (FVM) field. Code repurposed from DMPlexComputeGradientClementInterpolant */ 436*9371c9d4SSatish Balay static PetscErrorCode InterpolateGradient(DM dm, Vec locX, Vec locC) { 437c4762a1bSJed Brown DM_Plex *mesh = (DM_Plex *)dm->data; 438c4762a1bSJed Brown PetscInt debug = mesh->printFEM; 439c4762a1bSJed Brown DM dmC; 440c4762a1bSJed Brown PetscSection section; 441c4762a1bSJed Brown PetscQuadrature quad = NULL; 442c4762a1bSJed Brown PetscScalar *interpolant, *gradsum; 443c4762a1bSJed Brown PetscFEGeom fegeom; 444c4762a1bSJed Brown PetscReal *coords; 445c4762a1bSJed Brown const PetscReal *quadPoints, *quadWeights; 446c4762a1bSJed Brown PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, vStart, vEnd, v, field, fieldOffset; 447c4762a1bSJed Brown 448c4762a1bSJed Brown PetscFunctionBegin; 4499566063dSJacob Faibussowitsch PetscCall(VecGetDM(locC, &dmC)); 4509566063dSJacob Faibussowitsch PetscCall(VecSet(locC, 0.0)); 4519566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4529566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &coordDim)); 453c4762a1bSJed Brown fegeom.dimEmbed = coordDim; 4549566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 4559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &numFields)); 456c4762a1bSJed Brown for (field = 0; field < numFields; ++field) { 457c4762a1bSJed Brown PetscObject obj; 458c4762a1bSJed Brown PetscClassId id; 459c4762a1bSJed Brown PetscInt Nc; 460c4762a1bSJed Brown 4619566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 4629566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 463c4762a1bSJed Brown if (id == PETSCFE_CLASSID) { 464c4762a1bSJed Brown PetscFE fe = (PetscFE)obj; 465c4762a1bSJed Brown 4669566063dSJacob Faibussowitsch PetscCall(PetscFEGetQuadrature(fe, &quad)); 4679566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(fe, &Nc)); 468c4762a1bSJed Brown } else if (id == PETSCFV_CLASSID) { 469c4762a1bSJed Brown PetscFV fv = (PetscFV)obj; 470c4762a1bSJed Brown 4719566063dSJacob Faibussowitsch PetscCall(PetscFVGetQuadrature(fv, &quad)); 4729566063dSJacob Faibussowitsch PetscCall(PetscFVGetNumComponents(fv, &Nc)); 47363a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 474c4762a1bSJed Brown numComponents += Nc; 475c4762a1bSJed Brown } 4769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 47763a3b9bcSJacob Faibussowitsch PetscCheck(!(qNc != 1) || !(qNc != numComponents), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Quadrature components %" PetscInt_FMT " != %" PetscInt_FMT " field components", qNc, numComponents); 4789566063dSJacob 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)); 4799566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4809566063dSJacob Faibussowitsch PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd)); 481c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 482c4762a1bSJed Brown PetscInt *star = NULL; 483c4762a1bSJed Brown PetscInt starSize, st, d, fc; 484c4762a1bSJed Brown 4859566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gradsum, coordDim * numComponents)); 4869566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 487c4762a1bSJed Brown for (st = 0; st < starSize * 2; st += 2) { 488c4762a1bSJed Brown const PetscInt cell = star[st]; 489c4762a1bSJed Brown PetscScalar *grad = &gradsum[coordDim * numComponents]; 490c4762a1bSJed Brown PetscScalar *x = NULL; 491c4762a1bSJed Brown 492c4762a1bSJed Brown if ((cell < cStart) || (cell >= cEnd)) continue; 4939566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ)); 4949566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, locX, cell, NULL, &x)); 495c4762a1bSJed Brown for (field = 0, fieldOffset = 0; field < numFields; ++field) { 496c4762a1bSJed Brown PetscObject obj; 497c4762a1bSJed Brown PetscClassId id; 498c4762a1bSJed Brown PetscInt Nb, Nc, q, qc = 0; 499c4762a1bSJed Brown 5009566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(grad, coordDim * numComponents)); 5019566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, field, NULL, &obj)); 5029566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 503*9371c9d4SSatish Balay if (id == PETSCFE_CLASSID) { 504*9371c9d4SSatish Balay PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 505*9371c9d4SSatish Balay PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb)); 506*9371c9d4SSatish Balay } else if (id == PETSCFV_CLASSID) { 507*9371c9d4SSatish Balay PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 508*9371c9d4SSatish Balay Nb = 1; 509*9371c9d4SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 510c4762a1bSJed Brown for (q = 0; q < Nq; ++q) { 51163a3b9bcSJacob 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); 5129566063dSJacob Faibussowitsch if (id == PETSCFE_CLASSID) PetscCall(PetscFEInterpolateGradient_Static((PetscFE)obj, 1, &x[fieldOffset], &fegeom, q, interpolant)); 51363a3b9bcSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %" PetscInt_FMT, field); 514c4762a1bSJed Brown for (fc = 0; fc < Nc; ++fc) { 515c4762a1bSJed Brown const PetscReal wt = quadWeights[q * qNc + qc + fc]; 516c4762a1bSJed Brown 517c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) grad[fc * coordDim + d] += interpolant[fc * dim + d] * wt * fegeom.detJ[q]; 518c4762a1bSJed Brown } 519c4762a1bSJed Brown } 520c4762a1bSJed Brown fieldOffset += Nb; 521c4762a1bSJed Brown } 5229566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, locX, cell, NULL, &x)); 523c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 524*9371c9d4SSatish Balay for (d = 0; d < coordDim; ++d) { gradsum[fc * coordDim + d] += grad[fc * coordDim + d]; } 525c4762a1bSJed Brown } 526c4762a1bSJed Brown if (debug) { 52763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " gradient: [", cell)); 528c4762a1bSJed Brown for (fc = 0; fc < numComponents; ++fc) { 529c4762a1bSJed Brown for (d = 0; d < coordDim; ++d) { 5309566063dSJacob Faibussowitsch if (fc || d > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 5319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(grad[fc * coordDim + d]))); 532c4762a1bSJed Brown } 533c4762a1bSJed Brown } 5349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "]\n")); 535c4762a1bSJed Brown } 536c4762a1bSJed Brown } 5379566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, v, PETSC_FALSE, &starSize, &star)); 5389566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmC, NULL, locC, v, gradsum, INSERT_VALUES)); 539c4762a1bSJed Brown } 5409566063dSJacob Faibussowitsch PetscCall(PetscFree6(gradsum, interpolant, coords, fegeom.detJ, fegeom.J, fegeom.invJ)); 541c4762a1bSJed Brown PetscFunctionReturn(0); 542c4762a1bSJed Brown } 543c4762a1bSJed Brown 544*9371c9d4SSatish Balay static PetscErrorCode TestFieldGradientProjection(DM dm, DM sw, AppCtx *user) { 545c4762a1bSJed Brown MPI_Comm comm; 546c4762a1bSJed Brown KSP ksp; 547c4762a1bSJed Brown Mat M; /* FEM mass matrix */ 548c4762a1bSJed Brown Mat M_p; /* Particle mass matrix */ 549c4762a1bSJed Brown Vec f, rhs, fhat, grad; /* Particle field f, \int phi_i f, FEM field */ 550c4762a1bSJed Brown PetscReal pmoments[3]; /* \int f, \int x f, \int r^2 f */ 551c4762a1bSJed Brown PetscReal fmoments[3]; /* \int \hat f, \int x \hat f, \int r^2 \hat f */ 552c4762a1bSJed Brown PetscInt m; 553c4762a1bSJed Brown 554c4762a1bSJed Brown PetscFunctionBeginUser; 5559566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 5569566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 5579566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ptof_")); 5589566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &fhat)); 5599566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &rhs)); 5609566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &grad)); 561c4762a1bSJed Brown 5629566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 5639566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M_p, NULL, "-M_p_view")); 564c4762a1bSJed Brown 565c4762a1bSJed Brown /* make particle weight vector */ 5669566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 567c4762a1bSJed Brown 568c4762a1bSJed Brown /* create matrix RHS vector */ 5699566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rhs)); 5709566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 5719566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rhs, "rhs")); 5729566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view")); 573c4762a1bSJed Brown 5749566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 5759566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm, fhat, M, M, user)); 576c4762a1bSJed Brown 5779566063dSJacob Faibussowitsch PetscCall(InterpolateGradient(dm, fhat, grad)); 578c4762a1bSJed Brown 5799566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-M_view")); 5809566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M, M)); 5819566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 5829566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, grad)); 5839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fhat, "fhat")); 5849566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view")); 585c4762a1bSJed Brown 586c4762a1bSJed Brown /* Check moments of field */ 5879566063dSJacob Faibussowitsch PetscCall(computeParticleMoments(sw, pmoments, user)); 5889566063dSJacob Faibussowitsch PetscCall(computeFEMMoments(dm, grad, fmoments, user)); 5899566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2])); 590c4762a1bSJed Brown for (m = 0; m < 3; ++m) { 5911dca8a05SBarry 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); 592c4762a1bSJed Brown } 593c4762a1bSJed Brown 5949566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 5959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 5969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 5979566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &fhat)); 5989566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &rhs)); 5999566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &grad)); 600c4762a1bSJed Brown 601c4762a1bSJed Brown PetscFunctionReturn(0); 602c4762a1bSJed Brown } 603c4762a1bSJed Brown 604*9371c9d4SSatish Balay int main(int argc, char *argv[]) { 605c4762a1bSJed Brown MPI_Comm comm; 606c4762a1bSJed Brown DM dm, sw; 607c4762a1bSJed Brown AppCtx user; 608c4762a1bSJed Brown 609327415f7SBarry Smith PetscFunctionBeginUser; 6109566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 611c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 6129566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 6139566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm, &user)); 6149566063dSJacob Faibussowitsch PetscCall(CreateFEM(dm, &user)); 6159566063dSJacob Faibussowitsch PetscCall(CreateParticles(dm, &sw, &user)); 6169566063dSJacob Faibussowitsch PetscCall(TestL2ProjectionParticlesToField(dm, sw, &user)); 6179566063dSJacob Faibussowitsch PetscCall(TestL2ProjectionFieldToParticles(dm, sw, &user)); 6189566063dSJacob Faibussowitsch PetscCall(TestFieldGradientProjection(dm, sw, &user)); 6199566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 6219566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 622b122ec5aSJacob Faibussowitsch return 0; 623c4762a1bSJed Brown } 624c4762a1bSJed Brown 625c4762a1bSJed Brown /*TEST 626c4762a1bSJed Brown 627832e2b15SMatthew G. Knepley # Swarm does not handle complex or quad 628832e2b15SMatthew G. Knepley build: 629832e2b15SMatthew G. Knepley requires: !complex double 630832e2b15SMatthew G. Knepley 631c4762a1bSJed Brown test: 632c4762a1bSJed Brown suffix: proj_tri_0 633832e2b15SMatthew G. Knepley requires: triangle 634832e2b15SMatthew 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 635c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 636c4762a1bSJed Brown 637c4762a1bSJed Brown test: 638c4762a1bSJed Brown suffix: proj_tri_2_faces 639832e2b15SMatthew G. Knepley requires: triangle 640832e2b15SMatthew 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 641c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 642c4762a1bSJed Brown 643c4762a1bSJed Brown test: 644c4762a1bSJed Brown suffix: proj_quad_0 645832e2b15SMatthew G. Knepley requires: triangle 64630602db0SMatthew 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 647c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 648c4762a1bSJed Brown 649c4762a1bSJed Brown test: 650c4762a1bSJed Brown suffix: proj_quad_2_faces 651832e2b15SMatthew G. Knepley requires: triangle 65230602db0SMatthew 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 653c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 654c4762a1bSJed Brown 655c4762a1bSJed Brown test: 656c4762a1bSJed Brown suffix: proj_tri_5P 657832e2b15SMatthew G. Knepley requires: triangle 658832e2b15SMatthew 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 659c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 660c4762a1bSJed Brown 661c4762a1bSJed Brown test: 662c4762a1bSJed Brown suffix: proj_quad_5P 663832e2b15SMatthew G. Knepley requires: triangle 66430602db0SMatthew 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 665c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 666c4762a1bSJed Brown 667c4762a1bSJed Brown test: 668c4762a1bSJed Brown suffix: proj_tri_mdx 669832e2b15SMatthew G. Knepley requires: triangle 670832e2b15SMatthew 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 671c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 672c4762a1bSJed Brown 673c4762a1bSJed Brown test: 674c4762a1bSJed Brown suffix: proj_tri_mdx_5P 675832e2b15SMatthew G. Knepley requires: triangle 676832e2b15SMatthew 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 677c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 678c4762a1bSJed Brown 679c4762a1bSJed Brown test: 680c4762a1bSJed Brown suffix: proj_tri_3d 681832e2b15SMatthew G. Knepley requires: ctetgen 68230602db0SMatthew 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 683c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 684c4762a1bSJed Brown 685c4762a1bSJed Brown test: 686c4762a1bSJed Brown suffix: proj_tri_3d_2_faces 687832e2b15SMatthew G. Knepley requires: ctetgen 68830602db0SMatthew 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 689c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 690c4762a1bSJed Brown 691c4762a1bSJed Brown test: 692c4762a1bSJed Brown suffix: proj_tri_3d_5P 693832e2b15SMatthew G. Knepley requires: ctetgen 69430602db0SMatthew 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 695c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 696c4762a1bSJed Brown 697c4762a1bSJed Brown test: 698c4762a1bSJed Brown suffix: proj_tri_3d_mdx 699832e2b15SMatthew G. Knepley requires: ctetgen 70030602db0SMatthew 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 701c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 702c4762a1bSJed Brown 703c4762a1bSJed Brown test: 704c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P 705832e2b15SMatthew G. Knepley requires: ctetgen 70630602db0SMatthew 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 707c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 708c4762a1bSJed Brown 709c4762a1bSJed Brown test: 710c4762a1bSJed Brown suffix: proj_tri_3d_mdx_2_faces 711832e2b15SMatthew G. Knepley requires: ctetgen 71230602db0SMatthew 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 713c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 714c4762a1bSJed Brown 715c4762a1bSJed Brown test: 716c4762a1bSJed Brown suffix: proj_tri_3d_mdx_5P_2_faces 717832e2b15SMatthew G. Knepley requires: ctetgen 71830602db0SMatthew 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 719c4762a1bSJed Brown filter: grep -v marker | grep -v atomic | grep -v usage 720c4762a1bSJed Brown 721832e2b15SMatthew G. Knepley test: 722832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_scale 723832e2b15SMatthew G. Knepley requires: !complex 72430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 725832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 726832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 727832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 728832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type none 729832e2b15SMatthew G. Knepley 730832e2b15SMatthew G. Knepley test: 731832e2b15SMatthew G. Knepley suffix: proj_quad_lsqr_prec_scale 732832e2b15SMatthew G. Knepley requires: !complex 73330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 \ 734832e2b15SMatthew G. Knepley -petscspace_degree 2 -petscfe_default_quadrature_order 3 \ 735832e2b15SMatthew G. Knepley -particlesPerCell 200 \ 736832e2b15SMatthew G. Knepley -ptof_pc_type lu \ 737832e2b15SMatthew G. Knepley -ftop_ksp_rtol 1e-17 -ftop_ksp_type lsqr -ftop_pc_type lu -ftop_pc_factor_shift_type nonzero -block_diag_prec 738832e2b15SMatthew G. Knepley 739c4762a1bSJed Brown TEST*/ 740