xref: /petsc/src/snes/tutorials/ex64.c (revision 8f3c339631b3b5edcfac404dc77e5b2add68af7c)
1*b8cfa2afSStefano Zampini static char help[] = "Biot consolidation model discretized with finite elements,\n\
2*b8cfa2afSStefano Zampini using a parallel unstructured mesh (DMPLEX) to represent the domain.\n\
3*b8cfa2afSStefano Zampini We follow https://arxiv.org/pdf/1507.03199.\n\n\n";
4*b8cfa2afSStefano Zampini 
5*b8cfa2afSStefano Zampini #include <petscdmplex.h>
6*b8cfa2afSStefano Zampini #include <petscsnes.h>
7*b8cfa2afSStefano Zampini #include <petscds.h>
8*b8cfa2afSStefano Zampini 
9*b8cfa2afSStefano Zampini typedef struct {
10*b8cfa2afSStefano Zampini   PetscScalar mu;
11*b8cfa2afSStefano Zampini   PetscScalar lam;
12*b8cfa2afSStefano Zampini   PetscScalar alpha;
13*b8cfa2afSStefano Zampini   PetscScalar kappa;
14*b8cfa2afSStefano Zampini } Parameter;
15*b8cfa2afSStefano Zampini 
u_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 f[])16*b8cfa2afSStefano Zampini static void u_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 f[])
17*b8cfa2afSStefano Zampini {
18*b8cfa2afSStefano Zampini   const PetscReal mu = PetscRealPart(constants[0]);
19*b8cfa2afSStefano Zampini 
20*b8cfa2afSStefano Zampini   for (PetscInt c = 0; c < dim; ++c) {
21*b8cfa2afSStefano Zampini     for (PetscInt d = 0; d < dim; ++d) f[c * dim + d] = mu * (u_x[c * dim + d] + u_x[d * dim + c]);
22*b8cfa2afSStefano Zampini     f[c * dim + c] -= u[uOff[1]];
23*b8cfa2afSStefano Zampini   }
24*b8cfa2afSStefano Zampini }
25*b8cfa2afSStefano Zampini 
26*b8cfa2afSStefano Zampini /* Jfunction_testtrial */
Ju_1_u1p0(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[])27*b8cfa2afSStefano Zampini static void Ju_1_u1p0(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[])
28*b8cfa2afSStefano Zampini {
29*b8cfa2afSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) J[d * dim + d] = -1.0;
30*b8cfa2afSStefano Zampini }
31*b8cfa2afSStefano Zampini 
Ju_1_u1u1(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[])32*b8cfa2afSStefano Zampini static void Ju_1_u1u1(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[])
33*b8cfa2afSStefano Zampini {
34*b8cfa2afSStefano Zampini   const PetscReal mu = PetscRealPart(constants[0]);
35*b8cfa2afSStefano Zampini   const PetscInt  Nc = uOff[1] - uOff[0];
36*b8cfa2afSStefano Zampini 
37*b8cfa2afSStefano Zampini   for (PetscInt c = 0; c < Nc; ++c) {
38*b8cfa2afSStefano Zampini     for (PetscInt d = 0; d < dim; ++d) {
39*b8cfa2afSStefano Zampini       J[((c * Nc + c) * dim + d) * dim + d] += mu;
40*b8cfa2afSStefano Zampini       J[((c * Nc + d) * dim + d) * dim + c] += mu;
41*b8cfa2afSStefano Zampini     }
42*b8cfa2afSStefano Zampini   }
43*b8cfa2afSStefano Zampini }
44*b8cfa2afSStefano Zampini 
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 f[])45*b8cfa2afSStefano 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 f[])
46*b8cfa2afSStefano Zampini {
47*b8cfa2afSStefano Zampini   const PetscReal lam   = PetscRealPart(constants[1]);
48*b8cfa2afSStefano Zampini   const PetscReal alpha = PetscRealPart(constants[2]);
49*b8cfa2afSStefano Zampini 
50*b8cfa2afSStefano Zampini   f[0] = (-u[uOff[1]] + alpha * u[uOff[2]]) / lam;
51*b8cfa2afSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) f[0] -= u_x[d * dim + d];
52*b8cfa2afSStefano Zampini }
53*b8cfa2afSStefano Zampini 
Jp_0_p0u1(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[])54*b8cfa2afSStefano Zampini static void Jp_0_p0u1(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[])
55*b8cfa2afSStefano Zampini {
56*b8cfa2afSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) J[d * dim + d] = -1.0;
57*b8cfa2afSStefano Zampini }
58*b8cfa2afSStefano Zampini 
Jp_0_p0p0(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[])59*b8cfa2afSStefano Zampini static void Jp_0_p0p0(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[])
60*b8cfa2afSStefano Zampini {
61*b8cfa2afSStefano Zampini   const PetscReal lam = PetscRealPart(constants[1]);
62*b8cfa2afSStefano Zampini 
63*b8cfa2afSStefano Zampini   J[0] = -1.0 / lam;
64*b8cfa2afSStefano Zampini }
65*b8cfa2afSStefano Zampini 
pJp_0_p0p0(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[])66*b8cfa2afSStefano Zampini static void pJp_0_p0p0(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[])
67*b8cfa2afSStefano Zampini {
68*b8cfa2afSStefano Zampini   const PetscReal mu = PetscRealPart(constants[0]);
69*b8cfa2afSStefano Zampini 
70*b8cfa2afSStefano Zampini   J[0] = 1.0 / mu;
71*b8cfa2afSStefano Zampini }
72*b8cfa2afSStefano Zampini 
Jp_0_p0w0(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[])73*b8cfa2afSStefano Zampini static void Jp_0_p0w0(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[])
74*b8cfa2afSStefano Zampini {
75*b8cfa2afSStefano Zampini   const PetscReal lam   = PetscRealPart(constants[1]);
76*b8cfa2afSStefano Zampini   const PetscReal alpha = PetscRealPart(constants[2]);
77*b8cfa2afSStefano Zampini 
78*b8cfa2afSStefano Zampini   J[0] = alpha / lam;
79*b8cfa2afSStefano Zampini }
80*b8cfa2afSStefano Zampini 
w_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 f[])81*b8cfa2afSStefano Zampini static void w_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 f[])
82*b8cfa2afSStefano Zampini {
83*b8cfa2afSStefano Zampini   const PetscReal lam   = PetscRealPart(constants[1]);
84*b8cfa2afSStefano Zampini   const PetscReal alpha = PetscRealPart(constants[2]);
85*b8cfa2afSStefano Zampini 
86*b8cfa2afSStefano Zampini   f[0] = alpha / lam * (-2 * alpha * u[uOff[2]] + u[uOff[1]]);
87*b8cfa2afSStefano Zampini }
88*b8cfa2afSStefano Zampini 
Jw_0_w0p0(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[])89*b8cfa2afSStefano Zampini static void Jw_0_w0p0(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[])
90*b8cfa2afSStefano Zampini {
91*b8cfa2afSStefano Zampini   const PetscReal lam   = PetscRealPart(constants[1]);
92*b8cfa2afSStefano Zampini   const PetscReal alpha = PetscRealPart(constants[2]);
93*b8cfa2afSStefano Zampini 
94*b8cfa2afSStefano Zampini   J[0] = alpha / lam;
95*b8cfa2afSStefano Zampini }
96*b8cfa2afSStefano Zampini 
Jw_0_w0w0(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[])97*b8cfa2afSStefano Zampini static void Jw_0_w0w0(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*b8cfa2afSStefano Zampini {
99*b8cfa2afSStefano Zampini   const PetscReal lam   = PetscRealPart(constants[1]);
100*b8cfa2afSStefano Zampini   const PetscReal alpha = PetscRealPart(constants[2]);
101*b8cfa2afSStefano Zampini 
102*b8cfa2afSStefano Zampini   J[0] = -2 * PetscSqr(alpha) / lam;
103*b8cfa2afSStefano Zampini }
104*b8cfa2afSStefano Zampini 
pJw_0_w0w0(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[])105*b8cfa2afSStefano Zampini static void pJw_0_w0w0(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[])
106*b8cfa2afSStefano Zampini {
107*b8cfa2afSStefano Zampini   const PetscReal lam   = PetscRealPart(constants[1]);
108*b8cfa2afSStefano Zampini   const PetscReal alpha = PetscRealPart(constants[2]);
109*b8cfa2afSStefano Zampini 
110*b8cfa2afSStefano Zampini   J[0] = 2 * PetscSqr(alpha) / lam;
111*b8cfa2afSStefano Zampini }
112*b8cfa2afSStefano Zampini 
w_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 f[])113*b8cfa2afSStefano Zampini static void w_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 f[])
114*b8cfa2afSStefano Zampini {
115*b8cfa2afSStefano Zampini   const PetscReal kappa = PetscRealPart(constants[3]);
116*b8cfa2afSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) f[d] = -kappa * u_x[uOff_x[2] + d];
117*b8cfa2afSStefano Zampini }
118*b8cfa2afSStefano Zampini 
Jw_1_w1w1(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[])119*b8cfa2afSStefano Zampini static void Jw_1_w1w1(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*b8cfa2afSStefano Zampini {
121*b8cfa2afSStefano Zampini   const PetscReal kappa = PetscRealPart(constants[3]);
122*b8cfa2afSStefano Zampini 
123*b8cfa2afSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) J[d * dim + d] = -kappa;
124*b8cfa2afSStefano Zampini }
125*b8cfa2afSStefano Zampini 
pJw_1_w1w1(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[])126*b8cfa2afSStefano Zampini static void pJw_1_w1w1(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[])
127*b8cfa2afSStefano Zampini {
128*b8cfa2afSStefano Zampini   const PetscReal kappa = PetscRealPart(constants[3]);
129*b8cfa2afSStefano Zampini 
130*b8cfa2afSStefano Zampini   for (PetscInt d = 0; d < dim; ++d) J[d * dim + d] = kappa;
131*b8cfa2afSStefano Zampini }
132*b8cfa2afSStefano Zampini 
133*b8cfa2afSStefano Zampini typedef struct {
134*b8cfa2afSStefano Zampini   PetscScalar E;
135*b8cfa2afSStefano Zampini   PetscScalar nu;
136*b8cfa2afSStefano Zampini   PetscScalar alpha;
137*b8cfa2afSStefano Zampini   PetscScalar kappa;
138*b8cfa2afSStefano Zampini } AppCtx;
139*b8cfa2afSStefano Zampini 
ProcessOptions(MPI_Comm comm,AppCtx * user)140*b8cfa2afSStefano Zampini static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *user)
141*b8cfa2afSStefano Zampini {
142*b8cfa2afSStefano Zampini   PetscFunctionBeginUser;
143*b8cfa2afSStefano Zampini   user->E     = 1.0;
144*b8cfa2afSStefano Zampini   user->nu    = 0.3;
145*b8cfa2afSStefano Zampini   user->alpha = 0.5;
146*b8cfa2afSStefano Zampini   user->kappa = 1.0;
147*b8cfa2afSStefano Zampini   PetscOptionsBegin(comm, "", "Biot Problem Options", "DMPLEX");
148*b8cfa2afSStefano Zampini   PetscCall(PetscOptionsScalar("-E", "Young modulus", NULL, user->E, &user->E, NULL));
149*b8cfa2afSStefano Zampini   PetscCall(PetscOptionsScalar("-nu", "Poisson ratio", NULL, user->nu, &user->nu, NULL));
150*b8cfa2afSStefano Zampini   PetscCall(PetscOptionsScalar("-alpha", "Alpha", NULL, user->alpha, &user->alpha, NULL));
151*b8cfa2afSStefano Zampini   PetscCall(PetscOptionsScalar("-kappa", "kappa", NULL, user->kappa, &user->kappa, NULL));
152*b8cfa2afSStefano Zampini   PetscOptionsEnd();
153*b8cfa2afSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
154*b8cfa2afSStefano Zampini }
155*b8cfa2afSStefano Zampini 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)156*b8cfa2afSStefano Zampini static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
157*b8cfa2afSStefano Zampini {
158*b8cfa2afSStefano Zampini   PetscFunctionBeginUser;
159*b8cfa2afSStefano Zampini   PetscCall(DMCreate(comm, dm));
160*b8cfa2afSStefano Zampini   PetscCall(DMSetType(*dm, DMPLEX));
161*b8cfa2afSStefano Zampini   PetscCall(DMSetFromOptions(*dm));
162*b8cfa2afSStefano Zampini   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
163*b8cfa2afSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
164*b8cfa2afSStefano Zampini }
165*b8cfa2afSStefano Zampini 
SetupEqn(DM dm,AppCtx * user)166*b8cfa2afSStefano Zampini static PetscErrorCode SetupEqn(DM dm, AppCtx *user)
167*b8cfa2afSStefano Zampini {
168*b8cfa2afSStefano Zampini   PetscDS        ds;
169*b8cfa2afSStefano Zampini   DMLabel        label;
170*b8cfa2afSStefano Zampini   const PetscInt id = 1;
171*b8cfa2afSStefano Zampini   PetscScalar    constants[4];
172*b8cfa2afSStefano Zampini 
173*b8cfa2afSStefano Zampini   PetscFunctionBeginUser;
174*b8cfa2afSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
175*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetResidual(ds, 0, NULL, u_1));
176*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetResidual(ds, 1, p_0, NULL));
177*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetResidual(ds, 2, w_0, w_1));
178*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, Ju_1_u1u1));
179*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, Ju_1_u1p0, NULL));
180*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, Jp_0_p0u1, NULL, NULL));
181*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 1, 1, Jp_0_p0p0, NULL, NULL, NULL));
182*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 1, 2, Jp_0_p0w0, NULL, NULL, NULL));
183*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 2, 1, Jw_0_w0p0, NULL, NULL, NULL));
184*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobian(ds, 2, 2, Jw_0_w0w0, NULL, NULL, Jw_1_w1w1));
185*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, Ju_1_u1u1));
186*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobianPreconditioner(ds, 0, 1, NULL, NULL, Ju_1_u1p0, NULL));
187*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobianPreconditioner(ds, 1, 0, NULL, Jp_0_p0u1, NULL, NULL));
188*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobianPreconditioner(ds, 1, 1, pJp_0_p0p0, NULL, NULL, NULL));
189*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetJacobianPreconditioner(ds, 2, 2, pJw_0_w0w0, NULL, NULL, pJw_1_w1w1));
190*b8cfa2afSStefano Zampini 
191*b8cfa2afSStefano Zampini   PetscCall(DMGetLabel(dm, "marker", &label));
192*b8cfa2afSStefano Zampini   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall-d", label, 1, &id, 0, 0, NULL, NULL, NULL, NULL, NULL));
193*b8cfa2afSStefano Zampini   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall-p", label, 1, &id, 2, 0, NULL, NULL, NULL, NULL, NULL));
194*b8cfa2afSStefano Zampini 
195*b8cfa2afSStefano Zampini   constants[0] = user->E / (2.0 * (1.0 + user->nu));
196*b8cfa2afSStefano Zampini   constants[1] = user->nu * user->E / ((1.0 + user->nu) * (1.0 - 2.0 * user->nu));
197*b8cfa2afSStefano Zampini   constants[2] = user->alpha;
198*b8cfa2afSStefano Zampini   constants[3] = user->kappa;
199*b8cfa2afSStefano Zampini   PetscCall(PetscDSSetConstants(ds, 4, constants));
200*b8cfa2afSStefano Zampini 
201*b8cfa2afSStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "E = %g, nu = %g\n", (double)PetscRealPart(user->E), (double)PetscRealPart(user->nu)));
202*b8cfa2afSStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "mu = %g, lambda = %g\n", (double)PetscRealPart(constants[0]), (double)PetscRealPart(constants[1])));
203*b8cfa2afSStefano Zampini   PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "alpha = %g, kappa = %g\n", (double)PetscRealPart(constants[2]), (double)PetscRealPart(constants[3])));
204*b8cfa2afSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
205*b8cfa2afSStefano Zampini }
206*b8cfa2afSStefano Zampini 
SetupProblem(DM dm,PetscErrorCode (* setupEqn)(DM,AppCtx *),AppCtx * user)207*b8cfa2afSStefano Zampini static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user)
208*b8cfa2afSStefano Zampini {
209*b8cfa2afSStefano Zampini   DM              cdm = dm;
210*b8cfa2afSStefano Zampini   PetscQuadrature q   = NULL;
211*b8cfa2afSStefano Zampini   PetscBool       simplex;
212*b8cfa2afSStefano Zampini   PetscInt        dim, Nf = 3, f, Nc[3];
213*b8cfa2afSStefano Zampini   const char     *name[3]   = {"displacement", "totalpressure", "pressure"};
214*b8cfa2afSStefano Zampini   const char     *prefix[3] = {"displacement_", "totalpressure_", "pressure_"};
215*b8cfa2afSStefano Zampini 
216*b8cfa2afSStefano Zampini   PetscFunctionBegin;
217*b8cfa2afSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
218*b8cfa2afSStefano Zampini   PetscCall(DMPlexIsSimplex(dm, &simplex));
219*b8cfa2afSStefano Zampini   Nc[0] = dim;
220*b8cfa2afSStefano Zampini   Nc[1] = 1;
221*b8cfa2afSStefano Zampini   Nc[2] = 1;
222*b8cfa2afSStefano Zampini   for (f = 0; f < Nf; ++f) {
223*b8cfa2afSStefano Zampini     PetscFE fe;
224*b8cfa2afSStefano Zampini 
225*b8cfa2afSStefano Zampini     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe));
226*b8cfa2afSStefano Zampini     PetscCall(PetscObjectSetName((PetscObject)fe, name[f]));
227*b8cfa2afSStefano Zampini     if (!q) PetscCall(PetscFEGetQuadrature(fe, &q));
228*b8cfa2afSStefano Zampini     PetscCall(PetscFESetQuadrature(fe, q));
229*b8cfa2afSStefano Zampini     PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe));
230*b8cfa2afSStefano Zampini     PetscCall(PetscFEDestroy(&fe));
231*b8cfa2afSStefano Zampini   }
232*b8cfa2afSStefano Zampini   PetscCall(DMCreateDS(dm));
233*b8cfa2afSStefano Zampini   PetscCall((*setupEqn)(dm, user));
234*b8cfa2afSStefano Zampini   while (cdm) {
235*b8cfa2afSStefano Zampini     PetscCall(DMCopyDisc(dm, cdm));
236*b8cfa2afSStefano Zampini     PetscCall(DMGetCoarseDM(cdm, &cdm));
237*b8cfa2afSStefano Zampini   }
238*b8cfa2afSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
239*b8cfa2afSStefano Zampini }
240*b8cfa2afSStefano Zampini 
main(int argc,char ** argv)241*b8cfa2afSStefano Zampini int main(int argc, char **argv)
242*b8cfa2afSStefano Zampini {
243*b8cfa2afSStefano Zampini   SNES   snes;
244*b8cfa2afSStefano Zampini   DM     dm;
245*b8cfa2afSStefano Zampini   Vec    u;
246*b8cfa2afSStefano Zampini   AppCtx user;
247*b8cfa2afSStefano Zampini 
248*b8cfa2afSStefano Zampini   PetscFunctionBeginUser;
249*b8cfa2afSStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
250*b8cfa2afSStefano Zampini   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
251*b8cfa2afSStefano Zampini   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
252*b8cfa2afSStefano Zampini   PetscCall(SetupProblem(dm, SetupEqn, &user));
253*b8cfa2afSStefano Zampini   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
254*b8cfa2afSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
255*b8cfa2afSStefano Zampini   PetscCall(DMSetApplicationContext(dm, &user));
256*b8cfa2afSStefano Zampini 
257*b8cfa2afSStefano Zampini   PetscCall(SNESCreate(PetscObjectComm((PetscObject)dm), &snes));
258*b8cfa2afSStefano Zampini   PetscCall(SNESSetDM(snes, dm));
259*b8cfa2afSStefano Zampini   PetscCall(SNESSetType(snes, SNESKSPONLY));
260*b8cfa2afSStefano Zampini   PetscCall(SNESSetFromOptions(snes));
261*b8cfa2afSStefano Zampini 
262*b8cfa2afSStefano Zampini   /* Solve with random rhs */
263*b8cfa2afSStefano Zampini   PetscCall(DMCreateGlobalVector(dm, &u));
264*b8cfa2afSStefano Zampini   PetscCall(VecSetRandom(u, NULL));
265*b8cfa2afSStefano Zampini   PetscCall(SNESSolve(snes, NULL, u));
266*b8cfa2afSStefano Zampini 
267*b8cfa2afSStefano Zampini   PetscCall(VecDestroy(&u));
268*b8cfa2afSStefano Zampini   PetscCall(SNESDestroy(&snes));
269*b8cfa2afSStefano Zampini   PetscCall(DMDestroy(&dm));
270*b8cfa2afSStefano Zampini   PetscCall(PetscFinalize());
271*b8cfa2afSStefano Zampini   return 0;
272*b8cfa2afSStefano Zampini }
273*b8cfa2afSStefano Zampini 
274*b8cfa2afSStefano Zampini /*TEST
275*b8cfa2afSStefano Zampini 
276*b8cfa2afSStefano Zampini   test:
277*b8cfa2afSStefano Zampini     suffix: 2d_p2_p1_p2
278*b8cfa2afSStefano Zampini     args: -displacement_petscspace_degree 2 -totalpressure_petscspace_degree 1 -pressure_petscspace_degree 2 \
279*b8cfa2afSStefano Zampini           -dm_plex_box_faces 5,5 -dm_plex_separate_marker -dm_plex_simplex 0 \
280*b8cfa2afSStefano Zampini           -snes_error_if_not_converged -ksp_error_if_not_converged \
281*b8cfa2afSStefano Zampini           -ksp_type minres -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_totalpressure_pc_type jacobi -fieldsplit_displacement_pc_type gamg -fieldsplit_pressure_pc_type gamg
282*b8cfa2afSStefano Zampini 
283*b8cfa2afSStefano Zampini   test:
284*b8cfa2afSStefano Zampini     nsize: 4
285*b8cfa2afSStefano Zampini     suffix: 2d_p2_p1_p2_fetidp
286*b8cfa2afSStefano Zampini     args: -displacement_petscspace_degree 2 -totalpressure_petscspace_degree 1 -pressure_petscspace_degree 2 \
287*b8cfa2afSStefano Zampini           -petscpartitioner_type simple -dm_mat_type is -dm_plex_box_faces 3,3 -dm_refine 2 -dm_plex_separate_marker -dm_plex_simplex 0 \
288*b8cfa2afSStefano Zampini           -snes_error_if_not_converged -ksp_error_if_not_converged \
289*b8cfa2afSStefano Zampini           -ksp_type fetidp -ksp_fetidp_saddlepoint -ksp_fetidp_pressure_field 1,2 -ksp_fetidp_pressure_schur -fetidp_ksp_type cg -fetidp_ksp_norm_type natural \
290*b8cfa2afSStefano Zampini           -fetidp_fieldsplit_lag_ksp_type preonly \
291*b8cfa2afSStefano Zampini           -fetidp_bddc_pc_bddc_detect_disconnected -fetidp_bddc_pc_bddc_corner_selection -fetidp_bddc_pc_bddc_symmetric -fetidp_bddc_pc_bddc_coarse_redundant_pc_type cholesky \
292*b8cfa2afSStefano Zampini           -fetidp_pc_discrete_harmonic 1 -fetidp_harmonic_pc_type cholesky \
293*b8cfa2afSStefano Zampini           -fetidp_fieldsplit_p_pc_type bddc -fetidp_fieldsplit_p_ksp_type preonly
294*b8cfa2afSStefano Zampini 
295*b8cfa2afSStefano Zampini TEST*/
296