xref: /petsc/src/ts/tutorials/ex30.c (revision 204aa5235b0b5da05dfdf782f0d1b3760e61048e)
1811d88abSStefano Zampini static char help[] = "Biological network from https://link.springer.com/article/10.1007/s42967-023-00297-3\n\n\n";
2811d88abSStefano Zampini 
3811d88abSStefano Zampini #include <petscts.h>
4811d88abSStefano Zampini #include <petscsf.h>
5811d88abSStefano Zampini #include <petscdmplex.h>
6811d88abSStefano Zampini #include <petscdmplextransform.h>
7811d88abSStefano Zampini #include <petscdmforest.h>
8811d88abSStefano Zampini #include <petscviewerhdf5.h>
9811d88abSStefano Zampini #include <petscds.h>
10811d88abSStefano Zampini 
11811d88abSStefano Zampini /*
12811d88abSStefano Zampini     Here we solve the system of PDEs on \Omega \in R^2:
13811d88abSStefano Zampini 
14811d88abSStefano Zampini     * dC/dt - D^2 \Delta C - c^2 \nabla p \cross \nabla p + \alpha sqrt(||C||^2_F + eps)^(\gamma-2) C = 0
15811d88abSStefano Zampini     * - \nabla \cdot ((r + C) \nabla p) = S
16811d88abSStefano Zampini 
17811d88abSStefano Zampini     where:
18811d88abSStefano Zampini       C = symmetric 2x2 conductivity tensor
19811d88abSStefano Zampini       p = potential
20811d88abSStefano Zampini       S = source
21811d88abSStefano Zampini 
22811d88abSStefano Zampini     with natural boundary conditions on \partial\Omega:
23811d88abSStefano Zampini       \nabla C \cdot n  = 0
24811d88abSStefano Zampini       \nabla ((r + C)\nabla p) \cdot n  = 0
25811d88abSStefano Zampini 
26811d88abSStefano Zampini     Parameters:
27811d88abSStefano Zampini       D = diffusion constant
28811d88abSStefano Zampini       c = activation parameter
29811d88abSStefano Zampini       \alpha = metabolic coefficient
30811d88abSStefano Zampini       \gamma = metabolic exponent
31811d88abSStefano Zampini       r, eps are regularization parameters
32811d88abSStefano Zampini 
33811d88abSStefano Zampini     We use Lagrange elements for C_ij and P.
34811d88abSStefano Zampini */
35811d88abSStefano Zampini 
36811d88abSStefano Zampini typedef enum _fieldidx {
37811d88abSStefano Zampini   C_FIELD_ID = 0,
38811d88abSStefano Zampini   P_FIELD_ID,
39811d88abSStefano Zampini   NUM_FIELDS
40811d88abSStefano Zampini } FieldIdx;
41811d88abSStefano Zampini 
42811d88abSStefano Zampini typedef enum _constantidx {
43811d88abSStefano Zampini   R_ID = 0,
44811d88abSStefano Zampini   EPS_ID,
45811d88abSStefano Zampini   ALPHA_ID,
46811d88abSStefano Zampini   GAMMA_ID,
47811d88abSStefano Zampini   D_ID,
48811d88abSStefano Zampini   C2_ID,
49*204aa523SStefano Zampini   FIXC_ID,
50811d88abSStefano Zampini   NUM_CONSTANTS
51811d88abSStefano Zampini } ConstantIdx;
52811d88abSStefano Zampini 
53811d88abSStefano Zampini PetscLogStage SetupStage, SolveStage;
54811d88abSStefano Zampini 
55811d88abSStefano Zampini #define NORM2C(c00, c01, c11) PetscSqr(c00) + 2 * PetscSqr(c01) + PetscSqr(c11)
56811d88abSStefano Zampini 
57811d88abSStefano Zampini /* residual for C when tested against basis functions */
58811d88abSStefano Zampini static void C_0(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[])
59811d88abSStefano Zampini {
60811d88abSStefano Zampini   const PetscReal   c2       = PetscRealPart(constants[C2_ID]);
61811d88abSStefano Zampini   const PetscReal   alpha    = PetscRealPart(constants[ALPHA_ID]);
62811d88abSStefano Zampini   const PetscReal   gamma    = PetscRealPart(constants[GAMMA_ID]);
63811d88abSStefano Zampini   const PetscReal   eps      = PetscRealPart(constants[EPS_ID]);
64811d88abSStefano Zampini   const PetscScalar gradp[]  = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
65811d88abSStefano Zampini   const PetscScalar crossp[] = {gradp[0] * gradp[0], gradp[0] * gradp[1], gradp[1] * gradp[1]};
66811d88abSStefano Zampini   const PetscScalar C00      = u[uOff[C_FIELD_ID]];
67811d88abSStefano Zampini   const PetscScalar C01      = u[uOff[C_FIELD_ID] + 1];
68811d88abSStefano Zampini   const PetscScalar C11      = u[uOff[C_FIELD_ID] + 2];
69811d88abSStefano Zampini   const PetscScalar norm     = NORM2C(C00, C01, C11) + eps;
70811d88abSStefano Zampini   const PetscScalar nexp     = (gamma - 2.0) / 2.0;
71811d88abSStefano Zampini   const PetscScalar fnorm    = PetscPowScalar(norm, nexp);
72811d88abSStefano Zampini 
73811d88abSStefano Zampini   for (PetscInt k = 0; k < 3; k++) f0[k] = u_t[uOff[C_FIELD_ID] + k] - c2 * crossp[k] + alpha * fnorm * u[uOff[C_FIELD_ID] + k];
74811d88abSStefano Zampini }
75811d88abSStefano Zampini 
76811d88abSStefano Zampini /* Jacobian for C against C basis functions */
77811d88abSStefano Zampini static void JC_0_c0c0(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 J[])
78811d88abSStefano Zampini {
79811d88abSStefano Zampini   const PetscReal   alpha  = PetscRealPart(constants[ALPHA_ID]);
80811d88abSStefano Zampini   const PetscReal   gamma  = PetscRealPart(constants[GAMMA_ID]);
81811d88abSStefano Zampini   const PetscReal   eps    = PetscRealPart(constants[EPS_ID]);
82811d88abSStefano Zampini   const PetscScalar C00    = u[uOff[C_FIELD_ID]];
83811d88abSStefano Zampini   const PetscScalar C01    = u[uOff[C_FIELD_ID] + 1];
84811d88abSStefano Zampini   const PetscScalar C11    = u[uOff[C_FIELD_ID] + 2];
85811d88abSStefano Zampini   const PetscScalar norm   = NORM2C(C00, C01, C11) + eps;
86811d88abSStefano Zampini   const PetscScalar nexp   = (gamma - 2.0) / 2.0;
87811d88abSStefano Zampini   const PetscScalar fnorm  = PetscPowScalar(norm, nexp);
88811d88abSStefano Zampini   const PetscScalar dfnorm = nexp * PetscPowScalar(norm, nexp - 1.0);
89811d88abSStefano Zampini   const PetscScalar dC[]   = {2 * C00, 4 * C01, 2 * C11};
90811d88abSStefano Zampini 
91811d88abSStefano Zampini   for (PetscInt k = 0; k < 3; k++) {
92811d88abSStefano Zampini     for (PetscInt j = 0; j < 3; j++) J[k * 3 + j] = alpha * dfnorm * dC[j] * u[uOff[C_FIELD_ID] + k];
93811d88abSStefano Zampini     J[k * 3 + k] += alpha * fnorm + u_tShift;
94811d88abSStefano Zampini   }
95811d88abSStefano Zampini }
96811d88abSStefano Zampini 
97811d88abSStefano Zampini /* Jacobian for C against C basis functions and gradients of P basis functions */
98811d88abSStefano Zampini static void JC_0_c0p1(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 J[])
99811d88abSStefano Zampini {
100811d88abSStefano Zampini   const PetscReal   c2      = PetscRealPart(constants[C2_ID]);
101811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
102811d88abSStefano Zampini 
103811d88abSStefano Zampini   J[0] = -c2 * 2 * gradp[0];
104811d88abSStefano Zampini   J[1] = 0.0;
105811d88abSStefano Zampini   J[2] = -c2 * gradp[1];
106811d88abSStefano Zampini   J[3] = -c2 * gradp[0];
107811d88abSStefano Zampini   J[4] = 0.0;
108811d88abSStefano Zampini   J[5] = -c2 * 2 * gradp[1];
109811d88abSStefano Zampini }
110811d88abSStefano Zampini 
111811d88abSStefano Zampini /* residual for C when tested against gradients of basis functions */
112811d88abSStefano Zampini static void C_1(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[])
113811d88abSStefano Zampini {
114811d88abSStefano Zampini   const PetscReal D = PetscRealPart(constants[D_ID]);
115811d88abSStefano Zampini   for (PetscInt k = 0; k < 3; k++)
116811d88abSStefano Zampini     for (PetscInt d = 0; d < 2; d++) f1[k * 2 + d] = PetscSqr(D) * u_x[uOff_x[C_FIELD_ID] + k * 2 + d];
117811d88abSStefano Zampini }
118811d88abSStefano Zampini 
119811d88abSStefano Zampini /* Jacobian for C against gradients of C basis functions */
120811d88abSStefano Zampini static void JC_1_c1c1(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 J[])
121811d88abSStefano Zampini {
122811d88abSStefano Zampini   const PetscReal D = PetscRealPart(constants[D_ID]);
123811d88abSStefano Zampini   for (PetscInt k = 0; k < 3; k++)
124811d88abSStefano Zampini     for (PetscInt d = 0; d < 2; d++) J[k * (3 + 1) * 2 * 2 + d * 2 + d] = PetscSqr(D);
125811d88abSStefano Zampini }
126811d88abSStefano Zampini 
127811d88abSStefano Zampini /* residual for P when tested against basis functions.
128811d88abSStefano Zampini    The source term always comes from the auxiliary vec because it needs to have zero mean */
129811d88abSStefano Zampini static void P_0(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[])
130811d88abSStefano Zampini {
131811d88abSStefano Zampini   PetscScalar S = a[aOff[P_FIELD_ID]];
132811d88abSStefano Zampini 
133811d88abSStefano Zampini   f0[0] = -S;
134811d88abSStefano Zampini }
135811d88abSStefano Zampini 
136*204aa523SStefano Zampini static inline void QuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal x[2]);
137*204aa523SStefano Zampini 
138*204aa523SStefano Zampini /* compute shift to make C positive definite */
139*204aa523SStefano Zampini static inline PetscReal FIX_C(PetscScalar C00, PetscScalar C01, PetscScalar C11)
140*204aa523SStefano Zampini {
141*204aa523SStefano Zampini #if !PetscDefined(USE_COMPLEX)
142*204aa523SStefano Zampini   PetscReal eigs[2], s = 0.0;
143*204aa523SStefano Zampini 
144*204aa523SStefano Zampini   QuadraticRoots(1, -(C00 + C11), C00 * C11 - PetscSqr(C01), eigs);
145*204aa523SStefano Zampini   if (eigs[0] < 0 || eigs[1] < 0) s = -PetscMin(eigs[0], eigs[1]) + PETSC_SMALL;
146*204aa523SStefano Zampini   return s;
147*204aa523SStefano Zampini #else
148*204aa523SStefano Zampini   return 0.0;
149*204aa523SStefano Zampini #endif
150*204aa523SStefano Zampini }
151*204aa523SStefano Zampini 
152811d88abSStefano Zampini /* residual for P when tested against gradients of basis functions */
153811d88abSStefano Zampini static void P_1(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[])
154811d88abSStefano Zampini {
155811d88abSStefano Zampini   const PetscReal   r       = PetscRealPart(constants[R_ID]);
156*204aa523SStefano Zampini   const PetscScalar C00     = u[uOff[C_FIELD_ID]] + r;
157811d88abSStefano Zampini   const PetscScalar C01     = u[uOff[C_FIELD_ID] + 1];
158811d88abSStefano Zampini   const PetscScalar C10     = C01;
159*204aa523SStefano Zampini   const PetscScalar C11     = u[uOff[C_FIELD_ID] + 2] + r;
160811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
161*204aa523SStefano Zampini   const PetscBool   fix_c   = (PetscBool)(PetscRealPart(constants[FIXC_ID]) > 1.0);
162*204aa523SStefano Zampini   const PetscScalar s       = fix_c ? FIX_C(C00, C01, C11) : 0.0;
163811d88abSStefano Zampini 
164*204aa523SStefano Zampini   f1[0] = (C00 + s) * gradp[0] + C01 * gradp[1];
165*204aa523SStefano Zampini   f1[1] = C10 * gradp[0] + (C11 + s) * gradp[1];
166811d88abSStefano Zampini }
167811d88abSStefano Zampini 
168811d88abSStefano Zampini /* Same as above for the P-only subproblem for initial conditions: the conductivity values come from the auxiliary vec */
169811d88abSStefano Zampini static void P_1_aux(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[])
170811d88abSStefano Zampini {
171811d88abSStefano Zampini   const PetscReal   r       = PetscRealPart(constants[R_ID]);
172*204aa523SStefano Zampini   const PetscScalar C00     = a[aOff[C_FIELD_ID]] + r;
173811d88abSStefano Zampini   const PetscScalar C01     = a[aOff[C_FIELD_ID] + 1];
174811d88abSStefano Zampini   const PetscScalar C10     = C01;
175*204aa523SStefano Zampini   const PetscScalar C11     = a[aOff[C_FIELD_ID] + 2] + r;
176811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[0]], u_x[uOff_x[0] + 1]};
177*204aa523SStefano Zampini   const PetscBool   fix_c   = (PetscBool)(PetscRealPart(constants[FIXC_ID]) > 1.0);
178*204aa523SStefano Zampini   const PetscScalar s       = fix_c ? FIX_C(C00, C01, C11) : 0.0;
179811d88abSStefano Zampini 
180*204aa523SStefano Zampini   f1[0] = (C00 + s) * gradp[0] + C01 * gradp[1];
181*204aa523SStefano Zampini   f1[1] = C10 * gradp[0] + (C11 + s) * gradp[1];
182811d88abSStefano Zampini }
183811d88abSStefano Zampini 
184811d88abSStefano Zampini /* Jacobian for P against gradients of P basis functions */
185811d88abSStefano Zampini static void JP_1_p1p1(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 J[])
186811d88abSStefano Zampini {
187811d88abSStefano Zampini   const PetscReal   r     = PetscRealPart(constants[R_ID]);
188*204aa523SStefano Zampini   const PetscScalar C00   = u[uOff[C_FIELD_ID]] + r;
189811d88abSStefano Zampini   const PetscScalar C01   = u[uOff[C_FIELD_ID] + 1];
190811d88abSStefano Zampini   const PetscScalar C10   = C01;
191*204aa523SStefano Zampini   const PetscScalar C11   = u[uOff[C_FIELD_ID] + 2] + r;
192*204aa523SStefano Zampini   const PetscBool   fix_c = (PetscBool)(PetscRealPart(constants[FIXC_ID]) > 0.0);
193*204aa523SStefano Zampini   const PetscScalar s     = fix_c ? FIX_C(C00, C01, C11) : 0.0;
194811d88abSStefano Zampini 
195*204aa523SStefano Zampini   J[0] = C00 + s;
196811d88abSStefano Zampini   J[1] = C01;
197811d88abSStefano Zampini   J[2] = C10;
198*204aa523SStefano Zampini   J[3] = C11 + s;
199811d88abSStefano Zampini }
200811d88abSStefano Zampini 
201811d88abSStefano Zampini /* Same as above for the P-only subproblem for initial conditions */
202811d88abSStefano Zampini static void JP_1_p1p1_aux(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 J[])
203811d88abSStefano Zampini {
204811d88abSStefano Zampini   const PetscReal   r     = PetscRealPart(constants[R_ID]);
205*204aa523SStefano Zampini   const PetscScalar C00   = a[aOff[C_FIELD_ID]] + r;
206811d88abSStefano Zampini   const PetscScalar C01   = a[aOff[C_FIELD_ID] + 1];
207811d88abSStefano Zampini   const PetscScalar C10   = C01;
208*204aa523SStefano Zampini   const PetscScalar C11   = a[aOff[C_FIELD_ID] + 2] + r;
209*204aa523SStefano Zampini   const PetscBool   fix_c = (PetscBool)(PetscRealPart(constants[FIXC_ID]) > 0.0);
210*204aa523SStefano Zampini   const PetscScalar s     = fix_c ? FIX_C(C00, C01, C11) : 0.0;
211811d88abSStefano Zampini 
212*204aa523SStefano Zampini   J[0] = C00 + s;
213811d88abSStefano Zampini   J[1] = C01;
214811d88abSStefano Zampini   J[2] = C10;
215*204aa523SStefano Zampini   J[3] = C11 + s;
216811d88abSStefano Zampini }
217811d88abSStefano Zampini 
218811d88abSStefano Zampini /* Jacobian for P against gradients of P basis functions and C basis functions */
219811d88abSStefano Zampini static void JP_1_p1c0(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 J[])
220811d88abSStefano Zampini {
221811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
222811d88abSStefano Zampini 
223811d88abSStefano Zampini   J[0] = gradp[0];
224811d88abSStefano Zampini   J[1] = 0;
225811d88abSStefano Zampini   J[2] = gradp[1];
226811d88abSStefano Zampini   J[3] = gradp[0];
227811d88abSStefano Zampini   J[4] = 0;
228811d88abSStefano Zampini   J[5] = gradp[1];
229811d88abSStefano Zampini }
230811d88abSStefano Zampini 
231811d88abSStefano Zampini /* the source term S(x) = exp(-500*||x - x0||^2) */
232811d88abSStefano Zampini static PetscErrorCode source_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
233811d88abSStefano Zampini {
234811d88abSStefano Zampini   PetscReal *x0 = (PetscReal *)ctx;
235811d88abSStefano Zampini   PetscReal  n  = 0;
236811d88abSStefano Zampini 
237811d88abSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) n += (x[d] - x0[d]) * (x[d] - x0[d]);
238811d88abSStefano Zampini   u[0] = PetscExpReal(-500 * n);
239811d88abSStefano Zampini   return PETSC_SUCCESS;
240811d88abSStefano Zampini }
241811d88abSStefano Zampini 
242811d88abSStefano Zampini static PetscErrorCode source_1(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
243811d88abSStefano Zampini {
244811d88abSStefano Zampini   PetscScalar     ut[1];
245811d88abSStefano Zampini   const PetscReal x0[] = {0.25, 0.25};
246811d88abSStefano Zampini   const PetscReal x1[] = {0.75, 0.75};
247811d88abSStefano Zampini 
248811d88abSStefano Zampini   PetscCall(source_0(dim, time, x, Nf, ut, (void *)x0));
249811d88abSStefano Zampini   PetscCall(source_0(dim, time, x, Nf, u, (void *)x1));
250811d88abSStefano Zampini   u[0] += ut[0];
251811d88abSStefano Zampini   return PETSC_SUCCESS;
252811d88abSStefano Zampini }
253811d88abSStefano Zampini 
254811d88abSStefano Zampini /* functionals to be integrated: average -> \int_\Omega u dx */
255811d88abSStefano Zampini static void average(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 obj[])
256811d88abSStefano Zampini {
257811d88abSStefano Zampini   obj[0] = u[uOff[P_FIELD_ID]];
258811d88abSStefano Zampini }
259811d88abSStefano Zampini 
260811d88abSStefano Zampini /* stable implementation of roots of a*x^2 + b*x + c = 0 */
261811d88abSStefano Zampini static inline void QuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal x[2])
262811d88abSStefano Zampini {
263811d88abSStefano Zampini   PetscReal delta = PetscMax(b * b - 4 * a * c, 0); /* eigenvalues symmetric matrix */
264811d88abSStefano Zampini   PetscReal temp  = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(delta));
265811d88abSStefano Zampini 
266811d88abSStefano Zampini   x[0] = temp / a;
267811d88abSStefano Zampini   x[1] = c / temp;
268811d88abSStefano Zampini }
269811d88abSStefano Zampini 
270811d88abSStefano Zampini /* functionals to be integrated: energy -> D^2/2 * ||\nabla C||^2 + c^2\nabla p * (r + C) * \nabla p + \alpha/ \gamma * ||C||^\gamma */
271811d88abSStefano Zampini static void energy(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 obj[])
272811d88abSStefano Zampini {
273811d88abSStefano Zampini   const PetscReal   D         = PetscRealPart(constants[D_ID]);
274811d88abSStefano Zampini   const PetscReal   c2        = PetscRealPart(constants[C2_ID]);
275811d88abSStefano Zampini   const PetscReal   r         = PetscRealPart(constants[R_ID]);
276811d88abSStefano Zampini   const PetscReal   alpha     = PetscRealPart(constants[ALPHA_ID]);
277811d88abSStefano Zampini   const PetscReal   gamma     = PetscRealPart(constants[GAMMA_ID]);
278811d88abSStefano Zampini   const PetscScalar C00       = u[uOff[C_FIELD_ID]];
279811d88abSStefano Zampini   const PetscScalar C01       = u[uOff[C_FIELD_ID] + 1];
280811d88abSStefano Zampini   const PetscScalar C10       = C01;
281811d88abSStefano Zampini   const PetscScalar C11       = u[uOff[C_FIELD_ID] + 2];
282811d88abSStefano Zampini   const PetscScalar gradp[]   = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
283811d88abSStefano Zampini   const PetscScalar gradC00[] = {u_x[uOff_x[C_FIELD_ID] + 0], u_x[uOff_x[C_FIELD_ID] + 1]};
284811d88abSStefano Zampini   const PetscScalar gradC01[] = {u_x[uOff_x[C_FIELD_ID] + 2], u_x[uOff_x[C_FIELD_ID] + 3]};
285811d88abSStefano Zampini   const PetscScalar gradC11[] = {u_x[uOff_x[C_FIELD_ID] + 4], u_x[uOff_x[C_FIELD_ID] + 5]};
286811d88abSStefano Zampini   const PetscScalar normC     = NORM2C(C00, C01, C11);
287811d88abSStefano Zampini   const PetscScalar normgradC = NORM2C(gradC00[0], gradC01[0], gradC11[0]) + NORM2C(gradC00[1], gradC01[1], gradC11[1]);
288811d88abSStefano Zampini   const PetscScalar nexp      = gamma / 2.0;
289811d88abSStefano Zampini 
290811d88abSStefano Zampini   const PetscScalar t0 = PetscSqr(D) / 2.0 * normgradC;
291811d88abSStefano Zampini   const PetscScalar t1 = c2 * (gradp[0] * ((C00 + r) * gradp[0] + C01 * gradp[1]) + gradp[1] * (C10 * gradp[0] + (C11 + r) * gradp[1]));
292811d88abSStefano Zampini   const PetscScalar t2 = alpha / gamma * PetscPowScalar(normC, nexp);
293811d88abSStefano Zampini 
294811d88abSStefano Zampini   obj[0] = t0 + t1 + t2;
295811d88abSStefano Zampini }
296811d88abSStefano Zampini 
297811d88abSStefano Zampini /* functionals to be integrated: ellipticity_fail -> 0 means C+r is elliptic at quadrature point, otherwise it returns the absolute value of the most negative eigenvalue */
298811d88abSStefano Zampini static void ellipticity_fail(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 obj[])
299811d88abSStefano Zampini {
300811d88abSStefano Zampini   const PetscReal r   = PetscRealPart(constants[R_ID]);
301811d88abSStefano Zampini   const PetscReal C00 = PetscRealPart(u[uOff[C_FIELD_ID]] + r);
302811d88abSStefano Zampini   const PetscReal C01 = PetscRealPart(u[uOff[C_FIELD_ID] + 1]);
303811d88abSStefano Zampini   const PetscReal C11 = PetscRealPart(u[uOff[C_FIELD_ID] + 2] + r);
304811d88abSStefano Zampini 
305811d88abSStefano Zampini   PetscReal eigs[2];
306811d88abSStefano Zampini   QuadraticRoots(1, -(C00 + C11), C00 * C11 - PetscSqr(C01), eigs);
307811d88abSStefano Zampini   if (eigs[0] < 0 || eigs[1] < 0) obj[0] = -PetscMin(eigs[0], eigs[1]);
308811d88abSStefano Zampini   else obj[0] = 0.0;
309811d88abSStefano Zampini }
310811d88abSStefano Zampini 
311811d88abSStefano Zampini /* initial conditions for C: eq. 16 */
312811d88abSStefano Zampini static PetscErrorCode initial_conditions_C_0(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
313811d88abSStefano Zampini {
314811d88abSStefano Zampini   u[0] = 1;
315811d88abSStefano Zampini   u[1] = 0;
316811d88abSStefano Zampini   u[2] = 1;
317811d88abSStefano Zampini   return PETSC_SUCCESS;
318811d88abSStefano Zampini }
319811d88abSStefano Zampini 
320811d88abSStefano Zampini /* initial conditions for C: eq. 17 */
321811d88abSStefano Zampini static PetscErrorCode initial_conditions_C_1(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
322811d88abSStefano Zampini {
323811d88abSStefano Zampini   const PetscReal x = xx[0];
324811d88abSStefano Zampini   const PetscReal y = xx[1];
325811d88abSStefano Zampini 
326811d88abSStefano Zampini   u[0] = (2 - PetscAbsReal(x + y)) * PetscExpReal(-10 * PetscAbsReal(x - y));
327811d88abSStefano Zampini   u[1] = 0;
328811d88abSStefano Zampini   u[2] = (2 - PetscAbsReal(x + y)) * PetscExpReal(-10 * PetscAbsReal(x - y));
329811d88abSStefano Zampini   return PETSC_SUCCESS;
330811d88abSStefano Zampini }
331811d88abSStefano Zampini 
332811d88abSStefano Zampini /* initial conditions for C: eq. 18 */
333811d88abSStefano Zampini static PetscErrorCode initial_conditions_C_2(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
334811d88abSStefano Zampini {
335811d88abSStefano Zampini   u[0] = 0;
336811d88abSStefano Zampini   u[1] = 0;
337811d88abSStefano Zampini   u[2] = 0;
338811d88abSStefano Zampini   return PETSC_SUCCESS;
339811d88abSStefano Zampini }
340811d88abSStefano Zampini 
341811d88abSStefano Zampini /* functionals to be sampled: C * \grad p */
342811d88abSStefano Zampini static void flux(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 f[])
343811d88abSStefano Zampini {
344811d88abSStefano Zampini   const PetscScalar C00     = u[uOff[C_FIELD_ID]];
345811d88abSStefano Zampini   const PetscScalar C01     = u[uOff[C_FIELD_ID] + 1];
346811d88abSStefano Zampini   const PetscScalar C10     = C01;
347811d88abSStefano Zampini   const PetscScalar C11     = u[uOff[C_FIELD_ID] + 2];
348811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
349811d88abSStefano Zampini 
350811d88abSStefano Zampini   f[0] = C00 * gradp[0] + C01 * gradp[1];
351811d88abSStefano Zampini   f[1] = C10 * gradp[0] + C11 * gradp[1];
352811d88abSStefano Zampini }
353811d88abSStefano Zampini 
354811d88abSStefano Zampini /* functionals to be sampled: ||C|| */
355811d88abSStefano Zampini static void normc(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 f[])
356811d88abSStefano Zampini {
357811d88abSStefano Zampini   const PetscScalar C00 = u[uOff[C_FIELD_ID]];
358811d88abSStefano Zampini   const PetscScalar C01 = u[uOff[C_FIELD_ID] + 1];
359811d88abSStefano Zampini   const PetscScalar C11 = u[uOff[C_FIELD_ID] + 2];
360811d88abSStefano Zampini 
361811d88abSStefano Zampini   f[0] = PetscSqrtScalar(NORM2C(C00, C01, C11));
362811d88abSStefano Zampini }
363811d88abSStefano Zampini 
364811d88abSStefano Zampini /* functionals to be sampled: zero */
365811d88abSStefano Zampini static void zero(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 f[])
366811d88abSStefano Zampini {
367811d88abSStefano Zampini   f[0] = 0.0;
368811d88abSStefano Zampini }
369811d88abSStefano Zampini 
370811d88abSStefano Zampini /* functions to be sampled: zero function */
371811d88abSStefano Zampini static PetscErrorCode zerof(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
372811d88abSStefano Zampini {
373811d88abSStefano Zampini   for (PetscInt d = 0; d < Nc; ++d) u[d] = 0.0;
374811d88abSStefano Zampini   return PETSC_SUCCESS;
375811d88abSStefano Zampini }
376811d88abSStefano Zampini 
377811d88abSStefano Zampini /* functions to be sampled: constant function */
378811d88abSStefano Zampini static PetscErrorCode constantf(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
379811d88abSStefano Zampini {
380811d88abSStefano Zampini   PetscInt d;
381811d88abSStefano Zampini   for (d = 0; d < Nc; ++d) u[d] = 1.0;
382811d88abSStefano Zampini   return PETSC_SUCCESS;
383811d88abSStefano Zampini }
384811d88abSStefano Zampini 
385811d88abSStefano Zampini /* application context: customizable parameters */
386811d88abSStefano Zampini typedef struct {
387811d88abSStefano Zampini   PetscReal r;
388811d88abSStefano Zampini   PetscReal eps;
389811d88abSStefano Zampini   PetscReal alpha;
390811d88abSStefano Zampini   PetscReal gamma;
391811d88abSStefano Zampini   PetscReal D;
392811d88abSStefano Zampini   PetscReal c;
393811d88abSStefano Zampini   PetscInt  ic_num;
394811d88abSStefano Zampini   PetscInt  source_num;
395811d88abSStefano Zampini   PetscReal x0[2];
396*204aa523SStefano Zampini   PetscBool lump;
397811d88abSStefano Zampini   PetscBool amr;
398811d88abSStefano Zampini   PetscBool load;
399811d88abSStefano Zampini   char      load_filename[PETSC_MAX_PATH_LEN];
400811d88abSStefano Zampini   PetscBool save;
401811d88abSStefano Zampini   char      save_filename[PETSC_MAX_PATH_LEN];
402811d88abSStefano Zampini   PetscInt  save_every;
403811d88abSStefano Zampini   PetscBool test_restart;
404811d88abSStefano Zampini   PetscBool ellipticity;
405*204aa523SStefano Zampini   PetscInt  fix_c;
406811d88abSStefano Zampini } AppCtx;
407811d88abSStefano Zampini 
408811d88abSStefano Zampini /* process command line options */
409811d88abSStefano Zampini static PetscErrorCode ProcessOptions(AppCtx *options)
410811d88abSStefano Zampini {
411811d88abSStefano Zampini   PetscInt dim = PETSC_STATIC_ARRAY_LENGTH(options->x0);
412811d88abSStefano Zampini 
413811d88abSStefano Zampini   PetscFunctionBeginUser;
414811d88abSStefano Zampini   options->r            = 1.e-1;
415811d88abSStefano Zampini   options->eps          = 1.e-3;
416811d88abSStefano Zampini   options->alpha        = 0.75;
417811d88abSStefano Zampini   options->gamma        = 0.75;
418811d88abSStefano Zampini   options->c            = 5;
419811d88abSStefano Zampini   options->D            = 1.e-2;
420811d88abSStefano Zampini   options->ic_num       = 0;
421811d88abSStefano Zampini   options->source_num   = 0;
422811d88abSStefano Zampini   options->x0[0]        = 0.25;
423811d88abSStefano Zampini   options->x0[1]        = 0.25;
424*204aa523SStefano Zampini   options->lump         = PETSC_FALSE;
425811d88abSStefano Zampini   options->amr          = PETSC_FALSE;
426811d88abSStefano Zampini   options->load         = PETSC_FALSE;
427811d88abSStefano Zampini   options->save         = PETSC_FALSE;
428811d88abSStefano Zampini   options->save_every   = -1;
429811d88abSStefano Zampini   options->test_restart = PETSC_FALSE;
430811d88abSStefano Zampini   options->ellipticity  = PETSC_FALSE;
431*204aa523SStefano Zampini   options->fix_c        = 1; /* 1 means only Jac, 2 means function and Jac */
432811d88abSStefano Zampini 
433811d88abSStefano Zampini   PetscOptionsBegin(PETSC_COMM_WORLD, "", __FILE__, "DMPLEX");
434811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-alpha", "alpha", __FILE__, options->alpha, &options->alpha, NULL));
435811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-gamma", "gamma", __FILE__, options->gamma, &options->gamma, NULL));
436811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-c", "c", __FILE__, options->c, &options->c, NULL));
437811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-d", "D", __FILE__, options->D, &options->D, NULL));
438811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-eps", "eps", __FILE__, options->eps, &options->eps, NULL));
439811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-r", "r", __FILE__, options->r, &options->r, NULL));
440811d88abSStefano Zampini   PetscCall(PetscOptionsRealArray("-x0", "x0", __FILE__, options->x0, &dim, NULL));
441811d88abSStefano Zampini   PetscCall(PetscOptionsInt("-ic_num", "ic_num", __FILE__, options->ic_num, &options->ic_num, NULL));
442811d88abSStefano Zampini   PetscCall(PetscOptionsInt("-source_num", "source_num", __FILE__, options->source_num, &options->source_num, NULL));
443*204aa523SStefano Zampini   PetscCall(PetscOptionsBool("-lump", "use mass lumping", __FILE__, options->lump, &options->lump, NULL));
444*204aa523SStefano Zampini   PetscCall(PetscOptionsInt("-fix_c", "shift conductivity to always be positive semi-definite", __FILE__, options->fix_c, &options->fix_c, NULL));
445811d88abSStefano Zampini   PetscCall(PetscOptionsBool("-amr", "use adaptive mesh refinement", __FILE__, options->amr, &options->amr, NULL));
446811d88abSStefano Zampini   PetscCall(PetscOptionsBool("-test_restart", "test restarting files", __FILE__, options->test_restart, &options->test_restart, NULL));
447811d88abSStefano Zampini   if (!options->test_restart) {
448811d88abSStefano Zampini     PetscCall(PetscOptionsString("-load", "filename with data to be loaded for restarting", __FILE__, options->load_filename, options->load_filename, PETSC_MAX_PATH_LEN, &options->load));
449811d88abSStefano Zampini     PetscCall(PetscOptionsString("-save", "filename with data to be saved for restarting", __FILE__, options->save_filename, options->save_filename, PETSC_MAX_PATH_LEN, &options->save));
450811d88abSStefano Zampini     if (options->save) PetscCall(PetscOptionsInt("-save_every", "save every n timestep (-1 saves only the last)", __FILE__, options->save_every, &options->save_every, NULL));
451811d88abSStefano Zampini   }
452811d88abSStefano Zampini   PetscCall(PetscOptionsBool("-monitor_ellipticity", "Dump locations of ellipticity violation", __FILE__, options->ellipticity, &options->ellipticity, NULL));
453811d88abSStefano Zampini   PetscOptionsEnd();
454811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
455811d88abSStefano Zampini }
456811d88abSStefano Zampini 
457811d88abSStefano Zampini static PetscErrorCode SaveToFile(DM dm, Vec u, const char *filename)
458811d88abSStefano Zampini {
459811d88abSStefano Zampini #if defined(PETSC_HAVE_HDF5)
460811d88abSStefano Zampini   PetscViewerFormat format = PETSC_VIEWER_HDF5_PETSC;
461811d88abSStefano Zampini   PetscViewer       viewer;
462811d88abSStefano Zampini   DM                cdm       = dm;
463811d88abSStefano Zampini   PetscInt          numlevels = 0;
464811d88abSStefano Zampini 
465811d88abSStefano Zampini   PetscFunctionBeginUser;
466811d88abSStefano Zampini   while (cdm) {
467811d88abSStefano Zampini     numlevels++;
468811d88abSStefano Zampini     PetscCall(DMGetCoarseDM(cdm, &cdm));
469811d88abSStefano Zampini   }
470811d88abSStefano Zampini   /* Cannot be set programmatically */
471811d88abSStefano Zampini   PetscCall(PetscOptionsInsertString(NULL, "-dm_plex_view_hdf5_storage_version 3.0.0"));
472811d88abSStefano Zampini   PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject)dm), filename, FILE_MODE_WRITE, &viewer));
473811d88abSStefano Zampini   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "numlevels", PETSC_INT, &numlevels));
474811d88abSStefano Zampini   PetscCall(PetscViewerPushFormat(viewer, format));
475811d88abSStefano Zampini   for (PetscInt level = numlevels - 1; level >= 0; level--) {
476811d88abSStefano Zampini     PetscInt    cc, rr;
477811d88abSStefano Zampini     PetscBool   isRegular, isUniform;
478811d88abSStefano Zampini     const char *dmname;
479811d88abSStefano Zampini     char        groupname[PETSC_MAX_PATH_LEN];
480811d88abSStefano Zampini 
481811d88abSStefano Zampini     PetscCall(PetscSNPrintf(groupname, sizeof(groupname), "level_%" PetscInt_FMT, level));
482811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PushGroup(viewer, groupname));
483811d88abSStefano Zampini     PetscCall(PetscObjectGetName((PetscObject)dm, &dmname));
484811d88abSStefano Zampini     PetscCall(DMGetCoarsenLevel(dm, &cc));
485811d88abSStefano Zampini     PetscCall(DMGetRefineLevel(dm, &rr));
486811d88abSStefano Zampini     PetscCall(DMPlexGetRegularRefinement(dm, &isRegular));
487811d88abSStefano Zampini     PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
488811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "meshname", PETSC_STRING, dmname));
489811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "refinelevel", PETSC_INT, &rr));
490811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "coarsenlevel", PETSC_INT, &cc));
491811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "refRegular", PETSC_BOOL, &isRegular));
492811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "refUniform", PETSC_BOOL, &isUniform));
493811d88abSStefano Zampini     PetscCall(DMPlexTopologyView(dm, viewer));
494811d88abSStefano Zampini     PetscCall(DMPlexLabelsView(dm, viewer));
495811d88abSStefano Zampini     PetscCall(DMPlexCoordinatesView(dm, viewer));
496811d88abSStefano Zampini     PetscCall(DMPlexSectionView(dm, viewer, NULL));
497811d88abSStefano Zampini     if (level == numlevels - 1) {
498811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)u, "solution_"));
499811d88abSStefano Zampini       PetscCall(DMPlexGlobalVectorView(dm, viewer, NULL, u));
500811d88abSStefano Zampini     }
501811d88abSStefano Zampini     if (level) {
502811d88abSStefano Zampini       PetscInt        cStart, cEnd, ccStart, ccEnd, cpStart;
503811d88abSStefano Zampini       DMPolytopeType  ct;
504811d88abSStefano Zampini       DMPlexTransform tr;
505811d88abSStefano Zampini       DM              sdm;
506811d88abSStefano Zampini       PetscScalar    *array;
507811d88abSStefano Zampini       PetscSection    section;
508811d88abSStefano Zampini       Vec             map;
509811d88abSStefano Zampini       IS              gis;
510811d88abSStefano Zampini       const PetscInt *gidx;
511811d88abSStefano Zampini 
512811d88abSStefano Zampini       PetscCall(DMGetCoarseDM(dm, &cdm));
513811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
514811d88abSStefano Zampini       PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
515811d88abSStefano Zampini       PetscCall(PetscSectionSetChart(section, cStart, cEnd));
516811d88abSStefano Zampini       for (PetscInt c = cStart; c < cEnd; c++) PetscCall(PetscSectionSetDof(section, c, 1));
517811d88abSStefano Zampini       PetscCall(PetscSectionSetUp(section));
518811d88abSStefano Zampini 
519811d88abSStefano Zampini       PetscCall(DMClone(dm, &sdm));
520811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)sdm, "pdm"));
521811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)section, "pdm_section"));
522811d88abSStefano Zampini       PetscCall(DMSetLocalSection(sdm, section));
523811d88abSStefano Zampini       PetscCall(PetscSectionDestroy(&section));
524811d88abSStefano Zampini 
525811d88abSStefano Zampini       PetscCall(DMGetLocalVector(sdm, &map));
526811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)map, "pdm_map"));
527811d88abSStefano Zampini       PetscCall(VecGetArray(map, &array));
528811d88abSStefano Zampini       PetscCall(DMPlexTransformCreate(PETSC_COMM_SELF, &tr));
529811d88abSStefano Zampini       PetscCall(DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR));
530811d88abSStefano Zampini       PetscCall(DMPlexTransformSetDM(tr, cdm));
531811d88abSStefano Zampini       PetscCall(DMPlexTransformSetFromOptions(tr));
532811d88abSStefano Zampini       PetscCall(DMPlexTransformSetUp(tr));
533811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(cdm, 0, &ccStart, &ccEnd));
534811d88abSStefano Zampini       PetscCall(DMPlexGetChart(cdm, &cpStart, NULL));
535811d88abSStefano Zampini       PetscCall(DMPlexCreatePointNumbering(cdm, &gis));
536811d88abSStefano Zampini       PetscCall(ISGetIndices(gis, &gidx));
537811d88abSStefano Zampini       for (PetscInt c = ccStart; c < ccEnd; c++) {
538811d88abSStefano Zampini         PetscInt       *rsize, *rcone, *rornt, Nt;
539811d88abSStefano Zampini         DMPolytopeType *rct;
540811d88abSStefano Zampini         PetscInt        gnum = gidx[c - cpStart] >= 0 ? gidx[c - cpStart] : -(gidx[c - cpStart] + 1);
541811d88abSStefano Zampini 
542811d88abSStefano Zampini         PetscCall(DMPlexGetCellType(cdm, c, &ct));
543811d88abSStefano Zampini         PetscCall(DMPlexTransformCellTransform(tr, ct, c, NULL, &Nt, &rct, &rsize, &rcone, &rornt));
544811d88abSStefano Zampini         for (PetscInt r = 0; r < rsize[Nt - 1]; ++r) {
545811d88abSStefano Zampini           PetscInt pNew;
546811d88abSStefano Zampini 
547811d88abSStefano Zampini           PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[Nt - 1], c, r, &pNew));
548811d88abSStefano Zampini           array[pNew - cStart] = gnum;
549811d88abSStefano Zampini         }
550811d88abSStefano Zampini       }
551811d88abSStefano Zampini       PetscCall(ISRestoreIndices(gis, &gidx));
552811d88abSStefano Zampini       PetscCall(ISDestroy(&gis));
553811d88abSStefano Zampini       PetscCall(VecRestoreArray(map, &array));
554811d88abSStefano Zampini       PetscCall(DMPlexTransformDestroy(&tr));
555811d88abSStefano Zampini       PetscCall(DMPlexSectionView(dm, viewer, sdm));
556811d88abSStefano Zampini       PetscCall(DMPlexLocalVectorView(dm, viewer, sdm, map));
557811d88abSStefano Zampini       PetscCall(DMRestoreLocalVector(sdm, &map));
558811d88abSStefano Zampini       PetscCall(DMDestroy(&sdm));
559811d88abSStefano Zampini     }
560811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PopGroup(viewer));
561811d88abSStefano Zampini     PetscCall(DMGetCoarseDM(dm, &dm));
562811d88abSStefano Zampini   }
563811d88abSStefano Zampini   PetscCall(PetscViewerPopFormat(viewer));
564811d88abSStefano Zampini   PetscCall(PetscViewerDestroy(&viewer));
565811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
566811d88abSStefano Zampini #else
567811d88abSStefano Zampini   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Needs HDF5 support. Please reconfigure using --download-hdf5");
568811d88abSStefano Zampini #endif
569811d88abSStefano Zampini }
570811d88abSStefano Zampini 
571811d88abSStefano Zampini static PetscErrorCode LoadFromFile(MPI_Comm comm, const char *filename, DM *odm)
572811d88abSStefano Zampini {
573811d88abSStefano Zampini #if defined(PETSC_HAVE_HDF5)
574811d88abSStefano Zampini   PetscViewerFormat format = PETSC_VIEWER_HDF5_PETSC;
575811d88abSStefano Zampini   PetscViewer       viewer;
576811d88abSStefano Zampini   DM                dm, cdm = NULL;
577811d88abSStefano Zampini   PetscSF           sfXC      = NULL;
578811d88abSStefano Zampini   PetscInt          numlevels = -1;
579811d88abSStefano Zampini 
580811d88abSStefano Zampini   PetscFunctionBeginUser;
581811d88abSStefano Zampini   PetscCall(PetscViewerHDF5Open(comm, filename, FILE_MODE_READ, &viewer));
582811d88abSStefano Zampini   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "numlevels", PETSC_INT, NULL, &numlevels));
583811d88abSStefano Zampini   PetscCall(PetscViewerPushFormat(viewer, format));
584811d88abSStefano Zampini   for (PetscInt level = 0; level < numlevels; level++) {
585811d88abSStefano Zampini     char             groupname[PETSC_MAX_PATH_LEN], *dmname;
586811d88abSStefano Zampini     PetscSF          sfXB, sfBC, sfG;
587811d88abSStefano Zampini     PetscPartitioner part;
588811d88abSStefano Zampini     PetscInt         rr, cc;
589811d88abSStefano Zampini     PetscBool        isRegular, isUniform;
590811d88abSStefano Zampini 
591811d88abSStefano Zampini     PetscCall(DMCreate(comm, &dm));
592811d88abSStefano Zampini     PetscCall(DMSetType(dm, DMPLEX));
593811d88abSStefano Zampini     PetscCall(PetscSNPrintf(groupname, sizeof(groupname), "level_%" PetscInt_FMT, level));
594811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PushGroup(viewer, groupname));
595811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "meshname", PETSC_STRING, NULL, &dmname));
596811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "refinelevel", PETSC_INT, NULL, &rr));
597811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "coarsenlevel", PETSC_INT, NULL, &cc));
598811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "refRegular", PETSC_BOOL, NULL, &isRegular));
599811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "refUniform", PETSC_BOOL, NULL, &isUniform));
600811d88abSStefano Zampini     PetscCall(PetscObjectSetName((PetscObject)dm, dmname));
601811d88abSStefano Zampini     PetscCall(DMPlexTopologyLoad(dm, viewer, &sfXB));
602811d88abSStefano Zampini     PetscCall(DMPlexLabelsLoad(dm, viewer, sfXB));
603811d88abSStefano Zampini     PetscCall(DMPlexCoordinatesLoad(dm, viewer, sfXB));
604811d88abSStefano Zampini     PetscCall(DMPlexGetPartitioner(dm, &part));
605811d88abSStefano Zampini     if (!level) { /* partition the coarse level only */
606811d88abSStefano Zampini       PetscCall(PetscPartitionerSetFromOptions(part));
607811d88abSStefano Zampini     } else { /* propagate partitioning information from coarser to finer level */
608811d88abSStefano Zampini       DM           sdm;
609811d88abSStefano Zampini       Vec          map;
610811d88abSStefano Zampini       PetscSF      sf;
611811d88abSStefano Zampini       PetscLayout  clayout;
612811d88abSStefano Zampini       PetscScalar *array;
613811d88abSStefano Zampini       PetscInt    *cranks_leaf, *cranks_root, *npoints, *points, *ranks, *starts, *gidxs;
614811d88abSStefano Zampini       PetscInt     nparts, cStart, cEnd, nr, ccStart, ccEnd, cpStart, cpEnd;
615811d88abSStefano Zampini       PetscMPIInt  size, rank;
616811d88abSStefano Zampini 
617811d88abSStefano Zampini       PetscCall(DMClone(dm, &sdm));
618811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)sdm, "pdm"));
619811d88abSStefano Zampini       PetscCall(DMPlexSectionLoad(dm, viewer, sdm, sfXB, NULL, &sf));
620811d88abSStefano Zampini       PetscCall(DMGetLocalVector(sdm, &map));
621811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)map, "pdm_map"));
622811d88abSStefano Zampini       PetscCall(DMPlexLocalVectorLoad(dm, viewer, sdm, sf, map));
623811d88abSStefano Zampini 
624811d88abSStefano Zampini       PetscCallMPI(MPI_Comm_size(comm, &size));
625811d88abSStefano Zampini       PetscCallMPI(MPI_Comm_rank(comm, &rank));
626811d88abSStefano Zampini       nparts = size;
627811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
628811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(cdm, 0, &ccStart, &ccEnd));
629811d88abSStefano Zampini       PetscCall(DMPlexGetChart(cdm, &cpStart, &cpEnd));
630811d88abSStefano Zampini       PetscCall(PetscCalloc1(nparts, &npoints));
631811d88abSStefano Zampini       PetscCall(PetscMalloc4(cEnd - cStart, &points, cEnd - cStart, &ranks, nparts + 1, &starts, cEnd - cStart, &gidxs));
632811d88abSStefano Zampini       PetscCall(PetscSFGetGraph(sfXC, &nr, NULL, NULL, NULL));
633811d88abSStefano Zampini       PetscCall(PetscMalloc2(cpEnd - cpStart, &cranks_leaf, nr, &cranks_root));
634811d88abSStefano Zampini       for (PetscInt c = 0; c < cpEnd - cpStart; c++) cranks_leaf[c] = rank;
635811d88abSStefano Zampini       PetscCall(PetscSFReduceBegin(sfXC, MPIU_INT, cranks_leaf, cranks_root, MPI_REPLACE));
636811d88abSStefano Zampini       PetscCall(PetscSFReduceEnd(sfXC, MPIU_INT, cranks_leaf, cranks_root, MPI_REPLACE));
637811d88abSStefano Zampini 
638811d88abSStefano Zampini       PetscCall(VecGetArray(map, &array));
639811d88abSStefano Zampini       for (PetscInt c = 0; c < cEnd - cStart; c++) gidxs[c] = (PetscInt)PetscRealPart(array[c]);
640811d88abSStefano Zampini       PetscCall(VecRestoreArray(map, &array));
641811d88abSStefano Zampini 
642811d88abSStefano Zampini       PetscCall(PetscLayoutCreate(comm, &clayout));
643811d88abSStefano Zampini       PetscCall(PetscLayoutSetLocalSize(clayout, nr));
644811d88abSStefano Zampini       PetscCall(PetscSFSetGraphLayout(sf, clayout, cEnd - cStart, NULL, PETSC_OWN_POINTER, gidxs));
645811d88abSStefano Zampini       PetscCall(PetscLayoutDestroy(&clayout));
646811d88abSStefano Zampini 
647811d88abSStefano Zampini       PetscCall(PetscSFBcastBegin(sf, MPIU_INT, cranks_root, ranks, MPI_REPLACE));
648811d88abSStefano Zampini       PetscCall(PetscSFBcastEnd(sf, MPIU_INT, cranks_root, ranks, MPI_REPLACE));
649811d88abSStefano Zampini       PetscCall(PetscSFDestroy(&sf));
650811d88abSStefano Zampini       PetscCall(PetscFree2(cranks_leaf, cranks_root));
651811d88abSStefano Zampini       for (PetscInt c = 0; c < cEnd - cStart; c++) npoints[ranks[c]]++;
652811d88abSStefano Zampini 
653811d88abSStefano Zampini       starts[0] = 0;
654811d88abSStefano Zampini       for (PetscInt c = 0; c < nparts; c++) starts[c + 1] = starts[c] + npoints[c];
655811d88abSStefano Zampini       for (PetscInt c = 0; c < cEnd - cStart; c++) points[starts[ranks[c]]++] = c;
656811d88abSStefano Zampini       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
657811d88abSStefano Zampini       PetscCall(PetscPartitionerShellSetPartition(part, nparts, npoints, points));
658811d88abSStefano Zampini       PetscCall(PetscFree(npoints));
659811d88abSStefano Zampini       PetscCall(PetscFree4(points, ranks, starts, gidxs));
660811d88abSStefano Zampini       PetscCall(DMRestoreLocalVector(sdm, &map));
661811d88abSStefano Zampini       PetscCall(DMDestroy(&sdm));
662811d88abSStefano Zampini     }
663811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfXC));
664811d88abSStefano Zampini     PetscCall(DMPlexDistribute(dm, 0, &sfBC, odm));
665811d88abSStefano Zampini     if (*odm) {
666811d88abSStefano Zampini       PetscCall(DMDestroy(&dm));
667811d88abSStefano Zampini       dm   = *odm;
668811d88abSStefano Zampini       *odm = NULL;
669811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)dm, dmname));
670811d88abSStefano Zampini     }
671811d88abSStefano Zampini     if (sfBC) PetscCall(PetscSFCompose(sfXB, sfBC, &sfXC));
672811d88abSStefano Zampini     else {
673811d88abSStefano Zampini       PetscCall(PetscObjectReference((PetscObject)sfXB));
674811d88abSStefano Zampini       sfXC = sfXB;
675811d88abSStefano Zampini     }
676811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfXB));
677811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfBC));
678811d88abSStefano Zampini     PetscCall(DMSetCoarsenLevel(dm, cc));
679811d88abSStefano Zampini     PetscCall(DMSetRefineLevel(dm, rr));
680811d88abSStefano Zampini     PetscCall(DMPlexSetRegularRefinement(dm, isRegular));
681811d88abSStefano Zampini     PetscCall(DMPlexSetRefinementUniform(dm, isUniform));
682811d88abSStefano Zampini     PetscCall(DMPlexSectionLoad(dm, viewer, NULL, sfXC, &sfG, NULL));
683811d88abSStefano Zampini     if (level == numlevels - 1) {
684811d88abSStefano Zampini       Vec u;
685811d88abSStefano Zampini 
686811d88abSStefano Zampini       PetscCall(DMGetNamedGlobalVector(dm, "solution_", &u));
687811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)u, "solution_"));
688811d88abSStefano Zampini       PetscCall(DMPlexGlobalVectorLoad(dm, viewer, NULL, sfG, u));
689811d88abSStefano Zampini       PetscCall(DMRestoreNamedGlobalVector(dm, "solution_", &u));
690811d88abSStefano Zampini     }
691811d88abSStefano Zampini     PetscCall(PetscFree(dmname));
692811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfG));
693811d88abSStefano Zampini     PetscCall(DMSetCoarseDM(dm, cdm));
694811d88abSStefano Zampini     PetscCall(DMDestroy(&cdm));
695811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PopGroup(viewer));
696811d88abSStefano Zampini     cdm = dm;
697811d88abSStefano Zampini   }
698811d88abSStefano Zampini   *odm = dm;
699811d88abSStefano Zampini   PetscCall(PetscViewerPopFormat(viewer));
700811d88abSStefano Zampini   PetscCall(PetscViewerDestroy(&viewer));
701811d88abSStefano Zampini   PetscCall(PetscSFDestroy(&sfXC));
702811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
703811d88abSStefano Zampini #else
704811d88abSStefano Zampini   SETERRQ(comm, PETSC_ERR_SUP, "Needs HDF5 support. Please reconfigure using --download-hdf5");
705811d88abSStefano Zampini #endif
706811d88abSStefano Zampini }
707811d88abSStefano Zampini 
708811d88abSStefano Zampini /* Project source function and make it zero-mean */
709811d88abSStefano Zampini static PetscErrorCode ProjectSource(DM dm, PetscReal time, AppCtx *ctx)
710811d88abSStefano Zampini {
711811d88abSStefano Zampini   PetscInt    id = C_FIELD_ID;
712811d88abSStefano Zampini   DM          dmAux;
713811d88abSStefano Zampini   Vec         u, lu;
714811d88abSStefano Zampini   IS          is;
715811d88abSStefano Zampini   void       *ctxs[NUM_FIELDS];
716811d88abSStefano Zampini   PetscScalar vals[NUM_FIELDS];
717811d88abSStefano Zampini   PetscDS     ds;
718811d88abSStefano Zampini   PetscErrorCode (*funcs[NUM_FIELDS])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
719811d88abSStefano Zampini 
720811d88abSStefano Zampini   PetscFunctionBeginUser;
721811d88abSStefano Zampini   switch (ctx->source_num) {
722811d88abSStefano Zampini   case 0:
723811d88abSStefano Zampini     funcs[P_FIELD_ID] = source_0;
724811d88abSStefano Zampini     ctxs[P_FIELD_ID]  = ctx->x0;
725811d88abSStefano Zampini     break;
726811d88abSStefano Zampini   case 1:
727811d88abSStefano Zampini     funcs[P_FIELD_ID] = source_1;
728811d88abSStefano Zampini     ctxs[P_FIELD_ID]  = NULL;
729811d88abSStefano Zampini     break;
730811d88abSStefano Zampini   default:
731d8b4a066SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknown source");
732811d88abSStefano Zampini   }
733811d88abSStefano Zampini   funcs[C_FIELD_ID] = zerof;
734811d88abSStefano Zampini   ctxs[C_FIELD_ID]  = NULL;
735811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
736811d88abSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &u));
737811d88abSStefano Zampini   PetscCall(DMProjectFunction(dm, time, funcs, ctxs, INSERT_ALL_VALUES, u));
738811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, average));
739811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
740811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
741811d88abSStefano Zampini   PetscCall(VecShift(u, -vals[P_FIELD_ID]));
742811d88abSStefano Zampini   PetscCall(DMCreateSubDM(dm, 1, &id, &is, NULL));
743811d88abSStefano Zampini   PetscCall(VecISSet(u, is, 0));
744811d88abSStefano Zampini   PetscCall(ISDestroy(&is));
745811d88abSStefano Zampini 
746811d88abSStefano Zampini   /* Attach source vector as auxiliary vector:
747811d88abSStefano Zampini      Use a different DM to break ref cycles */
748811d88abSStefano Zampini   PetscCall(DMClone(dm, &dmAux));
749811d88abSStefano Zampini   PetscCall(DMCopyDisc(dm, dmAux));
750811d88abSStefano Zampini   PetscCall(DMCreateLocalVector(dmAux, &lu));
751811d88abSStefano Zampini   PetscCall(DMDestroy(&dmAux));
752811d88abSStefano Zampini   PetscCall(DMGlobalToLocal(dm, u, INSERT_VALUES, lu));
753811d88abSStefano Zampini   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, lu));
754811d88abSStefano Zampini   PetscCall(VecViewFromOptions(lu, NULL, "-aux_view"));
755811d88abSStefano Zampini   PetscCall(VecDestroy(&lu));
756811d88abSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &u));
757811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
758811d88abSStefano Zampini }
759811d88abSStefano Zampini 
760811d88abSStefano Zampini /* callback for the creation of the potential null space */
761811d88abSStefano Zampini static PetscErrorCode CreatePotentialNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
762811d88abSStefano Zampini {
763811d88abSStefano Zampini   Vec vec;
764811d88abSStefano Zampini   PetscErrorCode (*funcs[NUM_FIELDS])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zerof};
765811d88abSStefano Zampini 
766811d88abSStefano Zampini   PetscFunctionBeginUser;
767811d88abSStefano Zampini   funcs[nfield] = constantf;
768811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &vec));
769811d88abSStefano Zampini   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
770811d88abSStefano Zampini   PetscCall(VecNormalize(vec, NULL));
771811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)vec, "Potential Null Space"));
772811d88abSStefano Zampini   PetscCall(VecViewFromOptions(vec, NULL, "-potential_nullspace_view"));
773811d88abSStefano Zampini   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace));
774811d88abSStefano Zampini   /* break ref cycles */
775811d88abSStefano Zampini   PetscCall(VecSetDM(vec, NULL));
776811d88abSStefano Zampini   PetscCall(VecDestroy(&vec));
777811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
778811d88abSStefano Zampini }
779811d88abSStefano Zampini 
780*204aa523SStefano Zampini static PetscErrorCode DMGetLumpedMass(DM dm, PetscBool local, Vec *lumped_mass)
781*204aa523SStefano Zampini {
782*204aa523SStefano Zampini   PetscBool has;
783*204aa523SStefano Zampini 
784*204aa523SStefano Zampini   PetscFunctionBeginUser;
785*204aa523SStefano Zampini   if (local) {
786*204aa523SStefano Zampini     PetscCall(DMHasNamedLocalVector(dm, "lumped_mass", &has));
787*204aa523SStefano Zampini     PetscCall(DMGetNamedLocalVector(dm, "lumped_mass", lumped_mass));
788*204aa523SStefano Zampini   } else {
789*204aa523SStefano Zampini     PetscCall(DMHasNamedGlobalVector(dm, "lumped_mass", &has));
790*204aa523SStefano Zampini     PetscCall(DMGetNamedGlobalVector(dm, "lumped_mass", lumped_mass));
791*204aa523SStefano Zampini   }
792*204aa523SStefano Zampini   if (!has) {
793*204aa523SStefano Zampini     Vec w;
794*204aa523SStefano Zampini     IS  is;
795*204aa523SStefano Zampini 
796*204aa523SStefano Zampini     PetscCall(PetscObjectQuery((PetscObject)dm, "IS potential", (PetscObject *)&is));
797*204aa523SStefano Zampini     if (!is) {
798*204aa523SStefano Zampini       PetscInt fields[NUM_FIELDS] = {C_FIELD_ID, P_FIELD_ID};
799*204aa523SStefano Zampini 
800*204aa523SStefano Zampini       PetscCall(DMCreateSubDM(dm, NUM_FIELDS - 1, fields + 1, &is, NULL));
801*204aa523SStefano Zampini       PetscCall(PetscObjectCompose((PetscObject)dm, "IS potential", (PetscObject)is));
802*204aa523SStefano Zampini       PetscCall(PetscObjectDereference((PetscObject)is));
803*204aa523SStefano Zampini     }
804*204aa523SStefano Zampini     if (local) {
805*204aa523SStefano Zampini       Vec w2, wg;
806*204aa523SStefano Zampini 
807*204aa523SStefano Zampini       PetscCall(DMCreateMassMatrixLumped(dm, &w, NULL));
808*204aa523SStefano Zampini       PetscCall(DMGetGlobalVector(dm, &wg));
809*204aa523SStefano Zampini       PetscCall(DMGetLocalVector(dm, &w2));
810*204aa523SStefano Zampini       PetscCall(VecSet(w2, 0.0));
811*204aa523SStefano Zampini       PetscCall(VecSet(wg, 1.0));
812*204aa523SStefano Zampini       PetscCall(VecISSet(wg, is, 0.0));
813*204aa523SStefano Zampini       PetscCall(DMGlobalToLocal(dm, wg, INSERT_VALUES, w2));
814*204aa523SStefano Zampini       PetscCall(VecPointwiseMult(w, w, w2));
815*204aa523SStefano Zampini       PetscCall(DMRestoreGlobalVector(dm, &wg));
816*204aa523SStefano Zampini       PetscCall(DMRestoreLocalVector(dm, &w2));
817*204aa523SStefano Zampini     } else {
818*204aa523SStefano Zampini       PetscCall(DMCreateMassMatrixLumped(dm, NULL, &w));
819*204aa523SStefano Zampini       PetscCall(VecISSet(w, is, 0.0));
820*204aa523SStefano Zampini     }
821*204aa523SStefano Zampini     PetscCall(VecCopy(w, *lumped_mass));
822*204aa523SStefano Zampini     PetscCall(VecDestroy(&w));
823*204aa523SStefano Zampini   }
824*204aa523SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
825*204aa523SStefano Zampini }
826*204aa523SStefano Zampini 
827*204aa523SStefano Zampini static PetscErrorCode DMRestoreLumpedMass(DM dm, PetscBool local, Vec *lumped_mass)
828*204aa523SStefano Zampini {
829*204aa523SStefano Zampini   PetscFunctionBeginUser;
830*204aa523SStefano Zampini   if (local) PetscCall(DMRestoreNamedLocalVector(dm, "lumped_mass", lumped_mass));
831*204aa523SStefano Zampini   else PetscCall(DMRestoreNamedGlobalVector(dm, "lumped_mass", lumped_mass));
832*204aa523SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
833*204aa523SStefano Zampini }
834*204aa523SStefano Zampini 
835*204aa523SStefano Zampini /* callbacks for lumped mass matrix residual and Jacobian */
836*204aa523SStefano Zampini static PetscErrorCode DMPlexTSComputeIFunctionFEM_Lumped(DM dm, PetscReal time, Vec locX, Vec locX_t, Vec locF, void *user)
837*204aa523SStefano Zampini {
838*204aa523SStefano Zampini   Vec work, local_lumped_mass;
839*204aa523SStefano Zampini 
840*204aa523SStefano Zampini   PetscFunctionBeginUser;
841*204aa523SStefano Zampini   PetscCall(DMGetLumpedMass(dm, PETSC_TRUE, &local_lumped_mass));
842*204aa523SStefano Zampini   PetscCall(DMGetLocalVector(dm, &work));
843*204aa523SStefano Zampini   PetscCall(VecSet(work, 0.0));
844*204aa523SStefano Zampini   PetscCall(DMPlexTSComputeIFunctionFEM(dm, time, locX, work, locF, user));
845*204aa523SStefano Zampini   PetscCall(VecPointwiseMult(work, locX_t, local_lumped_mass));
846*204aa523SStefano Zampini   PetscCall(VecAXPY(locF, 1.0, work));
847*204aa523SStefano Zampini   PetscCall(DMRestoreLocalVector(dm, &work));
848*204aa523SStefano Zampini   PetscCall(DMRestoreLumpedMass(dm, PETSC_TRUE, &local_lumped_mass));
849*204aa523SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
850*204aa523SStefano Zampini }
851*204aa523SStefano Zampini 
852*204aa523SStefano Zampini static PetscErrorCode DMPlexTSComputeIJacobianFEM_Lumped(DM dm, PetscReal time, Vec locX, Vec locX_t, PetscReal X_tShift, Mat Jac, Mat JacP, void *user)
853*204aa523SStefano Zampini {
854*204aa523SStefano Zampini   Vec lumped_mass, work;
855*204aa523SStefano Zampini 
856*204aa523SStefano Zampini   PetscFunctionBeginUser;
857*204aa523SStefano Zampini   // XXX CHECK DIRK
858*204aa523SStefano Zampini   PetscCall(DMGetLumpedMass(dm, PETSC_FALSE, &lumped_mass));
859*204aa523SStefano Zampini   PetscCall(DMPlexTSComputeIJacobianFEM(dm, time, locX, locX_t, 0.0, Jac, JacP, user));
860*204aa523SStefano Zampini   PetscCall(DMGetGlobalVector(dm, &work));
861*204aa523SStefano Zampini   PetscCall(VecAXPBY(work, X_tShift, 0.0, lumped_mass));
862*204aa523SStefano Zampini   PetscCall(MatDiagonalSet(JacP, work, ADD_VALUES));
863*204aa523SStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &work));
864*204aa523SStefano Zampini   PetscCall(DMRestoreLumpedMass(dm, PETSC_FALSE, &lumped_mass));
865*204aa523SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
866*204aa523SStefano Zampini }
867*204aa523SStefano Zampini 
868811d88abSStefano Zampini /* customize residuals and Jacobians */
869811d88abSStefano Zampini static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
870811d88abSStefano Zampini {
871811d88abSStefano Zampini   PetscDS     ds;
872811d88abSStefano Zampini   PetscInt    cdim, dim;
873811d88abSStefano Zampini   PetscScalar constants[NUM_CONSTANTS];
874811d88abSStefano Zampini 
875811d88abSStefano Zampini   PetscFunctionBeginUser;
876811d88abSStefano Zampini   constants[R_ID]     = ctx->r;
877811d88abSStefano Zampini   constants[EPS_ID]   = ctx->eps;
878811d88abSStefano Zampini   constants[ALPHA_ID] = ctx->alpha;
879811d88abSStefano Zampini   constants[GAMMA_ID] = ctx->gamma;
880811d88abSStefano Zampini   constants[D_ID]     = ctx->D;
881811d88abSStefano Zampini   constants[C2_ID]    = ctx->c * ctx->c;
882*204aa523SStefano Zampini   constants[FIXC_ID]  = ctx->fix_c;
883811d88abSStefano Zampini 
884811d88abSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
885811d88abSStefano Zampini   PetscCall(DMGetCoordinateDim(dm, &cdim));
886811d88abSStefano Zampini   PetscCheck(dim == 2 && cdim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only for 2D meshes");
887811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
888811d88abSStefano Zampini   PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants));
889811d88abSStefano Zampini   PetscCall(PetscDSSetImplicit(ds, C_FIELD_ID, PETSC_TRUE));
890811d88abSStefano Zampini   PetscCall(PetscDSSetImplicit(ds, P_FIELD_ID, PETSC_TRUE));
891811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, energy));
892811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
893811d88abSStefano Zampini   PetscCall(PetscDSSetResidual(ds, C_FIELD_ID, C_0, C_1));
894811d88abSStefano Zampini   PetscCall(PetscDSSetResidual(ds, P_FIELD_ID, P_0, P_1));
895811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, C_FIELD_ID, C_FIELD_ID, JC_0_c0c0, NULL, NULL, JC_1_c1c1));
896811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, C_FIELD_ID, P_FIELD_ID, NULL, JC_0_c0p1, NULL, NULL));
897811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, P_FIELD_ID, C_FIELD_ID, NULL, NULL, JP_1_p1c0, NULL));
898811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, P_FIELD_ID, P_FIELD_ID, NULL, NULL, NULL, JP_1_p1p1));
899811d88abSStefano Zampini 
900811d88abSStefano Zampini   /* Attach potential nullspace */
901811d88abSStefano Zampini   PetscCall(DMSetNullSpaceConstructor(dm, P_FIELD_ID, CreatePotentialNullSpace));
902811d88abSStefano Zampini 
903811d88abSStefano Zampini   /* Attach source function as auxiliary vector */
904811d88abSStefano Zampini   PetscCall(ProjectSource(dm, 0, ctx));
905811d88abSStefano Zampini 
906811d88abSStefano Zampini   /* Add callbacks */
907*204aa523SStefano Zampini   if (ctx->lump) {
908*204aa523SStefano Zampini     PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM_Lumped, NULL));
909*204aa523SStefano Zampini     PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM_Lumped, NULL));
910*204aa523SStefano Zampini   } else {
911811d88abSStefano Zampini     PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, NULL));
912811d88abSStefano Zampini     PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, NULL));
913*204aa523SStefano Zampini   }
914*204aa523SStefano Zampini   /* This is not really needed because we use Neumann boundaries */
915*204aa523SStefano Zampini   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, NULL));
916811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
917811d88abSStefano Zampini }
918811d88abSStefano Zampini 
919811d88abSStefano Zampini /* setup discrete spaces and residuals */
920811d88abSStefano Zampini static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
921811d88abSStefano Zampini {
922811d88abSStefano Zampini   DM           plex, cdm = dm;
923811d88abSStefano Zampini   PetscFE      feC, feP;
924811d88abSStefano Zampini   PetscBool    simplex;
925811d88abSStefano Zampini   PetscInt     dim;
926811d88abSStefano Zampini   MPI_Comm     comm = PetscObjectComm((PetscObject)dm);
927811d88abSStefano Zampini   MatNullSpace nsp;
928811d88abSStefano Zampini 
929811d88abSStefano Zampini   PetscFunctionBeginUser;
930811d88abSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
931811d88abSStefano Zampini 
932811d88abSStefano Zampini   PetscCall(DMConvert(dm, DMPLEX, &plex));
933811d88abSStefano Zampini   PetscCall(DMPlexIsSimplex(plex, &simplex));
934811d88abSStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &simplex, 1, MPIU_BOOL, MPI_LOR, comm));
935811d88abSStefano Zampini   PetscCall(DMDestroy(&plex));
936811d88abSStefano Zampini 
937*204aa523SStefano Zampini   /* We model Cij with Cij = Cji -> dim*(dim+1)/2 components */
938811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(comm, dim, (dim * (dim + 1)) / 2, simplex, "c_", -1, &feC));
939811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)feC, "conductivity"));
940811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "p_", -1, &feP));
941811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)feP, "potential"));
942811d88abSStefano Zampini   PetscCall(MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, &nsp));
943811d88abSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)feP, "nullspace", (PetscObject)nsp));
944811d88abSStefano Zampini   PetscCall(MatNullSpaceDestroy(&nsp));
945*204aa523SStefano Zampini   PetscCall(PetscFECopyQuadrature(feP, feC));
946811d88abSStefano Zampini 
947811d88abSStefano Zampini   PetscCall(DMSetNumFields(dm, 2));
948811d88abSStefano Zampini   PetscCall(DMSetField(dm, C_FIELD_ID, NULL, (PetscObject)feC));
949811d88abSStefano Zampini   PetscCall(DMSetField(dm, P_FIELD_ID, NULL, (PetscObject)feP));
950811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feC));
951811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feP));
952811d88abSStefano Zampini   PetscCall(DMCreateDS(dm));
953811d88abSStefano Zampini   while (cdm) {
954811d88abSStefano Zampini     PetscCall(DMCopyDisc(dm, cdm));
955811d88abSStefano Zampini     PetscCall(SetupProblem(cdm, ctx));
956811d88abSStefano Zampini     PetscCall(DMGetCoarseDM(cdm, &cdm));
957811d88abSStefano Zampini   }
958811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
959811d88abSStefano Zampini }
960811d88abSStefano Zampini 
961811d88abSStefano Zampini /* Create mesh by command line options */
962811d88abSStefano Zampini static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
963811d88abSStefano Zampini {
964811d88abSStefano Zampini   PetscFunctionBeginUser;
965811d88abSStefano Zampini   if (ctx->load) {
966811d88abSStefano Zampini     PetscInt  refine = 0;
967811d88abSStefano Zampini     PetscBool isHierarchy;
968811d88abSStefano Zampini     DM       *dms;
969811d88abSStefano Zampini     char      typeName[256];
970811d88abSStefano Zampini     PetscBool flg;
971811d88abSStefano Zampini 
972811d88abSStefano Zampini     PetscCall(LoadFromFile(comm, ctx->load_filename, dm));
973811d88abSStefano Zampini     PetscOptionsBegin(comm, "", "Additional mesh options", "DMPLEX");
974811d88abSStefano Zampini     PetscCall(PetscOptionsFList("-dm_mat_type", "Matrix type used for created matrices", "DMSetMatType", MatList, MATAIJ, typeName, sizeof(typeName), &flg));
975811d88abSStefano Zampini     if (flg) PetscCall(DMSetMatType(*dm, typeName));
976811d88abSStefano Zampini     PetscCall(PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL, 0));
977811d88abSStefano Zampini     PetscCall(PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy, 0));
978811d88abSStefano Zampini     PetscOptionsEnd();
979811d88abSStefano Zampini     if (refine) {
980811d88abSStefano Zampini       PetscCall(SetupDiscretization(*dm, ctx));
981811d88abSStefano Zampini       PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE));
982811d88abSStefano Zampini     }
983811d88abSStefano Zampini     PetscCall(PetscCalloc1(refine, &dms));
984811d88abSStefano Zampini     if (isHierarchy) PetscCall(DMRefineHierarchy(*dm, refine, dms));
985811d88abSStefano Zampini     for (PetscInt r = 0; r < refine; r++) {
986811d88abSStefano Zampini       Mat M;
987811d88abSStefano Zampini       DM  dmr = dms[r];
988811d88abSStefano Zampini       Vec u, ur;
989811d88abSStefano Zampini 
990811d88abSStefano Zampini       if (!isHierarchy) {
99155bffe1bSPierre Jolivet         PetscCall(DMRefine(*dm, PetscObjectComm((PetscObject)*dm), &dmr));
992811d88abSStefano Zampini         PetscCall(DMSetCoarseDM(dmr, *dm));
993811d88abSStefano Zampini       }
994811d88abSStefano Zampini       PetscCall(DMCreateInterpolation(*dm, dmr, &M, NULL));
995811d88abSStefano Zampini       PetscCall(DMGetNamedGlobalVector(*dm, "solution_", &u));
996811d88abSStefano Zampini       PetscCall(DMGetNamedGlobalVector(dmr, "solution_", &ur));
997811d88abSStefano Zampini       PetscCall(MatInterpolate(M, u, ur));
998811d88abSStefano Zampini       PetscCall(DMRestoreNamedGlobalVector(*dm, "solution_", &u));
999811d88abSStefano Zampini       PetscCall(DMRestoreNamedGlobalVector(dmr, "solution_", &ur));
1000811d88abSStefano Zampini       PetscCall(MatDestroy(&M));
1001811d88abSStefano Zampini       if (!isHierarchy) PetscCall(DMSetCoarseDM(dmr, NULL));
1002811d88abSStefano Zampini       PetscCall(DMDestroy(dm));
1003811d88abSStefano Zampini       *dm = dmr;
1004811d88abSStefano Zampini     }
1005811d88abSStefano Zampini     if (refine && !isHierarchy) PetscCall(DMSetRefineLevel(*dm, 0));
1006811d88abSStefano Zampini     PetscCall(PetscFree(dms));
1007811d88abSStefano Zampini   } else {
1008811d88abSStefano Zampini     PetscCall(DMCreate(comm, dm));
1009811d88abSStefano Zampini     PetscCall(DMSetType(*dm, DMPLEX));
1010811d88abSStefano Zampini     PetscCall(DMSetFromOptions(*dm));
1011811d88abSStefano Zampini     PetscCall(DMLocalizeCoordinates(*dm));
1012811d88abSStefano Zampini     {
1013811d88abSStefano Zampini       char      convType[256];
1014811d88abSStefano Zampini       PetscBool flg;
1015811d88abSStefano Zampini       PetscOptionsBegin(comm, "", "Additional mesh options", "DMPLEX");
1016811d88abSStefano Zampini       PetscCall(PetscOptionsFList("-dm_plex_convert_type", "Convert DMPlex to another format", __FILE__, DMList, DMPLEX, convType, 256, &flg));
1017811d88abSStefano Zampini       PetscOptionsEnd();
1018811d88abSStefano Zampini       if (flg) {
1019811d88abSStefano Zampini         DM dmConv;
1020811d88abSStefano Zampini         PetscCall(DMConvert(*dm, convType, &dmConv));
1021811d88abSStefano Zampini         if (dmConv) {
1022811d88abSStefano Zampini           PetscCall(DMDestroy(dm));
1023811d88abSStefano Zampini           *dm = dmConv;
1024811d88abSStefano Zampini           PetscCall(DMSetFromOptions(*dm));
1025811d88abSStefano Zampini           PetscCall(DMSetUp(*dm));
1026811d88abSStefano Zampini         }
1027811d88abSStefano Zampini       }
1028811d88abSStefano Zampini     }
1029811d88abSStefano Zampini   }
1030811d88abSStefano Zampini   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1031811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1032811d88abSStefano Zampini }
1033811d88abSStefano Zampini 
1034811d88abSStefano Zampini /* Make potential field zero mean */
1035811d88abSStefano Zampini static PetscErrorCode ZeroMeanPotential(DM dm, Vec u)
1036811d88abSStefano Zampini {
1037811d88abSStefano Zampini   PetscScalar vals[NUM_FIELDS];
1038811d88abSStefano Zampini   PetscDS     ds;
1039811d88abSStefano Zampini   IS          is;
1040811d88abSStefano Zampini 
1041811d88abSStefano Zampini   PetscFunctionBeginUser;
1042811d88abSStefano Zampini   PetscCall(PetscObjectQuery((PetscObject)dm, "IS potential", (PetscObject *)&is));
1043811d88abSStefano Zampini   PetscCheck(is, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing potential IS");
1044811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
1045811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, average));
1046811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
1047811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
1048811d88abSStefano Zampini   PetscCall(VecISShift(u, is, -vals[P_FIELD_ID]));
1049811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
1050811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1051811d88abSStefano Zampini }
1052811d88abSStefano Zampini 
1053811d88abSStefano Zampini /* Compute initial conditions and exclude potential from local truncation error
1054811d88abSStefano Zampini    Since we are solving a DAE, once the initial conditions for the differential
1055811d88abSStefano Zampini    variables are set, we need to compute the corresponding value for the
1056811d88abSStefano Zampini    algebraic variables. We do so by creating a subDM for the potential only
1057811d88abSStefano Zampini    and solve a static problem with SNES */
1058811d88abSStefano Zampini static PetscErrorCode SetInitialConditionsAndTolerances(TS ts, PetscInt nv, Vec vecs[], PetscBool valid)
1059811d88abSStefano Zampini {
1060811d88abSStefano Zampini   DM         dm;
1061811d88abSStefano Zampini   Vec        tu, u, p, lsource, subaux, vatol, vrtol;
1062811d88abSStefano Zampini   PetscReal  t, atol, rtol;
1063811d88abSStefano Zampini   PetscInt   fields[NUM_FIELDS] = {C_FIELD_ID, P_FIELD_ID};
1064811d88abSStefano Zampini   IS         isp;
1065811d88abSStefano Zampini   DM         dmp;
1066811d88abSStefano Zampini   VecScatter sctp = NULL;
1067811d88abSStefano Zampini   PetscDS    ds;
1068811d88abSStefano Zampini   SNES       snes;
1069811d88abSStefano Zampini   KSP        ksp;
1070811d88abSStefano Zampini   PC         pc;
1071811d88abSStefano Zampini   AppCtx    *ctx;
1072811d88abSStefano Zampini 
1073811d88abSStefano Zampini   PetscFunctionBeginUser;
1074811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1075811d88abSStefano Zampini   PetscCall(TSGetApplicationContext(ts, &ctx));
1076*204aa523SStefano Zampini   if (valid) {
1077811d88abSStefano Zampini     PetscCall(DMCreateSubDM(dm, NUM_FIELDS - 1, fields + 1, &isp, NULL));
1078811d88abSStefano Zampini     PetscCall(PetscObjectCompose((PetscObject)dm, "IS potential", (PetscObject)isp));
1079811d88abSStefano Zampini     PetscCall(DMCreateGlobalVector(dm, &vatol));
1080811d88abSStefano Zampini     PetscCall(DMCreateGlobalVector(dm, &vrtol));
1081811d88abSStefano Zampini     PetscCall(TSGetTolerances(ts, &atol, NULL, &rtol, NULL));
1082811d88abSStefano Zampini     PetscCall(VecSet(vatol, atol));
1083811d88abSStefano Zampini     PetscCall(VecISSet(vatol, isp, -1));
1084811d88abSStefano Zampini     PetscCall(VecSet(vrtol, rtol));
1085811d88abSStefano Zampini     PetscCall(VecISSet(vrtol, isp, -1));
1086811d88abSStefano Zampini     PetscCall(TSSetTolerances(ts, atol, vatol, rtol, vrtol));
1087811d88abSStefano Zampini     PetscCall(VecDestroy(&vatol));
1088811d88abSStefano Zampini     PetscCall(VecDestroy(&vrtol));
1089811d88abSStefano Zampini     PetscCall(ISDestroy(&isp));
1090*204aa523SStefano Zampini     for (PetscInt i = 0; i < nv; i++) { PetscCall(ZeroMeanPotential(dm, vecs[i])); }
1091811d88abSStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
1092811d88abSStefano Zampini   }
1093811d88abSStefano Zampini   PetscCall(DMCreateSubDM(dm, NUM_FIELDS - 1, fields + 1, &isp, &dmp));
1094811d88abSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)dm, "IS potential", (PetscObject)isp));
1095811d88abSStefano Zampini   PetscCall(DMSetMatType(dmp, MATAIJ));
1096811d88abSStefano Zampini   PetscCall(DMGetDS(dmp, &ds));
1097811d88abSStefano Zampini   //PetscCall(PetscDSSetResidual(ds, 0, P_0, P_1_aux));
1098811d88abSStefano Zampini   PetscCall(PetscDSSetResidual(ds, 0, 0, P_1_aux));
1099811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, JP_1_p1p1_aux));
1100811d88abSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(dmp, PETSC_FALSE, NULL));
1101811d88abSStefano Zampini 
1102811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dmp, &p));
1103811d88abSStefano Zampini 
1104811d88abSStefano Zampini   PetscCall(SNESCreate(PetscObjectComm((PetscObject)dmp), &snes));
1105811d88abSStefano Zampini   PetscCall(SNESSetDM(snes, dmp));
1106811d88abSStefano Zampini   PetscCall(SNESSetOptionsPrefix(snes, "initial_"));
1107811d88abSStefano Zampini   PetscCall(SNESSetErrorIfNotConverged(snes, PETSC_TRUE));
1108811d88abSStefano Zampini   PetscCall(SNESGetKSP(snes, &ksp));
1109*204aa523SStefano Zampini   PetscCall(KSPSetType(ksp, KSPFGMRES));
1110811d88abSStefano Zampini   PetscCall(KSPGetPC(ksp, &pc));
1111811d88abSStefano Zampini   PetscCall(PCSetType(pc, PCGAMG));
1112811d88abSStefano Zampini   PetscCall(SNESSetFromOptions(snes));
1113811d88abSStefano Zampini   PetscCall(SNESSetUp(snes));
1114811d88abSStefano Zampini 
1115811d88abSStefano Zampini   /* Loop over input vectors and compute corresponding potential */
1116811d88abSStefano Zampini   for (PetscInt i = 0; i < nv; i++) {
1117811d88abSStefano Zampini     PetscErrorCode (*funcs[NUM_FIELDS])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
1118811d88abSStefano Zampini 
1119811d88abSStefano Zampini     u = vecs[i];
1120811d88abSStefano Zampini     if (!valid) { /* Assumes entries in u are not valid */
1121811d88abSStefano Zampini       PetscCall(TSGetTime(ts, &t));
1122811d88abSStefano Zampini       switch (ctx->ic_num) {
1123811d88abSStefano Zampini       case 0:
1124811d88abSStefano Zampini         funcs[C_FIELD_ID] = initial_conditions_C_0;
1125811d88abSStefano Zampini         break;
1126811d88abSStefano Zampini       case 1:
1127811d88abSStefano Zampini         funcs[C_FIELD_ID] = initial_conditions_C_1;
1128811d88abSStefano Zampini         break;
1129811d88abSStefano Zampini       case 2:
1130811d88abSStefano Zampini         funcs[C_FIELD_ID] = initial_conditions_C_2;
1131811d88abSStefano Zampini         break;
1132811d88abSStefano Zampini       default:
1133d8b4a066SPierre Jolivet         SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Unknown IC");
1134811d88abSStefano Zampini       }
1135811d88abSStefano Zampini       funcs[P_FIELD_ID] = zerof;
1136811d88abSStefano Zampini       PetscCall(DMProjectFunction(dm, t, funcs, NULL, INSERT_ALL_VALUES, u));
1137811d88abSStefano Zampini     }
1138811d88abSStefano Zampini 
1139811d88abSStefano Zampini     /* pass conductivity and source information via auxiliary data */
1140811d88abSStefano Zampini     PetscCall(DMGetGlobalVector(dm, &tu));
1141811d88abSStefano Zampini     PetscCall(VecCopy(u, tu));
1142811d88abSStefano Zampini     PetscCall(VecISSet(tu, isp, 0.0));
1143811d88abSStefano Zampini     PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &lsource));
1144811d88abSStefano Zampini     PetscCall(DMCreateLocalVector(dm, &subaux));
1145811d88abSStefano Zampini     PetscCall(DMGlobalToLocal(dm, tu, INSERT_VALUES, subaux));
1146811d88abSStefano Zampini     PetscCall(DMRestoreGlobalVector(dm, &tu));
1147811d88abSStefano Zampini     PetscCall(VecAXPY(subaux, 1.0, lsource));
1148811d88abSStefano Zampini     PetscCall(VecViewFromOptions(subaux, NULL, "-initial_aux_view"));
1149811d88abSStefano Zampini     PetscCall(DMSetAuxiliaryVec(dmp, NULL, 0, 0, subaux));
1150811d88abSStefano Zampini     PetscCall(VecDestroy(&subaux));
1151811d88abSStefano Zampini 
1152811d88abSStefano Zampini     /* solve the subproblem */
1153811d88abSStefano Zampini     if (!sctp) PetscCall(VecScatterCreate(u, isp, p, NULL, &sctp));
1154811d88abSStefano Zampini     PetscCall(VecScatterBegin(sctp, u, p, INSERT_VALUES, SCATTER_FORWARD));
1155811d88abSStefano Zampini     PetscCall(VecScatterEnd(sctp, u, p, INSERT_VALUES, SCATTER_FORWARD));
1156811d88abSStefano Zampini     PetscCall(SNESSolve(snes, NULL, p));
1157811d88abSStefano Zampini 
1158811d88abSStefano Zampini     /* scatter from potential only to full space */
1159811d88abSStefano Zampini     PetscCall(VecScatterBegin(sctp, p, u, INSERT_VALUES, SCATTER_REVERSE));
1160811d88abSStefano Zampini     PetscCall(VecScatterEnd(sctp, p, u, INSERT_VALUES, SCATTER_REVERSE));
1161811d88abSStefano Zampini     PetscCall(ZeroMeanPotential(dm, u));
1162811d88abSStefano Zampini   }
1163811d88abSStefano Zampini   PetscCall(VecDestroy(&p));
1164811d88abSStefano Zampini   PetscCall(DMDestroy(&dmp));
1165811d88abSStefano Zampini   PetscCall(SNESDestroy(&snes));
1166811d88abSStefano Zampini   PetscCall(VecScatterDestroy(&sctp));
1167811d88abSStefano Zampini 
1168811d88abSStefano Zampini   /* exclude potential from computation of the LTE */
1169811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &vatol));
1170811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &vrtol));
1171811d88abSStefano Zampini   PetscCall(TSGetTolerances(ts, &atol, NULL, &rtol, NULL));
1172811d88abSStefano Zampini   PetscCall(VecSet(vatol, atol));
1173811d88abSStefano Zampini   PetscCall(VecISSet(vatol, isp, -1));
1174811d88abSStefano Zampini   PetscCall(VecSet(vrtol, rtol));
1175811d88abSStefano Zampini   PetscCall(VecISSet(vrtol, isp, -1));
1176811d88abSStefano Zampini   PetscCall(TSSetTolerances(ts, atol, vatol, rtol, vrtol));
1177811d88abSStefano Zampini   PetscCall(VecDestroy(&vatol));
1178811d88abSStefano Zampini   PetscCall(VecDestroy(&vrtol));
1179811d88abSStefano Zampini   PetscCall(ISDestroy(&isp));
1180811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1181811d88abSStefano Zampini }
1182811d88abSStefano Zampini 
1183811d88abSStefano Zampini /* mesh adaption context */
1184811d88abSStefano Zampini typedef struct {
1185811d88abSStefano Zampini   VecTagger refineTag;
1186811d88abSStefano Zampini   DMLabel   adaptLabel;
1187811d88abSStefano Zampini   PetscInt  cnt;
1188811d88abSStefano Zampini } AdaptCtx;
1189811d88abSStefano Zampini 
1190811d88abSStefano Zampini static PetscErrorCode ResizeSetUp(TS ts, PetscInt nstep, PetscReal time, Vec u, PetscBool *resize, void *vctx)
1191811d88abSStefano Zampini {
1192811d88abSStefano Zampini   AdaptCtx *ctx = (AdaptCtx *)vctx;
1193811d88abSStefano Zampini   Vec       ellVecCells, ellVecCellsF;
1194811d88abSStefano Zampini   DM        dm, plex;
1195811d88abSStefano Zampini   PetscDS   ds;
1196811d88abSStefano Zampini   PetscReal norm;
1197811d88abSStefano Zampini   PetscInt  cStart, cEnd;
1198811d88abSStefano Zampini 
1199811d88abSStefano Zampini   PetscFunctionBeginUser;
1200811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1201811d88abSStefano Zampini   PetscCall(DMConvert(dm, DMPLEX, &plex));
1202811d88abSStefano Zampini   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
1203811d88abSStefano Zampini   PetscCall(DMDestroy(&plex));
1204811d88abSStefano Zampini   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)ts), NUM_FIELDS * (cEnd - cStart), PETSC_DECIDE, &ellVecCellsF));
1205811d88abSStefano Zampini   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)ts), cEnd - cStart, PETSC_DECIDE, &ellVecCells));
1206811d88abSStefano Zampini   PetscCall(VecSetBlockSize(ellVecCellsF, NUM_FIELDS));
1207811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
1208811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, ellipticity_fail));
1209811d88abSStefano Zampini   PetscCall(DMPlexComputeCellwiseIntegralFEM(dm, u, ellVecCellsF, NULL));
1210811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, energy));
1211811d88abSStefano Zampini   PetscCall(VecStrideGather(ellVecCellsF, C_FIELD_ID, ellVecCells, INSERT_VALUES));
1212811d88abSStefano Zampini   PetscCall(VecDestroy(&ellVecCellsF));
1213811d88abSStefano Zampini   PetscCall(VecNorm(ellVecCells, NORM_1, &norm));
1214811d88abSStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "STEP %d norm %g\n", (int)nstep, (double)norm));
1215811d88abSStefano Zampini   if (norm && !ctx->cnt) {
1216811d88abSStefano Zampini     IS refineIS;
1217811d88abSStefano Zampini 
1218811d88abSStefano Zampini     *resize = PETSC_TRUE;
1219811d88abSStefano Zampini     if (!ctx->refineTag) {
1220811d88abSStefano Zampini       VecTaggerBox refineBox;
1221811d88abSStefano Zampini       refineBox.min = PETSC_MACHINE_EPSILON;
1222811d88abSStefano Zampini       refineBox.max = PETSC_MAX_REAL;
1223811d88abSStefano Zampini 
1224811d88abSStefano Zampini       PetscCall(VecTaggerCreate(PETSC_COMM_SELF, &ctx->refineTag));
1225811d88abSStefano Zampini       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->refineTag, "refine_"));
1226811d88abSStefano Zampini       PetscCall(VecTaggerSetType(ctx->refineTag, VECTAGGERABSOLUTE));
1227811d88abSStefano Zampini       PetscCall(VecTaggerAbsoluteSetBox(ctx->refineTag, &refineBox));
1228811d88abSStefano Zampini       PetscCall(VecTaggerSetFromOptions(ctx->refineTag));
1229811d88abSStefano Zampini       PetscCall(VecTaggerSetUp(ctx->refineTag));
1230811d88abSStefano Zampini       PetscCall(PetscObjectViewFromOptions((PetscObject)ctx->refineTag, NULL, "-tag_view"));
1231811d88abSStefano Zampini     }
1232811d88abSStefano Zampini     PetscCall(DMLabelDestroy(&ctx->adaptLabel));
1233811d88abSStefano Zampini     PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)ts), "adapt", &ctx->adaptLabel));
1234811d88abSStefano Zampini     PetscCall(VecTaggerComputeIS(ctx->refineTag, ellVecCells, &refineIS, NULL));
1235811d88abSStefano Zampini     PetscCall(DMLabelSetStratumIS(ctx->adaptLabel, DM_ADAPT_REFINE, refineIS));
1236811d88abSStefano Zampini     PetscCall(ISDestroy(&refineIS));
1237811d88abSStefano Zampini #if 0
1238811d88abSStefano Zampini     void (*funcs[NUM_FIELDS])(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 f[]);
1239811d88abSStefano Zampini     Vec ellVec;
1240811d88abSStefano Zampini 
1241811d88abSStefano Zampini     funcs[P_FIELD_ID] = ellipticity_fail;
1242811d88abSStefano Zampini     funcs[C_FIELD_ID] = NULL;
1243811d88abSStefano Zampini 
1244811d88abSStefano Zampini     PetscCall(DMGetGlobalVector(dm, &ellVec));
1245811d88abSStefano Zampini     PetscCall(DMProjectField(dm, 0, u, funcs, INSERT_VALUES, ellVec));
1246811d88abSStefano Zampini     PetscCall(VecViewFromOptions(ellVec,NULL,"-view_amr_ell"));
1247811d88abSStefano Zampini     PetscCall(DMRestoreGlobalVector(dm, &ellVec));
1248811d88abSStefano Zampini #endif
1249811d88abSStefano Zampini     ctx->cnt++;
1250811d88abSStefano Zampini   } else {
1251811d88abSStefano Zampini     ctx->cnt = 0;
1252811d88abSStefano Zampini   }
1253811d88abSStefano Zampini   PetscCall(VecDestroy(&ellVecCells));
1254811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1255811d88abSStefano Zampini }
1256811d88abSStefano Zampini 
1257811d88abSStefano Zampini static PetscErrorCode ResizeTransfer(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], void *vctx)
1258811d88abSStefano Zampini {
1259811d88abSStefano Zampini   AdaptCtx *actx = (AdaptCtx *)vctx;
1260811d88abSStefano Zampini   AppCtx   *ctx;
1261811d88abSStefano Zampini   DM        dm, adm;
1262811d88abSStefano Zampini   PetscReal time;
1263811d88abSStefano Zampini 
1264811d88abSStefano Zampini   PetscFunctionBeginUser;
1265811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1266811d88abSStefano Zampini   PetscCheck(actx->adaptLabel, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Missing adaptLabel");
1267811d88abSStefano Zampini   PetscCall(DMAdaptLabel(dm, actx->adaptLabel, &adm));
1268811d88abSStefano Zampini   PetscCheck(adm, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Missing adapted DM");
1269811d88abSStefano Zampini   PetscCall(TSGetTime(ts, &time));
1270811d88abSStefano Zampini   PetscCall(DMLabelDestroy(&actx->adaptLabel));
1271811d88abSStefano Zampini   for (PetscInt i = 0; i < nv; i++) {
1272811d88abSStefano Zampini     PetscCall(DMCreateGlobalVector(adm, &vecsout[i]));
1273811d88abSStefano Zampini     PetscCall(DMForestTransferVec(dm, vecsin[i], adm, vecsout[i], PETSC_TRUE, time));
1274811d88abSStefano Zampini   }
1275811d88abSStefano Zampini   PetscCall(DMForestSetAdaptivityForest(adm, NULL));
1276811d88abSStefano Zampini   PetscCall(DMSetCoarseDM(adm, NULL));
1277811d88abSStefano Zampini   PetscCall(DMSetLocalSection(adm, NULL));
1278811d88abSStefano Zampini   PetscCall(TSSetDM(ts, adm));
1279811d88abSStefano Zampini   PetscCall(TSGetTime(ts, &time));
1280811d88abSStefano Zampini   PetscCall(TSGetApplicationContext(ts, &ctx));
1281*204aa523SStefano Zampini   PetscCall(DMSetNullSpaceConstructor(adm, P_FIELD_ID, CreatePotentialNullSpace));
1282811d88abSStefano Zampini   PetscCall(ProjectSource(adm, time, ctx));
1283811d88abSStefano Zampini   PetscCall(SetInitialConditionsAndTolerances(ts, nv, vecsout, PETSC_TRUE));
1284811d88abSStefano Zampini   PetscCall(DMDestroy(&adm));
1285811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1286811d88abSStefano Zampini }
1287811d88abSStefano Zampini 
1288811d88abSStefano Zampini static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1289811d88abSStefano Zampini {
1290811d88abSStefano Zampini   PetscFunctionBeginUser;
1291811d88abSStefano Zampini   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer));
1292811d88abSStefano Zampini   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERVTK));
1293811d88abSStefano Zampini   PetscCall(PetscViewerFileSetName(*viewer, filename));
1294811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1295811d88abSStefano Zampini }
1296811d88abSStefano Zampini 
1297811d88abSStefano Zampini /* Monitor relevant functionals */
1298811d88abSStefano Zampini static PetscErrorCode Monitor(TS ts, PetscInt stepnum, PetscReal time, Vec u, void *vctx)
1299811d88abSStefano Zampini {
1300811d88abSStefano Zampini   PetscScalar vals[2 * NUM_FIELDS];
1301811d88abSStefano Zampini   DM          dm;
1302811d88abSStefano Zampini   PetscDS     ds;
1303811d88abSStefano Zampini   AppCtx     *ctx;
1304811d88abSStefano Zampini 
1305811d88abSStefano Zampini   PetscFunctionBeginUser;
1306811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1307811d88abSStefano Zampini   PetscCall(TSGetApplicationContext(ts, &ctx));
1308811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
1309811d88abSStefano Zampini 
1310811d88abSStefano Zampini   /* monitor energy and potential average */
1311811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, average));
1312811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
1313811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
1314811d88abSStefano Zampini 
1315811d88abSStefano Zampini   /* monitor ellipticity_fail */
1316811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, ellipticity_fail));
1317811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals + NUM_FIELDS, NULL));
1318811d88abSStefano Zampini   if (ctx->ellipticity) {
1319811d88abSStefano Zampini     void (*funcs[NUM_FIELDS])(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 f[]);
1320811d88abSStefano Zampini     Vec         ellVec;
1321811d88abSStefano Zampini     PetscViewer viewer;
1322811d88abSStefano Zampini     char        filename[PETSC_MAX_PATH_LEN];
1323811d88abSStefano Zampini 
1324811d88abSStefano Zampini     funcs[P_FIELD_ID] = ellipticity_fail;
1325811d88abSStefano Zampini     funcs[C_FIELD_ID] = NULL;
1326811d88abSStefano Zampini 
1327811d88abSStefano Zampini     PetscCall(DMGetGlobalVector(dm, &ellVec));
1328811d88abSStefano Zampini     PetscCall(DMProjectField(dm, 0, u, funcs, INSERT_VALUES, ellVec));
1329811d88abSStefano Zampini     PetscCall(PetscSNPrintf(filename, sizeof filename, "ellipticity_fail-%03" PetscInt_FMT ".vtu", stepnum));
1330811d88abSStefano Zampini     PetscCall(OutputVTK(dm, filename, &viewer));
1331811d88abSStefano Zampini     PetscCall(VecView(ellVec, viewer));
1332811d88abSStefano Zampini     PetscCall(PetscViewerDestroy(&viewer));
1333811d88abSStefano Zampini     PetscCall(DMRestoreGlobalVector(dm, &ellVec));
1334811d88abSStefano Zampini   }
1335811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, energy));
1336811d88abSStefano Zampini 
1337811d88abSStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "%4" PetscInt_FMT " TS: time %g, energy %g, intp %g, ell %g\n", stepnum, (double)time, (double)PetscRealPart(vals[C_FIELD_ID]), (double)PetscRealPart(vals[P_FIELD_ID]), (double)PetscRealPart(vals[NUM_FIELDS + C_FIELD_ID])));
1338811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1339811d88abSStefano Zampini }
1340811d88abSStefano Zampini 
1341811d88abSStefano Zampini /* Save restart information */
1342811d88abSStefano Zampini static PetscErrorCode MonitorSave(TS ts, PetscInt steps, PetscReal time, Vec u, void *vctx)
1343811d88abSStefano Zampini {
1344811d88abSStefano Zampini   DM                dm;
1345811d88abSStefano Zampini   AppCtx           *ctx        = (AppCtx *)vctx;
1346811d88abSStefano Zampini   PetscInt          save_every = ctx->save_every;
1347811d88abSStefano Zampini   TSConvergedReason reason;
1348811d88abSStefano Zampini 
1349811d88abSStefano Zampini   PetscFunctionBeginUser;
1350811d88abSStefano Zampini   if (!ctx->save) PetscFunctionReturn(PETSC_SUCCESS);
1351811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1352811d88abSStefano Zampini   PetscCall(TSGetConvergedReason(ts, &reason));
1353811d88abSStefano Zampini   if ((save_every > 0 && steps % save_every == 0) || (save_every == -1 && reason) || save_every < -1) PetscCall(SaveToFile(dm, u, ctx->save_filename));
1354811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1355811d88abSStefano Zampini }
1356811d88abSStefano Zampini 
1357811d88abSStefano Zampini /* Make potential zero mean after SNES solve */
1358811d88abSStefano Zampini static PetscErrorCode PostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
1359811d88abSStefano Zampini {
1360811d88abSStefano Zampini   DM       dm;
1361811d88abSStefano Zampini   Vec      u = Y[stageindex];
1362*204aa523SStefano Zampini   SNES     snes;
1363*204aa523SStefano Zampini   PetscInt nits, lits, stepnum;
1364*204aa523SStefano Zampini   AppCtx  *ctx;
1365811d88abSStefano Zampini 
1366811d88abSStefano Zampini   PetscFunctionBeginUser;
1367811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1368811d88abSStefano Zampini   PetscCall(ZeroMeanPotential(dm, u));
1369*204aa523SStefano Zampini 
1370*204aa523SStefano Zampini   PetscCall(TSGetApplicationContext(ts, &ctx));
1371*204aa523SStefano Zampini   if (ctx->test_restart) PetscFunctionReturn(PETSC_SUCCESS);
1372*204aa523SStefano Zampini 
1373*204aa523SStefano Zampini   /* monitor linear and nonlinear iterations */
1374*204aa523SStefano Zampini   PetscCall(TSGetStepNumber(ts, &stepnum));
1375*204aa523SStefano Zampini   PetscCall(TSGetSNES(ts, &snes));
1376*204aa523SStefano Zampini   PetscCall(SNESGetIterationNumber(snes, &nits));
1377*204aa523SStefano Zampini   PetscCall(SNESGetLinearSolveIterations(snes, &lits));
1378*204aa523SStefano Zampini 
1379*204aa523SStefano Zampini   /* if function evals in TSDIRK are zero in the first stage, it is FSAL */
1380*204aa523SStefano Zampini   if (stageindex == 0) {
1381*204aa523SStefano Zampini     PetscBool dirk;
1382*204aa523SStefano Zampini     PetscInt  nf;
1383*204aa523SStefano Zampini 
1384*204aa523SStefano Zampini     PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSDIRK, &dirk));
1385*204aa523SStefano Zampini     PetscCall(SNESGetNumberFunctionEvals(snes, &nf));
1386*204aa523SStefano Zampini     if (dirk && nf == 0) nits = 0;
1387*204aa523SStefano Zampini   }
1388*204aa523SStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "         step %" PetscInt_FMT " stage %" PetscInt_FMT " nonlinear its %" PetscInt_FMT ", linear its %" PetscInt_FMT "\n", stepnum, stageindex, nits, lits));
1389811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1390811d88abSStefano Zampini }
1391811d88abSStefano Zampini 
1392811d88abSStefano Zampini static PetscErrorCode VecViewFlux(Vec u, const char *opts)
1393811d88abSStefano Zampini {
1394811d88abSStefano Zampini   Vec       fluxVec;
1395811d88abSStefano Zampini   DM        dmFlux, dm, plex;
1396811d88abSStefano Zampini   PetscInt  dim;
1397811d88abSStefano Zampini   PetscFE   feC, feFluxC, feNormC;
1398811d88abSStefano Zampini   PetscBool simplex, has;
1399811d88abSStefano Zampini 
1400811d88abSStefano Zampini   void (*funcs[])(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 f[]) = {normc, flux};
1401811d88abSStefano Zampini 
1402811d88abSStefano Zampini   PetscFunctionBeginUser;
1403811d88abSStefano Zampini   PetscCall(PetscOptionsHasName(NULL, NULL, opts, &has));
1404811d88abSStefano Zampini   if (!has) PetscFunctionReturn(PETSC_SUCCESS);
1405811d88abSStefano Zampini   PetscCall(VecGetDM(u, &dm));
1406811d88abSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
1407811d88abSStefano Zampini   PetscCall(DMGetField(dm, C_FIELD_ID, NULL, (PetscObject *)&feC));
1408811d88abSStefano Zampini   PetscCall(DMConvert(dm, DMPLEX, &plex));
1409811d88abSStefano Zampini   PetscCall(DMPlexIsSimplex(plex, &simplex));
1410811d88abSStefano Zampini   PetscCall(DMDestroy(&plex));
1411811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, simplex, "flux_", -1, &feFluxC));
1412811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, "normc_", -1, &feNormC));
1413811d88abSStefano Zampini   PetscCall(PetscFECopyQuadrature(feC, feFluxC));
1414811d88abSStefano Zampini   PetscCall(PetscFECopyQuadrature(feC, feNormC));
1415811d88abSStefano Zampini   PetscCall(DMClone(dm, &dmFlux));
1416811d88abSStefano Zampini   PetscCall(DMSetNumFields(dmFlux, 1));
1417811d88abSStefano Zampini   PetscCall(DMSetField(dmFlux, 0, NULL, (PetscObject)feNormC));
1418811d88abSStefano Zampini   /* paraview segfaults! */
1419811d88abSStefano Zampini   //PetscCall(DMSetField(dmFlux, 1, NULL, (PetscObject)feFluxC));
1420811d88abSStefano Zampini   PetscCall(DMCreateDS(dmFlux));
1421811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feFluxC));
1422811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feNormC));
1423811d88abSStefano Zampini 
1424811d88abSStefano Zampini   PetscCall(DMGetGlobalVector(dmFlux, &fluxVec));
1425811d88abSStefano Zampini   PetscCall(DMProjectField(dmFlux, 0.0, u, funcs, INSERT_VALUES, fluxVec));
1426811d88abSStefano Zampini   PetscCall(VecViewFromOptions(fluxVec, NULL, opts));
1427811d88abSStefano Zampini   PetscCall(DMRestoreGlobalVector(dmFlux, &fluxVec));
1428811d88abSStefano Zampini   PetscCall(DMDestroy(&dmFlux));
1429811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1430811d88abSStefano Zampini }
1431811d88abSStefano Zampini 
1432811d88abSStefano Zampini static PetscErrorCode Run(MPI_Comm comm, AppCtx *ctx)
1433811d88abSStefano Zampini {
1434811d88abSStefano Zampini   DM        dm;
1435811d88abSStefano Zampini   TS        ts;
1436811d88abSStefano Zampini   Vec       u;
1437811d88abSStefano Zampini   AdaptCtx *actx;
1438811d88abSStefano Zampini   PetscBool flg;
1439811d88abSStefano Zampini 
1440811d88abSStefano Zampini   PetscFunctionBeginUser;
1441811d88abSStefano Zampini   if (ctx->test_restart) {
1442811d88abSStefano Zampini     PetscViewer subviewer;
1443811d88abSStefano Zampini     PetscMPIInt rank;
1444811d88abSStefano Zampini 
1445811d88abSStefano Zampini     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1446811d88abSStefano Zampini     PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1447811d88abSStefano Zampini     if (ctx->load) PetscCall(PetscViewerASCIIPrintf(subviewer, "rank %d loading from %s\n", rank, ctx->load_filename));
1448811d88abSStefano Zampini     if (ctx->save) PetscCall(PetscViewerASCIIPrintf(subviewer, "rank %d saving to %s\n", rank, ctx->save_filename));
1449811d88abSStefano Zampini     PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1450811d88abSStefano Zampini     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1451811d88abSStefano Zampini   } else {
1452811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "----------------------------\n"));
1453811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "Simulation parameters:\n"));
1454811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  r    : %g\n", (double)ctx->r));
1455811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  eps  : %g\n", (double)ctx->eps));
1456811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  alpha: %g\n", (double)ctx->alpha));
1457811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  gamma: %g\n", (double)ctx->gamma));
1458811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  D    : %g\n", (double)ctx->D));
1459811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  c    : %g\n", (double)ctx->c));
1460811d88abSStefano Zampini     if (ctx->load) PetscCall(PetscPrintf(comm, "  load : %s\n", ctx->load_filename));
1461811d88abSStefano Zampini     else PetscCall(PetscPrintf(comm, "  IC   : %" PetscInt_FMT "\n", ctx->ic_num));
1462811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  S    : %" PetscInt_FMT "\n", ctx->source_num));
1463811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  x0   : (%g,%g)\n", (double)ctx->x0[0], (double)ctx->x0[1]));
1464811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "----------------------------\n"));
1465811d88abSStefano Zampini   }
1466811d88abSStefano Zampini 
1467811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePush(SetupStage));
1468811d88abSStefano Zampini   PetscCall(CreateMesh(comm, &dm, ctx));
1469811d88abSStefano Zampini   PetscCall(SetupDiscretization(dm, ctx));
1470811d88abSStefano Zampini 
1471811d88abSStefano Zampini   PetscCall(TSCreate(comm, &ts));
1472811d88abSStefano Zampini   PetscCall(TSSetApplicationContext(ts, ctx));
1473811d88abSStefano Zampini 
1474811d88abSStefano Zampini   PetscCall(TSSetDM(ts, dm));
1475811d88abSStefano Zampini   if (ctx->test_restart) {
1476811d88abSStefano Zampini     PetscViewer subviewer;
1477811d88abSStefano Zampini     PetscMPIInt rank;
1478811d88abSStefano Zampini     PetscInt    level;
1479811d88abSStefano Zampini 
1480811d88abSStefano Zampini     PetscCall(DMGetRefineLevel(dm, &level));
1481811d88abSStefano Zampini     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1482811d88abSStefano Zampini     PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1483811d88abSStefano Zampini     PetscCall(PetscViewerASCIIPrintf(subviewer, "rank %d DM refinement level %" PetscInt_FMT "\n", rank, level));
1484811d88abSStefano Zampini     PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1485811d88abSStefano Zampini     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1486811d88abSStefano Zampini   }
1487811d88abSStefano Zampini 
1488811d88abSStefano Zampini   if (ctx->test_restart) PetscCall(TSSetMaxSteps(ts, 1));
1489811d88abSStefano Zampini   PetscCall(TSSetMaxTime(ts, 10.0));
1490811d88abSStefano Zampini   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1491811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(TSMonitorSet(ts, Monitor, NULL, NULL));
1492811d88abSStefano Zampini   PetscCall(TSMonitorSet(ts, MonitorSave, ctx, NULL));
1493811d88abSStefano Zampini   PetscCall(PetscNew(&actx));
1494811d88abSStefano Zampini   if (ctx->amr) PetscCall(TSSetResize(ts, PETSC_TRUE, ResizeSetUp, ResizeTransfer, actx));
1495811d88abSStefano Zampini   PetscCall(TSSetPostStage(ts, PostStage));
1496811d88abSStefano Zampini   PetscCall(TSSetMaxSNESFailures(ts, -1));
1497811d88abSStefano Zampini   PetscCall(TSSetFromOptions(ts));
1498811d88abSStefano Zampini 
1499811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &u));
1500811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)u, "solution_"));
1501811d88abSStefano Zampini   PetscCall(DMHasNamedGlobalVector(dm, "solution_", &flg));
1502811d88abSStefano Zampini   if (flg) {
1503811d88abSStefano Zampini     Vec ru;
1504811d88abSStefano Zampini 
1505811d88abSStefano Zampini     PetscCall(DMGetNamedGlobalVector(dm, "solution_", &ru));
1506811d88abSStefano Zampini     PetscCall(VecCopy(ru, u));
1507811d88abSStefano Zampini     PetscCall(DMRestoreNamedGlobalVector(dm, "solution_", &ru));
1508811d88abSStefano Zampini   }
1509811d88abSStefano Zampini   PetscCall(SetInitialConditionsAndTolerances(ts, 1, &u, PETSC_FALSE));
1510811d88abSStefano Zampini   PetscCall(TSSetSolution(ts, u));
1511811d88abSStefano Zampini   PetscCall(VecDestroy(&u));
1512811d88abSStefano Zampini   PetscCall(DMDestroy(&dm));
1513811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePop());
1514811d88abSStefano Zampini 
1515811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePush(SolveStage));
1516811d88abSStefano Zampini   PetscCall(TSSolve(ts, NULL));
1517811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePop());
1518811d88abSStefano Zampini 
1519811d88abSStefano Zampini   PetscCall(TSGetSolution(ts, &u));
1520811d88abSStefano Zampini   PetscCall(VecViewFromOptions(u, NULL, "-final_view"));
1521811d88abSStefano Zampini   PetscCall(VecViewFlux(u, "-final_flux_view"));
1522811d88abSStefano Zampini 
1523811d88abSStefano Zampini   PetscCall(TSDestroy(&ts));
1524811d88abSStefano Zampini   PetscCall(VecTaggerDestroy(&actx->refineTag));
1525811d88abSStefano Zampini   PetscCall(PetscFree(actx));
1526811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1527811d88abSStefano Zampini }
1528811d88abSStefano Zampini 
1529811d88abSStefano Zampini int main(int argc, char **argv)
1530811d88abSStefano Zampini {
1531811d88abSStefano Zampini   AppCtx ctx;
1532811d88abSStefano Zampini 
1533811d88abSStefano Zampini   PetscFunctionBeginUser;
1534811d88abSStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1535811d88abSStefano Zampini   PetscCall(ProcessOptions(&ctx));
1536811d88abSStefano Zampini   PetscCall(PetscLogStageRegister("Setup", &SetupStage));
1537811d88abSStefano Zampini   PetscCall(PetscLogStageRegister("Solve", &SolveStage));
1538811d88abSStefano Zampini   if (ctx.test_restart) { /* Test sequences of save and loads */
1539811d88abSStefano Zampini     PetscMPIInt rank;
1540811d88abSStefano Zampini 
1541811d88abSStefano Zampini     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1542811d88abSStefano Zampini 
1543811d88abSStefano Zampini     /* test saving */
1544811d88abSStefano Zampini     ctx.load = PETSC_FALSE;
1545811d88abSStefano Zampini     ctx.save = PETSC_TRUE;
1546811d88abSStefano Zampini     /* sequential save */
1547811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test sequential save\n"));
1548811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.save_filename, sizeof(ctx.save_filename), "test_ex30_seq_%d.h5", rank));
1549811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_SELF, &ctx));
1550811d88abSStefano Zampini     /* parallel save */
1551811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test parallel save\n"));
1552811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.save_filename, sizeof(ctx.save_filename), "test_ex30_par.h5"));
1553811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1554811d88abSStefano Zampini 
1555811d88abSStefano Zampini     /* test loading */
1556811d88abSStefano Zampini     ctx.load = PETSC_TRUE;
1557811d88abSStefano Zampini     ctx.save = PETSC_FALSE;
1558811d88abSStefano Zampini     /* sequential and parallel runs from sequential save */
1559811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.load_filename, sizeof(ctx.load_filename), "test_ex30_seq_0.h5"));
1560811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test sequential load from sequential save\n"));
1561811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_SELF, &ctx));
1562811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test parallel load from sequential save\n"));
1563811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1564811d88abSStefano Zampini     /* sequential and parallel runs from parallel save */
1565811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.load_filename, sizeof(ctx.load_filename), "test_ex30_par.h5"));
1566811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test sequential load from parallel save\n"));
1567811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_SELF, &ctx));
1568811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test parallel load from parallel save\n"));
1569811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1570811d88abSStefano Zampini   } else { /* Run the simulation */
1571811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1572811d88abSStefano Zampini   }
1573811d88abSStefano Zampini   PetscCall(PetscFinalize());
1574811d88abSStefano Zampini   return 0;
1575811d88abSStefano Zampini }
1576811d88abSStefano Zampini 
1577811d88abSStefano Zampini /*TEST
1578811d88abSStefano Zampini 
1579811d88abSStefano Zampini   testset:
1580811d88abSStefano Zampini     args: -dm_plex_box_faces 3,3 -ksp_type preonly -pc_type svd -c_petscspace_degree 1 -p_petscspace_degree 1 -ts_max_steps 1 -initial_snes_test_jacobian -snes_test_jacobian -initial_snes_type ksponly -snes_type ksponly -petscpartitioner_type simple -dm_plex_simplex 0 -ts_adapt_type none
1581811d88abSStefano Zampini 
1582811d88abSStefano Zampini     test:
1583811d88abSStefano Zampini       suffix: 0
1584811d88abSStefano Zampini       nsize: {{1 2}}
1585*204aa523SStefano Zampini       args: -dm_refine 1 -lump {{0 1}}
1586811d88abSStefano Zampini 
1587811d88abSStefano Zampini     test:
1588811d88abSStefano Zampini       suffix: 0_dirk
1589811d88abSStefano Zampini       nsize: {{1 2}}
1590811d88abSStefano Zampini       args: -dm_refine 1 -ts_type dirk
1591811d88abSStefano Zampini 
1592811d88abSStefano Zampini     test:
1593811d88abSStefano Zampini       suffix: 0_dirk_mg
1594811d88abSStefano Zampini       nsize: {{1 2}}
1595*204aa523SStefano Zampini       args: -dm_refine_hierarchy 1 -ts_type dirk -pc_type mg -mg_levels_pc_type bjacobi -mg_levels_sub_pc_factor_levels 2 -mg_levels_sub_pc_factor_mat_ordering_type rcm -mg_levels_sub_pc_factor_reuse_ordering -mg_coarse_pc_type svd -lump {{0 1}}
1596*204aa523SStefano Zampini 
1597*204aa523SStefano Zampini     test:
1598*204aa523SStefano Zampini       suffix: 0_dirk_fieldsplit
1599*204aa523SStefano Zampini       nsize: {{1 2}}
1600*204aa523SStefano Zampini       args: -dm_refine 1 -ts_type dirk -pc_type fieldsplit -pc_fieldsplit_type multiplicative -lump {{0 1}}
1601811d88abSStefano Zampini 
1602811d88abSStefano Zampini     test:
1603811d88abSStefano Zampini       requires: p4est
1604811d88abSStefano Zampini       suffix: 0_p4est
1605811d88abSStefano Zampini       nsize: {{1 2}}
1606*204aa523SStefano Zampini       args: -dm_refine 1 -dm_plex_convert_type p4est -lump 0
1607811d88abSStefano Zampini 
1608811d88abSStefano Zampini     test:
1609811d88abSStefano Zampini       suffix: 0_periodic
1610811d88abSStefano Zampini       nsize: {{1 2}}
1611*204aa523SStefano Zampini       args: -dm_plex_box_bd periodic,periodic -dm_refine_pre 1 -lump {{0 1}}
1612811d88abSStefano Zampini 
1613811d88abSStefano Zampini     test:
1614811d88abSStefano Zampini       requires: p4est
1615811d88abSStefano Zampini       suffix: 0_p4est_periodic
1616811d88abSStefano Zampini       nsize: {{1 2}}
1617*204aa523SStefano Zampini       args: -dm_plex_box_bd periodic,periodic -dm_refine_pre 1 -dm_plex_convert_type p4est -lump 0
1618811d88abSStefano Zampini 
1619811d88abSStefano Zampini     test:
1620811d88abSStefano Zampini       requires: p4est
1621811d88abSStefano Zampini       suffix: 0_p4est_mg
1622811d88abSStefano Zampini       nsize: {{1 2}}
1623*204aa523SStefano Zampini       args: -dm_forest_minimum_refinement 0 -dm_forest_initial_refinement 2 -dm_plex_convert_type p4est -pc_type mg -mg_coarse_pc_type svd -mg_levels_pc_type svd -lump 0
1624811d88abSStefano Zampini 
1625811d88abSStefano Zampini   testset:
1626811d88abSStefano Zampini     requires: hdf5
1627811d88abSStefano Zampini     args: -test_restart -dm_plex_box_faces 3,3 -ksp_type preonly -pc_type mg -mg_levels_pc_type svd -c_petscspace_degree 1 -p_petscspace_degree 1 -petscpartitioner_type simple -test_restart
1628811d88abSStefano Zampini 
1629811d88abSStefano Zampini     test:
1630910b42cbSStefano Zampini       requires: !single
1631811d88abSStefano Zampini       suffix: restart
1632811d88abSStefano Zampini       nsize: {{1 2}separate output}
1633811d88abSStefano Zampini       args: -dm_refine_hierarchy {{0 1}separate output} -dm_plex_simplex 0
1634811d88abSStefano Zampini 
1635811d88abSStefano Zampini     test:
1636811d88abSStefano Zampini       requires: triangle
1637811d88abSStefano Zampini       suffix: restart_simplex
1638811d88abSStefano Zampini       nsize: {{1 2}separate output}
1639811d88abSStefano Zampini       args: -dm_refine_hierarchy {{0 1}separate output} -dm_plex_simplex 1
1640811d88abSStefano Zampini 
1641811d88abSStefano Zampini     test:
1642910b42cbSStefano Zampini       requires: !single
1643811d88abSStefano Zampini       suffix: restart_refonly
1644811d88abSStefano Zampini       nsize: {{1 2}separate output}
1645811d88abSStefano Zampini       args: -dm_refine 1 -dm_plex_simplex 0
1646811d88abSStefano Zampini 
1647811d88abSStefano Zampini     test:
1648811d88abSStefano Zampini       requires: triangle
1649811d88abSStefano Zampini       suffix: restart_simplex_refonly
1650811d88abSStefano Zampini       nsize: {{1 2}separate output}
1651811d88abSStefano Zampini       args: -dm_refine 1 -dm_plex_simplex 1
1652811d88abSStefano Zampini 
1653811d88abSStefano Zampini TEST*/
1654