xref: /petsc/src/dm/impls/swarm/tests/ex11.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
10ee167adSMatthew G. Knepley static char help[] = "Tests multifield and multicomponent L2 projection.\n";
20ee167adSMatthew G. Knepley 
30ee167adSMatthew G. Knepley #include <petscdmswarm.h>
40ee167adSMatthew G. Knepley #include <petscksp.h>
50ee167adSMatthew G. Knepley #include <petscdmplex.h>
60ee167adSMatthew G. Knepley #include <petscds.h>
70ee167adSMatthew G. Knepley 
80ee167adSMatthew G. Knepley typedef struct {
91898fd5cSMatthew G. Knepley   PetscBool           grad;  // Test gradient projection
100ee167adSMatthew G. Knepley   PetscBool           pass;  // Don't fail when moments are wrong
110ee167adSMatthew G. Knepley   PetscBool           fv;    // Use an FV discretization, instead of FE
120ee167adSMatthew G. Knepley   PetscInt            Npc;   // The number of partices per cell
130ee167adSMatthew G. Knepley   PetscInt            field; // The field to project
140ee167adSMatthew G. Knepley   PetscInt            Nm;    // The number of moments to match
150ee167adSMatthew G. Knepley   PetscReal           mtol;  // Tolerance for checking moment conservation
160ee167adSMatthew G. Knepley   PetscSimplePointFn *func;  // Function used to set particle weights
170ee167adSMatthew G. Knepley } AppCtx;
180ee167adSMatthew G. Knepley 
190ee167adSMatthew G. Knepley typedef enum {
200ee167adSMatthew G. Knepley   FUNCTION_CONSTANT,
210ee167adSMatthew G. Knepley   FUNCTION_LINEAR,
220ee167adSMatthew G. Knepley   FUNCTION_SIN,
230ee167adSMatthew G. Knepley   FUNCTION_X2_X4,
240ee167adSMatthew G. Knepley   FUNCTION_UNKNOWN,
250ee167adSMatthew G. Knepley   NUM_FUNCTIONS
260ee167adSMatthew G. Knepley } FunctionType;
270ee167adSMatthew G. Knepley const char *const FunctionTypes[] = {"constant", "linear", "sin", "x2_x4", "unknown", "FunctionType", "FUNCTION_", NULL};
280ee167adSMatthew G. Knepley 
constant(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)29*2a8381b2SBarry Smith static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
300ee167adSMatthew G. Knepley {
310ee167adSMatthew G. Knepley   u[0] = 1.0;
320ee167adSMatthew G. Knepley   return PETSC_SUCCESS;
330ee167adSMatthew G. Knepley }
340ee167adSMatthew G. Knepley 
linear(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)35*2a8381b2SBarry Smith static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
360ee167adSMatthew G. Knepley {
370ee167adSMatthew G. Knepley   u[0] = 0.0;
380ee167adSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) u[0] += x[d];
390ee167adSMatthew G. Knepley   return PETSC_SUCCESS;
400ee167adSMatthew G. Knepley }
410ee167adSMatthew G. Knepley 
sinx(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)42*2a8381b2SBarry Smith static PetscErrorCode sinx(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
430ee167adSMatthew G. Knepley {
440ee167adSMatthew G. Knepley   u[0] = 1;
450ee167adSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) u[0] *= PetscSinReal(2. * PETSC_PI * x[d]);
460ee167adSMatthew G. Knepley   return PETSC_SUCCESS;
470ee167adSMatthew G. Knepley }
480ee167adSMatthew G. Knepley 
x2_x4(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)49*2a8381b2SBarry Smith static PetscErrorCode x2_x4(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
500ee167adSMatthew G. Knepley {
510ee167adSMatthew G. Knepley   u[0] = 1;
520ee167adSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) u[0] *= PetscSqr(x[d]) - PetscSqr(PetscSqr(x[d]));
530ee167adSMatthew G. Knepley   return PETSC_SUCCESS;
540ee167adSMatthew G. Knepley }
550ee167adSMatthew G. Knepley 
ProcessOptions(MPI_Comm comm,AppCtx * options)560ee167adSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
570ee167adSMatthew G. Knepley {
580ee167adSMatthew G. Knepley   FunctionType func = FUNCTION_LINEAR;
590ee167adSMatthew G. Knepley   PetscBool    flg;
600ee167adSMatthew G. Knepley 
610ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
621898fd5cSMatthew G. Knepley   options->grad  = PETSC_FALSE;
630ee167adSMatthew G. Knepley   options->pass  = PETSC_FALSE;
640ee167adSMatthew G. Knepley   options->fv    = PETSC_FALSE;
650ee167adSMatthew G. Knepley   options->Npc   = 1;
660ee167adSMatthew G. Knepley   options->field = 0;
670ee167adSMatthew G. Knepley   options->Nm    = 1;
680ee167adSMatthew G. Knepley   options->mtol  = 100. * PETSC_MACHINE_EPSILON;
690ee167adSMatthew G. Knepley 
700ee167adSMatthew G. Knepley   PetscOptionsBegin(comm, "", "L2 Projection Options", "DMPLEX");
711898fd5cSMatthew G. Knepley   PetscCall(PetscOptionsBool("-grad", "Test gradient projection", __FILE__, options->grad, &options->grad, NULL));
720ee167adSMatthew G. Knepley   PetscCall(PetscOptionsBool("-pass", "Don't fail when moments are wrong", __FILE__, options->pass, &options->pass, NULL));
730ee167adSMatthew G. Knepley   PetscCall(PetscOptionsBool("-fv", "Use FV instead of FE", __FILE__, options->fv, &options->fv, NULL));
741898fd5cSMatthew G. Knepley   PetscCall(PetscOptionsBoundedInt("-npc", "Number of particles per cell", __FILE__, options->Npc, &options->Npc, NULL, 0));
750ee167adSMatthew G. Knepley   PetscCall(PetscOptionsBoundedInt("-field", "The field to project", __FILE__, options->field, &options->field, NULL, 0));
760ee167adSMatthew G. Knepley   PetscCall(PetscOptionsBoundedInt("-moments", "Number of moments to match", __FILE__, options->Nm, &options->Nm, NULL, 0));
770ee167adSMatthew G. Knepley   PetscCheck(options->Nm < 4, comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot match %" PetscInt_FMT " > 3 moments", options->Nm);
780ee167adSMatthew G. Knepley   PetscCall(PetscOptionsReal("-mtol", "Tolerance for moment checks", "ex2.c", options->mtol, &options->mtol, NULL));
790ee167adSMatthew G. Knepley   PetscCall(PetscOptionsEnum("-func", "Type of particle weight function", __FILE__, FunctionTypes, (PetscEnum)func, (PetscEnum *)&func, &flg));
800ee167adSMatthew G. Knepley   switch (func) {
810ee167adSMatthew G. Knepley   case FUNCTION_CONSTANT:
820ee167adSMatthew G. Knepley     options->func = constant;
830ee167adSMatthew G. Knepley     break;
840ee167adSMatthew G. Knepley   case FUNCTION_LINEAR:
850ee167adSMatthew G. Knepley     options->func = linear;
860ee167adSMatthew G. Knepley     break;
870ee167adSMatthew G. Knepley   case FUNCTION_SIN:
880ee167adSMatthew G. Knepley     options->func = sinx;
890ee167adSMatthew G. Knepley     break;
900ee167adSMatthew G. Knepley   case FUNCTION_X2_X4:
910ee167adSMatthew G. Knepley     options->func = x2_x4;
920ee167adSMatthew G. Knepley     break;
930ee167adSMatthew G. Knepley   default:
94d8b4a066SPierre Jolivet     PetscCheck(flg, comm, PETSC_ERR_ARG_WRONG, "Cannot handle function \"%s\"", FunctionTypes[func]);
950ee167adSMatthew G. Knepley   }
960ee167adSMatthew G. Knepley   PetscOptionsEnd();
970ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
980ee167adSMatthew G. Knepley }
990ee167adSMatthew G. Knepley 
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * user)1000ee167adSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *user)
1010ee167adSMatthew G. Knepley {
1020ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
1030ee167adSMatthew G. Knepley   PetscCall(DMCreate(comm, dm));
1040ee167adSMatthew G. Knepley   PetscCall(DMSetType(*dm, DMPLEX));
1050ee167adSMatthew G. Knepley   PetscCall(DMSetFromOptions(*dm));
1060ee167adSMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1070ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1080ee167adSMatthew G. Knepley }
1090ee167adSMatthew G. Knepley 
CreateDiscretization(DM dm,AppCtx * user)1100ee167adSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm, AppCtx *user)
1110ee167adSMatthew G. Knepley {
1120ee167adSMatthew G. Knepley   PetscFE        fe;
1130ee167adSMatthew G. Knepley   PetscFV        fv;
1140ee167adSMatthew G. Knepley   DMPolytopeType ct;
1150ee167adSMatthew G. Knepley   PetscInt       dim, cStart;
1160ee167adSMatthew G. Knepley 
1170ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
1180ee167adSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
1190ee167adSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1200ee167adSMatthew G. Knepley   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1210ee167adSMatthew G. Knepley   if (user->fv) {
1220ee167adSMatthew G. Knepley     PetscCall(PetscFVCreate(PETSC_COMM_SELF, &fv));
1230ee167adSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)fv, "fv"));
1240ee167adSMatthew G. Knepley     PetscCall(PetscFVSetNumComponents(fv, 1));
1250ee167adSMatthew G. Knepley     PetscCall(PetscFVSetSpatialDimension(fv, dim));
1260ee167adSMatthew G. Knepley     PetscCall(PetscFVCreateDualSpace(fv, ct));
1270ee167adSMatthew G. Knepley     PetscCall(PetscFVSetFromOptions(fv));
1280ee167adSMatthew G. Knepley     PetscCall(DMAddField(dm, NULL, (PetscObject)fv));
1290ee167adSMatthew G. Knepley     PetscCall(PetscFVDestroy(&fv));
1300ee167adSMatthew G. Knepley     PetscCall(PetscFVCreate(PETSC_COMM_SELF, &fv));
1310ee167adSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)fv, "fv2"));
1321898fd5cSMatthew G. Knepley     PetscCall(PetscFVSetNumComponents(fv, dim));
1330ee167adSMatthew G. Knepley     PetscCall(PetscFVSetSpatialDimension(fv, dim));
1340ee167adSMatthew G. Knepley     PetscCall(PetscFVCreateDualSpace(fv, ct));
1350ee167adSMatthew G. Knepley     PetscCall(PetscFVSetFromOptions(fv));
1360ee167adSMatthew G. Knepley     PetscCall(DMAddField(dm, NULL, (PetscObject)fv));
1370ee167adSMatthew G. Knepley     PetscCall(PetscFVDestroy(&fv));
1380ee167adSMatthew G. Knepley   } else {
1390ee167adSMatthew G. Knepley     PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, NULL, -1, &fe));
1400ee167adSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
1410ee167adSMatthew G. Knepley     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
1420ee167adSMatthew G. Knepley     PetscCall(PetscFEDestroy(&fe));
1431898fd5cSMatthew G. Knepley     PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, dim, ct, NULL, -1, &fe));
1440ee167adSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)fe, "fe2"));
1450ee167adSMatthew G. Knepley     PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
1460ee167adSMatthew G. Knepley     PetscCall(PetscFEDestroy(&fe));
1470ee167adSMatthew G. Knepley   }
1480ee167adSMatthew G. Knepley   PetscCall(DMCreateDS(dm));
1490ee167adSMatthew G. Knepley   if (user->fv) {
1500ee167adSMatthew G. Knepley     DMLabel  label;
1510ee167adSMatthew G. Knepley     PetscInt values[1] = {1};
1520ee167adSMatthew G. Knepley 
1530ee167adSMatthew G. Knepley     PetscCall(DMCreateLabel(dm, "ghost"));
1540ee167adSMatthew G. Knepley     PetscCall(DMGetLabel(dm, "marker", &label));
1550ee167adSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "dummy", label, 1, values, 0, 0, NULL, NULL, NULL, NULL, NULL));
1560ee167adSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "dummy", label, 1, values, 1, 0, NULL, NULL, NULL, NULL, NULL));
1570ee167adSMatthew G. Knepley   }
1580ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1590ee167adSMatthew G. Knepley }
1600ee167adSMatthew G. Knepley 
CreateGradDiscretization(DM dm,AppCtx * user)1611898fd5cSMatthew G. Knepley static PetscErrorCode CreateGradDiscretization(DM dm, AppCtx *user)
1621898fd5cSMatthew G. Knepley {
1631898fd5cSMatthew G. Knepley   PetscFE        fe;
1641898fd5cSMatthew G. Knepley   DMPolytopeType ct;
1651898fd5cSMatthew G. Knepley   PetscInt       dim, cStart;
1661898fd5cSMatthew G. Knepley 
1671898fd5cSMatthew G. Knepley   PetscFunctionBeginUser;
1681898fd5cSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
1691898fd5cSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
1701898fd5cSMatthew G. Knepley   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1711898fd5cSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, dim, ct, NULL, -1, &fe));
1721898fd5cSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
1731898fd5cSMatthew G. Knepley   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
1741898fd5cSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe));
1751898fd5cSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 2 * dim, ct, NULL, -1, &fe));
1761898fd5cSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fe, "fe2"));
1771898fd5cSMatthew G. Knepley   PetscCall(DMAddField(dm, NULL, (PetscObject)fe));
1781898fd5cSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe));
1791898fd5cSMatthew G. Knepley   PetscCall(DMCreateDS(dm));
1801898fd5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1811898fd5cSMatthew G. Knepley }
1821898fd5cSMatthew G. Knepley 
CreateSwarm(DM dm,DM * sw,AppCtx * user)1830ee167adSMatthew G. Knepley static PetscErrorCode CreateSwarm(DM dm, DM *sw, AppCtx *user)
1840ee167adSMatthew G. Knepley {
1850ee167adSMatthew G. Knepley   PetscScalar *coords, *wvals, *xvals;
1860ee167adSMatthew G. Knepley   PetscInt     Npc = user->Npc, dim, Np;
1870ee167adSMatthew G. Knepley 
1880ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
1890ee167adSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
1900ee167adSMatthew G. Knepley 
1910ee167adSMatthew G. Knepley   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
1920ee167adSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
1930ee167adSMatthew G. Knepley   PetscCall(DMSetType(*sw, DMSWARM));
1940ee167adSMatthew G. Knepley   PetscCall(DMSetDimension(*sw, dim));
1950ee167adSMatthew G. Knepley   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
1960ee167adSMatthew G. Knepley   PetscCall(DMSwarmSetCellDM(*sw, dm));
1970ee167adSMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
1980ee167adSMatthew G. Knepley   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "x_q", 2, PETSC_SCALAR));
1990ee167adSMatthew G. Knepley   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
2000ee167adSMatthew G. Knepley   PetscCall(DMSwarmInsertPointsUsingCellDM(*sw, DMSWARMPIC_LAYOUT_GAUSS, Npc));
2010ee167adSMatthew G. Knepley   PetscCall(DMSetFromOptions(*sw));
2020ee167adSMatthew G. Knepley 
2030ee167adSMatthew G. Knepley   PetscCall(DMSwarmGetLocalSize(*sw, &Np));
2040ee167adSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
2050ee167adSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, "w_q", NULL, NULL, (void **)&wvals));
2060ee167adSMatthew G. Knepley   PetscCall(DMSwarmGetField(*sw, "x_q", NULL, NULL, (void **)&xvals));
2070ee167adSMatthew G. Knepley   for (PetscInt p = 0; p < Np; ++p) {
2080ee167adSMatthew G. Knepley     PetscCall(user->func(dim, 0., &coords[p * dim], 1, &wvals[p], user));
2090ee167adSMatthew G. Knepley     for (PetscInt c = 0; c < 2; ++c) PetscCall(user->func(dim, 0., &coords[p * dim], 1, &xvals[p * 2 + c], user));
2100ee167adSMatthew G. Knepley   }
2110ee167adSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
2120ee167adSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, "w_q", NULL, NULL, (void **)&wvals));
2130ee167adSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(*sw, "x_q", NULL, NULL, (void **)&xvals));
2140ee167adSMatthew G. Knepley 
2150ee167adSMatthew G. Knepley   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
2160ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2170ee167adSMatthew G. Knepley }
2180ee167adSMatthew G. Knepley 
computeParticleMoments(DM sw,Vec u,PetscReal moments[3],AppCtx * user)2190ee167adSMatthew G. Knepley static PetscErrorCode computeParticleMoments(DM sw, Vec u, PetscReal moments[3], AppCtx *user)
2200ee167adSMatthew G. Knepley {
2210ee167adSMatthew G. Knepley   DM                 dm;
2220ee167adSMatthew G. Knepley   const PetscReal   *coords;
2230ee167adSMatthew G. Knepley   const PetscScalar *w;
2240ee167adSMatthew G. Knepley   PetscReal          mom[3] = {0.0, 0.0, 0.0};
2250ee167adSMatthew G. Knepley   PetscInt           dim, cStart, cEnd, Nc;
2260ee167adSMatthew G. Knepley 
2270ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
2280ee167adSMatthew G. Knepley   PetscCall(DMGetDimension(sw, &dim));
2290ee167adSMatthew G. Knepley   PetscCall(DMSwarmGetCellDM(sw, &dm));
2300ee167adSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
2310ee167adSMatthew G. Knepley   PetscCall(DMSwarmSortGetAccess(sw));
2320ee167adSMatthew G. Knepley   PetscCall(DMSwarmGetFieldInfo(sw, user->field ? "x_q" : "w_q", &Nc, NULL));
2330ee167adSMatthew G. Knepley   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
2340ee167adSMatthew G. Knepley   PetscCall(VecGetArrayRead(u, &w));
2350ee167adSMatthew G. Knepley   for (PetscInt cell = cStart; cell < cEnd; ++cell) {
2360ee167adSMatthew G. Knepley     PetscInt *pidx;
2370ee167adSMatthew G. Knepley     PetscInt  Np;
2380ee167adSMatthew G. Knepley 
2390ee167adSMatthew G. Knepley     PetscCall(DMSwarmSortGetPointsPerCell(sw, cell, &Np, &pidx));
2400ee167adSMatthew G. Knepley     for (PetscInt p = 0; p < Np; ++p) {
2410ee167adSMatthew G. Knepley       const PetscInt   idx = pidx[p];
2420ee167adSMatthew G. Knepley       const PetscReal *x   = &coords[idx * dim];
2430ee167adSMatthew G. Knepley 
2440ee167adSMatthew G. Knepley       for (PetscInt c = 0; c < Nc; ++c) {
2450ee167adSMatthew G. Knepley         mom[0] += PetscRealPart(w[idx * Nc + c]);
2460ee167adSMatthew G. Knepley         mom[1] += PetscRealPart(w[idx * Nc + c]) * x[0];
2470ee167adSMatthew G. Knepley         for (PetscInt d = 0; d < dim; ++d) mom[2] += PetscRealPart(w[idx * Nc + c]) * PetscSqr(x[d]);
2480ee167adSMatthew G. Knepley       }
2490ee167adSMatthew G. Knepley     }
250fe11354eSMatthew G. Knepley     PetscCall(DMSwarmSortRestorePointsPerCell(sw, cell, &Np, &pidx));
2510ee167adSMatthew G. Knepley   }
2520ee167adSMatthew G. Knepley   PetscCall(VecRestoreArrayRead(u, &w));
2530ee167adSMatthew G. Knepley   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
2540ee167adSMatthew G. Knepley   PetscCall(DMSwarmSortRestoreAccess(sw));
255462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(mom, moments, 3, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)sw)));
2560ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2570ee167adSMatthew G. Knepley }
2580ee167adSMatthew G. Knepley 
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[])2590ee167adSMatthew G. Knepley 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[])
2600ee167adSMatthew G. Knepley {
2610ee167adSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
2620ee167adSMatthew G. Knepley 
2630ee167adSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[0] += u[c];
2640ee167adSMatthew G. Knepley }
2650ee167adSMatthew G. Knepley 
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[])2660ee167adSMatthew G. Knepley 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[])
2670ee167adSMatthew G. Knepley {
2680ee167adSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
2690ee167adSMatthew G. Knepley 
2700ee167adSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[0] += x[0] * u[c];
2710ee167adSMatthew G. Knepley }
2720ee167adSMatthew G. Knepley 
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[])2730ee167adSMatthew G. Knepley 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[])
2740ee167adSMatthew G. Knepley {
2750ee167adSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
2760ee167adSMatthew G. Knepley 
2770ee167adSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c)
2780ee167adSMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) f0[0] += PetscSqr(x[d]) * u[c];
2790ee167adSMatthew G. Knepley }
2800ee167adSMatthew G. Knepley 
computeFieldMoments(DM dm,Vec u,PetscReal moments[3],AppCtx * user)2810ee167adSMatthew G. Knepley static PetscErrorCode computeFieldMoments(DM dm, Vec u, PetscReal moments[3], AppCtx *user)
2820ee167adSMatthew G. Knepley {
2830ee167adSMatthew G. Knepley   PetscDS     ds;
2840ee167adSMatthew G. Knepley   PetscScalar mom;
2850ee167adSMatthew G. Knepley 
2860ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
2870ee167adSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
2880ee167adSMatthew G. Knepley   PetscCall(PetscDSSetObjective(ds, 0, &f0_1));
2890ee167adSMatthew G. Knepley   mom = 0.;
2900ee167adSMatthew G. Knepley   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
2910ee167adSMatthew G. Knepley   moments[0] = PetscRealPart(mom);
2920ee167adSMatthew G. Knepley   PetscCall(PetscDSSetObjective(ds, 0, &f0_x));
2930ee167adSMatthew G. Knepley   mom = 0.;
2940ee167adSMatthew G. Knepley   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
2950ee167adSMatthew G. Knepley   moments[1] = PetscRealPart(mom);
2960ee167adSMatthew G. Knepley   PetscCall(PetscDSSetObjective(ds, 0, &f0_r2));
2970ee167adSMatthew G. Knepley   mom = 0.;
2980ee167adSMatthew G. Knepley   PetscCall(DMPlexComputeIntegralFEM(dm, u, &mom, user));
2990ee167adSMatthew G. Knepley   moments[2] = PetscRealPart(mom);
3000ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3010ee167adSMatthew G. Knepley }
3020ee167adSMatthew G. Knepley 
TestParticlesToField(DM sw,DM dm,Vec fhat,AppCtx * user)3030ee167adSMatthew G. Knepley static PetscErrorCode TestParticlesToField(DM sw, DM dm, Vec fhat, AppCtx *user)
3040ee167adSMatthew G. Knepley {
3050ee167adSMatthew G. Knepley   const char *fieldnames[1] = {user->field ? "x_q" : "w_q"};
3060ee167adSMatthew G. Knepley   Vec         fields[1]     = {fhat}, f;
3070ee167adSMatthew G. Knepley   PetscReal   pmoments[3]; // \int f, \int x f, \int r^2 f
3080ee167adSMatthew G. Knepley   PetscReal   fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f
3090ee167adSMatthew G. Knepley 
3100ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
3110ee167adSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, dm, 1, fieldnames, fields, SCATTER_FORWARD));
3120ee167adSMatthew G. Knepley 
3130ee167adSMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, fieldnames[0], &f));
3140ee167adSMatthew G. Knepley   PetscCall(computeParticleMoments(sw, f, pmoments, user));
3151898fd5cSMatthew G. Knepley   PetscCall(VecViewFromOptions(f, NULL, "-f_view"));
3160ee167adSMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, fieldnames[0], &f));
3170ee167adSMatthew G. Knepley   PetscCall(computeFieldMoments(dm, fhat, fmoments, user));
3180ee167adSMatthew G. Knepley   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
3190ee167adSMatthew G. Knepley   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
3200ee167adSMatthew G. Knepley   for (PetscInt m = 0; m < user->Nm; ++m) {
3210ee167adSMatthew G. Knepley     if (user->pass) {
3223a7d0413SPierre Jolivet       if (PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]) > user->mtol) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "p  projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", pmoments[0], pmoments[1], pmoments[2]));
3230ee167adSMatthew G. Knepley     } else {
3240ee167adSMatthew G. Knepley       PetscCheck(PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]) <= user->mtol, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Moment %" PetscInt_FMT " error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]),
3250ee167adSMatthew G. Knepley                  user->mtol);
3260ee167adSMatthew G. Knepley     }
3270ee167adSMatthew G. Knepley   }
3280ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3290ee167adSMatthew G. Knepley }
3300ee167adSMatthew G. Knepley 
TestFieldToParticles(DM sw,DM dm,Vec fhat,AppCtx * user)3310ee167adSMatthew G. Knepley static PetscErrorCode TestFieldToParticles(DM sw, DM dm, Vec fhat, AppCtx *user)
3320ee167adSMatthew G. Knepley {
3330ee167adSMatthew G. Knepley   const char *fieldnames[1] = {user->field ? "x_q" : "w_q"};
3340ee167adSMatthew G. Knepley   Vec         fields[1]     = {fhat}, f;
3350ee167adSMatthew G. Knepley   PetscReal   pmoments[3]; // \int f, \int x f, \int r^2 f
3360ee167adSMatthew G. Knepley   PetscReal   fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f
3370ee167adSMatthew G. Knepley 
3380ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
3390ee167adSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, dm, 1, fieldnames, fields, SCATTER_REVERSE));
3400ee167adSMatthew G. Knepley 
3410ee167adSMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, fieldnames[0], &f));
3420ee167adSMatthew G. Knepley   PetscCall(computeParticleMoments(sw, f, pmoments, user));
3430ee167adSMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, fieldnames[0], &f));
3440ee167adSMatthew G. Knepley   PetscCall(computeFieldMoments(dm, fhat, fmoments, user));
3450ee167adSMatthew G. Knepley   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "L2 projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", fmoments[0], fmoments[1], fmoments[2]));
3460ee167adSMatthew G. Knepley   for (PetscInt m = 0; m < user->Nm; ++m) {
3470ee167adSMatthew G. Knepley     if (user->pass) {
3483a7d0413SPierre Jolivet       if (PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]) > user->mtol) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "p  projection mass: %20.10e, x-momentum: %20.10e, energy: %20.10e.\n", pmoments[0], pmoments[1], pmoments[2]));
3490ee167adSMatthew G. Knepley     } else {
3500ee167adSMatthew G. Knepley       PetscCheck(PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]) <= user->mtol, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Moment %" PetscInt_FMT " error too large %g > %g", m, PetscAbsReal((fmoments[m] - pmoments[m]) / fmoments[m]),
3510ee167adSMatthew G. Knepley                  user->mtol);
3520ee167adSMatthew G. Knepley     }
3530ee167adSMatthew G. Knepley   }
3540ee167adSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3550ee167adSMatthew G. Knepley }
3560ee167adSMatthew G. Knepley 
TestParticlesToGradientField(DM sw,DM dm,Vec fhat,AppCtx * user)3571898fd5cSMatthew G. Knepley static PetscErrorCode TestParticlesToGradientField(DM sw, DM dm, Vec fhat, AppCtx *user)
3581898fd5cSMatthew G. Knepley {
3591898fd5cSMatthew G. Knepley   const char *fieldnames[1] = {"x_q"};
3601898fd5cSMatthew G. Knepley   Vec         fields[1]     = {fhat}, f;
3611898fd5cSMatthew G. Knepley   PetscReal   pmoments[3]; // \int f, \int x f, \int r^2 f
3621898fd5cSMatthew G. Knepley   PetscReal   fmoments[3]; // \int \hat f, \int x \hat f, \int r^2 \hat f
3631898fd5cSMatthew G. Knepley 
3641898fd5cSMatthew G. Knepley   PetscFunctionBeginUser;
3651898fd5cSMatthew G. Knepley   PetscCall(DMSwarmProjectGradientFields(sw, dm, 1, fieldnames, fields, SCATTER_FORWARD));
3661898fd5cSMatthew G. Knepley 
3671898fd5cSMatthew G. Knepley   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, fieldnames[0], &f));
3681898fd5cSMatthew G. Knepley   PetscCall(computeParticleMoments(sw, f, pmoments, user));
3691898fd5cSMatthew G. Knepley   PetscCall(VecViewFromOptions(f, NULL, "-f_view"));
3701898fd5cSMatthew G. Knepley   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, fieldnames[0], &f));
3711898fd5cSMatthew G. Knepley   PetscCall(computeFieldMoments(dm, fhat, fmoments, user));
3721898fd5cSMatthew G. Knepley   PetscCall(VecViewFromOptions(fhat, NULL, "-fhat_view"));
3731898fd5cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3741898fd5cSMatthew G. Knepley }
3751898fd5cSMatthew G. Knepley 
main(int argc,char * argv[])3760ee167adSMatthew G. Knepley int main(int argc, char *argv[])
3770ee167adSMatthew G. Knepley {
3780ee167adSMatthew G. Knepley   DM     dm, subdm, sw;
3790ee167adSMatthew G. Knepley   Vec    fhat;
3800ee167adSMatthew G. Knepley   IS     subis;
3810ee167adSMatthew G. Knepley   AppCtx user;
3820ee167adSMatthew G. Knepley 
3830ee167adSMatthew G. Knepley   PetscFunctionBeginUser;
3840ee167adSMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3850ee167adSMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3860ee167adSMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &user));
3870ee167adSMatthew G. Knepley   PetscCall(CreateDiscretization(dm, &user));
3880ee167adSMatthew G. Knepley   PetscCall(CreateSwarm(dm, &sw, &user));
3890ee167adSMatthew G. Knepley 
3900ee167adSMatthew G. Knepley   PetscCall(DMCreateSubDM(dm, 1, &user.field, &subis, &subdm));
3910ee167adSMatthew G. Knepley 
3920ee167adSMatthew G. Knepley   PetscCall(DMGetGlobalVector(subdm, &fhat));
3930ee167adSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fhat, "FEM f"));
3940ee167adSMatthew G. Knepley   PetscCall(TestParticlesToField(sw, subdm, fhat, &user));
3950ee167adSMatthew G. Knepley   PetscCall(TestFieldToParticles(sw, subdm, fhat, &user));
3960ee167adSMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(subdm, &fhat));
3970ee167adSMatthew G. Knepley 
3981898fd5cSMatthew G. Knepley   if (user.grad) {
3991898fd5cSMatthew G. Knepley     DM dmGrad, gsubdm;
4001898fd5cSMatthew G. Knepley     IS gsubis;
4011898fd5cSMatthew G. Knepley 
4021898fd5cSMatthew G. Knepley     PetscCall(DMClone(dm, &dmGrad));
4031898fd5cSMatthew G. Knepley     PetscCall(CreateGradDiscretization(dmGrad, &user));
4041898fd5cSMatthew G. Knepley     PetscCall(DMCreateSubDM(dmGrad, 1, &user.field, &gsubis, &gsubdm));
4051898fd5cSMatthew G. Knepley 
4061898fd5cSMatthew G. Knepley     PetscCall(DMGetGlobalVector(gsubdm, &fhat));
4071898fd5cSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)fhat, "FEM grad f"));
4081898fd5cSMatthew G. Knepley     PetscCall(TestParticlesToGradientField(sw, subdm, fhat, &user));
4091898fd5cSMatthew G. Knepley     //PetscCall(TestFieldToParticles(sw, subdm, fhat, &user));
4101898fd5cSMatthew G. Knepley     PetscCall(DMRestoreGlobalVector(gsubdm, &fhat));
4111898fd5cSMatthew G. Knepley     PetscCall(ISDestroy(&gsubis));
4121898fd5cSMatthew G. Knepley     PetscCall(DMDestroy(&gsubdm));
4131898fd5cSMatthew G. Knepley     PetscCall(DMDestroy(&dmGrad));
4141898fd5cSMatthew G. Knepley   }
4151898fd5cSMatthew G. Knepley 
4160ee167adSMatthew G. Knepley   PetscCall(ISDestroy(&subis));
4170ee167adSMatthew G. Knepley   PetscCall(DMDestroy(&subdm));
4180ee167adSMatthew G. Knepley   PetscCall(DMDestroy(&dm));
4190ee167adSMatthew G. Knepley   PetscCall(DMDestroy(&sw));
4200ee167adSMatthew G. Knepley   PetscCall(PetscFinalize());
4210ee167adSMatthew G. Knepley   return 0;
4220ee167adSMatthew G. Knepley }
4230ee167adSMatthew G. Knepley 
4240ee167adSMatthew G. Knepley /*TEST
4250ee167adSMatthew G. Knepley 
4260ee167adSMatthew G. Knepley   # Swarm does not handle complex or quad
4270ee167adSMatthew G. Knepley   build:
4280ee167adSMatthew G. Knepley     requires: !complex double
4290ee167adSMatthew G. Knepley 
4300ee167adSMatthew G. Knepley   testset:
4310ee167adSMatthew G. Knepley     requires: triangle
4320ee167adSMatthew G. Knepley     args: -dm_refine 1 -petscspace_degree 2 -moments 3 \
4330ee167adSMatthew G. Knepley           -ptof_pc_type lu  \
4340ee167adSMatthew G. Knepley           -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
4350ee167adSMatthew G. Knepley 
4360ee167adSMatthew G. Knepley     test:
4370ee167adSMatthew G. Knepley       suffix: tri_fe_f0
4380ee167adSMatthew G. Knepley       args: -field 0
4390ee167adSMatthew G. Knepley 
4400ee167adSMatthew G. Knepley     test:
4410ee167adSMatthew G. Knepley       suffix: tri_fe_f1
4420ee167adSMatthew G. Knepley       args: -field 1
4430ee167adSMatthew G. Knepley 
4440ee167adSMatthew G. Knepley     test:
4451898fd5cSMatthew G. Knepley       suffix: tri_fe_grad
4461898fd5cSMatthew G. Knepley       args: -field 0 -grad -gptof_ksp_type lsqr -gptof_pc_type none -gptof_ksp_rtol 1e-10
4471898fd5cSMatthew G. Knepley 
4481898fd5cSMatthew G. Knepley     # -gptof_ksp_converged_reason -gptof_ksp_lsqr_monitor to see the divergence solve
4491898fd5cSMatthew G. Knepley     test:
4500ee167adSMatthew G. Knepley       suffix: quad_fe_f0
4510ee167adSMatthew G. Knepley       args: -dm_plex_simplex 0 -field 0
4520ee167adSMatthew G. Knepley 
4530ee167adSMatthew G. Knepley     test:
4540ee167adSMatthew G. Knepley       suffix: quad_fe_f1
4550ee167adSMatthew G. Knepley       args: -dm_plex_simplex 0 -field 1
4560ee167adSMatthew G. Knepley 
4570ee167adSMatthew G. Knepley   testset:
4580ee167adSMatthew G. Knepley     requires: triangle
4590ee167adSMatthew G. Knepley     args: -dm_refine 1 -moments 1 -fv \
4600ee167adSMatthew G. Knepley           -ptof_pc_type lu \
4610ee167adSMatthew G. Knepley           -ftop_ksp_rtol 1e-15 -ftop_ksp_type lsqr -ftop_pc_type none
4620ee167adSMatthew G. Knepley 
4630ee167adSMatthew G. Knepley     test:
4640ee167adSMatthew G. Knepley       suffix: tri_fv_f0
4650ee167adSMatthew G. Knepley       args: -field 0
4660ee167adSMatthew G. Knepley 
4670ee167adSMatthew G. Knepley     test:
4680ee167adSMatthew G. Knepley       suffix: tri_fv_f1
4690ee167adSMatthew G. Knepley       args: -field 1
4700ee167adSMatthew G. Knepley 
4710ee167adSMatthew G. Knepley     test:
4720ee167adSMatthew G. Knepley       suffix: quad_fv_f0
4730ee167adSMatthew G. Knepley       args: -dm_plex_simplex 0 -field 0
4740ee167adSMatthew G. Knepley 
4750ee167adSMatthew G. Knepley     test:
4760ee167adSMatthew G. Knepley       suffix: quad_fv_f1
4770ee167adSMatthew G. Knepley       args: -dm_plex_simplex 0 -field 1
4780ee167adSMatthew G. Knepley 
4790ee167adSMatthew G. Knepley TEST*/
480