xref: /petsc/src/ts/tutorials/ex53.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
165876a83SMatthew G. Knepley static char help[] = "Time dependent Biot Poroelasticity problem with finite elements.\n\
265876a83SMatthew G. Knepley We solve three field, quasi-static poroelasticity in a rectangular\n\
365876a83SMatthew G. Knepley domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
465876a83SMatthew G. Knepley Contributed by: Robert Walker <rwalker6@buffalo.edu>\n\n\n";
565876a83SMatthew G. Knepley 
665876a83SMatthew G. Knepley #include <petscdmplex.h>
765876a83SMatthew G. Knepley #include <petscsnes.h>
865876a83SMatthew G. Knepley #include <petscts.h>
965876a83SMatthew G. Knepley #include <petscds.h>
1065876a83SMatthew G. Knepley #include <petscbag.h>
1165876a83SMatthew G. Knepley 
1265876a83SMatthew G. Knepley #include <petsc/private/tsimpl.h>
1365876a83SMatthew G. Knepley 
1465876a83SMatthew G. Knepley /* This presentation of poroelasticity is taken from
1565876a83SMatthew G. Knepley 
1665876a83SMatthew G. Knepley @book{Cheng2016,
1765876a83SMatthew G. Knepley   title={Poroelasticity},
1865876a83SMatthew G. Knepley   author={Cheng, Alexander H-D},
1965876a83SMatthew G. Knepley   volume={27},
2065876a83SMatthew G. Knepley   year={2016},
2165876a83SMatthew G. Knepley   publisher={Springer}
2265876a83SMatthew G. Knepley }
2365876a83SMatthew G. Knepley 
2465876a83SMatthew G. Knepley For visualization, use
2565876a83SMatthew G. Knepley 
2665876a83SMatthew G. Knepley   -dm_view hdf5:${PETSC_DIR}/sol.h5 -monitor_solution hdf5:${PETSC_DIR}/sol.h5::append
2765876a83SMatthew G. Knepley 
2865876a83SMatthew G. Knepley The weak form would then be, using test function $(v, q, \tau)$,
2965876a83SMatthew G. Knepley 
3065876a83SMatthew G. Knepley             (q, \frac{1}{M} \frac{dp}{dt}) + (q, \alpha \frac{d\varepsilon}{dt}) + (\nabla q, \kappa \nabla p) = (q, g)
3165876a83SMatthew G. Knepley  -(\nabla v, 2 G \epsilon) - (\nabla\cdot v, \frac{2 G \nu}{1 - 2\nu} \varepsilon) + \alpha (\nabla\cdot v, p) = (v, f)
3265876a83SMatthew G. Knepley                                                                           (\tau, \nabla \cdot u - \varepsilon) = 0
3365876a83SMatthew G. Knepley */
3465876a83SMatthew G. Knepley 
3565876a83SMatthew G. Knepley typedef enum {SOL_QUADRATIC_LINEAR, SOL_QUADRATIC_TRIG, SOL_TRIG_LINEAR, SOL_TERZAGHI, SOL_MANDEL, SOL_CRYER, NUM_SOLUTION_TYPES} SolutionType;
3665876a83SMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES+1] = {"quadratic_linear", "quadratic_trig", "trig_linear", "terzaghi", "mandel", "cryer", "unknown"};
3765876a83SMatthew G. Knepley 
3865876a83SMatthew G. Knepley typedef struct {
3965876a83SMatthew G. Knepley   PetscScalar mu;    /* shear modulus */
4065876a83SMatthew G. Knepley   PetscScalar K_u;   /* undrained bulk modulus */
4165876a83SMatthew G. Knepley   PetscScalar alpha; /* Biot effective stress coefficient */
4265876a83SMatthew G. Knepley   PetscScalar M;     /* Biot modulus */
4365876a83SMatthew G. Knepley   PetscScalar k;     /* (isotropic) permeability */
4465876a83SMatthew G. Knepley   PetscScalar mu_f;  /* fluid dynamic viscosity */
4565876a83SMatthew G. Knepley   PetscScalar P_0;   /* magnitude of vertical stress */
4665876a83SMatthew G. Knepley } Parameter;
4765876a83SMatthew G. Knepley 
4865876a83SMatthew G. Knepley typedef struct {
4965876a83SMatthew G. Knepley   /* Domain and mesh definition */
5030602db0SMatthew G. Knepley   PetscReal    xmin[3];     /* Lower left bottom corner of bounding box */
5130602db0SMatthew G. Knepley   PetscReal    xmax[3];     /* Upper right top corner of bounding box */
5265876a83SMatthew G. Knepley   /* Problem definition */
5365876a83SMatthew G. Knepley   SolutionType solType;     /* Type of exact solution */
5465876a83SMatthew G. Knepley   PetscBag     bag;         /* Problem parameters */
5565876a83SMatthew G. Knepley   PetscReal    t_r;         /* Relaxation time: 4 L^2 / c */
5624b15d09SMatthew G. Knepley   PetscReal    dtInitial;   /* Override the choice for first timestep */
5765876a83SMatthew G. Knepley   /* Exact solution terms */
5865876a83SMatthew G. Knepley   PetscInt     niter; /* Number of series term iterations in exact solutions */
5965876a83SMatthew G. Knepley   PetscReal    eps;   /* Precision value for root finding */
6065876a83SMatthew G. Knepley   PetscReal   *zeroArray; /* Array of root locations */
6165876a83SMatthew G. Knepley } AppCtx;
6265876a83SMatthew G. Knepley 
6365876a83SMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
6465876a83SMatthew G. Knepley {
6565876a83SMatthew G. Knepley   PetscInt c;
6665876a83SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
6765876a83SMatthew G. Knepley   return 0;
6865876a83SMatthew G. Knepley }
6965876a83SMatthew G. Knepley 
7065876a83SMatthew G. Knepley /* Quadratic space and linear time solution
7165876a83SMatthew G. Knepley 
7265876a83SMatthew G. Knepley   2D:
7365876a83SMatthew G. Knepley   u = x^2
7465876a83SMatthew G. Knepley   v = y^2 - 2xy
7565876a83SMatthew G. Knepley   p = (x + y) t
7665876a83SMatthew G. Knepley   e = 2y
7765876a83SMatthew G. Knepley   f = <2 G, 4 G + 2 \lambda > - <alpha t, alpha t>
7865876a83SMatthew G. Knepley   g = 0
7965876a83SMatthew G. Knepley   \epsilon = / 2x     -y    \
8065876a83SMatthew G. Knepley              \ -y   2y - 2x /
8165876a83SMatthew G. Knepley   Tr(\epsilon) = e = div u = 2y
8265876a83SMatthew G. Knepley   div \sigma = \partial_i 2 G \epsilon_{ij} + \partial_i \lambda \varepsilon \delta_{ij} - \partial_i \alpha p \delta_{ij}
8365876a83SMatthew G. Knepley     = 2 G < 2-1, 2 > + \lambda <0, 2> - alpha <t, t>
8465876a83SMatthew G. Knepley     = <2 G, 4 G + 2 \lambda> - <alpha t, alpha t>
8565876a83SMatthew G. Knepley   \frac{1}{M} \frac{dp}{dt} + \alpha \frac{d\varepsilon}{dt} - \nabla \cdot \kappa \nabla p
8665876a83SMatthew G. Knepley     = \frac{1}{M} \frac{dp}{dt} + \kappa \Delta p
8765876a83SMatthew G. Knepley     = (x + y)/M
8865876a83SMatthew G. Knepley 
8965876a83SMatthew G. Knepley   3D:
9065876a83SMatthew G. Knepley   u = x^2
9165876a83SMatthew G. Knepley   v = y^2 - 2xy
9265876a83SMatthew G. Knepley   w = z^2 - 2yz
9365876a83SMatthew G. Knepley   p = (x + y + z) t
9465876a83SMatthew G. Knepley   e = 2z
9565876a83SMatthew G. Knepley   f = <2 G, 4 G + 2 \lambda > - <alpha t, alpha t, alpha t>
9665876a83SMatthew G. Knepley   g = 0
9765876a83SMatthew G. Knepley   \varepsilon = / 2x     -y       0   \
9865876a83SMatthew G. Knepley                 | -y   2y - 2x   -z   |
9965876a83SMatthew G. Knepley                 \  0     -z    2z - 2y/
10065876a83SMatthew G. Knepley   Tr(\varepsilon) = div u = 2z
10165876a83SMatthew G. Knepley   div \sigma = \partial_i 2 G \epsilon_{ij} + \partial_i \lambda \varepsilon \delta_{ij} - \partial_i \alpha p \delta_{ij}
10265876a83SMatthew G. Knepley     = 2 G < 2-1, 2-1, 2 > + \lambda <0, 0, 2> - alpha <t, t, t>
10365876a83SMatthew G. Knepley     = <2 G, 2G, 4 G + 2 \lambda> - <alpha t, alpha t, alpha t>
10465876a83SMatthew G. Knepley */
10565876a83SMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
10665876a83SMatthew G. Knepley {
10765876a83SMatthew G. Knepley   PetscInt d;
10865876a83SMatthew G. Knepley 
10965876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
11065876a83SMatthew G. Knepley     u[d] = PetscSqr(x[d]) - (d > 0 ? 2.0 * x[d-1] * x[d] : 0.0);
11165876a83SMatthew G. Knepley   }
11265876a83SMatthew G. Knepley   return 0;
11365876a83SMatthew G. Knepley }
11465876a83SMatthew G. Knepley 
11565876a83SMatthew G. Knepley static PetscErrorCode linear_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
11665876a83SMatthew G. Knepley {
11765876a83SMatthew G. Knepley   u[0] = 2.0*x[dim-1];
11865876a83SMatthew G. Knepley   return 0;
11965876a83SMatthew G. Knepley }
12065876a83SMatthew G. Knepley 
12165876a83SMatthew G. Knepley static PetscErrorCode linear_linear_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
12265876a83SMatthew G. Knepley {
12365876a83SMatthew G. Knepley   PetscReal sum = 0.0;
12465876a83SMatthew G. Knepley   PetscInt  d;
12565876a83SMatthew G. Knepley 
12665876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += x[d];
12765876a83SMatthew G. Knepley   u[0] = sum*time;
12865876a83SMatthew G. Knepley   return 0;
12965876a83SMatthew G. Knepley }
13065876a83SMatthew G. Knepley 
13165876a83SMatthew G. Knepley static PetscErrorCode linear_linear_p_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
13265876a83SMatthew G. Knepley {
13365876a83SMatthew G. Knepley   PetscReal sum = 0.0;
13465876a83SMatthew G. Knepley   PetscInt  d;
13565876a83SMatthew G. Knepley 
13665876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += x[d];
13765876a83SMatthew G. Knepley   u[0] = sum;
13865876a83SMatthew G. Knepley   return 0;
13965876a83SMatthew G. Knepley }
14065876a83SMatthew G. Knepley 
14165876a83SMatthew G. Knepley static void f0_quadratic_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
14265876a83SMatthew G. Knepley                                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
14365876a83SMatthew G. Knepley                                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
14465876a83SMatthew G. Knepley                                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
14565876a83SMatthew G. Knepley {
14665876a83SMatthew G. Knepley   const PetscReal G      = PetscRealPart(constants[0]);
14765876a83SMatthew G. Knepley   const PetscReal K_u    = PetscRealPart(constants[1]);
14865876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
14965876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
15065876a83SMatthew G. Knepley   const PetscReal K_d    = K_u - alpha*alpha*M;
15165876a83SMatthew G. Knepley   const PetscReal lambda = K_d - (2.0 * G) / 3.0;
15265876a83SMatthew G. Knepley   PetscInt        d;
15365876a83SMatthew G. Knepley 
15465876a83SMatthew G. Knepley   for (d = 0; d < dim-1; ++d) {
15565876a83SMatthew G. Knepley     f0[d] -= 2.0*G - alpha*t;
15665876a83SMatthew G. Knepley   }
15765876a83SMatthew G. Knepley   f0[dim-1] -= 2.0*lambda + 4.0*G - alpha*t;
15865876a83SMatthew G. Knepley }
15965876a83SMatthew G. Knepley 
16065876a83SMatthew G. Knepley static void f0_quadratic_linear_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
16165876a83SMatthew G. Knepley                                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
16265876a83SMatthew G. Knepley                                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
16365876a83SMatthew G. Knepley                                   PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
16465876a83SMatthew G. Knepley {
16565876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
16665876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
16765876a83SMatthew G. Knepley   PetscReal       sum    = 0.0;
16865876a83SMatthew G. Knepley   PetscInt        d;
16965876a83SMatthew G. Knepley 
17065876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += x[d];
17165876a83SMatthew G. Knepley   f0[0] += u_t ? alpha*u_t[uOff[1]] : 0.0;
17265876a83SMatthew G. Knepley   f0[0] += u_t ? u_t[uOff[2]]/M     : 0.0;
17365876a83SMatthew G. Knepley   f0[0] -= sum/M;
17465876a83SMatthew G. Knepley }
17565876a83SMatthew G. Knepley 
17665876a83SMatthew G. Knepley /* Quadratic space and trigonometric time solution
17765876a83SMatthew G. Knepley 
17865876a83SMatthew G. Knepley   2D:
17965876a83SMatthew G. Knepley   u = x^2
18065876a83SMatthew G. Knepley   v = y^2 - 2xy
18165876a83SMatthew G. Knepley   p = (x + y) cos(t)
18265876a83SMatthew G. Knepley   e = 2y
18365876a83SMatthew G. Knepley   f = <2 G, 4 G + 2 \lambda > - <alpha cos(t), alpha cos(t)>
18465876a83SMatthew G. Knepley   g = 0
18565876a83SMatthew G. Knepley   \epsilon = / 2x     -y    \
18665876a83SMatthew G. Knepley              \ -y   2y - 2x /
18765876a83SMatthew G. Knepley   Tr(\epsilon) = e = div u = 2y
18865876a83SMatthew G. Knepley   div \sigma = \partial_i 2 G \epsilon_{ij} + \partial_i \lambda \varepsilon \delta_{ij} - \partial_i \alpha p \delta_{ij}
18965876a83SMatthew G. Knepley     = 2 G < 2-1, 2 > + \lambda <0, 2> - alpha <cos(t), cos(t)>
19065876a83SMatthew G. Knepley     = <2 G, 4 G + 2 \lambda> - <alpha cos(t), alpha cos(t)>
19165876a83SMatthew G. Knepley   \frac{1}{M} \frac{dp}{dt} + \alpha \frac{d\varepsilon}{dt} - \nabla \cdot \kappa \nabla p
19265876a83SMatthew G. Knepley     = \frac{1}{M} \frac{dp}{dt} + \kappa \Delta p
19365876a83SMatthew G. Knepley     = -(x + y)/M sin(t)
19465876a83SMatthew G. Knepley 
19565876a83SMatthew G. Knepley   3D:
19665876a83SMatthew G. Knepley   u = x^2
19765876a83SMatthew G. Knepley   v = y^2 - 2xy
19865876a83SMatthew G. Knepley   w = z^2 - 2yz
19965876a83SMatthew G. Knepley   p = (x + y + z) cos(t)
20065876a83SMatthew G. Knepley   e = 2z
20165876a83SMatthew G. Knepley   f = <2 G, 4 G + 2 \lambda > - <alpha cos(t), alpha cos(t), alpha cos(t)>
20265876a83SMatthew G. Knepley   g = 0
20365876a83SMatthew G. Knepley   \varepsilon = / 2x     -y       0   \
20465876a83SMatthew G. Knepley                 | -y   2y - 2x   -z   |
20565876a83SMatthew G. Knepley                 \  0     -z    2z - 2y/
20665876a83SMatthew G. Knepley   Tr(\varepsilon) = div u = 2z
20765876a83SMatthew G. Knepley   div \sigma = \partial_i 2 G \epsilon_{ij} + \partial_i \lambda \varepsilon \delta_{ij} - \partial_i \alpha p \delta_{ij}
20865876a83SMatthew G. Knepley     = 2 G < 2-1, 2-1, 2 > + \lambda <0, 0, 2> - alpha <cos(t), cos(t), cos(t)>
20965876a83SMatthew G. Knepley     = <2 G, 2G, 4 G + 2 \lambda> - <alpha cos(t), alpha cos(t), alpha cos(t)>
21065876a83SMatthew G. Knepley */
21165876a83SMatthew G. Knepley static PetscErrorCode linear_trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
21265876a83SMatthew G. Knepley {
21365876a83SMatthew G. Knepley   PetscReal sum = 0.0;
21465876a83SMatthew G. Knepley   PetscInt  d;
21565876a83SMatthew G. Knepley 
21665876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += x[d];
21765876a83SMatthew G. Knepley   u[0] = sum*PetscCosReal(time);
21865876a83SMatthew G. Knepley   return 0;
21965876a83SMatthew G. Knepley }
22065876a83SMatthew G. Knepley 
22165876a83SMatthew G. Knepley static PetscErrorCode linear_trig_p_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
22265876a83SMatthew G. Knepley {
22365876a83SMatthew G. Knepley   PetscReal sum = 0.0;
22465876a83SMatthew G. Knepley   PetscInt  d;
22565876a83SMatthew G. Knepley 
22665876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += x[d];
22765876a83SMatthew G. Knepley   u[0] = -sum*PetscSinReal(time);
22865876a83SMatthew G. Knepley   return 0;
22965876a83SMatthew G. Knepley }
23065876a83SMatthew G. Knepley 
23165876a83SMatthew G. Knepley static void f0_quadratic_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
23265876a83SMatthew G. Knepley                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
23365876a83SMatthew G. Knepley                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
23465876a83SMatthew G. Knepley                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
23565876a83SMatthew G. Knepley {
23665876a83SMatthew G. Knepley   const PetscReal G      = PetscRealPart(constants[0]);
23765876a83SMatthew G. Knepley   const PetscReal K_u    = PetscRealPart(constants[1]);
23865876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
23965876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
24065876a83SMatthew G. Knepley   const PetscReal K_d    = K_u - alpha*alpha*M;
24165876a83SMatthew G. Knepley   const PetscReal lambda = K_d - (2.0 * G) / 3.0;
24265876a83SMatthew G. Knepley   PetscInt        d;
24365876a83SMatthew G. Knepley 
24465876a83SMatthew G. Knepley   for (d = 0; d < dim-1; ++d) {
24565876a83SMatthew G. Knepley     f0[d] -= 2.0*G - alpha*PetscCosReal(t);
24665876a83SMatthew G. Knepley   }
24765876a83SMatthew G. Knepley   f0[dim-1] -= 2.0*lambda + 4.0*G - alpha*PetscCosReal(t);
24865876a83SMatthew G. Knepley }
24965876a83SMatthew G. Knepley 
25065876a83SMatthew G. Knepley static void f0_quadratic_trig_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
25165876a83SMatthew G. Knepley                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
25265876a83SMatthew G. Knepley                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
25365876a83SMatthew G. Knepley                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
25465876a83SMatthew G. Knepley {
25565876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
25665876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
25765876a83SMatthew G. Knepley   PetscReal       sum    = 0.0;
25865876a83SMatthew G. Knepley   PetscInt        d;
25965876a83SMatthew G. Knepley 
26065876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += x[d];
26165876a83SMatthew G. Knepley 
26265876a83SMatthew G. Knepley   f0[0] += u_t ? alpha*u_t[uOff[1]] : 0.0;
26365876a83SMatthew G. Knepley   f0[0] += u_t ? u_t[uOff[2]]/M     : 0.0;
26465876a83SMatthew G. Knepley   f0[0] += PetscSinReal(t)*sum/M;
26565876a83SMatthew G. Knepley }
26665876a83SMatthew G. Knepley 
26765876a83SMatthew G. Knepley /* Trigonometric space and linear time solution
26865876a83SMatthew G. Knepley 
26965876a83SMatthew G. Knepley u = sin(2 pi x)
27065876a83SMatthew G. Knepley v = sin(2 pi y) - 2xy
27165876a83SMatthew G. Knepley \varepsilon = / 2 pi cos(2 pi x)             -y        \
27265876a83SMatthew G. Knepley               \      -y          2 pi cos(2 pi y) - 2x /
27365876a83SMatthew G. Knepley Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x
27465876a83SMatthew G. Knepley div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
27565876a83SMatthew G. Knepley   = \lambda \partial_j 2 pi (cos(2 pi x) + cos(2 pi y)) + 2\mu < -4 pi^2 sin(2 pi x) - 1, -4 pi^2 sin(2 pi y) >
27665876a83SMatthew G. Knepley   = \lambda < -4 pi^2 sin(2 pi x) - 2, -4 pi^2 sin(2 pi y) > + \mu < -8 pi^2 sin(2 pi x) - 2, -8 pi^2 sin(2 pi y) >
27765876a83SMatthew G. Knepley 
27865876a83SMatthew G. Knepley   2D:
27965876a83SMatthew G. Knepley   u = sin(2 pi x)
28065876a83SMatthew G. Knepley   v = sin(2 pi y) - 2xy
28165876a83SMatthew G. Knepley   p = (cos(2 pi x) + cos(2 pi y)) t
28265876a83SMatthew G. Knepley   e = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x
28365876a83SMatthew G. Knepley   f = < -4 pi^2 sin(2 pi x) (2 G + lambda) - (2 G - 2 lambda), -4 pi^2 sin(2 pi y) (2G + lambda) > + 2 pi alpha t <sin(2 pi x), sin(2 pi y)>
28465876a83SMatthew G. Knepley   g = 0
28565876a83SMatthew G. Knepley   \varepsilon = / 2 pi cos(2 pi x)             -y        \
28665876a83SMatthew G. Knepley                 \      -y          2 pi cos(2 pi y) - 2x /
28765876a83SMatthew G. Knepley   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x
28865876a83SMatthew G. Knepley   div \sigma = \partial_i 2 G \epsilon_{ij} + \partial_i \lambda \varepsilon \delta_{ij} - \partial_i \alpha p \delta_{ij}
28965876a83SMatthew G. Knepley     = 2 G < -4 pi^2 sin(2 pi x) - 1, -4 pi^2 sin(2 pi y) > + \lambda <-4 pi^2 sin(2 pi x) - 2, -4 pi^2 sin(2 pi y)> - alpha <-2 pi sin(2 pi x) t, -2 pi sin(2 pi y) t>
29065876a83SMatthew G. Knepley     = < -4 pi^2 sin(2 pi x) (2 G + lambda) - (2 G + 2 lambda), -4 pi^2 sin(2 pi y) (2G + lambda) > + 2 pi alpha t <sin(2 pi x), sin(2 pi y)>
29165876a83SMatthew G. Knepley   \frac{1}{M} \frac{dp}{dt} + \alpha \frac{d\varepsilon}{dt} - \nabla \cdot \kappa \nabla p
29265876a83SMatthew G. Knepley     = \frac{1}{M} \frac{dp}{dt} + \kappa \Delta p
29365876a83SMatthew G. Knepley     = (cos(2 pi x) + cos(2 pi y))/M - 4 pi^2 \kappa (cos(2 pi x) + cos(2 pi y)) t
29465876a83SMatthew G. Knepley 
29565876a83SMatthew G. Knepley   3D:
29665876a83SMatthew G. Knepley   u = sin(2 pi x)
29765876a83SMatthew G. Knepley   v = sin(2 pi y) - 2xy
29865876a83SMatthew G. Knepley   v = sin(2 pi y) - 2yz
29965876a83SMatthew G. Knepley   p = (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) t
30065876a83SMatthew G. Knepley   e = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2y
30165876a83SMatthew G. Knepley   f = < -4 pi^2 sin(2 pi x) (2 G + lambda) - (2 G + 2 lambda),  -4 pi^2 sin(2 pi y) (2 G + lambda) - (2 G + 2 lambda), -4 pi^2 sin(2 pi z) (2G + lambda) > + 2 pi alpha t <sin(2 pi x), sin(2 pi y), , sin(2 pi z)>
30265876a83SMatthew G. Knepley   g = 0
30365876a83SMatthew G. Knepley   \varepsilon = / 2 pi cos(2 pi x)            -y                     0         \
30465876a83SMatthew G. Knepley                 |         -y       2 pi cos(2 pi y) - 2x            -z         |
30565876a83SMatthew G. Knepley                 \          0                  -z         2 pi cos(2 pi z) - 2y /
30665876a83SMatthew G. Knepley   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y
30765876a83SMatthew G. Knepley   div \sigma = \partial_i 2 G \epsilon_{ij} + \partial_i \lambda \varepsilon \delta_{ij} - \partial_i \alpha p \delta_{ij}
30865876a83SMatthew G. Knepley     = 2 G < -4 pi^2 sin(2 pi x) - 1, -4 pi^2 sin(2 pi y) - 1, -4 pi^2 sin(2 pi z) > + \lambda <-4 pi^2 sin(2 pi x) - 2, 4 pi^2 sin(2 pi y) - 2, -4 pi^2 sin(2 pi y)> - alpha <-2 pi sin(2 pi x) t, -2 pi sin(2 pi y) t, -2 pi sin(2 pi z) t>
30965876a83SMatthew G. Knepley     = < -4 pi^2 sin(2 pi x) (2 G + lambda) - (2 G + 2 lambda),  -4 pi^2 sin(2 pi y) (2 G + lambda) - (2 G + 2 lambda), -4 pi^2 sin(2 pi z) (2G + lambda) > + 2 pi alpha t <sin(2 pi x), sin(2 pi y), , sin(2 pi z)>
31065876a83SMatthew G. Knepley   \frac{1}{M} \frac{dp}{dt} + \alpha \frac{d\varepsilon}{dt} - \nabla \cdot \kappa \nabla p
31165876a83SMatthew G. Knepley     = \frac{1}{M} \frac{dp}{dt} + \kappa \Delta p
31265876a83SMatthew G. Knepley     = (cos(2 pi x) + cos(2 pi y) + cos(2 pi z))/M - 4 pi^2 \kappa (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) t
31365876a83SMatthew G. Knepley */
31465876a83SMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
31565876a83SMatthew G. Knepley {
31665876a83SMatthew G. Knepley   PetscInt d;
31765876a83SMatthew G. Knepley 
31865876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
31965876a83SMatthew G. Knepley     u[d] = PetscSinReal(2.*PETSC_PI*x[d]) - (d > 0 ? 2.0 * x[d-1] * x[d] : 0.0);
32065876a83SMatthew G. Knepley   }
32165876a83SMatthew G. Knepley   return 0;
32265876a83SMatthew G. Knepley }
32365876a83SMatthew G. Knepley 
32465876a83SMatthew G. Knepley static PetscErrorCode trig_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
32565876a83SMatthew G. Knepley {
32665876a83SMatthew G. Knepley   PetscReal sum = 0.0;
32765876a83SMatthew G. Knepley   PetscInt  d;
32865876a83SMatthew G. Knepley 
32965876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += 2.*PETSC_PI*PetscCosReal(2.*PETSC_PI*x[d]) - (d < dim-1 ? 2.*x[d] : 0.0);
33065876a83SMatthew G. Knepley   u[0] = sum;
33165876a83SMatthew G. Knepley   return 0;
33265876a83SMatthew G. Knepley }
33365876a83SMatthew G. Knepley 
33465876a83SMatthew G. Knepley static PetscErrorCode trig_linear_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
33565876a83SMatthew G. Knepley {
33665876a83SMatthew G. Knepley   PetscReal sum = 0.0;
33765876a83SMatthew G. Knepley   PetscInt  d;
33865876a83SMatthew G. Knepley 
33965876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += PetscCosReal(2.*PETSC_PI*x[d]);
34065876a83SMatthew G. Knepley   u[0] = sum*time;
34165876a83SMatthew G. Knepley   return 0;
34265876a83SMatthew G. Knepley }
34365876a83SMatthew G. Knepley 
34465876a83SMatthew G. Knepley static PetscErrorCode trig_linear_p_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
34565876a83SMatthew G. Knepley {
34665876a83SMatthew G. Knepley   PetscReal sum = 0.0;
34765876a83SMatthew G. Knepley   PetscInt  d;
34865876a83SMatthew G. Knepley 
34965876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += PetscCosReal(2.*PETSC_PI*x[d]);
35065876a83SMatthew G. Knepley   u[0] = sum;
35165876a83SMatthew G. Knepley   return 0;
35265876a83SMatthew G. Knepley }
35365876a83SMatthew G. Knepley 
35465876a83SMatthew G. Knepley static void f0_trig_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
35565876a83SMatthew G. Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
35665876a83SMatthew G. Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
35765876a83SMatthew G. Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
35865876a83SMatthew G. Knepley {
35965876a83SMatthew G. Knepley   const PetscReal G      = PetscRealPart(constants[0]);
36065876a83SMatthew G. Knepley   const PetscReal K_u    = PetscRealPart(constants[1]);
36165876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
36265876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
36365876a83SMatthew G. Knepley   const PetscReal K_d    = K_u - alpha*alpha*M;
36465876a83SMatthew G. Knepley   const PetscReal lambda = K_d - (2.0 * G) / 3.0;
36565876a83SMatthew G. Knepley   PetscInt        d;
36665876a83SMatthew G. Knepley 
36765876a83SMatthew G. Knepley   for (d = 0; d < dim-1; ++d) {
36865876a83SMatthew G. Knepley     f0[d] += PetscSqr(2.*PETSC_PI)*PetscSinReal(2.*PETSC_PI*x[d])*(2.*G + lambda) + 2.0*(G + lambda) - 2.*PETSC_PI*alpha*PetscSinReal(2.*PETSC_PI*x[d])*t;
36965876a83SMatthew G. Knepley   }
37065876a83SMatthew G. Knepley   f0[dim-1] += PetscSqr(2.*PETSC_PI)*PetscSinReal(2.*PETSC_PI*x[dim-1])*(2.*G + lambda) - 2.*PETSC_PI*alpha*PetscSinReal(2.*PETSC_PI*x[dim-1])*t;
37165876a83SMatthew G. Knepley }
37265876a83SMatthew G. Knepley 
37365876a83SMatthew G. Knepley static void f0_trig_linear_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
37465876a83SMatthew G. Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
37565876a83SMatthew G. Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
37665876a83SMatthew G. Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
37765876a83SMatthew G. Knepley {
37865876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
37965876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
38065876a83SMatthew G. Knepley   const PetscReal kappa  = PetscRealPart(constants[4]);
38165876a83SMatthew G. Knepley   PetscReal       sum    = 0.0;
38265876a83SMatthew G. Knepley   PetscInt        d;
38365876a83SMatthew G. Knepley 
38465876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) sum += PetscCosReal(2.*PETSC_PI*x[d]);
38565876a83SMatthew G. Knepley   f0[0] += u_t ? alpha*u_t[uOff[1]] : 0.0;
38665876a83SMatthew G. Knepley   f0[0] += u_t ? u_t[uOff[2]]/M     : 0.0;
38765876a83SMatthew G. Knepley   f0[0] -= sum/M - 4*PetscSqr(PETSC_PI)*kappa*sum*t;
38865876a83SMatthew G. Knepley }
38965876a83SMatthew G. Knepley 
39065876a83SMatthew G. Knepley /* Terzaghi Solutions */
39165876a83SMatthew G. Knepley /* The analytical solutions given here are drawn from chapter 7, section 3, */
39265876a83SMatthew G. Knepley /* "One-Dimensional Consolidation Problem," from Poroelasticity, by Cheng.  */
39365876a83SMatthew G. Knepley static PetscErrorCode terzaghi_drainage_pressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
39465876a83SMatthew G. Knepley {
39565876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
39665876a83SMatthew G. Knepley   Parameter     *param;
39765876a83SMatthew G. Knepley 
3985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
39965876a83SMatthew G. Knepley   if (time <= 0.0) {
40065876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
40165876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
40265876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
40365876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
40465876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
40565876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
40665876a83SMatthew G. Knepley     PetscScalar eta   = (3.0*alpha*G) / (3.0*K_d + 4.0*G);         /* -,       Cheng (B.11) */
40765876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
40865876a83SMatthew G. Knepley 
40965876a83SMatthew G. Knepley     u[0] = ((P_0*eta) / (G*S));
41065876a83SMatthew G. Knepley   } else {
41165876a83SMatthew G. Knepley     u[0] = 0.0;
41265876a83SMatthew G. Knepley   }
41365876a83SMatthew G. Knepley   return 0;
41465876a83SMatthew G. Knepley }
41565876a83SMatthew G. Knepley 
41665876a83SMatthew G. Knepley static PetscErrorCode terzaghi_initial_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
41765876a83SMatthew G. Knepley {
41865876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
41965876a83SMatthew G. Knepley   Parameter     *param;
42065876a83SMatthew G. Knepley 
4215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
42265876a83SMatthew G. Knepley   {
42365876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
42465876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
42565876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
42630602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
42765876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G)); /* -,       Cheng (B.9)  */
42865876a83SMatthew G. Knepley     PetscReal   zstar = x[1] / L;                                /* - */
42965876a83SMatthew G. Knepley 
43065876a83SMatthew G. Knepley     u[0] = 0.0;
43165876a83SMatthew G. Knepley     u[1] = ((P_0*L*(1.0 - 2.0*nu_u)) / (2.0*G*(1.0 - nu_u))) * (1.0 - zstar);
43265876a83SMatthew G. Knepley   }
43365876a83SMatthew G. Knepley   return 0;
43465876a83SMatthew G. Knepley }
43565876a83SMatthew G. Knepley 
43665876a83SMatthew G. Knepley static PetscErrorCode terzaghi_initial_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
43765876a83SMatthew G. Knepley {
43865876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
43965876a83SMatthew G. Knepley   Parameter     *param;
44065876a83SMatthew G. Knepley 
4415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
44265876a83SMatthew G. Knepley   {
44365876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
44465876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
44565876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
44665876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
44765876a83SMatthew G. Knepley 
44865876a83SMatthew G. Knepley     u[0] = -(P_0*(1.0 - 2.0*nu_u)) / (2.0*G*(1.0 - nu_u));
44965876a83SMatthew G. Knepley   }
45065876a83SMatthew G. Knepley   return 0;
45165876a83SMatthew G. Knepley }
45265876a83SMatthew G. Knepley 
45365876a83SMatthew G. Knepley static PetscErrorCode terzaghi_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
45465876a83SMatthew G. Knepley {
45565876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
45665876a83SMatthew G. Knepley   Parameter     *param;
45765876a83SMatthew G. Knepley 
4585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
45965876a83SMatthew G. Knepley   if (time < 0.0) {
4605f80ce2aSJacob Faibussowitsch     CHKERRQ(terzaghi_initial_u(dim, time, x, Nc, u, ctx));
46165876a83SMatthew G. Knepley   } else {
46265876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
46365876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
46465876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
46565876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
46665876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
46765876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
46830602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
46965876a83SMatthew G. Knepley     PetscInt    N     = user->niter, m;
47065876a83SMatthew G. Knepley 
47165876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
47265876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
47365876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
47465876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
47565876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
47665876a83SMatthew G. Knepley 
47765876a83SMatthew G. Knepley     PetscReal   zstar = x[1] / L;                                  /* - */
47865876a83SMatthew G. Knepley     PetscReal   tstar = PetscRealPart(c*time) / PetscSqr(2.0*L);   /* - */
47965876a83SMatthew G. Knepley     PetscScalar F2    = 0.0;
48065876a83SMatthew G. Knepley 
48165876a83SMatthew G. Knepley     for (m = 1; m < 2*N+1; ++m) {
48265876a83SMatthew G. Knepley       if (m%2 == 1) {
48365876a83SMatthew G. Knepley         F2 += (8.0 / PetscSqr(m*PETSC_PI)) * PetscCosReal(0.5*m*PETSC_PI*zstar) * (1.0 - PetscExpReal(-PetscSqr(m*PETSC_PI)*tstar));
48465876a83SMatthew G. Knepley       }
48565876a83SMatthew G. Knepley     }
48665876a83SMatthew G. Knepley     u[0] = 0.0;
48765876a83SMatthew G. Knepley     u[1] = ((P_0*L*(1.0 - 2.0*nu_u)) / (2.0*G*(1.0 - nu_u))) * (1.0 - zstar) + ((P_0*L*(nu_u - nu)) / (2.0*G*(1.0 - nu_u)*(1.0 - nu)))*F2; /* m */
48865876a83SMatthew G. Knepley   }
48965876a83SMatthew G. Knepley   return 0;
49065876a83SMatthew G. Knepley }
49165876a83SMatthew G. Knepley 
49265876a83SMatthew G. Knepley static PetscErrorCode terzaghi_2d_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
49365876a83SMatthew G. Knepley {
49465876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
49565876a83SMatthew G. Knepley   Parameter     *param;
49665876a83SMatthew G. Knepley 
4975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
49865876a83SMatthew G. Knepley   if (time < 0.0) {
4995f80ce2aSJacob Faibussowitsch     CHKERRQ(terzaghi_initial_eps(dim, time, x, Nc, u, ctx));
50065876a83SMatthew G. Knepley   } else {
50165876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
50265876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
50365876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
50465876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
50565876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
50665876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
50730602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
50865876a83SMatthew G. Knepley     PetscInt    N     = user->niter, m;
50965876a83SMatthew G. Knepley 
51065876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
51165876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
51265876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
51365876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
51465876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
51565876a83SMatthew G. Knepley 
51665876a83SMatthew G. Knepley     PetscReal   zstar = x[1] / L;                                  /* - */
51765876a83SMatthew G. Knepley     PetscReal   tstar = PetscRealPart(c*time) / PetscSqr(2.0*L);   /* - */
51865876a83SMatthew G. Knepley     PetscScalar F2_z  = 0.0;
51965876a83SMatthew G. Knepley 
52065876a83SMatthew G. Knepley     for (m = 1; m < 2*N+1; ++m) {
52165876a83SMatthew G. Knepley       if (m%2 == 1) {
52265876a83SMatthew G. Knepley         F2_z += (-4.0 / (m*PETSC_PI*L)) * PetscSinReal(0.5*m*PETSC_PI*zstar) * (1.0 - PetscExpReal(-PetscSqr(m*PETSC_PI)*tstar));
52365876a83SMatthew G. Knepley       }
52465876a83SMatthew G. Knepley     }
52565876a83SMatthew G. Knepley     u[0] = -((P_0*L*(1.0 - 2.0*nu_u)) / (2.0*G*(1.0 - nu_u)*L)) + ((P_0*L*(nu_u - nu)) / (2.0*G*(1.0 - nu_u)*(1.0 - nu)))*F2_z; /* - */
52665876a83SMatthew G. Knepley   }
52765876a83SMatthew G. Knepley   return 0;
52865876a83SMatthew G. Knepley }
52965876a83SMatthew G. Knepley 
53065876a83SMatthew G. Knepley // Pressure
53165876a83SMatthew G. Knepley static PetscErrorCode terzaghi_2d_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
53265876a83SMatthew G. Knepley {
53365876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
53465876a83SMatthew G. Knepley   Parameter     *param;
53565876a83SMatthew G. Knepley 
5365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
53765876a83SMatthew G. Knepley   if (time <= 0.0) {
5385f80ce2aSJacob Faibussowitsch     CHKERRQ(terzaghi_drainage_pressure(dim, time, x, Nc, u, ctx));
53965876a83SMatthew G. Knepley   } else {
54065876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
54165876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
54265876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
54365876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
54465876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
54565876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
54630602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
54765876a83SMatthew G. Knepley     PetscInt    N     = user->niter, m;
54865876a83SMatthew G. Knepley 
54965876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
55065876a83SMatthew G. Knepley     PetscScalar eta   = (3.0*alpha*G) / (3.0*K_d + 4.0*G);         /* -,       Cheng (B.11) */
55165876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
55265876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
55365876a83SMatthew G. Knepley 
55465876a83SMatthew G. Knepley     PetscReal   zstar = x[1] / L;                                  /* - */
55565876a83SMatthew G. Knepley     PetscReal   tstar = PetscRealPart(c*time) / PetscSqr(2.0*L);   /* - */
55665876a83SMatthew G. Knepley     PetscScalar F1    = 0.0;
55765876a83SMatthew G. Knepley 
5583c633725SBarry Smith     PetscCheck(PetscAbsScalar((1/M + (alpha*eta)/G) - S) <= 1.0e-10,PETSC_COMM_SELF, PETSC_ERR_PLIB, "S %g != check %g", S, (1/M + (alpha*eta)/G));
55965876a83SMatthew G. Knepley 
56065876a83SMatthew G. Knepley     for (m = 1; m < 2*N+1; ++m) {
56165876a83SMatthew G. Knepley       if (m%2 == 1) {
56265876a83SMatthew G. Knepley         F1 += (4.0 / (m*PETSC_PI)) * PetscSinReal(0.5*m*PETSC_PI*zstar) * PetscExpReal(-PetscSqr(m*PETSC_PI)*tstar);
56365876a83SMatthew G. Knepley       }
56465876a83SMatthew G. Knepley     }
56565876a83SMatthew G. Knepley     u[0] = ((P_0*eta) / (G*S)) * F1; /* Pa */
56665876a83SMatthew G. Knepley   }
56765876a83SMatthew G. Knepley   return 0;
56865876a83SMatthew G. Knepley }
56965876a83SMatthew G. Knepley 
57065876a83SMatthew G. Knepley static PetscErrorCode terzaghi_2d_u_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
57165876a83SMatthew G. Knepley {
57265876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
57365876a83SMatthew G. Knepley   Parameter     *param;
57465876a83SMatthew G. Knepley 
5755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
57665876a83SMatthew G. Knepley   if (time <= 0.0) {
57765876a83SMatthew G. Knepley     u[0] = 0.0;
57865876a83SMatthew G. Knepley     u[1] = 0.0;
57965876a83SMatthew G. Knepley   } else {
58065876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
58165876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
58265876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
58365876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
58465876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
58565876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
58630602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
58765876a83SMatthew G. Knepley     PetscInt    N     = user->niter, m;
58865876a83SMatthew G. Knepley 
58965876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
59065876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
59165876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
59265876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
59365876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
59465876a83SMatthew G. Knepley 
59565876a83SMatthew G. Knepley     PetscReal   zstar = x[1] / L;                                  /* - */
59665876a83SMatthew G. Knepley     PetscReal   tstar = PetscRealPart(c*time) / PetscSqr(2.0*L);   /* - */
59765876a83SMatthew G. Knepley     PetscScalar F2_t  = 0.0;
59865876a83SMatthew G. Knepley 
59965876a83SMatthew G. Knepley     for (m = 1; m < 2*N+1; ++m) {
60065876a83SMatthew G. Knepley       if (m%2 == 1) {
60165876a83SMatthew G. Knepley         F2_t += (2.0*c / PetscSqr(L)) * PetscCosReal(0.5*m*PETSC_PI*zstar) * PetscExpReal(-PetscSqr(m*PETSC_PI)*tstar);
60265876a83SMatthew G. Knepley       }
60365876a83SMatthew G. Knepley     }
60465876a83SMatthew G. Knepley     u[0] = 0.0;
60565876a83SMatthew G. Knepley     u[1] = ((P_0*L*(nu_u - nu)) / (2.0*G*(1.0 - nu_u)*(1.0 - nu)))*F2_t; /* m / s */
60665876a83SMatthew G. Knepley   }
60765876a83SMatthew G. Knepley   return 0;
60865876a83SMatthew G. Knepley }
60965876a83SMatthew G. Knepley 
61065876a83SMatthew G. Knepley static PetscErrorCode terzaghi_2d_eps_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
61165876a83SMatthew G. Knepley {
61265876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
61365876a83SMatthew G. Knepley   Parameter     *param;
61465876a83SMatthew G. Knepley 
6155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
61665876a83SMatthew G. Knepley   if (time <= 0.0) {
61765876a83SMatthew G. Knepley     u[0] = 0.0;
61865876a83SMatthew G. Knepley   } else {
61965876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
62065876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
62165876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
62265876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
62365876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
62465876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
62530602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
62665876a83SMatthew G. Knepley     PetscInt    N     = user->niter, m;
62765876a83SMatthew G. Knepley 
62865876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
62965876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
63065876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
63165876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
63265876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
63365876a83SMatthew G. Knepley 
63465876a83SMatthew G. Knepley     PetscReal   zstar = x[1] / L;                                  /* - */
63565876a83SMatthew G. Knepley     PetscReal   tstar = PetscRealPart(c*time) / PetscSqr(2.0*L);   /* - */
63665876a83SMatthew G. Knepley     PetscScalar F2_zt = 0.0;
63765876a83SMatthew G. Knepley 
63865876a83SMatthew G. Knepley     for (m = 1; m < 2*N+1; ++m) {
63965876a83SMatthew G. Knepley       if (m%2 == 1) {
64065876a83SMatthew G. Knepley         F2_zt += ((-m*PETSC_PI*c) / (L*L*L)) * PetscSinReal(0.5*m*PETSC_PI*zstar) * PetscExpReal(-PetscSqr(m*PETSC_PI)*tstar);
64165876a83SMatthew G. Knepley       }
64265876a83SMatthew G. Knepley     }
64365876a83SMatthew G. Knepley     u[0] = ((P_0*L*(nu_u - nu)) / (2.0*G*(1.0 - nu_u)*(1.0 - nu)))*F2_zt; /* 1 / s */
64465876a83SMatthew G. Knepley   }
64565876a83SMatthew G. Knepley   return 0;
64665876a83SMatthew G. Knepley }
64765876a83SMatthew G. Knepley 
64865876a83SMatthew G. Knepley static PetscErrorCode terzaghi_2d_p_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
64965876a83SMatthew G. Knepley {
65065876a83SMatthew G. Knepley 
65165876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
65265876a83SMatthew G. Knepley   Parameter     *param;
65365876a83SMatthew G. Knepley 
6545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
65565876a83SMatthew G. Knepley   if (time <= 0.0) {
65665876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
65765876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
65865876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
65965876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
66065876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
66165876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
66230602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
66365876a83SMatthew G. Knepley 
66465876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
66565876a83SMatthew G. Knepley     PetscScalar eta   = (3.0*alpha*G) / (3.0*K_d + 4.0*G);         /* -,       Cheng (B.11) */
66665876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
66765876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
66865876a83SMatthew G. Knepley 
66965876a83SMatthew G. Knepley     u[0] = -((P_0*eta) / (G*S)) * PetscSqr(0*PETSC_PI)*c / PetscSqr(2.0*L); /* Pa / s */
67065876a83SMatthew G. Knepley   } else {
67165876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
67265876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
67365876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
67465876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
67565876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
67665876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
67730602db0SMatthew G. Knepley     PetscReal   L     = user->xmax[1] - user->xmin[1]; /* m */
67865876a83SMatthew G. Knepley     PetscInt    N     = user->niter, m;
67965876a83SMatthew G. Knepley 
68065876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
68165876a83SMatthew G. Knepley     PetscScalar eta   = (3.0*alpha*G) / (3.0*K_d + 4.0*G);         /* -,       Cheng (B.11) */
68265876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
68365876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
68465876a83SMatthew G. Knepley 
68565876a83SMatthew G. Knepley     PetscReal   zstar = x[1] / L;                                  /* - */
68665876a83SMatthew G. Knepley     PetscReal   tstar = PetscRealPart(c*time) / PetscSqr(2.0*L);   /* - */
68765876a83SMatthew G. Knepley     PetscScalar F1_t  = 0.0;
68865876a83SMatthew G. Knepley 
6893c633725SBarry Smith     PetscCheck(PetscAbsScalar((1/M + (alpha*eta)/G) - S) <= 1.0e-10,PETSC_COMM_SELF, PETSC_ERR_PLIB, "S %g != check %g", S, (1/M + (alpha*eta)/G));
69065876a83SMatthew G. Knepley 
69165876a83SMatthew G. Knepley     for (m = 1; m < 2*N+1; ++m) {
69265876a83SMatthew G. Knepley       if (m%2 == 1) {
69365876a83SMatthew G. Knepley         F1_t += ((-m*PETSC_PI*c) / PetscSqr(L)) * PetscSinReal(0.5*m*PETSC_PI*zstar) * PetscExpReal(-PetscSqr(m*PETSC_PI)*tstar);
69465876a83SMatthew G. Knepley       }
69565876a83SMatthew G. Knepley     }
69665876a83SMatthew G. Knepley     u[0] = ((P_0*eta) / (G*S)) * F1_t; /* Pa / s */
69765876a83SMatthew G. Knepley   }
69865876a83SMatthew G. Knepley   return 0;
69965876a83SMatthew G. Knepley }
70065876a83SMatthew G. Knepley 
70165876a83SMatthew G. Knepley /* Mandel Solutions */
70265876a83SMatthew G. Knepley static PetscErrorCode mandel_drainage_pressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
70365876a83SMatthew G. Knepley {
70465876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
70565876a83SMatthew G. Knepley   Parameter     *param;
70665876a83SMatthew G. Knepley 
7075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
70865876a83SMatthew G. Knepley   if (time <= 0.0) {
70965876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
71065876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
71165876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
71265876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
71365876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
71465876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
71530602db0SMatthew G. Knepley     PetscReal   a     = 0.5*(user->xmax[0] - user->xmin[0]); /* m */
71665876a83SMatthew G. Knepley     PetscInt    N     = user->niter, n;
71765876a83SMatthew G. Knepley 
71865876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
71965876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
72065876a83SMatthew G. Knepley     PetscScalar B     = alpha*M / K_u;                             /* -,       Cheng (B.12) */
72165876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
72265876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
72365876a83SMatthew G. Knepley 
72465876a83SMatthew G. Knepley     PetscScalar A1    = 3.0 / (B * (1.0 + nu_u));
72565876a83SMatthew G. Knepley     PetscReal   aa    = 0.0;
72665876a83SMatthew G. Knepley     PetscReal   p     = 0.0;
72765876a83SMatthew G. Knepley     PetscReal   time  = 0.0;
72865876a83SMatthew G. Knepley 
72965876a83SMatthew G. Knepley     for (n = 1; n < N+1; ++n) {
73065876a83SMatthew G. Knepley       aa = user->zeroArray[n-1];
73165876a83SMatthew G. Knepley       p += (PetscSinReal(aa) / (aa - PetscSinReal(aa)*PetscCosReal(aa))) * (PetscCosReal( (aa*x[0]) / a) - PetscCosReal(aa)) * PetscExpReal(-1.0*(aa*aa * PetscRealPart(c) * time)/(a*a));
73265876a83SMatthew G. Knepley     }
73365876a83SMatthew G. Knepley     u[0] = ((2.0 * P_0) / (a*A1)) * p;
73465876a83SMatthew G. Knepley   } else {
73565876a83SMatthew G. Knepley     u[0] = 0.0;
73665876a83SMatthew G. Knepley   }
73765876a83SMatthew G. Knepley   return 0;
73865876a83SMatthew G. Knepley }
73965876a83SMatthew G. Knepley 
74065876a83SMatthew G. Knepley static PetscErrorCode mandel_initial_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
74165876a83SMatthew G. Knepley {
74265876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
74365876a83SMatthew G. Knepley   Parameter     *param;
74465876a83SMatthew G. Knepley 
7455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
74665876a83SMatthew G. Knepley   {
74765876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
74865876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
74965876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
75065876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
75165876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
75265876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
75330602db0SMatthew G. Knepley     PetscScalar a     = 0.5*(user->xmax[0] - user->xmin[0]); /* m */
75465876a83SMatthew G. Knepley     PetscInt    N     = user->niter, n;
75565876a83SMatthew G. Knepley 
75665876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
75765876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
75865876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
75965876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
76065876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
76165876a83SMatthew G. Knepley 
76265876a83SMatthew G. Knepley     PetscScalar A_s   = 0.0;
76365876a83SMatthew G. Knepley     PetscScalar B_s   = 0.0;
76465876a83SMatthew G. Knepley     PetscScalar time  = 0.0;
76565876a83SMatthew G. Knepley     PetscScalar alpha_n = 0.0;
76665876a83SMatthew G. Knepley 
76765876a83SMatthew G. Knepley     for (n = 1; n < N+1; ++n) {
76865876a83SMatthew G. Knepley       alpha_n = user->zeroArray[n-1];
76965876a83SMatthew G. Knepley       A_s += ((PetscSinReal(alpha_n) * PetscCosReal(alpha_n)) / (alpha_n - PetscSinReal(alpha_n) * PetscCosReal(alpha_n))) * PetscExpReal(-1*(alpha_n*alpha_n*c*time)/(a*a));
77065876a83SMatthew G. Knepley       B_s += (PetscCosReal(alpha_n) / (alpha_n - PetscSinReal(alpha_n)*PetscCosReal(alpha_n))) * PetscSinReal( (alpha_n * x[0])/a) * PetscExpReal(-1*(alpha_n*alpha_n*c*time)/(a*a));
77165876a83SMatthew G. Knepley     }
77265876a83SMatthew G. Knepley     u[0] = ((P_0*nu)/(2.0*G*a) - (P_0*nu_u)/(G*a) * A_s)* x[0] + P_0/G * B_s;
77365876a83SMatthew G. Knepley     u[1] = (-1*(P_0*(1.0-nu))/(2*G*a) + (P_0*(1-nu_u))/(G*a) * A_s)*x[1];
77465876a83SMatthew G. Knepley   }
77565876a83SMatthew G. Knepley   return 0;
77665876a83SMatthew G. Knepley }
77765876a83SMatthew G. Knepley 
77865876a83SMatthew G. Knepley static PetscErrorCode mandel_initial_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
77965876a83SMatthew G. Knepley {
78065876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
78165876a83SMatthew G. Knepley   Parameter     *param;
78265876a83SMatthew G. Knepley 
7835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
78465876a83SMatthew G. Knepley   {
78565876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
78665876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
78765876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
78865876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
78965876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
79065876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
79130602db0SMatthew G. Knepley     PetscReal   a     = 0.5*(user->xmax[0] - user->xmin[0]); /* m */
79265876a83SMatthew G. Knepley     PetscInt    N     = user->niter, n;
79365876a83SMatthew G. Knepley 
79465876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
79565876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
79665876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
79765876a83SMatthew G. Knepley     PetscReal   c     = PetscRealPart(kappa / S);                  /* m^2 / s, Cheng (B.16) */
79865876a83SMatthew G. Knepley 
79965876a83SMatthew G. Knepley     PetscReal   aa    = 0.0;
80065876a83SMatthew G. Knepley     PetscReal   eps_A = 0.0;
80165876a83SMatthew G. Knepley     PetscReal   eps_B = 0.0;
80265876a83SMatthew G. Knepley     PetscReal   eps_C = 0.0;
80365876a83SMatthew G. Knepley     PetscReal   time  = 0.0;
80465876a83SMatthew G. Knepley 
80565876a83SMatthew G. Knepley     for (n = 1; n < N+1; ++n) {
80665876a83SMatthew G. Knepley       aa     = user->zeroArray[n-1];
80765876a83SMatthew G. Knepley       eps_A += (aa * PetscExpReal( (-1.0*aa*aa*c*time)/(a*a))*PetscCosReal(aa)*PetscCosReal( (aa*x[0])/a)) / (a * (aa - PetscSinReal(aa)*PetscCosReal(aa)));
80865876a83SMatthew G. Knepley       eps_B += ( PetscExpReal( (-1.0*aa*aa*c*time)/(a*a))*PetscSinReal(aa)*PetscCosReal(aa)) / (aa - PetscSinReal(aa)*PetscCosReal(aa));
80965876a83SMatthew G. Knepley       eps_C += ( PetscExpReal( (-1.0*aa*aa*c*time)/(aa*aa))*PetscSinReal(aa)*PetscCosReal(aa)) / (aa - PetscSinReal(aa)*PetscCosReal(aa));
81065876a83SMatthew G. Knepley     }
81165876a83SMatthew G. Knepley     u[0] = (P_0/G)*eps_A + ( (P_0*nu)/(2.0*G*a)) - eps_B/(G*a) - (P_0*(1-nu))/(2*G*a) + eps_C/(G*a);
81265876a83SMatthew G. Knepley   }
81365876a83SMatthew G. Knepley   return 0;
81465876a83SMatthew G. Knepley }
81565876a83SMatthew G. Knepley 
81665876a83SMatthew G. Knepley // Displacement
81765876a83SMatthew G. Knepley static PetscErrorCode mandel_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
81865876a83SMatthew G. Knepley {
81965876a83SMatthew G. Knepley 
82065876a83SMatthew G. Knepley   Parameter  *param;
82165876a83SMatthew G. Knepley 
82265876a83SMatthew G. Knepley   AppCtx *user = (AppCtx *) ctx;
82365876a83SMatthew G. Knepley 
8245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
82565876a83SMatthew G. Knepley   if (time <= 0.0) {
8265f80ce2aSJacob Faibussowitsch     CHKERRQ(mandel_initial_u(dim, time, x, Nc, u, ctx));
82765876a83SMatthew G. Knepley   } else {
82865876a83SMatthew G. Knepley     PetscInt NITER = user->niter;
82965876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha;
83065876a83SMatthew G. Knepley     PetscScalar K_u = param->K_u;
83165876a83SMatthew G. Knepley     PetscScalar M = param->M;
83265876a83SMatthew G. Knepley     PetscScalar G = param->mu;
83365876a83SMatthew G. Knepley     PetscScalar k = param->k;
83465876a83SMatthew G. Knepley     PetscScalar mu_f = param->mu_f;
83565876a83SMatthew G. Knepley     PetscScalar F = param->P_0;
83665876a83SMatthew G. Knepley 
83765876a83SMatthew G. Knepley     PetscScalar K_d = K_u - alpha*alpha*M;
83865876a83SMatthew G. Knepley     PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
83965876a83SMatthew G. Knepley     PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
84065876a83SMatthew G. Knepley     PetscScalar kappa = k / mu_f;
84130602db0SMatthew G. Knepley     PetscReal   a = (user->xmax[0] - user->xmin[0]) / 2.0;
84265876a83SMatthew G. Knepley     PetscReal   c = PetscRealPart(((2.0*kappa*G) * (1.0 - nu) * (nu_u - nu)) / ( alpha*alpha * (1.0 - 2.0*nu) * (1.0 - nu_u)));
84365876a83SMatthew G. Knepley 
84465876a83SMatthew G. Knepley     // Series term
84565876a83SMatthew G. Knepley     PetscScalar A_x = 0.0;
84665876a83SMatthew G. Knepley     PetscScalar B_x = 0.0;
84765876a83SMatthew G. Knepley 
84865876a83SMatthew G. Knepley     for (PetscInt n=1; n < NITER+1; n++) {
84965876a83SMatthew G. Knepley       PetscReal alpha_n = user->zeroArray[n-1];
85065876a83SMatthew G. Knepley 
85165876a83SMatthew G. Knepley       A_x += ( (PetscSinReal(alpha_n) * PetscCosReal(alpha_n)) / (alpha_n - PetscSinReal(alpha_n) * PetscCosReal(alpha_n))) * PetscExpReal( -1*(alpha_n*alpha_n*c*time)/(a*a));
85265876a83SMatthew G. Knepley       B_x += ( PetscCosReal(alpha_n) / (alpha_n - PetscSinReal(alpha_n)*PetscCosReal(alpha_n))) * PetscSinReal( (alpha_n * x[0])/a) * PetscExpReal( -1*(alpha_n*alpha_n*c*time)/(a*a));
85365876a83SMatthew G. Knepley     }
85465876a83SMatthew G. Knepley     u[0] = ((F*nu)/(2.0*G*a) - (F*nu_u)/(G*a) * A_x)* x[0] + F/G * B_x;
85565876a83SMatthew G. Knepley     u[1] = (-1*(F*(1.0-nu))/(2*G*a) + (F*(1-nu_u))/(G*a) * A_x)*x[1];
85665876a83SMatthew G. Knepley   }
85765876a83SMatthew G. Knepley   return 0;
85865876a83SMatthew G. Knepley }
85965876a83SMatthew G. Knepley 
86065876a83SMatthew G. Knepley // Trace strain
86165876a83SMatthew G. Knepley static PetscErrorCode mandel_2d_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
86265876a83SMatthew G. Knepley {
86365876a83SMatthew G. Knepley 
86465876a83SMatthew G. Knepley   Parameter  *param;
86565876a83SMatthew G. Knepley 
86665876a83SMatthew G. Knepley   AppCtx *user = (AppCtx *) ctx;
86765876a83SMatthew G. Knepley 
8685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
86965876a83SMatthew G. Knepley   if (time <= 0.0) {
8705f80ce2aSJacob Faibussowitsch     CHKERRQ(mandel_initial_eps(dim, time, x, Nc, u, ctx));
87165876a83SMatthew G. Knepley   } else {
87265876a83SMatthew G. Knepley     PetscInt NITER = user->niter;
87365876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha;
87465876a83SMatthew G. Knepley     PetscScalar K_u = param->K_u;
87565876a83SMatthew G. Knepley     PetscScalar M = param->M;
87665876a83SMatthew G. Knepley     PetscScalar G = param->mu;
87765876a83SMatthew G. Knepley     PetscScalar k = param->k;
87865876a83SMatthew G. Knepley     PetscScalar mu_f = param->mu_f;
87965876a83SMatthew G. Knepley     PetscScalar F = param->P_0;
88065876a83SMatthew G. Knepley 
88165876a83SMatthew G. Knepley     PetscScalar K_d = K_u - alpha*alpha*M;
88265876a83SMatthew G. Knepley     PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
88365876a83SMatthew G. Knepley     PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
88465876a83SMatthew G. Knepley     PetscScalar kappa = k / mu_f;
88565876a83SMatthew G. Knepley     //const PetscScalar B = (alpha*M)/(K_d + alpha*alpha * M);
88665876a83SMatthew G. Knepley 
88765876a83SMatthew G. Knepley     //const PetscScalar b = (YMAX - YMIN) / 2.0;
88830602db0SMatthew G. Knepley     PetscScalar a = (user->xmax[0] - user->xmin[0]) / 2.0;
88965876a83SMatthew G. Knepley     PetscReal   c = PetscRealPart(((2.0*kappa*G) * (1.0 - nu) * (nu_u - nu)) / (alpha*alpha * (1.0 - 2.0*nu) * (1.0 - nu_u)));
89065876a83SMatthew G. Knepley 
89165876a83SMatthew G. Knepley     // Series term
89265876a83SMatthew G. Knepley     PetscScalar eps_A = 0.0;
89365876a83SMatthew G. Knepley     PetscScalar eps_B = 0.0;
89465876a83SMatthew G. Knepley     PetscScalar eps_C = 0.0;
89565876a83SMatthew G. Knepley 
89665876a83SMatthew G. Knepley     for (PetscInt n=1; n < NITER+1; n++)
89765876a83SMatthew G. Knepley     {
89865876a83SMatthew G. Knepley       PetscReal aa = user->zeroArray[n-1];
89965876a83SMatthew G. Knepley 
90065876a83SMatthew G. Knepley       eps_A += (aa * PetscExpReal( (-1.0*aa*aa*c*time)/(a*a))*PetscCosReal(aa)*PetscCosReal( (aa*x[0])/a)) / (a * (aa - PetscSinReal(aa)*PetscCosReal(aa)));
90165876a83SMatthew G. Knepley 
90265876a83SMatthew G. Knepley       eps_B += ( PetscExpReal( (-1.0*aa*aa*c*time)/(a*a))*PetscSinReal(aa)*PetscCosReal(aa)) / (aa - PetscSinReal(aa)*PetscCosReal(aa));
90365876a83SMatthew G. Knepley 
90465876a83SMatthew G. Knepley       eps_C += ( PetscExpReal( (-1.0*aa*aa*c*time)/(aa*aa))*PetscSinReal(aa)*PetscCosReal(aa)) / (aa - PetscSinReal(aa)*PetscCosReal(aa));
90565876a83SMatthew G. Knepley     }
90665876a83SMatthew G. Knepley 
90765876a83SMatthew G. Knepley     u[0] = (F/G)*eps_A + ( (F*nu)/(2.0*G*a)) - eps_B/(G*a) - (F*(1-nu))/(2*G*a) + eps_C/(G*a);
90865876a83SMatthew G. Knepley   }
90965876a83SMatthew G. Knepley   return 0;
91065876a83SMatthew G. Knepley 
91165876a83SMatthew G. Knepley }
91265876a83SMatthew G. Knepley 
91365876a83SMatthew G. Knepley // Pressure
91465876a83SMatthew G. Knepley static PetscErrorCode mandel_2d_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
91565876a83SMatthew G. Knepley {
91665876a83SMatthew G. Knepley 
91765876a83SMatthew G. Knepley   Parameter  *param;
91865876a83SMatthew G. Knepley 
91965876a83SMatthew G. Knepley   AppCtx *user = (AppCtx *) ctx;
92065876a83SMatthew G. Knepley 
9215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
92265876a83SMatthew G. Knepley   if (time <= 0.0) {
9235f80ce2aSJacob Faibussowitsch     CHKERRQ(mandel_drainage_pressure(dim, time, x, Nc, u, ctx));
92465876a83SMatthew G. Knepley   } else {
92565876a83SMatthew G. Knepley     PetscInt NITER = user->niter;
92665876a83SMatthew G. Knepley 
92765876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha;
92865876a83SMatthew G. Knepley     PetscScalar K_u = param->K_u;
92965876a83SMatthew G. Knepley     PetscScalar M = param->M;
93065876a83SMatthew G. Knepley     PetscScalar G = param->mu;
93165876a83SMatthew G. Knepley     PetscScalar k = param->k;
93265876a83SMatthew G. Knepley     PetscScalar mu_f = param->mu_f;
93365876a83SMatthew G. Knepley     PetscScalar F = param->P_0;
93465876a83SMatthew G. Knepley 
93565876a83SMatthew G. Knepley     PetscScalar K_d = K_u - alpha*alpha*M;
93665876a83SMatthew G. Knepley     PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
93765876a83SMatthew G. Knepley     PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
93865876a83SMatthew G. Knepley     PetscScalar kappa = k / mu_f;
93965876a83SMatthew G. Knepley     PetscScalar B = (alpha*M)/(K_d + alpha*alpha * M);
94065876a83SMatthew G. Knepley 
94130602db0SMatthew G. Knepley     PetscReal   a  = (user->xmax[0] - user->xmin[0]) / 2.0;
94265876a83SMatthew G. Knepley     PetscReal   c  = PetscRealPart(((2.0*kappa*G) * (1.0 - nu) * (nu_u - nu)) / (alpha*alpha * (1.0 - 2.0*nu) * (1.0 - nu_u)));
94365876a83SMatthew G. Knepley     PetscScalar A1 = 3.0 / (B * (1.0 + nu_u));
94465876a83SMatthew G. Knepley     //PetscScalar A2 = (alpha * (1.0 - 2.0*nu)) / (1.0 - nu);
94565876a83SMatthew G. Knepley 
94665876a83SMatthew G. Knepley     // Series term
94765876a83SMatthew G. Knepley     PetscScalar aa = 0.0;
94865876a83SMatthew G. Knepley     PetscScalar p  = 0.0;
94965876a83SMatthew G. Knepley 
95065876a83SMatthew G. Knepley     for (PetscInt n=1; n < NITER+1; n++)
95165876a83SMatthew G. Knepley     {
95265876a83SMatthew G. Knepley       aa = user->zeroArray[n-1];
95365876a83SMatthew G. Knepley       p += (PetscSinReal(aa)/ (aa - PetscSinReal(aa)*PetscCosReal(aa))) * (PetscCosReal( (aa*x[0]) / a) - PetscCosReal(aa)) * PetscExpReal(-1.0*(aa*aa * c * time)/(a*a));
95465876a83SMatthew G. Knepley     }
95565876a83SMatthew G. Knepley     u[0] = ((2.0 * F) / (a*A1)) * p;
95665876a83SMatthew G. Knepley   }
95765876a83SMatthew G. Knepley   return 0;
95865876a83SMatthew G. Knepley }
95965876a83SMatthew G. Knepley 
96065876a83SMatthew G. Knepley // Time derivative of displacement
96165876a83SMatthew G. Knepley static PetscErrorCode mandel_2d_u_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
96265876a83SMatthew G. Knepley {
96365876a83SMatthew G. Knepley 
96465876a83SMatthew G. Knepley   Parameter  *param;
96565876a83SMatthew G. Knepley 
96665876a83SMatthew G. Knepley   AppCtx *user = (AppCtx *) ctx;
96765876a83SMatthew G. Knepley 
9685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
96965876a83SMatthew G. Knepley 
97065876a83SMatthew G. Knepley   PetscInt NITER = user->niter;
97165876a83SMatthew G. Knepley   PetscScalar alpha = param->alpha;
97265876a83SMatthew G. Knepley   PetscScalar K_u = param->K_u;
97365876a83SMatthew G. Knepley   PetscScalar M = param->M;
97465876a83SMatthew G. Knepley   PetscScalar G = param->mu;
97565876a83SMatthew G. Knepley   PetscScalar F = param->P_0;
97665876a83SMatthew G. Knepley 
97765876a83SMatthew G. Knepley   PetscScalar K_d = K_u - alpha*alpha*M;
97865876a83SMatthew G. Knepley   PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
97965876a83SMatthew G. Knepley   PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
98065876a83SMatthew G. Knepley   PetscScalar kappa = param->k / param->mu_f;
98130602db0SMatthew G. Knepley   PetscReal   a = (user->xmax[0] - user->xmin[0]) / 2.0;
98265876a83SMatthew G. Knepley   PetscReal   c = PetscRealPart(((2.0*kappa*G) * (1.0 - nu) * (nu_u - nu)) / (alpha*alpha * (1.0 - 2.0*nu) * (1.0 - nu_u)));
98365876a83SMatthew G. Knepley 
98465876a83SMatthew G. Knepley   // Series term
98565876a83SMatthew G. Knepley   PetscScalar A_s_t = 0.0;
98665876a83SMatthew G. Knepley   PetscScalar B_s_t = 0.0;
98765876a83SMatthew G. Knepley 
98865876a83SMatthew G. Knepley   for (PetscInt n=1; n < NITER+1; n++)
98965876a83SMatthew G. Knepley   {
99065876a83SMatthew G. Knepley     PetscReal alpha_n = user->zeroArray[n-1];
99165876a83SMatthew G. Knepley 
99265876a83SMatthew G. Knepley     A_s_t += (-1.0*alpha_n*alpha_n*c*PetscExpReal( (-1.0*alpha_n*alpha_n*time)/(a*a))*PetscSinReal( (alpha_n*x[0])/a) * PetscCosReal(alpha_n)) / ( a*a*(alpha_n - PetscSinReal(alpha_n)*PetscCosReal(alpha_n)));
99365876a83SMatthew G. Knepley     B_s_t += (-1.0*alpha_n*alpha_n*c*PetscExpReal( (-1.0*alpha_n*alpha_n*time)/(a*a))*PetscSinReal(  alpha_n) * PetscCosReal(alpha_n)) / ( a*a*(alpha_n - PetscSinReal(alpha_n)*PetscCosReal(alpha_n)));
99465876a83SMatthew G. Knepley   }
99565876a83SMatthew G. Knepley 
99665876a83SMatthew G. Knepley   u[0] = (F/G)*A_s_t - ( (F*nu_u*x[0])/(G*a))*B_s_t;
99765876a83SMatthew G. Knepley   u[1] = ( (F*x[1]*(1 - nu_u)) / (G*a))*B_s_t;
99865876a83SMatthew G. Knepley 
99965876a83SMatthew G. Knepley   return 0;
100065876a83SMatthew G. Knepley }
100165876a83SMatthew G. Knepley 
100265876a83SMatthew G. Knepley // Time derivative of trace strain
100365876a83SMatthew G. Knepley static PetscErrorCode mandel_2d_eps_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
100465876a83SMatthew G. Knepley {
100565876a83SMatthew G. Knepley 
100665876a83SMatthew G. Knepley   Parameter  *param;
100765876a83SMatthew G. Knepley 
100865876a83SMatthew G. Knepley   AppCtx *user = (AppCtx *) ctx;
100965876a83SMatthew G. Knepley 
10105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
101165876a83SMatthew G. Knepley 
101265876a83SMatthew G. Knepley   PetscInt NITER = user->niter;
101365876a83SMatthew G. Knepley   PetscScalar alpha = param->alpha;
101465876a83SMatthew G. Knepley   PetscScalar K_u = param->K_u;
101565876a83SMatthew G. Knepley   PetscScalar M = param->M;
101665876a83SMatthew G. Knepley   PetscScalar G = param->mu;
101765876a83SMatthew G. Knepley   PetscScalar k = param->k;
101865876a83SMatthew G. Knepley   PetscScalar mu_f = param->mu_f;
101965876a83SMatthew G. Knepley   PetscScalar F = param->P_0;
102065876a83SMatthew G. Knepley 
102165876a83SMatthew G. Knepley   PetscScalar K_d = K_u - alpha*alpha*M;
102265876a83SMatthew G. Knepley   PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
102365876a83SMatthew G. Knepley   PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
102465876a83SMatthew G. Knepley   PetscScalar kappa = k / mu_f;
102565876a83SMatthew G. Knepley   //const PetscScalar B = (alpha*M)/(K_d + alpha*alpha * M);
102665876a83SMatthew G. Knepley 
102765876a83SMatthew G. Knepley   //const PetscScalar b = (YMAX - YMIN) / 2.0;
102830602db0SMatthew G. Knepley   PetscReal   a = (user->xmax[0] - user->xmin[0]) / 2.0;
102965876a83SMatthew G. Knepley   PetscReal   c = PetscRealPart(((2.0*kappa*G) * (1.0 - nu) * (nu_u - nu)) / (alpha*alpha * (1.0 - 2.0*nu) * (1.0 - nu_u)));
103065876a83SMatthew G. Knepley 
103165876a83SMatthew G. Knepley   // Series term
103265876a83SMatthew G. Knepley   PetscScalar eps_As = 0.0;
103365876a83SMatthew G. Knepley   PetscScalar eps_Bs = 0.0;
103465876a83SMatthew G. Knepley   PetscScalar eps_Cs = 0.0;
103565876a83SMatthew G. Knepley 
103665876a83SMatthew G. Knepley   for (PetscInt n=1; n < NITER+1; n++)
103765876a83SMatthew G. Knepley   {
103865876a83SMatthew G. Knepley     PetscReal alpha_n = user->zeroArray[n-1];
103965876a83SMatthew G. Knepley 
104065876a83SMatthew G. Knepley     eps_As += (-1.0*alpha_n*alpha_n*alpha_n*c*PetscExpReal( (-1.0*alpha_n*alpha_n*c*time)/(a*a))*PetscCosReal(alpha_n)*PetscCosReal( (alpha_n*x[0])/a)) / ( alpha_n*alpha_n*alpha_n*(alpha_n - PetscSinReal(alpha_n)*PetscCosReal(alpha_n)));
104165876a83SMatthew G. Knepley     eps_Bs += (-1.0*alpha_n*alpha_n*c*PetscExpReal( (-1.0*alpha_n*alpha_n*c*time)/(a*a))*PetscSinReal(alpha_n)*PetscCosReal(alpha_n)) / (alpha_n*alpha_n * (alpha_n - PetscSinReal(alpha_n)*PetscCosReal(alpha_n)));
104265876a83SMatthew G. Knepley     eps_Cs += (-1.0*alpha_n*alpha_n*c*PetscExpReal( (-1.0*alpha_n*alpha_n*c*time)/(a*a))*PetscSinReal(alpha_n)*PetscCosReal(alpha_n)) / (alpha_n*alpha_n * (alpha_n - PetscSinReal(alpha_n)*PetscCosReal(alpha_n)));
104365876a83SMatthew G. Knepley   }
104465876a83SMatthew G. Knepley 
104565876a83SMatthew G. Knepley   u[0] = (F/G)*eps_As - ( (F*nu_u)/(G*a))*eps_Bs + ( (F*(1-nu_u))/(G*a))*eps_Cs;
104665876a83SMatthew G. Knepley   return 0;
104765876a83SMatthew G. Knepley 
104865876a83SMatthew G. Knepley }
104965876a83SMatthew G. Knepley 
105065876a83SMatthew G. Knepley // Time derivative of pressure
105165876a83SMatthew G. Knepley static PetscErrorCode mandel_2d_p_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
105265876a83SMatthew G. Knepley {
105365876a83SMatthew G. Knepley 
105465876a83SMatthew G. Knepley   Parameter  *param;
105565876a83SMatthew G. Knepley 
105665876a83SMatthew G. Knepley   AppCtx *user = (AppCtx *) ctx;
105765876a83SMatthew G. Knepley 
10585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
105965876a83SMatthew G. Knepley 
106065876a83SMatthew G. Knepley   PetscScalar alpha = param->alpha;
106165876a83SMatthew G. Knepley   PetscScalar K_u = param->K_u;
106265876a83SMatthew G. Knepley   PetscScalar M = param->M;
106365876a83SMatthew G. Knepley   PetscScalar G = param->mu;
106465876a83SMatthew G. Knepley   PetscScalar F = param->P_0;
106565876a83SMatthew G. Knepley 
106665876a83SMatthew G. Knepley   PetscScalar K_d = K_u - alpha*alpha*M;
106765876a83SMatthew G. Knepley   PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
106865876a83SMatthew G. Knepley   PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
106965876a83SMatthew G. Knepley 
107030602db0SMatthew G. Knepley   PetscReal   a = (user->xmax[0] - user->xmin[0]) / 2.0;
107165876a83SMatthew G. Knepley   //PetscScalar A1 = 3.0 / (B * (1.0 + nu_u));
107265876a83SMatthew G. Knepley   //PetscScalar A2 = (alpha * (1.0 - 2.0*nu)) / (1.0 - nu);
107365876a83SMatthew G. Knepley 
107465876a83SMatthew G. Knepley   u[0] = ( (2.0*F*(-2.0*nu + 3.0*nu_u))/(3.0*a*alpha*(1.0 - 2.0*nu)));
107565876a83SMatthew G. Knepley 
107665876a83SMatthew G. Knepley   return 0;
107765876a83SMatthew G. Knepley }
107865876a83SMatthew G. Knepley 
107965876a83SMatthew G. Knepley /* Cryer Solutions */
108065876a83SMatthew G. Knepley static PetscErrorCode cryer_drainage_pressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
108165876a83SMatthew G. Knepley {
108265876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
108365876a83SMatthew G. Knepley   Parameter     *param;
108465876a83SMatthew G. Knepley 
10855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
108665876a83SMatthew G. Knepley   if (time <= 0.0) {
108765876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
108865876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
108965876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
109065876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
109165876a83SMatthew G. Knepley     PetscScalar B     = alpha*M / K_u; /* -, Cheng (B.12) */
109265876a83SMatthew G. Knepley 
109365876a83SMatthew G. Knepley     u[0] = P_0*B;
109465876a83SMatthew G. Knepley   } else {
109565876a83SMatthew G. Knepley     u[0] = 0.0;
109665876a83SMatthew G. Knepley   }
109765876a83SMatthew G. Knepley   return 0;
109865876a83SMatthew G. Knepley }
109965876a83SMatthew G. Knepley 
110065876a83SMatthew G. Knepley static PetscErrorCode cryer_initial_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
110165876a83SMatthew G. Knepley {
110265876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
110365876a83SMatthew G. Knepley   Parameter     *param;
110465876a83SMatthew G. Knepley 
11055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
110665876a83SMatthew G. Knepley   {
110765876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
110865876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
110965876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
111030602db0SMatthew G. Knepley     PetscReal   R_0   = user->xmax[1];  /* m */
111165876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
111265876a83SMatthew G. Knepley 
111365876a83SMatthew G. Knepley     PetscScalar u_0   = -P_0*R_0*(1. - 2.*nu_u) / (2.*G*(1. + nu_u)); /* Cheng (7.407) */
111465876a83SMatthew G. Knepley     PetscReal   u_sc  = PetscRealPart(u_0)/R_0;
111565876a83SMatthew G. Knepley 
111665876a83SMatthew G. Knepley     u[0] = u_sc * x[0];
111765876a83SMatthew G. Knepley     u[1] = u_sc * x[1];
111865876a83SMatthew G. Knepley     u[2] = u_sc * x[2];
111965876a83SMatthew G. Knepley   }
112065876a83SMatthew G. Knepley   return 0;
112165876a83SMatthew G. Knepley }
112265876a83SMatthew G. Knepley 
112365876a83SMatthew G. Knepley static PetscErrorCode cryer_initial_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
112465876a83SMatthew G. Knepley {
112565876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
112665876a83SMatthew G. Knepley   Parameter     *param;
112765876a83SMatthew G. Knepley 
11285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
112965876a83SMatthew G. Knepley   {
113065876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
113165876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
113265876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
113330602db0SMatthew G. Knepley     PetscReal   R_0   = user->xmax[1];  /* m */
113465876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
113565876a83SMatthew G. Knepley 
113665876a83SMatthew G. Knepley     PetscScalar u_0   = -P_0*R_0*(1. - 2.*nu_u) / (2.*G*(1. + nu_u)); /* Cheng (7.407) */
113765876a83SMatthew G. Knepley     PetscReal   u_sc  = PetscRealPart(u_0)/R_0;
113865876a83SMatthew G. Knepley 
113965876a83SMatthew G. Knepley     /* div R = 1/R^2 d/dR R^2 R = 3 */
114065876a83SMatthew G. Knepley     u[0] = 3.*u_sc;
114165876a83SMatthew G. Knepley     u[1] = 3.*u_sc;
114265876a83SMatthew G. Knepley     u[2] = 3.*u_sc;
114365876a83SMatthew G. Knepley   }
114465876a83SMatthew G. Knepley   return 0;
114565876a83SMatthew G. Knepley }
114665876a83SMatthew G. Knepley 
114765876a83SMatthew G. Knepley // Displacement
114865876a83SMatthew G. Knepley static PetscErrorCode cryer_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
114965876a83SMatthew G. Knepley {
115065876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
115165876a83SMatthew G. Knepley   Parameter     *param;
115265876a83SMatthew G. Knepley 
11535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
115465876a83SMatthew G. Knepley   if (time <= 0.0) {
11555f80ce2aSJacob Faibussowitsch     CHKERRQ(cryer_initial_u(dim, time, x, Nc, u, ctx));
115665876a83SMatthew G. Knepley   } else {
115765876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
115865876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
115965876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
116065876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
116165876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
116265876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
116330602db0SMatthew G. Knepley     PetscReal   R_0   = user->xmax[1];  /* m */
116465876a83SMatthew G. Knepley     PetscInt    N     = user->niter, n;
116565876a83SMatthew G. Knepley 
116665876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
116765876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
116865876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
116965876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
117065876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
117165876a83SMatthew G. Knepley     PetscScalar u_inf = -P_0*R_0*(1. - 2.*nu) / (2.*G*(1. + nu));  /* m,       Cheng (7.388) */
117265876a83SMatthew G. Knepley 
117365876a83SMatthew G. Knepley     PetscReal   R      = PetscSqrtReal(x[0]*x[0] + x[1]*x[1] + x[2]*x[2]);
117465876a83SMatthew G. Knepley     PetscReal   R_star = R/R_0;
117565876a83SMatthew G. Knepley     PetscReal   tstar  = PetscRealPart(c*time) / PetscSqr(R_0);    /* - */
117665876a83SMatthew G. Knepley     PetscReal   A_n    = 0.0;
117765876a83SMatthew G. Knepley     PetscScalar u_sc;
117865876a83SMatthew G. Knepley 
117965876a83SMatthew G. Knepley     for (n = 1; n < N+1; ++n) {
118065876a83SMatthew G. Knepley       const PetscReal x_n = user->zeroArray[n-1];
118165876a83SMatthew G. Knepley       const PetscReal E_n = PetscRealPart(PetscSqr(1 - nu)*PetscSqr(1 + nu_u)*x_n - 18.0*(1 + nu)*(nu_u - nu)*(1 - nu_u));
118265876a83SMatthew G. Knepley 
118365876a83SMatthew G. Knepley       /* m , Cheng (7.404) */
118465876a83SMatthew G. Knepley       A_n += PetscRealPart(
118565876a83SMatthew G. Knepley              (12.0*(1.0 + nu)*(nu_u - nu))/((1.0 - 2.0*nu)*E_n*PetscSqr(R_star)*x_n*PetscSinReal(PetscSqrtReal(x_n))) *
118665876a83SMatthew G. Knepley              (3.0*(nu_u - nu) * (PetscSinReal(R_star * PetscSqrtReal(x_n)) - R_star*PetscSqrtReal(x_n)*PetscCosReal(R_star * PetscSqrtReal(x_n)))
118765876a83SMatthew G. Knepley               + (1.0 - nu)*(1.0 - 2.0*nu)*PetscPowRealInt(R_star, 3)*x_n*PetscSinReal(PetscSqrtReal(x_n))) * PetscExpReal(-x_n * tstar));
118865876a83SMatthew G. Knepley     }
118965876a83SMatthew G. Knepley     u_sc = PetscRealPart(u_inf) * (R_star - A_n);
119065876a83SMatthew G. Knepley     u[0] = u_sc * x[0] / R;
119165876a83SMatthew G. Knepley     u[1] = u_sc * x[1] / R;
119265876a83SMatthew G. Knepley     u[2] = u_sc * x[2] / R;
119365876a83SMatthew G. Knepley   }
119465876a83SMatthew G. Knepley   return 0;
119565876a83SMatthew G. Knepley }
119665876a83SMatthew G. Knepley 
119765876a83SMatthew G. Knepley // Volumetric Strain
119865876a83SMatthew G. Knepley static PetscErrorCode cryer_3d_eps(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
119965876a83SMatthew G. Knepley {
120065876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
120165876a83SMatthew G. Knepley   Parameter     *param;
120265876a83SMatthew G. Knepley 
12035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
120465876a83SMatthew G. Knepley   if (time <= 0.0) {
12055f80ce2aSJacob Faibussowitsch     CHKERRQ(cryer_initial_eps(dim, time, x, Nc, u, ctx));
120665876a83SMatthew G. Knepley   } else {
120765876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
120865876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
120965876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
121065876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
121165876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
121265876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
121330602db0SMatthew G. Knepley     PetscReal   R_0   = user->xmax[1];  /* m */
121465876a83SMatthew G. Knepley     PetscInt    N     = user->niter, n;
121565876a83SMatthew G. Knepley 
121665876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
121765876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
121865876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
121965876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
122065876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
122165876a83SMatthew G. Knepley     PetscScalar u_inf = -P_0*R_0*(1. - 2.*nu) / (2.*G*(1. + nu));  /* m,       Cheng (7.388) */
122265876a83SMatthew G. Knepley 
122365876a83SMatthew G. Knepley     PetscReal   R      = PetscSqrtReal(x[0]*x[0] + x[1]*x[1] + x[2]*x[2]);
122465876a83SMatthew G. Knepley     PetscReal   R_star = R/R_0;
122565876a83SMatthew G. Knepley     PetscReal   tstar  = PetscRealPart(c*time) / PetscSqr(R_0);    /* - */
122665876a83SMatthew G. Knepley     PetscReal   divA_n = 0.0;
122765876a83SMatthew G. Knepley 
122865876a83SMatthew G. Knepley     if (R_star < PETSC_SMALL) {
122965876a83SMatthew G. Knepley       for (n = 1; n < N+1; ++n) {
123065876a83SMatthew G. Knepley         const PetscReal x_n = user->zeroArray[n-1];
123165876a83SMatthew G. Knepley         const PetscReal E_n = PetscRealPart(PetscSqr(1 - nu)*PetscSqr(1 + nu_u)*x_n - 18.0*(1 + nu)*(nu_u - nu)*(1 - nu_u));
123265876a83SMatthew G. Knepley 
123365876a83SMatthew G. Knepley         divA_n += PetscRealPart(
123465876a83SMatthew G. Knepley                   (12.0*(1.0 + nu)*(nu_u - nu))/((1.0 - 2.0*nu)*E_n*PetscSqr(R_star)*x_n*PetscSinReal(PetscSqrtReal(x_n))) *
123565876a83SMatthew G. Knepley                   (3.0*(nu_u - nu) * PetscSqrtReal(x_n) * ((2.0 + PetscSqr(R_star*PetscSqrtReal(x_n))) - 2.0*PetscCosReal(R_star * PetscSqrtReal(x_n)))
123665876a83SMatthew G. Knepley                   + 5.0 * (1.0 - nu)*(1.0 - 2.0*nu)*PetscPowRealInt(R_star, 2)*x_n*PetscSinReal(PetscSqrtReal(x_n))) * PetscExpReal(-x_n * tstar));
123765876a83SMatthew G. Knepley       }
123865876a83SMatthew G. Knepley     } else {
123965876a83SMatthew G. Knepley       for (n = 1; n < N+1; ++n) {
124065876a83SMatthew G. Knepley         const PetscReal x_n = user->zeroArray[n-1];
124165876a83SMatthew G. Knepley         const PetscReal E_n = PetscRealPart(PetscSqr(1 - nu)*PetscSqr(1 + nu_u)*x_n - 18.0*(1 + nu)*(nu_u - nu)*(1 - nu_u));
124265876a83SMatthew G. Knepley 
124365876a83SMatthew G. Knepley         divA_n += PetscRealPart(
124465876a83SMatthew G. Knepley                   (12.0*(1.0 + nu)*(nu_u - nu))/((1.0 - 2.0*nu)*E_n*PetscSqr(R_star)*x_n*PetscSinReal(PetscSqrtReal(x_n))) *
124565876a83SMatthew G. Knepley                   (3.0*(nu_u - nu) * PetscSqrtReal(x_n) * ((2.0/(R_star*PetscSqrtReal(x_n)) + R_star*PetscSqrtReal(x_n))*PetscSinReal(R_star * PetscSqrtReal(x_n)) - 2.0*PetscCosReal(R_star * PetscSqrtReal(x_n)))
124665876a83SMatthew G. Knepley                   + 5.0 * (1.0 - nu)*(1.0 - 2.0*nu)*PetscPowRealInt(R_star, 2)*x_n*PetscSinReal(PetscSqrtReal(x_n))) * PetscExpReal(-x_n * tstar));
124765876a83SMatthew G. Knepley       }
124865876a83SMatthew G. Knepley     }
124965876a83SMatthew G. Knepley     if (PetscAbsReal(divA_n) > 1e3) PetscPrintf(PETSC_COMM_SELF, "(%g, %g, %g) divA_n: %g\n", x[0], x[1], x[2], divA_n);
125065876a83SMatthew G. Knepley     u[0] = PetscRealPart(u_inf)/R_0 * (3.0 - divA_n);
125165876a83SMatthew G. Knepley   }
125265876a83SMatthew G. Knepley   return 0;
125365876a83SMatthew G. Knepley }
125465876a83SMatthew G. Knepley 
125565876a83SMatthew G. Knepley // Pressure
125665876a83SMatthew G. Knepley static PetscErrorCode cryer_3d_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
125765876a83SMatthew G. Knepley {
125865876a83SMatthew G. Knepley   AppCtx        *user = (AppCtx *) ctx;
125965876a83SMatthew G. Knepley   Parameter     *param;
126065876a83SMatthew G. Knepley 
12615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
126265876a83SMatthew G. Knepley   if (time <= 0.0) {
12635f80ce2aSJacob Faibussowitsch     CHKERRQ(cryer_drainage_pressure(dim, time, x, Nc, u, ctx));
126465876a83SMatthew G. Knepley   } else {
126565876a83SMatthew G. Knepley     PetscScalar alpha = param->alpha; /* -  */
126665876a83SMatthew G. Knepley     PetscScalar K_u   = param->K_u;   /* Pa */
126765876a83SMatthew G. Knepley     PetscScalar M     = param->M;     /* Pa */
126865876a83SMatthew G. Knepley     PetscScalar G     = param->mu;    /* Pa */
126965876a83SMatthew G. Knepley     PetscScalar P_0   = param->P_0;   /* Pa */
127030602db0SMatthew G. Knepley     PetscReal   R_0   = user->xmax[1];  /* m */
127165876a83SMatthew G. Knepley     PetscScalar kappa = param->k / param->mu_f;    /* m^2 / (Pa s) */
127265876a83SMatthew G. Knepley     PetscInt    N     = user->niter, n;
127365876a83SMatthew G. Knepley 
127465876a83SMatthew G. Knepley     PetscScalar K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
127565876a83SMatthew G. Knepley     PetscScalar eta   = (3.0*alpha*G) / (3.0*K_d + 4.0*G);         /* -,       Cheng (B.11) */
127665876a83SMatthew G. Knepley     PetscScalar nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
127765876a83SMatthew G. Knepley     PetscScalar nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
127865876a83SMatthew G. Knepley     PetscScalar S     = (3.0*K_u + 4.0*G) / (M*(3.0*K_d + 4.0*G)); /* Pa^{-1}, Cheng (B.14) */
127965876a83SMatthew G. Knepley     PetscScalar c     = kappa / S;                                 /* m^2 / s, Cheng (B.16) */
128065876a83SMatthew G. Knepley     PetscScalar R     = PetscSqrtReal(x[0]*x[0] + x[1]*x[1] + x[2]*x[2]);
128165876a83SMatthew G. Knepley 
128265876a83SMatthew G. Knepley     PetscScalar R_star = R / R_0;
128365876a83SMatthew G. Knepley     PetscScalar t_star = PetscRealPart(c * time) / PetscSqr(R_0);
128465876a83SMatthew G. Knepley     PetscReal   A_x    = 0.0;
128565876a83SMatthew G. Knepley 
128665876a83SMatthew G. Knepley     for (n = 1; n < N+1; ++n) {
128765876a83SMatthew G. Knepley       const PetscReal x_n = user->zeroArray[n-1];
128865876a83SMatthew G. Knepley       const PetscReal E_n = PetscRealPart(PetscSqr(1 - nu)*PetscSqr(1 + nu_u)*x_n - 18.0*(1 + nu)*(nu_u - nu)*(1 - nu_u));
128965876a83SMatthew G. Knepley 
129065876a83SMatthew G. Knepley       A_x += PetscRealPart(((18.0*PetscSqr(nu_u - nu)) / (eta * E_n)) * (PetscSinReal(R_star * PetscSqrtReal(x_n)) / (R_star * PetscSinReal(PetscSqrtReal(x_n))) - 1.0) * PetscExpReal(-x_n * t_star)); /* Cheng (7.395) */
129165876a83SMatthew G. Knepley     }
129265876a83SMatthew G. Knepley     u[0] = P_0 * A_x;
129365876a83SMatthew G. Knepley   }
129465876a83SMatthew G. Knepley   return 0;
129565876a83SMatthew G. Knepley }
129665876a83SMatthew G. Knepley 
129765876a83SMatthew G. Knepley /* Boundary Kernels */
129865876a83SMatthew G. Knepley static void f0_terzaghi_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
129965876a83SMatthew G. Knepley                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
130065876a83SMatthew G. Knepley                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
130165876a83SMatthew G. Knepley                                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
130265876a83SMatthew G. Knepley {
130365876a83SMatthew G. Knepley   const PetscReal P = PetscRealPart(constants[5]);
130465876a83SMatthew G. Knepley 
130565876a83SMatthew G. Knepley   f0[0] = 0.0;
130665876a83SMatthew G. Knepley   f0[1] = P;
130765876a83SMatthew G. Knepley }
130865876a83SMatthew G. Knepley 
130945480ffeSMatthew G. Knepley #if 0
131065876a83SMatthew G. Knepley static void f0_mandel_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
131165876a83SMatthew G. Knepley                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
131265876a83SMatthew G. Knepley                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
131365876a83SMatthew G. Knepley                                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
131465876a83SMatthew G. Knepley {
131565876a83SMatthew G. Knepley   // Uniform stress distribution
131665876a83SMatthew G. Knepley   /* PetscScalar xmax =  0.5;
131765876a83SMatthew G. Knepley   PetscScalar xmin = -0.5;
131865876a83SMatthew G. Knepley   PetscScalar ymax =  0.5;
131965876a83SMatthew G. Knepley   PetscScalar ymin = -0.5;
132065876a83SMatthew G. Knepley   PetscScalar P = constants[5];
132165876a83SMatthew G. Knepley   PetscScalar aL = (xmax - xmin) / 2.0;
132265876a83SMatthew G. Knepley   PetscScalar sigma_zz = -1.0*P / aL; */
132365876a83SMatthew G. Knepley 
132465876a83SMatthew G. Knepley   // Analytical (parabolic) stress distribution
132565876a83SMatthew G. Knepley   PetscReal a1, a2, am;
132665876a83SMatthew G. Knepley   PetscReal y1, y2, ym;
132765876a83SMatthew G. Knepley 
132865876a83SMatthew G. Knepley   PetscInt NITER = 500;
132965876a83SMatthew G. Knepley   PetscReal EPS = 0.000001;
133065876a83SMatthew G. Knepley   PetscReal zeroArray[500]; /* NITER */
133165876a83SMatthew G. Knepley   PetscReal xmax =  1.0;
133265876a83SMatthew G. Knepley   PetscReal xmin =  0.0;
133365876a83SMatthew G. Knepley   PetscReal ymax =  0.1;
133465876a83SMatthew G. Knepley   PetscReal ymin =  0.0;
133565876a83SMatthew G. Knepley   PetscReal lower[2], upper[2];
133665876a83SMatthew G. Knepley 
133765876a83SMatthew G. Knepley   lower[0] = xmin - (xmax - xmin) / 2.0;
133865876a83SMatthew G. Knepley   lower[1] = ymin - (ymax - ymin) / 2.0;
133965876a83SMatthew G. Knepley   upper[0] = xmax - (xmax - xmin) / 2.0;
134065876a83SMatthew G. Knepley   upper[1] = ymax - (ymax - ymin) / 2.0;
134165876a83SMatthew G. Knepley 
134265876a83SMatthew G. Knepley   xmin = lower[0];
134365876a83SMatthew G. Knepley   ymin = lower[1];
134465876a83SMatthew G. Knepley   xmax = upper[0];
134565876a83SMatthew G. Knepley   ymax = upper[1];
134665876a83SMatthew G. Knepley 
134765876a83SMatthew G. Knepley   PetscScalar G     = constants[0];
134865876a83SMatthew G. Knepley   PetscScalar K_u   = constants[1];
134965876a83SMatthew G. Knepley   PetscScalar alpha = constants[2];
135065876a83SMatthew G. Knepley   PetscScalar M     = constants[3];
135165876a83SMatthew G. Knepley   PetscScalar kappa = constants[4];
135265876a83SMatthew G. Knepley   PetscScalar F     = constants[5];
135365876a83SMatthew G. Knepley 
135465876a83SMatthew G. Knepley   PetscScalar K_d = K_u - alpha*alpha*M;
135565876a83SMatthew G. Knepley   PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
135665876a83SMatthew G. Knepley   PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
135765876a83SMatthew G. Knepley   PetscReal   aL = (xmax - xmin) / 2.0;
135865876a83SMatthew G. Knepley   PetscReal   c = PetscRealPart(((2.0*kappa*G) * (1.0 - nu) * (nu_u - nu)) / (alpha*alpha * (1.0 - 2.0*nu) * (1.0 - nu_u)));
135965876a83SMatthew G. Knepley   PetscScalar B = (3.0 * (nu_u - nu)) / ( alpha * (1.0 - 2.0*nu) * (1.0 + nu_u));
136065876a83SMatthew G. Knepley   PetscScalar A1 = 3.0 / (B * (1.0 + nu_u));
136165876a83SMatthew G. Knepley   PetscScalar A2 = (alpha * (1.0 - 2.0*nu)) / (1.0 - nu);
136265876a83SMatthew G. Knepley 
136365876a83SMatthew G. Knepley   // Generate zero values
136465876a83SMatthew G. Knepley   for (PetscInt i=1; i < NITER+1; i++)
136565876a83SMatthew G. Knepley   {
136665876a83SMatthew G. Knepley     a1 = ((PetscReal) i - 1.0) * PETSC_PI * PETSC_PI / 4.0 + EPS;
136765876a83SMatthew G. Knepley     a2 = a1 + PETSC_PI/2;
136865876a83SMatthew G. Knepley     for (PetscInt j=0; j<NITER; j++)
136965876a83SMatthew G. Knepley     {
137065876a83SMatthew G. Knepley       y1 = PetscTanReal(a1) - PetscRealPart(A1/A2)*a1;
137165876a83SMatthew G. Knepley       y2 = PetscTanReal(a2) - PetscRealPart(A1/A2)*a2;
137265876a83SMatthew G. Knepley       am = (a1 + a2)/2.0;
137365876a83SMatthew G. Knepley       ym = PetscTanReal(am) - PetscRealPart(A1/A2)*am;
137465876a83SMatthew G. Knepley       if ((ym*y1) > 0)
137565876a83SMatthew G. Knepley       {
137665876a83SMatthew G. Knepley         a1 = am;
137765876a83SMatthew G. Knepley       } else {
137865876a83SMatthew G. Knepley         a2 = am;
137965876a83SMatthew G. Knepley       }
138065876a83SMatthew G. Knepley       if (PetscAbsReal(y2) < EPS)
138165876a83SMatthew G. Knepley       {
138265876a83SMatthew G. Knepley         am = a2;
138365876a83SMatthew G. Knepley       }
138465876a83SMatthew G. Knepley     }
138565876a83SMatthew G. Knepley     zeroArray[i-1] = am;
138665876a83SMatthew G. Knepley   }
138765876a83SMatthew G. Knepley 
138865876a83SMatthew G. Knepley   // Solution for sigma_zz
138965876a83SMatthew G. Knepley   PetscScalar A_x = 0.0;
139065876a83SMatthew G. Knepley   PetscScalar B_x = 0.0;
139165876a83SMatthew G. Knepley 
139265876a83SMatthew G. Knepley   for (PetscInt n=1; n < NITER+1; n++)
139365876a83SMatthew G. Knepley   {
139465876a83SMatthew G. Knepley     PetscReal alpha_n = zeroArray[n-1];
139565876a83SMatthew G. Knepley 
139665876a83SMatthew G. Knepley     A_x += ( PetscSinReal(alpha_n) / (alpha_n - PetscSinReal(alpha_n) * PetscCosReal(alpha_n))) * PetscCosReal( (alpha_n * x[0]) / aL) * PetscExpReal( -1.0*( (alpha_n*alpha_n*c*t)/(aL*aL)));
139765876a83SMatthew G. Knepley     B_x += ( (PetscSinReal(alpha_n) * PetscCosReal(alpha_n))/(alpha_n - PetscSinReal(alpha_n) * PetscCosReal(alpha_n))) * PetscExpReal( -1.0*( (alpha_n*alpha_n*c*t)/(aL*aL)));
139865876a83SMatthew G. Knepley   }
139965876a83SMatthew G. Knepley 
140065876a83SMatthew G. Knepley   PetscScalar sigma_zz = -1.0*(F/aL) - ((2.0*F)/aL) * (A2/A1) * A_x + ((2.0*F)/aL) * B_x;
140165876a83SMatthew G. Knepley 
140265876a83SMatthew G. Knepley   if (x[1] == ymax) {
140365876a83SMatthew G. Knepley     f0[1] += sigma_zz;
140465876a83SMatthew G. Knepley   } else if (x[1] == ymin) {
140565876a83SMatthew G. Knepley     f0[1] -= sigma_zz;
140665876a83SMatthew G. Knepley   }
140765876a83SMatthew G. Knepley }
140845480ffeSMatthew G. Knepley #endif
140965876a83SMatthew G. Knepley 
141065876a83SMatthew G. Knepley static void f0_cryer_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
141165876a83SMatthew G. Knepley                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
141265876a83SMatthew G. Knepley                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
141365876a83SMatthew G. Knepley                                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
141465876a83SMatthew G. Knepley {
141565876a83SMatthew G. Knepley   const PetscReal P_0 = PetscRealPart(constants[5]);
141665876a83SMatthew G. Knepley   PetscInt        d;
141765876a83SMatthew G. Knepley 
141865876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[d] = -P_0*n[d];
141965876a83SMatthew G. Knepley }
142065876a83SMatthew G. Knepley 
142165876a83SMatthew G. Knepley /* Standard Kernels - Residual */
142265876a83SMatthew G. Knepley /* f0_e */
142365876a83SMatthew G. Knepley static void f0_epsilon(PetscInt dim, PetscInt Nf, PetscInt NfAux,
142465876a83SMatthew G. Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
142565876a83SMatthew G. Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
142665876a83SMatthew G. Knepley                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
142765876a83SMatthew G. Knepley {
142865876a83SMatthew G. Knepley   PetscInt d;
142965876a83SMatthew G. Knepley 
143065876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
143165876a83SMatthew G. Knepley     f0[0] += u_x[d*dim+d];
143265876a83SMatthew G. Knepley   }
143365876a83SMatthew G. Knepley   f0[0] -= u[uOff[1]];
143465876a83SMatthew G. Knepley }
143565876a83SMatthew G. Knepley 
143665876a83SMatthew G. Knepley /* f0_p */
143765876a83SMatthew G. Knepley static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
143865876a83SMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
143965876a83SMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
144065876a83SMatthew G. Knepley                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
144165876a83SMatthew G. Knepley {
144265876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
144365876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
144465876a83SMatthew G. Knepley 
144565876a83SMatthew G. Knepley   f0[0] += alpha*u_t[uOff[1]];
144665876a83SMatthew G. Knepley   f0[0] += u_t[uOff[2]]/M;
144730602db0SMatthew G. Knepley   if (f0[0] != f0[0]) abort();
144865876a83SMatthew G. Knepley }
144965876a83SMatthew G. Knepley 
145065876a83SMatthew G. Knepley /* f1_u */
145165876a83SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
145265876a83SMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
145365876a83SMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
145465876a83SMatthew G. Knepley                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
145565876a83SMatthew G. Knepley {
145665876a83SMatthew G. Knepley   const PetscInt  Nc     = dim;
145765876a83SMatthew G. Knepley   const PetscReal G      = PetscRealPart(constants[0]);
145865876a83SMatthew G. Knepley   const PetscReal K_u    = PetscRealPart(constants[1]);
145965876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
146065876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
146165876a83SMatthew G. Knepley   const PetscReal K_d    = K_u - alpha*alpha*M;
146265876a83SMatthew G. Knepley   const PetscReal lambda = K_d - (2.0 * G) / 3.0;
146365876a83SMatthew G. Knepley   PetscInt        c, d;
146465876a83SMatthew G. Knepley 
146565876a83SMatthew G. Knepley   for (c = 0; c < Nc; ++c)
146665876a83SMatthew G. Knepley   {
146765876a83SMatthew G. Knepley     for (d = 0; d < dim; ++d)
146865876a83SMatthew G. Knepley     {
146965876a83SMatthew G. Knepley       f1[c*dim+d] -= G*(u_x[c*dim+d] + u_x[d*dim+c]);
147065876a83SMatthew G. Knepley     }
147165876a83SMatthew G. Knepley     f1[c*dim+c] -= lambda*u[uOff[1]];
147265876a83SMatthew G. Knepley     f1[c*dim+c] += alpha*u[uOff[2]];
147365876a83SMatthew G. Knepley   }
147465876a83SMatthew G. Knepley }
147565876a83SMatthew G. Knepley 
147665876a83SMatthew G. Knepley /* f1_p */
147765876a83SMatthew G. Knepley static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
147865876a83SMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
147965876a83SMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
148065876a83SMatthew G. Knepley                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
148165876a83SMatthew G. Knepley {
148265876a83SMatthew G. Knepley   const PetscReal kappa = PetscRealPart(constants[4]);
148365876a83SMatthew G. Knepley   PetscInt        d;
148465876a83SMatthew G. Knepley 
148565876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
148665876a83SMatthew G. Knepley     f1[d] += kappa*u_x[uOff_x[2]+d];
148765876a83SMatthew G. Knepley   }
148865876a83SMatthew G. Knepley }
148965876a83SMatthew G. Knepley 
149065876a83SMatthew G. Knepley /*
149165876a83SMatthew G. Knepley   \partial_df \phi_fc g_{fc,gc,df,dg} \partial_dg \phi_gc
149265876a83SMatthew G. Knepley 
149365876a83SMatthew G. Knepley   \partial_df \phi_fc \lambda \delta_{fc,df} \sum_gc \partial_dg \phi_gc \delta_{gc,dg}
149465876a83SMatthew G. Knepley   = \partial_fc \phi_fc \sum_gc \partial_gc \phi_gc
149565876a83SMatthew G. Knepley */
149665876a83SMatthew G. Knepley 
149765876a83SMatthew G. Knepley /* Standard Kernels - Jacobian */
149865876a83SMatthew G. Knepley /* g0_ee */
149965876a83SMatthew G. Knepley static void g0_ee(PetscInt dim, PetscInt Nf, PetscInt NfAux,
150065876a83SMatthew G. Knepley            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
150165876a83SMatthew G. Knepley            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
150265876a83SMatthew G. Knepley            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
150365876a83SMatthew G. Knepley {
150465876a83SMatthew G. Knepley   g0[0] = -1.0;
150565876a83SMatthew G. Knepley }
150665876a83SMatthew G. Knepley 
150765876a83SMatthew G. Knepley /* g0_pe */
150865876a83SMatthew G. Knepley static void g0_pe(PetscInt dim, PetscInt Nf, PetscInt NfAux,
150965876a83SMatthew G. Knepley            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
151065876a83SMatthew G. Knepley            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
151165876a83SMatthew G. Knepley            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
151265876a83SMatthew G. Knepley {
151365876a83SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[2]);
151465876a83SMatthew G. Knepley 
151565876a83SMatthew G. Knepley   g0[0] = u_tShift*alpha;
151665876a83SMatthew G. Knepley }
151765876a83SMatthew G. Knepley 
151865876a83SMatthew G. Knepley /* g0_pp */
151965876a83SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
152065876a83SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
152165876a83SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
152265876a83SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
152365876a83SMatthew G. Knepley {
152465876a83SMatthew G. Knepley   const PetscReal M = PetscRealPart(constants[3]);
152565876a83SMatthew G. Knepley 
152665876a83SMatthew G. Knepley   g0[0] = u_tShift/M;
152765876a83SMatthew G. Knepley }
152865876a83SMatthew G. Knepley 
152965876a83SMatthew G. Knepley /* g1_eu */
153065876a83SMatthew G. Knepley static void g1_eu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
153165876a83SMatthew G. Knepley            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
153265876a83SMatthew G. Knepley            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
153365876a83SMatthew G. Knepley            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
153465876a83SMatthew G. Knepley {
153565876a83SMatthew G. Knepley   PetscInt d;
153665876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
153765876a83SMatthew G. Knepley }
153865876a83SMatthew G. Knepley 
153965876a83SMatthew G. Knepley /* g2_ue */
154065876a83SMatthew G. Knepley static void g2_ue(PetscInt dim, PetscInt Nf, PetscInt NfAux,
154165876a83SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
154265876a83SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
154365876a83SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
154465876a83SMatthew G. Knepley {
154565876a83SMatthew G. Knepley   const PetscReal G      = PetscRealPart(constants[0]);
154665876a83SMatthew G. Knepley   const PetscReal K_u    = PetscRealPart(constants[1]);
154765876a83SMatthew G. Knepley   const PetscReal alpha  = PetscRealPart(constants[2]);
154865876a83SMatthew G. Knepley   const PetscReal M      = PetscRealPart(constants[3]);
154965876a83SMatthew G. Knepley   const PetscReal K_d    = K_u - alpha*alpha*M;
155065876a83SMatthew G. Knepley   const PetscReal lambda = K_d - (2.0 * G) / 3.0;
155165876a83SMatthew G. Knepley   PetscInt        d;
155265876a83SMatthew G. Knepley 
155365876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
155465876a83SMatthew G. Knepley     g2[d*dim + d] -= lambda;
155565876a83SMatthew G. Knepley   }
155665876a83SMatthew G. Knepley }
155765876a83SMatthew G. Knepley /* g2_up */
155865876a83SMatthew G. Knepley static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
155965876a83SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
156065876a83SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
156165876a83SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
156265876a83SMatthew G. Knepley {
156365876a83SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[2]);
156465876a83SMatthew G. Knepley   PetscInt        d;
156565876a83SMatthew G. Knepley 
156665876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
156765876a83SMatthew G. Knepley     g2[d*dim + d] += alpha;
156865876a83SMatthew G. Knepley   }
156965876a83SMatthew G. Knepley }
157065876a83SMatthew G. Knepley 
157165876a83SMatthew G. Knepley /* g3_uu */
157265876a83SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
157365876a83SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
157465876a83SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
157565876a83SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
157665876a83SMatthew G. Knepley {
157765876a83SMatthew G. Knepley   const PetscInt  Nc = dim;
157865876a83SMatthew G. Knepley   const PetscReal G  = PetscRealPart(constants[0]);
157965876a83SMatthew G. Knepley   PetscInt        c, d;
158065876a83SMatthew G. Knepley 
158165876a83SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
158265876a83SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
158365876a83SMatthew G. Knepley       g3[((c*Nc + c)*dim + d)*dim + d] -= G;
158465876a83SMatthew G. Knepley       g3[((c*Nc + d)*dim + d)*dim + c] -= G;
158565876a83SMatthew G. Knepley     }
158665876a83SMatthew G. Knepley   }
158765876a83SMatthew G. Knepley }
158865876a83SMatthew G. Knepley 
158965876a83SMatthew G. Knepley /* g3_pp */
159065876a83SMatthew G. Knepley static void g3_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
159165876a83SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
159265876a83SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
159365876a83SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
159465876a83SMatthew G. Knepley {
159565876a83SMatthew G. Knepley   const PetscReal kappa = PetscRealPart(constants[4]);
159665876a83SMatthew G. Knepley   PetscInt        d;
159765876a83SMatthew G. Knepley 
159865876a83SMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] += kappa;
159965876a83SMatthew G. Knepley }
160065876a83SMatthew G. Knepley 
160165876a83SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
160265876a83SMatthew G. Knepley {
160365876a83SMatthew G. Knepley   PetscInt sol;
160465876a83SMatthew G. Knepley   PetscErrorCode ierr;
160565876a83SMatthew G. Knepley 
160665876a83SMatthew G. Knepley   PetscFunctionBeginUser;
160765876a83SMatthew G. Knepley   options->solType   = SOL_QUADRATIC_TRIG;
160865876a83SMatthew G. Knepley   options->niter     = 500;
160965876a83SMatthew G. Knepley   options->eps       = PETSC_SMALL;
161024b15d09SMatthew G. Knepley   options->dtInitial = -1.0;
161165876a83SMatthew G. Knepley 
161265876a83SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Biot Poroelasticity Options", "DMPLEX");CHKERRQ(ierr);
16135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-niter", "Number of series term iterations in exact solutions", "ex53.c", options->niter, &options->niter, NULL));
161465876a83SMatthew G. Knepley   sol  = options->solType;
16155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEList("-sol_type", "Type of exact solution", "ex53.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL));
161665876a83SMatthew G. Knepley   options->solType = (SolutionType) sol;
16175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-eps", "Precision value for root finding", "ex53.c", options->eps, &options->eps, NULL));
16185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-dt_initial", "Override the initial timestep", "ex53.c", options->dtInitial, &options->dtInitial, NULL));
161965876a83SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
162065876a83SMatthew G. Knepley   PetscFunctionReturn(0);
162165876a83SMatthew G. Knepley }
162265876a83SMatthew G. Knepley 
162365876a83SMatthew G. Knepley static PetscErrorCode mandelZeros(MPI_Comm comm, AppCtx *ctx, Parameter *param)
162465876a83SMatthew G. Knepley {
162565876a83SMatthew G. Knepley   //PetscBag       bag;
162665876a83SMatthew G. Knepley   PetscReal a1, a2, am;
162765876a83SMatthew G. Knepley   PetscReal y1, y2, ym;
162865876a83SMatthew G. Knepley 
162965876a83SMatthew G. Knepley   PetscFunctionBeginUser;
16305f80ce2aSJacob Faibussowitsch   //CHKERRQ(PetscBagGetData(ctx->bag, (void **) &param));
163165876a83SMatthew G. Knepley   PetscInt NITER = ctx->niter;
163265876a83SMatthew G. Knepley   PetscReal EPS = ctx->eps;
163365876a83SMatthew G. Knepley   //const PetscScalar YMAX = param->ymax;
163465876a83SMatthew G. Knepley   //const PetscScalar YMIN = param->ymin;
163565876a83SMatthew G. Knepley   PetscScalar alpha = param->alpha;
163665876a83SMatthew G. Knepley   PetscScalar K_u = param->K_u;
163765876a83SMatthew G. Knepley   PetscScalar M = param->M;
163865876a83SMatthew G. Knepley   PetscScalar G = param->mu;
163965876a83SMatthew G. Knepley   //const PetscScalar k = param->k;
164065876a83SMatthew G. Knepley   //const PetscScalar mu_f = param->mu_f;
164165876a83SMatthew G. Knepley   //const PetscScalar P_0 = param->P_0;
164265876a83SMatthew G. Knepley 
164365876a83SMatthew G. Knepley   PetscScalar K_d = K_u - alpha*alpha*M;
164465876a83SMatthew G. Knepley   PetscScalar nu = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));
164565876a83SMatthew G. Knepley   PetscScalar nu_u = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));
164665876a83SMatthew G. Knepley   //const PetscScalar kappa = k / mu_f;
164765876a83SMatthew G. Knepley 
164865876a83SMatthew G. Knepley   // Generate zero values
164965876a83SMatthew G. Knepley   for (PetscInt i=1; i < NITER+1; i++)
165065876a83SMatthew G. Knepley   {
165165876a83SMatthew G. Knepley     a1 = ((PetscReal) i - 1.0) * PETSC_PI * PETSC_PI / 4.0 + EPS;
165265876a83SMatthew G. Knepley     a2 = a1 + PETSC_PI/2;
165365876a83SMatthew G. Knepley     am = a1;
165465876a83SMatthew G. Knepley     for (PetscInt j=0; j<NITER; j++)
165565876a83SMatthew G. Knepley     {
165665876a83SMatthew G. Knepley       y1 = PetscTanReal(a1) - PetscRealPart((1.0 - nu)/(nu_u - nu))*a1;
165765876a83SMatthew G. Knepley       y2 = PetscTanReal(a2) - PetscRealPart((1.0 - nu)/(nu_u - nu))*a2;
165865876a83SMatthew G. Knepley       am = (a1 + a2)/2.0;
165965876a83SMatthew G. Knepley       ym = PetscTanReal(am) - PetscRealPart((1.0 - nu)/(nu_u - nu))*am;
166065876a83SMatthew G. Knepley       if ((ym*y1) > 0)
166165876a83SMatthew G. Knepley       {
166265876a83SMatthew G. Knepley         a1 = am;
166365876a83SMatthew G. Knepley       } else {
166465876a83SMatthew G. Knepley         a2 = am;
166565876a83SMatthew G. Knepley       }
166665876a83SMatthew G. Knepley       if (PetscAbsReal(y2) < EPS)
166765876a83SMatthew G. Knepley       {
166865876a83SMatthew G. Knepley         am = a2;
166965876a83SMatthew G. Knepley       }
167065876a83SMatthew G. Knepley     }
167165876a83SMatthew G. Knepley     ctx->zeroArray[i-1] = am;
167265876a83SMatthew G. Knepley   }
167365876a83SMatthew G. Knepley   PetscFunctionReturn(0);
167465876a83SMatthew G. Knepley }
167565876a83SMatthew G. Knepley 
167665876a83SMatthew G. Knepley static PetscReal CryerFunction(PetscReal nu_u, PetscReal nu, PetscReal x)
167765876a83SMatthew G. Knepley {
167865876a83SMatthew G. Knepley   return PetscTanReal(PetscSqrtReal(x))*(6.0*(nu_u - nu) - (1.0 - nu)*(1.0 + nu_u)*x) - (6.0*(nu_u - nu)*PetscSqrtReal(x));
167965876a83SMatthew G. Knepley }
168065876a83SMatthew G. Knepley 
168165876a83SMatthew G. Knepley static PetscErrorCode cryerZeros(MPI_Comm comm, AppCtx *ctx, Parameter *param)
168265876a83SMatthew G. Knepley {
168365876a83SMatthew G. Knepley   PetscReal   alpha = PetscRealPart(param->alpha); /* -  */
168465876a83SMatthew G. Knepley   PetscReal   K_u   = PetscRealPart(param->K_u);   /* Pa */
168565876a83SMatthew G. Knepley   PetscReal   M     = PetscRealPart(param->M);     /* Pa */
168665876a83SMatthew G. Knepley   PetscReal   G     = PetscRealPart(param->mu);    /* Pa */
168765876a83SMatthew G. Knepley   PetscInt    N     = ctx->niter, n;
168865876a83SMatthew G. Knepley 
168965876a83SMatthew G. Knepley   PetscReal   K_d   = K_u - alpha*alpha*M;                       /* Pa,      Cheng (B.5)  */
169065876a83SMatthew G. Knepley   PetscReal   nu    = (3.0*K_d - 2.0*G) / (2.0*(3.0*K_d + G));   /* -,       Cheng (B.8)  */
169165876a83SMatthew G. Knepley   PetscReal   nu_u  = (3.0*K_u - 2.0*G) / (2.0*(3.0*K_u + G));   /* -,       Cheng (B.9)  */
169265876a83SMatthew G. Knepley 
169365876a83SMatthew G. Knepley   PetscFunctionBeginUser;
169465876a83SMatthew G. Knepley   for (n = 1; n < N+1; ++n) {
169565876a83SMatthew G. Knepley     PetscReal tol = PetscPowReal(n, 1.5)*ctx->eps;
169665876a83SMatthew G. Knepley     PetscReal a1 = 0., a2 = 0., am = 0.;
169765876a83SMatthew G. Knepley     PetscReal y1, y2, ym;
169865876a83SMatthew G. Knepley     PetscInt  j, k = n-1;
169965876a83SMatthew G. Knepley 
170065876a83SMatthew G. Knepley     y1 = y2 = 1.;
170165876a83SMatthew G. Knepley     while (y1*y2 > 0) {
170265876a83SMatthew G. Knepley       ++k;
170365876a83SMatthew G. Knepley       a1 = PetscSqr(n*PETSC_PI) - k*PETSC_PI;
170465876a83SMatthew G. Knepley       a2 = PetscSqr(n*PETSC_PI) + k*PETSC_PI;
170565876a83SMatthew G. Knepley       y1 = CryerFunction(nu_u, nu, a1);
170665876a83SMatthew G. Knepley       y2 = CryerFunction(nu_u, nu, a2);
170765876a83SMatthew G. Knepley     }
170865876a83SMatthew G. Knepley     for (j = 0; j < 50000; ++j) {
170965876a83SMatthew G. Knepley       y1 = CryerFunction(nu_u, nu, a1);
171065876a83SMatthew G. Knepley       y2 = CryerFunction(nu_u, nu, a2);
17113c633725SBarry Smith       PetscCheck(y1*y2 <= 0,comm, PETSC_ERR_PLIB, "Invalid root finding initialization for root %D, (%g, %g)--(%g, %g)", n, a1, y1, a2, y2);
171265876a83SMatthew G. Knepley       am = (a1 + a2) / 2.0;
171365876a83SMatthew G. Knepley       ym = CryerFunction(nu_u, nu, am);
171465876a83SMatthew G. Knepley       if ((ym * y1) < 0) a2 = am;
171565876a83SMatthew G. Knepley       else               a1 = am;
171665876a83SMatthew G. Knepley       if (PetscAbsScalar(ym) < tol) break;
171765876a83SMatthew G. Knepley     }
17183c633725SBarry Smith     PetscCheck(PetscAbsScalar(ym) < tol,comm, PETSC_ERR_PLIB, "Root finding did not converge for root %D (%g)", n, PetscAbsScalar(ym));
171965876a83SMatthew G. Knepley     ctx->zeroArray[n-1] = am;
172065876a83SMatthew G. Knepley   }
172165876a83SMatthew G. Knepley   PetscFunctionReturn(0);
172265876a83SMatthew G. Knepley }
172365876a83SMatthew G. Knepley 
172465876a83SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
172565876a83SMatthew G. Knepley {
172665876a83SMatthew G. Knepley   PetscBag       bag;
172765876a83SMatthew G. Knepley   Parameter     *p;
172865876a83SMatthew G. Knepley 
172965876a83SMatthew G. Knepley   PetscFunctionBeginUser;
173065876a83SMatthew G. Knepley   /* setup PETSc parameter bag */
17315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(ctx->bag,(void**)&p));
17325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagSetName(ctx->bag,"par","Poroelastic Parameters"));
173365876a83SMatthew G. Knepley   bag  = ctx->bag;
173465876a83SMatthew G. Knepley   if (ctx->solType == SOL_TERZAGHI) {
173565876a83SMatthew G. Knepley     // Realistic values - Terzaghi
17365f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu,     3.0,                 "mu",    "Shear Modulus, Pa"));
17375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->K_u,    9.76,                "K_u",   "Undrained Bulk Modulus, Pa"));
17385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->alpha,  0.6,                 "alpha", "Biot Effective Stress Coefficient, -"));
17395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->M,      16.0,                "M",     "Biot Modulus, Pa"));
17405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->k,      1.5,                 "k",     "Isotropic Permeability, m**2"));
17415f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu_f,   1.0,                 "mu_f",  "Fluid Dynamic Viscosity, Pa*s"));
17425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->P_0,    1.0,                 "P_0",   "Magnitude of Vertical Stress, Pa"));
174365876a83SMatthew G. Knepley   } else if (ctx->solType == SOL_MANDEL) {
174465876a83SMatthew G. Knepley     // Realistic values - Mandel
17455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu,     0.75,                "mu",    "Shear Modulus, Pa"));
17465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->K_u,    2.6941176470588233,  "K_u",   "Undrained Bulk Modulus, Pa"));
17475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->alpha,  0.6,                 "alpha", "Biot Effective Stress Coefficient, -"));
17485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->M,      4.705882352941176,   "M",     "Biot Modulus, Pa"));
17495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->k,      1.5,                 "k",     "Isotropic Permeability, m**2"));
17505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu_f,   1.0,                 "mu_f",  "Fluid Dynamic Viscosity, Pa*s"));
17515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->P_0,    1.0,                 "P_0",   "Magnitude of Vertical Stress, Pa"));
175265876a83SMatthew G. Knepley   } else if (ctx->solType == SOL_CRYER) {
175365876a83SMatthew G. Knepley     // Realistic values - Mandel
17545f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu,     0.75,                "mu",    "Shear Modulus, Pa"));
17555f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->K_u,    2.6941176470588233,  "K_u",   "Undrained Bulk Modulus, Pa"));
17565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->alpha,  0.6,                 "alpha", "Biot Effective Stress Coefficient, -"));
17575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->M,      4.705882352941176,   "M",     "Biot Modulus, Pa"));
17585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->k,      1.5,                 "k",     "Isotropic Permeability, m**2"));
17595f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu_f,   1.0,                 "mu_f",  "Fluid Dynamic Viscosity, Pa*s"));
17605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->P_0,    1.0,                 "P_0",   "Magnitude of Vertical Stress, Pa"));
176165876a83SMatthew G. Knepley   } else {
176265876a83SMatthew G. Knepley     // Nonsense values
17635f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu,     1.0,                 "mu",    "Shear Modulus, Pa"));
17645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->K_u,    1.0,                 "K_u",   "Undrained Bulk Modulus, Pa"));
17655f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->alpha,  1.0,                 "alpha", "Biot Effective Stress Coefficient, -"));
17665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->M,      1.0,                 "M",     "Biot Modulus, Pa"));
17675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->k,      1.0,                 "k",     "Isotropic Permeability, m**2"));
17685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->mu_f,   1.0,                 "mu_f",  "Fluid Dynamic Viscosity, Pa*s"));
17695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscBagRegisterScalar(bag, &p->P_0,    1.0,                 "P_0",   "Magnitude of Vertical Stress, Pa"));
177065876a83SMatthew G. Knepley   }
17715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagSetFromOptions(bag));
177265876a83SMatthew G. Knepley   {
177365876a83SMatthew G. Knepley     PetscScalar K_d  = p->K_u - p->alpha*p->alpha*p->M;
177465876a83SMatthew G. Knepley     PetscScalar nu_u = (3.0*p->K_u - 2.0*p->mu) / (2.0*(3.0*p->K_u + p->mu));
177565876a83SMatthew G. Knepley     PetscScalar nu   = (3.0*K_d - 2.0*p->mu) / (2.0*(3.0*K_d + p->mu));
177665876a83SMatthew G. Knepley     PetscScalar S    = (3.0*p->K_u + 4.0*p->mu) / (p->M*(3.0*K_d + 4.0*p->mu));
177765876a83SMatthew G. Knepley     PetscReal   c    = PetscRealPart((p->k/p->mu_f) / S);
177865876a83SMatthew G. Knepley 
177965876a83SMatthew G. Knepley     PetscViewer       viewer;
178065876a83SMatthew G. Knepley     PetscViewerFormat format;
178165876a83SMatthew G. Knepley     PetscBool         flg;
178265876a83SMatthew G. Knepley 
178365876a83SMatthew G. Knepley     switch (ctx->solType) {
178465876a83SMatthew G. Knepley       case SOL_QUADRATIC_LINEAR:
178565876a83SMatthew G. Knepley       case SOL_QUADRATIC_TRIG:
178630602db0SMatthew G. Knepley       case SOL_TRIG_LINEAR: ctx->t_r = PetscSqr(ctx->xmax[0] - ctx->xmin[0])/c; break;
178730602db0SMatthew G. Knepley       case SOL_TERZAGHI:    ctx->t_r = PetscSqr(2.0*(ctx->xmax[1] - ctx->xmin[1]))/c; break;
178830602db0SMatthew G. Knepley       case SOL_MANDEL:      ctx->t_r = PetscSqr(2.0*(ctx->xmax[1] - ctx->xmin[1]))/c; break;
178930602db0SMatthew G. Knepley       case SOL_CRYER:       ctx->t_r = PetscSqr(ctx->xmax[1])/c; break;
179098921bdaSJacob Faibussowitsch       default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType);
179165876a83SMatthew G. Knepley     }
17925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
179365876a83SMatthew G. Knepley     if (flg) {
17945f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPushFormat(viewer, format));
17955f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscBagView(bag, viewer));
17965f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerFlush(viewer));
17975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerPopFormat(viewer));
17985f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscViewerDestroy(&viewer));
17995f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "  Max displacement: %g %g\n", p->P_0*(ctx->xmax[1] - ctx->xmin[1])*(1. - 2.*nu_u)/(2.*p->mu*(1. - nu_u)), p->P_0*(ctx->xmax[1] - ctx->xmin[1])*(1. - 2.*nu)/(2.*p->mu*(1. - nu))));
18005f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "  Relaxation time: %g\n", ctx->t_r));
180165876a83SMatthew G. Knepley     }
180265876a83SMatthew G. Knepley   }
180365876a83SMatthew G. Knepley   PetscFunctionReturn(0);
180465876a83SMatthew G. Knepley }
180565876a83SMatthew G. Knepley 
180665876a83SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
180765876a83SMatthew G. Knepley {
180865876a83SMatthew G. Knepley   PetscFunctionBeginUser;
18095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
18105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
18115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
18125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(*dm, user));
18135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
18145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetBoundingBox(*dm, user->xmin, user->xmax));
181565876a83SMatthew G. Knepley   PetscFunctionReturn(0);
181665876a83SMatthew G. Knepley }
181765876a83SMatthew G. Knepley 
181865876a83SMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
181965876a83SMatthew G. Knepley {
182065876a83SMatthew G. Knepley   PetscErrorCode (*exact[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
182165876a83SMatthew G. Knepley   PetscErrorCode (*exact_t[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
182245480ffeSMatthew G. Knepley   PetscDS          ds;
182345480ffeSMatthew G. Knepley   DMLabel          label;
182445480ffeSMatthew G. Knepley   PetscWeakForm    wf;
182565876a83SMatthew G. Knepley   Parameter       *param;
182665876a83SMatthew G. Knepley   PetscInt         id_mandel[2];
182765876a83SMatthew G. Knepley   PetscInt         comp[1];
182865876a83SMatthew G. Knepley   PetscInt         comp_mandel[2];
182945480ffeSMatthew G. Knepley   PetscInt         dim, id, bd, f;
183065876a83SMatthew G. Knepley 
183165876a83SMatthew G. Knepley   PetscFunctionBeginUser;
18325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
18335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
18345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetSpatialDimension(ds, &dim));
18355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagGetData(user->bag, (void **) &param));
183665876a83SMatthew G. Knepley   exact_t[0] = exact_t[1] = exact_t[2] = zero;
183765876a83SMatthew G. Knepley 
183865876a83SMatthew G. Knepley   /* Setup Problem Formulation and Boundary Conditions */
183965876a83SMatthew G. Knepley   switch (user->solType) {
184065876a83SMatthew G. Knepley   case SOL_QUADRATIC_LINEAR:
18415f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, f0_quadratic_linear_u, f1_u));
18425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 1, f0_epsilon,            NULL));
18435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 2, f0_quadratic_linear_p, f1_p));
18445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL,  NULL,  NULL,  g3_uu));
18455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL,  NULL,  g2_ue, NULL));
18465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 2, NULL,  NULL,  g2_up, NULL));
18475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL,  g1_eu, NULL,  NULL));
18485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 1, g0_ee, NULL,  NULL,  NULL));
18495f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 1, g0_pe, NULL,  NULL,  NULL));
18505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 2, g0_pp, NULL,  NULL,  g3_pp));
185165876a83SMatthew G. Knepley     exact[0]   = quadratic_u;
185265876a83SMatthew G. Knepley     exact[1]   = linear_eps;
185365876a83SMatthew G. Knepley     exact[2]   = linear_linear_p;
185465876a83SMatthew G. Knepley     exact_t[2] = linear_linear_p_t;
185565876a83SMatthew G. Knepley 
185665876a83SMatthew G. Knepley     id = 1;
18575f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall displacement", label, 1, &id, 0, 0, NULL, (void (*)(void)) exact[0], NULL, user, NULL));
18585f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall pressure",     label, 1, &id, 2, 0, NULL, (void (*)(void)) exact[2], (void (*)(void)) exact_t[2], user, NULL));
185965876a83SMatthew G. Knepley     break;
186065876a83SMatthew G. Knepley   case SOL_TRIG_LINEAR:
18615f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_linear_u, f1_u));
18625f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 1, f0_epsilon,       NULL));
18635f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 2, f0_trig_linear_p, f1_p));
18645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL,  NULL,  NULL,  g3_uu));
18655f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL,  NULL,  g2_ue, NULL));
18665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 2, NULL,  NULL,  g2_up, NULL));
18675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL,  g1_eu, NULL,  NULL));
18685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 1, g0_ee, NULL,  NULL,  NULL));
18695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 1, g0_pe, NULL,  NULL,  NULL));
18705f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 2, g0_pp, NULL,  NULL,  g3_pp));
187165876a83SMatthew G. Knepley     exact[0]   = trig_u;
187265876a83SMatthew G. Knepley     exact[1]   = trig_eps;
187365876a83SMatthew G. Knepley     exact[2]   = trig_linear_p;
187465876a83SMatthew G. Knepley     exact_t[2] = trig_linear_p_t;
187565876a83SMatthew G. Knepley 
187665876a83SMatthew G. Knepley     id = 1;
18775f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall displacement", label, 1, &id, 0, 0, NULL, (void (*)(void)) exact[0], NULL, user, NULL));
18785f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall pressure",     label, 1, &id, 2, 0, NULL, (void (*)(void)) exact[2], (void (*)(void)) exact_t[2], user, NULL));
187965876a83SMatthew G. Knepley     break;
188065876a83SMatthew G. Knepley   case SOL_QUADRATIC_TRIG:
18815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, f0_quadratic_trig_u, f1_u));
18825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 1, f0_epsilon,          NULL));
18835f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 2, f0_quadratic_trig_p, f1_p));
18845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL,  NULL,  NULL,  g3_uu));
18855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL,  NULL,  g2_ue, NULL));
18865f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 2, NULL,  NULL,  g2_up, NULL));
18875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL,  g1_eu, NULL,  NULL));
18885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 1, g0_ee, NULL,  NULL,  NULL));
18895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 1, g0_pe, NULL,  NULL,  NULL));
18905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 2, g0_pp, NULL,  NULL,  g3_pp));
189165876a83SMatthew G. Knepley     exact[0]   = quadratic_u;
189265876a83SMatthew G. Knepley     exact[1]   = linear_eps;
189365876a83SMatthew G. Knepley     exact[2]   = linear_trig_p;
189465876a83SMatthew G. Knepley     exact_t[2] = linear_trig_p_t;
189565876a83SMatthew G. Knepley 
189665876a83SMatthew G. Knepley     id = 1;
18975f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall displacement", label, 1, &id, 0, 0, NULL, (void (*)(void)) exact[0], NULL, user, NULL));
18985f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall pressure",     label, 1, &id, 2, 0, NULL, (void (*)(void)) exact[2], (void (*)(void)) exact_t[2], user, NULL));
189965876a83SMatthew G. Knepley     break;
190065876a83SMatthew G. Knepley   case SOL_TERZAGHI:
19015f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, NULL, f1_u));
19025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 1, f0_epsilon,     NULL));
19035f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 2, f0_p,           f1_p));
19045f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL,  NULL,  NULL,  g3_uu));
19055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL,  NULL,  g2_ue, NULL));
19065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 2, NULL,  NULL,  g2_up, NULL));
19075f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL,  g1_eu, NULL,  NULL));
19085f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 1, g0_ee, NULL,  NULL,  NULL));
19095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 1, g0_pe,  NULL,  NULL,  NULL));
19105f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 2, g0_pp,  NULL,  NULL,  g3_pp));
191165876a83SMatthew G. Knepley 
191265876a83SMatthew G. Knepley     exact[0] = terzaghi_2d_u;
191365876a83SMatthew G. Knepley     exact[1] = terzaghi_2d_eps;
191465876a83SMatthew G. Knepley     exact[2] = terzaghi_2d_p;
191565876a83SMatthew G. Knepley     exact_t[0] = terzaghi_2d_u_t;
191665876a83SMatthew G. Knepley     exact_t[1] = terzaghi_2d_eps_t;
191765876a83SMatthew G. Knepley     exact_t[2] = terzaghi_2d_p_t;
191865876a83SMatthew G. Knepley 
191965876a83SMatthew G. Knepley     id = 1;
19205f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL, "vertical stress",   label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
19215f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
19225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_terzaghi_bd_u, 0, NULL));
192345480ffeSMatthew G. Knepley 
192465876a83SMatthew G. Knepley     id = 3;
192565876a83SMatthew G. Knepley     comp[0] = 1;
19265f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed base",      label, 1, &id, 0, 1, comp, (void (*)(void)) zero, NULL, user, NULL));
192765876a83SMatthew G. Knepley     id = 2;
192865876a83SMatthew G. Knepley     comp[0] = 0;
19295f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed side",      label, 1, &id, 0, 1, comp, (void (*)(void)) zero, NULL, user, NULL));
193065876a83SMatthew G. Knepley     id = 4;
193165876a83SMatthew G. Knepley     comp[0] = 0;
19325f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed side",      label, 1, &id, 0, 1, comp, (void (*)(void)) zero, NULL, user, NULL));
193365876a83SMatthew G. Knepley     id = 1;
19345f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "drained surface", label, 1, &id, 2, 0, NULL, (void (*)(void)) terzaghi_drainage_pressure, NULL, user, NULL));
193565876a83SMatthew G. Knepley     break;
193665876a83SMatthew G. Knepley   case SOL_MANDEL:
19375f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, NULL, f1_u));
19385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 1, f0_epsilon,     NULL));
19395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 2, f0_p,           f1_p));
19405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL,  NULL,  NULL,  g3_uu));
19415f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL,  NULL,  g2_ue, NULL));
19425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 2, NULL,  NULL,  g2_up, NULL));
19435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL,  g1_eu, NULL,  NULL));
19445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 1, g0_ee, NULL,  NULL,  NULL));
19455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 1, g0_pe,  NULL,  NULL,  NULL));
19465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 2, g0_pp,  NULL,  NULL,  g3_pp));
194765876a83SMatthew G. Knepley 
19485f80ce2aSJacob Faibussowitsch     CHKERRQ(mandelZeros(PETSC_COMM_WORLD, user, param));
194965876a83SMatthew G. Knepley 
195065876a83SMatthew G. Knepley     exact[0] = mandel_2d_u;
195165876a83SMatthew G. Knepley     exact[1] = mandel_2d_eps;
195265876a83SMatthew G. Knepley     exact[2] = mandel_2d_p;
195365876a83SMatthew G. Knepley     exact_t[0] = mandel_2d_u_t;
195465876a83SMatthew G. Knepley     exact_t[1] = mandel_2d_eps_t;
195565876a83SMatthew G. Knepley     exact_t[2] = mandel_2d_p_t;
195665876a83SMatthew G. Knepley 
195765876a83SMatthew G. Knepley     id_mandel[0] = 3;
195865876a83SMatthew G. Knepley     id_mandel[1] = 1;
195965876a83SMatthew G. Knepley     //comp[0] = 1;
196065876a83SMatthew G. Knepley     comp_mandel[0] = 0;
196165876a83SMatthew G. Knepley     comp_mandel[1] = 1;
19625f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "vertical stress", label, 2, id_mandel, 0, 2, comp_mandel, (void (*)(void)) mandel_2d_u, NULL, user, NULL));
19635f80ce2aSJacob Faibussowitsch     //CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL, "vertical stress", "marker", 0, 1, comp, NULL, 2, id_mandel, user));
19645f80ce2aSJacob Faibussowitsch     //CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed base", "marker", 0, 1, comp, (void (*)(void)) zero, 2, id_mandel, user));
19655f80ce2aSJacob Faibussowitsch     //CHKERRQ(PetscDSSetBdResidual(ds, 0, f0_mandel_bd_u, NULL));
196665876a83SMatthew G. Knepley 
196765876a83SMatthew G. Knepley     id_mandel[0] = 2;
196865876a83SMatthew G. Knepley     id_mandel[1] = 4;
19695f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "drained surface", label, 2, id_mandel, 2, 0, NULL, (void (*)(void)) zero, NULL, user, NULL));
197065876a83SMatthew G. Knepley     break;
197165876a83SMatthew G. Knepley   case SOL_CRYER:
19725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 0, NULL, f1_u));
19735f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 1, f0_epsilon,     NULL));
19745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetResidual(ds, 2, f0_p,           f1_p));
19755f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL,  NULL,  NULL,  g3_uu));
19765f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL,  NULL,  g2_ue, NULL));
19775f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 0, 2, NULL,  NULL,  g2_up, NULL));
19785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL,  g1_eu, NULL,  NULL));
19795f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 1, 1, g0_ee, NULL,  NULL,  NULL));
19805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 1, g0_pe,  NULL,  NULL,  NULL));
19815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetJacobian(ds, 2, 2, g0_pp,  NULL,  NULL,  g3_pp));
198265876a83SMatthew G. Knepley 
19835f80ce2aSJacob Faibussowitsch     CHKERRQ(cryerZeros(PETSC_COMM_WORLD, user, param));
198465876a83SMatthew G. Knepley 
198565876a83SMatthew G. Knepley     exact[0] = cryer_3d_u;
198665876a83SMatthew G. Knepley     exact[1] = cryer_3d_eps;
198765876a83SMatthew G. Knepley     exact[2] = cryer_3d_p;
198865876a83SMatthew G. Knepley 
198965876a83SMatthew G. Knepley     id = 1;
19905f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL,   "normal stress",   label, 1, &id, 0, 0, NULL, NULL,                                     NULL, user, &bd));
19915f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
19925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_cryer_bd_u, 0, NULL));
199345480ffeSMatthew G. Knepley 
19945f80ce2aSJacob Faibussowitsch     CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "drained surface", label, 1, &id, 2, 0, NULL, (void (*)(void)) cryer_drainage_pressure, NULL, user, NULL));
199565876a83SMatthew G. Knepley     break;
199698921bdaSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType);
199765876a83SMatthew G. Knepley   }
199865876a83SMatthew G. Knepley   for (f = 0; f < 3; ++f) {
19995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetExactSolution(ds, f, exact[f], user));
20005f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, f, exact_t[f], user));
200165876a83SMatthew G. Knepley   }
200265876a83SMatthew G. Knepley 
200365876a83SMatthew G. Knepley   /* Setup constants */
200465876a83SMatthew G. Knepley   {
200565876a83SMatthew G. Knepley     PetscScalar constants[6];
200665876a83SMatthew G. Knepley     constants[0] = param->mu;            /* shear modulus, Pa */
200765876a83SMatthew G. Knepley     constants[1] = param->K_u;           /* undrained bulk modulus, Pa */
200865876a83SMatthew G. Knepley     constants[2] = param->alpha;         /* Biot effective stress coefficient, - */
200965876a83SMatthew G. Knepley     constants[3] = param->M;             /* Biot modulus, Pa */
201065876a83SMatthew G. Knepley     constants[4] = param->k/param->mu_f; /* Darcy coefficient, m**2 / Pa*s */
201165876a83SMatthew G. Knepley     constants[5] = param->P_0;           /* Magnitude of Vertical Stress, Pa */
20125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetConstants(ds, 6, constants));
201365876a83SMatthew G. Knepley   }
201465876a83SMatthew G. Knepley   PetscFunctionReturn(0);
201565876a83SMatthew G. Knepley }
201665876a83SMatthew G. Knepley 
20178cda7954SMatthew G. Knepley static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
201865876a83SMatthew G. Knepley {
201965876a83SMatthew G. Knepley   PetscFunctionBegin;
20205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateRigidBody(dm, origField, nullspace));
202165876a83SMatthew G. Knepley   PetscFunctionReturn(0);
202265876a83SMatthew G. Knepley }
202365876a83SMatthew G. Knepley 
202430602db0SMatthew G. Knepley static PetscErrorCode SetupFE(DM dm, PetscInt Nf, PetscInt Nc[], const char *name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx)
202565876a83SMatthew G. Knepley {
202665876a83SMatthew G. Knepley   AppCtx         *user = (AppCtx *) ctx;
202765876a83SMatthew G. Knepley   DM              cdm  = dm;
202865876a83SMatthew G. Knepley   PetscFE         fe;
202965876a83SMatthew G. Knepley   PetscQuadrature q = NULL;
203065876a83SMatthew G. Knepley   char            prefix[PETSC_MAX_PATH_LEN];
203165876a83SMatthew G. Knepley   PetscInt        dim, f;
203230602db0SMatthew G. Knepley   PetscBool       simplex;
203365876a83SMatthew G. Knepley 
203465876a83SMatthew G. Knepley   PetscFunctionBegin;
203565876a83SMatthew G. Knepley   /* Create finite element */
20365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
20375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsSimplex(dm, &simplex));
203865876a83SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
20395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name[f]));
20405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, name[f] ? prefix : NULL, -1, &fe));
20415f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) fe, name[f]));
20425f80ce2aSJacob Faibussowitsch     if (!q) CHKERRQ(PetscFEGetQuadrature(fe, &q));
20435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFESetQuadrature(fe, q));
20445f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetField(dm, f, NULL, (PetscObject) fe));
20455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEDestroy(&fe));
204665876a83SMatthew G. Knepley   }
20475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
20485f80ce2aSJacob Faibussowitsch   CHKERRQ((*setup)(dm, user));
204965876a83SMatthew G. Knepley   while (cdm) {
20505f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm, cdm));
20515f80ce2aSJacob Faibussowitsch     if (0) CHKERRQ(DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace));
205265876a83SMatthew G. Knepley     /* TODO: Check whether the boundary of coarse meshes is marked */
20535f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
205465876a83SMatthew G. Knepley   }
20555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
205665876a83SMatthew G. Knepley   PetscFunctionReturn(0);
205765876a83SMatthew G. Knepley }
205865876a83SMatthew G. Knepley 
205965876a83SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u)
206065876a83SMatthew G. Knepley {
206165876a83SMatthew G. Knepley   DM             dm;
206265876a83SMatthew G. Knepley   PetscReal      t;
206365876a83SMatthew G. Knepley 
206465876a83SMatthew G. Knepley   PetscFunctionBegin;
20655f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
20665f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts, &t));
206765876a83SMatthew G. Knepley   if (t <= 0.0) {
206865876a83SMatthew G. Knepley     PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
206965876a83SMatthew G. Knepley     void            *ctxs[3];
207065876a83SMatthew G. Knepley     AppCtx          *ctx;
207165876a83SMatthew G. Knepley 
20725f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetApplicationContext(dm, &ctx));
207365876a83SMatthew G. Knepley     switch (ctx->solType) {
207465876a83SMatthew G. Knepley       case SOL_TERZAGHI:
207565876a83SMatthew G. Knepley         funcs[0] = terzaghi_initial_u;         ctxs[0] = ctx;
207665876a83SMatthew G. Knepley         funcs[1] = terzaghi_initial_eps;       ctxs[1] = ctx;
207765876a83SMatthew G. Knepley         funcs[2] = terzaghi_drainage_pressure; ctxs[2] = ctx;
20785f80ce2aSJacob Faibussowitsch         CHKERRQ(DMProjectFunction(dm, t, funcs, ctxs, INSERT_VALUES, u));
207965876a83SMatthew G. Knepley         break;
208065876a83SMatthew G. Knepley       case SOL_MANDEL:
208165876a83SMatthew G. Knepley         funcs[0] = mandel_initial_u;         ctxs[0] = ctx;
208265876a83SMatthew G. Knepley         funcs[1] = mandel_initial_eps;       ctxs[1] = ctx;
208365876a83SMatthew G. Knepley         funcs[2] = mandel_drainage_pressure; ctxs[2] = ctx;
20845f80ce2aSJacob Faibussowitsch         CHKERRQ(DMProjectFunction(dm, t, funcs, ctxs, INSERT_VALUES, u));
208565876a83SMatthew G. Knepley         break;
208665876a83SMatthew G. Knepley       case SOL_CRYER:
208765876a83SMatthew G. Knepley         funcs[0] = cryer_initial_u;         ctxs[0] = ctx;
208865876a83SMatthew G. Knepley         funcs[1] = cryer_initial_eps;       ctxs[1] = ctx;
208965876a83SMatthew G. Knepley         funcs[2] = cryer_drainage_pressure; ctxs[2] = ctx;
20905f80ce2aSJacob Faibussowitsch         CHKERRQ(DMProjectFunction(dm, t, funcs, ctxs, INSERT_VALUES, u));
209165876a83SMatthew G. Knepley         break;
209265876a83SMatthew G. Knepley       default:
20935f80ce2aSJacob Faibussowitsch         CHKERRQ(DMComputeExactSolution(dm, t, u, NULL));
209465876a83SMatthew G. Knepley     }
209565876a83SMatthew G. Knepley   } else {
20965f80ce2aSJacob Faibussowitsch     CHKERRQ(DMComputeExactSolution(dm, t, u, NULL));
209765876a83SMatthew G. Knepley   }
209865876a83SMatthew G. Knepley   PetscFunctionReturn(0);
209965876a83SMatthew G. Knepley }
210065876a83SMatthew G. Knepley 
210165876a83SMatthew G. Knepley /* Need to create Viewer each time because HDF5 can get corrupted */
210265876a83SMatthew G. Knepley static PetscErrorCode SolutionMonitor(TS ts, PetscInt steps, PetscReal time, Vec u, void *mctx)
210365876a83SMatthew G. Knepley {
210465876a83SMatthew G. Knepley   DM                dm;
210565876a83SMatthew G. Knepley   Vec               exact;
210665876a83SMatthew G. Knepley   PetscViewer       viewer;
210765876a83SMatthew G. Knepley   PetscViewerFormat format;
210865876a83SMatthew G. Knepley   PetscOptions      options;
210965876a83SMatthew G. Knepley   const char       *prefix;
211065876a83SMatthew G. Knepley 
211165876a83SMatthew G. Knepley   PetscFunctionBegin;
21125f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
21135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetOptions((PetscObject) ts, &options));
21145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) ts, &prefix));
21155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetViewer(PetscObjectComm((PetscObject) ts), options, prefix, "-monitor_solution", &viewer, &format, NULL));
21165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &exact));
21175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMComputeExactSolution(dm, time, exact, NULL));
21185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOutputSequenceNumber(dm, steps, time));
21195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(exact, viewer));
21205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(u, viewer));
21215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &exact));
212265876a83SMatthew G. Knepley   {
212365876a83SMatthew G. Knepley     PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx);
212465876a83SMatthew G. Knepley     void            **ectxs;
212565876a83SMatthew G. Knepley     PetscReal        *err;
212665876a83SMatthew G. Knepley     PetscInt          Nf, f;
212765876a83SMatthew G. Knepley 
21285f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetNumFields(dm, &Nf));
21295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err));
213065876a83SMatthew G. Knepley     {
213165876a83SMatthew G. Knepley       PetscInt Nds, s;
213265876a83SMatthew G. Knepley 
21335f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetNumDS(dm, &Nds));
213465876a83SMatthew G. Knepley       for (s = 0; s < Nds; ++s) {
213565876a83SMatthew G. Knepley         PetscDS         ds;
213665876a83SMatthew G. Knepley         DMLabel         label;
213765876a83SMatthew G. Knepley         IS              fieldIS;
213865876a83SMatthew G. Knepley         const PetscInt *fields;
213965876a83SMatthew G. Knepley         PetscInt        dsNf, f;
214065876a83SMatthew G. Knepley 
21415f80ce2aSJacob Faibussowitsch         CHKERRQ(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds));
21425f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscDSGetNumFields(ds, &dsNf));
21435f80ce2aSJacob Faibussowitsch         CHKERRQ(ISGetIndices(fieldIS, &fields));
214465876a83SMatthew G. Knepley         for (f = 0; f < dsNf; ++f) {
214565876a83SMatthew G. Knepley           const PetscInt field = fields[f];
21465f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field]));
214765876a83SMatthew G. Knepley         }
21485f80ce2aSJacob Faibussowitsch         CHKERRQ(ISRestoreIndices(fieldIS, &fields));
214965876a83SMatthew G. Knepley       }
215065876a83SMatthew G. Knepley     }
21515f80ce2aSJacob Faibussowitsch     CHKERRQ(DMComputeL2FieldDiff(dm, time, exacts, ectxs, u, err));
21525f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), "Time: %g L_2 Error: [", time));
215365876a83SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
21545f80ce2aSJacob Faibussowitsch       if (f) CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), ", "));
21555f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), "%g", (double) err[f]));
215665876a83SMatthew G. Knepley     }
21575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) ts), "]\n"));
21585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree3(exacts, ectxs, err));
215965876a83SMatthew G. Knepley   }
21605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
216165876a83SMatthew G. Knepley   PetscFunctionReturn(0);
216265876a83SMatthew G. Knepley }
216365876a83SMatthew G. Knepley 
216465876a83SMatthew G. Knepley static PetscErrorCode SetupMonitor(TS ts, AppCtx *ctx)
216565876a83SMatthew G. Knepley {
216665876a83SMatthew G. Knepley   PetscViewer       viewer;
216765876a83SMatthew G. Knepley   PetscViewerFormat format;
216865876a83SMatthew G. Knepley   PetscOptions      options;
216965876a83SMatthew G. Knepley   const char       *prefix;
217065876a83SMatthew G. Knepley   PetscBool         flg;
217165876a83SMatthew G. Knepley 
217265876a83SMatthew G. Knepley   PetscFunctionBegin;
21735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetOptions((PetscObject) ts, &options));
21745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) ts, &prefix));
21755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetViewer(PetscObjectComm((PetscObject) ts), options, prefix, "-monitor_solution", &viewer, &format, &flg));
21765f80ce2aSJacob Faibussowitsch   if (flg) CHKERRQ(TSMonitorSet(ts, SolutionMonitor, ctx, NULL));
21775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
217865876a83SMatthew G. Knepley   PetscFunctionReturn(0);
217965876a83SMatthew G. Knepley }
218065876a83SMatthew G. Knepley 
218165876a83SMatthew G. Knepley static PetscErrorCode TSAdaptChoose_Terzaghi(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
218265876a83SMatthew G. Knepley {
218365876a83SMatthew G. Knepley   static PetscReal dtTarget = -1.0;
218465876a83SMatthew G. Knepley   PetscReal        dtInitial;
218565876a83SMatthew G. Knepley   DM               dm;
218665876a83SMatthew G. Knepley   AppCtx          *ctx;
218765876a83SMatthew G. Knepley   PetscInt         step;
218865876a83SMatthew G. Knepley 
218965876a83SMatthew G. Knepley   PetscFunctionBegin;
21905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
21915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dm, &ctx));
21925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts, &step));
219324b15d09SMatthew G. Knepley   dtInitial = ctx->dtInitial < 0.0 ? 1.0e-4*ctx->t_r : ctx->dtInitial;
219465876a83SMatthew G. Knepley   if (!step) {
219565876a83SMatthew G. Knepley     if (PetscAbsReal(dtInitial - h) > PETSC_SMALL) {
219665876a83SMatthew G. Knepley       *accept  = PETSC_FALSE;
219765876a83SMatthew G. Knepley       *next_h  = dtInitial;
219865876a83SMatthew G. Knepley       dtTarget = h;
219965876a83SMatthew G. Knepley     } else {
220065876a83SMatthew G. Knepley       *accept  = PETSC_TRUE;
220165876a83SMatthew G. Knepley       *next_h  = dtTarget < 0.0 ? dtInitial : dtTarget;
220265876a83SMatthew G. Knepley       dtTarget = -1.0;
220365876a83SMatthew G. Knepley     }
220465876a83SMatthew G. Knepley   } else {
220565876a83SMatthew G. Knepley     *accept = PETSC_TRUE;
220665876a83SMatthew G. Knepley     *next_h = h;
220765876a83SMatthew G. Knepley   }
220865876a83SMatthew G. Knepley   *next_sc = 0;  /* Reuse the same order scheme */
220965876a83SMatthew G. Knepley   *wlte    = -1; /* Weighted local truncation error was not evaluated */
221065876a83SMatthew G. Knepley   *wltea   = -1; /* Weighted absolute local truncation error was not evaluated */
221165876a83SMatthew G. Knepley   *wlter   = -1; /* Weighted relative local truncation error was not evaluated */
221265876a83SMatthew G. Knepley   PetscFunctionReturn(0);
221365876a83SMatthew G. Knepley }
221465876a83SMatthew G. Knepley 
221565876a83SMatthew G. Knepley int main(int argc, char **argv)
221665876a83SMatthew G. Knepley {
221765876a83SMatthew G. Knepley   AppCtx         ctx;       /* User-defined work context */
221865876a83SMatthew G. Knepley   DM             dm;        /* Problem specification */
221965876a83SMatthew G. Knepley   TS             ts;        /* Time Series / Nonlinear solver */
222065876a83SMatthew G. Knepley   Vec            u;         /* Solutions */
222165876a83SMatthew G. Knepley   const char    *name[3] = {"displacement", "tracestrain", "pressure"};
222265876a83SMatthew G. Knepley   PetscReal      t;
222330602db0SMatthew G. Knepley   PetscInt       dim, Nc[3];
222465876a83SMatthew G. Knepley 
2225*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
22265f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx));
22275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx.bag));
22285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(ctx.niter, &ctx.zeroArray));
22295f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
22305f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupParameters(PETSC_COMM_WORLD, &ctx));
223165876a83SMatthew G. Knepley   /* Primal System */
22325f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts));
22335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(dm, &ctx));
22345f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts, dm));
223565876a83SMatthew G. Knepley 
22365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
223730602db0SMatthew G. Knepley   Nc[0] = dim;
223865876a83SMatthew G. Knepley   Nc[1] = 1;
223965876a83SMatthew G. Knepley   Nc[2] = 1;
224065876a83SMatthew G. Knepley 
22415f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupFE(dm, 3, Nc, name, SetupPrimalProblem, &ctx));
22425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &u));
22435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
22445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
22455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
22465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
22475f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
22485f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetComputeInitialCondition(ts, SetInitialConditions));
22495f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupMonitor(ts, &ctx));
225065876a83SMatthew G. Knepley 
225165876a83SMatthew G. Knepley   if (ctx.solType != SOL_QUADRATIC_TRIG) {
225265876a83SMatthew G. Knepley     TSAdapt adapt;
225365876a83SMatthew G. Knepley 
22545f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetAdapt(ts, &adapt));
225565876a83SMatthew G. Knepley     adapt->ops->choose = TSAdaptChoose_Terzaghi;
225665876a83SMatthew G. Knepley   }
225765876a83SMatthew G. Knepley   if (ctx.solType == SOL_CRYER) {
225865876a83SMatthew G. Knepley     Mat          J;
225965876a83SMatthew G. Knepley     MatNullSpace sp;
226065876a83SMatthew G. Knepley 
22615f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetUp(ts));
22625f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetIJacobian(ts, &J, NULL, NULL, NULL));
22635f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateRigidBody(dm, 0, &sp));
22645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetNullSpace(J, sp));
22655f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNullSpaceDestroy(&sp));
226665876a83SMatthew G. Knepley   }
22675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts, &t));
22685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOutputSequenceNumber(dm, 0, t));
22695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSCheckFromOptions(ts, u));
22705f80ce2aSJacob Faibussowitsch   CHKERRQ(SetInitialConditions(ts, u));
22715f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) u, "solution"));
22725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts, u));
22735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSCheckFromOptions(ts, u));
22745f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts, &u));
22755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view"));
227665876a83SMatthew G. Knepley 
227765876a83SMatthew G. Knepley   /* Cleanup */
22785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
22795f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
22805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
22815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscBagDestroy(&ctx.bag));
22825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ctx.zeroArray));
2283*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
2284*b122ec5aSJacob Faibussowitsch   return 0;
228565876a83SMatthew G. Knepley }
228665876a83SMatthew G. Knepley 
228765876a83SMatthew G. Knepley /*TEST
228865876a83SMatthew G. Knepley 
228965876a83SMatthew G. Knepley   test:
229065876a83SMatthew G. Knepley     suffix: 2d_quad_linear
229165876a83SMatthew G. Knepley     requires: triangle
229265876a83SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 2 \
229365876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
229465876a83SMatthew G. Knepley       -dmts_check .0001 -ts_max_steps 5 -ts_monitor_extreme
229565876a83SMatthew G. Knepley 
229665876a83SMatthew G. Knepley   test:
229765876a83SMatthew G. Knepley     suffix: 3d_quad_linear
229865876a83SMatthew G. Knepley     requires: ctetgen
229930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -sol_type quadratic_linear -dm_refine 1 \
230065876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
230165876a83SMatthew G. Knepley       -dmts_check .0001 -ts_max_steps 5 -ts_monitor_extreme
230265876a83SMatthew G. Knepley 
230365876a83SMatthew G. Knepley   test:
230465876a83SMatthew G. Knepley     suffix: 2d_trig_linear
230565876a83SMatthew G. Knepley     requires: triangle
230665876a83SMatthew G. Knepley     args: -sol_type trig_linear -dm_refine 1 \
230765876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
230865876a83SMatthew G. Knepley       -dmts_check .0001 -ts_max_steps 5 -ts_dt 0.00001 -ts_monitor_extreme
230965876a83SMatthew G. Knepley 
231065876a83SMatthew G. Knepley   test:
231165876a83SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9, 2.1, 1.8]
231265876a83SMatthew G. Knepley     suffix: 2d_trig_linear_sconv
231365876a83SMatthew G. Knepley     requires: triangle
231465876a83SMatthew G. Knepley     args: -sol_type trig_linear -dm_refine 1 \
231565876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
231665876a83SMatthew G. Knepley       -convest_num_refine 1 -ts_convergence_estimate -ts_convergence_temporal 0 -ts_max_steps 1 -ts_dt 0.00001 -pc_type lu
231765876a83SMatthew G. Knepley 
231865876a83SMatthew G. Knepley   test:
231965876a83SMatthew G. Knepley     suffix: 3d_trig_linear
232065876a83SMatthew G. Knepley     requires: ctetgen
232130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -sol_type trig_linear -dm_refine 1 \
232265876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
232365876a83SMatthew G. Knepley       -dmts_check .0001 -ts_max_steps 2 -ts_monitor_extreme
232465876a83SMatthew G. Knepley 
232565876a83SMatthew G. Knepley   test:
232665876a83SMatthew G. Knepley     # -dm_refine 1 -convest_num_refine 2 gets L_2 convergence rate: [2.0, 2.1, 1.9]
232765876a83SMatthew G. Knepley     suffix: 3d_trig_linear_sconv
232865876a83SMatthew G. Knepley     requires: ctetgen
232930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -sol_type trig_linear -dm_refine 1 \
233065876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
233165876a83SMatthew G. Knepley       -convest_num_refine 1 -ts_convergence_estimate -ts_convergence_temporal 0 -ts_max_steps 1 -pc_type lu
233265876a83SMatthew G. Knepley 
233365876a83SMatthew G. Knepley   test:
233465876a83SMatthew G. Knepley     suffix: 2d_quad_trig
233565876a83SMatthew G. Knepley     requires: triangle
233665876a83SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_refine 2 \
233765876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
233865876a83SMatthew G. Knepley       -dmts_check .0001 -ts_max_steps 5 -ts_monitor_extreme
233965876a83SMatthew G. Knepley 
234065876a83SMatthew G. Knepley   test:
234165876a83SMatthew G. Knepley     # Using -dm_refine 4 gets the convergence rates to [0.95, 0.97, 0.90]
234265876a83SMatthew G. Knepley     suffix: 2d_quad_trig_tconv
234365876a83SMatthew G. Knepley     requires: triangle
234465876a83SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_refine 1 \
234565876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
234665876a83SMatthew G. Knepley       -convest_num_refine 3 -ts_convergence_estimate -ts_max_steps 5 -pc_type lu
234765876a83SMatthew G. Knepley 
234865876a83SMatthew G. Knepley   test:
234965876a83SMatthew G. Knepley     suffix: 3d_quad_trig
235065876a83SMatthew G. Knepley     requires: ctetgen
235130602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -sol_type quadratic_trig -dm_refine 1 \
235265876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
235365876a83SMatthew G. Knepley       -dmts_check .0001 -ts_max_steps 5 -ts_monitor_extreme
235465876a83SMatthew G. Knepley 
235565876a83SMatthew G. Knepley   test:
235665876a83SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 gets the convergence rates to [1.0, 1.0, 1.0]
235765876a83SMatthew G. Knepley     suffix: 3d_quad_trig_tconv
235865876a83SMatthew G. Knepley     requires: ctetgen
235930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -sol_type quadratic_trig -dm_refine 1 \
236065876a83SMatthew G. Knepley       -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
236165876a83SMatthew G. Knepley       -convest_num_refine 1 -ts_convergence_estimate -ts_max_steps 5 -pc_type lu
236265876a83SMatthew G. Knepley 
236330602db0SMatthew G. Knepley   testset:
236430602db0SMatthew G. Knepley     args: -sol_type terzaghi -dm_plex_simplex 0 -dm_plex_box_faces 1,8 -dm_plex_box_lower 0,0 -dm_plex_box_upper 10,10 -dm_plex_separate_marker \
236530602db0SMatthew G. Knepley           -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 -niter 16000 \
236630602db0SMatthew G. Knepley           -pc_type lu
236730602db0SMatthew G. Knepley 
236865876a83SMatthew G. Knepley     test:
236965876a83SMatthew G. Knepley       suffix: 2d_terzaghi
237030602db0SMatthew G. Knepley       requires: double
237130602db0SMatthew G. Knepley       args: -ts_dt 0.0028666667 -ts_max_steps 2 -ts_monitor -dmts_check .0001
237265876a83SMatthew G. Knepley 
237365876a83SMatthew G. Knepley     test:
237465876a83SMatthew G. Knepley       # -dm_plex_box_faces 1,64 -ts_max_steps 4 -convest_num_refine 3 gives L_2 convergence rate: [1.1, 1.1, 1.1]
237565876a83SMatthew G. Knepley       suffix: 2d_terzaghi_tconv
237630602db0SMatthew G. Knepley       args: -ts_dt 0.023 -ts_max_steps 2 -ts_convergence_estimate -convest_num_refine 1
237765876a83SMatthew G. Knepley 
237865876a83SMatthew G. Knepley     test:
237924b15d09SMatthew G. Knepley       # -dm_plex_box_faces 1,16 -convest_num_refine 4 gives L_2 convergence rate: [1.7, 1.2, 1.1]
238030602db0SMatthew G. Knepley       # if we add -displacement_petscspace_degree 3 -tracestrain_petscspace_degree 2 -pressure_petscspace_degree 2, we get [2.1, 1.6, 1.5]
238124b15d09SMatthew G. Knepley       suffix: 2d_terzaghi_sconv
238230602db0SMatthew G. Knepley       args: -ts_dt 1e-5 -dt_initial 1e-5 -ts_max_steps 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1
238330602db0SMatthew G. Knepley 
238430602db0SMatthew G. Knepley   testset:
238530602db0SMatthew G. Knepley     args: -sol_type mandel -dm_plex_simplex 0 -dm_plex_box_lower -0.5,-0.125 -dm_plex_box_upper 0.5,0.125 -dm_plex_separate_marker -dm_refine 1 \
238630602db0SMatthew G. Knepley           -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1 \
238730602db0SMatthew G. Knepley           -pc_type lu
238824b15d09SMatthew G. Knepley 
238924b15d09SMatthew G. Knepley     test:
239065876a83SMatthew G. Knepley       suffix: 2d_mandel
239130602db0SMatthew G. Knepley       requires: double
239230602db0SMatthew G. Knepley       args: -ts_dt 0.0028666667 -ts_max_steps 2 -ts_monitor -dmts_check .0001
239365876a83SMatthew G. Knepley 
239465876a83SMatthew G. Knepley     test:
239565876a83SMatthew G. Knepley       # -dm_refine 5 -ts_max_steps 4 -convest_num_refine 3 gives L_2 convergence rate: [0.26, -0.0058, 0.26]
239665876a83SMatthew G. Knepley       suffix: 2d_mandel_tconv
239730602db0SMatthew G. Knepley       args: -ts_dt 0.023 -ts_max_steps 2 -ts_convergence_estimate -convest_num_refine 1
239830602db0SMatthew G. Knepley 
239930602db0SMatthew G. Knepley   testset:
240030602db0SMatthew G. Knepley     requires: ctetgen !complex
240130602db0SMatthew G. Knepley     args: -sol_type cryer -dm_plex_dim 3 -dm_plex_shape ball \
240230602db0SMatthew G. Knepley           -displacement_petscspace_degree 2 -tracestrain_petscspace_degree 1 -pressure_petscspace_degree 1
240365876a83SMatthew G. Knepley 
240465876a83SMatthew G. Knepley     test:
240565876a83SMatthew G. Knepley       suffix: 3d_cryer
240630602db0SMatthew G. Knepley       args: -ts_dt 0.0028666667 -ts_max_time 0.014333 -ts_max_steps 2 -dmts_check .0001 \
240730602db0SMatthew G. Knepley             -pc_type svd
240865876a83SMatthew G. Knepley 
240965876a83SMatthew G. Knepley     test:
241065876a83SMatthew G. Knepley       # Displacement and Pressure converge. The analytic expression for trace strain is inaccurate at the origin
241165876a83SMatthew G. Knepley       # -bd_dm_refine 3 -ref_limit 0.00666667 -ts_max_steps 5 -convest_num_refine 2 gives L_2 convergence rate: [0.47, -0.43, 1.5]
241265876a83SMatthew G. Knepley       suffix: 3d_cryer_tconv
241330602db0SMatthew G. Knepley       args: -bd_dm_refine 1 -dm_refine_volume_limit_pre 0.00666667 \
241430602db0SMatthew G. Knepley             -ts_dt 0.023 -ts_max_time 0.092 -ts_max_steps 2 -ts_convergence_estimate -convest_num_refine 1 \
241530602db0SMatthew G. Knepley             -pc_type lu -pc_factor_shift_type nonzero
241665876a83SMatthew G. Knepley 
241765876a83SMatthew G. Knepley TEST*/
2418