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