xref: /petsc/src/ts/tutorials/ex30.c (revision 811d88ab63ac5868f4a3c762a0ed1431b7553237)
1*811d88abSStefano Zampini static char help[] = "Biological network from https://link.springer.com/article/10.1007/s42967-023-00297-3\n\n\n";
2*811d88abSStefano Zampini 
3*811d88abSStefano Zampini #include <petscts.h>
4*811d88abSStefano Zampini #include <petscsf.h>
5*811d88abSStefano Zampini #include <petscdmplex.h>
6*811d88abSStefano Zampini #include <petscdmplextransform.h>
7*811d88abSStefano Zampini #include <petscdmforest.h>
8*811d88abSStefano Zampini #include <petscviewerhdf5.h>
9*811d88abSStefano Zampini #include <petscds.h>
10*811d88abSStefano Zampini 
11*811d88abSStefano Zampini /*
12*811d88abSStefano Zampini     Here we solve the system of PDEs on \Omega \in R^2:
13*811d88abSStefano Zampini 
14*811d88abSStefano Zampini     * dC/dt - D^2 \Delta C - c^2 \nabla p \cross \nabla p + \alpha sqrt(||C||^2_F + eps)^(\gamma-2) C = 0
15*811d88abSStefano Zampini     * - \nabla \cdot ((r + C) \nabla p) = S
16*811d88abSStefano Zampini 
17*811d88abSStefano Zampini     where:
18*811d88abSStefano Zampini       C = symmetric 2x2 conductivity tensor
19*811d88abSStefano Zampini       p = potential
20*811d88abSStefano Zampini       S = source
21*811d88abSStefano Zampini 
22*811d88abSStefano Zampini     with natural boundary conditions on \partial\Omega:
23*811d88abSStefano Zampini       \nabla C \cdot n  = 0
24*811d88abSStefano Zampini       \nabla ((r + C)\nabla p) \cdot n  = 0
25*811d88abSStefano Zampini 
26*811d88abSStefano Zampini     Parameters:
27*811d88abSStefano Zampini       D = diffusion constant
28*811d88abSStefano Zampini       c = activation parameter
29*811d88abSStefano Zampini       \alpha = metabolic coefficient
30*811d88abSStefano Zampini       \gamma = metabolic exponent
31*811d88abSStefano Zampini       r, eps are regularization parameters
32*811d88abSStefano Zampini 
33*811d88abSStefano Zampini     We use Lagrange elements for C_ij and P.
34*811d88abSStefano Zampini */
35*811d88abSStefano Zampini 
36*811d88abSStefano Zampini typedef enum _fieldidx {
37*811d88abSStefano Zampini   C_FIELD_ID = 0,
38*811d88abSStefano Zampini   P_FIELD_ID,
39*811d88abSStefano Zampini   NUM_FIELDS
40*811d88abSStefano Zampini } FieldIdx;
41*811d88abSStefano Zampini 
42*811d88abSStefano Zampini typedef enum _constantidx {
43*811d88abSStefano Zampini   R_ID = 0,
44*811d88abSStefano Zampini   EPS_ID,
45*811d88abSStefano Zampini   ALPHA_ID,
46*811d88abSStefano Zampini   GAMMA_ID,
47*811d88abSStefano Zampini   D_ID,
48*811d88abSStefano Zampini   C2_ID,
49*811d88abSStefano Zampini   NUM_CONSTANTS
50*811d88abSStefano Zampini } ConstantIdx;
51*811d88abSStefano Zampini 
52*811d88abSStefano Zampini PetscLogStage SetupStage, SolveStage;
53*811d88abSStefano Zampini 
54*811d88abSStefano Zampini #define NORM2C(c00, c01, c11) PetscSqr(c00) + 2 * PetscSqr(c01) + PetscSqr(c11)
55*811d88abSStefano Zampini 
56*811d88abSStefano Zampini /* residual for C when tested against basis functions */
57*811d88abSStefano 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[])
58*811d88abSStefano Zampini {
59*811d88abSStefano Zampini   const PetscReal   c2       = PetscRealPart(constants[C2_ID]);
60*811d88abSStefano Zampini   const PetscReal   alpha    = PetscRealPart(constants[ALPHA_ID]);
61*811d88abSStefano Zampini   const PetscReal   gamma    = PetscRealPart(constants[GAMMA_ID]);
62*811d88abSStefano Zampini   const PetscReal   eps      = PetscRealPart(constants[EPS_ID]);
63*811d88abSStefano Zampini   const PetscScalar gradp[]  = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
64*811d88abSStefano Zampini   const PetscScalar crossp[] = {gradp[0] * gradp[0], gradp[0] * gradp[1], gradp[1] * gradp[1]};
65*811d88abSStefano Zampini   const PetscScalar C00      = u[uOff[C_FIELD_ID]];
66*811d88abSStefano Zampini   const PetscScalar C01      = u[uOff[C_FIELD_ID] + 1];
67*811d88abSStefano Zampini   const PetscScalar C11      = u[uOff[C_FIELD_ID] + 2];
68*811d88abSStefano Zampini   const PetscScalar norm     = NORM2C(C00, C01, C11) + eps;
69*811d88abSStefano Zampini   const PetscScalar nexp     = (gamma - 2.0) / 2.0;
70*811d88abSStefano Zampini   const PetscScalar fnorm    = PetscPowScalar(norm, nexp);
71*811d88abSStefano Zampini 
72*811d88abSStefano 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];
73*811d88abSStefano Zampini }
74*811d88abSStefano Zampini 
75*811d88abSStefano Zampini /* Jacobian for C against C basis functions */
76*811d88abSStefano 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[])
77*811d88abSStefano Zampini {
78*811d88abSStefano Zampini   const PetscReal   alpha  = PetscRealPart(constants[ALPHA_ID]);
79*811d88abSStefano Zampini   const PetscReal   gamma  = PetscRealPart(constants[GAMMA_ID]);
80*811d88abSStefano Zampini   const PetscReal   eps    = PetscRealPart(constants[EPS_ID]);
81*811d88abSStefano Zampini   const PetscScalar C00    = u[uOff[C_FIELD_ID]];
82*811d88abSStefano Zampini   const PetscScalar C01    = u[uOff[C_FIELD_ID] + 1];
83*811d88abSStefano Zampini   const PetscScalar C11    = u[uOff[C_FIELD_ID] + 2];
84*811d88abSStefano Zampini   const PetscScalar norm   = NORM2C(C00, C01, C11) + eps;
85*811d88abSStefano Zampini   const PetscScalar nexp   = (gamma - 2.0) / 2.0;
86*811d88abSStefano Zampini   const PetscScalar fnorm  = PetscPowScalar(norm, nexp);
87*811d88abSStefano Zampini   const PetscScalar dfnorm = nexp * PetscPowScalar(norm, nexp - 1.0);
88*811d88abSStefano Zampini   const PetscScalar dC[]   = {2 * C00, 4 * C01, 2 * C11};
89*811d88abSStefano Zampini 
90*811d88abSStefano Zampini   for (PetscInt k = 0; k < 3; k++) {
91*811d88abSStefano Zampini     for (PetscInt j = 0; j < 3; j++) J[k * 3 + j] = alpha * dfnorm * dC[j] * u[uOff[C_FIELD_ID] + k];
92*811d88abSStefano Zampini     J[k * 3 + k] += alpha * fnorm + u_tShift;
93*811d88abSStefano Zampini   }
94*811d88abSStefano Zampini }
95*811d88abSStefano Zampini 
96*811d88abSStefano Zampini /* Jacobian for C against C basis functions and gradients of P basis functions */
97*811d88abSStefano 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[])
98*811d88abSStefano Zampini {
99*811d88abSStefano Zampini   const PetscReal   c2      = PetscRealPart(constants[C2_ID]);
100*811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
101*811d88abSStefano Zampini 
102*811d88abSStefano Zampini   J[0] = -c2 * 2 * gradp[0];
103*811d88abSStefano Zampini   J[1] = 0.0;
104*811d88abSStefano Zampini   J[2] = -c2 * gradp[1];
105*811d88abSStefano Zampini   J[3] = -c2 * gradp[0];
106*811d88abSStefano Zampini   J[4] = 0.0;
107*811d88abSStefano Zampini   J[5] = -c2 * 2 * gradp[1];
108*811d88abSStefano Zampini }
109*811d88abSStefano Zampini 
110*811d88abSStefano Zampini /* residual for C when tested against gradients of basis functions */
111*811d88abSStefano 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[])
112*811d88abSStefano Zampini {
113*811d88abSStefano Zampini   const PetscReal D = PetscRealPart(constants[D_ID]);
114*811d88abSStefano Zampini   for (PetscInt k = 0; k < 3; k++)
115*811d88abSStefano Zampini     for (PetscInt d = 0; d < 2; d++) f1[k * 2 + d] = PetscSqr(D) * u_x[uOff_x[C_FIELD_ID] + k * 2 + d];
116*811d88abSStefano Zampini }
117*811d88abSStefano Zampini 
118*811d88abSStefano Zampini /* Jacobian for C against gradients of C basis functions */
119*811d88abSStefano 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[])
120*811d88abSStefano Zampini {
121*811d88abSStefano Zampini   const PetscReal D = PetscRealPart(constants[D_ID]);
122*811d88abSStefano Zampini   for (PetscInt k = 0; k < 3; k++)
123*811d88abSStefano Zampini     for (PetscInt d = 0; d < 2; d++) J[k * (3 + 1) * 2 * 2 + d * 2 + d] = PetscSqr(D);
124*811d88abSStefano Zampini }
125*811d88abSStefano Zampini 
126*811d88abSStefano Zampini /* residual for P when tested against basis functions.
127*811d88abSStefano Zampini    The source term always comes from the auxiliary vec because it needs to have zero mean */
128*811d88abSStefano 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[])
129*811d88abSStefano Zampini {
130*811d88abSStefano Zampini   PetscScalar S = a[aOff[P_FIELD_ID]];
131*811d88abSStefano Zampini 
132*811d88abSStefano Zampini   f0[0] = -S;
133*811d88abSStefano Zampini }
134*811d88abSStefano Zampini 
135*811d88abSStefano Zampini /* residual for P when tested against gradients of basis functions */
136*811d88abSStefano 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[])
137*811d88abSStefano Zampini {
138*811d88abSStefano Zampini   const PetscReal   r       = PetscRealPart(constants[R_ID]);
139*811d88abSStefano Zampini   const PetscScalar C00     = u[uOff[C_FIELD_ID]];
140*811d88abSStefano Zampini   const PetscScalar C01     = u[uOff[C_FIELD_ID] + 1];
141*811d88abSStefano Zampini   const PetscScalar C10     = C01;
142*811d88abSStefano Zampini   const PetscScalar C11     = u[uOff[C_FIELD_ID] + 2];
143*811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
144*811d88abSStefano Zampini 
145*811d88abSStefano Zampini   f1[0] = (C00 + r) * gradp[0] + C01 * gradp[1];
146*811d88abSStefano Zampini   f1[1] = C10 * gradp[0] + (C11 + r) * gradp[1];
147*811d88abSStefano Zampini }
148*811d88abSStefano Zampini 
149*811d88abSStefano Zampini /* Same as above for the P-only subproblem for initial conditions: the conductivity values come from the auxiliary vec */
150*811d88abSStefano 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[])
151*811d88abSStefano Zampini {
152*811d88abSStefano Zampini   const PetscReal   r       = PetscRealPart(constants[R_ID]);
153*811d88abSStefano Zampini   const PetscScalar C00     = a[aOff[C_FIELD_ID]];
154*811d88abSStefano Zampini   const PetscScalar C01     = a[aOff[C_FIELD_ID] + 1];
155*811d88abSStefano Zampini   const PetscScalar C10     = C01;
156*811d88abSStefano Zampini   const PetscScalar C11     = a[aOff[C_FIELD_ID] + 2];
157*811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[0]], u_x[uOff_x[0] + 1]};
158*811d88abSStefano Zampini 
159*811d88abSStefano Zampini   f1[0] = (C00 + r) * gradp[0] + C01 * gradp[1];
160*811d88abSStefano Zampini   f1[1] = C10 * gradp[0] + (C11 + r) * gradp[1];
161*811d88abSStefano Zampini }
162*811d88abSStefano Zampini 
163*811d88abSStefano Zampini /* Jacobian for P against gradients of P basis functions */
164*811d88abSStefano 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[])
165*811d88abSStefano Zampini {
166*811d88abSStefano Zampini   const PetscReal   r   = PetscRealPart(constants[R_ID]);
167*811d88abSStefano Zampini   const PetscScalar C00 = u[uOff[C_FIELD_ID]];
168*811d88abSStefano Zampini   const PetscScalar C01 = u[uOff[C_FIELD_ID] + 1];
169*811d88abSStefano Zampini   const PetscScalar C10 = C01;
170*811d88abSStefano Zampini   const PetscScalar C11 = u[uOff[C_FIELD_ID] + 2];
171*811d88abSStefano Zampini 
172*811d88abSStefano Zampini   J[0] = C00 + r;
173*811d88abSStefano Zampini   J[1] = C01;
174*811d88abSStefano Zampini   J[2] = C10;
175*811d88abSStefano Zampini   J[3] = C11 + r;
176*811d88abSStefano Zampini }
177*811d88abSStefano Zampini 
178*811d88abSStefano Zampini /* Same as above for the P-only subproblem for initial conditions */
179*811d88abSStefano 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[])
180*811d88abSStefano Zampini {
181*811d88abSStefano Zampini   const PetscReal   r   = PetscRealPart(constants[R_ID]);
182*811d88abSStefano Zampini   const PetscScalar C00 = a[aOff[C_FIELD_ID]];
183*811d88abSStefano Zampini   const PetscScalar C01 = a[aOff[C_FIELD_ID] + 1];
184*811d88abSStefano Zampini   const PetscScalar C10 = C01;
185*811d88abSStefano Zampini   const PetscScalar C11 = a[aOff[C_FIELD_ID] + 2];
186*811d88abSStefano Zampini 
187*811d88abSStefano Zampini   J[0] = C00 + r;
188*811d88abSStefano Zampini   J[1] = C01;
189*811d88abSStefano Zampini   J[2] = C10;
190*811d88abSStefano Zampini   J[3] = C11 + r;
191*811d88abSStefano Zampini }
192*811d88abSStefano Zampini 
193*811d88abSStefano Zampini /* Jacobian for P against gradients of P basis functions and C basis functions */
194*811d88abSStefano 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[])
195*811d88abSStefano Zampini {
196*811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
197*811d88abSStefano Zampini 
198*811d88abSStefano Zampini   J[0] = gradp[0];
199*811d88abSStefano Zampini   J[1] = 0;
200*811d88abSStefano Zampini   J[2] = gradp[1];
201*811d88abSStefano Zampini   J[3] = gradp[0];
202*811d88abSStefano Zampini   J[4] = 0;
203*811d88abSStefano Zampini   J[5] = gradp[1];
204*811d88abSStefano Zampini }
205*811d88abSStefano Zampini 
206*811d88abSStefano Zampini /* the source term S(x) = exp(-500*||x - x0||^2) */
207*811d88abSStefano Zampini static PetscErrorCode source_0(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
208*811d88abSStefano Zampini {
209*811d88abSStefano Zampini   PetscReal *x0 = (PetscReal *)ctx;
210*811d88abSStefano Zampini   PetscReal  n  = 0;
211*811d88abSStefano Zampini 
212*811d88abSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) n += (x[d] - x0[d]) * (x[d] - x0[d]);
213*811d88abSStefano Zampini   u[0] = PetscExpReal(-500 * n);
214*811d88abSStefano Zampini   return PETSC_SUCCESS;
215*811d88abSStefano Zampini }
216*811d88abSStefano Zampini 
217*811d88abSStefano Zampini static PetscErrorCode source_1(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
218*811d88abSStefano Zampini {
219*811d88abSStefano Zampini   PetscScalar     ut[1];
220*811d88abSStefano Zampini   const PetscReal x0[] = {0.25, 0.25};
221*811d88abSStefano Zampini   const PetscReal x1[] = {0.75, 0.75};
222*811d88abSStefano Zampini 
223*811d88abSStefano Zampini   PetscCall(source_0(dim, time, x, Nf, ut, (void *)x0));
224*811d88abSStefano Zampini   PetscCall(source_0(dim, time, x, Nf, u, (void *)x1));
225*811d88abSStefano Zampini   u[0] += ut[0];
226*811d88abSStefano Zampini   return PETSC_SUCCESS;
227*811d88abSStefano Zampini }
228*811d88abSStefano Zampini 
229*811d88abSStefano Zampini /* functionals to be integrated: average -> \int_\Omega u dx */
230*811d88abSStefano 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[])
231*811d88abSStefano Zampini {
232*811d88abSStefano Zampini   obj[0] = u[uOff[P_FIELD_ID]];
233*811d88abSStefano Zampini }
234*811d88abSStefano Zampini 
235*811d88abSStefano Zampini /* stable implementation of roots of a*x^2 + b*x + c = 0 */
236*811d88abSStefano Zampini static inline void QuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal x[2])
237*811d88abSStefano Zampini {
238*811d88abSStefano Zampini   PetscReal delta = PetscMax(b * b - 4 * a * c, 0); /* eigenvalues symmetric matrix */
239*811d88abSStefano Zampini   PetscReal temp  = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(delta));
240*811d88abSStefano Zampini 
241*811d88abSStefano Zampini   x[0] = temp / a;
242*811d88abSStefano Zampini   x[1] = c / temp;
243*811d88abSStefano Zampini }
244*811d88abSStefano Zampini 
245*811d88abSStefano Zampini /* functionals to be integrated: energy -> D^2/2 * ||\nabla C||^2 + c^2\nabla p * (r + C) * \nabla p + \alpha/ \gamma * ||C||^\gamma */
246*811d88abSStefano 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[])
247*811d88abSStefano Zampini {
248*811d88abSStefano Zampini   const PetscReal   D         = PetscRealPart(constants[D_ID]);
249*811d88abSStefano Zampini   const PetscReal   c2        = PetscRealPart(constants[C2_ID]);
250*811d88abSStefano Zampini   const PetscReal   r         = PetscRealPart(constants[R_ID]);
251*811d88abSStefano Zampini   const PetscReal   alpha     = PetscRealPart(constants[ALPHA_ID]);
252*811d88abSStefano Zampini   const PetscReal   gamma     = PetscRealPart(constants[GAMMA_ID]);
253*811d88abSStefano Zampini   const PetscScalar C00       = u[uOff[C_FIELD_ID]];
254*811d88abSStefano Zampini   const PetscScalar C01       = u[uOff[C_FIELD_ID] + 1];
255*811d88abSStefano Zampini   const PetscScalar C10       = C01;
256*811d88abSStefano Zampini   const PetscScalar C11       = u[uOff[C_FIELD_ID] + 2];
257*811d88abSStefano Zampini   const PetscScalar gradp[]   = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
258*811d88abSStefano Zampini   const PetscScalar gradC00[] = {u_x[uOff_x[C_FIELD_ID] + 0], u_x[uOff_x[C_FIELD_ID] + 1]};
259*811d88abSStefano Zampini   const PetscScalar gradC01[] = {u_x[uOff_x[C_FIELD_ID] + 2], u_x[uOff_x[C_FIELD_ID] + 3]};
260*811d88abSStefano Zampini   const PetscScalar gradC11[] = {u_x[uOff_x[C_FIELD_ID] + 4], u_x[uOff_x[C_FIELD_ID] + 5]};
261*811d88abSStefano Zampini   const PetscScalar normC     = NORM2C(C00, C01, C11);
262*811d88abSStefano Zampini   const PetscScalar normgradC = NORM2C(gradC00[0], gradC01[0], gradC11[0]) + NORM2C(gradC00[1], gradC01[1], gradC11[1]);
263*811d88abSStefano Zampini   const PetscScalar nexp      = gamma / 2.0;
264*811d88abSStefano Zampini 
265*811d88abSStefano Zampini   const PetscScalar t0 = PetscSqr(D) / 2.0 * normgradC;
266*811d88abSStefano Zampini   const PetscScalar t1 = c2 * (gradp[0] * ((C00 + r) * gradp[0] + C01 * gradp[1]) + gradp[1] * (C10 * gradp[0] + (C11 + r) * gradp[1]));
267*811d88abSStefano Zampini   const PetscScalar t2 = alpha / gamma * PetscPowScalar(normC, nexp);
268*811d88abSStefano Zampini 
269*811d88abSStefano Zampini   obj[0] = t0 + t1 + t2;
270*811d88abSStefano Zampini }
271*811d88abSStefano Zampini 
272*811d88abSStefano 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 */
273*811d88abSStefano 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[])
274*811d88abSStefano Zampini {
275*811d88abSStefano Zampini   const PetscReal r   = PetscRealPart(constants[R_ID]);
276*811d88abSStefano Zampini   const PetscReal C00 = PetscRealPart(u[uOff[C_FIELD_ID]] + r);
277*811d88abSStefano Zampini   const PetscReal C01 = PetscRealPart(u[uOff[C_FIELD_ID] + 1]);
278*811d88abSStefano Zampini   const PetscReal C11 = PetscRealPart(u[uOff[C_FIELD_ID] + 2] + r);
279*811d88abSStefano Zampini 
280*811d88abSStefano Zampini   PetscReal eigs[2];
281*811d88abSStefano Zampini   QuadraticRoots(1, -(C00 + C11), C00 * C11 - PetscSqr(C01), eigs);
282*811d88abSStefano Zampini   if (eigs[0] < 0 || eigs[1] < 0) obj[0] = -PetscMin(eigs[0], eigs[1]);
283*811d88abSStefano Zampini   else obj[0] = 0.0;
284*811d88abSStefano Zampini }
285*811d88abSStefano Zampini 
286*811d88abSStefano Zampini /* initial conditions for C: eq. 16 */
287*811d88abSStefano Zampini static PetscErrorCode initial_conditions_C_0(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
288*811d88abSStefano Zampini {
289*811d88abSStefano Zampini   u[0] = 1;
290*811d88abSStefano Zampini   u[1] = 0;
291*811d88abSStefano Zampini   u[2] = 1;
292*811d88abSStefano Zampini   return PETSC_SUCCESS;
293*811d88abSStefano Zampini }
294*811d88abSStefano Zampini 
295*811d88abSStefano Zampini /* initial conditions for C: eq. 17 */
296*811d88abSStefano Zampini static PetscErrorCode initial_conditions_C_1(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
297*811d88abSStefano Zampini {
298*811d88abSStefano Zampini   const PetscReal x = xx[0];
299*811d88abSStefano Zampini   const PetscReal y = xx[1];
300*811d88abSStefano Zampini 
301*811d88abSStefano Zampini   u[0] = (2 - PetscAbsReal(x + y)) * PetscExpReal(-10 * PetscAbsReal(x - y));
302*811d88abSStefano Zampini   u[1] = 0;
303*811d88abSStefano Zampini   u[2] = (2 - PetscAbsReal(x + y)) * PetscExpReal(-10 * PetscAbsReal(x - y));
304*811d88abSStefano Zampini   return PETSC_SUCCESS;
305*811d88abSStefano Zampini }
306*811d88abSStefano Zampini 
307*811d88abSStefano Zampini /* initial conditions for C: eq. 18 */
308*811d88abSStefano Zampini static PetscErrorCode initial_conditions_C_2(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
309*811d88abSStefano Zampini {
310*811d88abSStefano Zampini   u[0] = 0;
311*811d88abSStefano Zampini   u[1] = 0;
312*811d88abSStefano Zampini   u[2] = 0;
313*811d88abSStefano Zampini   return PETSC_SUCCESS;
314*811d88abSStefano Zampini }
315*811d88abSStefano Zampini 
316*811d88abSStefano Zampini /* functionals to be sampled: C * \grad p */
317*811d88abSStefano 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[])
318*811d88abSStefano Zampini {
319*811d88abSStefano Zampini   const PetscScalar C00     = u[uOff[C_FIELD_ID]];
320*811d88abSStefano Zampini   const PetscScalar C01     = u[uOff[C_FIELD_ID] + 1];
321*811d88abSStefano Zampini   const PetscScalar C10     = C01;
322*811d88abSStefano Zampini   const PetscScalar C11     = u[uOff[C_FIELD_ID] + 2];
323*811d88abSStefano Zampini   const PetscScalar gradp[] = {u_x[uOff_x[P_FIELD_ID]], u_x[uOff_x[P_FIELD_ID] + 1]};
324*811d88abSStefano Zampini 
325*811d88abSStefano Zampini   f[0] = C00 * gradp[0] + C01 * gradp[1];
326*811d88abSStefano Zampini   f[1] = C10 * gradp[0] + C11 * gradp[1];
327*811d88abSStefano Zampini }
328*811d88abSStefano Zampini 
329*811d88abSStefano Zampini /* functionals to be sampled: ||C|| */
330*811d88abSStefano 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[])
331*811d88abSStefano Zampini {
332*811d88abSStefano Zampini   const PetscScalar C00 = u[uOff[C_FIELD_ID]];
333*811d88abSStefano Zampini   const PetscScalar C01 = u[uOff[C_FIELD_ID] + 1];
334*811d88abSStefano Zampini   const PetscScalar C11 = u[uOff[C_FIELD_ID] + 2];
335*811d88abSStefano Zampini 
336*811d88abSStefano Zampini   f[0] = PetscSqrtScalar(NORM2C(C00, C01, C11));
337*811d88abSStefano Zampini }
338*811d88abSStefano Zampini 
339*811d88abSStefano Zampini /* functionals to be sampled: zero */
340*811d88abSStefano 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[])
341*811d88abSStefano Zampini {
342*811d88abSStefano Zampini   f[0] = 0.0;
343*811d88abSStefano Zampini }
344*811d88abSStefano Zampini 
345*811d88abSStefano Zampini /* functions to be sampled: zero function */
346*811d88abSStefano Zampini static PetscErrorCode zerof(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
347*811d88abSStefano Zampini {
348*811d88abSStefano Zampini   for (PetscInt d = 0; d < Nc; ++d) u[d] = 0.0;
349*811d88abSStefano Zampini   return PETSC_SUCCESS;
350*811d88abSStefano Zampini }
351*811d88abSStefano Zampini 
352*811d88abSStefano Zampini /* functions to be sampled: constant function */
353*811d88abSStefano Zampini static PetscErrorCode constantf(PetscInt dim, PetscReal time, const PetscReal xx[], PetscInt Nc, PetscScalar *u, void *ctx)
354*811d88abSStefano Zampini {
355*811d88abSStefano Zampini   PetscInt d;
356*811d88abSStefano Zampini   for (d = 0; d < Nc; ++d) u[d] = 1.0;
357*811d88abSStefano Zampini   return PETSC_SUCCESS;
358*811d88abSStefano Zampini }
359*811d88abSStefano Zampini 
360*811d88abSStefano Zampini /* application context: customizable parameters */
361*811d88abSStefano Zampini typedef struct {
362*811d88abSStefano Zampini   PetscReal r;
363*811d88abSStefano Zampini   PetscReal eps;
364*811d88abSStefano Zampini   PetscReal alpha;
365*811d88abSStefano Zampini   PetscReal gamma;
366*811d88abSStefano Zampini   PetscReal D;
367*811d88abSStefano Zampini   PetscReal c;
368*811d88abSStefano Zampini   PetscInt  ic_num;
369*811d88abSStefano Zampini   PetscInt  source_num;
370*811d88abSStefano Zampini   PetscReal x0[2];
371*811d88abSStefano Zampini   PetscBool amr;
372*811d88abSStefano Zampini   PetscBool load;
373*811d88abSStefano Zampini   char      load_filename[PETSC_MAX_PATH_LEN];
374*811d88abSStefano Zampini   PetscBool save;
375*811d88abSStefano Zampini   char      save_filename[PETSC_MAX_PATH_LEN];
376*811d88abSStefano Zampini   PetscInt  save_every;
377*811d88abSStefano Zampini   PetscBool test_restart;
378*811d88abSStefano Zampini   PetscBool ellipticity;
379*811d88abSStefano Zampini } AppCtx;
380*811d88abSStefano Zampini 
381*811d88abSStefano Zampini /* process command line options */
382*811d88abSStefano Zampini static PetscErrorCode ProcessOptions(AppCtx *options)
383*811d88abSStefano Zampini {
384*811d88abSStefano Zampini   PetscInt dim = PETSC_STATIC_ARRAY_LENGTH(options->x0);
385*811d88abSStefano Zampini 
386*811d88abSStefano Zampini   PetscFunctionBeginUser;
387*811d88abSStefano Zampini   options->r            = 1.e-1;
388*811d88abSStefano Zampini   options->eps          = 1.e-3;
389*811d88abSStefano Zampini   options->alpha        = 0.75;
390*811d88abSStefano Zampini   options->gamma        = 0.75;
391*811d88abSStefano Zampini   options->c            = 5;
392*811d88abSStefano Zampini   options->D            = 1.e-2;
393*811d88abSStefano Zampini   options->ic_num       = 0;
394*811d88abSStefano Zampini   options->source_num   = 0;
395*811d88abSStefano Zampini   options->x0[0]        = 0.25;
396*811d88abSStefano Zampini   options->x0[1]        = 0.25;
397*811d88abSStefano Zampini   options->amr          = PETSC_FALSE;
398*811d88abSStefano Zampini   options->load         = PETSC_FALSE;
399*811d88abSStefano Zampini   options->save         = PETSC_FALSE;
400*811d88abSStefano Zampini   options->save_every   = -1;
401*811d88abSStefano Zampini   options->test_restart = PETSC_FALSE;
402*811d88abSStefano Zampini   options->ellipticity  = PETSC_FALSE;
403*811d88abSStefano Zampini 
404*811d88abSStefano Zampini   PetscOptionsBegin(PETSC_COMM_WORLD, "", __FILE__, "DMPLEX");
405*811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-alpha", "alpha", __FILE__, options->alpha, &options->alpha, NULL));
406*811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-gamma", "gamma", __FILE__, options->gamma, &options->gamma, NULL));
407*811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-c", "c", __FILE__, options->c, &options->c, NULL));
408*811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-d", "D", __FILE__, options->D, &options->D, NULL));
409*811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-eps", "eps", __FILE__, options->eps, &options->eps, NULL));
410*811d88abSStefano Zampini   PetscCall(PetscOptionsReal("-r", "r", __FILE__, options->r, &options->r, NULL));
411*811d88abSStefano Zampini   PetscCall(PetscOptionsRealArray("-x0", "x0", __FILE__, options->x0, &dim, NULL));
412*811d88abSStefano Zampini   PetscCall(PetscOptionsInt("-ic_num", "ic_num", __FILE__, options->ic_num, &options->ic_num, NULL));
413*811d88abSStefano Zampini   PetscCall(PetscOptionsInt("-source_num", "source_num", __FILE__, options->source_num, &options->source_num, NULL));
414*811d88abSStefano Zampini   PetscCall(PetscOptionsBool("-amr", "use adaptive mesh refinement", __FILE__, options->amr, &options->amr, NULL));
415*811d88abSStefano Zampini   PetscCall(PetscOptionsBool("-test_restart", "test restarting files", __FILE__, options->test_restart, &options->test_restart, NULL));
416*811d88abSStefano Zampini   if (!options->test_restart) {
417*811d88abSStefano 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));
418*811d88abSStefano 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));
419*811d88abSStefano 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));
420*811d88abSStefano Zampini   }
421*811d88abSStefano Zampini   PetscCall(PetscOptionsBool("-monitor_ellipticity", "Dump locations of ellipticity violation", __FILE__, options->ellipticity, &options->ellipticity, NULL));
422*811d88abSStefano Zampini   PetscOptionsEnd();
423*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
424*811d88abSStefano Zampini }
425*811d88abSStefano Zampini 
426*811d88abSStefano Zampini static PetscErrorCode SaveToFile(DM dm, Vec u, const char *filename)
427*811d88abSStefano Zampini {
428*811d88abSStefano Zampini #if defined(PETSC_HAVE_HDF5)
429*811d88abSStefano Zampini   PetscViewerFormat format = PETSC_VIEWER_HDF5_PETSC;
430*811d88abSStefano Zampini   PetscViewer       viewer;
431*811d88abSStefano Zampini   DM                cdm       = dm;
432*811d88abSStefano Zampini   PetscInt          numlevels = 0;
433*811d88abSStefano Zampini 
434*811d88abSStefano Zampini   PetscFunctionBeginUser;
435*811d88abSStefano Zampini   while (cdm) {
436*811d88abSStefano Zampini     numlevels++;
437*811d88abSStefano Zampini     PetscCall(DMGetCoarseDM(cdm, &cdm));
438*811d88abSStefano Zampini   }
439*811d88abSStefano Zampini   /* Cannot be set programmatically */
440*811d88abSStefano Zampini   PetscCall(PetscOptionsInsertString(NULL, "-dm_plex_view_hdf5_storage_version 3.0.0"));
441*811d88abSStefano Zampini   PetscCall(PetscViewerHDF5Open(PetscObjectComm((PetscObject)dm), filename, FILE_MODE_WRITE, &viewer));
442*811d88abSStefano Zampini   PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "numlevels", PETSC_INT, &numlevels));
443*811d88abSStefano Zampini   PetscCall(PetscViewerPushFormat(viewer, format));
444*811d88abSStefano Zampini   for (PetscInt level = numlevels - 1; level >= 0; level--) {
445*811d88abSStefano Zampini     PetscInt    cc, rr;
446*811d88abSStefano Zampini     PetscBool   isRegular, isUniform;
447*811d88abSStefano Zampini     const char *dmname;
448*811d88abSStefano Zampini     char        groupname[PETSC_MAX_PATH_LEN];
449*811d88abSStefano Zampini 
450*811d88abSStefano Zampini     PetscCall(PetscSNPrintf(groupname, sizeof(groupname), "level_%" PetscInt_FMT, level));
451*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PushGroup(viewer, groupname));
452*811d88abSStefano Zampini     PetscCall(PetscObjectGetName((PetscObject)dm, &dmname));
453*811d88abSStefano Zampini     PetscCall(DMGetCoarsenLevel(dm, &cc));
454*811d88abSStefano Zampini     PetscCall(DMGetRefineLevel(dm, &rr));
455*811d88abSStefano Zampini     PetscCall(DMPlexGetRegularRefinement(dm, &isRegular));
456*811d88abSStefano Zampini     PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
457*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "meshname", PETSC_STRING, dmname));
458*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "refinelevel", PETSC_INT, &rr));
459*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "coarsenlevel", PETSC_INT, &cc));
460*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "refRegular", PETSC_BOOL, &isRegular));
461*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5WriteAttribute(viewer, NULL, "refUniform", PETSC_BOOL, &isUniform));
462*811d88abSStefano Zampini     PetscCall(DMPlexTopologyView(dm, viewer));
463*811d88abSStefano Zampini     PetscCall(DMPlexLabelsView(dm, viewer));
464*811d88abSStefano Zampini     PetscCall(DMPlexCoordinatesView(dm, viewer));
465*811d88abSStefano Zampini     PetscCall(DMPlexSectionView(dm, viewer, NULL));
466*811d88abSStefano Zampini     if (level == numlevels - 1) {
467*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)u, "solution_"));
468*811d88abSStefano Zampini       PetscCall(DMPlexGlobalVectorView(dm, viewer, NULL, u));
469*811d88abSStefano Zampini     }
470*811d88abSStefano Zampini     if (level) {
471*811d88abSStefano Zampini       PetscInt        cStart, cEnd, ccStart, ccEnd, cpStart;
472*811d88abSStefano Zampini       DMPolytopeType  ct;
473*811d88abSStefano Zampini       DMPlexTransform tr;
474*811d88abSStefano Zampini       DM              sdm;
475*811d88abSStefano Zampini       PetscScalar    *array;
476*811d88abSStefano Zampini       PetscSection    section;
477*811d88abSStefano Zampini       Vec             map;
478*811d88abSStefano Zampini       IS              gis;
479*811d88abSStefano Zampini       const PetscInt *gidx;
480*811d88abSStefano Zampini 
481*811d88abSStefano Zampini       PetscCall(DMGetCoarseDM(dm, &cdm));
482*811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
483*811d88abSStefano Zampini       PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
484*811d88abSStefano Zampini       PetscCall(PetscSectionSetChart(section, cStart, cEnd));
485*811d88abSStefano Zampini       for (PetscInt c = cStart; c < cEnd; c++) PetscCall(PetscSectionSetDof(section, c, 1));
486*811d88abSStefano Zampini       PetscCall(PetscSectionSetUp(section));
487*811d88abSStefano Zampini 
488*811d88abSStefano Zampini       PetscCall(DMClone(dm, &sdm));
489*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)sdm, "pdm"));
490*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)section, "pdm_section"));
491*811d88abSStefano Zampini       PetscCall(DMSetLocalSection(sdm, section));
492*811d88abSStefano Zampini       PetscCall(PetscSectionDestroy(&section));
493*811d88abSStefano Zampini 
494*811d88abSStefano Zampini       PetscCall(DMGetLocalVector(sdm, &map));
495*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)map, "pdm_map"));
496*811d88abSStefano Zampini       PetscCall(VecGetArray(map, &array));
497*811d88abSStefano Zampini       PetscCall(DMPlexTransformCreate(PETSC_COMM_SELF, &tr));
498*811d88abSStefano Zampini       PetscCall(DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR));
499*811d88abSStefano Zampini       PetscCall(DMPlexTransformSetDM(tr, cdm));
500*811d88abSStefano Zampini       PetscCall(DMPlexTransformSetFromOptions(tr));
501*811d88abSStefano Zampini       PetscCall(DMPlexTransformSetUp(tr));
502*811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(cdm, 0, &ccStart, &ccEnd));
503*811d88abSStefano Zampini       PetscCall(DMPlexGetChart(cdm, &cpStart, NULL));
504*811d88abSStefano Zampini       PetscCall(DMPlexCreatePointNumbering(cdm, &gis));
505*811d88abSStefano Zampini       PetscCall(ISGetIndices(gis, &gidx));
506*811d88abSStefano Zampini       for (PetscInt c = ccStart; c < ccEnd; c++) {
507*811d88abSStefano Zampini         PetscInt       *rsize, *rcone, *rornt, Nt;
508*811d88abSStefano Zampini         DMPolytopeType *rct;
509*811d88abSStefano Zampini         PetscInt        gnum = gidx[c - cpStart] >= 0 ? gidx[c - cpStart] : -(gidx[c - cpStart] + 1);
510*811d88abSStefano Zampini 
511*811d88abSStefano Zampini         PetscCall(DMPlexGetCellType(cdm, c, &ct));
512*811d88abSStefano Zampini         PetscCall(DMPlexTransformCellTransform(tr, ct, c, NULL, &Nt, &rct, &rsize, &rcone, &rornt));
513*811d88abSStefano Zampini         for (PetscInt r = 0; r < rsize[Nt - 1]; ++r) {
514*811d88abSStefano Zampini           PetscInt pNew;
515*811d88abSStefano Zampini 
516*811d88abSStefano Zampini           PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[Nt - 1], c, r, &pNew));
517*811d88abSStefano Zampini           array[pNew - cStart] = gnum;
518*811d88abSStefano Zampini         }
519*811d88abSStefano Zampini       }
520*811d88abSStefano Zampini       PetscCall(ISRestoreIndices(gis, &gidx));
521*811d88abSStefano Zampini       PetscCall(ISDestroy(&gis));
522*811d88abSStefano Zampini       PetscCall(VecRestoreArray(map, &array));
523*811d88abSStefano Zampini       PetscCall(DMPlexTransformDestroy(&tr));
524*811d88abSStefano Zampini       PetscCall(DMPlexSectionView(dm, viewer, sdm));
525*811d88abSStefano Zampini       PetscCall(DMPlexLocalVectorView(dm, viewer, sdm, map));
526*811d88abSStefano Zampini       PetscCall(DMRestoreLocalVector(sdm, &map));
527*811d88abSStefano Zampini       PetscCall(DMDestroy(&sdm));
528*811d88abSStefano Zampini     }
529*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PopGroup(viewer));
530*811d88abSStefano Zampini     PetscCall(DMGetCoarseDM(dm, &dm));
531*811d88abSStefano Zampini   }
532*811d88abSStefano Zampini   PetscCall(PetscViewerPopFormat(viewer));
533*811d88abSStefano Zampini   PetscCall(PetscViewerDestroy(&viewer));
534*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
535*811d88abSStefano Zampini #else
536*811d88abSStefano Zampini   SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Needs HDF5 support. Please reconfigure using --download-hdf5");
537*811d88abSStefano Zampini #endif
538*811d88abSStefano Zampini }
539*811d88abSStefano Zampini 
540*811d88abSStefano Zampini static PetscErrorCode LoadFromFile(MPI_Comm comm, const char *filename, DM *odm)
541*811d88abSStefano Zampini {
542*811d88abSStefano Zampini #if defined(PETSC_HAVE_HDF5)
543*811d88abSStefano Zampini   PetscViewerFormat format = PETSC_VIEWER_HDF5_PETSC;
544*811d88abSStefano Zampini   PetscViewer       viewer;
545*811d88abSStefano Zampini   DM                dm, cdm = NULL;
546*811d88abSStefano Zampini   PetscSF           sfXC      = NULL;
547*811d88abSStefano Zampini   PetscInt          numlevels = -1;
548*811d88abSStefano Zampini 
549*811d88abSStefano Zampini   PetscFunctionBeginUser;
550*811d88abSStefano Zampini   PetscCall(PetscViewerHDF5Open(comm, filename, FILE_MODE_READ, &viewer));
551*811d88abSStefano Zampini   PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "numlevels", PETSC_INT, NULL, &numlevels));
552*811d88abSStefano Zampini   PetscCall(PetscViewerPushFormat(viewer, format));
553*811d88abSStefano Zampini   for (PetscInt level = 0; level < numlevels; level++) {
554*811d88abSStefano Zampini     char             groupname[PETSC_MAX_PATH_LEN], *dmname;
555*811d88abSStefano Zampini     PetscSF          sfXB, sfBC, sfG;
556*811d88abSStefano Zampini     PetscPartitioner part;
557*811d88abSStefano Zampini     PetscInt         rr, cc;
558*811d88abSStefano Zampini     PetscBool        isRegular, isUniform;
559*811d88abSStefano Zampini 
560*811d88abSStefano Zampini     PetscCall(DMCreate(comm, &dm));
561*811d88abSStefano Zampini     PetscCall(DMSetType(dm, DMPLEX));
562*811d88abSStefano Zampini     PetscCall(PetscSNPrintf(groupname, sizeof(groupname), "level_%" PetscInt_FMT, level));
563*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PushGroup(viewer, groupname));
564*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "meshname", PETSC_STRING, NULL, &dmname));
565*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "refinelevel", PETSC_INT, NULL, &rr));
566*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "coarsenlevel", PETSC_INT, NULL, &cc));
567*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "refRegular", PETSC_BOOL, NULL, &isRegular));
568*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5ReadAttribute(viewer, NULL, "refUniform", PETSC_BOOL, NULL, &isUniform));
569*811d88abSStefano Zampini     PetscCall(PetscObjectSetName((PetscObject)dm, dmname));
570*811d88abSStefano Zampini     PetscCall(DMPlexTopologyLoad(dm, viewer, &sfXB));
571*811d88abSStefano Zampini     PetscCall(DMPlexLabelsLoad(dm, viewer, sfXB));
572*811d88abSStefano Zampini     PetscCall(DMPlexCoordinatesLoad(dm, viewer, sfXB));
573*811d88abSStefano Zampini     PetscCall(DMPlexGetPartitioner(dm, &part));
574*811d88abSStefano Zampini     if (!level) { /* partition the coarse level only */
575*811d88abSStefano Zampini       PetscCall(PetscPartitionerSetFromOptions(part));
576*811d88abSStefano Zampini     } else { /* propagate partitioning information from coarser to finer level */
577*811d88abSStefano Zampini       DM           sdm;
578*811d88abSStefano Zampini       Vec          map;
579*811d88abSStefano Zampini       PetscSF      sf;
580*811d88abSStefano Zampini       PetscLayout  clayout;
581*811d88abSStefano Zampini       PetscScalar *array;
582*811d88abSStefano Zampini       PetscInt    *cranks_leaf, *cranks_root, *npoints, *points, *ranks, *starts, *gidxs;
583*811d88abSStefano Zampini       PetscInt     nparts, cStart, cEnd, nr, ccStart, ccEnd, cpStart, cpEnd;
584*811d88abSStefano Zampini       PetscMPIInt  size, rank;
585*811d88abSStefano Zampini 
586*811d88abSStefano Zampini       PetscCall(DMClone(dm, &sdm));
587*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)sdm, "pdm"));
588*811d88abSStefano Zampini       PetscCall(DMPlexSectionLoad(dm, viewer, sdm, sfXB, NULL, &sf));
589*811d88abSStefano Zampini       PetscCall(DMGetLocalVector(sdm, &map));
590*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)map, "pdm_map"));
591*811d88abSStefano Zampini       PetscCall(DMPlexLocalVectorLoad(dm, viewer, sdm, sf, map));
592*811d88abSStefano Zampini 
593*811d88abSStefano Zampini       PetscCallMPI(MPI_Comm_size(comm, &size));
594*811d88abSStefano Zampini       PetscCallMPI(MPI_Comm_rank(comm, &rank));
595*811d88abSStefano Zampini       nparts = size;
596*811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
597*811d88abSStefano Zampini       PetscCall(DMPlexGetHeightStratum(cdm, 0, &ccStart, &ccEnd));
598*811d88abSStefano Zampini       PetscCall(DMPlexGetChart(cdm, &cpStart, &cpEnd));
599*811d88abSStefano Zampini       PetscCall(PetscCalloc1(nparts, &npoints));
600*811d88abSStefano Zampini       PetscCall(PetscMalloc4(cEnd - cStart, &points, cEnd - cStart, &ranks, nparts + 1, &starts, cEnd - cStart, &gidxs));
601*811d88abSStefano Zampini       PetscCall(PetscSFGetGraph(sfXC, &nr, NULL, NULL, NULL));
602*811d88abSStefano Zampini       PetscCall(PetscMalloc2(cpEnd - cpStart, &cranks_leaf, nr, &cranks_root));
603*811d88abSStefano Zampini       for (PetscInt c = 0; c < cpEnd - cpStart; c++) cranks_leaf[c] = rank;
604*811d88abSStefano Zampini       PetscCall(PetscSFReduceBegin(sfXC, MPIU_INT, cranks_leaf, cranks_root, MPI_REPLACE));
605*811d88abSStefano Zampini       PetscCall(PetscSFReduceEnd(sfXC, MPIU_INT, cranks_leaf, cranks_root, MPI_REPLACE));
606*811d88abSStefano Zampini 
607*811d88abSStefano Zampini       PetscCall(VecGetArray(map, &array));
608*811d88abSStefano Zampini       for (PetscInt c = 0; c < cEnd - cStart; c++) gidxs[c] = (PetscInt)PetscRealPart(array[c]);
609*811d88abSStefano Zampini       PetscCall(VecRestoreArray(map, &array));
610*811d88abSStefano Zampini 
611*811d88abSStefano Zampini       PetscCall(PetscLayoutCreate(comm, &clayout));
612*811d88abSStefano Zampini       PetscCall(PetscLayoutSetLocalSize(clayout, nr));
613*811d88abSStefano Zampini       PetscCall(PetscSFSetGraphLayout(sf, clayout, cEnd - cStart, NULL, PETSC_OWN_POINTER, gidxs));
614*811d88abSStefano Zampini       PetscCall(PetscLayoutDestroy(&clayout));
615*811d88abSStefano Zampini 
616*811d88abSStefano Zampini       PetscCall(PetscSFBcastBegin(sf, MPIU_INT, cranks_root, ranks, MPI_REPLACE));
617*811d88abSStefano Zampini       PetscCall(PetscSFBcastEnd(sf, MPIU_INT, cranks_root, ranks, MPI_REPLACE));
618*811d88abSStefano Zampini       PetscCall(PetscSFDestroy(&sf));
619*811d88abSStefano Zampini       PetscCall(PetscFree2(cranks_leaf, cranks_root));
620*811d88abSStefano Zampini       for (PetscInt c = 0; c < cEnd - cStart; c++) npoints[ranks[c]]++;
621*811d88abSStefano Zampini 
622*811d88abSStefano Zampini       starts[0] = 0;
623*811d88abSStefano Zampini       for (PetscInt c = 0; c < nparts; c++) starts[c + 1] = starts[c] + npoints[c];
624*811d88abSStefano Zampini       for (PetscInt c = 0; c < cEnd - cStart; c++) points[starts[ranks[c]]++] = c;
625*811d88abSStefano Zampini       PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
626*811d88abSStefano Zampini       PetscCall(PetscPartitionerShellSetPartition(part, nparts, npoints, points));
627*811d88abSStefano Zampini       PetscCall(PetscFree(npoints));
628*811d88abSStefano Zampini       PetscCall(PetscFree4(points, ranks, starts, gidxs));
629*811d88abSStefano Zampini       PetscCall(DMRestoreLocalVector(sdm, &map));
630*811d88abSStefano Zampini       PetscCall(DMDestroy(&sdm));
631*811d88abSStefano Zampini     }
632*811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfXC));
633*811d88abSStefano Zampini     PetscCall(DMPlexDistribute(dm, 0, &sfBC, odm));
634*811d88abSStefano Zampini     if (*odm) {
635*811d88abSStefano Zampini       PetscCall(DMDestroy(&dm));
636*811d88abSStefano Zampini       dm   = *odm;
637*811d88abSStefano Zampini       *odm = NULL;
638*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)dm, dmname));
639*811d88abSStefano Zampini     }
640*811d88abSStefano Zampini     if (sfBC) PetscCall(PetscSFCompose(sfXB, sfBC, &sfXC));
641*811d88abSStefano Zampini     else {
642*811d88abSStefano Zampini       PetscCall(PetscObjectReference((PetscObject)sfXB));
643*811d88abSStefano Zampini       sfXC = sfXB;
644*811d88abSStefano Zampini     }
645*811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfXB));
646*811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfBC));
647*811d88abSStefano Zampini     PetscCall(DMSetCoarsenLevel(dm, cc));
648*811d88abSStefano Zampini     PetscCall(DMSetRefineLevel(dm, rr));
649*811d88abSStefano Zampini     PetscCall(DMPlexSetRegularRefinement(dm, isRegular));
650*811d88abSStefano Zampini     PetscCall(DMPlexSetRefinementUniform(dm, isUniform));
651*811d88abSStefano Zampini     PetscCall(DMPlexSectionLoad(dm, viewer, NULL, sfXC, &sfG, NULL));
652*811d88abSStefano Zampini     if (level == numlevels - 1) {
653*811d88abSStefano Zampini       Vec u;
654*811d88abSStefano Zampini 
655*811d88abSStefano Zampini       PetscCall(DMGetNamedGlobalVector(dm, "solution_", &u));
656*811d88abSStefano Zampini       PetscCall(PetscObjectSetName((PetscObject)u, "solution_"));
657*811d88abSStefano Zampini       PetscCall(DMPlexGlobalVectorLoad(dm, viewer, NULL, sfG, u));
658*811d88abSStefano Zampini       PetscCall(DMRestoreNamedGlobalVector(dm, "solution_", &u));
659*811d88abSStefano Zampini     }
660*811d88abSStefano Zampini     PetscCall(PetscFree(dmname));
661*811d88abSStefano Zampini     PetscCall(PetscSFDestroy(&sfG));
662*811d88abSStefano Zampini     PetscCall(DMSetCoarseDM(dm, cdm));
663*811d88abSStefano Zampini     PetscCall(DMDestroy(&cdm));
664*811d88abSStefano Zampini     PetscCall(PetscViewerHDF5PopGroup(viewer));
665*811d88abSStefano Zampini     cdm = dm;
666*811d88abSStefano Zampini   }
667*811d88abSStefano Zampini   *odm = dm;
668*811d88abSStefano Zampini   PetscCall(PetscViewerPopFormat(viewer));
669*811d88abSStefano Zampini   PetscCall(PetscViewerDestroy(&viewer));
670*811d88abSStefano Zampini   PetscCall(PetscSFDestroy(&sfXC));
671*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
672*811d88abSStefano Zampini #else
673*811d88abSStefano Zampini   SETERRQ(comm, PETSC_ERR_SUP, "Needs HDF5 support. Please reconfigure using --download-hdf5");
674*811d88abSStefano Zampini #endif
675*811d88abSStefano Zampini }
676*811d88abSStefano Zampini 
677*811d88abSStefano Zampini /* Project source function and make it zero-mean */
678*811d88abSStefano Zampini static PetscErrorCode ProjectSource(DM dm, PetscReal time, AppCtx *ctx)
679*811d88abSStefano Zampini {
680*811d88abSStefano Zampini   PetscInt    id = C_FIELD_ID;
681*811d88abSStefano Zampini   DM          dmAux;
682*811d88abSStefano Zampini   Vec         u, lu;
683*811d88abSStefano Zampini   IS          is;
684*811d88abSStefano Zampini   void       *ctxs[NUM_FIELDS];
685*811d88abSStefano Zampini   PetscScalar vals[NUM_FIELDS];
686*811d88abSStefano Zampini   PetscDS     ds;
687*811d88abSStefano Zampini   PetscErrorCode (*funcs[NUM_FIELDS])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
688*811d88abSStefano Zampini 
689*811d88abSStefano Zampini   PetscFunctionBeginUser;
690*811d88abSStefano Zampini   switch (ctx->source_num) {
691*811d88abSStefano Zampini   case 0:
692*811d88abSStefano Zampini     funcs[P_FIELD_ID] = source_0;
693*811d88abSStefano Zampini     ctxs[P_FIELD_ID]  = ctx->x0;
694*811d88abSStefano Zampini     break;
695*811d88abSStefano Zampini   case 1:
696*811d88abSStefano Zampini     funcs[P_FIELD_ID] = source_1;
697*811d88abSStefano Zampini     ctxs[P_FIELD_ID]  = NULL;
698*811d88abSStefano Zampini     break;
699*811d88abSStefano Zampini   default:
700*811d88abSStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unknwon source");
701*811d88abSStefano Zampini   }
702*811d88abSStefano Zampini   funcs[C_FIELD_ID] = zerof;
703*811d88abSStefano Zampini   ctxs[C_FIELD_ID]  = NULL;
704*811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
705*811d88abSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &u));
706*811d88abSStefano Zampini   PetscCall(DMProjectFunction(dm, time, funcs, ctxs, INSERT_ALL_VALUES, u));
707*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, average));
708*811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
709*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
710*811d88abSStefano Zampini   PetscCall(VecShift(u, -vals[P_FIELD_ID]));
711*811d88abSStefano Zampini   PetscCall(DMCreateSubDM(dm, 1, &id, &is, NULL));
712*811d88abSStefano Zampini   PetscCall(VecISSet(u, is, 0));
713*811d88abSStefano Zampini   PetscCall(ISDestroy(&is));
714*811d88abSStefano Zampini 
715*811d88abSStefano Zampini   /* Attach source vector as auxiliary vector:
716*811d88abSStefano Zampini      Use a different DM to break ref cycles */
717*811d88abSStefano Zampini   PetscCall(DMClone(dm, &dmAux));
718*811d88abSStefano Zampini   PetscCall(DMCopyDisc(dm, dmAux));
719*811d88abSStefano Zampini   PetscCall(DMCreateLocalVector(dmAux, &lu));
720*811d88abSStefano Zampini   PetscCall(DMDestroy(&dmAux));
721*811d88abSStefano Zampini   PetscCall(DMGlobalToLocal(dm, u, INSERT_VALUES, lu));
722*811d88abSStefano Zampini   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, lu));
723*811d88abSStefano Zampini   PetscCall(VecViewFromOptions(lu, NULL, "-aux_view"));
724*811d88abSStefano Zampini   PetscCall(VecDestroy(&lu));
725*811d88abSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &u));
726*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
727*811d88abSStefano Zampini }
728*811d88abSStefano Zampini 
729*811d88abSStefano Zampini /* callback for the creation of the potential null space */
730*811d88abSStefano Zampini static PetscErrorCode CreatePotentialNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
731*811d88abSStefano Zampini {
732*811d88abSStefano Zampini   Vec vec;
733*811d88abSStefano Zampini   PetscErrorCode (*funcs[NUM_FIELDS])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zerof};
734*811d88abSStefano Zampini 
735*811d88abSStefano Zampini   PetscFunctionBeginUser;
736*811d88abSStefano Zampini   funcs[nfield] = constantf;
737*811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &vec));
738*811d88abSStefano Zampini   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
739*811d88abSStefano Zampini   PetscCall(VecNormalize(vec, NULL));
740*811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)vec, "Potential Null Space"));
741*811d88abSStefano Zampini   PetscCall(VecViewFromOptions(vec, NULL, "-potential_nullspace_view"));
742*811d88abSStefano Zampini   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace));
743*811d88abSStefano Zampini   /* break ref cycles */
744*811d88abSStefano Zampini   PetscCall(VecSetDM(vec, NULL));
745*811d88abSStefano Zampini   PetscCall(VecDestroy(&vec));
746*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
747*811d88abSStefano Zampini }
748*811d88abSStefano Zampini 
749*811d88abSStefano Zampini /* customize residuals and Jacobians */
750*811d88abSStefano Zampini static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
751*811d88abSStefano Zampini {
752*811d88abSStefano Zampini   PetscDS     ds;
753*811d88abSStefano Zampini   PetscInt    cdim, dim;
754*811d88abSStefano Zampini   PetscScalar constants[NUM_CONSTANTS];
755*811d88abSStefano Zampini 
756*811d88abSStefano Zampini   PetscFunctionBeginUser;
757*811d88abSStefano Zampini   constants[R_ID]     = ctx->r;
758*811d88abSStefano Zampini   constants[EPS_ID]   = ctx->eps;
759*811d88abSStefano Zampini   constants[ALPHA_ID] = ctx->alpha;
760*811d88abSStefano Zampini   constants[GAMMA_ID] = ctx->gamma;
761*811d88abSStefano Zampini   constants[D_ID]     = ctx->D;
762*811d88abSStefano Zampini   constants[C2_ID]    = ctx->c * ctx->c;
763*811d88abSStefano Zampini 
764*811d88abSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
765*811d88abSStefano Zampini   PetscCall(DMGetCoordinateDim(dm, &cdim));
766*811d88abSStefano Zampini   PetscCheck(dim == 2 && cdim == 2, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Only for 2D meshes");
767*811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
768*811d88abSStefano Zampini   PetscCall(PetscDSSetConstants(ds, NUM_CONSTANTS, constants));
769*811d88abSStefano Zampini   PetscCall(PetscDSSetImplicit(ds, C_FIELD_ID, PETSC_TRUE));
770*811d88abSStefano Zampini   PetscCall(PetscDSSetImplicit(ds, P_FIELD_ID, PETSC_TRUE));
771*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, energy));
772*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
773*811d88abSStefano Zampini   PetscCall(PetscDSSetResidual(ds, C_FIELD_ID, C_0, C_1));
774*811d88abSStefano Zampini   PetscCall(PetscDSSetResidual(ds, P_FIELD_ID, P_0, P_1));
775*811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, C_FIELD_ID, C_FIELD_ID, JC_0_c0c0, NULL, NULL, JC_1_c1c1));
776*811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, C_FIELD_ID, P_FIELD_ID, NULL, JC_0_c0p1, NULL, NULL));
777*811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, P_FIELD_ID, C_FIELD_ID, NULL, NULL, JP_1_p1c0, NULL));
778*811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, P_FIELD_ID, P_FIELD_ID, NULL, NULL, NULL, JP_1_p1p1));
779*811d88abSStefano Zampini 
780*811d88abSStefano Zampini   /* Attach potential nullspace */
781*811d88abSStefano Zampini   PetscCall(DMSetNullSpaceConstructor(dm, P_FIELD_ID, CreatePotentialNullSpace));
782*811d88abSStefano Zampini 
783*811d88abSStefano Zampini   /* Attach source function as auxiliary vector */
784*811d88abSStefano Zampini   PetscCall(ProjectSource(dm, 0, ctx));
785*811d88abSStefano Zampini 
786*811d88abSStefano Zampini   /* Add callbacks */
787*811d88abSStefano Zampini   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, NULL));
788*811d88abSStefano Zampini   PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, NULL));
789*811d88abSStefano Zampini   PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, NULL));
790*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
791*811d88abSStefano Zampini }
792*811d88abSStefano Zampini 
793*811d88abSStefano Zampini /* setup discrete spaces and residuals */
794*811d88abSStefano Zampini static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
795*811d88abSStefano Zampini {
796*811d88abSStefano Zampini   DM           plex, cdm = dm;
797*811d88abSStefano Zampini   PetscFE      feC, feP;
798*811d88abSStefano Zampini   PetscBool    simplex;
799*811d88abSStefano Zampini   PetscInt     dim;
800*811d88abSStefano Zampini   MPI_Comm     comm = PetscObjectComm((PetscObject)dm);
801*811d88abSStefano Zampini   MatNullSpace nsp;
802*811d88abSStefano Zampini 
803*811d88abSStefano Zampini   PetscFunctionBeginUser;
804*811d88abSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
805*811d88abSStefano Zampini 
806*811d88abSStefano Zampini   PetscCall(DMConvert(dm, DMPLEX, &plex));
807*811d88abSStefano Zampini   PetscCall(DMPlexIsSimplex(plex, &simplex));
808*811d88abSStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &simplex, 1, MPIU_BOOL, MPI_LOR, comm));
809*811d88abSStefano Zampini   PetscCall(DMDestroy(&plex));
810*811d88abSStefano Zampini 
811*811d88abSStefano Zampini   /* We model Cij in H^1 with Cij = Cji -> dim*(dim+1)/2 components */
812*811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(comm, dim, (dim * (dim + 1)) / 2, simplex, "c_", -1, &feC));
813*811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)feC, "conductivity"));
814*811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "p_", -1, &feP));
815*811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)feP, "potential"));
816*811d88abSStefano Zampini   PetscCall(MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, &nsp));
817*811d88abSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)feP, "nullspace", (PetscObject)nsp));
818*811d88abSStefano Zampini   PetscCall(MatNullSpaceDestroy(&nsp));
819*811d88abSStefano Zampini   PetscCall(PetscFECopyQuadrature(feC, feP));
820*811d88abSStefano Zampini 
821*811d88abSStefano Zampini   PetscCall(DMSetNumFields(dm, 2));
822*811d88abSStefano Zampini   PetscCall(DMSetField(dm, C_FIELD_ID, NULL, (PetscObject)feC));
823*811d88abSStefano Zampini   PetscCall(DMSetField(dm, P_FIELD_ID, NULL, (PetscObject)feP));
824*811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feC));
825*811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feP));
826*811d88abSStefano Zampini   PetscCall(DMCreateDS(dm));
827*811d88abSStefano Zampini   while (cdm) {
828*811d88abSStefano Zampini     PetscCall(DMCopyDisc(dm, cdm));
829*811d88abSStefano Zampini     PetscCall(SetupProblem(cdm, ctx));
830*811d88abSStefano Zampini     PetscCall(DMGetCoarseDM(cdm, &cdm));
831*811d88abSStefano Zampini   }
832*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
833*811d88abSStefano Zampini }
834*811d88abSStefano Zampini 
835*811d88abSStefano Zampini /* Create mesh by command line options */
836*811d88abSStefano Zampini static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
837*811d88abSStefano Zampini {
838*811d88abSStefano Zampini   PetscFunctionBeginUser;
839*811d88abSStefano Zampini   if (ctx->load) {
840*811d88abSStefano Zampini     PetscInt  refine = 0;
841*811d88abSStefano Zampini     PetscBool isHierarchy;
842*811d88abSStefano Zampini     DM       *dms;
843*811d88abSStefano Zampini     char      typeName[256];
844*811d88abSStefano Zampini     PetscBool flg;
845*811d88abSStefano Zampini 
846*811d88abSStefano Zampini     PetscCall(LoadFromFile(comm, ctx->load_filename, dm));
847*811d88abSStefano Zampini     PetscOptionsBegin(comm, "", "Additional mesh options", "DMPLEX");
848*811d88abSStefano Zampini     PetscCall(PetscOptionsFList("-dm_mat_type", "Matrix type used for created matrices", "DMSetMatType", MatList, MATAIJ, typeName, sizeof(typeName), &flg));
849*811d88abSStefano Zampini     if (flg) PetscCall(DMSetMatType(*dm, typeName));
850*811d88abSStefano Zampini     PetscCall(PetscOptionsBoundedInt("-dm_refine", "The number of uniform refinements", "DMCreate", refine, &refine, NULL, 0));
851*811d88abSStefano Zampini     PetscCall(PetscOptionsBoundedInt("-dm_refine_hierarchy", "The number of uniform refinements", "DMCreate", refine, &refine, &isHierarchy, 0));
852*811d88abSStefano Zampini     PetscOptionsEnd();
853*811d88abSStefano Zampini     if (refine) {
854*811d88abSStefano Zampini       PetscCall(SetupDiscretization(*dm, ctx));
855*811d88abSStefano Zampini       PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE));
856*811d88abSStefano Zampini     }
857*811d88abSStefano Zampini     PetscCall(PetscCalloc1(refine, &dms));
858*811d88abSStefano Zampini     if (isHierarchy) PetscCall(DMRefineHierarchy(*dm, refine, dms));
859*811d88abSStefano Zampini     for (PetscInt r = 0; r < refine; r++) {
860*811d88abSStefano Zampini       Mat M;
861*811d88abSStefano Zampini       DM  dmr = dms[r];
862*811d88abSStefano Zampini       Vec u, ur;
863*811d88abSStefano Zampini 
864*811d88abSStefano Zampini       if (!isHierarchy) {
865*811d88abSStefano Zampini         PetscCall(DMRefine(*dm, PetscObjectComm((PetscObject)dm), &dmr));
866*811d88abSStefano Zampini         PetscCall(DMSetCoarseDM(dmr, *dm));
867*811d88abSStefano Zampini       }
868*811d88abSStefano Zampini       PetscCall(DMCreateInterpolation(*dm, dmr, &M, NULL));
869*811d88abSStefano Zampini       PetscCall(DMGetNamedGlobalVector(*dm, "solution_", &u));
870*811d88abSStefano Zampini       PetscCall(DMGetNamedGlobalVector(dmr, "solution_", &ur));
871*811d88abSStefano Zampini       PetscCall(MatInterpolate(M, u, ur));
872*811d88abSStefano Zampini       PetscCall(DMRestoreNamedGlobalVector(*dm, "solution_", &u));
873*811d88abSStefano Zampini       PetscCall(DMRestoreNamedGlobalVector(dmr, "solution_", &ur));
874*811d88abSStefano Zampini       PetscCall(MatDestroy(&M));
875*811d88abSStefano Zampini       if (!isHierarchy) PetscCall(DMSetCoarseDM(dmr, NULL));
876*811d88abSStefano Zampini       PetscCall(DMDestroy(dm));
877*811d88abSStefano Zampini       *dm = dmr;
878*811d88abSStefano Zampini     }
879*811d88abSStefano Zampini     if (refine && !isHierarchy) PetscCall(DMSetRefineLevel(*dm, 0));
880*811d88abSStefano Zampini     PetscCall(PetscFree(dms));
881*811d88abSStefano Zampini   } else {
882*811d88abSStefano Zampini     PetscCall(DMCreate(comm, dm));
883*811d88abSStefano Zampini     PetscCall(DMSetType(*dm, DMPLEX));
884*811d88abSStefano Zampini     PetscCall(DMSetFromOptions(*dm));
885*811d88abSStefano Zampini     PetscCall(DMLocalizeCoordinates(*dm));
886*811d88abSStefano Zampini     {
887*811d88abSStefano Zampini       char      convType[256];
888*811d88abSStefano Zampini       PetscBool flg;
889*811d88abSStefano Zampini       PetscOptionsBegin(comm, "", "Additional mesh options", "DMPLEX");
890*811d88abSStefano Zampini       PetscCall(PetscOptionsFList("-dm_plex_convert_type", "Convert DMPlex to another format", __FILE__, DMList, DMPLEX, convType, 256, &flg));
891*811d88abSStefano Zampini       PetscOptionsEnd();
892*811d88abSStefano Zampini       if (flg) {
893*811d88abSStefano Zampini         DM dmConv;
894*811d88abSStefano Zampini         PetscCall(DMConvert(*dm, convType, &dmConv));
895*811d88abSStefano Zampini         if (dmConv) {
896*811d88abSStefano Zampini           PetscCall(DMDestroy(dm));
897*811d88abSStefano Zampini           *dm = dmConv;
898*811d88abSStefano Zampini           PetscCall(DMSetFromOptions(*dm));
899*811d88abSStefano Zampini           PetscCall(DMSetUp(*dm));
900*811d88abSStefano Zampini         }
901*811d88abSStefano Zampini       }
902*811d88abSStefano Zampini     }
903*811d88abSStefano Zampini   }
904*811d88abSStefano Zampini   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
905*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
906*811d88abSStefano Zampini }
907*811d88abSStefano Zampini 
908*811d88abSStefano Zampini /* Make potential field zero mean */
909*811d88abSStefano Zampini static PetscErrorCode ZeroMeanPotential(DM dm, Vec u)
910*811d88abSStefano Zampini {
911*811d88abSStefano Zampini   PetscScalar vals[NUM_FIELDS];
912*811d88abSStefano Zampini   PetscDS     ds;
913*811d88abSStefano Zampini   IS          is;
914*811d88abSStefano Zampini 
915*811d88abSStefano Zampini   PetscFunctionBeginUser;
916*811d88abSStefano Zampini   PetscCall(PetscObjectQuery((PetscObject)dm, "IS potential", (PetscObject *)&is));
917*811d88abSStefano Zampini   PetscCheck(is, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing potential IS");
918*811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
919*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, average));
920*811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
921*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
922*811d88abSStefano Zampini   PetscCall(VecISShift(u, is, -vals[P_FIELD_ID]));
923*811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
924*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
925*811d88abSStefano Zampini }
926*811d88abSStefano Zampini 
927*811d88abSStefano Zampini /* Compute initial conditions and exclude potential from local truncation error
928*811d88abSStefano Zampini    Since we are solving a DAE, once the initial conditions for the differential
929*811d88abSStefano Zampini    variables are set, we need to compute the corresponding value for the
930*811d88abSStefano Zampini    algebraic variables. We do so by creating a subDM for the potential only
931*811d88abSStefano Zampini    and solve a static problem with SNES */
932*811d88abSStefano Zampini static PetscErrorCode SetInitialConditionsAndTolerances(TS ts, PetscInt nv, Vec vecs[], PetscBool valid)
933*811d88abSStefano Zampini {
934*811d88abSStefano Zampini   DM         dm;
935*811d88abSStefano Zampini   Vec        tu, u, p, lsource, subaux, vatol, vrtol;
936*811d88abSStefano Zampini   PetscReal  t, atol, rtol;
937*811d88abSStefano Zampini   PetscInt   fields[NUM_FIELDS] = {C_FIELD_ID, P_FIELD_ID};
938*811d88abSStefano Zampini   IS         isp;
939*811d88abSStefano Zampini   DM         dmp;
940*811d88abSStefano Zampini   VecScatter sctp = NULL;
941*811d88abSStefano Zampini   PetscDS    ds;
942*811d88abSStefano Zampini   SNES       snes;
943*811d88abSStefano Zampini   KSP        ksp;
944*811d88abSStefano Zampini   PC         pc;
945*811d88abSStefano Zampini   AppCtx    *ctx;
946*811d88abSStefano Zampini 
947*811d88abSStefano Zampini   PetscFunctionBeginUser;
948*811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
949*811d88abSStefano Zampini   PetscCall(TSGetApplicationContext(ts, &ctx));
950*811d88abSStefano Zampini   if (nv > 1) {
951*811d88abSStefano Zampini     PetscCall(DMCreateSubDM(dm, NUM_FIELDS - 1, fields + 1, &isp, NULL));
952*811d88abSStefano Zampini     PetscCall(PetscObjectCompose((PetscObject)dm, "IS potential", (PetscObject)isp));
953*811d88abSStefano Zampini     PetscCall(DMCreateGlobalVector(dm, &vatol));
954*811d88abSStefano Zampini     PetscCall(DMCreateGlobalVector(dm, &vrtol));
955*811d88abSStefano Zampini     PetscCall(TSGetTolerances(ts, &atol, NULL, &rtol, NULL));
956*811d88abSStefano Zampini     PetscCall(VecSet(vatol, atol));
957*811d88abSStefano Zampini     PetscCall(VecISSet(vatol, isp, -1));
958*811d88abSStefano Zampini     PetscCall(VecSet(vrtol, rtol));
959*811d88abSStefano Zampini     PetscCall(VecISSet(vrtol, isp, -1));
960*811d88abSStefano Zampini     PetscCall(TSSetTolerances(ts, atol, vatol, rtol, vrtol));
961*811d88abSStefano Zampini     PetscCall(VecDestroy(&vatol));
962*811d88abSStefano Zampini     PetscCall(VecDestroy(&vrtol));
963*811d88abSStefano Zampini     PetscCall(ISDestroy(&isp));
964*811d88abSStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
965*811d88abSStefano Zampini   }
966*811d88abSStefano Zampini   PetscCall(DMCreateSubDM(dm, NUM_FIELDS - 1, fields + 1, &isp, &dmp));
967*811d88abSStefano Zampini   PetscCall(PetscObjectCompose((PetscObject)dm, "IS potential", (PetscObject)isp));
968*811d88abSStefano Zampini   PetscCall(DMSetMatType(dmp, MATAIJ));
969*811d88abSStefano Zampini   PetscCall(DMGetDS(dmp, &ds));
970*811d88abSStefano Zampini   //PetscCall(PetscDSSetResidual(ds, 0, P_0, P_1_aux));
971*811d88abSStefano Zampini   PetscCall(PetscDSSetResidual(ds, 0, 0, P_1_aux));
972*811d88abSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, JP_1_p1p1_aux));
973*811d88abSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(dmp, PETSC_FALSE, NULL));
974*811d88abSStefano Zampini 
975*811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dmp, &p));
976*811d88abSStefano Zampini 
977*811d88abSStefano Zampini   PetscCall(SNESCreate(PetscObjectComm((PetscObject)dmp), &snes));
978*811d88abSStefano Zampini   PetscCall(SNESSetDM(snes, dmp));
979*811d88abSStefano Zampini   PetscCall(SNESSetOptionsPrefix(snes, "initial_"));
980*811d88abSStefano Zampini   PetscCall(SNESSetErrorIfNotConverged(snes, PETSC_TRUE));
981*811d88abSStefano Zampini   PetscCall(SNESGetKSP(snes, &ksp));
982*811d88abSStefano Zampini   PetscCall(KSPGetPC(ksp, &pc));
983*811d88abSStefano Zampini   PetscCall(PCSetType(pc, PCGAMG));
984*811d88abSStefano Zampini   PetscCall(SNESSetFromOptions(snes));
985*811d88abSStefano Zampini   PetscCall(SNESSetUp(snes));
986*811d88abSStefano Zampini 
987*811d88abSStefano Zampini   /* Loop over input vectors and compute corresponding potential */
988*811d88abSStefano Zampini   for (PetscInt i = 0; i < nv; i++) {
989*811d88abSStefano Zampini     PetscErrorCode (*funcs[NUM_FIELDS])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
990*811d88abSStefano Zampini 
991*811d88abSStefano Zampini     u = vecs[i];
992*811d88abSStefano Zampini     if (!valid) { /* Assumes entries in u are not valid */
993*811d88abSStefano Zampini       PetscCall(TSGetTime(ts, &t));
994*811d88abSStefano Zampini       switch (ctx->ic_num) {
995*811d88abSStefano Zampini       case 0:
996*811d88abSStefano Zampini         funcs[C_FIELD_ID] = initial_conditions_C_0;
997*811d88abSStefano Zampini         break;
998*811d88abSStefano Zampini       case 1:
999*811d88abSStefano Zampini         funcs[C_FIELD_ID] = initial_conditions_C_1;
1000*811d88abSStefano Zampini         break;
1001*811d88abSStefano Zampini       case 2:
1002*811d88abSStefano Zampini         funcs[C_FIELD_ID] = initial_conditions_C_2;
1003*811d88abSStefano Zampini         break;
1004*811d88abSStefano Zampini       default:
1005*811d88abSStefano Zampini         SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Unknwon IC");
1006*811d88abSStefano Zampini       }
1007*811d88abSStefano Zampini       funcs[P_FIELD_ID] = zerof;
1008*811d88abSStefano Zampini       PetscCall(DMProjectFunction(dm, t, funcs, NULL, INSERT_ALL_VALUES, u));
1009*811d88abSStefano Zampini     }
1010*811d88abSStefano Zampini 
1011*811d88abSStefano Zampini     /* pass conductivity and source information via auxiliary data */
1012*811d88abSStefano Zampini     PetscCall(DMGetGlobalVector(dm, &tu));
1013*811d88abSStefano Zampini     PetscCall(VecCopy(u, tu));
1014*811d88abSStefano Zampini     PetscCall(VecISSet(tu, isp, 0.0));
1015*811d88abSStefano Zampini     PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &lsource));
1016*811d88abSStefano Zampini     PetscCall(DMCreateLocalVector(dm, &subaux));
1017*811d88abSStefano Zampini     PetscCall(DMGlobalToLocal(dm, tu, INSERT_VALUES, subaux));
1018*811d88abSStefano Zampini     PetscCall(DMRestoreGlobalVector(dm, &tu));
1019*811d88abSStefano Zampini     PetscCall(VecAXPY(subaux, 1.0, lsource));
1020*811d88abSStefano Zampini     PetscCall(VecViewFromOptions(subaux, NULL, "-initial_aux_view"));
1021*811d88abSStefano Zampini     PetscCall(DMSetAuxiliaryVec(dmp, NULL, 0, 0, subaux));
1022*811d88abSStefano Zampini     PetscCall(VecDestroy(&subaux));
1023*811d88abSStefano Zampini 
1024*811d88abSStefano Zampini     /* solve the subproblem */
1025*811d88abSStefano Zampini     if (!sctp) PetscCall(VecScatterCreate(u, isp, p, NULL, &sctp));
1026*811d88abSStefano Zampini     PetscCall(VecScatterBegin(sctp, u, p, INSERT_VALUES, SCATTER_FORWARD));
1027*811d88abSStefano Zampini     PetscCall(VecScatterEnd(sctp, u, p, INSERT_VALUES, SCATTER_FORWARD));
1028*811d88abSStefano Zampini     PetscCall(SNESSolve(snes, NULL, p));
1029*811d88abSStefano Zampini 
1030*811d88abSStefano Zampini     /* scatter from potential only to full space */
1031*811d88abSStefano Zampini     PetscCall(VecScatterBegin(sctp, p, u, INSERT_VALUES, SCATTER_REVERSE));
1032*811d88abSStefano Zampini     PetscCall(VecScatterEnd(sctp, p, u, INSERT_VALUES, SCATTER_REVERSE));
1033*811d88abSStefano Zampini     PetscCall(ZeroMeanPotential(dm, u));
1034*811d88abSStefano Zampini   }
1035*811d88abSStefano Zampini   PetscCall(VecDestroy(&p));
1036*811d88abSStefano Zampini   PetscCall(DMDestroy(&dmp));
1037*811d88abSStefano Zampini   PetscCall(SNESDestroy(&snes));
1038*811d88abSStefano Zampini   PetscCall(VecScatterDestroy(&sctp));
1039*811d88abSStefano Zampini 
1040*811d88abSStefano Zampini   /* exclude potential from computation of the LTE */
1041*811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &vatol));
1042*811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &vrtol));
1043*811d88abSStefano Zampini   PetscCall(TSGetTolerances(ts, &atol, NULL, &rtol, NULL));
1044*811d88abSStefano Zampini   PetscCall(VecSet(vatol, atol));
1045*811d88abSStefano Zampini   PetscCall(VecISSet(vatol, isp, -1));
1046*811d88abSStefano Zampini   PetscCall(VecSet(vrtol, rtol));
1047*811d88abSStefano Zampini   PetscCall(VecISSet(vrtol, isp, -1));
1048*811d88abSStefano Zampini   PetscCall(TSSetTolerances(ts, atol, vatol, rtol, vrtol));
1049*811d88abSStefano Zampini   PetscCall(VecDestroy(&vatol));
1050*811d88abSStefano Zampini   PetscCall(VecDestroy(&vrtol));
1051*811d88abSStefano Zampini   PetscCall(ISDestroy(&isp));
1052*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1053*811d88abSStefano Zampini }
1054*811d88abSStefano Zampini 
1055*811d88abSStefano Zampini /* mesh adaption context */
1056*811d88abSStefano Zampini typedef struct {
1057*811d88abSStefano Zampini   VecTagger refineTag;
1058*811d88abSStefano Zampini   DMLabel   adaptLabel;
1059*811d88abSStefano Zampini   PetscInt  cnt;
1060*811d88abSStefano Zampini } AdaptCtx;
1061*811d88abSStefano Zampini 
1062*811d88abSStefano Zampini static PetscErrorCode ResizeSetUp(TS ts, PetscInt nstep, PetscReal time, Vec u, PetscBool *resize, void *vctx)
1063*811d88abSStefano Zampini {
1064*811d88abSStefano Zampini   AdaptCtx *ctx = (AdaptCtx *)vctx;
1065*811d88abSStefano Zampini   Vec       ellVecCells, ellVecCellsF;
1066*811d88abSStefano Zampini   DM        dm, plex;
1067*811d88abSStefano Zampini   PetscDS   ds;
1068*811d88abSStefano Zampini   PetscReal norm;
1069*811d88abSStefano Zampini   PetscInt  cStart, cEnd;
1070*811d88abSStefano Zampini 
1071*811d88abSStefano Zampini   PetscFunctionBeginUser;
1072*811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1073*811d88abSStefano Zampini   PetscCall(DMConvert(dm, DMPLEX, &plex));
1074*811d88abSStefano Zampini   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
1075*811d88abSStefano Zampini   PetscCall(DMDestroy(&plex));
1076*811d88abSStefano Zampini   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)ts), NUM_FIELDS * (cEnd - cStart), PETSC_DECIDE, &ellVecCellsF));
1077*811d88abSStefano Zampini   PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)ts), cEnd - cStart, PETSC_DECIDE, &ellVecCells));
1078*811d88abSStefano Zampini   PetscCall(VecSetBlockSize(ellVecCellsF, NUM_FIELDS));
1079*811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
1080*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, ellipticity_fail));
1081*811d88abSStefano Zampini   PetscCall(DMPlexComputeCellwiseIntegralFEM(dm, u, ellVecCellsF, NULL));
1082*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, energy));
1083*811d88abSStefano Zampini   PetscCall(VecStrideGather(ellVecCellsF, C_FIELD_ID, ellVecCells, INSERT_VALUES));
1084*811d88abSStefano Zampini   PetscCall(VecDestroy(&ellVecCellsF));
1085*811d88abSStefano Zampini   PetscCall(VecNorm(ellVecCells, NORM_1, &norm));
1086*811d88abSStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "STEP %d norm %g\n", (int)nstep, (double)norm));
1087*811d88abSStefano Zampini   if (norm && !ctx->cnt) {
1088*811d88abSStefano Zampini     IS refineIS;
1089*811d88abSStefano Zampini 
1090*811d88abSStefano Zampini     *resize = PETSC_TRUE;
1091*811d88abSStefano Zampini     if (!ctx->refineTag) {
1092*811d88abSStefano Zampini       VecTaggerBox refineBox;
1093*811d88abSStefano Zampini       refineBox.min = PETSC_MACHINE_EPSILON;
1094*811d88abSStefano Zampini       refineBox.max = PETSC_MAX_REAL;
1095*811d88abSStefano Zampini 
1096*811d88abSStefano Zampini       PetscCall(VecTaggerCreate(PETSC_COMM_SELF, &ctx->refineTag));
1097*811d88abSStefano Zampini       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)ctx->refineTag, "refine_"));
1098*811d88abSStefano Zampini       PetscCall(VecTaggerSetType(ctx->refineTag, VECTAGGERABSOLUTE));
1099*811d88abSStefano Zampini       PetscCall(VecTaggerAbsoluteSetBox(ctx->refineTag, &refineBox));
1100*811d88abSStefano Zampini       PetscCall(VecTaggerSetFromOptions(ctx->refineTag));
1101*811d88abSStefano Zampini       PetscCall(VecTaggerSetUp(ctx->refineTag));
1102*811d88abSStefano Zampini       PetscCall(PetscObjectViewFromOptions((PetscObject)ctx->refineTag, NULL, "-tag_view"));
1103*811d88abSStefano Zampini     }
1104*811d88abSStefano Zampini     PetscCall(DMLabelDestroy(&ctx->adaptLabel));
1105*811d88abSStefano Zampini     PetscCall(DMLabelCreate(PetscObjectComm((PetscObject)ts), "adapt", &ctx->adaptLabel));
1106*811d88abSStefano Zampini     PetscCall(VecTaggerComputeIS(ctx->refineTag, ellVecCells, &refineIS, NULL));
1107*811d88abSStefano Zampini     PetscCall(DMLabelSetStratumIS(ctx->adaptLabel, DM_ADAPT_REFINE, refineIS));
1108*811d88abSStefano Zampini     PetscCall(ISDestroy(&refineIS));
1109*811d88abSStefano Zampini #if 0
1110*811d88abSStefano 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[]);
1111*811d88abSStefano Zampini     Vec ellVec;
1112*811d88abSStefano Zampini 
1113*811d88abSStefano Zampini     funcs[P_FIELD_ID] = ellipticity_fail;
1114*811d88abSStefano Zampini     funcs[C_FIELD_ID] = NULL;
1115*811d88abSStefano Zampini 
1116*811d88abSStefano Zampini     PetscCall(DMGetGlobalVector(dm, &ellVec));
1117*811d88abSStefano Zampini     PetscCall(DMProjectField(dm, 0, u, funcs, INSERT_VALUES, ellVec));
1118*811d88abSStefano Zampini     PetscCall(VecViewFromOptions(ellVec,NULL,"-view_amr_ell"));
1119*811d88abSStefano Zampini     PetscCall(DMRestoreGlobalVector(dm, &ellVec));
1120*811d88abSStefano Zampini #endif
1121*811d88abSStefano Zampini     ctx->cnt++;
1122*811d88abSStefano Zampini   } else {
1123*811d88abSStefano Zampini     ctx->cnt = 0;
1124*811d88abSStefano Zampini   }
1125*811d88abSStefano Zampini   PetscCall(VecDestroy(&ellVecCells));
1126*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1127*811d88abSStefano Zampini }
1128*811d88abSStefano Zampini 
1129*811d88abSStefano Zampini static PetscErrorCode ResizeTransfer(TS ts, PetscInt nv, Vec vecsin[], Vec vecsout[], void *vctx)
1130*811d88abSStefano Zampini {
1131*811d88abSStefano Zampini   AdaptCtx *actx = (AdaptCtx *)vctx;
1132*811d88abSStefano Zampini   AppCtx   *ctx;
1133*811d88abSStefano Zampini   DM        dm, adm;
1134*811d88abSStefano Zampini   PetscReal time;
1135*811d88abSStefano Zampini 
1136*811d88abSStefano Zampini   PetscFunctionBeginUser;
1137*811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1138*811d88abSStefano Zampini   PetscCheck(actx->adaptLabel, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Missing adaptLabel");
1139*811d88abSStefano Zampini   PetscCall(DMAdaptLabel(dm, actx->adaptLabel, &adm));
1140*811d88abSStefano Zampini   PetscCheck(adm, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Missing adapted DM");
1141*811d88abSStefano Zampini   PetscCall(TSGetTime(ts, &time));
1142*811d88abSStefano Zampini   PetscCall(DMLabelDestroy(&actx->adaptLabel));
1143*811d88abSStefano Zampini   for (PetscInt i = 0; i < nv; i++) {
1144*811d88abSStefano Zampini     PetscCall(DMCreateGlobalVector(adm, &vecsout[i]));
1145*811d88abSStefano Zampini     PetscCall(DMForestTransferVec(dm, vecsin[i], adm, vecsout[i], PETSC_TRUE, time));
1146*811d88abSStefano Zampini   }
1147*811d88abSStefano Zampini   PetscCall(DMForestSetAdaptivityForest(adm, NULL));
1148*811d88abSStefano Zampini   PetscCall(DMSetCoarseDM(adm, NULL));
1149*811d88abSStefano Zampini   PetscCall(DMSetLocalSection(adm, NULL));
1150*811d88abSStefano Zampini   PetscCall(TSSetDM(ts, adm));
1151*811d88abSStefano Zampini   PetscCall(TSGetTime(ts, &time));
1152*811d88abSStefano Zampini   PetscCall(TSGetApplicationContext(ts, &ctx));
1153*811d88abSStefano Zampini   PetscCall(ProjectSource(adm, time, ctx));
1154*811d88abSStefano Zampini   PetscCall(SetInitialConditionsAndTolerances(ts, nv, vecsout, PETSC_TRUE));
1155*811d88abSStefano Zampini   PetscCall(DMDestroy(&adm));
1156*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1157*811d88abSStefano Zampini }
1158*811d88abSStefano Zampini 
1159*811d88abSStefano Zampini static PetscErrorCode OutputVTK(DM dm, const char *filename, PetscViewer *viewer)
1160*811d88abSStefano Zampini {
1161*811d88abSStefano Zampini   PetscFunctionBeginUser;
1162*811d88abSStefano Zampini   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)dm), viewer));
1163*811d88abSStefano Zampini   PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERVTK));
1164*811d88abSStefano Zampini   PetscCall(PetscViewerFileSetName(*viewer, filename));
1165*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1166*811d88abSStefano Zampini }
1167*811d88abSStefano Zampini 
1168*811d88abSStefano Zampini /* Monitor relevant functionals */
1169*811d88abSStefano Zampini static PetscErrorCode Monitor(TS ts, PetscInt stepnum, PetscReal time, Vec u, void *vctx)
1170*811d88abSStefano Zampini {
1171*811d88abSStefano Zampini   PetscScalar vals[2 * NUM_FIELDS];
1172*811d88abSStefano Zampini   DM          dm;
1173*811d88abSStefano Zampini   PetscDS     ds;
1174*811d88abSStefano Zampini   SNES        snes;
1175*811d88abSStefano Zampini   PetscInt    nits, lits;
1176*811d88abSStefano Zampini   AppCtx     *ctx;
1177*811d88abSStefano Zampini 
1178*811d88abSStefano Zampini   PetscFunctionBeginUser;
1179*811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1180*811d88abSStefano Zampini   PetscCall(TSGetApplicationContext(ts, &ctx));
1181*811d88abSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
1182*811d88abSStefano Zampini 
1183*811d88abSStefano Zampini   /* monitor energy and potential average */
1184*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, average));
1185*811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals, NULL));
1186*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, P_FIELD_ID, zero));
1187*811d88abSStefano Zampini 
1188*811d88abSStefano Zampini   /* monitor ellipticity_fail */
1189*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, ellipticity_fail));
1190*811d88abSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, u, vals + NUM_FIELDS, NULL));
1191*811d88abSStefano Zampini   if (ctx->ellipticity) {
1192*811d88abSStefano 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[]);
1193*811d88abSStefano Zampini     Vec         ellVec;
1194*811d88abSStefano Zampini     PetscViewer viewer;
1195*811d88abSStefano Zampini     char        filename[PETSC_MAX_PATH_LEN];
1196*811d88abSStefano Zampini 
1197*811d88abSStefano Zampini     funcs[P_FIELD_ID] = ellipticity_fail;
1198*811d88abSStefano Zampini     funcs[C_FIELD_ID] = NULL;
1199*811d88abSStefano Zampini 
1200*811d88abSStefano Zampini     PetscCall(DMGetGlobalVector(dm, &ellVec));
1201*811d88abSStefano Zampini     PetscCall(DMProjectField(dm, 0, u, funcs, INSERT_VALUES, ellVec));
1202*811d88abSStefano Zampini     PetscCall(PetscSNPrintf(filename, sizeof filename, "ellipticity_fail-%03" PetscInt_FMT ".vtu", stepnum));
1203*811d88abSStefano Zampini     PetscCall(OutputVTK(dm, filename, &viewer));
1204*811d88abSStefano Zampini     PetscCall(VecView(ellVec, viewer));
1205*811d88abSStefano Zampini     PetscCall(PetscViewerDestroy(&viewer));
1206*811d88abSStefano Zampini     PetscCall(DMRestoreGlobalVector(dm, &ellVec));
1207*811d88abSStefano Zampini   }
1208*811d88abSStefano Zampini   PetscCall(PetscDSSetObjective(ds, C_FIELD_ID, energy));
1209*811d88abSStefano Zampini 
1210*811d88abSStefano Zampini   /* monitor linear and nonlinear iterations */
1211*811d88abSStefano Zampini   PetscCall(TSGetSNES(ts, &snes));
1212*811d88abSStefano Zampini   PetscCall(SNESGetIterationNumber(snes, &nits));
1213*811d88abSStefano Zampini   PetscCall(SNESGetLinearSolveIterations(snes, &lits));
1214*811d88abSStefano Zampini 
1215*811d88abSStefano 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])));
1216*811d88abSStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)ts), "         nonlinear its %" PetscInt_FMT ", linear its %" PetscInt_FMT "\n", nits, lits));
1217*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1218*811d88abSStefano Zampini }
1219*811d88abSStefano Zampini 
1220*811d88abSStefano Zampini /* Save restart information */
1221*811d88abSStefano Zampini static PetscErrorCode MonitorSave(TS ts, PetscInt steps, PetscReal time, Vec u, void *vctx)
1222*811d88abSStefano Zampini {
1223*811d88abSStefano Zampini   DM                dm;
1224*811d88abSStefano Zampini   AppCtx           *ctx        = (AppCtx *)vctx;
1225*811d88abSStefano Zampini   PetscInt          save_every = ctx->save_every;
1226*811d88abSStefano Zampini   TSConvergedReason reason;
1227*811d88abSStefano Zampini 
1228*811d88abSStefano Zampini   PetscFunctionBeginUser;
1229*811d88abSStefano Zampini   if (!ctx->save) PetscFunctionReturn(PETSC_SUCCESS);
1230*811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1231*811d88abSStefano Zampini   PetscCall(TSGetConvergedReason(ts, &reason));
1232*811d88abSStefano Zampini   if ((save_every > 0 && steps % save_every == 0) || (save_every == -1 && reason) || save_every < -1) PetscCall(SaveToFile(dm, u, ctx->save_filename));
1233*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1234*811d88abSStefano Zampini }
1235*811d88abSStefano Zampini 
1236*811d88abSStefano Zampini /* Make potential zero mean after SNES solve */
1237*811d88abSStefano Zampini static PetscErrorCode PostStage(TS ts, PetscReal stagetime, PetscInt stageindex, Vec *Y)
1238*811d88abSStefano Zampini {
1239*811d88abSStefano Zampini   DM  dm;
1240*811d88abSStefano Zampini   Vec u = Y[stageindex];
1241*811d88abSStefano Zampini 
1242*811d88abSStefano Zampini   PetscFunctionBeginUser;
1243*811d88abSStefano Zampini   PetscCall(TSGetDM(ts, &dm));
1244*811d88abSStefano Zampini   PetscCall(ZeroMeanPotential(dm, u));
1245*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1246*811d88abSStefano Zampini }
1247*811d88abSStefano Zampini 
1248*811d88abSStefano Zampini static PetscErrorCode VecViewFlux(Vec u, const char *opts)
1249*811d88abSStefano Zampini {
1250*811d88abSStefano Zampini   Vec       fluxVec;
1251*811d88abSStefano Zampini   DM        dmFlux, dm, plex;
1252*811d88abSStefano Zampini   PetscInt  dim;
1253*811d88abSStefano Zampini   PetscFE   feC, feFluxC, feNormC;
1254*811d88abSStefano Zampini   PetscBool simplex, has;
1255*811d88abSStefano Zampini 
1256*811d88abSStefano 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};
1257*811d88abSStefano Zampini 
1258*811d88abSStefano Zampini   PetscFunctionBeginUser;
1259*811d88abSStefano Zampini   PetscCall(PetscOptionsHasName(NULL, NULL, opts, &has));
1260*811d88abSStefano Zampini   if (!has) PetscFunctionReturn(PETSC_SUCCESS);
1261*811d88abSStefano Zampini   PetscCall(VecGetDM(u, &dm));
1262*811d88abSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
1263*811d88abSStefano Zampini   PetscCall(DMGetField(dm, C_FIELD_ID, NULL, (PetscObject *)&feC));
1264*811d88abSStefano Zampini   PetscCall(DMConvert(dm, DMPLEX, &plex));
1265*811d88abSStefano Zampini   PetscCall(DMPlexIsSimplex(plex, &simplex));
1266*811d88abSStefano Zampini   PetscCall(DMDestroy(&plex));
1267*811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, simplex, "flux_", -1, &feFluxC));
1268*811d88abSStefano Zampini   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, "normc_", -1, &feNormC));
1269*811d88abSStefano Zampini   PetscCall(PetscFECopyQuadrature(feC, feFluxC));
1270*811d88abSStefano Zampini   PetscCall(PetscFECopyQuadrature(feC, feNormC));
1271*811d88abSStefano Zampini   PetscCall(DMClone(dm, &dmFlux));
1272*811d88abSStefano Zampini   PetscCall(DMSetNumFields(dmFlux, 1));
1273*811d88abSStefano Zampini   PetscCall(DMSetField(dmFlux, 0, NULL, (PetscObject)feNormC));
1274*811d88abSStefano Zampini   /* paraview segfaults! */
1275*811d88abSStefano Zampini   //PetscCall(DMSetField(dmFlux, 1, NULL, (PetscObject)feFluxC));
1276*811d88abSStefano Zampini   PetscCall(DMCreateDS(dmFlux));
1277*811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feFluxC));
1278*811d88abSStefano Zampini   PetscCall(PetscFEDestroy(&feNormC));
1279*811d88abSStefano Zampini 
1280*811d88abSStefano Zampini   PetscCall(DMGetGlobalVector(dmFlux, &fluxVec));
1281*811d88abSStefano Zampini   PetscCall(DMProjectField(dmFlux, 0.0, u, funcs, INSERT_VALUES, fluxVec));
1282*811d88abSStefano Zampini   PetscCall(VecViewFromOptions(fluxVec, NULL, opts));
1283*811d88abSStefano Zampini   PetscCall(DMRestoreGlobalVector(dmFlux, &fluxVec));
1284*811d88abSStefano Zampini   PetscCall(DMDestroy(&dmFlux));
1285*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1286*811d88abSStefano Zampini }
1287*811d88abSStefano Zampini 
1288*811d88abSStefano Zampini static PetscErrorCode Run(MPI_Comm comm, AppCtx *ctx)
1289*811d88abSStefano Zampini {
1290*811d88abSStefano Zampini   DM        dm;
1291*811d88abSStefano Zampini   TS        ts;
1292*811d88abSStefano Zampini   Vec       u;
1293*811d88abSStefano Zampini   AdaptCtx *actx;
1294*811d88abSStefano Zampini   PetscBool flg;
1295*811d88abSStefano Zampini 
1296*811d88abSStefano Zampini   PetscFunctionBeginUser;
1297*811d88abSStefano Zampini   if (ctx->test_restart) {
1298*811d88abSStefano Zampini     PetscViewer subviewer;
1299*811d88abSStefano Zampini     PetscMPIInt rank;
1300*811d88abSStefano Zampini 
1301*811d88abSStefano Zampini     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1302*811d88abSStefano Zampini     PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1303*811d88abSStefano Zampini     if (ctx->load) PetscCall(PetscViewerASCIIPrintf(subviewer, "rank %d loading from %s\n", rank, ctx->load_filename));
1304*811d88abSStefano Zampini     if (ctx->save) PetscCall(PetscViewerASCIIPrintf(subviewer, "rank %d saving to %s\n", rank, ctx->save_filename));
1305*811d88abSStefano Zampini     PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1306*811d88abSStefano Zampini     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1307*811d88abSStefano Zampini   } else {
1308*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "----------------------------\n"));
1309*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "Simulation parameters:\n"));
1310*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  r    : %g\n", (double)ctx->r));
1311*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  eps  : %g\n", (double)ctx->eps));
1312*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  alpha: %g\n", (double)ctx->alpha));
1313*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  gamma: %g\n", (double)ctx->gamma));
1314*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  D    : %g\n", (double)ctx->D));
1315*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  c    : %g\n", (double)ctx->c));
1316*811d88abSStefano Zampini     if (ctx->load) PetscCall(PetscPrintf(comm, "  load : %s\n", ctx->load_filename));
1317*811d88abSStefano Zampini     else PetscCall(PetscPrintf(comm, "  IC   : %" PetscInt_FMT "\n", ctx->ic_num));
1318*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  S    : %" PetscInt_FMT "\n", ctx->source_num));
1319*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "  x0   : (%g,%g)\n", (double)ctx->x0[0], (double)ctx->x0[1]));
1320*811d88abSStefano Zampini     PetscCall(PetscPrintf(comm, "----------------------------\n"));
1321*811d88abSStefano Zampini   }
1322*811d88abSStefano Zampini 
1323*811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePush(SetupStage));
1324*811d88abSStefano Zampini   PetscCall(CreateMesh(comm, &dm, ctx));
1325*811d88abSStefano Zampini   PetscCall(SetupDiscretization(dm, ctx));
1326*811d88abSStefano Zampini 
1327*811d88abSStefano Zampini   PetscCall(TSCreate(comm, &ts));
1328*811d88abSStefano Zampini   PetscCall(TSSetApplicationContext(ts, ctx));
1329*811d88abSStefano Zampini 
1330*811d88abSStefano Zampini   PetscCall(TSSetDM(ts, dm));
1331*811d88abSStefano Zampini   if (ctx->test_restart) {
1332*811d88abSStefano Zampini     PetscViewer subviewer;
1333*811d88abSStefano Zampini     PetscMPIInt rank;
1334*811d88abSStefano Zampini     PetscInt    level;
1335*811d88abSStefano Zampini 
1336*811d88abSStefano Zampini     PetscCall(DMGetRefineLevel(dm, &level));
1337*811d88abSStefano Zampini     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1338*811d88abSStefano Zampini     PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1339*811d88abSStefano Zampini     PetscCall(PetscViewerASCIIPrintf(subviewer, "rank %d DM refinement level %" PetscInt_FMT "\n", rank, level));
1340*811d88abSStefano Zampini     PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, comm, &subviewer));
1341*811d88abSStefano Zampini     PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD));
1342*811d88abSStefano Zampini   }
1343*811d88abSStefano Zampini 
1344*811d88abSStefano Zampini   if (ctx->test_restart) PetscCall(TSSetMaxSteps(ts, 1));
1345*811d88abSStefano Zampini   PetscCall(TSSetMaxTime(ts, 10.0));
1346*811d88abSStefano Zampini   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1347*811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(TSMonitorSet(ts, Monitor, NULL, NULL));
1348*811d88abSStefano Zampini   PetscCall(TSMonitorSet(ts, MonitorSave, ctx, NULL));
1349*811d88abSStefano Zampini   PetscCall(PetscNew(&actx));
1350*811d88abSStefano Zampini   if (ctx->amr) PetscCall(TSSetResize(ts, PETSC_TRUE, ResizeSetUp, ResizeTransfer, actx));
1351*811d88abSStefano Zampini   PetscCall(TSSetPostStage(ts, PostStage));
1352*811d88abSStefano Zampini   PetscCall(TSSetMaxSNESFailures(ts, -1));
1353*811d88abSStefano Zampini   PetscCall(TSSetFromOptions(ts));
1354*811d88abSStefano Zampini 
1355*811d88abSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &u));
1356*811d88abSStefano Zampini   PetscCall(PetscObjectSetName((PetscObject)u, "solution_"));
1357*811d88abSStefano Zampini   PetscCall(DMHasNamedGlobalVector(dm, "solution_", &flg));
1358*811d88abSStefano Zampini   if (flg) {
1359*811d88abSStefano Zampini     Vec ru;
1360*811d88abSStefano Zampini 
1361*811d88abSStefano Zampini     PetscCall(DMGetNamedGlobalVector(dm, "solution_", &ru));
1362*811d88abSStefano Zampini     PetscCall(VecCopy(ru, u));
1363*811d88abSStefano Zampini     PetscCall(DMRestoreNamedGlobalVector(dm, "solution_", &ru));
1364*811d88abSStefano Zampini   }
1365*811d88abSStefano Zampini   PetscCall(SetInitialConditionsAndTolerances(ts, 1, &u, PETSC_FALSE));
1366*811d88abSStefano Zampini   PetscCall(TSSetSolution(ts, u));
1367*811d88abSStefano Zampini   PetscCall(VecDestroy(&u));
1368*811d88abSStefano Zampini   PetscCall(DMDestroy(&dm));
1369*811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePop());
1370*811d88abSStefano Zampini 
1371*811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePush(SolveStage));
1372*811d88abSStefano Zampini   PetscCall(TSSolve(ts, NULL));
1373*811d88abSStefano Zampini   if (!ctx->test_restart) PetscCall(PetscLogStagePop());
1374*811d88abSStefano Zampini 
1375*811d88abSStefano Zampini   PetscCall(TSGetSolution(ts, &u));
1376*811d88abSStefano Zampini   PetscCall(VecViewFromOptions(u, NULL, "-final_view"));
1377*811d88abSStefano Zampini   PetscCall(VecViewFlux(u, "-final_flux_view"));
1378*811d88abSStefano Zampini 
1379*811d88abSStefano Zampini   PetscCall(TSDestroy(&ts));
1380*811d88abSStefano Zampini   PetscCall(VecTaggerDestroy(&actx->refineTag));
1381*811d88abSStefano Zampini   PetscCall(PetscFree(actx));
1382*811d88abSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1383*811d88abSStefano Zampini }
1384*811d88abSStefano Zampini 
1385*811d88abSStefano Zampini int main(int argc, char **argv)
1386*811d88abSStefano Zampini {
1387*811d88abSStefano Zampini   AppCtx ctx;
1388*811d88abSStefano Zampini 
1389*811d88abSStefano Zampini   PetscFunctionBeginUser;
1390*811d88abSStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1391*811d88abSStefano Zampini   PetscCall(ProcessOptions(&ctx));
1392*811d88abSStefano Zampini   PetscCall(PetscLogStageRegister("Setup", &SetupStage));
1393*811d88abSStefano Zampini   PetscCall(PetscLogStageRegister("Solve", &SolveStage));
1394*811d88abSStefano Zampini   if (ctx.test_restart) { /* Test sequences of save and loads */
1395*811d88abSStefano Zampini     PetscMPIInt rank;
1396*811d88abSStefano Zampini 
1397*811d88abSStefano Zampini     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1398*811d88abSStefano Zampini 
1399*811d88abSStefano Zampini     /* test saving */
1400*811d88abSStefano Zampini     ctx.load = PETSC_FALSE;
1401*811d88abSStefano Zampini     ctx.save = PETSC_TRUE;
1402*811d88abSStefano Zampini     /* sequential save */
1403*811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test sequential save\n"));
1404*811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.save_filename, sizeof(ctx.save_filename), "test_ex30_seq_%d.h5", rank));
1405*811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_SELF, &ctx));
1406*811d88abSStefano Zampini     /* parallel save */
1407*811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test parallel save\n"));
1408*811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.save_filename, sizeof(ctx.save_filename), "test_ex30_par.h5"));
1409*811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1410*811d88abSStefano Zampini 
1411*811d88abSStefano Zampini     /* test loading */
1412*811d88abSStefano Zampini     ctx.load = PETSC_TRUE;
1413*811d88abSStefano Zampini     ctx.save = PETSC_FALSE;
1414*811d88abSStefano Zampini     /* sequential and parallel runs from sequential save */
1415*811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.load_filename, sizeof(ctx.load_filename), "test_ex30_seq_0.h5"));
1416*811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test sequential load from sequential save\n"));
1417*811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_SELF, &ctx));
1418*811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test parallel load from sequential save\n"));
1419*811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1420*811d88abSStefano Zampini     /* sequential and parallel runs from parallel save */
1421*811d88abSStefano Zampini     PetscCall(PetscSNPrintf(ctx.load_filename, sizeof(ctx.load_filename), "test_ex30_par.h5"));
1422*811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test sequential load from parallel save\n"));
1423*811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_SELF, &ctx));
1424*811d88abSStefano Zampini     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test parallel load from parallel save\n"));
1425*811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1426*811d88abSStefano Zampini   } else { /* Run the simulation */
1427*811d88abSStefano Zampini     PetscCall(Run(PETSC_COMM_WORLD, &ctx));
1428*811d88abSStefano Zampini   }
1429*811d88abSStefano Zampini   PetscCall(PetscFinalize());
1430*811d88abSStefano Zampini   return 0;
1431*811d88abSStefano Zampini }
1432*811d88abSStefano Zampini 
1433*811d88abSStefano Zampini /*TEST
1434*811d88abSStefano Zampini 
1435*811d88abSStefano Zampini   testset:
1436*811d88abSStefano 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
1437*811d88abSStefano Zampini 
1438*811d88abSStefano Zampini     test:
1439*811d88abSStefano Zampini       suffix: 0
1440*811d88abSStefano Zampini       nsize: {{1 2}}
1441*811d88abSStefano Zampini       args: -dm_refine 1
1442*811d88abSStefano Zampini 
1443*811d88abSStefano Zampini     test:
1444*811d88abSStefano Zampini       suffix: 0_dirk
1445*811d88abSStefano Zampini       nsize: {{1 2}}
1446*811d88abSStefano Zampini       args: -dm_refine 1 -ts_type dirk
1447*811d88abSStefano Zampini 
1448*811d88abSStefano Zampini     test:
1449*811d88abSStefano Zampini       suffix: 0_dirk_mg
1450*811d88abSStefano Zampini       nsize: {{1 2}}
1451*811d88abSStefano 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
1452*811d88abSStefano Zampini 
1453*811d88abSStefano Zampini     test:
1454*811d88abSStefano Zampini       requires: p4est
1455*811d88abSStefano Zampini       suffix: 0_p4est
1456*811d88abSStefano Zampini       nsize: {{1 2}}
1457*811d88abSStefano Zampini       args: -dm_refine 1 -dm_plex_convert_type p4est
1458*811d88abSStefano Zampini 
1459*811d88abSStefano Zampini     test:
1460*811d88abSStefano Zampini       suffix: 0_periodic
1461*811d88abSStefano Zampini       nsize: {{1 2}}
1462*811d88abSStefano Zampini       args: -dm_plex_box_bd periodic,periodic -dm_refine_pre 1
1463*811d88abSStefano Zampini 
1464*811d88abSStefano Zampini     test:
1465*811d88abSStefano Zampini       requires: p4est
1466*811d88abSStefano Zampini       suffix: 0_p4est_periodic
1467*811d88abSStefano Zampini       nsize: {{1 2}}
1468*811d88abSStefano Zampini       args: -dm_plex_box_bd periodic,periodic -dm_refine_pre 1 -dm_plex_convert_type p4est
1469*811d88abSStefano Zampini 
1470*811d88abSStefano Zampini     test:
1471*811d88abSStefano Zampini       requires: p4est
1472*811d88abSStefano Zampini       suffix: 0_p4est_mg
1473*811d88abSStefano Zampini       nsize: {{1 2}}
1474*811d88abSStefano 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
1475*811d88abSStefano Zampini 
1476*811d88abSStefano Zampini   testset:
1477*811d88abSStefano Zampini     requires: hdf5
1478*811d88abSStefano 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
1479*811d88abSStefano Zampini 
1480*811d88abSStefano Zampini     test:
1481*811d88abSStefano Zampini       suffix: restart
1482*811d88abSStefano Zampini       nsize: {{1 2}separate output}
1483*811d88abSStefano Zampini       args: -dm_refine_hierarchy {{0 1}separate output} -dm_plex_simplex 0
1484*811d88abSStefano Zampini 
1485*811d88abSStefano Zampini     test:
1486*811d88abSStefano Zampini       requires: triangle
1487*811d88abSStefano Zampini       suffix: restart_simplex
1488*811d88abSStefano Zampini       nsize: {{1 2}separate output}
1489*811d88abSStefano Zampini       args: -dm_refine_hierarchy {{0 1}separate output} -dm_plex_simplex 1
1490*811d88abSStefano Zampini 
1491*811d88abSStefano Zampini     test:
1492*811d88abSStefano Zampini       suffix: restart_refonly
1493*811d88abSStefano Zampini       nsize: {{1 2}separate output}
1494*811d88abSStefano Zampini       args: -dm_refine 1 -dm_plex_simplex 0
1495*811d88abSStefano Zampini 
1496*811d88abSStefano Zampini     test:
1497*811d88abSStefano Zampini       requires: triangle
1498*811d88abSStefano Zampini       suffix: restart_simplex_refonly
1499*811d88abSStefano Zampini       nsize: {{1 2}separate output}
1500*811d88abSStefano Zampini       args: -dm_refine 1 -dm_plex_simplex 1
1501*811d88abSStefano Zampini 
1502*811d88abSStefano Zampini TEST*/
1503