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
10188af4bfSBarry Smith -convest_num_refine 0 -ts_time_step 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
ProcessOptions(MPI_Comm comm,AppCtx * options)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
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)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
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[])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
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[])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 */
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[])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
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[])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> */
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[])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
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[])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
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[])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
CreateFEM(DM dm,AppCtx * user)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
CreatePoisson(DM dm,AppCtx * user)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
CreateSwarm(DM dm,AppCtx * user,DM * sw)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 }
22484467f86SMatthew G. Knepley {
22584467f86SMatthew G. Knepley const char *fieldnames[2] = {DMSwarmPICField_coor, "velocity"};
22684467f86SMatthew G. Knepley
22784467f86SMatthew G. Knepley PetscCall(DMSwarmVectorDefineFields(*sw, 2, fieldnames));
22884467f86SMatthew G. Knepley }
229c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
230c78c1075SJoseph Pusztay }
231c78c1075SJoseph Pusztay
ComputeFieldAtParticles_Coulomb(SNES snes,DM sw,PetscReal E[])232c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles_Coulomb(SNES snes, DM sw, PetscReal E[])
233c78c1075SJoseph Pusztay {
234c78c1075SJoseph Pusztay PetscReal *coords;
235c78c1075SJoseph Pusztay PetscInt dim, d, Np, p, q;
236c78c1075SJoseph Pusztay PetscMPIInt size;
237c78c1075SJoseph Pusztay
238c78c1075SJoseph Pusztay PetscFunctionBegin;
239c78c1075SJoseph Pusztay PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)snes), &size));
240c78c1075SJoseph Pusztay PetscCheck(size == 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Coulomb code only works in serial");
241c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
242c78c1075SJoseph Pusztay PetscCall(DMSwarmGetLocalSize(sw, &Np));
243c78c1075SJoseph Pusztay
244c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
245c78c1075SJoseph Pusztay for (p = 0; p < Np; ++p) {
246c78c1075SJoseph Pusztay PetscReal *pcoord = &coords[p * dim];
247c78c1075SJoseph Pusztay PetscReal *pE = &E[p * dim];
248c78c1075SJoseph Pusztay /* Calculate field at particle p due to particle q */
249c78c1075SJoseph Pusztay for (q = 0; q < Np; ++q) {
250c78c1075SJoseph Pusztay PetscReal *qcoord = &coords[q * dim];
251c78c1075SJoseph Pusztay PetscReal rpq[3], r;
252c78c1075SJoseph Pusztay
253c78c1075SJoseph Pusztay if (p == q) continue;
254c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) rpq[d] = pcoord[d] - qcoord[d];
255c78c1075SJoseph Pusztay r = DMPlex_NormD_Internal(dim, rpq);
256c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) pE[d] += rpq[d] / PetscPowRealInt(r, 3);
257c78c1075SJoseph Pusztay }
258c78c1075SJoseph Pusztay }
259c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
260c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
261c78c1075SJoseph Pusztay }
262c78c1075SJoseph Pusztay
ComputeFieldAtParticles_Primal(SNES snes,DM sw,PetscReal E[])263c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles_Primal(SNES snes, DM sw, PetscReal E[])
264c78c1075SJoseph Pusztay {
265c78c1075SJoseph Pusztay DM dm;
266c78c1075SJoseph Pusztay PetscDS ds;
267c78c1075SJoseph Pusztay PetscFE fe;
268c78c1075SJoseph Pusztay Mat M_p;
269c78c1075SJoseph Pusztay Vec phi, locPhi, rho, f;
270c78c1075SJoseph Pusztay PetscReal *coords, chargeTol = 1e-13;
271c78c1075SJoseph Pusztay PetscInt dim, d, cStart, cEnd, c, Np;
27284467f86SMatthew G. Knepley const char **oldFields;
27384467f86SMatthew G. Knepley PetscInt Nf;
27484467f86SMatthew G. Knepley const char **tmp;
275c78c1075SJoseph Pusztay
276c78c1075SJoseph Pusztay PetscFunctionBegin;
277c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
278c78c1075SJoseph Pusztay PetscCall(DMSwarmGetLocalSize(sw, &Np));
2790bf7c1c5SMatthew G. Knepley PetscCall(SNESGetDM(snes, &dm));
2800bf7c1c5SMatthew G. Knepley
28184467f86SMatthew G. Knepley PetscCall(DMSwarmVectorGetField(sw, &Nf, &tmp));
28284467f86SMatthew G. Knepley PetscCall(PetscMalloc1(Nf, &oldFields));
28384467f86SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscStrallocpy(tmp[f], (char **)&oldFields[f]));
2840bf7c1c5SMatthew G. Knepley PetscCall(DMSwarmVectorDefineField(sw, "w_q"));
2850bf7c1c5SMatthew G. Knepley PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
28684467f86SMatthew G. Knepley PetscCall(DMSwarmVectorDefineFields(sw, Nf, oldFields));
28784467f86SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(oldFields[f]));
28884467f86SMatthew G. Knepley PetscCall(PetscFree(oldFields));
289c78c1075SJoseph Pusztay
290c78c1075SJoseph Pusztay /* Create the charges rho */
291c78c1075SJoseph Pusztay PetscCall(DMGetGlobalVector(dm, &rho));
292c78c1075SJoseph Pusztay PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
293c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
294c78c1075SJoseph Pusztay PetscCall(PetscObjectSetName((PetscObject)f, "particle weight"));
295c78c1075SJoseph Pusztay PetscCall(MatMultTranspose(M_p, f, rho));
296c78c1075SJoseph Pusztay PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
297c78c1075SJoseph Pusztay PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
298c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
299c78c1075SJoseph Pusztay PetscCall(MatDestroy(&M_p));
300c78c1075SJoseph Pusztay {
301c78c1075SJoseph Pusztay PetscScalar sum;
302c78c1075SJoseph Pusztay PetscInt n;
303c78c1075SJoseph Pusztay PetscReal phi_0 = 1.; /* (sigma*sigma*sigma)*(timeScale*timeScale)/(m_e*q_e*epsi_0)*/
304c78c1075SJoseph Pusztay
305c78c1075SJoseph Pusztay /* Remove constant from rho */
306c78c1075SJoseph Pusztay PetscCall(VecGetSize(rho, &n));
307c78c1075SJoseph Pusztay PetscCall(VecSum(rho, &sum));
308c78c1075SJoseph Pusztay PetscCall(VecShift(rho, -sum / n));
309c78c1075SJoseph Pusztay PetscCall(VecSum(rho, &sum));
310c78c1075SJoseph Pusztay PetscCheck(PetscAbsScalar(sum) < chargeTol, PetscObjectComm((PetscObject)sw), PETSC_ERR_PLIB, "Charge should have no DC component: %g", (double)PetscRealPart(sum));
311c78c1075SJoseph Pusztay /* Nondimensionalize rho */
312c78c1075SJoseph Pusztay PetscCall(VecScale(rho, phi_0));
313c78c1075SJoseph Pusztay }
314c78c1075SJoseph Pusztay PetscCall(VecViewFromOptions(rho, NULL, "-poisson_rho_view"));
315c78c1075SJoseph Pusztay
316c78c1075SJoseph Pusztay PetscCall(DMGetGlobalVector(dm, &phi));
317c78c1075SJoseph Pusztay PetscCall(PetscObjectSetName((PetscObject)phi, "potential"));
318c78c1075SJoseph Pusztay PetscCall(VecSet(phi, 0.0));
319c78c1075SJoseph Pusztay PetscCall(SNESSolve(snes, rho, phi));
320c78c1075SJoseph Pusztay PetscCall(DMRestoreGlobalVector(dm, &rho));
321c78c1075SJoseph Pusztay PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
322c78c1075SJoseph Pusztay
323c78c1075SJoseph Pusztay PetscCall(DMGetLocalVector(dm, &locPhi));
324c78c1075SJoseph Pusztay PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi));
325c78c1075SJoseph Pusztay PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi));
326c78c1075SJoseph Pusztay PetscCall(DMRestoreGlobalVector(dm, &phi));
327c78c1075SJoseph Pusztay
328c78c1075SJoseph Pusztay PetscCall(DMGetDS(dm, &ds));
329c78c1075SJoseph Pusztay PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
330c78c1075SJoseph Pusztay PetscCall(DMSwarmSortGetAccess(sw));
331c78c1075SJoseph Pusztay PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
332c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
333c78c1075SJoseph Pusztay for (c = cStart; c < cEnd; ++c) {
334c78c1075SJoseph Pusztay PetscTabulation tab;
335c78c1075SJoseph Pusztay PetscScalar *clPhi = NULL;
336c78c1075SJoseph Pusztay PetscReal *pcoord, *refcoord;
337c78c1075SJoseph Pusztay PetscReal v[3], J[9], invJ[9], detJ;
338c78c1075SJoseph Pusztay PetscInt *points;
339c78c1075SJoseph Pusztay PetscInt Ncp, cp;
340c78c1075SJoseph Pusztay
341c78c1075SJoseph Pusztay PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points));
342c78c1075SJoseph Pusztay PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
343c78c1075SJoseph Pusztay PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
344c78c1075SJoseph Pusztay for (cp = 0; cp < Ncp; ++cp)
345c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d];
346c78c1075SJoseph Pusztay PetscCall(DMPlexCoordinatesToReference(dm, c, Ncp, pcoord, refcoord));
347c78c1075SJoseph Pusztay PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
348c78c1075SJoseph Pusztay PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v, J, invJ, &detJ));
349c78c1075SJoseph Pusztay PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi));
350c78c1075SJoseph Pusztay for (cp = 0; cp < Ncp; ++cp) {
351c78c1075SJoseph Pusztay const PetscReal *basisDer = tab->T[1];
352c78c1075SJoseph Pusztay const PetscInt p = points[cp];
353c78c1075SJoseph Pusztay
354c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) E[p * dim + d] = 0.;
355c78c1075SJoseph Pusztay PetscCall(PetscFEFreeInterpolateGradient_Static(fe, basisDer, clPhi, dim, invJ, NULL, cp, &E[p * dim]));
356c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) E[p * dim + d] *= -1.0;
357c78c1075SJoseph Pusztay }
358c78c1075SJoseph Pusztay PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi));
359c78c1075SJoseph Pusztay PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
360c78c1075SJoseph Pusztay PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
361c78c1075SJoseph Pusztay PetscCall(PetscTabulationDestroy(&tab));
36284467f86SMatthew G. Knepley PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points));
363c78c1075SJoseph Pusztay }
364c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
365c78c1075SJoseph Pusztay PetscCall(DMSwarmSortRestoreAccess(sw));
366c78c1075SJoseph Pusztay PetscCall(DMRestoreLocalVector(dm, &locPhi));
367c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
368c78c1075SJoseph Pusztay }
369c78c1075SJoseph Pusztay
ComputeFieldAtParticles_Mixed(SNES snes,DM sw,PetscReal E[])370c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles_Mixed(SNES snes, DM sw, PetscReal E[])
371c78c1075SJoseph Pusztay {
372c78c1075SJoseph Pusztay DM dm, potential_dm;
373c78c1075SJoseph Pusztay IS potential_IS;
374c78c1075SJoseph Pusztay PetscDS ds;
375c78c1075SJoseph Pusztay PetscFE fe;
376c78c1075SJoseph Pusztay PetscFEGeom feGeometry;
377c78c1075SJoseph Pusztay Mat M_p;
378c78c1075SJoseph Pusztay Vec phi, locPhi, rho, f, temp_rho;
379c78c1075SJoseph Pusztay PetscQuadrature q;
380c78c1075SJoseph Pusztay PetscReal *coords, chargeTol = 1e-13;
38184467f86SMatthew G. Knepley PetscInt dim, d, cStart, cEnd, c, Np, pot_field = 1;
38284467f86SMatthew G. Knepley const char **oldFields;
38384467f86SMatthew G. Knepley PetscInt Nf;
38484467f86SMatthew G. Knepley const char **tmp;
385c78c1075SJoseph Pusztay
386c78c1075SJoseph Pusztay PetscFunctionBegin;
387c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
388c78c1075SJoseph Pusztay PetscCall(DMSwarmGetLocalSize(sw, &Np));
389c78c1075SJoseph Pusztay
390c78c1075SJoseph Pusztay /* Create the charges rho */
391c78c1075SJoseph Pusztay PetscCall(SNESGetDM(snes, &dm));
392c78c1075SJoseph Pusztay PetscCall(DMGetGlobalVector(dm, &rho));
393c78c1075SJoseph Pusztay PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
39484467f86SMatthew G. Knepley PetscCall(DMCreateSubDM(dm, 1, &pot_field, &potential_IS, &potential_dm));
3950bf7c1c5SMatthew G. Knepley
39684467f86SMatthew G. Knepley PetscCall(DMSwarmVectorGetField(sw, &Nf, &tmp));
39784467f86SMatthew G. Knepley PetscCall(PetscMalloc1(Nf, &oldFields));
39884467f86SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscStrallocpy(tmp[f], (char **)&oldFields[f]));
3990bf7c1c5SMatthew G. Knepley PetscCall(DMSwarmVectorDefineField(sw, "w_q"));
400c78c1075SJoseph Pusztay PetscCall(DMCreateMassMatrix(sw, potential_dm, &M_p));
40184467f86SMatthew G. Knepley PetscCall(DMSwarmVectorDefineFields(sw, Nf, oldFields));
40284467f86SMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscFree(oldFields[f]));
40384467f86SMatthew G. Knepley PetscCall(PetscFree(oldFields));
4040bf7c1c5SMatthew G. Knepley
405c78c1075SJoseph Pusztay PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
406c78c1075SJoseph Pusztay PetscCall(DMGetGlobalVector(potential_dm, &temp_rho));
407c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
408c78c1075SJoseph Pusztay PetscCall(PetscObjectSetName((PetscObject)f, "particle weight"));
409c78c1075SJoseph Pusztay PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
410c78c1075SJoseph Pusztay PetscCall(MatMultTranspose(M_p, f, temp_rho));
411c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
412c78c1075SJoseph Pusztay PetscCall(MatDestroy(&M_p));
413c78c1075SJoseph Pusztay PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
414c78c1075SJoseph Pusztay PetscCall(VecViewFromOptions(rho, NULL, "-poisson_rho_view"));
415c78c1075SJoseph Pusztay PetscCall(VecISCopy(rho, potential_IS, SCATTER_FORWARD, temp_rho));
416c78c1075SJoseph Pusztay PetscCall(DMRestoreGlobalVector(potential_dm, &temp_rho));
417c78c1075SJoseph Pusztay PetscCall(DMDestroy(&potential_dm));
418c78c1075SJoseph Pusztay PetscCall(ISDestroy(&potential_IS));
419c78c1075SJoseph Pusztay {
420c78c1075SJoseph Pusztay PetscScalar sum;
421c78c1075SJoseph Pusztay PetscInt n;
422c78c1075SJoseph Pusztay PetscReal phi_0 = 1.; /*(sigma*sigma*sigma)*(timeScale*timeScale)/(m_e*q_e*epsi_0);*/
423c78c1075SJoseph Pusztay
424c78c1075SJoseph Pusztay /*Remove constant from rho*/
425c78c1075SJoseph Pusztay PetscCall(VecGetSize(rho, &n));
426c78c1075SJoseph Pusztay PetscCall(VecSum(rho, &sum));
427c78c1075SJoseph Pusztay PetscCall(VecShift(rho, -sum / n));
428c78c1075SJoseph Pusztay PetscCall(VecSum(rho, &sum));
429c78c1075SJoseph Pusztay PetscCheck(PetscAbsScalar(sum) < chargeTol, PetscObjectComm((PetscObject)sw), PETSC_ERR_PLIB, "Charge should have no DC component: %g", (double)PetscRealPart(sum));
430c78c1075SJoseph Pusztay /* Nondimensionalize rho */
431c78c1075SJoseph Pusztay PetscCall(VecScale(rho, phi_0));
432c78c1075SJoseph Pusztay }
433c78c1075SJoseph Pusztay PetscCall(DMGetGlobalVector(dm, &phi));
434c78c1075SJoseph Pusztay PetscCall(PetscObjectSetName((PetscObject)phi, "potential"));
435c78c1075SJoseph Pusztay PetscCall(VecSet(phi, 0.0));
436c78c1075SJoseph Pusztay PetscCall(SNESSolve(snes, NULL, phi));
437c78c1075SJoseph Pusztay PetscCall(DMRestoreGlobalVector(dm, &rho));
438c78c1075SJoseph Pusztay PetscCall(VecViewFromOptions(phi, NULL, "-phi_view"));
439c78c1075SJoseph Pusztay
440c78c1075SJoseph Pusztay PetscCall(DMGetLocalVector(dm, &locPhi));
441c78c1075SJoseph Pusztay PetscCall(DMGlobalToLocalBegin(dm, phi, INSERT_VALUES, locPhi));
442c78c1075SJoseph Pusztay PetscCall(DMGlobalToLocalEnd(dm, phi, INSERT_VALUES, locPhi));
443c78c1075SJoseph Pusztay PetscCall(DMRestoreGlobalVector(dm, &phi));
444c78c1075SJoseph Pusztay
445c78c1075SJoseph Pusztay PetscCall(DMGetDS(dm, &ds));
446c78c1075SJoseph Pusztay PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
447c78c1075SJoseph Pusztay PetscCall(DMSwarmSortGetAccess(sw));
448c78c1075SJoseph Pusztay PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
449c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
450c78c1075SJoseph Pusztay for (c = cStart; c < cEnd; ++c) {
451c78c1075SJoseph Pusztay PetscTabulation tab;
452c78c1075SJoseph Pusztay PetscScalar *clPhi = NULL;
453c78c1075SJoseph Pusztay PetscReal *pcoord, *refcoord;
454c78c1075SJoseph Pusztay PetscReal v[3], J[9], invJ[9], detJ;
455c78c1075SJoseph Pusztay PetscInt *points;
456c78c1075SJoseph Pusztay PetscInt Ncp, cp;
457c78c1075SJoseph Pusztay
458c78c1075SJoseph Pusztay PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points));
459c78c1075SJoseph Pusztay PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
460c78c1075SJoseph Pusztay PetscCall(DMGetWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
461c78c1075SJoseph Pusztay for (cp = 0; cp < Ncp; ++cp)
462c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) pcoord[cp * dim + d] = coords[points[cp] * dim + d];
463c78c1075SJoseph Pusztay PetscCall(DMPlexCoordinatesToReference(dm, c, Ncp, pcoord, refcoord));
464c78c1075SJoseph Pusztay PetscCall(PetscFECreateTabulation(fe, 1, Ncp, refcoord, 1, &tab));
465c78c1075SJoseph Pusztay PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v, J, invJ, &detJ));
466c78c1075SJoseph Pusztay PetscCall(DMPlexVecGetClosure(dm, NULL, locPhi, c, NULL, &clPhi));
467c78c1075SJoseph Pusztay for (cp = 0; cp < Ncp; ++cp) {
468c78c1075SJoseph Pusztay const PetscInt p = points[cp];
469c78c1075SJoseph Pusztay
470c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) E[p * dim + d] = 0.;
471c78c1075SJoseph Pusztay PetscCall(PetscFEGetQuadrature(fe, &q));
472c78c1075SJoseph Pusztay PetscCall(PetscFECreateCellGeometry(fe, q, &feGeometry));
473c78c1075SJoseph Pusztay PetscCall(PetscFEInterpolateAtPoints_Static(fe, tab, clPhi, &feGeometry, cp, &E[p * dim]));
474c78c1075SJoseph Pusztay PetscCall(PetscFEDestroyCellGeometry(fe, &feGeometry));
475c78c1075SJoseph Pusztay }
476c78c1075SJoseph Pusztay PetscCall(DMPlexVecRestoreClosure(dm, NULL, locPhi, c, NULL, &clPhi));
477c78c1075SJoseph Pusztay PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &pcoord));
478c78c1075SJoseph Pusztay PetscCall(DMRestoreWorkArray(dm, Ncp * dim, MPIU_REAL, &refcoord));
479c78c1075SJoseph Pusztay PetscCall(PetscTabulationDestroy(&tab));
48084467f86SMatthew G. Knepley PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points));
481c78c1075SJoseph Pusztay }
482c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
483c78c1075SJoseph Pusztay PetscCall(DMSwarmSortRestoreAccess(sw));
484c78c1075SJoseph Pusztay PetscCall(DMRestoreLocalVector(dm, &locPhi));
485c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
486c78c1075SJoseph Pusztay }
487c78c1075SJoseph Pusztay
ComputeFieldAtParticles(SNES snes,DM sw,PetscReal E[])488c78c1075SJoseph Pusztay static PetscErrorCode ComputeFieldAtParticles(SNES snes, DM sw, PetscReal E[])
489c78c1075SJoseph Pusztay {
490c78c1075SJoseph Pusztay AppCtx *ctx;
491c78c1075SJoseph Pusztay PetscInt dim, Np;
492c78c1075SJoseph Pusztay
493c78c1075SJoseph Pusztay PetscFunctionBegin;
494c78c1075SJoseph Pusztay PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
495c78c1075SJoseph Pusztay PetscValidHeaderSpecific(sw, DM_CLASSID, 2);
4964f572ea9SToby Isaac PetscAssertPointer(E, 3);
497c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
498c78c1075SJoseph Pusztay PetscCall(DMSwarmGetLocalSize(sw, &Np));
499c78c1075SJoseph Pusztay PetscCall(DMGetApplicationContext(sw, &ctx));
500c78c1075SJoseph Pusztay PetscCall(PetscArrayzero(E, Np * dim));
501c78c1075SJoseph Pusztay
502c78c1075SJoseph Pusztay switch (ctx->em) {
503c78c1075SJoseph Pusztay case EM_PRIMAL:
504c78c1075SJoseph Pusztay PetscCall(ComputeFieldAtParticles_Primal(snes, sw, E));
505c78c1075SJoseph Pusztay break;
506c78c1075SJoseph Pusztay case EM_COULOMB:
507c78c1075SJoseph Pusztay PetscCall(ComputeFieldAtParticles_Coulomb(snes, sw, E));
508c78c1075SJoseph Pusztay break;
509c78c1075SJoseph Pusztay case EM_MIXED:
510c78c1075SJoseph Pusztay PetscCall(ComputeFieldAtParticles_Mixed(snes, sw, E));
511c78c1075SJoseph Pusztay break;
512c78c1075SJoseph Pusztay case EM_NONE:
513c78c1075SJoseph Pusztay break;
514c78c1075SJoseph Pusztay default:
515c78c1075SJoseph Pusztay SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No solver for electrostatic model %s", EMTypes[ctx->em]);
516c78c1075SJoseph Pusztay }
517c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
518c78c1075SJoseph Pusztay }
519c78c1075SJoseph Pusztay
RHSFunction(TS ts,PetscReal t,Vec U,Vec G,PetscCtx ctx)520*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec G, PetscCtx ctx)
521c78c1075SJoseph Pusztay {
522c78c1075SJoseph Pusztay DM sw;
523c78c1075SJoseph Pusztay SNES snes = ((AppCtx *)ctx)->snes;
524c78c1075SJoseph Pusztay const PetscReal *coords, *vel;
525c78c1075SJoseph Pusztay const PetscScalar *u;
526c78c1075SJoseph Pusztay PetscScalar *g;
527c78c1075SJoseph Pusztay PetscReal *E;
528c78c1075SJoseph Pusztay PetscInt dim, d, Np, p;
529c78c1075SJoseph Pusztay
530c78c1075SJoseph Pusztay PetscFunctionBeginUser;
531c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
532c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
533c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
534c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
535c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
536c78c1075SJoseph Pusztay PetscCall(VecGetLocalSize(U, &Np));
537c78c1075SJoseph Pusztay PetscCall(VecGetArrayRead(U, &u));
538c78c1075SJoseph Pusztay PetscCall(VecGetArray(G, &g));
539c78c1075SJoseph Pusztay
540c78c1075SJoseph Pusztay PetscCall(ComputeFieldAtParticles(snes, sw, E));
541c78c1075SJoseph Pusztay
542c78c1075SJoseph Pusztay Np /= 2 * dim;
543c78c1075SJoseph Pusztay for (p = 0; p < Np; ++p) {
544c78c1075SJoseph Pusztay const PetscReal x0 = coords[p * dim + 0];
545c78c1075SJoseph Pusztay const PetscReal vy0 = vel[p * dim + 1];
546c78c1075SJoseph Pusztay const PetscReal omega = vy0 / x0;
547c78c1075SJoseph Pusztay
548c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) {
549c78c1075SJoseph Pusztay g[(p * 2 + 0) * dim + d] = u[(p * 2 + 1) * dim + d];
550c78c1075SJoseph Pusztay g[(p * 2 + 1) * dim + d] = E[p * dim + d] - PetscSqr(omega) * u[(p * 2 + 0) * dim + d];
551c78c1075SJoseph Pusztay }
552c78c1075SJoseph Pusztay }
553c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
554c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
555c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
556c78c1075SJoseph Pusztay PetscCall(VecRestoreArrayRead(U, &u));
557c78c1075SJoseph Pusztay PetscCall(VecRestoreArray(G, &g));
558c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
559c78c1075SJoseph Pusztay }
560c78c1075SJoseph Pusztay
561c78c1075SJoseph Pusztay /* J_{ij} = dF_i/dx_j
562c78c1075SJoseph Pusztay J_p = ( 0 1)
563c78c1075SJoseph Pusztay (-w^2 0)
564c78c1075SJoseph Pusztay TODO Now there is another term with w^2 from the electric field. I think we will need to invert the operator.
565c78c1075SJoseph Pusztay Perhaps we can approximate the Jacobian using only the cellwise P-P gradient from Coulomb
566c78c1075SJoseph Pusztay */
RHSJacobian(TS ts,PetscReal t,Vec U,Mat J,Mat P,PetscCtx ctx)567*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat J, Mat P, PetscCtx ctx)
568c78c1075SJoseph Pusztay {
569c78c1075SJoseph Pusztay DM sw;
570c78c1075SJoseph Pusztay const PetscReal *coords, *vel;
571c78c1075SJoseph Pusztay PetscInt dim, d, Np, p, rStart;
572c78c1075SJoseph Pusztay
573c78c1075SJoseph Pusztay PetscFunctionBeginUser;
574c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
575c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
576c78c1075SJoseph Pusztay PetscCall(VecGetLocalSize(U, &Np));
577c78c1075SJoseph Pusztay PetscCall(MatGetOwnershipRange(J, &rStart, NULL));
578c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
579c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
580c78c1075SJoseph Pusztay Np /= 2 * dim;
581c78c1075SJoseph Pusztay for (p = 0; p < Np; ++p) {
582c78c1075SJoseph Pusztay const PetscReal x0 = coords[p * dim + 0];
583c78c1075SJoseph Pusztay const PetscReal vy0 = vel[p * dim + 1];
584c78c1075SJoseph Pusztay const PetscReal omega = vy0 / x0;
585c78c1075SJoseph Pusztay PetscScalar vals[4] = {0., 1., -PetscSqr(omega), 0.};
586c78c1075SJoseph Pusztay
587c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) {
588c78c1075SJoseph Pusztay const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart};
589c78c1075SJoseph Pusztay PetscCall(MatSetValues(J, 2, rows, 2, rows, vals, INSERT_VALUES));
590c78c1075SJoseph Pusztay }
591c78c1075SJoseph Pusztay }
592c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
593c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
594c78c1075SJoseph Pusztay PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
595c78c1075SJoseph Pusztay PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
596c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
597c78c1075SJoseph Pusztay }
598c78c1075SJoseph Pusztay
RHSFunctionX(TS ts,PetscReal t,Vec V,Vec Xres,PetscCtx ctx)599*2a8381b2SBarry Smith static PetscErrorCode RHSFunctionX(TS ts, PetscReal t, Vec V, Vec Xres, PetscCtx ctx)
600c78c1075SJoseph Pusztay {
601c78c1075SJoseph Pusztay DM sw;
602c78c1075SJoseph Pusztay const PetscScalar *v;
603c78c1075SJoseph Pusztay PetscScalar *xres;
604c78c1075SJoseph Pusztay PetscInt Np, p, dim, d;
605c78c1075SJoseph Pusztay
606c78c1075SJoseph Pusztay PetscFunctionBeginUser;
607c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
608c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
609c78c1075SJoseph Pusztay PetscCall(VecGetLocalSize(Xres, &Np));
610c78c1075SJoseph Pusztay Np /= dim;
611c78c1075SJoseph Pusztay PetscCall(VecGetArrayRead(V, &v));
612c78c1075SJoseph Pusztay PetscCall(VecGetArray(Xres, &xres));
613c78c1075SJoseph Pusztay for (p = 0; p < Np; ++p) {
614aa624791SPierre Jolivet for (d = 0; d < dim; ++d) xres[p * dim + d] = v[p * dim + d];
615c78c1075SJoseph Pusztay }
616c78c1075SJoseph Pusztay PetscCall(VecRestoreArrayRead(V, &v));
617c78c1075SJoseph Pusztay PetscCall(VecRestoreArray(Xres, &xres));
618c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
619c78c1075SJoseph Pusztay }
620c78c1075SJoseph Pusztay
RHSFunctionV(TS ts,PetscReal t,Vec X,Vec Vres,PetscCtx ctx)621*2a8381b2SBarry Smith static PetscErrorCode RHSFunctionV(TS ts, PetscReal t, Vec X, Vec Vres, PetscCtx ctx)
622c78c1075SJoseph Pusztay {
623c78c1075SJoseph Pusztay DM sw;
624c78c1075SJoseph Pusztay SNES snes = ((AppCtx *)ctx)->snes;
625c78c1075SJoseph Pusztay const PetscScalar *x;
626c78c1075SJoseph Pusztay const PetscReal *coords, *vel;
627c78c1075SJoseph Pusztay PetscScalar *vres;
628c78c1075SJoseph Pusztay PetscReal *E;
629c78c1075SJoseph Pusztay PetscInt Np, p, dim, d;
630c78c1075SJoseph Pusztay
631c78c1075SJoseph Pusztay PetscFunctionBeginUser;
632c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
633c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
634c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
635c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
636c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
637c78c1075SJoseph Pusztay PetscCall(VecGetLocalSize(Vres, &Np));
638c78c1075SJoseph Pusztay PetscCall(VecGetArrayRead(X, &x));
639c78c1075SJoseph Pusztay PetscCall(VecGetArray(Vres, &vres));
640c78c1075SJoseph Pusztay PetscCheck(dim == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension must be 2");
641c78c1075SJoseph Pusztay
642c78c1075SJoseph Pusztay PetscCall(ComputeFieldAtParticles(snes, sw, E));
643c78c1075SJoseph Pusztay
644c78c1075SJoseph Pusztay Np /= dim;
645c78c1075SJoseph Pusztay for (p = 0; p < Np; ++p) {
646c78c1075SJoseph Pusztay const PetscReal x0 = coords[p * dim + 0];
647c78c1075SJoseph Pusztay const PetscReal vy0 = vel[p * dim + 1];
648c78c1075SJoseph Pusztay const PetscReal omega = vy0 / x0;
649c78c1075SJoseph Pusztay
650c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) vres[p * dim + d] = E[p * dim + d] - PetscSqr(omega) * x[p * dim + d];
651c78c1075SJoseph Pusztay }
652c78c1075SJoseph Pusztay PetscCall(VecRestoreArrayRead(X, &x));
653c78c1075SJoseph Pusztay PetscCall(VecRestoreArray(Vres, &vres));
654c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
655c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
656c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
657c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
658c78c1075SJoseph Pusztay }
659c78c1075SJoseph Pusztay
660f940b0e3Sdanofinn /* Discrete Gradients Formulation: S, F, gradF (G) */
RHSJacobianS(TS ts,PetscReal t,Vec U,Mat S,PetscCtx ctx)661*2a8381b2SBarry Smith PetscErrorCode RHSJacobianS(TS ts, PetscReal t, Vec U, Mat S, PetscCtx ctx)
662f940b0e3Sdanofinn {
663f940b0e3Sdanofinn PetscScalar vals[4] = {0., 1., -1., 0.};
664f940b0e3Sdanofinn DM sw;
665f940b0e3Sdanofinn PetscInt dim, d, Np, p, rStart;
666f940b0e3Sdanofinn
667f940b0e3Sdanofinn PetscFunctionBeginUser;
668f940b0e3Sdanofinn PetscCall(TSGetDM(ts, &sw));
669f940b0e3Sdanofinn PetscCall(DMGetDimension(sw, &dim));
670f940b0e3Sdanofinn PetscCall(VecGetLocalSize(U, &Np));
671f940b0e3Sdanofinn PetscCall(MatGetOwnershipRange(S, &rStart, NULL));
672f940b0e3Sdanofinn Np /= 2 * dim;
673f940b0e3Sdanofinn for (p = 0; p < Np; ++p) {
674f940b0e3Sdanofinn for (d = 0; d < dim; ++d) {
675f940b0e3Sdanofinn const PetscInt rows[2] = {(p * 2 + 0) * dim + d + rStart, (p * 2 + 1) * dim + d + rStart};
676f940b0e3Sdanofinn PetscCall(MatSetValues(S, 2, rows, 2, rows, vals, INSERT_VALUES));
677f940b0e3Sdanofinn }
678f940b0e3Sdanofinn }
679f940b0e3Sdanofinn PetscCall(MatAssemblyBegin(S, MAT_FINAL_ASSEMBLY));
680f940b0e3Sdanofinn PetscCall(MatAssemblyEnd(S, MAT_FINAL_ASSEMBLY));
681f940b0e3Sdanofinn PetscFunctionReturn(PETSC_SUCCESS);
682f940b0e3Sdanofinn }
683f940b0e3Sdanofinn
RHSObjectiveF(TS ts,PetscReal t,Vec U,PetscScalar * F,PetscCtx ctx)684*2a8381b2SBarry Smith PetscErrorCode RHSObjectiveF(TS ts, PetscReal t, Vec U, PetscScalar *F, PetscCtx ctx)
685f940b0e3Sdanofinn {
686f940b0e3Sdanofinn SNES snes = ((AppCtx *)ctx)->snes;
687f940b0e3Sdanofinn DM dm, sw;
688f940b0e3Sdanofinn const PetscScalar *u, *phi_vals;
689f940b0e3Sdanofinn PetscInt dim, Np, cStart, cEnd;
690f940b0e3Sdanofinn PetscReal *vel, *coords, m_p = 1., q_p = -1.;
691f940b0e3Sdanofinn Vec phi;
692f940b0e3Sdanofinn
693f940b0e3Sdanofinn PetscFunctionBeginUser;
694f940b0e3Sdanofinn PetscCall(TSGetDM(ts, &sw));
695f940b0e3Sdanofinn PetscCall(DMGetDimension(sw, &dim));
696f940b0e3Sdanofinn PetscCall(SNESGetDM(snes, &dm));
697f940b0e3Sdanofinn PetscCall(VecGetArrayRead(U, &u));
698f940b0e3Sdanofinn PetscCall(VecGetLocalSize(U, &Np));
699f940b0e3Sdanofinn PetscCall(DMGetGlobalVector(dm, &phi));
700f940b0e3Sdanofinn PetscCall(VecViewFromOptions(phi, NULL, "-phi_view_dg"));
701f940b0e3Sdanofinn PetscCall(PetscObjectSetName((PetscObject)phi, "potential"));
702f940b0e3Sdanofinn PetscInt phi_size;
703f940b0e3Sdanofinn PetscCall(VecGetSize(phi, &phi_size));
704f940b0e3Sdanofinn PetscCall(VecGetArrayRead(phi, &phi_vals));
705f940b0e3Sdanofinn PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
706f940b0e3Sdanofinn PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
707f940b0e3Sdanofinn
708f940b0e3Sdanofinn PetscCall(DMSwarmSortGetAccess(sw));
709f940b0e3Sdanofinn PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
710f940b0e3Sdanofinn Np /= 2 * dim;
711f940b0e3Sdanofinn for (PetscInt c = cStart; c < cEnd; ++c) {
712f940b0e3Sdanofinn PetscInt *points;
713f940b0e3Sdanofinn PetscInt Ncp;
714f940b0e3Sdanofinn PetscReal E = phi_vals[c];
715f940b0e3Sdanofinn
716f940b0e3Sdanofinn PetscCall(DMSwarmSortGetPointsPerCell(sw, c, &Ncp, &points));
717f940b0e3Sdanofinn for (PetscInt cp = 0; cp < Ncp; ++cp) {
718f940b0e3Sdanofinn const PetscInt p = points[cp];
719f940b0e3Sdanofinn const PetscReal x0 = coords[p * dim + 0];
720f940b0e3Sdanofinn const PetscReal vy0 = vel[p * dim + 1];
721f940b0e3Sdanofinn const PetscReal omega = vy0 / x0;
722f940b0e3Sdanofinn const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &u[(p * 2 + 1) * dim]);
723f940b0e3Sdanofinn const PetscReal x2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 0) * dim], &u[(p * 2 + 0) * dim]);
724f940b0e3Sdanofinn E += 0.5 * q_p * m_p * (v2) + 0.5 * PetscSqr(omega) * (x2);
725f940b0e3Sdanofinn
726f940b0e3Sdanofinn *F += E;
727f940b0e3Sdanofinn }
728f940b0e3Sdanofinn PetscCall(DMSwarmSortRestorePointsPerCell(sw, c, &Ncp, &points));
729f940b0e3Sdanofinn }
730f940b0e3Sdanofinn PetscCall(DMSwarmSortRestoreAccess(sw));
731f940b0e3Sdanofinn PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
732f940b0e3Sdanofinn PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
733f940b0e3Sdanofinn PetscCall(VecRestoreArrayRead(phi, &phi_vals));
734f940b0e3Sdanofinn PetscCall(DMRestoreGlobalVector(dm, &phi));
735f940b0e3Sdanofinn // PetscCall(DMSwarmRestoreField(sw, "potential", NULL, NULL, (void **)&pot));
736f940b0e3Sdanofinn PetscCall(VecRestoreArrayRead(U, &u));
737f940b0e3Sdanofinn PetscFunctionReturn(PETSC_SUCCESS);
738f940b0e3Sdanofinn }
739f940b0e3Sdanofinn
740f940b0e3Sdanofinn /* dF/dx = q E dF/dv = v */
RHSFunctionG(TS ts,PetscReal t,Vec U,Vec G,PetscCtx ctx)741*2a8381b2SBarry Smith PetscErrorCode RHSFunctionG(TS ts, PetscReal t, Vec U, Vec G, PetscCtx ctx)
742f940b0e3Sdanofinn {
743f940b0e3Sdanofinn DM sw;
744f940b0e3Sdanofinn SNES snes = ((AppCtx *)ctx)->snes;
745f940b0e3Sdanofinn const PetscReal *coords, *vel;
746f940b0e3Sdanofinn const PetscScalar *u;
747f940b0e3Sdanofinn PetscScalar *g;
748f940b0e3Sdanofinn PetscReal *E, m_p = 1., q_p = -1.;
749f940b0e3Sdanofinn PetscInt dim, d, Np, p;
750f940b0e3Sdanofinn
751f940b0e3Sdanofinn PetscFunctionBeginUser;
752f940b0e3Sdanofinn PetscCall(TSGetDM(ts, &sw));
753f940b0e3Sdanofinn PetscCall(DMGetDimension(sw, &dim));
754f940b0e3Sdanofinn PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
755f940b0e3Sdanofinn PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
756f940b0e3Sdanofinn PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
757f940b0e3Sdanofinn PetscCall(DMSwarmGetLocalSize(sw, &Np));
758f940b0e3Sdanofinn PetscCall(VecGetArrayRead(U, &u));
759f940b0e3Sdanofinn PetscCall(VecGetArray(G, &g));
760f940b0e3Sdanofinn
761f940b0e3Sdanofinn int COMPUTEFIELD;
762f940b0e3Sdanofinn PetscCall(PetscLogEventRegister("COMPFIELDATPART", TS_CLASSID, &COMPUTEFIELD));
763f940b0e3Sdanofinn PetscCall(PetscLogEventBegin(COMPUTEFIELD, 0, 0, 0, 0));
764f940b0e3Sdanofinn PetscCall(ComputeFieldAtParticles(snes, sw, E));
765f940b0e3Sdanofinn PetscCall(PetscLogEventEnd(COMPUTEFIELD, 0, 0, 0, 0));
766f940b0e3Sdanofinn for (p = 0; p < Np; ++p) {
767f940b0e3Sdanofinn const PetscReal x0 = coords[p * dim + 0];
768f940b0e3Sdanofinn const PetscReal vy0 = vel[p * dim + 1];
769f940b0e3Sdanofinn const PetscReal omega = vy0 / x0;
770f940b0e3Sdanofinn for (d = 0; d < dim; ++d) {
771f940b0e3Sdanofinn g[(p * 2 + 0) * dim + d] = -(q_p / m_p) * E[p * dim + d] + PetscSqr(omega) * u[(p * 2 + 0) * dim + d];
772f940b0e3Sdanofinn g[(p * 2 + 1) * dim + d] = u[(p * 2 + 1) * dim + d];
773f940b0e3Sdanofinn }
774f940b0e3Sdanofinn }
775f940b0e3Sdanofinn PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
776f940b0e3Sdanofinn PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
777f940b0e3Sdanofinn PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
778f940b0e3Sdanofinn PetscCall(VecRestoreArrayRead(U, &u));
779f940b0e3Sdanofinn PetscCall(VecRestoreArray(G, &g));
780f940b0e3Sdanofinn PetscFunctionReturn(PETSC_SUCCESS);
781f940b0e3Sdanofinn }
782f940b0e3Sdanofinn
CreateSolution(TS ts)783c78c1075SJoseph Pusztay static PetscErrorCode CreateSolution(TS ts)
784c78c1075SJoseph Pusztay {
785c78c1075SJoseph Pusztay DM sw;
786c78c1075SJoseph Pusztay Vec u;
787c78c1075SJoseph Pusztay PetscInt dim, Np;
788c78c1075SJoseph Pusztay
789c78c1075SJoseph Pusztay PetscFunctionBegin;
790c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
791c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
792c78c1075SJoseph Pusztay PetscCall(DMSwarmGetLocalSize(sw, &Np));
793c78c1075SJoseph Pusztay PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
794c78c1075SJoseph Pusztay PetscCall(VecSetBlockSize(u, dim));
795c78c1075SJoseph Pusztay PetscCall(VecSetSizes(u, 2 * Np * dim, PETSC_DECIDE));
796c78c1075SJoseph Pusztay PetscCall(VecSetUp(u));
797c78c1075SJoseph Pusztay PetscCall(TSSetSolution(ts, u));
798c78c1075SJoseph Pusztay PetscCall(VecDestroy(&u));
799c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
800c78c1075SJoseph Pusztay }
801c78c1075SJoseph Pusztay
SetProblem(TS ts)802c78c1075SJoseph Pusztay static PetscErrorCode SetProblem(TS ts)
803c78c1075SJoseph Pusztay {
804c78c1075SJoseph Pusztay AppCtx *user;
805c78c1075SJoseph Pusztay DM sw;
806c78c1075SJoseph Pusztay
807c78c1075SJoseph Pusztay PetscFunctionBegin;
808c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
809*2a8381b2SBarry Smith PetscCall(DMGetApplicationContext(sw, &user));
810c78c1075SJoseph Pusztay // Define unified system for (X, V)
811c78c1075SJoseph Pusztay {
812c78c1075SJoseph Pusztay Mat J;
813c78c1075SJoseph Pusztay PetscInt dim, Np;
814c78c1075SJoseph Pusztay
815c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
816c78c1075SJoseph Pusztay PetscCall(DMSwarmGetLocalSize(sw, &Np));
817c78c1075SJoseph Pusztay PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
818c78c1075SJoseph Pusztay PetscCall(MatSetSizes(J, 2 * Np * dim, 2 * Np * dim, PETSC_DECIDE, PETSC_DECIDE));
819c78c1075SJoseph Pusztay PetscCall(MatSetBlockSize(J, 2 * dim));
820c78c1075SJoseph Pusztay PetscCall(MatSetFromOptions(J));
821c78c1075SJoseph Pusztay PetscCall(MatSetUp(J));
822c78c1075SJoseph Pusztay PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, user));
823c78c1075SJoseph Pusztay PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, user));
824c78c1075SJoseph Pusztay PetscCall(MatDestroy(&J));
825c78c1075SJoseph Pusztay }
826c78c1075SJoseph Pusztay /* Define split system for X and V */
827c78c1075SJoseph Pusztay {
828c78c1075SJoseph Pusztay Vec u;
829c78c1075SJoseph Pusztay IS isx, isv, istmp;
830c78c1075SJoseph Pusztay const PetscInt *idx;
831c78c1075SJoseph Pusztay PetscInt dim, Np, rstart;
832c78c1075SJoseph Pusztay
833c78c1075SJoseph Pusztay PetscCall(TSGetSolution(ts, &u));
834c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
835c78c1075SJoseph Pusztay PetscCall(DMSwarmGetLocalSize(sw, &Np));
836c78c1075SJoseph Pusztay PetscCall(VecGetOwnershipRange(u, &rstart, NULL));
837c78c1075SJoseph Pusztay PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 0, 2, &istmp));
838c78c1075SJoseph Pusztay PetscCall(ISGetIndices(istmp, &idx));
839c78c1075SJoseph Pusztay PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isx));
840c78c1075SJoseph Pusztay PetscCall(ISRestoreIndices(istmp, &idx));
841c78c1075SJoseph Pusztay PetscCall(ISDestroy(&istmp));
842c78c1075SJoseph Pusztay PetscCall(ISCreateStride(PETSC_COMM_WORLD, Np, (rstart / dim) + 1, 2, &istmp));
843c78c1075SJoseph Pusztay PetscCall(ISGetIndices(istmp, &idx));
844c78c1075SJoseph Pusztay PetscCall(ISCreateBlock(PETSC_COMM_WORLD, dim, Np, idx, PETSC_COPY_VALUES, &isv));
845c78c1075SJoseph Pusztay PetscCall(ISRestoreIndices(istmp, &idx));
846c78c1075SJoseph Pusztay PetscCall(ISDestroy(&istmp));
847c78c1075SJoseph Pusztay PetscCall(TSRHSSplitSetIS(ts, "position", isx));
848c78c1075SJoseph Pusztay PetscCall(TSRHSSplitSetIS(ts, "momentum", isv));
849c78c1075SJoseph Pusztay PetscCall(ISDestroy(&isx));
850c78c1075SJoseph Pusztay PetscCall(ISDestroy(&isv));
851c78c1075SJoseph Pusztay PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunctionX, user));
852c78c1075SJoseph Pusztay PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunctionV, user));
853c78c1075SJoseph Pusztay }
854f940b0e3Sdanofinn // Define symplectic formulation U_t = S . G, where G = grad F
855f940b0e3Sdanofinn {
856f940b0e3Sdanofinn PetscCall(TSDiscGradSetFormulation(ts, RHSJacobianS, RHSObjectiveF, RHSFunctionG, user));
857f940b0e3Sdanofinn }
858c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
859c78c1075SJoseph Pusztay }
860c78c1075SJoseph Pusztay
circleSingleX(PetscInt dim,PetscReal time,const PetscReal unused[],PetscInt p,PetscScalar x[],PetscCtx ctx)861*2a8381b2SBarry Smith PetscErrorCode circleSingleX(PetscInt dim, PetscReal time, const PetscReal unused[], PetscInt p, PetscScalar x[], PetscCtx ctx)
862c78c1075SJoseph Pusztay {
863c78c1075SJoseph Pusztay x[0] = p + 1.;
864c78c1075SJoseph Pusztay x[1] = 0.;
865c78c1075SJoseph Pusztay return PETSC_SUCCESS;
866c78c1075SJoseph Pusztay }
867c78c1075SJoseph Pusztay
circleSingleV(PetscInt dim,PetscReal time,const PetscReal unused[],PetscInt p,PetscScalar v[],PetscCtx ctx)868*2a8381b2SBarry Smith PetscErrorCode circleSingleV(PetscInt dim, PetscReal time, const PetscReal unused[], PetscInt p, PetscScalar v[], PetscCtx ctx)
869c78c1075SJoseph Pusztay {
870c78c1075SJoseph Pusztay v[0] = 0.;
871c78c1075SJoseph Pusztay v[1] = PetscSqrtReal(1000. / (p + 1.));
872c78c1075SJoseph Pusztay return PETSC_SUCCESS;
873c78c1075SJoseph Pusztay }
874c78c1075SJoseph Pusztay
875c78c1075SJoseph Pusztay /* Put 5 particles into each circle */
circleMultipleX(PetscInt dim,PetscReal time,const PetscReal unused[],PetscInt p,PetscScalar x[],PetscCtx ctx)876*2a8381b2SBarry Smith PetscErrorCode circleMultipleX(PetscInt dim, PetscReal time, const PetscReal unused[], PetscInt p, PetscScalar x[], PetscCtx ctx)
877c78c1075SJoseph Pusztay {
878c78c1075SJoseph Pusztay const PetscInt n = 5;
879c78c1075SJoseph Pusztay const PetscReal r0 = (p / n) + 1.;
880c78c1075SJoseph Pusztay const PetscReal th0 = (2. * PETSC_PI * (p % n)) / n;
881c78c1075SJoseph Pusztay
882c78c1075SJoseph Pusztay x[0] = r0 * PetscCosReal(th0);
883c78c1075SJoseph Pusztay x[1] = r0 * PetscSinReal(th0);
884c78c1075SJoseph Pusztay return PETSC_SUCCESS;
885c78c1075SJoseph Pusztay }
886c78c1075SJoseph Pusztay
887c78c1075SJoseph Pusztay /* Put 5 particles into each circle */
circleMultipleV(PetscInt dim,PetscReal time,const PetscReal unused[],PetscInt p,PetscScalar v[],PetscCtx ctx)888*2a8381b2SBarry Smith PetscErrorCode circleMultipleV(PetscInt dim, PetscReal time, const PetscReal unused[], PetscInt p, PetscScalar v[], PetscCtx ctx)
889c78c1075SJoseph Pusztay {
890c78c1075SJoseph Pusztay const PetscInt n = 5;
891c78c1075SJoseph Pusztay const PetscReal r0 = (p / n) + 1.;
892c78c1075SJoseph Pusztay const PetscReal th0 = (2. * PETSC_PI * (p % n)) / n;
893c78c1075SJoseph Pusztay const PetscReal omega = PetscSqrtReal(1000. / r0) / r0;
894c78c1075SJoseph Pusztay
895c78c1075SJoseph Pusztay v[0] = -r0 * omega * PetscSinReal(th0);
896c78c1075SJoseph Pusztay v[1] = r0 * omega * PetscCosReal(th0);
897c78c1075SJoseph Pusztay return PETSC_SUCCESS;
898c78c1075SJoseph Pusztay }
899c78c1075SJoseph Pusztay
900c78c1075SJoseph Pusztay /*
901c78c1075SJoseph Pusztay InitializeSolveAndSwarm - Set the solution values to the swarm coordinates and velocities, and also possibly set the initial values.
902c78c1075SJoseph Pusztay
903c78c1075SJoseph Pusztay Input Parameters:
904c78c1075SJoseph Pusztay + ts - The TS
905aaa8cc7dSPierre Jolivet - useInitial - Flag to also set the initial conditions to the current coordinates and velocities and setup the problem
906c78c1075SJoseph Pusztay
9072fe279fdSBarry Smith Output Parameter:
908c78c1075SJoseph Pusztay . u - The initialized solution vector
909c78c1075SJoseph Pusztay
910c78c1075SJoseph Pusztay Level: advanced
911c78c1075SJoseph Pusztay
912c78c1075SJoseph Pusztay .seealso: InitializeSolve()
913c78c1075SJoseph Pusztay */
InitializeSolveAndSwarm(TS ts,PetscBool useInitial)914c78c1075SJoseph Pusztay static PetscErrorCode InitializeSolveAndSwarm(TS ts, PetscBool useInitial)
915c78c1075SJoseph Pusztay {
916c78c1075SJoseph Pusztay DM sw;
917c78c1075SJoseph Pusztay Vec u, gc, gv, gc0, gv0;
918c78c1075SJoseph Pusztay IS isx, isv;
919c78c1075SJoseph Pusztay AppCtx *user;
920c78c1075SJoseph Pusztay
921c78c1075SJoseph Pusztay PetscFunctionBeginUser;
922c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
923c78c1075SJoseph Pusztay PetscCall(DMGetApplicationContext(sw, &user));
924c78c1075SJoseph Pusztay if (useInitial) {
925c78c1075SJoseph Pusztay PetscReal v0[1] = {1.};
926c78c1075SJoseph Pusztay
927c78c1075SJoseph Pusztay PetscCall(DMSwarmInitializeCoordinates(sw));
928c78c1075SJoseph Pusztay PetscCall(DMSwarmInitializeVelocitiesFromOptions(sw, v0));
929c78c1075SJoseph Pusztay PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
9304a2752e6SStefano Zampini PetscCall(TSReset(ts));
9313dd8e6d7SStefano Zampini PetscCall(CreateSolution(ts));
9323dd8e6d7SStefano Zampini PetscCall(SetProblem(ts));
933c78c1075SJoseph Pusztay }
934c78c1075SJoseph Pusztay PetscCall(TSGetSolution(ts, &u));
935c78c1075SJoseph Pusztay PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
936c78c1075SJoseph Pusztay PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
937c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
938c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initCoordinates", &gc0));
939c78c1075SJoseph Pusztay if (useInitial) PetscCall(VecCopy(gc, gc0));
940c78c1075SJoseph Pusztay PetscCall(VecISCopy(u, isx, SCATTER_FORWARD, gc));
941c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
942c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initCoordinates", &gc0));
943c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
944c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "initVelocity", &gv0));
945c78c1075SJoseph Pusztay if (useInitial) PetscCall(VecCopy(gv, gv0));
946c78c1075SJoseph Pusztay PetscCall(VecISCopy(u, isv, SCATTER_FORWARD, gv));
947c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
948c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "initVelocity", &gv0));
949c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
950c78c1075SJoseph Pusztay }
951c78c1075SJoseph Pusztay
InitializeSolve(TS ts,Vec u)952c78c1075SJoseph Pusztay static PetscErrorCode InitializeSolve(TS ts, Vec u)
953c78c1075SJoseph Pusztay {
954c78c1075SJoseph Pusztay PetscFunctionBegin;
955c78c1075SJoseph Pusztay PetscCall(TSSetSolution(ts, u));
956c78c1075SJoseph Pusztay PetscCall(InitializeSolveAndSwarm(ts, PETSC_TRUE));
957c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
958c78c1075SJoseph Pusztay }
959c78c1075SJoseph Pusztay
ComputeError(TS ts,Vec U,Vec E)960c78c1075SJoseph Pusztay static PetscErrorCode ComputeError(TS ts, Vec U, Vec E)
961c78c1075SJoseph Pusztay {
962c78c1075SJoseph Pusztay MPI_Comm comm;
963c78c1075SJoseph Pusztay DM sw;
964c78c1075SJoseph Pusztay AppCtx *user;
965c78c1075SJoseph Pusztay const PetscScalar *u;
966c78c1075SJoseph Pusztay const PetscReal *coords, *vel;
967c78c1075SJoseph Pusztay PetscScalar *e;
968c78c1075SJoseph Pusztay PetscReal t;
969c78c1075SJoseph Pusztay PetscInt dim, Np, p;
970c78c1075SJoseph Pusztay
971c78c1075SJoseph Pusztay PetscFunctionBeginUser;
972c78c1075SJoseph Pusztay PetscCall(PetscObjectGetComm((PetscObject)ts, &comm));
973c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
974c78c1075SJoseph Pusztay PetscCall(DMGetApplicationContext(sw, &user));
975c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
976c78c1075SJoseph Pusztay PetscCall(TSGetSolveTime(ts, &t));
977c78c1075SJoseph Pusztay PetscCall(VecGetArray(E, &e));
978c78c1075SJoseph Pusztay PetscCall(VecGetArrayRead(U, &u));
979c78c1075SJoseph Pusztay PetscCall(VecGetLocalSize(U, &Np));
980c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
981c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "initVelocity", NULL, NULL, (void **)&vel));
982c78c1075SJoseph Pusztay Np /= 2 * dim;
983c78c1075SJoseph Pusztay for (p = 0; p < Np; ++p) {
984c78c1075SJoseph Pusztay /* TODO generalize initial conditions and project into plane instead of assuming x-y */
985c78c1075SJoseph Pusztay const PetscReal r0 = DMPlex_NormD_Internal(dim, &coords[p * dim]);
986c78c1075SJoseph Pusztay const PetscReal th0 = PetscAtan2Real(coords[p * dim + 1], coords[p * dim + 0]);
987c78c1075SJoseph Pusztay const PetscReal v0 = DMPlex_NormD_Internal(dim, &vel[p * dim]);
988c78c1075SJoseph Pusztay const PetscReal omega = v0 / r0;
989c78c1075SJoseph Pusztay const PetscReal ct = PetscCosReal(omega * t + th0);
990c78c1075SJoseph Pusztay const PetscReal st = PetscSinReal(omega * t + th0);
991c78c1075SJoseph Pusztay const PetscScalar *x = &u[(p * 2 + 0) * dim];
992c78c1075SJoseph Pusztay const PetscScalar *v = &u[(p * 2 + 1) * dim];
993c78c1075SJoseph Pusztay const PetscReal xe[3] = {r0 * ct, r0 * st, 0.0};
994c78c1075SJoseph Pusztay const PetscReal ve[3] = {-v0 * st, v0 * ct, 0.0};
995c78c1075SJoseph Pusztay PetscInt d;
996c78c1075SJoseph Pusztay
997c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) {
998c78c1075SJoseph Pusztay e[(p * 2 + 0) * dim + d] = x[d] - xe[d];
999c78c1075SJoseph Pusztay e[(p * 2 + 1) * dim + d] = v[d] - ve[d];
1000c78c1075SJoseph Pusztay }
1001c78c1075SJoseph Pusztay if (user->error) {
1002c78c1075SJoseph Pusztay const PetscReal en = 0.5 * DMPlex_DotRealD_Internal(dim, v, v);
1003c78c1075SJoseph Pusztay const PetscReal exen = 0.5 * PetscSqr(v0);
10043771a240SStefano 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)));
1005c78c1075SJoseph Pusztay }
1006c78c1075SJoseph Pusztay }
1007c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initCoordinates", NULL, NULL, (void **)&coords));
1008c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "initVelocity", NULL, NULL, (void **)&vel));
1009c78c1075SJoseph Pusztay PetscCall(VecRestoreArrayRead(U, &u));
1010c78c1075SJoseph Pusztay PetscCall(VecRestoreArray(E, &e));
1011c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
1012c78c1075SJoseph Pusztay }
1013c78c1075SJoseph Pusztay
EnergyMonitor(TS ts,PetscInt step,PetscReal t,Vec U,PetscCtx ctx)1014*2a8381b2SBarry Smith static PetscErrorCode EnergyMonitor(TS ts, PetscInt step, PetscReal t, Vec U, PetscCtx ctx)
1015c78c1075SJoseph Pusztay {
1016c78c1075SJoseph Pusztay const PetscInt ostep = ((AppCtx *)ctx)->ostep;
1017c78c1075SJoseph Pusztay const EMType em = ((AppCtx *)ctx)->em;
1018c78c1075SJoseph Pusztay DM sw;
1019c78c1075SJoseph Pusztay const PetscScalar *u;
1020c78c1075SJoseph Pusztay PetscReal *coords, *E;
1021c78c1075SJoseph Pusztay PetscReal enKin = 0., enEM = 0.;
1022c78c1075SJoseph Pusztay PetscInt dim, d, Np, p, q;
1023c78c1075SJoseph Pusztay
1024c78c1075SJoseph Pusztay PetscFunctionBeginUser;
1025c78c1075SJoseph Pusztay if (step % ostep == 0) {
1026c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
1027c78c1075SJoseph Pusztay PetscCall(DMGetDimension(sw, &dim));
1028c78c1075SJoseph Pusztay PetscCall(VecGetArrayRead(U, &u));
1029c78c1075SJoseph Pusztay PetscCall(VecGetLocalSize(U, &Np));
1030c78c1075SJoseph Pusztay Np /= 2 * dim;
1031c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1032c78c1075SJoseph Pusztay PetscCall(DMSwarmGetField(sw, "E_field", NULL, NULL, (void **)&E));
1033c78c1075SJoseph Pusztay if (!step) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "Time Step Part Energy\n"));
1034c78c1075SJoseph Pusztay for (p = 0; p < Np; ++p) {
1035c78c1075SJoseph Pusztay const PetscReal v2 = DMPlex_DotRealD_Internal(dim, &u[(p * 2 + 1) * dim], &u[(p * 2 + 1) * dim]);
1036c78c1075SJoseph Pusztay PetscReal *pcoord = &coords[p * dim];
1037c78c1075SJoseph Pusztay
1038c78c1075SJoseph Pusztay PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " %5" PetscInt_FMT " %10.4lf\n", (double)t, step, p, (double)(0.5 * v2)));
1039c78c1075SJoseph Pusztay enKin += 0.5 * v2;
1040c78c1075SJoseph Pusztay if (em == EM_NONE) {
1041c78c1075SJoseph Pusztay continue;
1042c78c1075SJoseph Pusztay } else if (em == EM_COULOMB) {
1043c78c1075SJoseph Pusztay for (q = p + 1; q < Np; ++q) {
1044c78c1075SJoseph Pusztay PetscReal *qcoord = &coords[q * dim];
1045c78c1075SJoseph Pusztay PetscReal rpq[3], r;
1046c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) rpq[d] = pcoord[d] - qcoord[d];
1047c78c1075SJoseph Pusztay r = DMPlex_NormD_Internal(dim, rpq);
1048c78c1075SJoseph Pusztay enEM += 1. / r;
1049c78c1075SJoseph Pusztay }
1050c78c1075SJoseph Pusztay } else if (em == EM_PRIMAL || em == EM_MIXED) {
1051c78c1075SJoseph Pusztay for (d = 0; d < dim; ++d) enEM += E[p * dim + d];
1052c78c1075SJoseph Pusztay }
1053c78c1075SJoseph Pusztay }
1054c78c1075SJoseph Pusztay PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " KE\t %10.4lf\n", (double)t, step, (double)enKin));
1055c78c1075SJoseph Pusztay PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " PE\t %1.10g\n", (double)t, step, (double)enEM));
1056c78c1075SJoseph Pusztay PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject)ts), "%.6lf %4" PetscInt_FMT " E\t %10.4lf\n", (double)t, step, (double)(enKin + enEM)));
1057c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, DMSwarmPICField_coor, NULL, NULL, (void **)&coords));
1058c78c1075SJoseph Pusztay PetscCall(DMSwarmRestoreField(sw, "E_field", NULL, NULL, (void **)&E));
1059c78c1075SJoseph Pusztay PetscCall(PetscSynchronizedFlush(PetscObjectComm((PetscObject)ts), NULL));
1060c78c1075SJoseph Pusztay PetscCall(VecRestoreArrayRead(U, &u));
1061c78c1075SJoseph Pusztay }
1062c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
1063c78c1075SJoseph Pusztay }
1064c78c1075SJoseph Pusztay
SetUpMigrateParticles(TS ts,PetscInt n,PetscReal t,Vec x,PetscBool * flg,PetscCtx ctx)1065*2a8381b2SBarry Smith static PetscErrorCode SetUpMigrateParticles(TS ts, PetscInt n, PetscReal t, Vec x, PetscBool *flg, PetscCtx ctx)
1066c78c1075SJoseph Pusztay {
1067c78c1075SJoseph Pusztay DM sw;
1068c78c1075SJoseph Pusztay
1069c78c1075SJoseph Pusztay PetscFunctionBeginUser;
10703dd8e6d7SStefano Zampini *flg = PETSC_TRUE;
1071c78c1075SJoseph Pusztay PetscCall(TSGetDM(ts, &sw));
1072c78c1075SJoseph Pusztay PetscCall(DMViewFromOptions(sw, NULL, "-migrate_view_pre"));
1073c78c1075SJoseph Pusztay {
1074c78c1075SJoseph Pusztay Vec u, gc, gv;
1075c78c1075SJoseph Pusztay IS isx, isv;
1076c78c1075SJoseph Pusztay
1077c78c1075SJoseph Pusztay PetscCall(TSGetSolution(ts, &u));
1078c78c1075SJoseph Pusztay PetscCall(TSRHSSplitGetIS(ts, "position", &isx));
1079c78c1075SJoseph Pusztay PetscCall(TSRHSSplitGetIS(ts, "momentum", &isv));
1080c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
1081c78c1075SJoseph Pusztay PetscCall(VecISCopy(u, isx, SCATTER_REVERSE, gc));
1082c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, DMSwarmPICField_coor, &gc));
1083c78c1075SJoseph Pusztay PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "velocity", &gv));
1084c78c1075SJoseph Pusztay PetscCall(VecISCopy(u, isv, SCATTER_REVERSE, gv));
1085c78c1075SJoseph Pusztay PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "velocity", &gv));
1086c78c1075SJoseph Pusztay }
10873dd8e6d7SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
10883dd8e6d7SStefano Zampini }
10893dd8e6d7SStefano Zampini
MigrateParticles(TS ts,PetscInt nv,Vec vecsin[],Vec vecsout[],PetscCtx ctx)1090*2a8381b2SBarry Smith static PetscErrorCode MigrateParticles(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], PetscCtx ctx)
10913dd8e6d7SStefano Zampini {
10923dd8e6d7SStefano Zampini DM sw;
10933dd8e6d7SStefano Zampini
10943dd8e6d7SStefano Zampini PetscFunctionBeginUser;
10953dd8e6d7SStefano Zampini PetscCall(TSGetDM(ts, &sw));
1096c78c1075SJoseph Pusztay PetscCall(DMSwarmMigrate(sw, PETSC_TRUE));
10973dd8e6d7SStefano Zampini PetscCall(CreateSolution(ts));
10983dd8e6d7SStefano Zampini PetscCall(SetProblem(ts));
1099c78c1075SJoseph Pusztay PetscCall(InitializeSolveAndSwarm(ts, PETSC_FALSE));
1100c78c1075SJoseph Pusztay PetscFunctionReturn(PETSC_SUCCESS);
1101c78c1075SJoseph Pusztay }
1102c78c1075SJoseph Pusztay
main(int argc,char ** argv)1103c78c1075SJoseph Pusztay int main(int argc, char **argv)
1104c78c1075SJoseph Pusztay {
1105c78c1075SJoseph Pusztay DM dm, sw;
1106c78c1075SJoseph Pusztay TS ts;
1107c78c1075SJoseph Pusztay Vec u;
1108c78c1075SJoseph Pusztay AppCtx user;
1109c78c1075SJoseph Pusztay
1110c78c1075SJoseph Pusztay PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1111c78c1075SJoseph Pusztay PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1112c78c1075SJoseph Pusztay PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1113c78c1075SJoseph Pusztay PetscCall(CreatePoisson(dm, &user));
1114c78c1075SJoseph Pusztay PetscCall(CreateSwarm(dm, &user, &sw));
1115c78c1075SJoseph Pusztay PetscCall(DMSetApplicationContext(sw, &user));
1116c78c1075SJoseph Pusztay
1117c78c1075SJoseph Pusztay PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1118c78c1075SJoseph Pusztay PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1119c78c1075SJoseph Pusztay PetscCall(TSSetDM(ts, sw));
1120c78c1075SJoseph Pusztay PetscCall(TSSetMaxTime(ts, 0.1));
1121c78c1075SJoseph Pusztay PetscCall(TSSetTimeStep(ts, 0.00001));
1122c78c1075SJoseph Pusztay PetscCall(TSSetMaxSteps(ts, 100));
1123c78c1075SJoseph Pusztay PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1124c78c1075SJoseph Pusztay PetscCall(TSMonitorSet(ts, EnergyMonitor, &user, NULL));
1125c78c1075SJoseph Pusztay PetscCall(TSSetFromOptions(ts));
1126c78c1075SJoseph Pusztay PetscCall(TSSetComputeInitialCondition(ts, InitializeSolve));
1127c78c1075SJoseph Pusztay PetscCall(TSSetComputeExactError(ts, ComputeError));
1128ecc87898SStefano Zampini PetscCall(TSSetResize(ts, PETSC_FALSE, SetUpMigrateParticles, MigrateParticles, NULL));
1129c78c1075SJoseph Pusztay
1130c78c1075SJoseph Pusztay PetscCall(CreateSolution(ts));
1131c78c1075SJoseph Pusztay PetscCall(TSGetSolution(ts, &u));
1132c78c1075SJoseph Pusztay PetscCall(TSComputeInitialCondition(ts, u));
1133c78c1075SJoseph Pusztay PetscCall(TSSolve(ts, NULL));
1134c78c1075SJoseph Pusztay
1135c78c1075SJoseph Pusztay PetscCall(SNESDestroy(&user.snes));
1136c78c1075SJoseph Pusztay PetscCall(TSDestroy(&ts));
1137c78c1075SJoseph Pusztay PetscCall(DMDestroy(&sw));
1138c78c1075SJoseph Pusztay PetscCall(DMDestroy(&dm));
1139c78c1075SJoseph Pusztay PetscCall(PetscFinalize());
1140c78c1075SJoseph Pusztay return 0;
1141c78c1075SJoseph Pusztay }
1142c78c1075SJoseph Pusztay
1143c78c1075SJoseph Pusztay /*TEST
1144c78c1075SJoseph Pusztay
1145c78c1075SJoseph Pusztay build:
1146c78c1075SJoseph Pusztay requires: double !complex
1147c78c1075SJoseph Pusztay
1148c78c1075SJoseph Pusztay testset:
1149c78c1075SJoseph Pusztay requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1150c78c1075SJoseph 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 \
1151c78c1075SJoseph Pusztay -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1152c78c1075SJoseph Pusztay -ts_type basicsymplectic\
1153188af4bfSBarry Smith -dm_view -output_step 50 -ts_time_step 0.01 -ts_max_time 10.0 -ts_max_steps 10
1154c78c1075SJoseph Pusztay test:
1155c78c1075SJoseph Pusztay suffix: none_bsi_2d_1
1156c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 1 -em_type none -error
1157c78c1075SJoseph Pusztay test:
1158c78c1075SJoseph Pusztay suffix: none_bsi_2d_2
1159c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 2 -em_type none -error
1160c78c1075SJoseph Pusztay test:
1161c78c1075SJoseph Pusztay suffix: none_bsi_2d_3
1162c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 3 -em_type none -error
1163c78c1075SJoseph Pusztay test:
1164c78c1075SJoseph Pusztay suffix: none_bsi_2d_4
1165c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 4 -em_type none -error
1166c78c1075SJoseph Pusztay test:
1167c78c1075SJoseph Pusztay suffix: coulomb_bsi_2d_1
1168c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 1
1169c78c1075SJoseph Pusztay test:
1170c78c1075SJoseph Pusztay suffix: coulomb_bsi_2d_2
1171c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 2
1172c78c1075SJoseph Pusztay test:
1173c78c1075SJoseph Pusztay suffix: coulomb_bsi_2d_3
1174c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 3
1175c78c1075SJoseph Pusztay test:
1176c78c1075SJoseph Pusztay suffix: coulomb_bsi_2d_4
1177c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 4
1178c78c1075SJoseph Pusztay
1179c78c1075SJoseph Pusztay testset:
1180c78c1075SJoseph Pusztay requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1181c78c1075SJoseph 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 \
1182c78c1075SJoseph Pusztay -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1183c78c1075SJoseph Pusztay -ts_type basicsymplectic\
1184c78c1075SJoseph Pusztay -em_type primal -em_pc_type svd\
1185188af4bfSBarry Smith -dm_view -output_step 50 -error -ts_time_step 0.01 -ts_max_time 10.0 -ts_max_steps 10\
1186c78c1075SJoseph Pusztay -petscspace_degree 2 -petscfe_default_quadrature_order 3 -sigma 1.0e-8 -timeScale 2.0e-14
1187c78c1075SJoseph Pusztay test:
1188c78c1075SJoseph Pusztay suffix: poisson_bsi_2d_1
1189c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 1
1190c78c1075SJoseph Pusztay test:
1191c78c1075SJoseph Pusztay suffix: poisson_bsi_2d_2
1192c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 2
1193c78c1075SJoseph Pusztay test:
1194c78c1075SJoseph Pusztay suffix: poisson_bsi_2d_3
1195c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 3
1196c78c1075SJoseph Pusztay test:
1197c78c1075SJoseph Pusztay suffix: poisson_bsi_2d_4
1198c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 4
1199c78c1075SJoseph Pusztay
1200c78c1075SJoseph Pusztay testset:
1201c78c1075SJoseph Pusztay requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1202f940b0e3Sdanofinn 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 \
1203f940b0e3Sdanofinn -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1204f940b0e3Sdanofinn -ts_convergence_estimate -convest_num_refine 2 -em_type primal \
1205f940b0e3Sdanofinn -mat_type baij -em_ksp_error_if_not_converged -em_pc_type svd \
1206188af4bfSBarry Smith -dm_view -output_step 50 -error -ts_time_step 0.01 -ts_max_time 10.0 -ts_max_steps 10 \
1207f940b0e3Sdanofinn -sigma 1.0e-8 -timeScale 2.0e-14
1208c78c1075SJoseph Pusztay test:
1209c78c1075SJoseph Pusztay suffix: im_2d_0
1210f940b0e3Sdanofinn args: -ts_type theta -ts_theta_theta 0.5
1211f940b0e3Sdanofinn test:
1212f940b0e3Sdanofinn suffix: dg_2d_none
1213f940b0e3Sdanofinn args: -ts_type discgrad -ts_discgrad_type none -snes_type qn
1214f940b0e3Sdanofinn test:
1215f940b0e3Sdanofinn suffix: dg_2d_average
1216f940b0e3Sdanofinn args: -ts_type discgrad -ts_discgrad_type average -snes_type qn
1217f940b0e3Sdanofinn test:
1218f940b0e3Sdanofinn suffix: dg_2d_gonzalez
1219f940b0e3Sdanofinn args: -ts_type discgrad -ts_discgrad_type gonzalez -snes_fd -snes_type newtonls -snes_fd -pc_type lu
1220c78c1075SJoseph Pusztay
1221c78c1075SJoseph Pusztay testset:
1222c78c1075SJoseph Pusztay requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1223c78c1075SJoseph 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 \
1224c78c1075SJoseph Pusztay -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV -dm_swarm_num_species 1\
1225c78c1075SJoseph Pusztay -ts_type basicsymplectic -ts_convergence_estimate -convest_num_refine 2 \
1226c78c1075SJoseph Pusztay -em_snes_type ksponly -em_pc_type svd -em_type primal -petscspace_degree 1\
1227c78c1075SJoseph Pusztay -dm_view -output_step 50\
1228188af4bfSBarry Smith -pc_type svd -sigma 1.0e-8 -timeScale 2.0e-14 -ts_time_step 0.01 -ts_max_time 1.0 -ts_max_steps 10
1229c78c1075SJoseph Pusztay test:
1230c78c1075SJoseph Pusztay suffix: bsi_2d_mesh_1
1231c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 4
1232c78c1075SJoseph Pusztay test:
1233c78c1075SJoseph Pusztay suffix: bsi_2d_mesh_1_par_2
1234c78c1075SJoseph Pusztay nsize: 2
1235c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 4
1236c78c1075SJoseph Pusztay test:
1237c78c1075SJoseph Pusztay suffix: bsi_2d_mesh_1_par_3
1238c78c1075SJoseph Pusztay nsize: 3
1239c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 4
1240c78c1075SJoseph Pusztay test:
1241c78c1075SJoseph Pusztay suffix: bsi_2d_mesh_1_par_4
1242c78c1075SJoseph Pusztay nsize: 4
1243c78c1075SJoseph Pusztay args: -ts_basicsymplectic_type 4 -dm_swarm_num_particles 0,0,2,0
1244c78c1075SJoseph Pusztay
1245c78c1075SJoseph Pusztay testset:
1246c78c1075SJoseph Pusztay requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1247c78c1075SJoseph 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 \
1248c78c1075SJoseph Pusztay -dm_swarm_num_particles 10 -dm_swarm_coordinate_function circleMultipleX -dm_swarm_velocity_function circleMultipleV \
1249c78c1075SJoseph Pusztay -ts_convergence_estimate -convest_num_refine 2 \
1250c78c1075SJoseph Pusztay -em_pc_type lu\
1251c78c1075SJoseph Pusztay -dm_view -output_step 50 -error\
1252188af4bfSBarry Smith -sigma 1.0e-8 -timeScale 2.0e-14 -ts_time_step 0.01 -ts_max_time 10.0 -ts_max_steps 10
1253c78c1075SJoseph Pusztay test:
1254c78c1075SJoseph Pusztay suffix: bsi_2d_multiple_1
1255c78c1075SJoseph Pusztay args: -ts_type basicsymplectic -ts_basicsymplectic_type 1
1256c78c1075SJoseph Pusztay test:
1257c78c1075SJoseph Pusztay suffix: bsi_2d_multiple_2
1258c78c1075SJoseph Pusztay args: -ts_type basicsymplectic -ts_basicsymplectic_type 2
1259c78c1075SJoseph Pusztay test:
1260c78c1075SJoseph Pusztay suffix: bsi_2d_multiple_3
1261188af4bfSBarry Smith args: -ts_type basicsymplectic -ts_basicsymplectic_type 3 -ts_time_step 0.001
1262c78c1075SJoseph Pusztay test:
1263c78c1075SJoseph Pusztay suffix: im_2d_multiple_0
1264c78c1075SJoseph Pusztay args: -ts_type theta -ts_theta_theta 0.5 \
1265f940b0e3Sdanofinn -mat_type baij -em_ksp_error_if_not_converged -em_pc_type lu
1266c78c1075SJoseph Pusztay
1267c78c1075SJoseph Pusztay testset:
1268c78c1075SJoseph Pusztay requires: defined(PETSC_HAVE_EXECUTABLE_EXPORT)
1269c78c1075SJoseph 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 \
1270c78c1075SJoseph Pusztay -dm_swarm_num_particles 2 -dm_swarm_coordinate_function circleSingleX -dm_swarm_velocity_function circleSingleV \
1271c78c1075SJoseph Pusztay -em_pc_type fieldsplit -ksp_rtol 1e-10 -em_ksp_type preonly -em_type mixed -em_ksp_error_if_not_converged\
1272c78c1075SJoseph Pusztay -dm_view -output_step 50 -error -dm_refine 0\
1273188af4bfSBarry Smith -pc_type svd -sigma 1.0e-8 -timeScale 2.0e-14 -ts_time_step 0.01 -ts_max_time 10.0 -ts_max_steps 10
1274c78c1075SJoseph Pusztay test:
1275c78c1075SJoseph Pusztay suffix: bsi_4_rt_1
1276c78c1075SJoseph Pusztay args: -ts_type basicsymplectic -ts_basicsymplectic_type 4\
1277c78c1075SJoseph Pusztay -pc_fieldsplit_detect_saddle_point\
1278c78c1075SJoseph Pusztay -pc_type fieldsplit\
1279c78c1075SJoseph Pusztay -pc_fieldsplit_type schur\
1280c78c1075SJoseph Pusztay -pc_fieldsplit_schur_precondition full \
1281c78c1075SJoseph Pusztay -field_petscspace_degree 2\
1282c78c1075SJoseph Pusztay -field_petscfe_default_quadrature_order 1\
1283c78c1075SJoseph Pusztay -field_petscspace_type sum \
1284c78c1075SJoseph Pusztay -field_petscspace_variables 2 \
1285c78c1075SJoseph Pusztay -field_petscspace_components 2 \
1286c78c1075SJoseph Pusztay -field_petscspace_sum_spaces 2 \
1287c78c1075SJoseph Pusztay -field_petscspace_sum_concatenate true \
1288c78c1075SJoseph Pusztay -field_sumcomp_0_petscspace_variables 2 \
1289c78c1075SJoseph Pusztay -field_sumcomp_0_petscspace_type tensor \
1290c78c1075SJoseph Pusztay -field_sumcomp_0_petscspace_tensor_spaces 2 \
1291c78c1075SJoseph Pusztay -field_sumcomp_0_petscspace_tensor_uniform false \
1292c78c1075SJoseph Pusztay -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
1293c78c1075SJoseph Pusztay -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
1294c78c1075SJoseph Pusztay -field_sumcomp_1_petscspace_variables 2 \
1295c78c1075SJoseph Pusztay -field_sumcomp_1_petscspace_type tensor \
1296c78c1075SJoseph Pusztay -field_sumcomp_1_petscspace_tensor_spaces 2 \
1297c78c1075SJoseph Pusztay -field_sumcomp_1_petscspace_tensor_uniform false \
1298c78c1075SJoseph Pusztay -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
1299c78c1075SJoseph Pusztay -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
1300c78c1075SJoseph Pusztay -field_petscdualspace_form_degree -1 \
1301c78c1075SJoseph Pusztay -field_petscdualspace_order 1 \
1302c78c1075SJoseph Pusztay -field_petscdualspace_lagrange_trimmed true\
1303c78c1075SJoseph Pusztay -ksp_gmres_restart 500
1304c78c1075SJoseph Pusztay
1305c78c1075SJoseph Pusztay TEST*/
1306