xref: /petsc/src/dm/impls/swarm/tutorials/ex1.c (revision 4cd81d5996809478a54114bcbeca34ead748de02)
1*4cd81d59SMatthew G. Knepley static char help[] = "Example program demonstrating projection between particle and finite element spaces\n\n";
2*4cd81d59SMatthew G. Knepley 
3*4cd81d59SMatthew G. Knepley #include "petscdmplex.h"
4*4cd81d59SMatthew G. Knepley #include "petscds.h"
5*4cd81d59SMatthew G. Knepley #include "petscdmswarm.h"
6*4cd81d59SMatthew G. Knepley #include "petscksp.h"
7*4cd81d59SMatthew G. Knepley 
8*4cd81d59SMatthew G. Knepley static PetscErrorCode crd_func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
9*4cd81d59SMatthew G. Knepley {
10*4cd81d59SMatthew G. Knepley   PetscInt i;
11*4cd81d59SMatthew G. Knepley   PetscFunctionBeginUser;
12*4cd81d59SMatthew G. Knepley   for (i = 0; i < dim; ++i) u[i] = x[i];
13*4cd81d59SMatthew G. Knepley   PetscFunctionReturn(0);
14*4cd81d59SMatthew G. Knepley }
15*4cd81d59SMatthew G. Knepley 
16*4cd81d59SMatthew G. Knepley int main(int argc, char **argv)
17*4cd81d59SMatthew G. Knepley {
18*4cd81d59SMatthew G. Knepley   DM               dm, crddm, sw;
19*4cd81d59SMatthew G. Knepley   PetscFE          fe;
20*4cd81d59SMatthew G. Knepley   KSP              ksp;
21*4cd81d59SMatthew G. Knepley   Mat              M_p, M;
22*4cd81d59SMatthew G. Knepley   Vec              f, rho, rhs, crd_vec;
23*4cd81d59SMatthew G. Knepley   PetscInt         dim = 2, Nc = 1, timestep = 0, N, i, idx[3];
24*4cd81d59SMatthew G. Knepley   PetscInt         Np = 10, p, field = 0, zero = 0, bs, cells[] = {40, 20, 20};
25*4cd81d59SMatthew G. Knepley   PetscReal        time = 0.0,  norm;
26*4cd81d59SMatthew G. Knepley   PetscBool        removePoints = PETSC_TRUE, flg;
27*4cd81d59SMatthew G. Knepley   const PetscReal  *xx, *vv;
28*4cd81d59SMatthew G. Knepley   PetscReal        *wq, *coords, lo[] = {-1,-1,-1}, hi[]={1,1,1}, h[3];
29*4cd81d59SMatthew G. Knepley   PetscDataType    dtype;
30*4cd81d59SMatthew G. Knepley   PetscBool        interpolate = PETSC_TRUE;
31*4cd81d59SMatthew G. Knepley   PetscErrorCode   ierr;
32*4cd81d59SMatthew G. Knepley   PetscErrorCode   (*initu[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *);
33*4cd81d59SMatthew G. Knepley 
34*4cd81d59SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
35*4cd81d59SMatthew G. Knepley   /* get options */
36*4cd81d59SMatthew G. Knepley   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Options for Fokker-Plank collision operator", "none");CHKERRQ(ierr);
37*4cd81d59SMatthew G. Knepley   ierr = PetscOptionsRangeInt("-dim", "dim (2 or 3)", "ex1.c", dim, &dim, NULL,2,3);CHKERRQ(ierr);
38*4cd81d59SMatthew G. Knepley   i = dim;  ierr = PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells of grid", "ex1.c", cells, &i, &flg);CHKERRQ(ierr);
39*4cd81d59SMatthew G. Knepley   i = dim;  ierr = PetscOptionsRealArray("-dm_plex_box_lower", "Low corner of grid", "ex1.c", lo, &i, &flg);CHKERRQ(ierr);
40*4cd81d59SMatthew G. Knepley   i = dim;  ierr = PetscOptionsRealArray("-dm_plex_box_upper", "High corner of grid", "ex1.c", hi, &i, &flg);CHKERRQ(ierr);
41*4cd81d59SMatthew G. Knepley   for (i=0;i<dim;i++) {
42*4cd81d59SMatthew G. Knepley     h[i] = (hi[i] - lo[i])/cells[i];
43*4cd81d59SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF," lo = %g hi = %g n = %D h = %g\n",lo[i],hi[i],cells[i],h[i]);CHKERRQ(ierr);
44*4cd81d59SMatthew G. Knepley   }
45*4cd81d59SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
46*4cd81d59SMatthew G. Knepley   /* Create a mesh */
47*4cd81d59SMatthew G. Knepley   ierr = DMPlexCreateBoxMesh(PETSC_COMM_SELF, dim, PETSC_FALSE, NULL, NULL, NULL, NULL, interpolate, &dm);CHKERRQ(ierr);
48*4cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)dm, "Potential Grid");CHKERRQ(ierr);
49*4cd81d59SMatthew G. Knepley   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
50*4cd81d59SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe);CHKERRQ(ierr);
51*4cd81d59SMatthew G. Knepley   ierr = PetscFESetFromOptions(fe);CHKERRQ(ierr);
52*4cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)fe, "fe");CHKERRQ(ierr);
53*4cd81d59SMatthew G. Knepley   ierr = DMSetField(dm, field, NULL, (PetscObject)fe);CHKERRQ(ierr);
54*4cd81d59SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
55*4cd81d59SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
56*4cd81d59SMatthew G. Knepley   /* Create particle swarm */
57*4cd81d59SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_SELF, &sw);CHKERRQ(ierr);
58*4cd81d59SMatthew G. Knepley   ierr = DMSetType(sw, DMSWARM);CHKERRQ(ierr);
59*4cd81d59SMatthew G. Knepley   ierr = DMSetDimension(sw, dim);CHKERRQ(ierr);
60*4cd81d59SMatthew G. Knepley   ierr = DMSwarmSetType(sw, DMSWARM_PIC);CHKERRQ(ierr);
61*4cd81d59SMatthew G. Knepley   ierr = DMSwarmSetCellDM(sw, dm);CHKERRQ(ierr);
62*4cd81d59SMatthew G. Knepley   ierr = DMSwarmRegisterPetscDatatypeField(sw, "w_q", Nc, PETSC_SCALAR);CHKERRQ(ierr);
63*4cd81d59SMatthew G. Knepley   ierr = DMSwarmFinalizeFieldRegister(sw);CHKERRQ(ierr);
64*4cd81d59SMatthew G. Knepley   ierr = DMSwarmSetLocalSizes(sw, Np, zero);CHKERRQ(ierr);
65*4cd81d59SMatthew G. Knepley   ierr = DMSetFromOptions(sw);CHKERRQ(ierr);
66*4cd81d59SMatthew G. Knepley   ierr = DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
67*4cd81d59SMatthew G. Knepley   ierr = DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
68*4cd81d59SMatthew G. Knepley   for (p=0;p<Np;p++) {
69*4cd81d59SMatthew G. Knepley     coords[p*2+0]  = -PetscCosReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI);
70*4cd81d59SMatthew G. Knepley     coords[p*2+1] =   PetscSinReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI);
71*4cd81d59SMatthew G. Knepley     wq[p]          = 1.0;
72*4cd81d59SMatthew G. Knepley   }
73*4cd81d59SMatthew G. Knepley   ierr = DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
74*4cd81d59SMatthew G. Knepley   ierr = DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
75*4cd81d59SMatthew G. Knepley   ierr = DMSwarmMigrate(sw, removePoints);CHKERRQ(ierr);
76*4cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)sw, "Particle Grid");CHKERRQ(ierr);
77*4cd81d59SMatthew G. Knepley   ierr = DMViewFromOptions(sw, NULL, "-swarm_view");CHKERRQ(ierr);
78*4cd81d59SMatthew G. Knepley   /* Project particles to field */
79*4cd81d59SMatthew G. Knepley   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
80*4cd81d59SMatthew G. Knepley   ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr);
81*4cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &rho);CHKERRQ(ierr);
82*4cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)rho, "rho");CHKERRQ(ierr);
83*4cd81d59SMatthew G. Knepley   ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr);
84*4cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)f, "weights");CHKERRQ(ierr);
85*4cd81d59SMatthew G. Knepley   ierr = MatMultTranspose(M_p, f, rho);CHKERRQ(ierr);
86*4cd81d59SMatthew G. Knepley 
87*4cd81d59SMatthew G. Knepley   /* Visualize mesh field */
88*4cd81d59SMatthew G. Knepley   ierr = DMSetOutputSequenceNumber(dm, timestep, time);CHKERRQ(ierr);
89*4cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(rho, NULL, "-rho_view");CHKERRQ(ierr);
90*4cd81d59SMatthew G. Knepley 
91*4cd81d59SMatthew G. Knepley   /* create coordinate DM */
92*4cd81d59SMatthew G. Knepley   ierr = DMClone(dm, &crddm);CHKERRQ(ierr);
93*4cd81d59SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, PETSC_FALSE, "", PETSC_DECIDE, &fe);CHKERRQ(ierr);
94*4cd81d59SMatthew G. Knepley   ierr = PetscFESetFromOptions(fe);CHKERRQ(ierr);
95*4cd81d59SMatthew G. Knepley   ierr = DMSetField(crddm, field, NULL, (PetscObject)fe);CHKERRQ(ierr);
96*4cd81d59SMatthew G. Knepley   ierr = DMCreateDS(crddm);CHKERRQ(ierr);
97*4cd81d59SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
98*4cd81d59SMatthew G. Knepley   /* project coordiantes to vertices */
99*4cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(crddm, &crd_vec);CHKERRQ(ierr);
100*4cd81d59SMatthew G. Knepley   initu[0] = crd_func;
101*4cd81d59SMatthew G. Knepley   ierr = DMProjectFunction(crddm, 0.0, initu, NULL, INSERT_ALL_VALUES, crd_vec);CHKERRQ(ierr);
102*4cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(crd_vec, NULL, "-coord_view");CHKERRQ(ierr);
103*4cd81d59SMatthew G. Knepley   /* iterate over mesh data and get indices */
104*4cd81d59SMatthew G. Knepley   ierr = VecGetArrayRead(crd_vec,&xx);CHKERRQ(ierr);
105*4cd81d59SMatthew G. Knepley   ierr = VecGetArrayRead(rho,&vv);CHKERRQ(ierr);
106*4cd81d59SMatthew G. Knepley   ierr = VecGetLocalSize(rho,&N);CHKERRQ(ierr);
107*4cd81d59SMatthew G. Knepley   for (p=0;p<N;p++) {
108*4cd81d59SMatthew G. Knepley     for (i=0;i<dim;i++) idx[i] = (PetscInt)((xx[p*dim+i] - lo[i])/h[i] + 1.e-8);
109*4cd81d59SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF,"(%D,%D) = %g\n",idx[0],idx[1],vv[p]);CHKERRQ(ierr);
110*4cd81d59SMatthew G. Knepley     /* access grid data here */
111*4cd81d59SMatthew G. Knepley   }
112*4cd81d59SMatthew G. Knepley   ierr = VecRestoreArrayRead(crd_vec,&xx);CHKERRQ(ierr);
113*4cd81d59SMatthew G. Knepley   ierr = VecRestoreArrayRead(rho,&vv);CHKERRQ(ierr);
114*4cd81d59SMatthew G. Knepley   ierr = VecDestroy(&crd_vec);CHKERRQ(ierr);
115*4cd81d59SMatthew G. Knepley   /* Project field to particles */
116*4cd81d59SMatthew G. Knepley   /*   This gives f_p = M_p^+ M f */
117*4cd81d59SMatthew G. Knepley   ierr = DMCreateMassMatrix(dm, dm, &M);CHKERRQ(ierr);
118*4cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &rhs);CHKERRQ(ierr);
119*4cd81d59SMatthew G. Knepley   if (0) {
120*4cd81d59SMatthew G. Knepley     ierr = MatMult(M, rho, rhs);CHKERRQ(ierr);  /* this is what you would do for an FE solve */
121*4cd81d59SMatthew G. Knepley   } else {
122*4cd81d59SMatthew G. Knepley     ierr = VecCopy(rho, rhs);CHKERRQ(ierr); /* Indentity: M^1 M rho */
123*4cd81d59SMatthew G. Knepley   }
124*4cd81d59SMatthew G. Knepley   ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr);
125*4cd81d59SMatthew G. Knepley   ierr = KSPSetOptionsPrefix(ksp, "ftop_");CHKERRQ(ierr);
126*4cd81d59SMatthew G. Knepley   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
127*4cd81d59SMatthew G. Knepley   ierr = KSPSetOperators(ksp, M_p, M_p);CHKERRQ(ierr);
128*4cd81d59SMatthew G. Knepley   ierr = KSPSolveTranspose(ksp, rhs, f);CHKERRQ(ierr);
129*4cd81d59SMatthew G. Knepley   ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
130*4cd81d59SMatthew G. Knepley   ierr = VecDestroy(&rhs);CHKERRQ(ierr);
131*4cd81d59SMatthew G. Knepley 
132*4cd81d59SMatthew G. Knepley   /* Visualize particle field */
133*4cd81d59SMatthew G. Knepley   ierr = DMSetOutputSequenceNumber(sw, timestep, time);CHKERRQ(ierr);
134*4cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(f, NULL, "-weights_view");CHKERRQ(ierr);
135*4cd81d59SMatthew G. Knepley   ierr = VecNorm(f,NORM_1,&norm);CHKERRQ(ierr);
136*4cd81d59SMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_SELF,"Total number density = %g\n", norm);CHKERRQ(ierr);
137*4cd81d59SMatthew G. Knepley   /* Cleanup */
138*4cd81d59SMatthew G. Knepley   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr);
139*4cd81d59SMatthew G. Knepley   ierr = MatDestroy(&M);CHKERRQ(ierr);
140*4cd81d59SMatthew G. Knepley   ierr = MatDestroy(&M_p);CHKERRQ(ierr);
141*4cd81d59SMatthew G. Knepley   ierr = VecDestroy(&rho);CHKERRQ(ierr);
142*4cd81d59SMatthew G. Knepley   ierr = DMDestroy(&sw);CHKERRQ(ierr);
143*4cd81d59SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
144*4cd81d59SMatthew G. Knepley   ierr = DMDestroy(&crddm);CHKERRQ(ierr);
145*4cd81d59SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
146*4cd81d59SMatthew G. Knepley   ierr = PetscFinalize();
147*4cd81d59SMatthew G. Knepley   /* PetscFunctionReturn(0); */
148*4cd81d59SMatthew G. Knepley   return ierr;
149*4cd81d59SMatthew G. Knepley }
150*4cd81d59SMatthew G. Knepley 
151*4cd81d59SMatthew G. Knepley /*TEST
152*4cd81d59SMatthew G. Knepley 
153*4cd81d59SMatthew G. Knepley   build:
154*4cd81d59SMatthew G. Knepley     requires: !complex
155*4cd81d59SMatthew G. Knepley 
156*4cd81d59SMatthew G. Knepley   test:
157*4cd81d59SMatthew G. Knepley     suffix: 0
158*4cd81d59SMatthew G. Knepley     requires: double
159*4cd81d59SMatthew G. Knepley     args: -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 -swarm_view
160*4cd81d59SMatthew G. Knepley     filter: grep -v DM_ | grep -v atomic
161*4cd81d59SMatthew G. Knepley     filter_output: grep -v atomic
162*4cd81d59SMatthew G. Knepley 
163*4cd81d59SMatthew G. Knepley TEST*/
164