xref: /petsc/src/snes/tests/ex15.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1 static char help[] = "Mixed element discretization of the Poisson equation.\n\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscdmswarm.h>
5 #include <petscds.h>
6 #include <petscsnes.h>
7 #include <petscconvest.h>
8 #include <petscbag.h>
9 
10 /*
11 The Poisson equation
12 
13   -\Delta\phi = f
14 
15 can be rewritten in first order form
16 
17   q - \nabla\phi  &= 0
18   -\nabla \cdot q &= f
19 */
20 
21 typedef enum {
22   SIGMA,
23   NUM_CONSTANTS
24 } ConstantType;
25 typedef struct {
26   PetscReal sigma; /* Nondimensional charge per length in x */
27 } Parameter;
28 
29 typedef enum {
30   SOL_CONST,
31   SOL_LINEAR,
32   SOL_QUADRATIC,
33   SOL_TRIG,
34   SOL_TRIGX,
35   SOL_PARTICLES,
36   NUM_SOL_TYPES
37 } SolType;
38 static const char *solTypes[] = {"const", "linear", "quadratic", "trig", "trigx", "particles"};
39 
40 typedef struct {
41   SolType   solType; /* MMS solution type */
42   PetscBag  bag;     /* Problem parameters */
43   PetscBool particleRHS;
44   PetscInt  Np;
45 } AppCtx;
46 
47 /* SOLUTION CONST: \phi = 1, q = 0, f = 0 */
const_phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)48 static PetscErrorCode const_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
49 {
50   *u = 1.0;
51   return PETSC_SUCCESS;
52 }
53 
const_q(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)54 static PetscErrorCode const_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
55 {
56   for (PetscInt d = 0; d < dim; ++d) u[d] = 0.0;
57   return PETSC_SUCCESS;
58 }
59 
60 /* SOLUTION LINEAR: \phi = 2y, q = <0, 2>, f = 0 */
linear_phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)61 static PetscErrorCode linear_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
62 {
63   u[0] = 2. * x[1];
64   return PETSC_SUCCESS;
65 }
66 
linear_q(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)67 static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
68 {
69   u[0] = 0.;
70   u[1] = 2.;
71   return PETSC_SUCCESS;
72 }
73 
74 /* SOLUTION QUADRATIC: \phi = x (2\pi - x) + (1 + y) (1 - y), q = <2\pi - 2 x, - 2 y> = <2\pi, 0> - 2 x, f = -4 */
quadratic_phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)75 static PetscErrorCode quadratic_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
76 {
77   u[0] = x[0] * (6.283185307179586 - x[0]) + (1. + x[1]) * (1. - x[1]);
78   return PETSC_SUCCESS;
79 }
80 
quadratic_q(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)81 static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
82 {
83   u[0] = 6.283185307179586 - 2. * x[0];
84   u[1] = -2. * x[1];
85   return PETSC_SUCCESS;
86 }
87 
quadratic_q_bc(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)88 static PetscErrorCode quadratic_q_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
89 {
90   u[0] = x[1] > 0. ? -2. * x[1] : 2. * x[1];
91   return PETSC_SUCCESS;
92 }
93 
f0_quadratic_phi(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[])94 static void f0_quadratic_phi(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[])
95 {
96   for (PetscInt d = 0; d < dim; ++d) f0[0] -= -2.0;
97 }
98 
99 /* SOLUTION TRIG: \phi = sin(x) + (1/3 - y^2), q = <cos(x), -2 y>, f = sin(x) + 2 */
trig_phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)100 static PetscErrorCode trig_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
101 {
102   u[0] = PetscSinReal(x[0]) + (1. / 3. - x[1] * x[1]);
103   return PETSC_SUCCESS;
104 }
105 
trig_q(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)106 static PetscErrorCode trig_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
107 {
108   u[0] = PetscCosReal(x[0]);
109   u[1] = -2. * x[1];
110   return PETSC_SUCCESS;
111 }
112 
trig_q_bc(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)113 static PetscErrorCode trig_q_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
114 {
115   u[0] = x[1] > 0. ? -2. * x[1] : 2. * x[1];
116   return PETSC_SUCCESS;
117 }
118 
f0_trig_phi(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[])119 static void f0_trig_phi(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[])
120 {
121   f0[0] += PetscSinReal(x[0]) + 2.;
122 }
123 
124 /* SOLUTION TRIGX: \phi = sin(x), q = <cos(x), 0>, f = sin(x) */
trigx_phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)125 static PetscErrorCode trigx_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
126 {
127   u[0] = PetscSinReal(x[0]);
128   return PETSC_SUCCESS;
129 }
130 
trigx_q(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)131 static PetscErrorCode trigx_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
132 {
133   u[0] = PetscCosReal(x[0]);
134   u[1] = 0.;
135   return PETSC_SUCCESS;
136 }
137 
trigx_q_bc(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)138 static PetscErrorCode trigx_q_bc(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
139 {
140   u[0] = x[1] > 0. ? -2. * x[1] : 2. * x[1];
141   return PETSC_SUCCESS;
142 }
143 
f0_trigx_phi(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[])144 static void f0_trigx_phi(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[])
145 {
146   f0[0] += PetscSinReal(x[0]);
147 }
148 
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[])149 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[])
150 {
151   for (PetscInt d = 0; d < dim; ++d) f0[d] += u[uOff[0] + d];
152 }
153 
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[])154 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[])
155 {
156   for (PetscInt d = 0; d < dim; ++d) f1[d * dim + d] = u[uOff[1]];
157 }
158 
f0_phi(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[])159 static void f0_phi(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[])
160 {
161   for (PetscInt d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0] + d * dim + d];
162 }
163 
f0_phi_backgroundCharge(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[])164 static void f0_phi_backgroundCharge(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[])
165 {
166   f0[0] += constants[SIGMA];
167   for (PetscInt d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0] + d * dim + d];
168 }
169 
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[])170 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[])
171 {
172   for (PetscInt d = 0; d < dim; ++d) g0[d * dim + d] = 1.0;
173 }
174 
g2_qphi(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[])175 static void g2_qphi(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[])
176 {
177   for (PetscInt d = 0; d < dim; ++d) g2[d * dim + d] = 1.0;
178 }
179 
g1_phiq(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[])180 static void g1_phiq(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[])
181 {
182   for (PetscInt d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
183 }
184 
185 /* SOLUTION PARTICLES: \phi = sigma, q = <cos(x), 0>, f = sin(x) */
particles_phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)186 static PetscErrorCode particles_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
187 {
188   u[0] = 0.0795775;
189   return PETSC_SUCCESS;
190 }
191 
particles_q(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)192 static PetscErrorCode particles_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
193 {
194   u[0] = 0.;
195   u[1] = 0.;
196   return PETSC_SUCCESS;
197 }
198 
ProcessOptions(MPI_Comm comm,AppCtx * options)199 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
200 {
201   PetscInt sol;
202 
203   PetscFunctionBeginUser;
204   options->solType     = SOL_CONST;
205   options->particleRHS = PETSC_FALSE;
206   options->Np          = 100;
207 
208   PetscOptionsBegin(comm, "", "Mixed Poisson Options", "DMPLEX");
209   PetscCall(PetscOptionsBool("-particleRHS", "Flag to user particle RHS and background charge", "ex9.c", options->particleRHS, &options->particleRHS, NULL));
210   sol = options->solType;
211   PetscCall(PetscOptionsEList("-sol_type", "The MMS solution type", "ex12.c", solTypes, NUM_SOL_TYPES, solTypes[sol], &sol, NULL));
212   options->solType = (SolType)sol;
213   PetscOptionsEnd();
214   PetscFunctionReturn(PETSC_SUCCESS);
215 }
216 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)217 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
218 {
219   PetscFunctionBeginUser;
220   PetscCall(DMCreate(comm, dm));
221   PetscCall(DMSetType(*dm, DMPLEX));
222   PetscCall(DMSetFromOptions(*dm));
223   PetscCall(DMSetApplicationContext(*dm, user));
224   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
225   PetscFunctionReturn(PETSC_SUCCESS);
226 }
227 
SetupPrimalProblem(DM dm,AppCtx * user)228 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
229 {
230   PetscDS        ds;
231   PetscWeakForm  wf;
232   DMLabel        label;
233   const PetscInt id = 1;
234 
235   PetscFunctionBeginUser;
236   PetscCall(DMGetDS(dm, &ds));
237   PetscCall(PetscDSGetWeakForm(ds, &wf));
238   PetscCall(DMGetLabel(dm, "marker", &label));
239   PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q));
240   if (user->particleRHS) {
241     PetscCall(PetscDSSetResidual(ds, 1, f0_phi_backgroundCharge, NULL));
242   } else {
243     PetscCall(PetscDSSetResidual(ds, 1, f0_phi, NULL));
244   }
245   PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL));
246   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qphi, NULL));
247   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_phiq, NULL, NULL));
248   switch (user->solType) {
249   case SOL_CONST:
250     PetscCall(PetscDSSetExactSolution(ds, 0, const_q, user));
251     PetscCall(PetscDSSetExactSolution(ds, 1, const_phi, user));
252     break;
253   case SOL_LINEAR:
254     PetscCall(PetscDSSetExactSolution(ds, 0, linear_q, user));
255     PetscCall(PetscDSSetExactSolution(ds, 1, linear_phi, user));
256     break;
257   case SOL_QUADRATIC:
258     PetscCall(PetscWeakFormAddResidual(wf, NULL, 0, 1, 0, f0_quadratic_phi, NULL));
259     PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user));
260     PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_phi, user));
261     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)quadratic_q_bc, NULL, user, NULL));
262     break;
263   case SOL_TRIG:
264     PetscCall(PetscWeakFormAddResidual(wf, NULL, 0, 1, 0, f0_trig_phi, NULL));
265     PetscCall(PetscDSSetExactSolution(ds, 0, trig_q, user));
266     PetscCall(PetscDSSetExactSolution(ds, 1, trig_phi, user));
267     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)trig_q_bc, NULL, user, NULL));
268     break;
269   case SOL_TRIGX:
270     PetscCall(PetscWeakFormAddResidual(wf, NULL, 0, 1, 0, f0_trigx_phi, NULL));
271     PetscCall(PetscDSSetExactSolution(ds, 0, trigx_q, user));
272     PetscCall(PetscDSSetExactSolution(ds, 1, trigx_phi, user));
273     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)trigx_q_bc, NULL, user, NULL));
274     break;
275   case SOL_PARTICLES:
276     PetscCall(PetscDSSetExactSolution(ds, 0, particles_q, user));
277     PetscCall(PetscDSSetExactSolution(ds, 1, particles_phi, user));
278     break;
279   default:
280     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid solution type: %d", user->solType);
281   }
282   PetscFunctionReturn(PETSC_SUCCESS);
283 }
284 
SetupDiscretization(DM dm,PetscInt Nf,const char * names[],PetscErrorCode (* setup)(DM,AppCtx *),AppCtx * user)285 static PetscErrorCode SetupDiscretization(DM dm, PetscInt Nf, const char *names[], PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
286 {
287   DM             cdm = dm;
288   PetscFE        fe;
289   DMPolytopeType ct;
290   PetscInt       dim, cStart;
291   char           prefix[PETSC_MAX_PATH_LEN];
292 
293   PetscFunctionBeginUser;
294   PetscCall(DMGetDimension(dm, &dim));
295   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
296   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
297   for (PetscInt f = 0; f < Nf; ++f) {
298     PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", names[f]));
299     PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, prefix, -1, &fe));
300     PetscCall(PetscObjectSetName((PetscObject)fe, names[f]));
301     if (f > 0) {
302       PetscFE fe0;
303 
304       PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe0));
305       PetscCall(PetscFECopyQuadrature(fe0, fe));
306     }
307     PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
308     PetscCall(PetscFEDestroy(&fe));
309   }
310   PetscCall(DMCreateDS(dm));
311   PetscCall((*setup)(dm, user));
312   while (cdm) {
313     PetscCall(DMCopyDisc(dm, cdm));
314     PetscCall(DMGetCoarseDM(cdm, &cdm));
315   }
316   PetscFunctionReturn(PETSC_SUCCESS);
317 }
318 
InitializeWeights(DM sw,AppCtx * user)319 static PetscErrorCode InitializeWeights(DM sw, AppCtx *user)
320 {
321   PetscScalar *weight;
322   PetscInt     Np;
323   PetscReal    weightsum = 0.0;
324 
325   PetscFunctionBegin;
326   PetscCall(DMSwarmGetLocalSize(sw, &Np));
327   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
328   PetscCall(DMSwarmSortGetAccess(sw));
329   for (PetscInt p = 0; p < Np; ++p) {
330     weight[p] = 1.0 / Np;
331     weightsum += PetscRealPart(weight[p]);
332   }
333   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Np = %" PetscInt_FMT "\n", Np));
334   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "weightsum = %1.10f\n", (double)weightsum));
335   PetscCall(DMSwarmSortRestoreAccess(sw));
336   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
337   PetscFunctionReturn(PETSC_SUCCESS);
338 }
339 
CreateSwarm(DM dm,AppCtx * user,DM * sw)340 static PetscErrorCode CreateSwarm(DM dm, AppCtx *user, DM *sw)
341 {
342   PetscInt dim;
343 
344   PetscFunctionBeginUser;
345   PetscCall(DMGetDimension(dm, &dim));
346   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), sw));
347   PetscCall(DMSetType(*sw, DMSWARM));
348   PetscCall(DMSetDimension(*sw, dim));
349   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
350   PetscCall(DMSwarmSetCellDM(*sw, dm));
351   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", 1, PETSC_SCALAR));
352   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "species", 1, PETSC_INT));
353   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
354   PetscCall(DMSwarmComputeLocalSizeFromOptions(*sw));
355   PetscCall(DMSwarmInitializeCoordinates(*sw));
356   PetscCall(InitializeWeights(*sw, user));
357   PetscCall(DMSetFromOptions(*sw));
358   PetscCall(DMSetApplicationContext(*sw, user));
359   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particles"));
360   PetscCall(DMViewFromOptions(*sw, NULL, "-sw_view"));
361   PetscCall(DMSwarmVectorDefineField(*sw, "w_q"));
362   PetscFunctionReturn(PETSC_SUCCESS);
363 }
364 
SetupParameters(MPI_Comm comm,AppCtx * ctx)365 static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
366 {
367   PetscBag   bag;
368   Parameter *p;
369 
370   PetscFunctionBeginUser;
371   /* setup PETSc parameter bag */
372   PetscCall(PetscBagGetData(ctx->bag, &p));
373   PetscCall(PetscBagSetName(ctx->bag, "par", "Parameters"));
374   bag = ctx->bag;
375   PetscCall(PetscBagRegisterScalar(bag, &p->sigma, 1.0, "sigma", "Charge per unit area, C/m^3"));
376   PetscCall(PetscBagSetFromOptions(bag));
377   {
378     PetscViewer       viewer;
379     PetscViewerFormat format;
380     PetscBool         flg;
381 
382     PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
383     if (flg) {
384       PetscCall(PetscViewerPushFormat(viewer, format));
385       PetscCall(PetscBagView(bag, viewer));
386       PetscCall(PetscViewerFlush(viewer));
387       PetscCall(PetscViewerPopFormat(viewer));
388       PetscCall(PetscViewerDestroy(&viewer));
389     }
390   }
391   PetscFunctionReturn(PETSC_SUCCESS);
392 }
393 
InitializeConstants(DM sw,AppCtx * user)394 static PetscErrorCode InitializeConstants(DM sw, AppCtx *user)
395 {
396   DM         dm;
397   PetscReal *weight, totalCharge, totalWeight = 0., gmin[3], gmax[3];
398   PetscInt   Np, p, dim;
399 
400   PetscFunctionBegin;
401   PetscCall(DMSwarmGetCellDM(sw, &dm));
402   PetscCall(DMGetDimension(sw, &dim));
403   PetscCall(DMSwarmGetLocalSize(sw, &Np));
404   PetscCall(DMGetBoundingBox(dm, gmin, gmax));
405   PetscCall(DMSwarmGetField(sw, "w_q", NULL, NULL, (void **)&weight));
406   for (p = 0; p < Np; ++p) totalWeight += weight[p];
407   totalCharge = -1.0 * totalWeight;
408   PetscCall(DMSwarmRestoreField(sw, "w_q", NULL, NULL, (void **)&weight));
409   {
410     Parameter *param;
411     PetscReal  Area;
412 
413     PetscCall(PetscBagGetData(user->bag, &param));
414     switch (dim) {
415     case 1:
416       Area = (gmax[0] - gmin[0]);
417       break;
418     case 2:
419       Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]);
420       break;
421     case 3:
422       Area = (gmax[0] - gmin[0]) * (gmax[1] - gmin[1]) * (gmax[2] - gmin[2]);
423       break;
424     default:
425       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Dimension %" PetscInt_FMT " not supported", dim);
426     }
427     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "dim = %" PetscInt_FMT "\ttotalWeight = %f\ttotalCharge = %f, Total Area = %f\n", dim, (double)totalWeight, (double)totalCharge, (double)Area));
428     param->sigma = PetscAbsReal(totalCharge / (Area));
429 
430     PetscCall(PetscPrintf(PETSC_COMM_SELF, "sigma: %g\n", (double)param->sigma));
431   }
432   /* Setup Constants */
433   {
434     PetscDS    ds;
435     Parameter *param;
436     PetscCall(PetscBagGetData(user->bag, &param));
437     PetscScalar constants[NUM_CONSTANTS];
438     constants[SIGMA] = param->sigma;
439     PetscCall(DMGetDS(dm, &ds));
440     PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants));
441   }
442   PetscFunctionReturn(PETSC_SUCCESS);
443 }
444 
main(int argc,char ** argv)445 int main(int argc, char **argv)
446 {
447   DM          dm, sw;
448   SNES        snes;
449   Vec         u;
450   AppCtx      user;
451   const char *names[] = {"q", "phi"};
452 
453   PetscFunctionBeginUser;
454   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
455   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
456   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
457   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
458   PetscCall(SNESSetDM(snes, dm));
459   PetscCall(SetupDiscretization(dm, 2, names, SetupPrimalProblem, &user));
460   if (user.particleRHS) {
461     PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag));
462     PetscCall(CreateSwarm(dm, &user, &sw));
463     PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
464     PetscCall(InitializeConstants(sw, &user));
465   }
466   PetscCall(DMCreateGlobalVector(dm, &u));
467   PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
468   PetscCall(SNESSetFromOptions(snes));
469   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
470   PetscCall(DMSNESCheckFromOptions(snes, u));
471   if (user.particleRHS) {
472     DM       potential_dm;
473     IS       potential_IS;
474     Mat      M_p;
475     Vec      rho, f, temp_rho;
476     PetscInt fields = 1;
477 
478     PetscCall(DMGetGlobalVector(dm, &rho));
479     PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
480     PetscCall(DMCreateSubDM(dm, 1, &fields, &potential_IS, &potential_dm));
481     PetscCall(DMCreateMassMatrix(sw, potential_dm, &M_p));
482     PetscCall(MatViewFromOptions(M_p, NULL, "-mp_view"));
483     PetscCall(DMGetGlobalVector(potential_dm, &temp_rho));
484     PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
485     PetscCall(PetscObjectSetName((PetscObject)f, "particle weight"));
486     PetscCall(VecViewFromOptions(f, NULL, "-weights_view"));
487     PetscCall(MatMultTranspose(M_p, f, temp_rho));
488     PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
489     PetscCall(MatDestroy(&M_p));
490     PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
491     PetscCall(VecViewFromOptions(rho, NULL, "-poisson_rho_view"));
492     PetscCall(VecISCopy(rho, potential_IS, SCATTER_FORWARD, temp_rho));
493     PetscCall(VecViewFromOptions(temp_rho, NULL, "-rho_view"));
494     PetscCall(DMRestoreGlobalVector(potential_dm, &temp_rho));
495     PetscCall(DMDestroy(&potential_dm));
496     PetscCall(ISDestroy(&potential_IS));
497 
498     PetscCall(SNESSolve(snes, rho, u));
499     PetscCall(DMRestoreGlobalVector(dm, &rho));
500   } else {
501     PetscCall(SNESSolve(snes, NULL, u));
502   }
503   PetscCall(VecDestroy(&u));
504   PetscCall(SNESDestroy(&snes));
505   PetscCall(DMDestroy(&dm));
506   if (user.particleRHS) {
507     PetscCall(DMDestroy(&sw));
508     PetscCall(PetscBagDestroy(&user.bag));
509   }
510   PetscCall(PetscFinalize());
511   return PETSC_SUCCESS;
512 }
513 
514 /*TEST
515 
516   # RT1-P0 on quads
517   testset:
518     args: -dm_plex_simplex 0 -dm_plex_box_bd periodic,none -dm_plex_box_faces 3,1 \
519           -dm_plex_box_lower 0,-1 -dm_plex_box_upper 6.283185307179586,1\
520           -phi_petscspace_degree 0 \
521           -phi_petscdualspace_lagrange_use_moments \
522           -phi_petscdualspace_lagrange_moment_order 2 \
523           -q_petscfe_default_quadrature_order 1 \
524           -q_petscspace_type sum \
525           -q_petscspace_variables 2 \
526           -q_petscspace_components 2 \
527           -q_petscspace_sum_spaces 2 \
528           -q_petscspace_sum_concatenate true \
529           -q_sumcomp_0_petscspace_variables 2 \
530           -q_sumcomp_0_petscspace_type tensor \
531           -q_sumcomp_0_petscspace_tensor_spaces 2 \
532           -q_sumcomp_0_petscspace_tensor_uniform false \
533           -q_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
534           -q_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
535           -q_sumcomp_1_petscspace_variables 2 \
536           -q_sumcomp_1_petscspace_type tensor \
537           -q_sumcomp_1_petscspace_tensor_spaces 2 \
538           -q_sumcomp_1_petscspace_tensor_uniform false \
539           -q_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
540           -q_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
541           -q_petscdualspace_form_degree -1 \
542           -q_petscdualspace_order 1 \
543           -q_petscdualspace_lagrange_trimmed true \
544           -ksp_error_if_not_converged \
545           -pc_type fieldsplit -pc_fieldsplit_type schur \
546           -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full
547 
548     # The Jacobian test is meaningless here
549     test:
550           suffix: quad_hdiv_0
551           args: -dmsnes_check
552           filter: sed -e "s/Taylor approximation converging at order.*''//"
553 
554     # The Jacobian test is meaningless here
555     test:
556           suffix: quad_hdiv_1
557           args: -sol_type linear -dmsnes_check
558           filter: sed -e "s/Taylor approximation converging at order.*''//"
559 
560     test:
561           suffix: quad_hdiv_2
562           args: -sol_type quadratic -dmsnes_check \
563                 -fieldsplit_q_pc_type lu -fieldsplit_phi_pc_type svd
564 
565     test:
566           suffix: quad_hdiv_3
567           args: -sol_type trig \
568                 -fieldsplit_q_pc_type lu -fieldsplit_phi_pc_type svd
569           output_file: output/empty.out
570 
571     test:
572           suffix: quad_hdiv_4
573           requires: !single
574           args: -sol_type trigx \
575                 -fieldsplit_q_pc_type lu -fieldsplit_phi_pc_type svd
576           output_file: output/empty.out
577 
578     test:
579           suffix: particle_hdiv_5
580           requires: !complex double
581           args: -dm_swarm_num_particles 100 -dm_swarm_coordinate_density constant \
582                 -particleRHS -sol_type particles \
583                 -fieldsplit_q_pc_type lu -fieldsplit_phi_pc_type svd
584 
585 TEST*/
586