xref: /petsc/src/dm/impls/swarm/tutorials/ex1.c (revision 7f71fc7e399543653e1829a2685bc041125727a8)
14cd81d59SMatthew G. Knepley static char help[] = "Example program demonstrating projection between particle and finite element spaces\n\n";
24cd81d59SMatthew G. Knepley 
346233b44SBarry Smith #include <petscdmplex.h>
446233b44SBarry Smith #include <petscds.h>
546233b44SBarry Smith #include <petscdmswarm.h>
646233b44SBarry Smith #include <petscksp.h>
74cd81d59SMatthew G. Knepley 
main(int argc,char ** argv)8d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
9d71ae5a4SJacob Faibussowitsch {
103c611800SMark Adams   DM             dm, sw;
114cd81d59SMatthew G. Knepley   PetscFE        fe;
12fc3eb83cSMatthew G. Knepley   Vec            u_f;
13fc3eb83cSMatthew G. Knepley   DMPolytopeType ct;
14fc3eb83cSMatthew G. Knepley   PetscInt       dim, Nc = 1, faces[3];
15fc3eb83cSMatthew G. Knepley   PetscInt       Np = 10, field = 0, zero = 0, bs, cStart;
16fc3eb83cSMatthew G. Knepley   PetscReal      energy_0 = 0, energy_1 = 0;
1730602db0SMatthew G. Knepley   PetscReal      lo[3], hi[3], h[3];
1830602db0SMatthew G. Knepley   PetscBool      removePoints = PETSC_TRUE;
1930602db0SMatthew G. Knepley   PetscReal     *wq, *coords;
204cd81d59SMatthew G. Knepley   PetscDataType  dtype;
214cd81d59SMatthew G. Knepley 
22327415f7SBarry Smith   PetscFunctionBeginUser;
239566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
244cd81d59SMatthew G. Knepley   /* Create a mesh */
259566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
269566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
279566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
289566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
2930602db0SMatthew G. Knepley 
309566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
31fc3eb83cSMatthew G. Knepley   bs = dim;
32fc3eb83cSMatthew G. Knepley   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &bs, NULL));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-np", &Np, NULL));
349566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dm, lo, hi));
35fc3eb83cSMatthew G. Knepley   for (PetscInt i = 0; i < dim; ++i) {
3630602db0SMatthew G. Knepley     h[i] = (hi[i] - lo[i]) / faces[i];
3763a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, " lo = %g hi = %g n = %" PetscInt_FMT " h = %g\n", (double)lo[i], (double)hi[i], faces[i], (double)h[i]));
3830602db0SMatthew G. Knepley   }
39fc3eb83cSMatthew G. Knepley   // Create FE space
40fc3eb83cSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
41fc3eb83cSMatthew G. Knepley   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
42fc3eb83cSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, Nc, ct, NULL, PETSC_DECIDE, &fe));
439566063dSJacob Faibussowitsch   PetscCall(PetscFESetFromOptions(fe));
449566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
459566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, field, NULL, (PetscObject)fe));
469566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
479566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
48fc3eb83cSMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &u_f));
49fc3eb83cSMatthew G. Knepley   // Create particle swarm
509566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_SELF, &sw));
519566063dSJacob Faibussowitsch   PetscCall(DMSetType(sw, DMSWARM));
529566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(sw, dim));
539566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(sw, DMSWARM_PIC));
549566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(sw, dm));
559566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(sw, "w_q", Nc, PETSC_SCALAR));
569566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(sw));
579566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
589566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(sw));
599566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
609566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
61fc3eb83cSMatthew G. Knepley   for (PetscInt p = 0; p < Np; ++p) {
624cd81d59SMatthew G. Knepley     coords[p * 2 + 0] = -PetscCosReal((PetscReal)(p + 1) / (PetscReal)(Np + 1) * PETSC_PI);
634cd81d59SMatthew G. Knepley     coords[p * 2 + 1] = PetscSinReal((PetscReal)(p + 1) / (PetscReal)(Np + 1) * PETSC_PI);
644cd81d59SMatthew G. Knepley     wq[p]             = 1.0;
653c611800SMark Adams     energy_0 += wq[p] * (PetscSqr(coords[p * 2 + 0]) + PetscSqr(coords[p * 2 + 1]));
664cd81d59SMatthew G. Knepley   }
679566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
689566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
699566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate(sw, removePoints));
709566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
719566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(sw, NULL, "-swarm_view"));
723c611800SMark Adams 
73fc3eb83cSMatthew G. Knepley   // Project between particles and continuum field
74fc3eb83cSMatthew G. Knepley   const char *fieldnames[1] = {"w_q"};
75fc3eb83cSMatthew G. Knepley   Vec         fields[1]     = {u_f};
76*953494dbSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, NULL, 1, fieldnames, fields, SCATTER_FORWARD));
77*953494dbSMatthew G. Knepley   PetscCall(DMSwarmProjectFields(sw, NULL, 1, fieldnames, fields, SCATTER_REVERSE));
783c611800SMark Adams 
79fc3eb83cSMatthew G. Knepley   // Compute energy
809566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
819566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
82fc3eb83cSMatthew G. Knepley   for (PetscInt p = 0; p < Np; ++p) energy_1 += wq[p] * (PetscSqr(coords[p * 2 + 0]) + PetscSqr(coords[p * 2 + 1]));
839566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
849566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
85fc3eb83cSMatthew G. Knepley   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Energy = %20.12e error = %20.12e\n", (double)energy_0, (double)((energy_1 - energy_0) / energy_0)));
86fc3eb83cSMatthew G. Knepley 
87fc3eb83cSMatthew G. Knepley   // Cleanup
88fc3eb83cSMatthew G. Knepley   PetscCall(VecDestroy(&u_f));
899566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&sw));
909566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
919566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
92b122ec5aSJacob Faibussowitsch   return 0;
934cd81d59SMatthew G. Knepley }
944cd81d59SMatthew G. Knepley 
954cd81d59SMatthew G. Knepley /*TEST
964cd81d59SMatthew G. Knepley 
974cd81d59SMatthew G. Knepley   build:
984cd81d59SMatthew G. Knepley     requires: !complex
994cd81d59SMatthew G. Knepley 
1004cd81d59SMatthew G. Knepley   test:
1014cd81d59SMatthew G. Knepley     suffix: 0
1023c611800SMark Adams     requires: double triangle
103fc3eb83cSMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 \
104fc3eb83cSMatthew G. Knepley            -np 50 -petscspace_degree 2 \
105fc3eb83cSMatthew G. Knepley            -ptof_ksp_type cg -ptof_pc_type ilu -ptof_ksp_rtol 1.e-14 \
106fc3eb83cSMatthew G. Knepley            -ftop_ksp_type lsqr -ftop_pc_type none -ftop_ksp_rtol 1.e-14 \
107fc3eb83cSMatthew G. Knepley            -dm_view -swarm_view
1083c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
1093c611800SMark Adams 
1103c611800SMark Adams   test:
1113c611800SMark Adams     suffix: bjacobi
1123c611800SMark Adams     requires: double triangle
113fc3eb83cSMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 \
114fc3eb83cSMatthew G. Knepley           -np 50 -petscspace_degree 2 -dm_plex_hash_location \
115fc3eb83cSMatthew G. Knepley           -ptof_ksp_type cg -ptof_pc_type ilu -ptof_ksp_rtol 1.e-14 \
116fc3eb83cSMatthew G. Knepley           -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero \
117fc3eb83cSMatthew G. Knepley           -dm_view -swarm_view -ftop_ksp_rtol 1.e-14
1184cd81d59SMatthew G. Knepley     filter: grep -v DM_ | grep -v atomic
1194cd81d59SMatthew G. Knepley 
1204cd81d59SMatthew G. Knepley TEST*/
121