xref: /petsc/src/snes/tutorials/ex77.c (revision e600fa544e2bb197ca2af9b6e65ea465976dec56)
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 
44c4762a1bSJed Brown typedef enum {RUN_FULL, RUN_TEST} RunType;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown typedef struct {
47c4762a1bSJed Brown   RunType   runType; /* Whether to run tests, or solve the full problem */
48c4762a1bSJed Brown   PetscReal mu;      /* The shear modulus */
49c4762a1bSJed Brown   PetscReal p_wall;  /* The wall pressure */
50c4762a1bSJed Brown } AppCtx;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown #if 0
539fbee547SJacob Faibussowitsch static inline void Det2D(PetscReal *detJ, const PetscReal J[])
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   *detJ = J[0]*J[3] - J[1]*J[2];
56c4762a1bSJed Brown }
57c4762a1bSJed Brown #endif
58c4762a1bSJed Brown 
599fbee547SJacob Faibussowitsch static inline void Det3D(PetscReal *detJ, const PetscScalar J[])
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   *detJ = PetscRealPart(J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
62c4762a1bSJed Brown                         J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
63c4762a1bSJed Brown                         J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown #if 0
679fbee547SJacob Faibussowitsch static inline void Cof2D(PetscReal C[], const PetscReal A[])
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   C[0] =  A[3];
70c4762a1bSJed Brown   C[1] = -A[2];
71c4762a1bSJed Brown   C[2] = -A[1];
72c4762a1bSJed Brown   C[3] =  A[0];
73c4762a1bSJed Brown }
74c4762a1bSJed Brown #endif
75c4762a1bSJed Brown 
769fbee547SJacob Faibussowitsch static inline void Cof3D(PetscReal C[], const PetscScalar A[])
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   C[0*3+0] = PetscRealPart(A[1*3+1]*A[2*3+2] - A[1*3+2]*A[2*3+1]);
79c4762a1bSJed Brown   C[0*3+1] = PetscRealPart(A[1*3+2]*A[2*3+0] - A[1*3+0]*A[2*3+2]);
80c4762a1bSJed Brown   C[0*3+2] = PetscRealPart(A[1*3+0]*A[2*3+1] - A[1*3+1]*A[2*3+0]);
81c4762a1bSJed Brown   C[1*3+0] = PetscRealPart(A[0*3+2]*A[2*3+1] - A[0*3+1]*A[2*3+2]);
82c4762a1bSJed Brown   C[1*3+1] = PetscRealPart(A[0*3+0]*A[2*3+2] - A[0*3+2]*A[2*3+0]);
83c4762a1bSJed Brown   C[1*3+2] = PetscRealPart(A[0*3+1]*A[2*3+0] - A[0*3+0]*A[2*3+1]);
84c4762a1bSJed Brown   C[2*3+0] = PetscRealPart(A[0*3+1]*A[1*3+2] - A[0*3+2]*A[1*3+1]);
85c4762a1bSJed Brown   C[2*3+1] = PetscRealPart(A[0*3+2]*A[1*3+0] - A[0*3+0]*A[1*3+2]);
86c4762a1bSJed Brown   C[2*3+2] = PetscRealPart(A[0*3+0]*A[1*3+1] - A[0*3+1]*A[1*3+0]);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89c4762a1bSJed Brown PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
90c4762a1bSJed Brown {
91c4762a1bSJed Brown   u[0] = 0.0;
92c4762a1bSJed Brown   return 0;
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
96c4762a1bSJed Brown {
97c4762a1bSJed Brown   const PetscInt Ncomp = dim;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   PetscInt       comp;
100c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0;
101c4762a1bSJed Brown   return 0;
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   const PetscInt Ncomp = dim;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   PetscInt       comp;
109c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp];
110c4762a1bSJed Brown   return 0;
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
114c4762a1bSJed Brown {
115c4762a1bSJed Brown   AppCtx *user = (AppCtx *) ctx;
116c4762a1bSJed Brown   u[0] = user->mu;
117c4762a1bSJed Brown   return 0;
118c4762a1bSJed Brown }
119c4762a1bSJed Brown 
120c4762a1bSJed Brown PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
121c4762a1bSJed Brown {
122c4762a1bSJed Brown   AppCtx *user = (AppCtx *) ctx;
123c4762a1bSJed Brown   u[0] = user->p_wall;
124c4762a1bSJed Brown   return 0;
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
128c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
129c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
130c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
131c4762a1bSJed Brown {
132c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
133c4762a1bSJed Brown   const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
134c4762a1bSJed Brown   PetscReal       cofu_x[9/*Ncomp*dim*/], detu_x, p = PetscRealPart(u[Ncomp]);
135c4762a1bSJed Brown   PetscInt        comp, d;
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   Cof3D(cofu_x, u_x);
138c4762a1bSJed Brown   Det3D(&detu_x, u_x);
139c4762a1bSJed Brown   p += kappa * (detu_x - 1.0);
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* f1 is the first Piola-Kirchhoff tensor */
142c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
143c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
144c4762a1bSJed Brown       f1[comp*dim+d] = mu * u_x[comp*dim+d] + p * cofu_x[comp*dim+d];
145c4762a1bSJed Brown     }
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
150c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
151c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
152c4762a1bSJed Brown           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
155c4762a1bSJed Brown   const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
156c4762a1bSJed Brown   PetscReal       cofu_x[9/*Ncomp*dim*/], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]);
157c4762a1bSJed Brown   PetscInt        compI, compJ, d1, d2;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   Cof3D(cofu_x, u_x);
160c4762a1bSJed Brown   Det3D(&detu_x, u_x);
161c4762a1bSJed Brown   p += kappa * (detu_x - 1.0);
162c4762a1bSJed Brown   pp = p/detu_x + kappa;
163c4762a1bSJed Brown   pm = p/detu_x;
164c4762a1bSJed Brown 
165c4762a1bSJed 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} */
166c4762a1bSJed Brown   for (compI = 0; compI < Ncomp; ++compI) {
167c4762a1bSJed Brown     for (compJ = 0; compJ < Ncomp; ++compJ) {
168c4762a1bSJed Brown       const PetscReal G = (compI == compJ) ? 1.0 : 0.0;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown       for (d1 = 0; d1 < dim; ++d1) {
171c4762a1bSJed Brown         for (d2 = 0; d2 < dim; ++d2) {
172c4762a1bSJed Brown           const PetscReal g = (d1 == d2) ? 1.0 : 0.0;
173c4762a1bSJed Brown 
174c4762a1bSJed 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];
175c4762a1bSJed Brown         }
176c4762a1bSJed Brown       }
177c4762a1bSJed Brown     }
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
182c4762a1bSJed Brown     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
183c4762a1bSJed Brown     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
184c4762a1bSJed Brown     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
185c4762a1bSJed Brown {
186c4762a1bSJed Brown   const PetscInt    Ncomp = dim;
187c4762a1bSJed Brown   const PetscScalar p = a[aOff[1]];
188c4762a1bSJed Brown   PetscReal         cofu_x[9/*Ncomp*dim*/];
189c4762a1bSJed Brown   PetscInt          comp, d;
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   Cof3D(cofu_x, u_x);
192c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
193c4762a1bSJed Brown     for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp*dim+d] * n[d];
194c4762a1bSJed Brown     f0[comp] *= p;
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown void g1_bd_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
199c4762a1bSJed Brown     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
200c4762a1bSJed Brown     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
201c4762a1bSJed Brown     PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   const PetscInt Ncomp = dim;
204c4762a1bSJed Brown   PetscScalar    p = a[aOff[1]];
205c4762a1bSJed Brown   PetscReal      cofu_x[9/*Ncomp*dim*/], m[3/*Ncomp*/], detu_x;
206c4762a1bSJed Brown   PetscInt       comp, compI, compJ, d;
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   Cof3D(cofu_x, u_x);
209c4762a1bSJed Brown   Det3D(&detu_x, u_x);
210c4762a1bSJed Brown   p /= detu_x;
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) for (d = 0, m[comp] = 0.0; d < dim; ++d) m[comp] += cofu_x[comp*dim+d] * n[d];
213c4762a1bSJed Brown   for (compI = 0; compI < Ncomp; ++compI) {
214c4762a1bSJed Brown     for (compJ = 0; compJ < Ncomp; ++compJ) {
215c4762a1bSJed Brown       for (d = 0; d < dim; ++d) {
216c4762a1bSJed Brown         g1[(compI*Ncomp+compJ)*dim+d] = p * (m[compI] * cofu_x[compJ*dim+d] - cofu_x[compI*dim+d] * m[compJ]);
217c4762a1bSJed Brown       }
218c4762a1bSJed Brown     }
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown void f0_p_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
223c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
224c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
225c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
226c4762a1bSJed Brown {
227c4762a1bSJed Brown   PetscReal detu_x;
228c4762a1bSJed Brown   Det3D(&detu_x, u_x);
229c4762a1bSJed Brown   f0[0] = detu_x - 1.0;
230c4762a1bSJed Brown }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown void g1_pu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
233c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
234c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
235c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
236c4762a1bSJed Brown {
237c4762a1bSJed Brown   PetscReal cofu_x[9/*Ncomp*dim*/];
238c4762a1bSJed Brown   PetscInt  compI, d;
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   Cof3D(cofu_x, u_x);
241c4762a1bSJed Brown   for (compI = 0; compI < dim; ++compI)
242c4762a1bSJed Brown     for (d = 0; d < dim; ++d) g1[compI*dim+d] = cofu_x[compI*dim+d];
243c4762a1bSJed Brown }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown void g2_up_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
246c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
247c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
248c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
249c4762a1bSJed Brown {
250c4762a1bSJed Brown   PetscReal cofu_x[9/*Ncomp*dim*/];
251c4762a1bSJed Brown   PetscInt  compI, d;
252c4762a1bSJed Brown 
253c4762a1bSJed Brown   Cof3D(cofu_x, u_x);
254c4762a1bSJed Brown   for (compI = 0; compI < dim; ++compI)
255c4762a1bSJed Brown     for (d = 0; d < dim; ++d) g2[compI*dim+d] = cofu_x[compI*dim+d];
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
259c4762a1bSJed Brown {
260c4762a1bSJed Brown   const char    *runTypes[2] = {"full", "test"};
261c4762a1bSJed Brown   PetscInt       run;
262c4762a1bSJed Brown   PetscErrorCode ierr;
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   PetscFunctionBeginUser;
265c4762a1bSJed Brown   options->runType      = RUN_FULL;
266c4762a1bSJed Brown   options->mu           = 1.0;
267c4762a1bSJed Brown   options->p_wall       = 0.4;
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");CHKERRQ(ierr);
270c4762a1bSJed Brown   run  = options->runType;
271c4762a1bSJed Brown   ierr = PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr);
272c4762a1bSJed Brown   options->runType = (RunType) run;
273c4762a1bSJed Brown   ierr = PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL);CHKERRQ(ierr);
274c4762a1bSJed Brown   ierr = PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL);CHKERRQ(ierr);
2751e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
276c4762a1bSJed Brown   PetscFunctionReturn(0);
277c4762a1bSJed Brown }
278c4762a1bSJed Brown 
279c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
280c4762a1bSJed Brown {
281c4762a1bSJed Brown   PetscErrorCode ierr;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   PetscFunctionBeginUser;
28430602db0SMatthew G. Knepley   /* TODO The P1 coordinate space gives wrong results when compared to the affine version. Track this down */
28530602db0SMatthew G. Knepley   if (0) {
28630602db0SMatthew G. Knepley     ierr = DMPlexCreateBoxMesh(comm, 3, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
28730602db0SMatthew G. Knepley   } else {
28830602db0SMatthew G. Knepley     ierr = DMCreate(comm, dm);CHKERRQ(ierr);
28930602db0SMatthew G. Knepley     ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
29030602db0SMatthew G. Knepley   }
29130602db0SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
292c4762a1bSJed Brown   /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */
29330602db0SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr);
294c4762a1bSJed Brown   {
295c4762a1bSJed Brown     DM              cdm;
296c4762a1bSJed Brown     DMLabel         label;
297c4762a1bSJed Brown     IS              is;
29830602db0SMatthew G. Knepley     PetscInt        d, dim, b, f, Nf;
299c4762a1bSJed Brown     const PetscInt *faces;
300c4762a1bSJed Brown     PetscInt        csize;
301c4762a1bSJed Brown     PetscScalar    *coords = NULL;
302c4762a1bSJed Brown     PetscSection    cs;
303c4762a1bSJed Brown     Vec             coordinates ;
304c4762a1bSJed Brown 
30530602db0SMatthew G. Knepley     ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
306c4762a1bSJed Brown     ierr = DMCreateLabel(*dm, "boundary");CHKERRQ(ierr);
307c4762a1bSJed Brown     ierr = DMGetLabel(*dm, "boundary", &label);CHKERRQ(ierr);
308c4762a1bSJed Brown     ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr);
309c4762a1bSJed Brown     ierr = DMGetStratumIS(*dm, "boundary", 1,  &is);CHKERRQ(ierr);
310c4762a1bSJed Brown     if (is) {
311c4762a1bSJed Brown       PetscReal faceCoord;
312c4762a1bSJed Brown       PetscInt  v;
313c4762a1bSJed Brown 
314c4762a1bSJed Brown       ierr = ISGetLocalSize(is, &Nf);CHKERRQ(ierr);
315c4762a1bSJed Brown       ierr = ISGetIndices(is, &faces);CHKERRQ(ierr);
316c4762a1bSJed Brown 
317c4762a1bSJed Brown       ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr);
318c4762a1bSJed Brown       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
319c4762a1bSJed Brown       ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr);
320c4762a1bSJed Brown 
321c4762a1bSJed Brown       /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
322c4762a1bSJed Brown       for (f = 0; f < Nf; ++f) {
323c4762a1bSJed Brown         ierr = DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr);
324c4762a1bSJed Brown         /* Calculate mean coordinate vector */
325c4762a1bSJed Brown         for (d = 0; d < dim; ++d) {
326c4762a1bSJed Brown           const PetscInt Nv = csize/dim;
327c4762a1bSJed Brown           faceCoord = 0.0;
328c4762a1bSJed Brown           for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]);
329c4762a1bSJed Brown           faceCoord /= Nv;
330c4762a1bSJed Brown           for (b = 0; b < 2; ++b) {
331c4762a1bSJed Brown             if (PetscAbs(faceCoord - b*1.0) < PETSC_SMALL) {
332c4762a1bSJed Brown               ierr = DMSetLabelValue(*dm, "Faces", faces[f], d*2+b+1);CHKERRQ(ierr);
333c4762a1bSJed Brown             }
334c4762a1bSJed Brown           }
335c4762a1bSJed Brown         }
336c4762a1bSJed Brown         ierr = DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr);
337c4762a1bSJed Brown       }
338c4762a1bSJed Brown       ierr = ISRestoreIndices(is, &faces);CHKERRQ(ierr);
339c4762a1bSJed Brown     }
340c4762a1bSJed Brown     ierr = ISDestroy(&is);CHKERRQ(ierr);
341c4762a1bSJed Brown   }
342c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
343c4762a1bSJed Brown   PetscFunctionReturn(0);
344c4762a1bSJed Brown }
345c4762a1bSJed Brown 
346c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user)
347c4762a1bSJed Brown {
34845480ffeSMatthew G. Knepley   PetscDS        ds;
34945480ffeSMatthew G. Knepley   PetscWeakForm  wf;
35045480ffeSMatthew G. Knepley   DMLabel        label;
35145480ffeSMatthew G. Knepley   PetscInt       bd;
352c4762a1bSJed Brown   PetscErrorCode ierr;
353c4762a1bSJed Brown 
354c4762a1bSJed Brown   PetscFunctionBeginUser;
35545480ffeSMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
35645480ffeSMatthew G. Knepley   ierr = PetscDSSetResidual(ds, 0, NULL, f1_u_3d);CHKERRQ(ierr);
35745480ffeSMatthew G. Knepley   ierr = PetscDSSetResidual(ds, 1, f0_p_3d, NULL);CHKERRQ(ierr);
35845480ffeSMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu_3d);CHKERRQ(ierr);
35945480ffeSMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up_3d, NULL);CHKERRQ(ierr);
36045480ffeSMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL,  NULL);CHKERRQ(ierr);
361c4762a1bSJed Brown 
36245480ffeSMatthew G. Knepley   ierr = DMGetLabel(dm, "Faces", &label);CHKERRQ(ierr);
36345480ffeSMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr);
36445480ffeSMatthew G. Knepley   ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
36506ad1575SMatthew G. Knepley   ierr = PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_3d, 0, NULL);CHKERRQ(ierr);
36606ad1575SMatthew G. Knepley   ierr = PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL);CHKERRQ(ierr);
367c4762a1bSJed Brown 
36845480ffeSMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (void (*)(void)) coordinates, NULL, user, NULL);CHKERRQ(ierr);
369c4762a1bSJed Brown   PetscFunctionReturn(0);
370c4762a1bSJed Brown }
371c4762a1bSJed Brown 
372c4762a1bSJed Brown PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
373c4762a1bSJed Brown {
374c4762a1bSJed Brown   PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure};
375c4762a1bSJed Brown   Vec            A;
376c4762a1bSJed Brown   void          *ctxs[2];
377c4762a1bSJed Brown   PetscErrorCode ierr;
378c4762a1bSJed Brown 
379c4762a1bSJed Brown   PetscFunctionBegin;
380c4762a1bSJed Brown   ctxs[0] = user; ctxs[1] = user;
381c4762a1bSJed Brown   ierr = DMCreateLocalVector(dmAux, &A);CHKERRQ(ierr);
382c4762a1bSJed Brown   ierr = DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A);CHKERRQ(ierr);
3839a2a23afSMatthew G. Knepley   ierr = DMSetAuxiliaryVec(dm, NULL, 0, A);CHKERRQ(ierr);
384c4762a1bSJed Brown   ierr = VecDestroy(&A);CHKERRQ(ierr);
385c4762a1bSJed Brown   PetscFunctionReturn(0);
386c4762a1bSJed Brown }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user)
389c4762a1bSJed Brown {
390c4762a1bSJed Brown   /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
391c4762a1bSJed Brown   DM             subdm;
392c4762a1bSJed Brown   MatNullSpace   nearNullSpace;
393c4762a1bSJed Brown   PetscInt       fields = 0;
394c4762a1bSJed Brown   PetscObject    deformation;
395c4762a1bSJed Brown   PetscErrorCode ierr;
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   PetscFunctionBeginUser;
398c4762a1bSJed Brown   ierr = DMCreateSubDM(dm, 1, &fields, NULL, &subdm);CHKERRQ(ierr);
39956cf3b9cSMatthew G. Knepley   ierr = DMPlexCreateRigidBody(subdm, 0, &nearNullSpace);CHKERRQ(ierr);
400c4762a1bSJed Brown   ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr);
401c4762a1bSJed Brown   ierr = PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace);CHKERRQ(ierr);
402c4762a1bSJed Brown   ierr = DMDestroy(&subdm);CHKERRQ(ierr);
403c4762a1bSJed Brown   ierr = MatNullSpaceDestroy(&nearNullSpace);CHKERRQ(ierr);
404c4762a1bSJed Brown   PetscFunctionReturn(0);
405c4762a1bSJed Brown }
406c4762a1bSJed Brown 
407c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user)
408c4762a1bSJed Brown {
409c4762a1bSJed Brown   DM             dmAux, coordDM;
410c4762a1bSJed Brown   PetscInt       f;
411c4762a1bSJed Brown   PetscErrorCode ierr;
412c4762a1bSJed Brown 
413c4762a1bSJed Brown   PetscFunctionBegin;
414c4762a1bSJed Brown   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
415c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
416c4762a1bSJed Brown   if (!feAux) PetscFunctionReturn(0);
417c4762a1bSJed Brown   ierr = DMClone(dm, &dmAux);CHKERRQ(ierr);
418c4762a1bSJed Brown   ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr);
419c4762a1bSJed Brown   for (f = 0; f < NfAux; ++f) {ierr = DMSetField(dmAux, f, NULL, (PetscObject) feAux[f]);CHKERRQ(ierr);}
420c4762a1bSJed Brown   ierr = DMCreateDS(dmAux);CHKERRQ(ierr);
421c4762a1bSJed Brown   ierr = SetupMaterial(dm, dmAux, user);CHKERRQ(ierr);
422c4762a1bSJed Brown   ierr = DMDestroy(&dmAux);CHKERRQ(ierr);
423c4762a1bSJed Brown   PetscFunctionReturn(0);
424c4762a1bSJed Brown }
425c4762a1bSJed Brown 
426c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
427c4762a1bSJed Brown {
428c4762a1bSJed Brown   DM              cdm = dm;
429c4762a1bSJed Brown   PetscFE         fe[2], feAux[2];
43030602db0SMatthew G. Knepley   PetscBool       simplex;
43130602db0SMatthew G. Knepley   PetscInt        dim;
432c4762a1bSJed Brown   MPI_Comm        comm;
433c4762a1bSJed Brown   PetscErrorCode  ierr;
434c4762a1bSJed Brown 
435c4762a1bSJed Brown   PetscFunctionBeginUser;
436c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
43730602db0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
43830602db0SMatthew G. Knepley   ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr);
439c4762a1bSJed Brown   /* Create finite element */
440c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
441c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[0], "deformation");CHKERRQ(ierr);
442c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
443c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
446c4762a1bSJed Brown 
447c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]);CHKERRQ(ierr);
448c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) feAux[0], "elasticityMaterial");CHKERRQ(ierr);
449c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], feAux[0]);CHKERRQ(ierr);
450c4762a1bSJed 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 */
451c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]);CHKERRQ(ierr);
452c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) feAux[1], "wall_pressure");CHKERRQ(ierr);
453c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], feAux[1]);CHKERRQ(ierr);
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
456c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
457c4762a1bSJed Brown   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
458c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
459c4762a1bSJed Brown   ierr = SetupProblem(dm, dim, user);CHKERRQ(ierr);
460c4762a1bSJed Brown   while (cdm) {
461c4762a1bSJed Brown     ierr = SetupAuxDM(cdm, 2, feAux, user);CHKERRQ(ierr);
462408cafa0SMatthew G. Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
463c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
464c4762a1bSJed Brown   }
465c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
466c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
467c4762a1bSJed Brown   ierr = PetscFEDestroy(&feAux[0]);CHKERRQ(ierr);
468c4762a1bSJed Brown   ierr = PetscFEDestroy(&feAux[1]);CHKERRQ(ierr);
469c4762a1bSJed Brown   PetscFunctionReturn(0);
470c4762a1bSJed Brown }
471c4762a1bSJed Brown 
472c4762a1bSJed Brown int main(int argc, char **argv)
473c4762a1bSJed Brown {
474c4762a1bSJed Brown   SNES           snes;                 /* nonlinear solver */
475c4762a1bSJed Brown   DM             dm;                   /* problem definition */
476c4762a1bSJed Brown   Vec            u,r;                  /* solution, residual vectors */
477c4762a1bSJed Brown   Mat            A,J;                  /* Jacobian matrix */
478c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
479c4762a1bSJed Brown   PetscInt       its;                  /* iterations for convergence */
480c4762a1bSJed Brown   PetscErrorCode ierr;
481c4762a1bSJed Brown 
482c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
483c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
484c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
485c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
486c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
487c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
490c4762a1bSJed Brown   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
491c4762a1bSJed Brown   ierr = SetupNearNullSpace(dm, &user);CHKERRQ(ierr);
492c4762a1bSJed Brown 
493c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
494c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   ierr = DMSetMatType(dm,MATAIJ);CHKERRQ(ierr);
497c4762a1bSJed Brown   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
498c4762a1bSJed Brown   A = J;
499c4762a1bSJed Brown 
500c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
501c4762a1bSJed Brown   ierr = SNESSetJacobian(snes, A, J, NULL, NULL);CHKERRQ(ierr);
502c4762a1bSJed Brown 
503c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
504c4762a1bSJed Brown 
505c4762a1bSJed Brown   {
506c4762a1bSJed Brown     PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx);
507c4762a1bSJed Brown     initialGuess[0] = coordinates; initialGuess[1] = zero_scalar;
508c4762a1bSJed Brown     ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
509c4762a1bSJed Brown   }
510c4762a1bSJed Brown 
511c4762a1bSJed Brown   if (user.runType == RUN_FULL) {
512c4762a1bSJed Brown     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
513c4762a1bSJed Brown     ierr = SNESGetIterationNumber(snes, &its);CHKERRQ(ierr);
514c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);CHKERRQ(ierr);
515c4762a1bSJed Brown   } else {
516c4762a1bSJed Brown     PetscReal res = 0.0;
517c4762a1bSJed Brown 
518c4762a1bSJed Brown     /* Check initial guess */
519c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr);
520c4762a1bSJed Brown     ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
521c4762a1bSJed Brown     /* Check residual */
522c4762a1bSJed Brown     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
523c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr);
524c4762a1bSJed Brown     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
525c4762a1bSJed Brown     ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
526c4762a1bSJed Brown     ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
527c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
528c4762a1bSJed Brown     /* Check Jacobian */
529c4762a1bSJed Brown     {
530c4762a1bSJed Brown       Vec b;
531c4762a1bSJed Brown 
532c4762a1bSJed Brown       ierr = SNESComputeJacobian(snes, u, A, A);CHKERRQ(ierr);
533c4762a1bSJed Brown       ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
534c4762a1bSJed Brown       ierr = VecSet(r, 0.0);CHKERRQ(ierr);
535c4762a1bSJed Brown       ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr);
536c4762a1bSJed Brown       ierr = MatMult(A, u, r);CHKERRQ(ierr);
537c4762a1bSJed Brown       ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr);
538c4762a1bSJed Brown       ierr = VecDestroy(&b);CHKERRQ(ierr);
539c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr);
540c4762a1bSJed Brown       ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
541c4762a1bSJed Brown       ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
542c4762a1bSJed Brown       ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
543c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
544c4762a1bSJed Brown     }
545c4762a1bSJed Brown   }
546c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
547c4762a1bSJed Brown 
548c4762a1bSJed Brown   if (A != J) {ierr = MatDestroy(&A);CHKERRQ(ierr);}
549c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
550c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
551c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
552c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
553c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
554c4762a1bSJed Brown   ierr = PetscFinalize();
555c4762a1bSJed Brown   return ierr;
556c4762a1bSJed Brown }
557c4762a1bSJed Brown 
558c4762a1bSJed Brown /*TEST
559c4762a1bSJed Brown 
560c4762a1bSJed Brown   build:
561c4762a1bSJed Brown     requires: !complex
562c4762a1bSJed Brown 
56330602db0SMatthew G. Knepley   testset:
56430602db0SMatthew G. Knepley     requires: ctetgen
56530602db0SMatthew G. Knepley     args: -run_type full -dm_plex_dim 3 \
56630602db0SMatthew G. Knepley           -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 \
56730602db0SMatthew G. Knepley           -snes_rtol 1e-05 -snes_monitor_short -snes_converged_reason \
56830602db0SMatthew G. Knepley           -ksp_type fgmres -ksp_rtol 1e-10 -ksp_monitor_short -ksp_converged_reason \
56930602db0SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper \
57030602db0SMatthew G. Knepley             -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu \
57130602db0SMatthew G. Knepley             -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
57230602db0SMatthew G. Knepley 
573c4762a1bSJed Brown     test:
574c4762a1bSJed Brown       suffix: 0
57530602db0SMatthew G. Knepley       requires: !single
57630602db0SMatthew G. Knepley       args: -dm_refine 2 \
57730602db0SMatthew G. Knepley             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
578c4762a1bSJed Brown 
579c4762a1bSJed Brown     test:
580c4762a1bSJed Brown       suffix: 1
58130602db0SMatthew G. Knepley       requires: superlu_dist
582c4762a1bSJed Brown       nsize: 2
583*e600fa54SMatthew G. Knepley       args: -dm_refine 0 -petscpartitioner_type simple \
58430602db0SMatthew G. Knepley             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
585c4762a1bSJed Brown       timeoutfactor: 2
586c4762a1bSJed Brown 
587c4762a1bSJed Brown     test:
588c4762a1bSJed Brown       suffix: 4
58930602db0SMatthew G. Knepley       requires: superlu_dist
590c4762a1bSJed Brown       nsize: 2
591*e600fa54SMatthew G. Knepley       args: -dm_refine 0 -petscpartitioner_type simple \
59230602db0SMatthew G. Knepley             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
593c4762a1bSJed Brown       output_file: output/ex77_1.out
594c4762a1bSJed Brown 
595c4762a1bSJed Brown     test:
59630602db0SMatthew G. Knepley       suffix: 2
59730602db0SMatthew G. Knepley       requires: !single
59830602db0SMatthew G. Knepley       args: -dm_refine 2 \
59930602db0SMatthew G. Knepley             -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
600c4762a1bSJed Brown 
601c4762a1bSJed Brown     test:
602c4762a1bSJed Brown       suffix: 2_par
60330602db0SMatthew G. Knepley       requires: superlu_dist !single
604c4762a1bSJed Brown       nsize: 4
605*e600fa54SMatthew G. Knepley       args: -dm_refine 2 -petscpartitioner_type simple \
60630602db0SMatthew G. Knepley             -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
607c4762a1bSJed Brown       output_file: output/ex77_2.out
608c4762a1bSJed Brown 
609c4762a1bSJed Brown TEST*/
610