xref: /petsc/src/dm/impls/swarm/tutorials/ex1.c (revision 30602db00db74b7e41a0c75e517aefe6711423f0)
14cd81d59SMatthew G. Knepley static char help[] = "Example program demonstrating projection between particle and finite element spaces\n\n";
24cd81d59SMatthew G. Knepley 
34cd81d59SMatthew G. Knepley #include "petscdmplex.h"
44cd81d59SMatthew G. Knepley #include "petscds.h"
54cd81d59SMatthew G. Knepley #include "petscdmswarm.h"
64cd81d59SMatthew G. Knepley #include "petscksp.h"
74cd81d59SMatthew G. Knepley 
84cd81d59SMatthew G. Knepley static PetscErrorCode crd_func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
94cd81d59SMatthew G. Knepley {
104cd81d59SMatthew G. Knepley   PetscInt i;
114cd81d59SMatthew G. Knepley   PetscFunctionBeginUser;
124cd81d59SMatthew G. Knepley   for (i = 0; i < dim; ++i) u[i] = x[i];
134cd81d59SMatthew G. Knepley   PetscFunctionReturn(0);
144cd81d59SMatthew G. Knepley }
154cd81d59SMatthew G. Knepley 
164cd81d59SMatthew G. Knepley int main(int argc, char **argv)
174cd81d59SMatthew G. Knepley {
184cd81d59SMatthew G. Knepley   DM              dm, crddm, sw;
194cd81d59SMatthew G. Knepley   PetscFE         fe;
204cd81d59SMatthew G. Knepley   KSP             ksp;
214cd81d59SMatthew G. Knepley   Mat             M_p, M;
224cd81d59SMatthew G. Knepley   Vec             f, rho, rhs, crd_vec;
23*30602db0SMatthew G. Knepley   PetscInt        dim, Nc = 1, timestep = 0, N, i, idx[3], faces[3];
24*30602db0SMatthew G. Knepley   PetscInt        Np = 10, p, field = 0, zero = 0, bs;
254cd81d59SMatthew G. Knepley   PetscReal       time = 0.0,  norm;
26*30602db0SMatthew G. Knepley   PetscReal       lo[3], hi[3], h[3];
27*30602db0SMatthew G. Knepley   PetscBool       removePoints = PETSC_TRUE;
284cd81d59SMatthew G. Knepley   const PetscReal *xx, *vv;
29*30602db0SMatthew G. Knepley   PetscReal       *wq, *coords;
304cd81d59SMatthew G. Knepley   PetscDataType   dtype;
314cd81d59SMatthew G. Knepley   PetscErrorCode  (*initu[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *);
32*30602db0SMatthew G. Knepley   PetscErrorCode  ierr;
334cd81d59SMatthew G. Knepley 
344cd81d59SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
354cd81d59SMatthew G. Knepley   /* Create a mesh */
36*30602db0SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);
37*30602db0SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
38*30602db0SMatthew G. Knepley   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
394cd81d59SMatthew G. Knepley   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
40*30602db0SMatthew G. Knepley 
41*30602db0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
42*30602db0SMatthew G. Knepley   i    = dim;
43*30602db0SMatthew G. Knepley   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL);CHKERRQ(ierr);
44*30602db0SMatthew G. Knepley   ierr = DMGetBoundingBox(dm, lo, hi);CHKERRQ(ierr);
45*30602db0SMatthew G. Knepley   for (i=0;i<dim;i++) {
46*30602db0SMatthew G. Knepley     h[i] = (hi[i] - lo[i])/faces[i];
47*30602db0SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF," lo = %g hi = %g n = %D h = %g\n",lo[i],hi[i],faces[i],h[i]);CHKERRQ(ierr);
48*30602db0SMatthew G. Knepley   }
49*30602db0SMatthew G. Knepley 
504cd81d59SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe);CHKERRQ(ierr);
514cd81d59SMatthew G. Knepley   ierr = PetscFESetFromOptions(fe);CHKERRQ(ierr);
524cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)fe, "fe");CHKERRQ(ierr);
534cd81d59SMatthew G. Knepley   ierr = DMSetField(dm, field, NULL, (PetscObject)fe);CHKERRQ(ierr);
544cd81d59SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
554cd81d59SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
564cd81d59SMatthew G. Knepley   /* Create particle swarm */
574cd81d59SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_SELF, &sw);CHKERRQ(ierr);
584cd81d59SMatthew G. Knepley   ierr = DMSetType(sw, DMSWARM);CHKERRQ(ierr);
594cd81d59SMatthew G. Knepley   ierr = DMSetDimension(sw, dim);CHKERRQ(ierr);
604cd81d59SMatthew G. Knepley   ierr = DMSwarmSetType(sw, DMSWARM_PIC);CHKERRQ(ierr);
614cd81d59SMatthew G. Knepley   ierr = DMSwarmSetCellDM(sw, dm);CHKERRQ(ierr);
624cd81d59SMatthew G. Knepley   ierr = DMSwarmRegisterPetscDatatypeField(sw, "w_q", Nc, PETSC_SCALAR);CHKERRQ(ierr);
634cd81d59SMatthew G. Knepley   ierr = DMSwarmFinalizeFieldRegister(sw);CHKERRQ(ierr);
644cd81d59SMatthew G. Knepley   ierr = DMSwarmSetLocalSizes(sw, Np, zero);CHKERRQ(ierr);
654cd81d59SMatthew G. Knepley   ierr = DMSetFromOptions(sw);CHKERRQ(ierr);
664cd81d59SMatthew G. Knepley   ierr = DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
674cd81d59SMatthew G. Knepley   ierr = DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
684cd81d59SMatthew G. Knepley   for (p=0;p<Np;p++) {
694cd81d59SMatthew G. Knepley     coords[p*2+0]  = -PetscCosReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI);
704cd81d59SMatthew G. Knepley     coords[p*2+1] =   PetscSinReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI);
714cd81d59SMatthew G. Knepley     wq[p]          = 1.0;
724cd81d59SMatthew G. Knepley   }
734cd81d59SMatthew G. Knepley   ierr = DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
744cd81d59SMatthew G. Knepley   ierr = DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
754cd81d59SMatthew G. Knepley   ierr = DMSwarmMigrate(sw, removePoints);CHKERRQ(ierr);
764cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)sw, "Particle Grid");CHKERRQ(ierr);
774cd81d59SMatthew G. Knepley   ierr = DMViewFromOptions(sw, NULL, "-swarm_view");CHKERRQ(ierr);
784cd81d59SMatthew G. Knepley   /* Project particles to field */
794cd81d59SMatthew G. Knepley   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
804cd81d59SMatthew G. Knepley   ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr);
814cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &rho);CHKERRQ(ierr);
824cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)rho, "rho");CHKERRQ(ierr);
834cd81d59SMatthew G. Knepley   ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr);
844cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)f, "weights");CHKERRQ(ierr);
854cd81d59SMatthew G. Knepley   ierr = MatMultTranspose(M_p, f, rho);CHKERRQ(ierr);
864cd81d59SMatthew G. Knepley 
874cd81d59SMatthew G. Knepley   /* Visualize mesh field */
884cd81d59SMatthew G. Knepley   ierr = DMSetOutputSequenceNumber(dm, timestep, time);CHKERRQ(ierr);
894cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(rho, NULL, "-rho_view");CHKERRQ(ierr);
904cd81d59SMatthew G. Knepley 
914cd81d59SMatthew G. Knepley   /* create coordinate DM */
924cd81d59SMatthew G. Knepley   ierr = DMClone(dm, &crddm);CHKERRQ(ierr);
934cd81d59SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, PETSC_FALSE, "", PETSC_DECIDE, &fe);CHKERRQ(ierr);
944cd81d59SMatthew G. Knepley   ierr = PetscFESetFromOptions(fe);CHKERRQ(ierr);
954cd81d59SMatthew G. Knepley   ierr = DMSetField(crddm, field, NULL, (PetscObject)fe);CHKERRQ(ierr);
964cd81d59SMatthew G. Knepley   ierr = DMCreateDS(crddm);CHKERRQ(ierr);
974cd81d59SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
984cd81d59SMatthew G. Knepley   /* project coordiantes to vertices */
994cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(crddm, &crd_vec);CHKERRQ(ierr);
1004cd81d59SMatthew G. Knepley   initu[0] = crd_func;
1014cd81d59SMatthew G. Knepley   ierr = DMProjectFunction(crddm, 0.0, initu, NULL, INSERT_ALL_VALUES, crd_vec);CHKERRQ(ierr);
1024cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(crd_vec, NULL, "-coord_view");CHKERRQ(ierr);
1034cd81d59SMatthew G. Knepley   /* iterate over mesh data and get indices */
1044cd81d59SMatthew G. Knepley   ierr = VecGetArrayRead(crd_vec,&xx);CHKERRQ(ierr);
1054cd81d59SMatthew G. Knepley   ierr = VecGetArrayRead(rho,&vv);CHKERRQ(ierr);
1064cd81d59SMatthew G. Knepley   ierr = VecGetLocalSize(rho,&N);CHKERRQ(ierr);
1074cd81d59SMatthew G. Knepley   for (p=0;p<N;p++) {
1084cd81d59SMatthew G. Knepley     for (i=0;i<dim;i++) idx[i] = (PetscInt)((xx[p*dim+i] - lo[i])/h[i] + 1.e-8);
1094cd81d59SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"(%D,%D) = %g\n",idx[0],idx[1],vv[p]);CHKERRQ(ierr);
1104cd81d59SMatthew G. Knepley     /* access grid data here */
1114cd81d59SMatthew G. Knepley   }
1124cd81d59SMatthew G. Knepley   ierr = VecRestoreArrayRead(crd_vec,&xx);CHKERRQ(ierr);
1134cd81d59SMatthew G. Knepley   ierr = VecRestoreArrayRead(rho,&vv);CHKERRQ(ierr);
1144cd81d59SMatthew G. Knepley   ierr = VecDestroy(&crd_vec);CHKERRQ(ierr);
1154cd81d59SMatthew G. Knepley   /* Project field to particles */
1164cd81d59SMatthew G. Knepley   /*   This gives f_p = M_p^+ M f */
1174cd81d59SMatthew G. Knepley   ierr = DMCreateMassMatrix(dm, dm, &M);CHKERRQ(ierr);
1184cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &rhs);CHKERRQ(ierr);
1194cd81d59SMatthew G. Knepley   if (0) {
1204cd81d59SMatthew G. Knepley     ierr = MatMult(M, rho, rhs);CHKERRQ(ierr);  /* this is what you would do for an FE solve */
1214cd81d59SMatthew G. Knepley   } else {
122*30602db0SMatthew G. Knepley     ierr = VecCopy(rho, rhs);CHKERRQ(ierr); /* Identity: M^1 M rho */
1234cd81d59SMatthew G. Knepley   }
1244cd81d59SMatthew G. Knepley   ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr);
1254cd81d59SMatthew G. Knepley   ierr = KSPSetOptionsPrefix(ksp, "ftop_");CHKERRQ(ierr);
1264cd81d59SMatthew G. Knepley   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
1274cd81d59SMatthew G. Knepley   ierr = KSPSetOperators(ksp, M_p, M_p);CHKERRQ(ierr);
1284cd81d59SMatthew G. Knepley   ierr = KSPSolveTranspose(ksp, rhs, f);CHKERRQ(ierr);
1294cd81d59SMatthew G. Knepley   ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
1304cd81d59SMatthew G. Knepley   ierr = VecDestroy(&rhs);CHKERRQ(ierr);
1314cd81d59SMatthew G. Knepley 
1324cd81d59SMatthew G. Knepley   /* Visualize particle field */
1334cd81d59SMatthew G. Knepley   ierr = DMSetOutputSequenceNumber(sw, timestep, time);CHKERRQ(ierr);
1344cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(f, NULL, "-weights_view");CHKERRQ(ierr);
1354cd81d59SMatthew G. Knepley   ierr = VecNorm(f,NORM_1,&norm);CHKERRQ(ierr);
1364cd81d59SMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_SELF,"Total number density = %g\n", norm);CHKERRQ(ierr);
1374cd81d59SMatthew G. Knepley   /* Cleanup */
1384cd81d59SMatthew G. Knepley   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr);
1394cd81d59SMatthew G. Knepley   ierr = MatDestroy(&M);CHKERRQ(ierr);
1404cd81d59SMatthew G. Knepley   ierr = MatDestroy(&M_p);CHKERRQ(ierr);
1414cd81d59SMatthew G. Knepley   ierr = VecDestroy(&rho);CHKERRQ(ierr);
1424cd81d59SMatthew G. Knepley   ierr = DMDestroy(&sw);CHKERRQ(ierr);
1434cd81d59SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
1444cd81d59SMatthew G. Knepley   ierr = DMDestroy(&crddm);CHKERRQ(ierr);
1454cd81d59SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
1464cd81d59SMatthew G. Knepley   ierr = PetscFinalize();
1474cd81d59SMatthew G. Knepley   return ierr;
1484cd81d59SMatthew G. Knepley }
1494cd81d59SMatthew G. Knepley 
1504cd81d59SMatthew G. Knepley /*TEST
1514cd81d59SMatthew G. Knepley 
1524cd81d59SMatthew G. Knepley   build:
1534cd81d59SMatthew G. Knepley     requires: !complex
1544cd81d59SMatthew G. Knepley 
1554cd81d59SMatthew G. Knepley   test:
1564cd81d59SMatthew G. Knepley     suffix: 0
1574cd81d59SMatthew G. Knepley     requires: double
158*30602db0SMatthew 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 -petscspace_degree 2 -ftop_ksp_type lsqr -ftop_pc_type none -dm_view
1594cd81d59SMatthew G. Knepley     filter: grep -v DM_ | grep -v atomic
1604cd81d59SMatthew G. Knepley 
1614cd81d59SMatthew G. Knepley TEST*/
162