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