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