1c4762a1bSJed Brown static char help[] = "Nonlinear elasticity problem in 3d with simplicial finite elements.\n\
2c4762a1bSJed Brown We solve a nonlinear elasticity problem, modelled as an incompressible Neo-Hookean solid, \n\
3c4762a1bSJed Brown with pressure loading in a rectangular domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4c4762a1bSJed Brown
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown Nonlinear elasticity problem, which we discretize using the finite
7c4762a1bSJed Brown element method on an unstructured mesh. This uses both Dirichlet boundary conditions (fixed faces)
8c4762a1bSJed Brown and nonlinear Neumann boundary conditions (pressure loading).
9c4762a1bSJed Brown The Lagrangian density (modulo boundary conditions) for this problem is given by
10c4762a1bSJed Brown \begin{equation}
11c4762a1bSJed Brown \frac{\mu}{2} (\mathrm{Tr}{C}-3) + J p + \frac{\kappa}{2} (J-1).
12c4762a1bSJed Brown \end{equation}
13c4762a1bSJed Brown
14c4762a1bSJed Brown Discretization:
15c4762a1bSJed Brown
16c4762a1bSJed Brown We use PetscFE to generate a tabulation of the finite element basis functions
17c4762a1bSJed Brown at quadrature points. We can currently generate an arbitrary order Lagrange
18c4762a1bSJed Brown element.
19c4762a1bSJed Brown
20c4762a1bSJed Brown Field Data:
21c4762a1bSJed Brown
22c4762a1bSJed Brown DMPLEX data is organized by point, and the closure operation just stacks up the
23c4762a1bSJed Brown data from each sieve point in the closure. Thus, for a P_2-P_1 Stokes element, we
24c4762a1bSJed Brown have
25c4762a1bSJed Brown
26c4762a1bSJed Brown cl{e} = {f e_0 e_1 e_2 v_0 v_1 v_2}
27c4762a1bSJed Brown x = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} p_{v_0} u_{v_1} v_{v_1} p_{v_1} u_{v_2} v_{v_2} p_{v_2}]
28c4762a1bSJed Brown
29c4762a1bSJed Brown The problem here is that we would like to loop over each field separately for
30c4762a1bSJed Brown integration. Therefore, the closure visitor in DMPlexVecGetClosure() reorders
31c4762a1bSJed Brown the data so that each field is contiguous
32c4762a1bSJed Brown
33c4762a1bSJed Brown x' = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} u_{v_1} v_{v_1} u_{v_2} v_{v_2} p_{v_0} p_{v_1} p_{v_2}]
34c4762a1bSJed Brown
35c4762a1bSJed Brown Likewise, DMPlexVecSetClosure() takes data partitioned by field, and correctly
36c4762a1bSJed Brown puts it into the Sieve ordering.
37c4762a1bSJed Brown
38c4762a1bSJed Brown */
39c4762a1bSJed Brown
40c4762a1bSJed Brown #include <petscdmplex.h>
41c4762a1bSJed Brown #include <petscsnes.h>
42c4762a1bSJed Brown #include <petscds.h>
43c4762a1bSJed Brown
449371c9d4SSatish Balay typedef enum {
459371c9d4SSatish Balay RUN_FULL,
469371c9d4SSatish Balay RUN_TEST
479371c9d4SSatish Balay } RunType;
48c4762a1bSJed Brown
49c4762a1bSJed Brown typedef struct {
50c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */
51c4762a1bSJed Brown PetscReal mu; /* The shear modulus */
52c4762a1bSJed Brown PetscReal p_wall; /* The wall pressure */
53c4762a1bSJed Brown } AppCtx;
54c4762a1bSJed Brown
55c4762a1bSJed Brown #if 0
569fbee547SJacob Faibussowitsch static inline void Det2D(PetscReal *detJ, const PetscReal J[])
57c4762a1bSJed Brown {
58c4762a1bSJed Brown *detJ = J[0]*J[3] - J[1]*J[2];
59c4762a1bSJed Brown }
60c4762a1bSJed Brown #endif
61c4762a1bSJed Brown
Det3D(PetscReal * detJ,const PetscScalar J[])62d71ae5a4SJacob Faibussowitsch static inline void Det3D(PetscReal *detJ, const PetscScalar J[])
63d71ae5a4SJacob Faibussowitsch {
649371c9d4SSatish Balay *detJ = PetscRealPart(J[0 * 3 + 0] * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]) + J[0 * 3 + 1] * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]) + J[0 * 3 + 2] * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]));
65c4762a1bSJed Brown }
66c4762a1bSJed Brown
67c4762a1bSJed Brown #if 0
689fbee547SJacob Faibussowitsch static inline void Cof2D(PetscReal C[], const PetscReal A[])
69c4762a1bSJed Brown {
70c4762a1bSJed Brown C[0] = A[3];
71c4762a1bSJed Brown C[1] = -A[2];
72c4762a1bSJed Brown C[2] = -A[1];
73c4762a1bSJed Brown C[3] = A[0];
74c4762a1bSJed Brown }
75c4762a1bSJed Brown #endif
76c4762a1bSJed Brown
Cof3D(PetscReal C[],const PetscScalar A[])77d71ae5a4SJacob Faibussowitsch static inline void Cof3D(PetscReal C[], const PetscScalar A[])
78d71ae5a4SJacob Faibussowitsch {
79c4762a1bSJed Brown C[0 * 3 + 0] = PetscRealPart(A[1 * 3 + 1] * A[2 * 3 + 2] - A[1 * 3 + 2] * A[2 * 3 + 1]);
80c4762a1bSJed Brown C[0 * 3 + 1] = PetscRealPart(A[1 * 3 + 2] * A[2 * 3 + 0] - A[1 * 3 + 0] * A[2 * 3 + 2]);
81c4762a1bSJed Brown C[0 * 3 + 2] = PetscRealPart(A[1 * 3 + 0] * A[2 * 3 + 1] - A[1 * 3 + 1] * A[2 * 3 + 0]);
82c4762a1bSJed Brown C[1 * 3 + 0] = PetscRealPart(A[0 * 3 + 2] * A[2 * 3 + 1] - A[0 * 3 + 1] * A[2 * 3 + 2]);
83c4762a1bSJed Brown C[1 * 3 + 1] = PetscRealPart(A[0 * 3 + 0] * A[2 * 3 + 2] - A[0 * 3 + 2] * A[2 * 3 + 0]);
84c4762a1bSJed Brown C[1 * 3 + 2] = PetscRealPart(A[0 * 3 + 1] * A[2 * 3 + 0] - A[0 * 3 + 0] * A[2 * 3 + 1]);
85c4762a1bSJed Brown C[2 * 3 + 0] = PetscRealPart(A[0 * 3 + 1] * A[1 * 3 + 2] - A[0 * 3 + 2] * A[1 * 3 + 1]);
86c4762a1bSJed Brown C[2 * 3 + 1] = PetscRealPart(A[0 * 3 + 2] * A[1 * 3 + 0] - A[0 * 3 + 0] * A[1 * 3 + 2]);
87c4762a1bSJed Brown C[2 * 3 + 2] = PetscRealPart(A[0 * 3 + 0] * A[1 * 3 + 1] - A[0 * 3 + 1] * A[1 * 3 + 0]);
88c4762a1bSJed Brown }
89c4762a1bSJed Brown
zero_scalar(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)90*2a8381b2SBarry Smith PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
91d71ae5a4SJacob Faibussowitsch {
92c4762a1bSJed Brown u[0] = 0.0;
933ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
94c4762a1bSJed Brown }
95c4762a1bSJed Brown
zero_vector(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)96*2a8381b2SBarry Smith PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
97d71ae5a4SJacob Faibussowitsch {
98c4762a1bSJed Brown const PetscInt Ncomp = dim;
99c4762a1bSJed Brown
100c4762a1bSJed Brown PetscInt comp;
101c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0;
1023ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
103c4762a1bSJed Brown }
104c4762a1bSJed Brown
coordinates(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)105*2a8381b2SBarry Smith PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
106d71ae5a4SJacob Faibussowitsch {
107c4762a1bSJed Brown const PetscInt Ncomp = dim;
108c4762a1bSJed Brown
109c4762a1bSJed Brown PetscInt comp;
110c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp];
1113ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
112c4762a1bSJed Brown }
113c4762a1bSJed Brown
elasticityMaterial(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)114*2a8381b2SBarry Smith PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
115d71ae5a4SJacob Faibussowitsch {
116c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
117c4762a1bSJed Brown u[0] = user->mu;
1183ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
119c4762a1bSJed Brown }
120c4762a1bSJed Brown
wallPressure(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)121*2a8381b2SBarry Smith PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
122d71ae5a4SJacob Faibussowitsch {
123c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx;
124c4762a1bSJed Brown u[0] = user->p_wall;
1253ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
126c4762a1bSJed Brown }
127c4762a1bSJed Brown
f1_u_3d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f1[])128d71ae5a4SJacob Faibussowitsch void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
129d71ae5a4SJacob Faibussowitsch {
130c4762a1bSJed Brown const PetscInt Ncomp = dim;
131c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
132c4762a1bSJed Brown PetscReal cofu_x[9 /* Ncomp*dim */], detu_x, p = PetscRealPart(u[Ncomp]);
133c4762a1bSJed Brown PetscInt comp, d;
134c4762a1bSJed Brown
135c4762a1bSJed Brown Cof3D(cofu_x, u_x);
136c4762a1bSJed Brown Det3D(&detu_x, u_x);
137c4762a1bSJed Brown p += kappa * (detu_x - 1.0);
138c4762a1bSJed Brown
139c4762a1bSJed Brown /* f1 is the first Piola-Kirchhoff tensor */
140c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) {
141ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f1[comp * dim + d] = mu * u_x[comp * dim + d] + p * cofu_x[comp * dim + d];
142c4762a1bSJed Brown }
143c4762a1bSJed Brown }
144c4762a1bSJed Brown
g3_uu_3d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g3[])145d71ae5a4SJacob Faibussowitsch void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
146d71ae5a4SJacob Faibussowitsch {
147c4762a1bSJed Brown const PetscInt Ncomp = dim;
148c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
149c4762a1bSJed Brown PetscReal cofu_x[9 /* Ncomp*dim */], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]);
150c4762a1bSJed Brown PetscInt compI, compJ, d1, d2;
151c4762a1bSJed Brown
152c4762a1bSJed Brown Cof3D(cofu_x, u_x);
153c4762a1bSJed Brown Det3D(&detu_x, u_x);
154c4762a1bSJed Brown p += kappa * (detu_x - 1.0);
155c4762a1bSJed Brown pp = p / detu_x + kappa;
156c4762a1bSJed Brown pm = p / detu_x;
157c4762a1bSJed Brown
158c4762a1bSJed Brown /* g3 is the first elasticity tensor, i.e. A_i^I_j^J = S^{IJ} g_{ij} + C^{KILJ} F^k_K F^l_L g_{kj} g_{li} */
159c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) {
160c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) {
161c4762a1bSJed Brown const PetscReal G = (compI == compJ) ? 1.0 : 0.0;
162c4762a1bSJed Brown
163c4762a1bSJed Brown for (d1 = 0; d1 < dim; ++d1) {
164c4762a1bSJed Brown for (d2 = 0; d2 < dim; ++d2) {
165c4762a1bSJed Brown const PetscReal g = (d1 == d2) ? 1.0 : 0.0;
166c4762a1bSJed Brown
167c4762a1bSJed Brown g3[((compI * Ncomp + compJ) * dim + d1) * dim + d2] = g * G * mu + pp * cofu_x[compI * dim + d1] * cofu_x[compJ * dim + d2] - pm * cofu_x[compI * dim + d2] * cofu_x[compJ * dim + d1];
168c4762a1bSJed Brown }
169c4762a1bSJed Brown }
170c4762a1bSJed Brown }
171c4762a1bSJed Brown }
172c4762a1bSJed Brown }
173c4762a1bSJed Brown
f0_bd_u_3d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])174d71ae5a4SJacob Faibussowitsch void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
175d71ae5a4SJacob Faibussowitsch {
176c4762a1bSJed Brown const PetscInt Ncomp = dim;
177c4762a1bSJed Brown const PetscScalar p = a[aOff[1]];
178c4762a1bSJed Brown PetscReal cofu_x[9 /* Ncomp*dim */];
179c4762a1bSJed Brown PetscInt comp, d;
180c4762a1bSJed Brown
181c4762a1bSJed Brown Cof3D(cofu_x, u_x);
182c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) {
183c4762a1bSJed Brown for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp * dim + d] * n[d];
184c4762a1bSJed Brown f0[comp] *= p;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown }
187c4762a1bSJed Brown
g1_bd_uu_3d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],const PetscReal n[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g1[])188d71ae5a4SJacob Faibussowitsch void g1_bd_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
189d71ae5a4SJacob Faibussowitsch {
190c4762a1bSJed Brown const PetscInt Ncomp = dim;
191c4762a1bSJed Brown PetscScalar p = a[aOff[1]];
192c4762a1bSJed Brown PetscReal cofu_x[9 /* Ncomp*dim */], m[3 /* Ncomp */], detu_x;
193c4762a1bSJed Brown PetscInt comp, compI, compJ, d;
194c4762a1bSJed Brown
195c4762a1bSJed Brown Cof3D(cofu_x, u_x);
196c4762a1bSJed Brown Det3D(&detu_x, u_x);
197c4762a1bSJed Brown p /= detu_x;
198c4762a1bSJed Brown
1999371c9d4SSatish Balay for (comp = 0; comp < Ncomp; ++comp)
2009371c9d4SSatish Balay for (d = 0, m[comp] = 0.0; d < dim; ++d) m[comp] += cofu_x[comp * dim + d] * n[d];
201c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) {
202c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) {
203ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g1[(compI * Ncomp + compJ) * dim + d] = p * (m[compI] * cofu_x[compJ * dim + d] - cofu_x[compI * dim + d] * m[compJ]);
204c4762a1bSJed Brown }
205c4762a1bSJed Brown }
206c4762a1bSJed Brown }
207c4762a1bSJed Brown
f0_p_3d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])208d71ae5a4SJacob Faibussowitsch void f0_p_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
209d71ae5a4SJacob Faibussowitsch {
210c4762a1bSJed Brown PetscReal detu_x;
211c4762a1bSJed Brown Det3D(&detu_x, u_x);
212c4762a1bSJed Brown f0[0] = detu_x - 1.0;
213c4762a1bSJed Brown }
214c4762a1bSJed Brown
g1_pu_3d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g1[])215d71ae5a4SJacob Faibussowitsch void g1_pu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
216d71ae5a4SJacob Faibussowitsch {
217c4762a1bSJed Brown PetscReal cofu_x[9 /* Ncomp*dim */];
218c4762a1bSJed Brown PetscInt compI, d;
219c4762a1bSJed Brown
220c4762a1bSJed Brown Cof3D(cofu_x, u_x);
221c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI)
222c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[compI * dim + d] = cofu_x[compI * dim + d];
223c4762a1bSJed Brown }
224c4762a1bSJed Brown
g2_up_3d(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g2[])225d71ae5a4SJacob Faibussowitsch void g2_up_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
226d71ae5a4SJacob Faibussowitsch {
227c4762a1bSJed Brown PetscReal cofu_x[9 /* Ncomp*dim */];
228c4762a1bSJed Brown PetscInt compI, d;
229c4762a1bSJed Brown
230c4762a1bSJed Brown Cof3D(cofu_x, u_x);
231c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI)
232c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[compI * dim + d] = cofu_x[compI * dim + d];
233c4762a1bSJed Brown }
234c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)235d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
236d71ae5a4SJacob Faibussowitsch {
237c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"};
238c4762a1bSJed Brown PetscInt run;
239c4762a1bSJed Brown
240c4762a1bSJed Brown PetscFunctionBeginUser;
241c4762a1bSJed Brown options->runType = RUN_FULL;
242c4762a1bSJed Brown options->mu = 1.0;
243c4762a1bSJed Brown options->p_wall = 0.4;
244d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");
245c4762a1bSJed Brown run = options->runType;
2469566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL));
247c4762a1bSJed Brown options->runType = (RunType)run;
2489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL));
2499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL));
250d0609cedSBarry Smith PetscOptionsEnd();
2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)254d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
255d71ae5a4SJacob Faibussowitsch {
256c4762a1bSJed Brown PetscFunctionBeginUser;
2579566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
2589566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
2599566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
260c4762a1bSJed Brown /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */
2619566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
262c4762a1bSJed Brown {
263c4762a1bSJed Brown DM cdm;
264c4762a1bSJed Brown DMLabel label;
265c4762a1bSJed Brown IS is;
26630602db0SMatthew G. Knepley PetscInt d, dim, b, f, Nf;
267c4762a1bSJed Brown const PetscInt *faces;
268c4762a1bSJed Brown PetscInt csize;
269c4762a1bSJed Brown PetscScalar *coords = NULL;
270c4762a1bSJed Brown PetscSection cs;
271c4762a1bSJed Brown Vec coordinates;
272c4762a1bSJed Brown
2739566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim));
2749566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "boundary"));
2759566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "boundary", &label));
2769566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(*dm, 1, label));
2779566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(*dm, "boundary", 1, &is));
278c4762a1bSJed Brown if (is) {
279c4762a1bSJed Brown PetscReal faceCoord;
280c4762a1bSJed Brown PetscInt v;
281c4762a1bSJed Brown
2829566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &Nf));
2839566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &faces));
284c4762a1bSJed Brown
2859566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(*dm, &coordinates));
2869566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm));
2879566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cdm, &cs));
288c4762a1bSJed Brown
289c4762a1bSJed Brown /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
290c4762a1bSJed Brown for (f = 0; f < Nf; ++f) {
2919566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
292c4762a1bSJed Brown /* Calculate mean coordinate vector */
293c4762a1bSJed Brown for (d = 0; d < dim; ++d) {
294c4762a1bSJed Brown const PetscInt Nv = csize / dim;
295c4762a1bSJed Brown faceCoord = 0.0;
296c4762a1bSJed Brown for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v * dim + d]);
297c4762a1bSJed Brown faceCoord /= Nv;
298c4762a1bSJed Brown for (b = 0; b < 2; ++b) {
29948a46eb9SPierre Jolivet if (PetscAbs(faceCoord - b * 1.0) < PETSC_SMALL) PetscCall(DMSetLabelValue(*dm, "Faces", faces[f], d * 2 + b + 1));
300c4762a1bSJed Brown }
301c4762a1bSJed Brown }
3029566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
303c4762a1bSJed Brown }
3049566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &faces));
305c4762a1bSJed Brown }
3069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is));
307c4762a1bSJed Brown }
3089566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
310c4762a1bSJed Brown }
311c4762a1bSJed Brown
SetupProblem(DM dm,PetscInt dim,AppCtx * user)312d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user)
313d71ae5a4SJacob Faibussowitsch {
31445480ffeSMatthew G. Knepley PetscDS ds;
31545480ffeSMatthew G. Knepley PetscWeakForm wf;
31645480ffeSMatthew G. Knepley DMLabel label;
31745480ffeSMatthew G. Knepley PetscInt bd;
318c4762a1bSJed Brown
319c4762a1bSJed Brown PetscFunctionBeginUser;
3209566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
3219566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u_3d));
3229566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p_3d, NULL));
3239566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu_3d));
3249566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up_3d, NULL));
3259566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL, NULL));
326c4762a1bSJed Brown
3279566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label));
3289566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd));
3299566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3309566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_3d, 0, NULL));
3319566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL));
332c4762a1bSJed Brown
33357d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (PetscVoidFn *)coordinates, NULL, user, NULL));
3343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
335c4762a1bSJed Brown }
336c4762a1bSJed Brown
SetupMaterial(DM dm,DM dmAux,AppCtx * user)337d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
338d71ae5a4SJacob Faibussowitsch {
339*2a8381b2SBarry Smith PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], PetscCtx ctx) = {elasticityMaterial, wallPressure};
340c4762a1bSJed Brown Vec A;
341*2a8381b2SBarry Smith PetscCtx ctxs[2];
342c4762a1bSJed Brown
343c4762a1bSJed Brown PetscFunctionBegin;
3449371c9d4SSatish Balay ctxs[0] = user;
3459371c9d4SSatish Balay ctxs[1] = user;
3469566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmAux, &A));
3479566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A));
3489566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, A));
3499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&A));
3503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown
SetupNearNullSpace(DM dm,AppCtx * user)353d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user)
354d71ae5a4SJacob Faibussowitsch {
355c4762a1bSJed Brown /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
356c4762a1bSJed Brown DM subdm;
357c4762a1bSJed Brown MatNullSpace nearNullSpace;
358c4762a1bSJed Brown PetscInt fields = 0;
359c4762a1bSJed Brown PetscObject deformation;
360c4762a1bSJed Brown
361c4762a1bSJed Brown PetscFunctionBeginUser;
3629566063dSJacob Faibussowitsch PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm));
3639566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace));
3649566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation));
3659566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject)nearNullSpace));
3669566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm));
3679566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nearNullSpace));
3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
369c4762a1bSJed Brown }
370c4762a1bSJed Brown
SetupAuxDM(DM dm,PetscInt NfAux,PetscFE feAux[],AppCtx * user)371d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user)
372d71ae5a4SJacob Faibussowitsch {
373c4762a1bSJed Brown DM dmAux, coordDM;
374c4762a1bSJed Brown PetscInt f;
375c4762a1bSJed Brown
376c4762a1bSJed Brown PetscFunctionBegin;
377c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
3789566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM));
3793ba16761SJacob Faibussowitsch if (!feAux) PetscFunctionReturn(PETSC_SUCCESS);
3809566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAux));
3819566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDM(dmAux, coordDM));
3829566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(DMSetField(dmAux, f, NULL, (PetscObject)feAux[f]));
3839566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmAux));
3849566063dSJacob Faibussowitsch PetscCall(SetupMaterial(dm, dmAux, user));
3859566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAux));
3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
387c4762a1bSJed Brown }
388c4762a1bSJed Brown
SetupDiscretization(DM dm,AppCtx * user)389d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
390d71ae5a4SJacob Faibussowitsch {
391c4762a1bSJed Brown DM cdm = dm;
392c4762a1bSJed Brown PetscFE fe[2], feAux[2];
39330602db0SMatthew G. Knepley PetscBool simplex;
39430602db0SMatthew G. Knepley PetscInt dim;
395c4762a1bSJed Brown MPI_Comm comm;
396c4762a1bSJed Brown
397c4762a1bSJed Brown PetscFunctionBeginUser;
3989566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3999566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
4009566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex));
401c4762a1bSJed Brown /* Create finite element */
4029566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]));
4039566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "deformation"));
4049566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
4059566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
406c4762a1bSJed Brown
4079566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
408c4762a1bSJed Brown
4099566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]));
4109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feAux[0], "elasticityMaterial"));
4119566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[0]));
412c4762a1bSJed Brown /* It is not yet possible to define a field on a submesh (e.g. a boundary), so we will use a normal finite element */
4139566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]));
4149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feAux[1], "wall_pressure"));
4159566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[1]));
416c4762a1bSJed Brown
417c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */
4189566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
4199566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
4209566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
4219566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, dim, user));
422c4762a1bSJed Brown while (cdm) {
4239566063dSJacob Faibussowitsch PetscCall(SetupAuxDM(cdm, 2, feAux, user));
4249566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
4259566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
426c4762a1bSJed Brown }
4279566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0]));
4289566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1]));
4299566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feAux[0]));
4309566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feAux[1]));
4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
432c4762a1bSJed Brown }
433c4762a1bSJed Brown
main(int argc,char ** argv)434d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
435d71ae5a4SJacob Faibussowitsch {
436c4762a1bSJed Brown SNES snes; /* nonlinear solver */
437c4762a1bSJed Brown DM dm; /* problem definition */
438c4762a1bSJed Brown Vec u, r; /* solution, residual vectors */
439c4762a1bSJed Brown Mat A, J; /* Jacobian matrix */
440c4762a1bSJed Brown AppCtx user; /* user-defined work context */
441c4762a1bSJed Brown PetscInt its; /* iterations for convergence */
442c4762a1bSJed Brown
443327415f7SBarry Smith PetscFunctionBeginUser;
4449566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4459566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
4469566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
4479566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
4489566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm));
4499566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user));
450c4762a1bSJed Brown
4519566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user));
4529566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL));
4539566063dSJacob Faibussowitsch PetscCall(SetupNearNullSpace(dm, &user));
454c4762a1bSJed Brown
4559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u));
4564f25160eSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)u, "u"));
4579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r));
458c4762a1bSJed Brown
4599566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm, MATAIJ));
4609566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J));
461c4762a1bSJed Brown A = J;
462c4762a1bSJed Brown
4636493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
4649566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, A, J, NULL, NULL));
465c4762a1bSJed Brown
4669566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
467c4762a1bSJed Brown
468c4762a1bSJed Brown {
469*2a8381b2SBarry Smith PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
4709371c9d4SSatish Balay initialGuess[0] = coordinates;
4719371c9d4SSatish Balay initialGuess[1] = zero_scalar;
4729566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
473c4762a1bSJed Brown }
474c4762a1bSJed Brown
475c4762a1bSJed Brown if (user.runType == RUN_FULL) {
4769566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u));
4779566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its));
47863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
479c4762a1bSJed Brown } else {
480c4762a1bSJed Brown PetscReal res = 0.0;
481c4762a1bSJed Brown
482c4762a1bSJed Brown /* Check initial guess */
4839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n"));
4849566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
485c4762a1bSJed Brown /* Check residual */
4869566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, u, r));
4879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n"));
48893ec0da9SPierre Jolivet PetscCall(VecFilter(r, 1.0e-10));
4899566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
4909566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res));
4919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res));
492c4762a1bSJed Brown /* Check Jacobian */
493c4762a1bSJed Brown {
494c4762a1bSJed Brown Vec b;
495c4762a1bSJed Brown
4969566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, u, A, A));
4979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b));
4989566063dSJacob Faibussowitsch PetscCall(VecSet(r, 0.0));
4999566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, r, b));
5009566063dSJacob Faibussowitsch PetscCall(MatMult(A, u, r));
5019566063dSJacob Faibussowitsch PetscCall(VecAXPY(r, 1.0, b));
5029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
5039566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n"));
50493ec0da9SPierre Jolivet PetscCall(VecFilter(r, 1.0e-10));
5059566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
5069566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res));
5079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res));
508c4762a1bSJed Brown }
509c4762a1bSJed Brown }
5109566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
511c4762a1bSJed Brown
5129566063dSJacob Faibussowitsch if (A != J) PetscCall(MatDestroy(&A));
5139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
5149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
5159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
5169566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
5179566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
5189566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
519b122ec5aSJacob Faibussowitsch return 0;
520c4762a1bSJed Brown }
521c4762a1bSJed Brown
522c4762a1bSJed Brown /*TEST
523c4762a1bSJed Brown
524c4762a1bSJed Brown build:
525c4762a1bSJed Brown requires: !complex
526c4762a1bSJed Brown
52730602db0SMatthew G. Knepley testset:
52830602db0SMatthew G. Knepley requires: ctetgen
52930602db0SMatthew G. Knepley args: -run_type full -dm_plex_dim 3 \
53030602db0SMatthew G. Knepley -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 \
53130602db0SMatthew G. Knepley -snes_rtol 1e-05 -snes_monitor_short -snes_converged_reason \
53230602db0SMatthew G. Knepley -ksp_type fgmres -ksp_rtol 1e-10 -ksp_monitor_short -ksp_converged_reason \
53330602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper \
53430602db0SMatthew G. Knepley -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu \
53530602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
53630602db0SMatthew G. Knepley
537c4762a1bSJed Brown test:
538c4762a1bSJed Brown suffix: 0
53930602db0SMatthew G. Knepley requires: !single
54030602db0SMatthew G. Knepley args: -dm_refine 2 \
54130602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
542c4762a1bSJed Brown
543c4762a1bSJed Brown test:
544c4762a1bSJed Brown suffix: 1
54530602db0SMatthew G. Knepley requires: superlu_dist
546c4762a1bSJed Brown nsize: 2
547e600fa54SMatthew G. Knepley args: -dm_refine 0 -petscpartitioner_type simple \
54830602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
549c4762a1bSJed Brown timeoutfactor: 2
550c4762a1bSJed Brown
551c4762a1bSJed Brown test:
552c4762a1bSJed Brown suffix: 4
55330602db0SMatthew G. Knepley requires: superlu_dist
554c4762a1bSJed Brown nsize: 2
555e600fa54SMatthew G. Knepley args: -dm_refine 0 -petscpartitioner_type simple \
55630602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
557c4762a1bSJed Brown output_file: output/ex77_1.out
558c4762a1bSJed Brown
559c4762a1bSJed Brown test:
56030602db0SMatthew G. Knepley suffix: 2
56130602db0SMatthew G. Knepley requires: !single
56230602db0SMatthew G. Knepley args: -dm_refine 2 \
56330602db0SMatthew G. Knepley -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
564c4762a1bSJed Brown
565c4762a1bSJed Brown test:
566c4762a1bSJed Brown suffix: 2_par
56730602db0SMatthew G. Knepley requires: superlu_dist !single
568c4762a1bSJed Brown nsize: 4
569e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple \
57030602db0SMatthew G. Knepley -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
571c4762a1bSJed Brown output_file: output/ex77_2.out
572c4762a1bSJed Brown
573c4762a1bSJed Brown TEST*/
574