xref: /petsc/src/dm/impls/swarm/tests/ex6.c (revision c78c1075938cde1a2216af5a6ea940b358ba7f27)
1*c78c1075SJoseph Pusztay static char help[] = "Vlasov-Poisson example of central orbits\n";
2*c78c1075SJoseph Pusztay 
3*c78c1075SJoseph Pusztay /*
4*c78c1075SJoseph Pusztay   To visualize the orbit, we can used
5*c78c1075SJoseph Pusztay 
6*c78c1075SJoseph Pusztay     -ts_monitor_sp_swarm -ts_monitor_sp_swarm_retain -1 -ts_monitor_sp_swarm_phase 0 -draw_size 500,500
7*c78c1075SJoseph Pusztay 
8*c78c1075SJoseph Pusztay   and we probably want it to run fast and not check convergence
9*c78c1075SJoseph Pusztay 
10*c78c1075SJoseph Pusztay     -convest_num_refine 0 -ts_dt 0.01 -ts_max_steps 100 -ts_max_time 100 -output_step 10
11*c78c1075SJoseph Pusztay */
12*c78c1075SJoseph Pusztay 
13*c78c1075SJoseph Pusztay #include <petscts.h>
14*c78c1075SJoseph Pusztay #include <petscdmplex.h>
15*c78c1075SJoseph Pusztay #include <petscdmswarm.h>
16*c78c1075SJoseph Pusztay #include <petsc/private/dmpleximpl.h> /* For norm and dot */
17*c78c1075SJoseph Pusztay #include <petscfe.h>
18*c78c1075SJoseph Pusztay #include <petscds.h>
19*c78c1075SJoseph Pusztay #include <petsc/private/petscfeimpl.h> /* For interpolation */
20*c78c1075SJoseph Pusztay #include <petscksp.h>
21*c78c1075SJoseph Pusztay #include "petscsnes.h"
22*c78c1075SJoseph Pusztay 
23*c78c1075SJoseph Pusztay PETSC_EXTERN PetscErrorCode circleSingleX(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
24*c78c1075SJoseph Pusztay PETSC_EXTERN PetscErrorCode circleSingleV(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
25*c78c1075SJoseph Pusztay PETSC_EXTERN PetscErrorCode circleMultipleX(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
26*c78c1075SJoseph Pusztay PETSC_EXTERN PetscErrorCode circleMultipleV(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
27*c78c1075SJoseph Pusztay 
28*c78c1075SJoseph Pusztay const char *EMTypes[] = {"primal", "mixed", "coulomb", "none", "EMType", "EM_", NULL};
29*c78c1075SJoseph Pusztay typedef enum {
30*c78c1075SJoseph Pusztay   EM_PRIMAL,
31*c78c1075SJoseph Pusztay   EM_MIXED,
32*c78c1075SJoseph Pusztay   EM_COULOMB,
33*c78c1075SJoseph Pusztay   EM_NONE
34*c78c1075SJoseph Pusztay } EMType;
35*c78c1075SJoseph Pusztay 
36*c78c1075SJoseph Pusztay typedef struct {
37*c78c1075SJoseph Pusztay   PetscBool error;     /* Flag for printing the error */
38*c78c1075SJoseph Pusztay   PetscInt  ostep;     /* print the energy at each ostep time steps */
39*c78c1075SJoseph Pusztay   PetscReal timeScale; /* Nondimensionalizing time scale */
40*c78c1075SJoseph Pusztay   PetscReal sigma;     /* Linear charge per box length */
41*c78c1075SJoseph Pusztay   EMType    em;        /* Type of electrostatic model */
42*c78c1075SJoseph Pusztay   SNES      snes;
43*c78c1075SJoseph Pusztay } AppCtx;
44*c78c1075SJoseph Pusztay 
45*c78c1075SJoseph Pusztay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
46*c78c1075SJoseph Pusztay {
47*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
48*c78c1075SJoseph Pusztay   options->error     = PETSC_FALSE;
49*c78c1075SJoseph Pusztay   options->ostep     = 100;
50*c78c1075SJoseph Pusztay   options->timeScale = 1.0e-6;
51*c78c1075SJoseph Pusztay   options->sigma     = 1.;
52*c78c1075SJoseph Pusztay   options->em        = EM_COULOMB;
53*c78c1075SJoseph Pusztay 
54*c78c1075SJoseph Pusztay   PetscOptionsBegin(comm, "", "Central Orbit Options", "DMSWARM");
55*c78c1075SJoseph Pusztay   PetscCall(PetscOptionsBool("-error", "Flag to print the error", "ex6.c", options->error, &options->error, NULL));
56*c78c1075SJoseph Pusztay   PetscCall(PetscOptionsInt("-output_step", "Number of time steps between output", "ex6.c", options->ostep, &options->ostep, NULL));
57*c78c1075SJoseph Pusztay   PetscCall(PetscOptionsReal("-sigma", "Linear charge per box length", "ex6.c", options->sigma, &options->sigma, NULL));
58*c78c1075SJoseph Pusztay   PetscCall(PetscOptionsReal("-timeScale", "Nondimensionalizing time scale", "ex6.c", options->timeScale, &options->timeScale, NULL));
59*c78c1075SJoseph Pusztay   PetscCall(PetscOptionsEnum("-em_type", "Type of electrostatic solver", "ex6.c", EMTypes, (PetscEnum)options->em, (PetscEnum *)&options->em, NULL));
60*c78c1075SJoseph Pusztay   PetscOptionsEnd();
61*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
62*c78c1075SJoseph Pusztay }
63*c78c1075SJoseph Pusztay 
64*c78c1075SJoseph Pusztay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
65*c78c1075SJoseph Pusztay {
66*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
67*c78c1075SJoseph Pusztay   PetscCall(DMCreate(comm, dm));
68*c78c1075SJoseph Pusztay   PetscCall(DMSetType(*dm, DMPLEX));
69*c78c1075SJoseph Pusztay   PetscCall(DMSetFromOptions(*dm));
70*c78c1075SJoseph Pusztay   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
71*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
72*c78c1075SJoseph Pusztay }
73*c78c1075SJoseph Pusztay 
74*c78c1075SJoseph Pusztay static void laplacian_f1(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 f1[])
75*c78c1075SJoseph Pusztay {
76*c78c1075SJoseph Pusztay   PetscInt d;
77*c78c1075SJoseph Pusztay   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
78*c78c1075SJoseph Pusztay }
79*c78c1075SJoseph Pusztay 
80*c78c1075SJoseph Pusztay static void laplacian_g3(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 g3[])
81*c78c1075SJoseph Pusztay {
82*c78c1075SJoseph Pusztay   PetscInt d;
83*c78c1075SJoseph Pusztay   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
84*c78c1075SJoseph Pusztay }
85*c78c1075SJoseph Pusztay 
86*c78c1075SJoseph Pusztay /*
87*c78c1075SJoseph Pusztay    /  I   grad\ /q\ = /0\
88*c78c1075SJoseph Pusztay    \-div    0 / \u/   \f/
89*c78c1075SJoseph Pusztay */
90*c78c1075SJoseph Pusztay static void f0_q(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[])
91*c78c1075SJoseph Pusztay {
92*c78c1075SJoseph Pusztay   for (PetscInt c = 0; c < dim; ++c) f0[c] += u[uOff[0] + c];
93*c78c1075SJoseph Pusztay }
94*c78c1075SJoseph Pusztay 
95*c78c1075SJoseph Pusztay static void f1_q(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 f1[])
96*c78c1075SJoseph Pusztay {
97*c78c1075SJoseph Pusztay   for (PetscInt d = 0; d < dim; ++d) f1[d * dim + d] += u[uOff[1]];
98*c78c1075SJoseph Pusztay }
99*c78c1075SJoseph Pusztay 
100*c78c1075SJoseph Pusztay /* <t, q> */
101*c78c1075SJoseph Pusztay static void g0_qq(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[])
102*c78c1075SJoseph Pusztay {
103*c78c1075SJoseph Pusztay   PetscInt c;
104*c78c1075SJoseph Pusztay   for (c = 0; c < dim; ++c) g0[c * dim + c] += 1.0;
105*c78c1075SJoseph Pusztay }
106*c78c1075SJoseph Pusztay 
107*c78c1075SJoseph Pusztay static void g2_qu(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 g2[])
108*c78c1075SJoseph Pusztay {
109*c78c1075SJoseph Pusztay   for (PetscInt d = 0; d < dim; ++d) g2[d * dim + d] += 1.0;
110*c78c1075SJoseph Pusztay }
111*c78c1075SJoseph Pusztay 
112*c78c1075SJoseph Pusztay static void g1_uq(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 g1[])
113*c78c1075SJoseph Pusztay {
114*c78c1075SJoseph Pusztay   for (PetscInt d = 0; d < dim; ++d) g1[d * dim + d] += 1.0;
115*c78c1075SJoseph Pusztay }
116*c78c1075SJoseph Pusztay 
117*c78c1075SJoseph Pusztay static PetscErrorCode CreateFEM(DM dm, AppCtx *user)
118*c78c1075SJoseph Pusztay {
119*c78c1075SJoseph Pusztay   PetscFE   feu, feq;
120*c78c1075SJoseph Pusztay   PetscDS   ds;
121*c78c1075SJoseph Pusztay   PetscBool simplex;
122*c78c1075SJoseph Pusztay   PetscInt  dim;
123*c78c1075SJoseph Pusztay 
124*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
125*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(dm, &dim));
126*c78c1075SJoseph Pusztay   PetscCall(DMPlexIsSimplex(dm, &simplex));
127*c78c1075SJoseph Pusztay   if (user->em == EM_MIXED) {
128*c78c1075SJoseph Pusztay     DMLabel label;
129*c78c1075SJoseph Pusztay 
130*c78c1075SJoseph Pusztay     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", PETSC_DETERMINE, &feq));
131*c78c1075SJoseph Pusztay     PetscCall(PetscObjectSetName((PetscObject)feq, "field"));
132*c78c1075SJoseph Pusztay     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", PETSC_DETERMINE, &feu));
133*c78c1075SJoseph Pusztay     PetscCall(PetscObjectSetName((PetscObject)feu, "potential"));
134*c78c1075SJoseph Pusztay     PetscCall(PetscFECopyQuadrature(feq, feu));
135*c78c1075SJoseph Pusztay     PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq));
136*c78c1075SJoseph Pusztay     PetscCall(DMSetField(dm, 1, NULL, (PetscObject)feu));
137*c78c1075SJoseph Pusztay     PetscCall(DMCreateDS(dm));
138*c78c1075SJoseph Pusztay     PetscCall(PetscFEDestroy(&feu));
139*c78c1075SJoseph Pusztay     PetscCall(PetscFEDestroy(&feq));
140*c78c1075SJoseph Pusztay 
141*c78c1075SJoseph Pusztay     PetscCall(DMGetLabel(dm, "marker", &label));
142*c78c1075SJoseph Pusztay     PetscCall(DMGetDS(dm, &ds));
143*c78c1075SJoseph Pusztay     PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q));
144*c78c1075SJoseph Pusztay     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL));
145*c78c1075SJoseph Pusztay     PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qu, NULL));
146*c78c1075SJoseph Pusztay     PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_uq, NULL, NULL));
147*c78c1075SJoseph Pusztay   } else if (user->em == EM_PRIMAL) {
148*c78c1075SJoseph Pusztay     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, PETSC_DETERMINE, &feu));
149*c78c1075SJoseph Pusztay     PetscCall(PetscObjectSetName((PetscObject)feu, "potential"));
150*c78c1075SJoseph Pusztay     PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feu));
151*c78c1075SJoseph Pusztay     PetscCall(DMCreateDS(dm));
152*c78c1075SJoseph Pusztay     PetscCall(PetscFEDestroy(&feu));
153*c78c1075SJoseph Pusztay     PetscCall(DMGetDS(dm, &ds));
154*c78c1075SJoseph Pusztay     PetscCall(PetscDSSetResidual(ds, 0, NULL, laplacian_f1));
155*c78c1075SJoseph Pusztay     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, laplacian_g3));
156*c78c1075SJoseph Pusztay   }
157*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
158*c78c1075SJoseph Pusztay }
159*c78c1075SJoseph Pusztay 
160*c78c1075SJoseph Pusztay static PetscErrorCode CreatePoisson(DM dm, AppCtx *user)
161*c78c1075SJoseph Pusztay {
162*c78c1075SJoseph Pusztay   SNES         snes;
163*c78c1075SJoseph Pusztay   Mat          J;
164*c78c1075SJoseph Pusztay   MatNullSpace nullSpace;
165*c78c1075SJoseph Pusztay 
166*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
167*c78c1075SJoseph Pusztay   PetscCall(CreateFEM(dm, user));
168*c78c1075SJoseph Pusztay   PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes));
169*c78c1075SJoseph Pusztay   PetscCall(SNESSetOptionsPrefix(snes, "em_"));
170*c78c1075SJoseph Pusztay   PetscCall(SNESSetDM(snes, dm));
171*c78c1075SJoseph Pusztay   PetscCall(DMPlexSetSNESLocalFEM(dm, user, user, user));
172*c78c1075SJoseph Pusztay   PetscCall(SNESSetFromOptions(snes));
173*c78c1075SJoseph Pusztay 
174*c78c1075SJoseph Pusztay   PetscCall(DMCreateMatrix(dm, &J));
175*c78c1075SJoseph Pusztay   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_TRUE, 0, NULL, &nullSpace));
176*c78c1075SJoseph Pusztay   PetscCall(MatSetNullSpace(J, nullSpace));
177*c78c1075SJoseph Pusztay   PetscCall(MatNullSpaceDestroy(&nullSpace));
178*c78c1075SJoseph Pusztay   PetscCall(SNESSetJacobian(snes, J, J, NULL, NULL));
179*c78c1075SJoseph Pusztay   PetscCall(MatDestroy(&J));
180*c78c1075SJoseph Pusztay   user->snes = snes;
181*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
182*c78c1075SJoseph Pusztay }
183*c78c1075SJoseph Pusztay 
184*c78c1075SJoseph Pusztay static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw)
185*c78c1075SJoseph Pusztay {
186*c78c1075SJoseph Pusztay   PetscReal v0[1] = {1.};
187*c78c1075SJoseph Pusztay   PetscInt  dim;
188*c78c1075SJoseph Pusztay 
189*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
190*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(dm, &dim));
191*c78c1075SJoseph Pusztay   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
192*c78c1075SJoseph Pusztay   PetscCall(DMSetType(*sw, DMSWARM));
193*c78c1075SJoseph Pusztay   PetscCall(DMSetDimension(*sw, dim));
194*c78c1075SJoseph Pusztay   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
195*c78c1075SJoseph Pusztay   PetscCall(DMSwarmSetCellDM(*sw, dm));
196*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
197*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "velocity", dim, PETSC_REAL));
198*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
199*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initCoordinates", dim, PETSC_REAL));
200*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "initVelocity", dim, PETSC_REAL));
201*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "E_field", dim, PETSC_REAL));
202*c78c1075SJoseph Pusztay   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
203*c78c1075SJoseph Pusztay   PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
204*c78c1075SJoseph Pusztay   PetscCall(DMSwarmInitializeCoordinates(*sw));
205*c78c1075SJoseph Pusztay   PetscCall(DMSwarmInitializeVelocitiesFromOptions(*sw, v0));
206*c78c1075SJoseph Pusztay   PetscCall(DMSetFromOptions(*sw));
207*c78c1075SJoseph Pusztay   PetscCall(DMSetApplicationContext(*sw, user));
208*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
209*c78c1075SJoseph Pusztay   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
210*c78c1075SJoseph Pusztay   {
211*c78c1075SJoseph Pusztay     Vec gc, gc0, gv, gv0;
212*c78c1075SJoseph Pusztay 
213*c78c1075SJoseph Pusztay     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
214*c78c1075SJoseph Pusztay     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initCoordinates", &gc0));
215*c78c1075SJoseph Pusztay     PetscCall(VecCopy(gc, gc0));
216*c78c1075SJoseph Pusztay     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, DMSwarmPICField_coor, &gc));
217*c78c1075SJoseph Pusztay     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initCoordinates", &gc0));
218*c78c1075SJoseph Pusztay     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "velocity", &gv));
219*c78c1075SJoseph Pusztay     PetscCall(DMSwarmCreateGlobalVectorFromField(*sw, "initVelocity", &gv0));
220*c78c1075SJoseph Pusztay     PetscCall(VecCopy(gv, gv0));
221*c78c1075SJoseph Pusztay     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "velocity", &gv));
222*c78c1075SJoseph Pusztay     PetscCall(DMSwarmDestroyGlobalVectorFromField(*sw, "initVelocity", &gv0));
223*c78c1075SJoseph Pusztay   }
224*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
225*c78c1075SJoseph Pusztay }
226*c78c1075SJoseph Pusztay 
227*c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles_Coulomb(SNES snes, DM sw, PetscReal E[])
228*c78c1075SJoseph Pusztay {
229*c78c1075SJoseph Pusztay   PetscReal  *coords;
230*c78c1075SJoseph Pusztay   PetscInt    dim, d, Np, p, q;
231*c78c1075SJoseph Pusztay   PetscMPIInt size;
232*c78c1075SJoseph Pusztay 
233*c78c1075SJoseph Pusztay   PetscFunctionBegin;
234*c78c1075SJoseph Pusztay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)snes), &size));
235*c78c1075SJoseph Pusztay   PetscCheck(size == 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Coulomb code only works in serial");
236*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
237*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetLocalSize(sw, &Np));
238*c78c1075SJoseph Pusztay 
239*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
240*c78c1075SJoseph Pusztay   for (p = 0; p < Np; ++p) {
241*c78c1075SJoseph Pusztay     PetscReal *pcoord = &coords[p * dim];
242*c78c1075SJoseph Pusztay     PetscReal *pE     = &E[p * dim];
243*c78c1075SJoseph Pusztay     /* Calculate field at particle p due to particle q */
244*c78c1075SJoseph Pusztay     for (q = 0; q < Np; ++q) {
245*c78c1075SJoseph Pusztay       PetscReal *qcoord = &coords[q * dim];
246*c78c1075SJoseph Pusztay       PetscReal  rpq[3], r;
247*c78c1075SJoseph Pusztay 
248*c78c1075SJoseph Pusztay       if (p == q) continue;
249*c78c1075SJoseph Pusztay       for (d = 0; d < dim; ++d) rpq[d] = pcoord[d] - qcoord[d];
250*c78c1075SJoseph Pusztay       r = DMPlex_NormD_Internal(dim, rpq);
251*c78c1075SJoseph Pusztay       for (d = 0; d < dim; ++d) pE[d] += rpq[d] / PetscPowRealInt(r, 3);
252*c78c1075SJoseph Pusztay     }
253*c78c1075SJoseph Pusztay   }
254*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
255*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
256*c78c1075SJoseph Pusztay }
257*c78c1075SJoseph Pusztay 
258*c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles_Primal(SNES snes, DM sw, PetscReal E[])
259*c78c1075SJoseph Pusztay {
260*c78c1075SJoseph Pusztay   DM         dm;
261*c78c1075SJoseph Pusztay   PetscDS    ds;
262*c78c1075SJoseph Pusztay   PetscFE    fe;
263*c78c1075SJoseph Pusztay   Mat        M_p;
264*c78c1075SJoseph Pusztay   Vec        phi, locPhi, rho, f;
265*c78c1075SJoseph Pusztay   PetscReal *coords, chargeTol = 1e-13;
266*c78c1075SJoseph Pusztay   PetscInt   dim, d, cStart, cEnd, c, Np;
267*c78c1075SJoseph Pusztay 
268*c78c1075SJoseph Pusztay   PetscFunctionBegin;
269*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
270*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetLocalSize(sw, &Np));
271*c78c1075SJoseph Pusztay 
272*c78c1075SJoseph Pusztay   /* Create the charges rho */
273*c78c1075SJoseph Pusztay   PetscCall(SNESGetDM(snes, &dm));
274*c78c1075SJoseph Pusztay   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
275*c78c1075SJoseph Pusztay   PetscCall(DMGetGlobalVector(dm, &rho));
276*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
277*c78c1075SJoseph Pusztay   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
278*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)f, "particle weight"));
279*c78c1075SJoseph Pusztay   PetscCall(MatMultTranspose(M_p, f, rho));
280*c78c1075SJoseph Pusztay   PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
281*c78c1075SJoseph Pusztay   PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
282*c78c1075SJoseph Pusztay   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
283*c78c1075SJoseph Pusztay   PetscCall(MatDestroy(&M_p));
284*c78c1075SJoseph Pusztay   {
285*c78c1075SJoseph Pusztay     PetscScalar sum;
286*c78c1075SJoseph Pusztay     PetscInt    n;
287*c78c1075SJoseph Pusztay     PetscReal   phi_0 = 1.; /* (sigma*sigma*sigma)*(timeScale*timeScale)/(m_e*q_e*epsi_0)*/
288*c78c1075SJoseph Pusztay 
289*c78c1075SJoseph Pusztay     /* Remove constant from rho */
290*c78c1075SJoseph Pusztay     PetscCall(VecGetSize(rho, &n));
291*c78c1075SJoseph Pusztay     PetscCall(VecSum(rho, &sum));
292*c78c1075SJoseph Pusztay     PetscCall(VecShift(rho, -sum / n));
293*c78c1075SJoseph Pusztay     PetscCall(VecSum(rho, &sum));
294*c78c1075SJoseph Pusztay     PetscCheck(PetscAbsScalar(sum) < chargeTol, PetscObjectComm((PetscObject)sw), PETSC_ERR_PLIB, "Charge should have no DC component: %g", (double)PetscRealPart(sum));
295*c78c1075SJoseph Pusztay     /* Nondimensionalize rho */
296*c78c1075SJoseph Pusztay     PetscCall(VecScale(rho, phi_0));
297*c78c1075SJoseph Pusztay   }
298*c78c1075SJoseph Pusztay   PetscCall(VecViewFromOptions(rho, NULL, "-poisson_rho_view"));
299*c78c1075SJoseph Pusztay 
300*c78c1075SJoseph Pusztay   PetscCall(DMGetGlobalVector(dm, &phi));
301*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)phi, "potential"));
302*c78c1075SJoseph Pusztay   PetscCall(VecSet(phi, 0.0));
303*c78c1075SJoseph Pusztay   PetscCall(SNESSolve(snes, rho, phi));
304*c78c1075SJoseph Pusztay   PetscCall(DMRestoreGlobalVector(dm, &rho));
305*c78c1075SJoseph Pusztay   PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
306*c78c1075SJoseph Pusztay 
307*c78c1075SJoseph Pusztay   PetscCall(DMGetLocalVector(dm, &locPhi));
308*c78c1075SJoseph Pusztay   PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi));
309*c78c1075SJoseph Pusztay   PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi));
310*c78c1075SJoseph Pusztay   PetscCall(DMRestoreGlobalVector(dm, &phi));
311*c78c1075SJoseph Pusztay 
312*c78c1075SJoseph Pusztay   PetscCall(DMGetDS(dm, &ds));
313*c78c1075SJoseph Pusztay   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
314*c78c1075SJoseph Pusztay   PetscCall(DMSwarmSortGetAccess(sw));
315*c78c1075SJoseph Pusztay   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
316*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
317*c78c1075SJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
318*c78c1075SJoseph Pusztay     PetscTabulation tab;
319*c78c1075SJoseph Pusztay     PetscScalar    *clPhi = NULL;
320*c78c1075SJoseph Pusztay     PetscReal      *pcoord, *refcoord;
321*c78c1075SJoseph Pusztay     PetscReal       v[3], J[9], invJ[9], detJ;
322*c78c1075SJoseph Pusztay     PetscInt       *points;
323*c78c1075SJoseph Pusztay     PetscInt        Ncp, cp;
324*c78c1075SJoseph Pusztay 
325*c78c1075SJoseph Pusztay     PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points));
326*c78c1075SJoseph Pusztay     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
327*c78c1075SJoseph Pusztay     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
328*c78c1075SJoseph Pusztay     for (cp = 0; cp < Ncp; ++cp)
329*c78c1075SJoseph Pusztay       for (d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d];
330*c78c1075SJoseph Pusztay     PetscCall(DMPlexCoordinatesToReference(dm, c, Ncp, pcoord, refcoord));
331*c78c1075SJoseph Pusztay     PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
332*c78c1075SJoseph Pusztay     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v, J, invJ, &detJ));
333*c78c1075SJoseph Pusztay     PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi));
334*c78c1075SJoseph Pusztay     for (cp = 0; cp < Ncp; ++cp) {
335*c78c1075SJoseph Pusztay       const PetscReal *basisDer = tab->T[1];
336*c78c1075SJoseph Pusztay       const PetscInt   p        = points[cp];
337*c78c1075SJoseph Pusztay 
338*c78c1075SJoseph Pusztay       for (d = 0; d < dim; ++d) E[p * dim + d] = 0.;
339*c78c1075SJoseph Pusztay       PetscCall(PetscFEFreeInterpolateGradient_Static(fe, basisDer, clPhi, dim, invJ, NULL, cp, &E[p * dim]));
340*c78c1075SJoseph Pusztay       for (d = 0; d < dim; ++d) E[p * dim + d] *= -1.0;
341*c78c1075SJoseph Pusztay     }
342*c78c1075SJoseph Pusztay     PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi));
343*c78c1075SJoseph Pusztay     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
344*c78c1075SJoseph Pusztay     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
345*c78c1075SJoseph Pusztay     PetscCall(PetscTabulationDestroy(&tab));
346*c78c1075SJoseph Pusztay     PetscCall(PetscFree(points));
347*c78c1075SJoseph Pusztay   }
348*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
349*c78c1075SJoseph Pusztay   PetscCall(DMSwarmSortRestoreAccess(sw));
350*c78c1075SJoseph Pusztay   PetscCall(DMRestoreLocalVector(dm, &locPhi));
351*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
352*c78c1075SJoseph Pusztay }
353*c78c1075SJoseph Pusztay 
354*c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles_Mixed(SNES snes, DM sw, PetscReal E[])
355*c78c1075SJoseph Pusztay {
356*c78c1075SJoseph Pusztay   DM              dm, potential_dm;
357*c78c1075SJoseph Pusztay   IS              potential_IS;
358*c78c1075SJoseph Pusztay   PetscDS         ds;
359*c78c1075SJoseph Pusztay   PetscFE         fe;
360*c78c1075SJoseph Pusztay   PetscFEGeom     feGeometry;
361*c78c1075SJoseph Pusztay   Mat             M_p;
362*c78c1075SJoseph Pusztay   Vec             phi, locPhi, rho, f, temp_rho;
363*c78c1075SJoseph Pusztay   PetscQuadrature q;
364*c78c1075SJoseph Pusztay   PetscReal      *coords, chargeTol = 1e-13;
365*c78c1075SJoseph Pusztay   PetscInt        dim, d, cStart, cEnd, c, Np, fields = 1;
366*c78c1075SJoseph Pusztay 
367*c78c1075SJoseph Pusztay   PetscFunctionBegin;
368*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
369*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetLocalSize(sw, &Np));
370*c78c1075SJoseph Pusztay 
371*c78c1075SJoseph Pusztay   /* Create the charges rho */
372*c78c1075SJoseph Pusztay   PetscCall(SNESGetDM(snes, &dm));
373*c78c1075SJoseph Pusztay   PetscCall(DMGetGlobalVector(dm, &rho));
374*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
375*c78c1075SJoseph Pusztay 
376*c78c1075SJoseph Pusztay   PetscCall(DMCreateSubDM(dm, 1, &fields, &potential_IS, &potential_dm));
377*c78c1075SJoseph Pusztay   PetscCall(DMCreateMassMatrix(sw, potential_dm, &M_p));
378*c78c1075SJoseph Pusztay   PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
379*c78c1075SJoseph Pusztay   PetscCall(DMGetGlobalVector(potential_dm, &temp_rho));
380*c78c1075SJoseph Pusztay   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
381*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)f, "particle weight"));
382*c78c1075SJoseph Pusztay   PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
383*c78c1075SJoseph Pusztay   PetscCall(MatMultTranspose(M_p, f, temp_rho));
384*c78c1075SJoseph Pusztay   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
385*c78c1075SJoseph Pusztay   PetscCall(MatDestroy(&M_p));
386*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
387*c78c1075SJoseph Pusztay   PetscCall(VecViewFromOptions(rho, NULL, "-poisson_rho_view"));
388*c78c1075SJoseph Pusztay   PetscCall(VecISCopy(rho, potential_IS, SCATTER_FORWARD, temp_rho));
389*c78c1075SJoseph Pusztay   PetscCall(DMRestoreGlobalVector(potential_dm, &temp_rho));
390*c78c1075SJoseph Pusztay   PetscCall(DMDestroy(&potential_dm));
391*c78c1075SJoseph Pusztay   PetscCall(ISDestroy(&potential_IS));
392*c78c1075SJoseph Pusztay   {
393*c78c1075SJoseph Pusztay     PetscScalar sum;
394*c78c1075SJoseph Pusztay     PetscInt    n;
395*c78c1075SJoseph Pusztay     PetscReal   phi_0 = 1.; /*(sigma*sigma*sigma)*(timeScale*timeScale)/(m_e*q_e*epsi_0);*/
396*c78c1075SJoseph Pusztay 
397*c78c1075SJoseph Pusztay     /*Remove constant from rho*/
398*c78c1075SJoseph Pusztay     PetscCall(VecGetSize(rho, &n));
399*c78c1075SJoseph Pusztay     PetscCall(VecSum(rho, &sum));
400*c78c1075SJoseph Pusztay     PetscCall(VecShift(rho, -sum / n));
401*c78c1075SJoseph Pusztay     PetscCall(VecSum(rho, &sum));
402*c78c1075SJoseph Pusztay     PetscCheck(PetscAbsScalar(sum) < chargeTol, PetscObjectComm((PetscObject)sw), PETSC_ERR_PLIB, "Charge should have no DC component: %g", (double)PetscRealPart(sum));
403*c78c1075SJoseph Pusztay     /* Nondimensionalize rho */
404*c78c1075SJoseph Pusztay     PetscCall(VecScale(rho, phi_0));
405*c78c1075SJoseph Pusztay   }
406*c78c1075SJoseph Pusztay   PetscCall(DMGetGlobalVector(dm, &phi));
407*c78c1075SJoseph Pusztay   PetscCall(PetscObjectSetName((PetscObject)phi, "potential"));
408*c78c1075SJoseph Pusztay   PetscCall(VecSet(phi, 0.0));
409*c78c1075SJoseph Pusztay   PetscCall(SNESSolve(snes, NULL, phi));
410*c78c1075SJoseph Pusztay   PetscCall(DMRestoreGlobalVector(dm, &rho));
411*c78c1075SJoseph Pusztay   PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
412*c78c1075SJoseph Pusztay 
413*c78c1075SJoseph Pusztay   PetscCall(DMGetLocalVector(dm, &locPhi));
414*c78c1075SJoseph Pusztay   PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi));
415*c78c1075SJoseph Pusztay   PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi));
416*c78c1075SJoseph Pusztay   PetscCall(DMRestoreGlobalVector(dm, &phi));
417*c78c1075SJoseph Pusztay 
418*c78c1075SJoseph Pusztay   PetscCall(DMGetDS(dm, &ds));
419*c78c1075SJoseph Pusztay   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
420*c78c1075SJoseph Pusztay   PetscCall(DMSwarmSortGetAccess(sw));
421*c78c1075SJoseph Pusztay   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
422*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
423*c78c1075SJoseph Pusztay   for (c = cStart; c < cEnd; ++c) {
424*c78c1075SJoseph Pusztay     PetscTabulation tab;
425*c78c1075SJoseph Pusztay     PetscScalar    *clPhi = NULL;
426*c78c1075SJoseph Pusztay     PetscReal      *pcoord, *refcoord;
427*c78c1075SJoseph Pusztay     PetscReal       v[3], J[9], invJ[9], detJ;
428*c78c1075SJoseph Pusztay     PetscInt       *points;
429*c78c1075SJoseph Pusztay     PetscInt        Ncp, cp;
430*c78c1075SJoseph Pusztay 
431*c78c1075SJoseph Pusztay     PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points));
432*c78c1075SJoseph Pusztay     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
433*c78c1075SJoseph Pusztay     PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
434*c78c1075SJoseph Pusztay     for (cp = 0; cp < Ncp; ++cp)
435*c78c1075SJoseph Pusztay       for (d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d];
436*c78c1075SJoseph Pusztay     PetscCall(DMPlexCoordinatesToReference(dm, c, Ncp, pcoord, refcoord));
437*c78c1075SJoseph Pusztay     PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
438*c78c1075SJoseph Pusztay     PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v, J, invJ, &detJ));
439*c78c1075SJoseph Pusztay     PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi));
440*c78c1075SJoseph Pusztay     for (cp = 0; cp < Ncp; ++cp) {
441*c78c1075SJoseph Pusztay       const PetscInt p = points[cp];
442*c78c1075SJoseph Pusztay 
443*c78c1075SJoseph Pusztay       for (d = 0; d < dim; ++d) E[p * dim + d] = 0.;
444*c78c1075SJoseph Pusztay       PetscCall(PetscFEGetQuadrature(fe, &q));
445*c78c1075SJoseph Pusztay       PetscCall(PetscFECreateCellGeometry(fe, q, &feGeometry));
446*c78c1075SJoseph Pusztay       PetscCall(PetscFEInterpolateAtPoints_Static(fe, tab, clPhi, &feGeometry, cp, &E[p * dim]));
447*c78c1075SJoseph Pusztay       PetscCall(PetscFEDestroyCellGeometry(fe, &feGeometry));
448*c78c1075SJoseph Pusztay     }
449*c78c1075SJoseph Pusztay     PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi));
450*c78c1075SJoseph Pusztay     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
451*c78c1075SJoseph Pusztay     PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
452*c78c1075SJoseph Pusztay     PetscCall(PetscTabulationDestroy(&tab));
453*c78c1075SJoseph Pusztay     PetscCall(PetscFree(points));
454*c78c1075SJoseph Pusztay   }
455*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
456*c78c1075SJoseph Pusztay   PetscCall(DMSwarmSortRestoreAccess(sw));
457*c78c1075SJoseph Pusztay   PetscCall(DMRestoreLocalVector(dm, &locPhi));
458*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
459*c78c1075SJoseph Pusztay }
460*c78c1075SJoseph Pusztay 
461*c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles(SNES snes, DM sw, PetscReal E[])
462*c78c1075SJoseph Pusztay {
463*c78c1075SJoseph Pusztay   AppCtx  *ctx;
464*c78c1075SJoseph Pusztay   PetscInt dim, Np;
465*c78c1075SJoseph Pusztay 
466*c78c1075SJoseph Pusztay   PetscFunctionBegin;
467*c78c1075SJoseph Pusztay   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
468*c78c1075SJoseph Pusztay   PetscValidHeaderSpecific(sw, DM_CLASSID, 2);
469*c78c1075SJoseph Pusztay   PetscValidRealPointer(E, 3);
470*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
471*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetLocalSize(sw, &Np));
472*c78c1075SJoseph Pusztay   PetscCall(DMGetApplicationContext(sw, &ctx));
473*c78c1075SJoseph Pusztay   PetscCall(PetscArrayzero(E, Np * dim));
474*c78c1075SJoseph Pusztay 
475*c78c1075SJoseph Pusztay   switch (ctx->em) {
476*c78c1075SJoseph Pusztay   case EM_PRIMAL:
477*c78c1075SJoseph Pusztay     PetscCall(ComputeFieldAtParticles_Primal(snes, sw, E));
478*c78c1075SJoseph Pusztay     break;
479*c78c1075SJoseph Pusztay   case EM_COULOMB:
480*c78c1075SJoseph Pusztay     PetscCall(ComputeFieldAtParticles_Coulomb(snes, sw, E));
481*c78c1075SJoseph Pusztay     break;
482*c78c1075SJoseph Pusztay   case EM_MIXED:
483*c78c1075SJoseph Pusztay     PetscCall(ComputeFieldAtParticles_Mixed(snes, sw, E));
484*c78c1075SJoseph Pusztay     break;
485*c78c1075SJoseph Pusztay   case EM_NONE:
486*c78c1075SJoseph Pusztay     break;
487*c78c1075SJoseph Pusztay   default:
488*c78c1075SJoseph Pusztay     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No solver for electrostatic model %s", EMTypes[ctx->em]);
489*c78c1075SJoseph Pusztay   }
490*c78c1075SJoseph Pusztay 
491*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
492*c78c1075SJoseph Pusztay }
493*c78c1075SJoseph Pusztay 
494*c78c1075SJoseph Pusztay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, void *ctx)
495*c78c1075SJoseph Pusztay {
496*c78c1075SJoseph Pusztay   DM                 sw;
497*c78c1075SJoseph Pusztay   SNES               snes = ((AppCtx *)ctx)->snes;
498*c78c1075SJoseph Pusztay   const PetscReal   *coords, *vel;
499*c78c1075SJoseph Pusztay   const PetscScalar *u;
500*c78c1075SJoseph Pusztay   PetscScalar       *g;
501*c78c1075SJoseph Pusztay   PetscReal         *E;
502*c78c1075SJoseph Pusztay   PetscInt           dim, d, Np, p;
503*c78c1075SJoseph Pusztay 
504*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
505*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
506*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
507*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
508*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
509*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
510*c78c1075SJoseph Pusztay   PetscCall(VecGetLocalSize(U, &Np));
511*c78c1075SJoseph Pusztay   PetscCall(VecGetArrayRead(U, &u));
512*c78c1075SJoseph Pusztay   PetscCall(VecGetArray(G, &g));
513*c78c1075SJoseph Pusztay 
514*c78c1075SJoseph Pusztay   PetscCall(ComputeFieldAtParticles(snes, sw, E));
515*c78c1075SJoseph Pusztay 
516*c78c1075SJoseph Pusztay   Np /= 2 * dim;
517*c78c1075SJoseph Pusztay   for (p = 0; p < Np; ++p) {
518*c78c1075SJoseph Pusztay     const PetscReal x0    = coords[p * dim + 0];
519*c78c1075SJoseph Pusztay     const PetscReal vy0   = vel[p * dim + 1];
520*c78c1075SJoseph Pusztay     const PetscReal omega = vy0 / x0;
521*c78c1075SJoseph Pusztay 
522*c78c1075SJoseph Pusztay     for (d = 0; d < dim; ++d) {
523*c78c1075SJoseph Pusztay       g[(p * 2 + 0) * dim + d] = u[(p * 2 + 1) * dim + d];
524*c78c1075SJoseph Pusztay       g[(p * 2 + 1) * dim + d] = E[p * dim + d] - PetscSqr(omega) * u[(p * 2 + 0) * dim + d];
525*c78c1075SJoseph Pusztay     }
526*c78c1075SJoseph Pusztay   }
527*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
528*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
529*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
530*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArrayRead(U, &u));
531*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArray(G, &g));
532*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
533*c78c1075SJoseph Pusztay }
534*c78c1075SJoseph Pusztay 
535*c78c1075SJoseph Pusztay /* J_{ij} = dF_i/dx_j
536*c78c1075SJoseph Pusztay    J_p = (  0   1)
537*c78c1075SJoseph Pusztay          (-w^2  0)
538*c78c1075SJoseph Pusztay    TODO Now there is another term with w^2 from the electric field. I think we will need to invert the operator.
539*c78c1075SJoseph Pusztay         Perhaps we can approximate the Jacobian using only the cellwise P-P gradient from Coulomb
540*c78c1075SJoseph Pusztay */
541*c78c1075SJoseph Pusztay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat P, void *ctx)
542*c78c1075SJoseph Pusztay {
543*c78c1075SJoseph Pusztay   DM               sw;
544*c78c1075SJoseph Pusztay   const PetscReal *coords, *vel;
545*c78c1075SJoseph Pusztay   PetscInt         dim, d, Np, p, rStart;
546*c78c1075SJoseph Pusztay 
547*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
548*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
549*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
550*c78c1075SJoseph Pusztay   PetscCall(VecGetLocalSize(U, &Np));
551*c78c1075SJoseph Pusztay   PetscCall(MatGetOwnershipRange(J, &rStart, NULL));
552*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
553*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
554*c78c1075SJoseph Pusztay   Np /= 2 * dim;
555*c78c1075SJoseph Pusztay   for (p = 0; p < Np; ++p) {
556*c78c1075SJoseph Pusztay     const PetscReal x0      = coords[p * dim + 0];
557*c78c1075SJoseph Pusztay     const PetscReal vy0     = vel[p * dim + 1];
558*c78c1075SJoseph Pusztay     const PetscReal omega   = vy0 / x0;
559*c78c1075SJoseph Pusztay     PetscScalar     vals[4] = {0., 1., -PetscSqr(omega), 0.};
560*c78c1075SJoseph Pusztay 
561*c78c1075SJoseph Pusztay     for (d = 0; d < dim; ++d) {
562*c78c1075SJoseph Pusztay       const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart};
563*c78c1075SJoseph Pusztay       PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES));
564*c78c1075SJoseph Pusztay     }
565*c78c1075SJoseph Pusztay   }
566*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
567*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
568*c78c1075SJoseph Pusztay   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
569*c78c1075SJoseph Pusztay   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
570*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
571*c78c1075SJoseph Pusztay }
572*c78c1075SJoseph Pusztay 
573*c78c1075SJoseph Pusztay static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx)
574*c78c1075SJoseph Pusztay {
575*c78c1075SJoseph Pusztay   DM                 sw;
576*c78c1075SJoseph Pusztay   const PetscScalar *v;
577*c78c1075SJoseph Pusztay   PetscScalar       *xres;
578*c78c1075SJoseph Pusztay   PetscInt           Np, p, dim, d;
579*c78c1075SJoseph Pusztay 
580*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
581*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
582*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
583*c78c1075SJoseph Pusztay   PetscCall(VecGetLocalSize(Xres, &Np));
584*c78c1075SJoseph Pusztay   Np /= dim;
585*c78c1075SJoseph Pusztay   PetscCall(VecGetArrayRead(V, &v));
586*c78c1075SJoseph Pusztay   PetscCall(VecGetArray(Xres, &xres));
587*c78c1075SJoseph Pusztay   for (p = 0; p < Np; ++p) {
588*c78c1075SJoseph Pusztay     for (d = 0; d < dim; ++d) { xres[p * dim + d] = v[p * dim + d]; }
589*c78c1075SJoseph Pusztay   }
590*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArrayRead(V, &v));
591*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArray(Xres, &xres));
592*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
593*c78c1075SJoseph Pusztay }
594*c78c1075SJoseph Pusztay 
595*c78c1075SJoseph Pusztay static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx)
596*c78c1075SJoseph Pusztay {
597*c78c1075SJoseph Pusztay   DM                 sw;
598*c78c1075SJoseph Pusztay   SNES               snes = ((AppCtx *)ctx)->snes;
599*c78c1075SJoseph Pusztay   const PetscScalar *x;
600*c78c1075SJoseph Pusztay   const PetscReal   *coords, *vel;
601*c78c1075SJoseph Pusztay   PetscScalar       *vres;
602*c78c1075SJoseph Pusztay   PetscReal         *E;
603*c78c1075SJoseph Pusztay   PetscInt           Np, p, dim, d;
604*c78c1075SJoseph Pusztay 
605*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
606*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
607*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
608*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
609*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
610*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
611*c78c1075SJoseph Pusztay   PetscCall(VecGetLocalSize(Vres, &Np));
612*c78c1075SJoseph Pusztay   PetscCall(VecGetArrayRead(X, &x));
613*c78c1075SJoseph Pusztay   PetscCall(VecGetArray(Vres, &vres));
614*c78c1075SJoseph Pusztay   PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension must be 2");
615*c78c1075SJoseph Pusztay 
616*c78c1075SJoseph Pusztay   PetscCall(ComputeFieldAtParticles(snes, sw, E));
617*c78c1075SJoseph Pusztay 
618*c78c1075SJoseph Pusztay   Np /= dim;
619*c78c1075SJoseph Pusztay   for (p = 0; p < Np; ++p) {
620*c78c1075SJoseph Pusztay     const PetscReal x0    = coords[p * dim + 0];
621*c78c1075SJoseph Pusztay     const PetscReal vy0   = vel[p * dim + 1];
622*c78c1075SJoseph Pusztay     const PetscReal omega = vy0 / x0;
623*c78c1075SJoseph Pusztay 
624*c78c1075SJoseph Pusztay     for (d = 0; d < dim; ++d) vres[p * dim + d] = E[p * dim + d] - PetscSqr(omega) * x[p * dim + d];
625*c78c1075SJoseph Pusztay   }
626*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArrayRead(X, &x));
627*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArray(Vres, &vres));
628*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
629*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
630*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
631*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
632*c78c1075SJoseph Pusztay }
633*c78c1075SJoseph Pusztay 
634*c78c1075SJoseph Pusztay static PetscErrorCode CreateSolution(TS ts)
635*c78c1075SJoseph Pusztay {
636*c78c1075SJoseph Pusztay   DM       sw;
637*c78c1075SJoseph Pusztay   Vec      u;
638*c78c1075SJoseph Pusztay   PetscInt dim, Np;
639*c78c1075SJoseph Pusztay 
640*c78c1075SJoseph Pusztay   PetscFunctionBegin;
641*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
642*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
643*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetLocalSize(sw, &Np));
644*c78c1075SJoseph Pusztay   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
645*c78c1075SJoseph Pusztay   PetscCall(VecSetBlockSize(u, dim));
646*c78c1075SJoseph Pusztay   PetscCall(VecSetSizes(u, 2 * Np * dim, PETSC_DECIDE));
647*c78c1075SJoseph Pusztay   PetscCall(VecSetUp(u));
648*c78c1075SJoseph Pusztay   PetscCall(TSSetSolution(ts, u));
649*c78c1075SJoseph Pusztay   PetscCall(VecDestroy(&u));
650*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
651*c78c1075SJoseph Pusztay }
652*c78c1075SJoseph Pusztay 
653*c78c1075SJoseph Pusztay static PetscErrorCode SetProblem(TS ts)
654*c78c1075SJoseph Pusztay {
655*c78c1075SJoseph Pusztay   AppCtx *user;
656*c78c1075SJoseph Pusztay   DM      sw;
657*c78c1075SJoseph Pusztay 
658*c78c1075SJoseph Pusztay   PetscFunctionBegin;
659*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
660*c78c1075SJoseph Pusztay   PetscCall(DMGetApplicationContext(sw, (void **)&user));
661*c78c1075SJoseph Pusztay   // Define unified system for (X, V)
662*c78c1075SJoseph Pusztay   {
663*c78c1075SJoseph Pusztay     Mat      J;
664*c78c1075SJoseph Pusztay     PetscInt dim, Np;
665*c78c1075SJoseph Pusztay 
666*c78c1075SJoseph Pusztay     PetscCall(DMGetDimension(sw, &dim));
667*c78c1075SJoseph Pusztay     PetscCall(DMSwarmGetLocalSize(sw, &Np));
668*c78c1075SJoseph Pusztay     PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
669*c78c1075SJoseph Pusztay     PetscCall(MatSetSizes(J, 2 * Np * dim, 2 * Np * dim, PETSC_DECIDE, PETSC_DECIDE));
670*c78c1075SJoseph Pusztay     PetscCall(MatSetBlockSize(J, 2 * dim));
671*c78c1075SJoseph Pusztay     PetscCall(MatSetFromOptions(J));
672*c78c1075SJoseph Pusztay     PetscCall(MatSetUp(J));
673*c78c1075SJoseph Pusztay     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, user));
674*c78c1075SJoseph Pusztay     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, user));
675*c78c1075SJoseph Pusztay     PetscCall(MatDestroy(&J));
676*c78c1075SJoseph Pusztay   }
677*c78c1075SJoseph Pusztay   /* Define split system for X and V */
678*c78c1075SJoseph Pusztay   {
679*c78c1075SJoseph Pusztay     Vec             u;
680*c78c1075SJoseph Pusztay     IS              isx, isv, istmp;
681*c78c1075SJoseph Pusztay     const PetscInt *idx;
682*c78c1075SJoseph Pusztay     PetscInt        dim, Np, rstart;
683*c78c1075SJoseph Pusztay 
684*c78c1075SJoseph Pusztay     PetscCall(TSGetSolution(ts, &u));
685*c78c1075SJoseph Pusztay     PetscCall(DMGetDimension(sw, &dim));
686*c78c1075SJoseph Pusztay     PetscCall(DMSwarmGetLocalSize(sw, &Np));
687*c78c1075SJoseph Pusztay     PetscCall(VecGetOwnershipRange(u, &rstart, NULL));
688*c78c1075SJoseph Pusztay     PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 0, 2, &istmp));
689*c78c1075SJoseph Pusztay     PetscCall(ISGetIndices(istmp, &idx));
690*c78c1075SJoseph Pusztay     PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isx));
691*c78c1075SJoseph Pusztay     PetscCall(ISRestoreIndices(istmp, &idx));
692*c78c1075SJoseph Pusztay     PetscCall(ISDestroy(&istmp));
693*c78c1075SJoseph Pusztay     PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 1, 2, &istmp));
694*c78c1075SJoseph Pusztay     PetscCall(ISGetIndices(istmp, &idx));
695*c78c1075SJoseph Pusztay     PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isv));
696*c78c1075SJoseph Pusztay     PetscCall(ISRestoreIndices(istmp, &idx));
697*c78c1075SJoseph Pusztay     PetscCall(ISDestroy(&istmp));
698*c78c1075SJoseph Pusztay     PetscCall(TSRHSSplitSetIS(ts, "position", isx));
699*c78c1075SJoseph Pusztay     PetscCall(TSRHSSplitSetIS(ts, "momentum", isv));
700*c78c1075SJoseph Pusztay     PetscCall(ISDestroy(&isx));
701*c78c1075SJoseph Pusztay     PetscCall(ISDestroy(&isv));
702*c78c1075SJoseph Pusztay     PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, user));
703*c78c1075SJoseph Pusztay     PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, user));
704*c78c1075SJoseph Pusztay   }
705*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
706*c78c1075SJoseph Pusztay }
707*c78c1075SJoseph Pusztay 
708*c78c1075SJoseph Pusztay static PetscErrorCode DMSwarmTSRedistribute(TS ts)
709*c78c1075SJoseph Pusztay {
710*c78c1075SJoseph Pusztay   DM        sw;
711*c78c1075SJoseph Pusztay   Vec       u;
712*c78c1075SJoseph Pusztay   PetscReal t, maxt, dt;
713*c78c1075SJoseph Pusztay   PetscInt  n, maxn;
714*c78c1075SJoseph Pusztay 
715*c78c1075SJoseph Pusztay   PetscFunctionBegin;
716*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
717*c78c1075SJoseph Pusztay   PetscCall(TSGetTime(ts, &t));
718*c78c1075SJoseph Pusztay   PetscCall(TSGetMaxTime(ts, &maxt));
719*c78c1075SJoseph Pusztay   PetscCall(TSGetTimeStep(ts, &dt));
720*c78c1075SJoseph Pusztay   PetscCall(TSGetStepNumber(ts, &n));
721*c78c1075SJoseph Pusztay   PetscCall(TSGetMaxSteps(ts, &maxn));
722*c78c1075SJoseph Pusztay 
723*c78c1075SJoseph Pusztay   PetscCall(TSReset(ts));
724*c78c1075SJoseph Pusztay   PetscCall(TSSetDM(ts, sw));
725*c78c1075SJoseph Pusztay   /* TODO Check whether TS was set from options */
726*c78c1075SJoseph Pusztay   PetscCall(TSSetFromOptions(ts));
727*c78c1075SJoseph Pusztay   PetscCall(TSSetTime(ts, t));
728*c78c1075SJoseph Pusztay   PetscCall(TSSetMaxTime(ts, maxt));
729*c78c1075SJoseph Pusztay   PetscCall(TSSetTimeStep(ts, dt));
730*c78c1075SJoseph Pusztay   PetscCall(TSSetStepNumber(ts, n));
731*c78c1075SJoseph Pusztay   PetscCall(TSSetMaxSteps(ts, maxn));
732*c78c1075SJoseph Pusztay 
733*c78c1075SJoseph Pusztay   PetscCall(CreateSolution(ts));
734*c78c1075SJoseph Pusztay   PetscCall(SetProblem(ts));
735*c78c1075SJoseph Pusztay   PetscCall(TSGetSolution(ts, &u));
736*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
737*c78c1075SJoseph Pusztay }
738*c78c1075SJoseph Pusztay 
739*c78c1075SJoseph Pusztay PetscErrorCode circleSingleX(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx)
740*c78c1075SJoseph Pusztay {
741*c78c1075SJoseph Pusztay   x[0] = p + 1.;
742*c78c1075SJoseph Pusztay   x[1] = 0.;
743*c78c1075SJoseph Pusztay   return PETSC_SUCCESS;
744*c78c1075SJoseph Pusztay }
745*c78c1075SJoseph Pusztay 
746*c78c1075SJoseph Pusztay PetscErrorCode circleSingleV(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar v[], void *ctx)
747*c78c1075SJoseph Pusztay {
748*c78c1075SJoseph Pusztay   v[0] = 0.;
749*c78c1075SJoseph Pusztay   v[1] = PetscSqrtReal(1000. / (p + 1.));
750*c78c1075SJoseph Pusztay   return PETSC_SUCCESS;
751*c78c1075SJoseph Pusztay }
752*c78c1075SJoseph Pusztay 
753*c78c1075SJoseph Pusztay /* Put 5 particles into each circle */
754*c78c1075SJoseph Pusztay PetscErrorCode circleMultipleX(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar x[], void *ctx)
755*c78c1075SJoseph Pusztay {
756*c78c1075SJoseph Pusztay   const PetscInt  n   = 5;
757*c78c1075SJoseph Pusztay   const PetscReal r0  = (p / n) + 1.;
758*c78c1075SJoseph Pusztay   const PetscReal th0 = (2. * PETSC_PI * (p % n)) / n;
759*c78c1075SJoseph Pusztay 
760*c78c1075SJoseph Pusztay   x[0] = r0 * PetscCosReal(th0);
761*c78c1075SJoseph Pusztay   x[1] = r0 * PetscSinReal(th0);
762*c78c1075SJoseph Pusztay   return PETSC_SUCCESS;
763*c78c1075SJoseph Pusztay }
764*c78c1075SJoseph Pusztay 
765*c78c1075SJoseph Pusztay /* Put 5 particles into each circle */
766*c78c1075SJoseph Pusztay PetscErrorCode circleMultipleV(PetscInt dim, PetscReal time, const PetscReal dummy[], PetscInt p, PetscScalar v[], void *ctx)
767*c78c1075SJoseph Pusztay {
768*c78c1075SJoseph Pusztay   const PetscInt  n     = 5;
769*c78c1075SJoseph Pusztay   const PetscReal r0    = (p / n) + 1.;
770*c78c1075SJoseph Pusztay   const PetscReal th0   = (2. * PETSC_PI * (p % n)) / n;
771*c78c1075SJoseph Pusztay   const PetscReal omega = PetscSqrtReal(1000. / r0) / r0;
772*c78c1075SJoseph Pusztay 
773*c78c1075SJoseph Pusztay   v[0] = -r0 * omega * PetscSinReal(th0);
774*c78c1075SJoseph Pusztay   v[1] = r0 * omega * PetscCosReal(th0);
775*c78c1075SJoseph Pusztay   return PETSC_SUCCESS;
776*c78c1075SJoseph Pusztay }
777*c78c1075SJoseph Pusztay 
778*c78c1075SJoseph Pusztay /*
779*c78c1075SJoseph Pusztay   InitializeSolveAndSwarm - Set the solution values to the swarm coordinates and velocities, and also possibly set the initial values.
780*c78c1075SJoseph Pusztay 
781*c78c1075SJoseph Pusztay   Input Parameters:
782*c78c1075SJoseph Pusztay + ts         - The TS
783*c78c1075SJoseph Pusztay - useInitial - Flag to also set the initial conditions to the current coodinates and velocities and setup the problem
784*c78c1075SJoseph Pusztay 
785*c78c1075SJoseph Pusztay   Output Parameters:
786*c78c1075SJoseph Pusztay . u - The initialized solution vector
787*c78c1075SJoseph Pusztay 
788*c78c1075SJoseph Pusztay   Level: advanced
789*c78c1075SJoseph Pusztay 
790*c78c1075SJoseph Pusztay .seealso: InitializeSolve()
791*c78c1075SJoseph Pusztay */
792*c78c1075SJoseph Pusztay static PetscErrorCode InitializeSolveAndSwarm(TS ts, PetscBool useInitial)
793*c78c1075SJoseph Pusztay {
794*c78c1075SJoseph Pusztay   DM      sw;
795*c78c1075SJoseph Pusztay   Vec     u, gc, gv, gc0, gv0;
796*c78c1075SJoseph Pusztay   IS      isx, isv;
797*c78c1075SJoseph Pusztay   AppCtx *user;
798*c78c1075SJoseph Pusztay 
799*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
800*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
801*c78c1075SJoseph Pusztay   PetscCall(DMGetApplicationContext(sw, &user));
802*c78c1075SJoseph Pusztay   if (useInitial) {
803*c78c1075SJoseph Pusztay     PetscReal v0[1] = {1.};
804*c78c1075SJoseph Pusztay 
805*c78c1075SJoseph Pusztay     PetscCall(DMSwarmInitializeCoordinates(sw));
806*c78c1075SJoseph Pusztay     PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0));
807*c78c1075SJoseph Pusztay     PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
808*c78c1075SJoseph Pusztay     PetscCall(DMSwarmTSRedistribute(ts));
809*c78c1075SJoseph Pusztay   }
810*c78c1075SJoseph Pusztay   PetscCall(TSGetSolution(ts, &u));
811*c78c1075SJoseph Pusztay   PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
812*c78c1075SJoseph Pusztay   PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
813*c78c1075SJoseph Pusztay   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
814*c78c1075SJoseph Pusztay   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates", &gc0));
815*c78c1075SJoseph Pusztay   if (useInitial) PetscCall(VecCopy(gc, gc0));
816*c78c1075SJoseph Pusztay   PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc));
817*c78c1075SJoseph Pusztay   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
818*c78c1075SJoseph Pusztay   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates", &gc0));
819*c78c1075SJoseph Pusztay   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
820*c78c1075SJoseph Pusztay   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initVelocity", &gv0));
821*c78c1075SJoseph Pusztay   if (useInitial) PetscCall(VecCopy(gv, gv0));
822*c78c1075SJoseph Pusztay   PetscCall(VecISCopy(u, isv, SCATTER_FORWARD, gv));
823*c78c1075SJoseph Pusztay   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
824*c78c1075SJoseph Pusztay   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initVelocity", &gv0));
825*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
826*c78c1075SJoseph Pusztay }
827*c78c1075SJoseph Pusztay 
828*c78c1075SJoseph Pusztay static PetscErrorCode InitializeSolve(TS ts, Vec u)
829*c78c1075SJoseph Pusztay {
830*c78c1075SJoseph Pusztay   PetscFunctionBegin;
831*c78c1075SJoseph Pusztay   PetscCall(TSSetSolution(ts, u));
832*c78c1075SJoseph Pusztay   PetscCall(InitializeSolveAndSwarm(ts, PETSC_TRUE));
833*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
834*c78c1075SJoseph Pusztay }
835*c78c1075SJoseph Pusztay 
836*c78c1075SJoseph Pusztay static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
837*c78c1075SJoseph Pusztay {
838*c78c1075SJoseph Pusztay   MPI_Comm           comm;
839*c78c1075SJoseph Pusztay   DM                 sw;
840*c78c1075SJoseph Pusztay   AppCtx            *user;
841*c78c1075SJoseph Pusztay   const PetscScalar *u;
842*c78c1075SJoseph Pusztay   const PetscReal   *coords, *vel;
843*c78c1075SJoseph Pusztay   PetscScalar       *e;
844*c78c1075SJoseph Pusztay   PetscReal          t;
845*c78c1075SJoseph Pusztay   PetscInt           dim, Np, p;
846*c78c1075SJoseph Pusztay 
847*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
848*c78c1075SJoseph Pusztay   PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
849*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
850*c78c1075SJoseph Pusztay   PetscCall(DMGetApplicationContext(sw, &user));
851*c78c1075SJoseph Pusztay   PetscCall(DMGetDimension(sw, &dim));
852*c78c1075SJoseph Pusztay   PetscCall(TSGetSolveTime(ts, &t));
853*c78c1075SJoseph Pusztay   PetscCall(VecGetArray(E, &e));
854*c78c1075SJoseph Pusztay   PetscCall(VecGetArrayRead(U, &u));
855*c78c1075SJoseph Pusztay   PetscCall(VecGetLocalSize(U, &Np));
856*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
857*c78c1075SJoseph Pusztay   PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
858*c78c1075SJoseph Pusztay   Np /= 2 * dim;
859*c78c1075SJoseph Pusztay   for (p = 0; p < Np; ++p) {
860*c78c1075SJoseph Pusztay     /* TODO generalize initial conditions and project into plane instead of assuming x-y */
861*c78c1075SJoseph Pusztay     const PetscReal    r0    = DMPlex_NormD_Internal(dim, &coords[p * dim]);
862*c78c1075SJoseph Pusztay     const PetscReal    th0   = PetscAtan2Real(coords[p * dim + 1], coords[p * dim + 0]);
863*c78c1075SJoseph Pusztay     const PetscReal    v0    = DMPlex_NormD_Internal(dim, &vel[p * dim]);
864*c78c1075SJoseph Pusztay     const PetscReal    omega = v0 / r0;
865*c78c1075SJoseph Pusztay     const PetscReal    ct    = PetscCosReal(omega * t + th0);
866*c78c1075SJoseph Pusztay     const PetscReal    st    = PetscSinReal(omega * t + th0);
867*c78c1075SJoseph Pusztay     const PetscScalar *x     = &u[(p * 2 + 0) * dim];
868*c78c1075SJoseph Pusztay     const PetscScalar *v     = &u[(p * 2 + 1) * dim];
869*c78c1075SJoseph Pusztay     const PetscReal    xe[3] = {r0 * ct, r0 * st, 0.0};
870*c78c1075SJoseph Pusztay     const PetscReal    ve[3] = {-v0 * st, v0 * ct, 0.0};
871*c78c1075SJoseph Pusztay     PetscInt           d;
872*c78c1075SJoseph Pusztay 
873*c78c1075SJoseph Pusztay     for (d = 0; d < dim; ++d) {
874*c78c1075SJoseph Pusztay       e[(p * 2 + 0) * dim + d] = x[d] - xe[d];
875*c78c1075SJoseph Pusztay       e[(p * 2 + 1) * dim + d] = v[d] - ve[d];
876*c78c1075SJoseph Pusztay     }
877*c78c1075SJoseph Pusztay     if (user->error) {
878*c78c1075SJoseph Pusztay       const PetscReal en   = 0.5 * DMPlex_DotRealD_Internal(dim, v, v);
879*c78c1075SJoseph Pusztay       const PetscReal exen = 0.5 * PetscSqr(v0);
880*c78c1075SJoseph Pusztay       PetscCall(PetscPrintf(comm, "t %.4g: p%" PetscInt_FMT " error [%.2g %.2g] sol [(%.6lf %.6lf) (%.6lf %.6lf)] exact [(%.6lf %.6lf) (%.6lf %.6lf)] energy/exact energy %g / %g (%.10lf%%)\n", (double)t, p, (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 0) * dim]), (double)DMPlex_NormD_Internal(dim, &e[(p * 2 + 1) * dim]), (double)x[0], (double)x[1], (double)v[0], (double)v[1], (double)xe[0], (double)xe[1], (double)ve[0], (double)ve[1], (double)en, (double)exen, (double)(PetscAbsReal(exen - en) * 100. / exen)));
881*c78c1075SJoseph Pusztay     }
882*c78c1075SJoseph Pusztay   }
883*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
884*c78c1075SJoseph Pusztay   PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
885*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArrayRead(U, &u));
886*c78c1075SJoseph Pusztay   PetscCall(VecRestoreArray(E, &e));
887*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
888*c78c1075SJoseph Pusztay }
889*c78c1075SJoseph Pusztay 
890*c78c1075SJoseph Pusztay static PetscErrorCode EnergyMonitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx)
891*c78c1075SJoseph Pusztay {
892*c78c1075SJoseph Pusztay   const PetscInt     ostep = ((AppCtx *)ctx)->ostep;
893*c78c1075SJoseph Pusztay   const EMType       em    = ((AppCtx *)ctx)->em;
894*c78c1075SJoseph Pusztay   DM                 sw;
895*c78c1075SJoseph Pusztay   const PetscScalar *u;
896*c78c1075SJoseph Pusztay   PetscReal         *coords, *E;
897*c78c1075SJoseph Pusztay   PetscReal          enKin = 0., enEM = 0.;
898*c78c1075SJoseph Pusztay   PetscInt           dim, d, Np, p, q;
899*c78c1075SJoseph Pusztay 
900*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
901*c78c1075SJoseph Pusztay   if (step % ostep == 0) {
902*c78c1075SJoseph Pusztay     PetscCall(TSGetDM(ts, &sw));
903*c78c1075SJoseph Pusztay     PetscCall(DMGetDimension(sw, &dim));
904*c78c1075SJoseph Pusztay     PetscCall(VecGetArrayRead(U, &u));
905*c78c1075SJoseph Pusztay     PetscCall(VecGetLocalSize(U, &Np));
906*c78c1075SJoseph Pusztay     Np /= 2 * dim;
907*c78c1075SJoseph Pusztay     PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
908*c78c1075SJoseph Pusztay     PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
909*c78c1075SJoseph Pusztay     if (!step) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Time     Step Part     Energy\n"));
910*c78c1075SJoseph Pusztay     for (p = 0; p < Np; ++p) {
911*c78c1075SJoseph Pusztay       const PetscReal v2     = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &u[(p * 2 + 1) * dim]);
912*c78c1075SJoseph Pusztay       PetscReal      *pcoord = &coords[p * dim];
913*c78c1075SJoseph Pusztay 
914*c78c1075SJoseph Pusztay       PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " %5" PetscInt_FMT " %10.4lf\n", (double)t, step, p, (double)(0.5 * v2)));
915*c78c1075SJoseph Pusztay       enKin += 0.5 * v2;
916*c78c1075SJoseph Pusztay       if (em == EM_NONE) {
917*c78c1075SJoseph Pusztay         continue;
918*c78c1075SJoseph Pusztay       } else if (em == EM_COULOMB) {
919*c78c1075SJoseph Pusztay         for (q = p + 1; q < Np; ++q) {
920*c78c1075SJoseph Pusztay           PetscReal *qcoord = &coords[q * dim];
921*c78c1075SJoseph Pusztay           PetscReal  rpq[3], r;
922*c78c1075SJoseph Pusztay           for (d = 0; d < dim; ++d) rpq[d] = pcoord[d] - qcoord[d];
923*c78c1075SJoseph Pusztay           r = DMPlex_NormD_Internal(dim, rpq);
924*c78c1075SJoseph Pusztay           enEM += 1. / r;
925*c78c1075SJoseph Pusztay         }
926*c78c1075SJoseph Pusztay       } else if (em == EM_PRIMAL || em == EM_MIXED) {
927*c78c1075SJoseph Pusztay         for (d = 0; d < dim; ++d) enEM += E[p * dim + d];
928*c78c1075SJoseph Pusztay       }
929*c78c1075SJoseph Pusztay     }
930*c78c1075SJoseph Pusztay     PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " KE\t    %10.4lf\n", (double)t, step, (double)enKin));
931*c78c1075SJoseph Pusztay     PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " PE\t    %1.10g\n", (double)t, step, (double)enEM));
932*c78c1075SJoseph Pusztay     PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " E\t    %10.4lf\n", (double)t, step, (double)(enKin + enEM)));
933*c78c1075SJoseph Pusztay     PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
934*c78c1075SJoseph Pusztay     PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
935*c78c1075SJoseph Pusztay     PetscCall(PetscSynchronizedFlush(PetscObjectComm((PetscObject)ts), NULL));
936*c78c1075SJoseph Pusztay     PetscCall(VecRestoreArrayRead(U, &u));
937*c78c1075SJoseph Pusztay   }
938*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
939*c78c1075SJoseph Pusztay }
940*c78c1075SJoseph Pusztay 
941*c78c1075SJoseph Pusztay static PetscErrorCode MigrateParticles(TS ts)
942*c78c1075SJoseph Pusztay {
943*c78c1075SJoseph Pusztay   DM sw;
944*c78c1075SJoseph Pusztay 
945*c78c1075SJoseph Pusztay   PetscFunctionBeginUser;
946*c78c1075SJoseph Pusztay   PetscCall(TSGetDM(ts, &sw));
947*c78c1075SJoseph Pusztay   PetscCall(DMViewFromOptions(sw, NULL, "-migrate_view_pre"));
948*c78c1075SJoseph Pusztay   {
949*c78c1075SJoseph Pusztay     Vec u, gc, gv;
950*c78c1075SJoseph Pusztay     IS  isx, isv;
951*c78c1075SJoseph Pusztay 
952*c78c1075SJoseph Pusztay     PetscCall(TSGetSolution(ts, &u));
953*c78c1075SJoseph Pusztay     PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
954*c78c1075SJoseph Pusztay     PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
955*c78c1075SJoseph Pusztay     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
956*c78c1075SJoseph Pusztay     PetscCall(VecISCopy(u, isx, SCATTER_REVERSE, gc));
957*c78c1075SJoseph Pusztay     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
958*c78c1075SJoseph Pusztay     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
959*c78c1075SJoseph Pusztay     PetscCall(VecISCopy(u, isv, SCATTER_REVERSE, gv));
960*c78c1075SJoseph Pusztay     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
961*c78c1075SJoseph Pusztay   }
962*c78c1075SJoseph Pusztay   PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
963*c78c1075SJoseph Pusztay   PetscCall(DMSwarmTSRedistribute(ts));
964*c78c1075SJoseph Pusztay   PetscCall(InitializeSolveAndSwarm(ts, PETSC_FALSE));
965*c78c1075SJoseph Pusztay   PetscFunctionReturn(PETSC_SUCCESS);
966*c78c1075SJoseph Pusztay }
967*c78c1075SJoseph Pusztay 
968*c78c1075SJoseph Pusztay int main(int argc, char **argv)
969*c78c1075SJoseph Pusztay {
970*c78c1075SJoseph Pusztay   DM     dm, sw;
971*c78c1075SJoseph Pusztay   TS     ts;
972*c78c1075SJoseph Pusztay   Vec    u;
973*c78c1075SJoseph Pusztay   AppCtx user;
974*c78c1075SJoseph Pusztay 
975*c78c1075SJoseph Pusztay   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
976*c78c1075SJoseph Pusztay   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
977*c78c1075SJoseph Pusztay   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
978*c78c1075SJoseph Pusztay   PetscCall(CreatePoisson(dm, &user));
979*c78c1075SJoseph Pusztay   PetscCall(CreateSwarm(dm, &user, &sw));
980*c78c1075SJoseph Pusztay   PetscCall(DMSetApplicationContext(sw, &user));
981*c78c1075SJoseph Pusztay 
982*c78c1075SJoseph Pusztay   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
983*c78c1075SJoseph Pusztay   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
984*c78c1075SJoseph Pusztay   PetscCall(TSSetDM(ts, sw));
985*c78c1075SJoseph Pusztay   PetscCall(TSSetMaxTime(ts, 0.1));
986*c78c1075SJoseph Pusztay   PetscCall(TSSetTimeStep(ts, 0.00001));
987*c78c1075SJoseph Pusztay   PetscCall(TSSetMaxSteps(ts, 100));
988*c78c1075SJoseph Pusztay   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
989*c78c1075SJoseph Pusztay   PetscCall(TSMonitorSet(ts, EnergyMonitor, &user, NULL));
990*c78c1075SJoseph Pusztay   PetscCall(TSSetFromOptions(ts));
991*c78c1075SJoseph Pusztay   PetscCall(DMSwarmVectorDefineField(sw, "velocity"));
992*c78c1075SJoseph Pusztay   PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
993*c78c1075SJoseph Pusztay   PetscCall(TSSetComputeExactError(ts, ComputeError));
994*c78c1075SJoseph Pusztay   PetscCall(TSSetPostStep(ts, MigrateParticles));
995*c78c1075SJoseph Pusztay 
996*c78c1075SJoseph Pusztay   PetscCall(CreateSolution(ts));
997*c78c1075SJoseph Pusztay   PetscCall(TSGetSolution(ts, &u));
998*c78c1075SJoseph Pusztay   PetscCall(TSComputeInitialCondition(ts, u));
999*c78c1075SJoseph Pusztay   PetscCall(TSSolve(ts, NULL));
1000*c78c1075SJoseph Pusztay 
1001*c78c1075SJoseph Pusztay   PetscCall(SNESDestroy(&user.snes));
1002*c78c1075SJoseph Pusztay   PetscCall(TSDestroy(&ts));
1003*c78c1075SJoseph Pusztay   PetscCall(DMDestroy(&sw));
1004*c78c1075SJoseph Pusztay   PetscCall(DMDestroy(&dm));
1005*c78c1075SJoseph Pusztay   PetscCall(PetscFinalize());
1006*c78c1075SJoseph Pusztay   return 0;
1007*c78c1075SJoseph Pusztay }
1008*c78c1075SJoseph Pusztay 
1009*c78c1075SJoseph Pusztay /*TEST
1010*c78c1075SJoseph Pusztay 
1011*c78c1075SJoseph Pusztay    build:
1012*c78c1075SJoseph Pusztay      requires: double !complex
1013*c78c1075SJoseph Pusztay 
1014*c78c1075SJoseph Pusztay    testset:
1015*c78c1075SJoseph Pusztay      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1016*c78c1075SJoseph Pusztay      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
1017*c78c1075SJoseph Pusztay            -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1018*c78c1075SJoseph Pusztay            -ts_type basicsymplectic\
1019*c78c1075SJoseph Pusztay            -dm_view -output_step 50 -ts_dt 0.01 -ts_max_time 10.0 -ts_max_steps 10
1020*c78c1075SJoseph Pusztay      test:
1021*c78c1075SJoseph Pusztay        suffix: none_bsi_2d_1
1022*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 1 -em_type none -error
1023*c78c1075SJoseph Pusztay      test:
1024*c78c1075SJoseph Pusztay        suffix: none_bsi_2d_2
1025*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 2 -em_type none -error
1026*c78c1075SJoseph Pusztay      test:
1027*c78c1075SJoseph Pusztay        suffix: none_bsi_2d_3
1028*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 3 -em_type none -error
1029*c78c1075SJoseph Pusztay      test:
1030*c78c1075SJoseph Pusztay        suffix: none_bsi_2d_4
1031*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 4 -em_type none -error
1032*c78c1075SJoseph Pusztay      test:
1033*c78c1075SJoseph Pusztay        suffix: coulomb_bsi_2d_1
1034*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 1
1035*c78c1075SJoseph Pusztay      test:
1036*c78c1075SJoseph Pusztay        suffix: coulomb_bsi_2d_2
1037*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 2
1038*c78c1075SJoseph Pusztay      test:
1039*c78c1075SJoseph Pusztay        suffix: coulomb_bsi_2d_3
1040*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 3
1041*c78c1075SJoseph Pusztay      test:
1042*c78c1075SJoseph Pusztay        suffix: coulomb_bsi_2d_4
1043*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 4
1044*c78c1075SJoseph Pusztay 
1045*c78c1075SJoseph Pusztay    testset:
1046*c78c1075SJoseph Pusztay      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1047*c78c1075SJoseph Pusztay      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
1048*c78c1075SJoseph Pusztay            -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1049*c78c1075SJoseph Pusztay            -ts_type basicsymplectic\
1050*c78c1075SJoseph Pusztay            -em_type primal -em_pc_type svd\
1051*c78c1075SJoseph Pusztay            -dm_view -output_step 50 -error -ts_dt 0.01 -ts_max_time 10.0 -ts_max_steps 10\
1052*c78c1075SJoseph Pusztay            -petscspace_degree 2 -petscfe_default_quadrature_order 3 -sigma 1.0e-8 -timeScale 2.0e-14
1053*c78c1075SJoseph Pusztay      test:
1054*c78c1075SJoseph Pusztay        suffix: poisson_bsi_2d_1
1055*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 1
1056*c78c1075SJoseph Pusztay      test:
1057*c78c1075SJoseph Pusztay        suffix: poisson_bsi_2d_2
1058*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 2
1059*c78c1075SJoseph Pusztay      test:
1060*c78c1075SJoseph Pusztay        suffix: poisson_bsi_2d_3
1061*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 3
1062*c78c1075SJoseph Pusztay      test:
1063*c78c1075SJoseph Pusztay        suffix: poisson_bsi_2d_4
1064*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 4
1065*c78c1075SJoseph Pusztay 
1066*c78c1075SJoseph Pusztay    testset:
1067*c78c1075SJoseph Pusztay      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1068*c78c1075SJoseph Pusztay      args: -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1069*c78c1075SJoseph Pusztay            -ts_type theta -ts_theta_theta 0.5 -ts_convergence_estimate -convest_num_refine 2 \
1070*c78c1075SJoseph Pusztay              -mat_type baij -ksp_error_if_not_converged -em_pc_type svd\
1071*c78c1075SJoseph Pusztay            -dm_view -output_step 50 -error -ts_dt 0.01 -ts_max_time 10.0 -ts_max_steps 10\
1072*c78c1075SJoseph Pusztay            -pc_type svd -sigma 1.0e-8 -timeScale 2.0e-14
1073*c78c1075SJoseph Pusztay      test:
1074*c78c1075SJoseph Pusztay        suffix: im_2d_0
1075*c78c1075SJoseph Pusztay        args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 1,1 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5
1076*c78c1075SJoseph Pusztay 
1077*c78c1075SJoseph Pusztay    testset:
1078*c78c1075SJoseph Pusztay      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1079*c78c1075SJoseph Pusztay      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 10,10 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 -petscpartitioner_type simple \
1080*c78c1075SJoseph Pusztay            -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV -dm_swarm_num_species 1\
1081*c78c1075SJoseph Pusztay            -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
1082*c78c1075SJoseph Pusztay            -em_snes_type ksponly -em_pc_type svd -em_type primal -petscspace_degree 1\
1083*c78c1075SJoseph Pusztay            -dm_view -output_step 50\
1084*c78c1075SJoseph Pusztay            -pc_type svd -sigma 1.0e-8 -timeScale 2.0e-14 -ts_dt 0.01 -ts_max_time 1.0 -ts_max_steps 10
1085*c78c1075SJoseph Pusztay      test:
1086*c78c1075SJoseph Pusztay        suffix: bsi_2d_mesh_1
1087*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 4
1088*c78c1075SJoseph Pusztay      test:
1089*c78c1075SJoseph Pusztay        suffix: bsi_2d_mesh_1_par_2
1090*c78c1075SJoseph Pusztay        nsize: 2
1091*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 4
1092*c78c1075SJoseph Pusztay      test:
1093*c78c1075SJoseph Pusztay        suffix: bsi_2d_mesh_1_par_3
1094*c78c1075SJoseph Pusztay        nsize: 3
1095*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 4
1096*c78c1075SJoseph Pusztay      test:
1097*c78c1075SJoseph Pusztay        suffix: bsi_2d_mesh_1_par_4
1098*c78c1075SJoseph Pusztay        nsize: 4
1099*c78c1075SJoseph Pusztay        args: -ts_basicsymplectic_type 4 -dm_swarm_num_particles 0,0,2,0
1100*c78c1075SJoseph Pusztay 
1101*c78c1075SJoseph Pusztay    testset:
1102*c78c1075SJoseph Pusztay      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1103*c78c1075SJoseph Pusztay      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 10,10 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
1104*c78c1075SJoseph Pusztay            -dm_swarm_num_particles 10 -dm_swarm_coordinate_function circleMultipleX -dm_swarm_velocity_function circleMultipleV \
1105*c78c1075SJoseph Pusztay            -ts_convergence_estimate -convest_num_refine 2 \
1106*c78c1075SJoseph Pusztay              -em_pc_type lu\
1107*c78c1075SJoseph Pusztay            -dm_view -output_step 50 -error\
1108*c78c1075SJoseph Pusztay            -pc_type svd -sigma 1.0e-8 -timeScale 2.0e-14 -ts_dt 0.01 -ts_max_time 10.0 -ts_max_steps 10
1109*c78c1075SJoseph Pusztay      test:
1110*c78c1075SJoseph Pusztay        suffix: bsi_2d_multiple_1
1111*c78c1075SJoseph Pusztay        args: -ts_type basicsymplectic -ts_basicsymplectic_type 1
1112*c78c1075SJoseph Pusztay      test:
1113*c78c1075SJoseph Pusztay        suffix: bsi_2d_multiple_2
1114*c78c1075SJoseph Pusztay        args: -ts_type basicsymplectic -ts_basicsymplectic_type 2
1115*c78c1075SJoseph Pusztay      test:
1116*c78c1075SJoseph Pusztay        suffix: bsi_2d_multiple_3
1117*c78c1075SJoseph Pusztay        args: -ts_type basicsymplectic -ts_basicsymplectic_type 3 -ts_dt 0.001
1118*c78c1075SJoseph Pusztay      test:
1119*c78c1075SJoseph Pusztay        suffix: im_2d_multiple_0
1120*c78c1075SJoseph Pusztay        args: -ts_type theta -ts_theta_theta 0.5 \
1121*c78c1075SJoseph Pusztay                -mat_type baij -ksp_error_if_not_converged -em_pc_type lu
1122*c78c1075SJoseph Pusztay 
1123*c78c1075SJoseph Pusztay    testset:
1124*c78c1075SJoseph Pusztay      requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1125*c78c1075SJoseph Pusztay      args: -dm_plex_dim 2 -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -dm_plex_box_lower -5,-5 -dm_plex_box_upper 5,5 \
1126*c78c1075SJoseph Pusztay            -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1127*c78c1075SJoseph Pusztay            -em_pc_type fieldsplit -ksp_rtol 1e-10 -em_ksp_type preonly -em_type mixed -em_ksp_error_if_not_converged\
1128*c78c1075SJoseph Pusztay            -dm_view -output_step 50 -error -dm_refine 0\
1129*c78c1075SJoseph Pusztay            -pc_type svd -sigma 1.0e-8 -timeScale 2.0e-14 -ts_dt 0.01 -ts_max_time 10.0 -ts_max_steps 10
1130*c78c1075SJoseph Pusztay      test:
1131*c78c1075SJoseph Pusztay        suffix: bsi_4_rt_1
1132*c78c1075SJoseph Pusztay        args: -ts_type basicsymplectic -ts_basicsymplectic_type 4\
1133*c78c1075SJoseph Pusztay              -pc_fieldsplit_detect_saddle_point\
1134*c78c1075SJoseph Pusztay              -pc_type fieldsplit\
1135*c78c1075SJoseph Pusztay              -pc_fieldsplit_type schur\
1136*c78c1075SJoseph Pusztay              -pc_fieldsplit_schur_precondition full \
1137*c78c1075SJoseph Pusztay              -field_petscspace_degree 2\
1138*c78c1075SJoseph Pusztay              -field_petscfe_default_quadrature_order 1\
1139*c78c1075SJoseph Pusztay              -field_petscspace_type sum \
1140*c78c1075SJoseph Pusztay              -field_petscspace_variables 2 \
1141*c78c1075SJoseph Pusztay              -field_petscspace_components 2 \
1142*c78c1075SJoseph Pusztay              -field_petscspace_sum_spaces 2 \
1143*c78c1075SJoseph Pusztay              -field_petscspace_sum_concatenate true \
1144*c78c1075SJoseph Pusztay              -field_sumcomp_0_petscspace_variables 2 \
1145*c78c1075SJoseph Pusztay              -field_sumcomp_0_petscspace_type tensor \
1146*c78c1075SJoseph Pusztay              -field_sumcomp_0_petscspace_tensor_spaces 2 \
1147*c78c1075SJoseph Pusztay              -field_sumcomp_0_petscspace_tensor_uniform false \
1148*c78c1075SJoseph Pusztay              -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
1149*c78c1075SJoseph Pusztay              -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
1150*c78c1075SJoseph Pusztay              -field_sumcomp_1_petscspace_variables 2 \
1151*c78c1075SJoseph Pusztay              -field_sumcomp_1_petscspace_type tensor \
1152*c78c1075SJoseph Pusztay              -field_sumcomp_1_petscspace_tensor_spaces 2 \
1153*c78c1075SJoseph Pusztay              -field_sumcomp_1_petscspace_tensor_uniform false \
1154*c78c1075SJoseph Pusztay              -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
1155*c78c1075SJoseph Pusztay              -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
1156*c78c1075SJoseph Pusztay              -field_petscdualspace_form_degree -1 \
1157*c78c1075SJoseph Pusztay              -field_petscdualspace_order 1 \
1158*c78c1075SJoseph Pusztay              -field_petscdualspace_lagrange_trimmed true\
1159*c78c1075SJoseph Pusztay              -ksp_gmres_restart 500
1160*c78c1075SJoseph Pusztay 
1161*c78c1075SJoseph Pusztay TEST*/
1162