xref: /petsc/src/snes/tutorials/ex16.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
1c4762a1bSJed Brown static char help[] = "Large-deformation Elasticity Buckling Example";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*F-----------------------------------------------------------------------
4c4762a1bSJed Brown 
5c4762a1bSJed Brown     This example solves the 3D large deformation elasticity problem
6c4762a1bSJed Brown 
7c4762a1bSJed Brown \begin{equation}
8c4762a1bSJed Brown  \int_{\Omega}F \cdot S : \nabla v d\Omega + \int_{\Omega} (loading)\mathbf{e}_y \cdot v d\Omega = 0
9c4762a1bSJed Brown \end{equation}
10c4762a1bSJed Brown 
11c4762a1bSJed Brown     F is the deformation gradient, and S is the second Piola-Kirchhoff tensor from the Saint Venant-Kirchhoff model of
12c4762a1bSJed Brown     hyperelasticity.  \Omega is a (arc) angle subsection of a cylindrical shell of thickness (height), inner radius
13c4762a1bSJed Brown     (rad) and width (width).  The problem is discretized using Q1 finite elements on a logically structured grid.
14da81f932SPierre Jolivet     Homogeneous Dirichlet boundary conditions are applied at the centers of the ends of the sphere.
15c4762a1bSJed Brown 
16c4762a1bSJed Brown     This example is tunable with the following options:
17c4762a1bSJed Brown     -rad : the radius of the circle
18c4762a1bSJed Brown     -arc : set the angle of the arch represented
19c4762a1bSJed Brown     -loading : set the bulk loading (the mass)
20c4762a1bSJed Brown     -ploading : set the point loading at the top
21c4762a1bSJed Brown     -height : set the height of the arch
22c4762a1bSJed Brown     -width : set the width of the arch
23c4762a1bSJed Brown     -view_line : print initial and final offsets of the centerline of the
24c4762a1bSJed Brown                  beam along the x direction
25c4762a1bSJed Brown 
26c4762a1bSJed Brown     The material properties may be modified using either:
27c4762a1bSJed Brown     -mu      : the first lame' parameter
28c4762a1bSJed Brown     -lambda  : the second lame' parameter
29c4762a1bSJed Brown 
30c4762a1bSJed Brown     Or:
31c4762a1bSJed Brown     -young   : the Young's modulus
32c4762a1bSJed Brown     -poisson : the poisson ratio
33c4762a1bSJed Brown 
34c4762a1bSJed Brown     This example is meant to show the strain placed upon the nonlinear solvers when trying to "snap through" the arch
35c4762a1bSJed Brown     using the loading.  Under certain parameter regimes, the arch will invert under the load, and the number of Newton
3697276fddSZach Atkins     steps will jump considerably.  Composed nonlinear solvers may be used to mitigate this difficulty.  This example
3797276fddSZach Atkins     also demonstrates the use of the arc length continuation method NEWTONAL, which avoids the numerical difficulties
3897276fddSZach Atkins     of the snap-through via tracing the equilibrium path through load increments.
39c4762a1bSJed Brown 
40c4762a1bSJed Brown     The initial setup follows the example in pg. 268 of "Nonlinear Finite Element Methods" by Peter Wriggers, but is a
41c4762a1bSJed Brown     3D extension.
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   ------------------------------------------------------------------------F*/
44c4762a1bSJed Brown 
45c4762a1bSJed Brown #include <petscsnes.h>
46c4762a1bSJed Brown #include <petscdm.h>
47c4762a1bSJed Brown #include <petscdmda.h>
48c4762a1bSJed Brown 
49c4762a1bSJed Brown #define QP0 0.2113248654051871
50c4762a1bSJed Brown #define QP1 0.7886751345948129
51c4762a1bSJed Brown #define NQ  2
52c4762a1bSJed Brown #define NB  2
53c4762a1bSJed Brown #define NEB 8
54c4762a1bSJed Brown #define NEQ 8
55c4762a1bSJed Brown #define NPB 24
56c4762a1bSJed Brown 
57c4762a1bSJed Brown #define NVALS NEB *NEQ
58c4762a1bSJed Brown const PetscReal pts[NQ] = {QP0, QP1};
59c4762a1bSJed Brown const PetscReal wts[NQ] = {0.5, 0.5};
60c4762a1bSJed Brown 
61c4762a1bSJed Brown PetscScalar vals[NVALS];
62c4762a1bSJed Brown PetscScalar grad[3 * NVALS];
63c4762a1bSJed Brown 
64c4762a1bSJed Brown typedef PetscScalar Field[3];
65c4762a1bSJed Brown typedef PetscScalar CoordField[3];
66c4762a1bSJed Brown 
67c4762a1bSJed Brown typedef PetscScalar JacField[9];
68c4762a1bSJed Brown 
69c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *, Field ***, Field ***, void *);
70c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *, Field ***, Mat, Mat, void *);
71c4762a1bSJed Brown PetscErrorCode DisplayLine(SNES, Vec);
723274405dSPierre Jolivet PetscErrorCode FormElements(void);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown typedef struct {
75c4762a1bSJed Brown   PetscReal loading;
76c4762a1bSJed Brown   PetscReal mu;
77c4762a1bSJed Brown   PetscReal lambda;
78c4762a1bSJed Brown   PetscReal rad;
79c4762a1bSJed Brown   PetscReal height;
80c4762a1bSJed Brown   PetscReal width;
81c4762a1bSJed Brown   PetscReal arc;
82c4762a1bSJed Brown   PetscReal ploading;
8397276fddSZach Atkins   PetscReal load_factor;
84c4762a1bSJed Brown } AppCtx;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown PetscErrorCode InitialGuess(DM, AppCtx *, Vec);
87c4762a1bSJed Brown PetscErrorCode FormRHS(DM, AppCtx *, Vec);
88c4762a1bSJed Brown PetscErrorCode FormCoordinates(DM, AppCtx *);
8997276fddSZach Atkins PetscErrorCode TangentLoad(SNES, Vec, Vec, void *);
90c4762a1bSJed Brown 
main(int argc,char ** argv)91d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
92d71ae5a4SJacob Faibussowitsch {
93c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
94c4762a1bSJed Brown   PetscInt  mx, my, its;
95c4762a1bSJed Brown   MPI_Comm  comm;
96c4762a1bSJed Brown   SNES      snes;
97c4762a1bSJed Brown   DM        da;
98c4762a1bSJed Brown   Vec       x, X, b;
9997276fddSZach Atkins   PetscBool youngflg, poissonflg, muflg, lambdaflg, alflg, view = PETSC_FALSE, viewline = PETSC_FALSE;
100c4762a1bSJed Brown   PetscReal poisson = 0.2, young = 4e4;
101c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN]     = "ex16.vts";
102c4762a1bSJed Brown   char      filename_def[PETSC_MAX_PATH_LEN] = "ex16_def.vts";
103c4762a1bSJed Brown 
104327415f7SBarry Smith   PetscFunctionBeginUser;
105c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1069566063dSJacob Faibussowitsch   PetscCall(FormElements());
107c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1089566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &snes));
1099566063dSJacob Faibussowitsch   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 11, 2, 2, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 3, 1, NULL, NULL, NULL, &da));
1109566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1119566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1129566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, (DM)da));
113c4762a1bSJed Brown 
1149566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
115c4762a1bSJed Brown   user.loading     = 0.0;
116c4762a1bSJed Brown   user.arc         = PETSC_PI / 3.;
117c4762a1bSJed Brown   user.mu          = 4.0;
118c4762a1bSJed Brown   user.lambda      = 1.0;
119c4762a1bSJed Brown   user.rad         = 100.0;
120c4762a1bSJed Brown   user.height      = 3.;
121c4762a1bSJed Brown   user.width       = 1.;
122c4762a1bSJed Brown   user.ploading    = -5e3;
12397276fddSZach Atkins   user.load_factor = 1.0;
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-arc", &user.arc, NULL));
1269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, &muflg));
1279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lambda", &user.lambda, &lambdaflg));
1289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-rad", &user.rad, NULL));
1299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-height", &user.height, NULL));
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-width", &user.width, NULL));
1319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-loading", &user.loading, NULL));
1329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-ploading", &user.ploading, NULL));
1339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-poisson", &poisson, &poissonflg));
1349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-young", &young, &youngflg));
135*b6555650SPierre Jolivet   if (youngflg || poissonflg || !(muflg || lambdaflg)) {
136c4762a1bSJed Brown     /* set the lame' parameters based upon the poisson ratio and young's modulus */
137c4762a1bSJed Brown     user.lambda = poisson * young / ((1. + poisson) * (1. - 2. * poisson));
138c4762a1bSJed Brown     user.mu     = young / (2. * (1. + poisson));
139c4762a1bSJed Brown   }
1409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view", &view, NULL));
1419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_line", &viewline, NULL));
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "x_disp"));
1449566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "y_disp"));
1459566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 2, "z_disp"));
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
1489566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode (*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user));
1498434afd1SBarry Smith   PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobianFn *)FormJacobianLocal, &user));
1509566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
15197276fddSZach Atkins   PetscCall(SNESNewtonALSetFunction(snes, TangentLoad, &user));
15297276fddSZach Atkins   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONAL, &alflg));
15397276fddSZach Atkins   if (alflg) user.load_factor = 0.0;
15497276fddSZach Atkins 
1559566063dSJacob Faibussowitsch   PetscCall(FormCoordinates(da, &user));
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
1589566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &b));
1599566063dSJacob Faibussowitsch   PetscCall(InitialGuess(da, &user, x));
1609566063dSJacob Faibussowitsch   PetscCall(FormRHS(da, &user, b));
161c4762a1bSJed Brown 
1629566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "lambda: %f mu: %f\n", (double)user.lambda, (double)user.mu));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* show a cross-section of the initial state */
1651baa6e33SBarry Smith   if (viewline) PetscCall(DisplayLine(snes, x));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /* get the loaded configuration */
1689566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, b, x));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
17163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
1729566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &X));
173c4762a1bSJed Brown   /* show a cross-section of the final state */
1741baa6e33SBarry Smith   if (viewline) PetscCall(DisplayLine(snes, X));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   if (view) {
177c4762a1bSJed Brown     PetscViewer viewer;
178c4762a1bSJed Brown     Vec         coords;
1799566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD, filename, FILE_MODE_WRITE, &viewer));
1809566063dSJacob Faibussowitsch     PetscCall(VecView(x, viewer));
1819566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1829566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da, &coords));
1839566063dSJacob Faibussowitsch     PetscCall(VecAXPY(coords, 1.0, x));
1849566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD, filename_def, FILE_MODE_WRITE, &viewer));
1859566063dSJacob Faibussowitsch     PetscCall(VecView(x, viewer));
1869566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
187c4762a1bSJed Brown   }
188c4762a1bSJed Brown 
1899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1919566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1929566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1939566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
194b122ec5aSJacob Faibussowitsch   return 0;
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
OnBoundary(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz)197d71ae5a4SJacob Faibussowitsch PetscInt OnBoundary(PetscInt i, PetscInt j, PetscInt k, PetscInt mx, PetscInt my, PetscInt mz)
198d71ae5a4SJacob Faibussowitsch {
199c4762a1bSJed Brown   if ((i == 0 || i == mx - 1) && j == my - 1) return 1;
200c4762a1bSJed Brown   return 0;
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
BoundaryValue(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz,PetscScalar * val,AppCtx * user)203d71ae5a4SJacob Faibussowitsch void BoundaryValue(PetscInt i, PetscInt j, PetscInt k, PetscInt mx, PetscInt my, PetscInt mz, PetscScalar *val, AppCtx *user)
204d71ae5a4SJacob Faibussowitsch {
205c4762a1bSJed Brown   /* reference coordinates */
206f4f49eeaSPierre Jolivet   PetscReal p_x = ((PetscReal)i) / ((PetscReal)(mx - 1));
207f4f49eeaSPierre Jolivet   PetscReal p_y = ((PetscReal)j) / ((PetscReal)(my - 1));
208f4f49eeaSPierre Jolivet   PetscReal p_z = ((PetscReal)k) / ((PetscReal)(mz - 1));
209c4762a1bSJed Brown   PetscReal o_x = p_x;
210c4762a1bSJed Brown   PetscReal o_y = p_y;
211c4762a1bSJed Brown   PetscReal o_z = p_z;
212c4762a1bSJed Brown   val[0]        = o_x - p_x;
213c4762a1bSJed Brown   val[1]        = o_y - p_y;
214c4762a1bSJed Brown   val[2]        = o_z - p_z;
215c4762a1bSJed Brown }
216c4762a1bSJed Brown 
InvertTensor(PetscScalar * t,PetscScalar * ti,PetscReal * dett)217d71ae5a4SJacob Faibussowitsch void InvertTensor(PetscScalar *t, PetscScalar *ti, PetscReal *dett)
218d71ae5a4SJacob Faibussowitsch {
219c4762a1bSJed Brown   const PetscScalar a   = t[0];
220c4762a1bSJed Brown   const PetscScalar b   = t[1];
221c4762a1bSJed Brown   const PetscScalar c   = t[2];
222c4762a1bSJed Brown   const PetscScalar d   = t[3];
223c4762a1bSJed Brown   const PetscScalar e   = t[4];
224c4762a1bSJed Brown   const PetscScalar f   = t[5];
225c4762a1bSJed Brown   const PetscScalar g   = t[6];
226c4762a1bSJed Brown   const PetscScalar h   = t[7];
227c4762a1bSJed Brown   const PetscScalar i   = t[8];
228c4762a1bSJed Brown   const PetscReal   det = PetscRealPart(a * (e * i - f * h) - b * (i * d - f * g) + c * (d * h - e * g));
229c4762a1bSJed Brown   const PetscReal   di  = 1. / det;
230c4762a1bSJed Brown   if (ti) {
231c4762a1bSJed Brown     const PetscScalar A  = (e * i - f * h);
232c4762a1bSJed Brown     const PetscScalar B  = -(d * i - f * g);
233c4762a1bSJed Brown     const PetscScalar C  = (d * h - e * g);
234c4762a1bSJed Brown     const PetscScalar D  = -(b * i - c * h);
235c4762a1bSJed Brown     const PetscScalar E  = (a * i - c * g);
236c4762a1bSJed Brown     const PetscScalar F  = -(a * h - b * g);
237c4762a1bSJed Brown     const PetscScalar G  = (b * f - c * e);
238c4762a1bSJed Brown     const PetscScalar H  = -(a * f - c * d);
239c4762a1bSJed Brown     const PetscScalar II = (a * e - b * d);
240c4762a1bSJed Brown     ti[0]                = di * A;
241c4762a1bSJed Brown     ti[1]                = di * D;
242c4762a1bSJed Brown     ti[2]                = di * G;
243c4762a1bSJed Brown     ti[3]                = di * B;
244c4762a1bSJed Brown     ti[4]                = di * E;
245c4762a1bSJed Brown     ti[5]                = di * H;
246c4762a1bSJed Brown     ti[6]                = di * C;
247c4762a1bSJed Brown     ti[7]                = di * F;
248c4762a1bSJed Brown     ti[8]                = di * II;
249c4762a1bSJed Brown   }
250c4762a1bSJed Brown   if (dett) *dett = det;
251c4762a1bSJed Brown }
252c4762a1bSJed Brown 
TensorTensor(PetscScalar * a,PetscScalar * b,PetscScalar * c)253d71ae5a4SJacob Faibussowitsch void TensorTensor(PetscScalar *a, PetscScalar *b, PetscScalar *c)
254d71ae5a4SJacob Faibussowitsch {
255c4762a1bSJed Brown   PetscInt i, j, m;
256c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
257c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
258c4762a1bSJed Brown       c[i + 3 * j] = 0;
2599371c9d4SSatish Balay       for (m = 0; m < 3; m++) c[i + 3 * j] += a[m + 3 * j] * b[i + 3 * m];
260c4762a1bSJed Brown     }
261c4762a1bSJed Brown   }
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
TensorTransposeTensor(PetscScalar * a,PetscScalar * b,PetscScalar * c)264d71ae5a4SJacob Faibussowitsch void TensorTransposeTensor(PetscScalar *a, PetscScalar *b, PetscScalar *c)
265d71ae5a4SJacob Faibussowitsch {
266c4762a1bSJed Brown   PetscInt i, j, m;
267c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
268c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
269c4762a1bSJed Brown       c[i + 3 * j] = 0;
2709371c9d4SSatish Balay       for (m = 0; m < 3; m++) c[i + 3 * j] += a[3 * m + j] * b[i + 3 * m];
271c4762a1bSJed Brown     }
272c4762a1bSJed Brown   }
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
TensorVector(PetscScalar * rot,PetscScalar * vec,PetscScalar * tvec)275d71ae5a4SJacob Faibussowitsch void TensorVector(PetscScalar *rot, PetscScalar *vec, PetscScalar *tvec)
276d71ae5a4SJacob Faibussowitsch {
277c4762a1bSJed Brown   tvec[0] = rot[0] * vec[0] + rot[1] * vec[1] + rot[2] * vec[2];
278c4762a1bSJed Brown   tvec[1] = rot[3] * vec[0] + rot[4] * vec[1] + rot[5] * vec[2];
279c4762a1bSJed Brown   tvec[2] = rot[6] * vec[0] + rot[7] * vec[1] + rot[8] * vec[2];
280c4762a1bSJed Brown }
281c4762a1bSJed Brown 
DeformationGradient(Field * ex,PetscInt qi,PetscInt qj,PetscInt qk,PetscScalar * invJ,PetscScalar * F)282d71ae5a4SJacob Faibussowitsch void DeformationGradient(Field *ex, PetscInt qi, PetscInt qj, PetscInt qk, PetscScalar *invJ, PetscScalar *F)
283d71ae5a4SJacob Faibussowitsch {
284c4762a1bSJed Brown   PetscInt ii, jj, kk, l;
285ad540459SPierre Jolivet   for (l = 0; l < 9; l++) F[l] = 0.;
286c4762a1bSJed Brown   F[0] = 1.;
287c4762a1bSJed Brown   F[4] = 1.;
288c4762a1bSJed Brown   F[8] = 1.;
289c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
290c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
291c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
292c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
293c4762a1bSJed Brown         PetscInt    idx  = ii + jj * NB + kk * NB * NB;
294c4762a1bSJed Brown         PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
295c4762a1bSJed Brown         PetscScalar lgrad[3];
296c4762a1bSJed Brown         TensorVector(invJ, &grad[3 * bidx], lgrad);
2979371c9d4SSatish Balay         F[0] += lgrad[0] * ex[idx][0];
2989371c9d4SSatish Balay         F[1] += lgrad[1] * ex[idx][0];
2999371c9d4SSatish Balay         F[2] += lgrad[2] * ex[idx][0];
3009371c9d4SSatish Balay         F[3] += lgrad[0] * ex[idx][1];
3019371c9d4SSatish Balay         F[4] += lgrad[1] * ex[idx][1];
3029371c9d4SSatish Balay         F[5] += lgrad[2] * ex[idx][1];
3039371c9d4SSatish Balay         F[6] += lgrad[0] * ex[idx][2];
3049371c9d4SSatish Balay         F[7] += lgrad[1] * ex[idx][2];
3059371c9d4SSatish Balay         F[8] += lgrad[2] * ex[idx][2];
306c4762a1bSJed Brown       }
307c4762a1bSJed Brown     }
308c4762a1bSJed Brown   }
309c4762a1bSJed Brown }
310c4762a1bSJed Brown 
DeformationGradientJacobian(PetscInt qi,PetscInt qj,PetscInt qk,PetscInt ii,PetscInt jj,PetscInt kk,PetscInt fld,PetscScalar * invJ,PetscScalar * dF)311d71ae5a4SJacob Faibussowitsch void DeformationGradientJacobian(PetscInt qi, PetscInt qj, PetscInt qk, PetscInt ii, PetscInt jj, PetscInt kk, PetscInt fld, PetscScalar *invJ, PetscScalar *dF)
312d71ae5a4SJacob Faibussowitsch {
313c4762a1bSJed Brown   PetscInt    l;
314c4762a1bSJed Brown   PetscScalar lgrad[3];
315c4762a1bSJed Brown   PetscInt    idx  = ii + jj * NB + kk * NB * NB;
316c4762a1bSJed Brown   PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
317ad540459SPierre Jolivet   for (l = 0; l < 9; l++) dF[l] = 0.;
318c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
319c4762a1bSJed Brown   TensorVector(invJ, &grad[3 * bidx], lgrad);
3209371c9d4SSatish Balay   dF[3 * fld]     = lgrad[0];
3219371c9d4SSatish Balay   dF[3 * fld + 1] = lgrad[1];
3229371c9d4SSatish Balay   dF[3 * fld + 2] = lgrad[2];
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
LagrangeGreenStrain(PetscScalar * F,PetscScalar * E)325d71ae5a4SJacob Faibussowitsch void LagrangeGreenStrain(PetscScalar *F, PetscScalar *E)
326d71ae5a4SJacob Faibussowitsch {
327c4762a1bSJed Brown   PetscInt i, j, m;
328c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
329c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
330c4762a1bSJed Brown       E[i + 3 * j] = 0;
3319371c9d4SSatish Balay       for (m = 0; m < 3; m++) E[i + 3 * j] += 0.5 * F[3 * m + j] * F[i + 3 * m];
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown   }
334ad540459SPierre Jolivet   for (i = 0; i < 3; i++) E[i + 3 * i] -= 0.5;
335c4762a1bSJed Brown }
336c4762a1bSJed Brown 
SaintVenantKirchoff(PetscReal lambda,PetscReal mu,PetscScalar * F,PetscScalar * S)337d71ae5a4SJacob Faibussowitsch void SaintVenantKirchoff(PetscReal lambda, PetscReal mu, PetscScalar *F, PetscScalar *S)
338d71ae5a4SJacob Faibussowitsch {
339c4762a1bSJed Brown   PetscInt    i, j;
340c4762a1bSJed Brown   PetscScalar E[9];
341c4762a1bSJed Brown   PetscScalar trE = 0;
342c4762a1bSJed Brown   LagrangeGreenStrain(F, E);
343ad540459SPierre Jolivet   for (i = 0; i < 3; i++) trE += E[i + 3 * i];
344c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
345c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
346c4762a1bSJed Brown       S[i + 3 * j] = 2. * mu * E[i + 3 * j];
347ad540459SPierre Jolivet       if (i == j) S[i + 3 * j] += trE * lambda;
348c4762a1bSJed Brown     }
349c4762a1bSJed Brown   }
350c4762a1bSJed Brown }
351c4762a1bSJed Brown 
SaintVenantKirchoffJacobian(PetscReal lambda,PetscReal mu,PetscScalar * F,PetscScalar * dF,PetscScalar * dS)352d71ae5a4SJacob Faibussowitsch void SaintVenantKirchoffJacobian(PetscReal lambda, PetscReal mu, PetscScalar *F, PetscScalar *dF, PetscScalar *dS)
353d71ae5a4SJacob Faibussowitsch {
354c4762a1bSJed Brown   PetscScalar FtdF[9], dE[9];
355c4762a1bSJed Brown   PetscInt    i, j;
356c4762a1bSJed Brown   PetscScalar dtrE = 0.;
357c4762a1bSJed Brown   TensorTransposeTensor(dF, F, dE);
358c4762a1bSJed Brown   TensorTransposeTensor(F, dF, FtdF);
359c4762a1bSJed Brown   for (i = 0; i < 9; i++) dE[i] += FtdF[i];
360c4762a1bSJed Brown   for (i = 0; i < 9; i++) dE[i] *= 0.5;
361ad540459SPierre Jolivet   for (i = 0; i < 3; i++) dtrE += dE[i + 3 * i];
362c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
363c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
364c4762a1bSJed Brown       dS[i + 3 * j] = 2. * mu * dE[i + 3 * j];
365ad540459SPierre Jolivet       if (i == j) dS[i + 3 * j] += dtrE * lambda;
366c4762a1bSJed Brown     }
367c4762a1bSJed Brown   }
368c4762a1bSJed Brown }
369c4762a1bSJed Brown 
FormElements()370d71ae5a4SJacob Faibussowitsch PetscErrorCode FormElements()
371d71ae5a4SJacob Faibussowitsch {
372c4762a1bSJed Brown   PetscInt  i, j, k, ii, jj, kk;
373c4762a1bSJed Brown   PetscReal bx, by, bz, dbx, dby, dbz;
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   PetscFunctionBegin;
376c4762a1bSJed Brown   /* construct the basis function values and derivatives */
377c4762a1bSJed Brown   for (k = 0; k < NB; k++) {
378c4762a1bSJed Brown     for (j = 0; j < NB; j++) {
379c4762a1bSJed Brown       for (i = 0; i < NB; i++) {
380c4762a1bSJed Brown         /* loop over the quadrature points */
381c4762a1bSJed Brown         for (kk = 0; kk < NQ; kk++) {
382c4762a1bSJed Brown           for (jj = 0; jj < NQ; jj++) {
383c4762a1bSJed Brown             for (ii = 0; ii < NQ; ii++) {
384c4762a1bSJed Brown               PetscInt idx = ii + NQ * jj + NQ * NQ * kk + NEQ * i + NEQ * NB * j + NEQ * NB * NB * k;
385c4762a1bSJed Brown               bx           = pts[ii];
386c4762a1bSJed Brown               by           = pts[jj];
387c4762a1bSJed Brown               bz           = pts[kk];
388c4762a1bSJed Brown               dbx          = 1.;
389c4762a1bSJed Brown               dby          = 1.;
390c4762a1bSJed Brown               dbz          = 1.;
3919371c9d4SSatish Balay               if (i == 0) {
3929371c9d4SSatish Balay                 bx  = 1. - bx;
3939371c9d4SSatish Balay                 dbx = -1;
3949371c9d4SSatish Balay               }
3959371c9d4SSatish Balay               if (j == 0) {
3969371c9d4SSatish Balay                 by  = 1. - by;
3979371c9d4SSatish Balay                 dby = -1;
3989371c9d4SSatish Balay               }
3999371c9d4SSatish Balay               if (k == 0) {
4009371c9d4SSatish Balay                 bz  = 1. - bz;
4019371c9d4SSatish Balay                 dbz = -1;
4029371c9d4SSatish Balay               }
403c4762a1bSJed Brown               vals[idx]         = bx * by * bz;
404c4762a1bSJed Brown               grad[3 * idx + 0] = dbx * by * bz;
405c4762a1bSJed Brown               grad[3 * idx + 1] = dby * bx * bz;
406c4762a1bSJed Brown               grad[3 * idx + 2] = dbz * bx * by;
407c4762a1bSJed Brown             }
408c4762a1bSJed Brown           }
409c4762a1bSJed Brown         }
410c4762a1bSJed Brown       }
411c4762a1bSJed Brown     }
412c4762a1bSJed Brown   }
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
414c4762a1bSJed Brown }
415c4762a1bSJed Brown 
GatherElementData(PetscInt mx,PetscInt my,PetscInt mz,Field *** x,CoordField *** c,PetscInt i,PetscInt j,PetscInt k,Field * ex,CoordField * ec,AppCtx * user)416d71ae5a4SJacob Faibussowitsch void GatherElementData(PetscInt mx, PetscInt my, PetscInt mz, Field ***x, CoordField ***c, PetscInt i, PetscInt j, PetscInt k, Field *ex, CoordField *ec, AppCtx *user)
417d71ae5a4SJacob Faibussowitsch {
418c4762a1bSJed Brown   PetscInt m;
419c4762a1bSJed Brown   PetscInt ii, jj, kk;
420c4762a1bSJed Brown   /* gather the data -- loop over element unknowns */
421c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
422c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
423c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
424c4762a1bSJed Brown         PetscInt idx = ii + jj * NB + kk * NB * NB;
425c4762a1bSJed Brown         /* decouple the boundary nodes for the displacement variables */
426c4762a1bSJed Brown         if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz)) {
427c4762a1bSJed Brown           BoundaryValue(i + ii, j + jj, k + kk, mx, my, mz, ex[idx], user);
428c4762a1bSJed Brown         } else {
429ad540459SPierre Jolivet           for (m = 0; m < 3; m++) ex[idx][m] = x[k + kk][j + jj][i + ii][m];
430c4762a1bSJed Brown         }
431ad540459SPierre Jolivet         for (m = 0; m < 3; m++) ec[idx][m] = c[k + kk][j + jj][i + ii][m];
432c4762a1bSJed Brown       }
433c4762a1bSJed Brown     }
434c4762a1bSJed Brown   }
435c4762a1bSJed Brown }
436c4762a1bSJed Brown 
QuadraturePointGeometricJacobian(CoordField * ec,PetscInt qi,PetscInt qj,PetscInt qk,PetscScalar * J)437d71ae5a4SJacob Faibussowitsch void QuadraturePointGeometricJacobian(CoordField *ec, PetscInt qi, PetscInt qj, PetscInt qk, PetscScalar *J)
438d71ae5a4SJacob Faibussowitsch {
439c4762a1bSJed Brown   PetscInt ii, jj, kk;
440c4762a1bSJed Brown   /* construct the gradient at the given quadrature point named by i,j,k */
441ad540459SPierre Jolivet   for (ii = 0; ii < 9; ii++) J[ii] = 0;
442c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
443c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
444c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
445c4762a1bSJed Brown         PetscInt idx  = ii + jj * NB + kk * NB * NB;
446c4762a1bSJed Brown         PetscInt bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
4479371c9d4SSatish Balay         J[0] += grad[3 * bidx + 0] * ec[idx][0];
4489371c9d4SSatish Balay         J[1] += grad[3 * bidx + 1] * ec[idx][0];
4499371c9d4SSatish Balay         J[2] += grad[3 * bidx + 2] * ec[idx][0];
4509371c9d4SSatish Balay         J[3] += grad[3 * bidx + 0] * ec[idx][1];
4519371c9d4SSatish Balay         J[4] += grad[3 * bidx + 1] * ec[idx][1];
4529371c9d4SSatish Balay         J[5] += grad[3 * bidx + 2] * ec[idx][1];
4539371c9d4SSatish Balay         J[6] += grad[3 * bidx + 0] * ec[idx][2];
4549371c9d4SSatish Balay         J[7] += grad[3 * bidx + 1] * ec[idx][2];
4559371c9d4SSatish Balay         J[8] += grad[3 * bidx + 2] * ec[idx][2];
456c4762a1bSJed Brown       }
457c4762a1bSJed Brown     }
458c4762a1bSJed Brown   }
459c4762a1bSJed Brown }
460c4762a1bSJed Brown 
FormElementJacobian(Field * ex,CoordField * ec,Field * ef,Field * eq,PetscScalar * ej,AppCtx * user)46197276fddSZach Atkins void FormElementJacobian(Field *ex, CoordField *ec, Field *ef, Field *eq, PetscScalar *ej, AppCtx *user)
462d71ae5a4SJacob Faibussowitsch {
463c4762a1bSJed Brown   PetscReal   vol;
464c4762a1bSJed Brown   PetscScalar J[9];
465c4762a1bSJed Brown   PetscScalar invJ[9];
466c4762a1bSJed Brown   PetscScalar F[9], S[9], dF[9], dS[9], dFS[9], FdS[9], FS[9];
467c4762a1bSJed Brown   PetscReal   scl;
468c4762a1bSJed Brown   PetscInt    i, j, k, l, ii, jj, kk, ll, qi, qj, qk, m;
469c4762a1bSJed Brown 
4709371c9d4SSatish Balay   if (ej)
4719371c9d4SSatish Balay     for (i = 0; i < NPB * NPB; i++) ej[i] = 0.;
4729371c9d4SSatish Balay   if (ef)
4739371c9d4SSatish Balay     for (i = 0; i < NEB; i++) {
4749371c9d4SSatish Balay       ef[i][0] = 0.;
4759371c9d4SSatish Balay       ef[i][1] = 0.;
4769371c9d4SSatish Balay       ef[i][2] = 0.;
4779371c9d4SSatish Balay     }
47897276fddSZach Atkins   if (eq)
47997276fddSZach Atkins     for (i = 0; i < NEB; i++) {
48097276fddSZach Atkins       eq[i][0] = 0.;
48197276fddSZach Atkins       eq[i][1] = 0.;
48297276fddSZach Atkins       eq[i][2] = 0.;
48397276fddSZach Atkins     }
484c4762a1bSJed Brown   /* loop over quadrature */
485c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
486c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
487c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
488c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec, qi, qj, qk, J);
489c4762a1bSJed Brown         InvertTensor(J, invJ, &vol);
490c4762a1bSJed Brown         scl = vol * wts[qi] * wts[qj] * wts[qk];
491c4762a1bSJed Brown         DeformationGradient(ex, qi, qj, qk, invJ, F);
492c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda, user->mu, F, S);
493c4762a1bSJed Brown         /* form the function */
494c4762a1bSJed Brown         if (ef) {
495c4762a1bSJed Brown           TensorTensor(F, S, FS);
496c4762a1bSJed Brown           for (kk = 0; kk < NB; kk++) {
497c4762a1bSJed Brown             for (jj = 0; jj < NB; jj++) {
498c4762a1bSJed Brown               for (ii = 0; ii < NB; ii++) {
499c4762a1bSJed Brown                 PetscInt    idx  = ii + jj * NB + kk * NB * NB;
500c4762a1bSJed Brown                 PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
501c4762a1bSJed Brown                 PetscScalar lgrad[3];
502c4762a1bSJed Brown                 TensorVector(invJ, &grad[3 * bidx], lgrad);
503c4762a1bSJed Brown                 /* mu*F : grad phi_{u,v,w} */
504ad540459SPierre Jolivet                 for (m = 0; m < 3; m++) ef[idx][m] += scl * (lgrad[0] * FS[3 * m + 0] + lgrad[1] * FS[3 * m + 1] + lgrad[2] * FS[3 * m + 2]);
50597276fddSZach Atkins                 ef[idx][1] -= user->load_factor * scl * user->loading * vals[bidx];
50697276fddSZach Atkins               }
50797276fddSZach Atkins             }
50897276fddSZach Atkins           }
50997276fddSZach Atkins         }
51097276fddSZach Atkins         if (eq) {
51197276fddSZach Atkins           for (kk = 0; kk < NB; kk++) {
51297276fddSZach Atkins             for (jj = 0; jj < NB; jj++) {
51397276fddSZach Atkins               for (ii = 0; ii < NB; ii++) {
51497276fddSZach Atkins                 PetscInt idx  = ii + jj * NB + kk * NB * NB;
51597276fddSZach Atkins                 PetscInt bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
51697276fddSZach Atkins                 /* external force vector */
51797276fddSZach Atkins                 eq[idx][1] += scl * user->loading * vals[bidx];
518c4762a1bSJed Brown               }
519c4762a1bSJed Brown             }
520c4762a1bSJed Brown           }
521c4762a1bSJed Brown         }
522c4762a1bSJed Brown         /* form the jacobian */
523c4762a1bSJed Brown         if (ej) {
524c4762a1bSJed Brown           /* loop over trialfunctions */
525c4762a1bSJed Brown           for (k = 0; k < NB; k++) {
526c4762a1bSJed Brown             for (j = 0; j < NB; j++) {
527c4762a1bSJed Brown               for (i = 0; i < NB; i++) {
528c4762a1bSJed Brown                 for (l = 0; l < 3; l++) {
529c4762a1bSJed Brown                   PetscInt tridx = l + 3 * (i + j * NB + k * NB * NB);
530c4762a1bSJed Brown                   DeformationGradientJacobian(qi, qj, qk, i, j, k, l, invJ, dF);
531c4762a1bSJed Brown                   SaintVenantKirchoffJacobian(user->lambda, user->mu, F, dF, dS);
532c4762a1bSJed Brown                   TensorTensor(dF, S, dFS);
533c4762a1bSJed Brown                   TensorTensor(F, dS, FdS);
534c4762a1bSJed Brown                   for (m = 0; m < 9; m++) dFS[m] += FdS[m];
535c4762a1bSJed Brown                   /* loop over testfunctions */
536c4762a1bSJed Brown                   for (kk = 0; kk < NB; kk++) {
537c4762a1bSJed Brown                     for (jj = 0; jj < NB; jj++) {
538c4762a1bSJed Brown                       for (ii = 0; ii < NB; ii++) {
539c4762a1bSJed Brown                         PetscInt    idx  = ii + jj * NB + kk * NB * NB;
540c4762a1bSJed Brown                         PetscInt    bidx = 8 * idx + qi + NQ * qj + NQ * NQ * qk;
541c4762a1bSJed Brown                         PetscScalar lgrad[3];
542c4762a1bSJed Brown                         TensorVector(invJ, &grad[3 * bidx], lgrad);
543c4762a1bSJed Brown                         for (ll = 0; ll < 3; ll++) {
544c4762a1bSJed Brown                           PetscInt teidx = ll + 3 * (ii + jj * NB + kk * NB * NB);
5459371c9d4SSatish Balay                           ej[teidx + NPB * tridx] += scl * (lgrad[0] * dFS[3 * ll + 0] + lgrad[1] * dFS[3 * ll + 1] + lgrad[2] * dFS[3 * ll + 2]);
546c4762a1bSJed Brown                         }
547c4762a1bSJed Brown                       }
548c4762a1bSJed Brown                     }
549c4762a1bSJed Brown                   } /* end of testfunctions */
550c4762a1bSJed Brown                 }
551c4762a1bSJed Brown               }
552c4762a1bSJed Brown             }
553c4762a1bSJed Brown           } /* end of trialfunctions */
554c4762a1bSJed Brown         }
555c4762a1bSJed Brown       }
556c4762a1bSJed Brown     }
557c4762a1bSJed Brown   } /* end of quadrature points */
558c4762a1bSJed Brown }
559c4762a1bSJed Brown 
ApplyBCsElement(PetscInt mx,PetscInt my,PetscInt mz,PetscInt i,PetscInt j,PetscInt k,PetscScalar * jacobian)560d71ae5a4SJacob Faibussowitsch void ApplyBCsElement(PetscInt mx, PetscInt my, PetscInt mz, PetscInt i, PetscInt j, PetscInt k, PetscScalar *jacobian)
561d71ae5a4SJacob Faibussowitsch {
562c4762a1bSJed Brown   PetscInt ii, jj, kk, ll, ei, ej, ek, el;
563c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
564c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
565c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
566c4762a1bSJed Brown         for (ll = 0; ll < 3; ll++) {
567c4762a1bSJed Brown           PetscInt tridx = ll + 3 * (ii + jj * NB + kk * NB * NB);
568c4762a1bSJed Brown           for (ek = 0; ek < NB; ek++) {
569c4762a1bSJed Brown             for (ej = 0; ej < NB; ej++) {
570c4762a1bSJed Brown               for (ei = 0; ei < NB; ei++) {
571c4762a1bSJed Brown                 for (el = 0; el < 3; el++) {
572c4762a1bSJed Brown                   if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz) || OnBoundary(i + ei, j + ej, k + ek, mx, my, mz)) {
573c4762a1bSJed Brown                     PetscInt teidx = el + 3 * (ei + ej * NB + ek * NB * NB);
574c4762a1bSJed Brown                     if (teidx == tridx) {
575c4762a1bSJed Brown                       jacobian[tridx + NPB * teidx] = 1.;
576c4762a1bSJed Brown                     } else {
577c4762a1bSJed Brown                       jacobian[tridx + NPB * teidx] = 0.;
578c4762a1bSJed Brown                     }
579c4762a1bSJed Brown                   }
580c4762a1bSJed Brown                 }
581c4762a1bSJed Brown               }
582c4762a1bSJed Brown             }
583c4762a1bSJed Brown           }
584c4762a1bSJed Brown         }
585c4762a1bSJed Brown       }
586c4762a1bSJed Brown     }
587c4762a1bSJed Brown   }
588c4762a1bSJed Brown }
589c4762a1bSJed Brown 
FormJacobianLocal(DMDALocalInfo * info,Field *** x,Mat jacpre,Mat jac,void * ptr)590d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, Field ***x, Mat jacpre, Mat jac, void *ptr)
591d71ae5a4SJacob Faibussowitsch {
592c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
593c4762a1bSJed Brown   AppCtx     *user = (AppCtx *)ptr;
594c4762a1bSJed Brown   PetscInt    i, j, k, m, l;
595c4762a1bSJed Brown   PetscInt    ii, jj, kk;
596c4762a1bSJed Brown   PetscScalar ej[NPB * NPB];
597c4762a1bSJed Brown   PetscScalar vals[NPB * NPB];
598c4762a1bSJed Brown   Field       ex[NEB];
599c4762a1bSJed Brown   CoordField  ec[NEB];
600c4762a1bSJed Brown 
601c4762a1bSJed Brown   PetscInt      xs = info->xs, ys = info->ys, zs = info->zs;
602c4762a1bSJed Brown   PetscInt      xm = info->xm, ym = info->ym, zm = info->zm;
603c4762a1bSJed Brown   PetscInt      xes, yes, zes, xee, yee, zee;
604c4762a1bSJed Brown   PetscInt      mx = info->mx, my = info->my, mz = info->mz;
605c4762a1bSJed Brown   DM            cda;
606c4762a1bSJed Brown   CoordField ***c;
607c4762a1bSJed Brown   Vec           C;
608c4762a1bSJed Brown   PetscInt      nrows;
609c4762a1bSJed Brown   MatStencil    col[NPB], row[NPB];
610c4762a1bSJed Brown   PetscScalar   v[9];
611c4762a1bSJed Brown 
612c4762a1bSJed Brown   PetscFunctionBegin;
6139566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(info->da, &cda));
6149566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(info->da, &C));
6159566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
6169566063dSJacob Faibussowitsch   PetscCall(MatScale(jac, 0.0));
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   xes = xs;
619c4762a1bSJed Brown   yes = ys;
620c4762a1bSJed Brown   zes = zs;
621c4762a1bSJed Brown   xee = xs + xm;
622c4762a1bSJed Brown   yee = ys + ym;
623c4762a1bSJed Brown   zee = zs + zm;
624c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
625c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
626c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
627c4762a1bSJed Brown   if (xs + xm == mx) xee = xs + xm - 1;
628c4762a1bSJed Brown   if (ys + ym == my) yee = ys + ym - 1;
629c4762a1bSJed Brown   if (zs + zm == mz) zee = zs + zm - 1;
630c4762a1bSJed Brown   for (k = zes; k < zee; k++) {
631c4762a1bSJed Brown     for (j = yes; j < yee; j++) {
632c4762a1bSJed Brown       for (i = xes; i < xee; i++) {
633c4762a1bSJed Brown         GatherElementData(mx, my, mz, x, c, i, j, k, ex, ec, user);
63497276fddSZach Atkins         FormElementJacobian(ex, ec, NULL, NULL, ej, user);
635c4762a1bSJed Brown         ApplyBCsElement(mx, my, mz, i, j, k, ej);
636c4762a1bSJed Brown         nrows = 0.;
637c4762a1bSJed Brown         for (kk = 0; kk < NB; kk++) {
638c4762a1bSJed Brown           for (jj = 0; jj < NB; jj++) {
639c4762a1bSJed Brown             for (ii = 0; ii < NB; ii++) {
640c4762a1bSJed Brown               PetscInt idx = ii + jj * 2 + kk * 4;
641c4762a1bSJed Brown               for (m = 0; m < 3; m++) {
642c4762a1bSJed Brown                 col[3 * idx + m].i = i + ii;
643c4762a1bSJed Brown                 col[3 * idx + m].j = j + jj;
644c4762a1bSJed Brown                 col[3 * idx + m].k = k + kk;
645c4762a1bSJed Brown                 col[3 * idx + m].c = m;
646c4762a1bSJed Brown                 if (i + ii >= xs && i + ii < xm + xs && j + jj >= ys && j + jj < ys + ym && k + kk >= zs && k + kk < zs + zm) {
647c4762a1bSJed Brown                   row[nrows].i = i + ii;
648c4762a1bSJed Brown                   row[nrows].j = j + jj;
649c4762a1bSJed Brown                   row[nrows].k = k + kk;
650c4762a1bSJed Brown                   row[nrows].c = m;
651c4762a1bSJed Brown                   for (l = 0; l < NPB; l++) vals[NPB * nrows + l] = ej[NPB * (3 * idx + m) + l];
652c4762a1bSJed Brown                   nrows++;
653c4762a1bSJed Brown                 }
654c4762a1bSJed Brown               }
655c4762a1bSJed Brown             }
656c4762a1bSJed Brown           }
657c4762a1bSJed Brown         }
6589566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jac, nrows, row, NPB, col, vals, ADD_VALUES));
659c4762a1bSJed Brown       }
660c4762a1bSJed Brown     }
661c4762a1bSJed Brown   }
662c4762a1bSJed Brown 
6639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FLUSH_ASSEMBLY));
6649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FLUSH_ASSEMBLY));
665c4762a1bSJed Brown 
666c4762a1bSJed Brown   /* set the diagonal */
6679371c9d4SSatish Balay   v[0] = 1.;
6689371c9d4SSatish Balay   v[1] = 0.;
6699371c9d4SSatish Balay   v[2] = 0.;
6709371c9d4SSatish Balay   v[3] = 0.;
6719371c9d4SSatish Balay   v[4] = 1.;
6729371c9d4SSatish Balay   v[5] = 0.;
6739371c9d4SSatish Balay   v[6] = 0.;
6749371c9d4SSatish Balay   v[7] = 0.;
6759371c9d4SSatish Balay   v[8] = 1.;
676c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
677c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
678c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
679c4762a1bSJed Brown         if (OnBoundary(i, j, k, mx, my, mz)) {
680c4762a1bSJed Brown           for (m = 0; m < 3; m++) {
681c4762a1bSJed Brown             col[m].i = i;
682c4762a1bSJed Brown             col[m].j = j;
683c4762a1bSJed Brown             col[m].k = k;
684c4762a1bSJed Brown             col[m].c = m;
685c4762a1bSJed Brown           }
6869566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(jac, 3, col, 3, col, v, INSERT_VALUES));
687c4762a1bSJed Brown         }
688c4762a1bSJed Brown       }
689c4762a1bSJed Brown     }
690c4762a1bSJed Brown   }
691c4762a1bSJed Brown 
6929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
6939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
694c4762a1bSJed Brown 
6959566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
6963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
697c4762a1bSJed Brown }
698c4762a1bSJed Brown 
FormFunctionLocal(DMDALocalInfo * info,Field *** x,Field *** f,void * ptr)699d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field ***x, Field ***f, void *ptr)
700d71ae5a4SJacob Faibussowitsch {
701c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
702c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
703c4762a1bSJed Brown   PetscInt i, j, k, l;
704c4762a1bSJed Brown   PetscInt ii, jj, kk;
705c4762a1bSJed Brown 
706c4762a1bSJed Brown   Field      ef[NEB];
707c4762a1bSJed Brown   Field      ex[NEB];
708c4762a1bSJed Brown   CoordField ec[NEB];
709c4762a1bSJed Brown 
710c4762a1bSJed Brown   PetscInt      xs = info->xs, ys = info->ys, zs = info->zs;
711c4762a1bSJed Brown   PetscInt      xm = info->xm, ym = info->ym, zm = info->zm;
712c4762a1bSJed Brown   PetscInt      xes, yes, zes, xee, yee, zee;
713c4762a1bSJed Brown   PetscInt      mx = info->mx, my = info->my, mz = info->mz;
714c4762a1bSJed Brown   DM            cda;
715c4762a1bSJed Brown   CoordField ***c;
716c4762a1bSJed Brown   Vec           C;
717c4762a1bSJed Brown 
718c4762a1bSJed Brown   PetscFunctionBegin;
7199566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(info->da, &cda));
7209566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(info->da, &C));
7219566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
7229566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(info->da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
7239566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(info->da, &xs, &ys, &zs, &xm, &ym, &zm));
724c4762a1bSJed Brown 
725c4762a1bSJed Brown   /* loop over elements */
726c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
727c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
728c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
729ad540459SPierre Jolivet         for (l = 0; l < 3; l++) f[k][j][i][l] = 0.;
730c4762a1bSJed Brown       }
731c4762a1bSJed Brown     }
732c4762a1bSJed Brown   }
733c4762a1bSJed Brown   /* element starts and ends */
734c4762a1bSJed Brown   xes = xs;
735c4762a1bSJed Brown   yes = ys;
736c4762a1bSJed Brown   zes = zs;
737c4762a1bSJed Brown   xee = xs + xm;
738c4762a1bSJed Brown   yee = ys + ym;
739c4762a1bSJed Brown   zee = zs + zm;
740c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
741c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
742c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
743c4762a1bSJed Brown   if (xs + xm == mx) xee = xs + xm - 1;
744c4762a1bSJed Brown   if (ys + ym == my) yee = ys + ym - 1;
745c4762a1bSJed Brown   if (zs + zm == mz) zee = zs + zm - 1;
746c4762a1bSJed Brown   for (k = zes; k < zee; k++) {
747c4762a1bSJed Brown     for (j = yes; j < yee; j++) {
748c4762a1bSJed Brown       for (i = xes; i < xee; i++) {
749c4762a1bSJed Brown         GatherElementData(mx, my, mz, x, c, i, j, k, ex, ec, user);
75097276fddSZach Atkins         FormElementJacobian(ex, ec, ef, NULL, NULL, user);
751c4762a1bSJed Brown         /* put this element's additions into the residuals */
752c4762a1bSJed Brown         for (kk = 0; kk < NB; kk++) {
753c4762a1bSJed Brown           for (jj = 0; jj < NB; jj++) {
754c4762a1bSJed Brown             for (ii = 0; ii < NB; ii++) {
755c4762a1bSJed Brown               PetscInt idx = ii + jj * NB + kk * NB * NB;
756c4762a1bSJed Brown               if (k + kk >= zs && j + jj >= ys && i + ii >= xs && k + kk < zs + zm && j + jj < ys + ym && i + ii < xs + xm) {
757c4762a1bSJed Brown                 if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz)) {
7589371c9d4SSatish Balay                   for (l = 0; l < 3; l++) f[k + kk][j + jj][i + ii][l] = x[k + kk][j + jj][i + ii][l] - ex[idx][l];
759c4762a1bSJed Brown                 } else {
7609371c9d4SSatish Balay                   for (l = 0; l < 3; l++) f[k + kk][j + jj][i + ii][l] += ef[idx][l];
761c4762a1bSJed Brown                 }
762c4762a1bSJed Brown               }
763c4762a1bSJed Brown             }
764c4762a1bSJed Brown           }
765c4762a1bSJed Brown         }
766c4762a1bSJed Brown       }
767c4762a1bSJed Brown     }
768c4762a1bSJed Brown   }
7699566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
7703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
771c4762a1bSJed Brown }
772c4762a1bSJed Brown 
TangentLoad(SNES snes,Vec X,Vec Q,void * ptr)77397276fddSZach Atkins PetscErrorCode TangentLoad(SNES snes, Vec X, Vec Q, void *ptr)
774d71ae5a4SJacob Faibussowitsch {
775c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
776c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
777c4762a1bSJed Brown   PetscInt xs, ys, zs;
778c4762a1bSJed Brown   PetscInt xm, ym, zm;
779c4762a1bSJed Brown   PetscInt mx, my, mz;
78097276fddSZach Atkins   DM       da;
78197276fddSZach Atkins   Vec      Xl, Ql;
78297276fddSZach Atkins   Field ***x, ***q;
78397276fddSZach Atkins   PetscInt i, j, k, l;
78497276fddSZach Atkins   PetscInt ii, jj, kk;
78597276fddSZach Atkins 
78697276fddSZach Atkins   Field      eq[NEB];
78797276fddSZach Atkins   Field      ex[NEB];
78897276fddSZach Atkins   CoordField ec[NEB];
78997276fddSZach Atkins 
79097276fddSZach Atkins   PetscInt      xes, yes, zes, xee, yee, zee;
791c4762a1bSJed Brown   DM            cda;
792c4762a1bSJed Brown   CoordField ***c;
793c4762a1bSJed Brown   Vec           C;
794c4762a1bSJed Brown 
795c4762a1bSJed Brown   PetscFunctionBegin;
79697276fddSZach Atkins   /* update user context with current load parameter */
79797276fddSZach Atkins   PetscCall(SNESNewtonALGetLoadParameter(snes, &user->load_factor));
798c4762a1bSJed Brown 
7999566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
8009566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xl));
80197276fddSZach Atkins   PetscCall(DMGetLocalVector(da, &Ql));
80297276fddSZach Atkins   PetscCall(DMGlobalToLocal(da, X, INSERT_VALUES, Xl));
80397276fddSZach Atkins 
8049566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xl, &x));
80597276fddSZach Atkins   PetscCall(DMDAVecGetArray(da, Ql, &q));
806c4762a1bSJed Brown 
8079566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
8089566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &C));
8099566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
8109566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
8119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
812c4762a1bSJed Brown 
81397276fddSZach Atkins   /* loop over elements */
814c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
815c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
816c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
81797276fddSZach Atkins         for (l = 0; l < 3; l++) q[k][j][i][l] = 0.;
818c4762a1bSJed Brown       }
819c4762a1bSJed Brown     }
820c4762a1bSJed Brown   }
82197276fddSZach Atkins   /* element starts and ends */
82297276fddSZach Atkins   xes = xs;
82397276fddSZach Atkins   yes = ys;
82497276fddSZach Atkins   zes = zs;
82597276fddSZach Atkins   xee = xs + xm;
82697276fddSZach Atkins   yee = ys + ym;
82797276fddSZach Atkins   zee = zs + zm;
82897276fddSZach Atkins   if (xs > 0) xes = xs - 1;
82997276fddSZach Atkins   if (ys > 0) yes = ys - 1;
83097276fddSZach Atkins   if (zs > 0) zes = zs - 1;
83197276fddSZach Atkins   if (xs + xm == mx) xee = xs + xm - 1;
83297276fddSZach Atkins   if (ys + ym == my) yee = ys + ym - 1;
83397276fddSZach Atkins   if (zs + zm == mz) zee = zs + zm - 1;
83497276fddSZach Atkins   for (k = zes; k < zee; k++) {
83597276fddSZach Atkins     for (j = yes; j < yee; j++) {
83697276fddSZach Atkins       for (i = xes; i < xee; i++) {
83797276fddSZach Atkins         GatherElementData(mx, my, mz, x, c, i, j, k, ex, ec, user);
83897276fddSZach Atkins         FormElementJacobian(ex, ec, NULL, eq, NULL, user);
83997276fddSZach Atkins         /* put this element's additions into the residuals */
84097276fddSZach Atkins         for (kk = 0; kk < NB; kk++) {
84197276fddSZach Atkins           for (jj = 0; jj < NB; jj++) {
84297276fddSZach Atkins             for (ii = 0; ii < NB; ii++) {
84397276fddSZach Atkins               PetscInt idx = ii + jj * NB + kk * NB * NB;
84497276fddSZach Atkins               if (k + kk >= zs && j + jj >= ys && i + ii >= xs && k + kk < zs + zm && j + jj < ys + ym && i + ii < xs + xm) {
84597276fddSZach Atkins                 if (!OnBoundary(i + ii, j + jj, k + kk, mx, my, mz)) {
84697276fddSZach Atkins                   for (l = 0; l < 3; l++) q[k + kk][j + jj][i + ii][l] += eq[idx][l];
847c4762a1bSJed Brown                 }
848c4762a1bSJed Brown               }
849c4762a1bSJed Brown             }
850c4762a1bSJed Brown           }
851c4762a1bSJed Brown         }
852c4762a1bSJed Brown       }
853c4762a1bSJed Brown     }
854c4762a1bSJed Brown   }
8559566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xl, &x));
85697276fddSZach Atkins   PetscCall(DMDAVecRestoreArray(da, Ql, &q));
85797276fddSZach Atkins   PetscCall(VecZeroEntries(Q));
85897276fddSZach Atkins   PetscCall(DMLocalToGlobal(da, Ql, INSERT_VALUES, Q));
85997276fddSZach Atkins   PetscCall(DMRestoreLocalVector(da, &Ql));
8609566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xl));
8619566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
8623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
863c4762a1bSJed Brown }
864c4762a1bSJed Brown 
FormCoordinates(DM da,AppCtx * user)865d71ae5a4SJacob Faibussowitsch PetscErrorCode FormCoordinates(DM da, AppCtx *user)
866d71ae5a4SJacob Faibussowitsch {
867c4762a1bSJed Brown   Vec           coords;
868c4762a1bSJed Brown   DM            cda;
869c4762a1bSJed Brown   PetscInt      mx, my, mz;
870c4762a1bSJed Brown   PetscInt      i, j, k, xs, ys, zs, xm, ym, zm;
871c4762a1bSJed Brown   CoordField ***x;
872c4762a1bSJed Brown 
873c4762a1bSJed Brown   PetscFunctionBegin;
8749566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
8759566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(cda, &coords));
8769566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
8779566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
8789566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, coords, &x));
879c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
880c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
881c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
882f4f49eeaSPierre Jolivet         PetscReal cx  = ((PetscReal)i) / ((PetscReal)(mx - 1));
883f4f49eeaSPierre Jolivet         PetscReal cy  = ((PetscReal)j) / ((PetscReal)(my - 1));
884f4f49eeaSPierre Jolivet         PetscReal cz  = ((PetscReal)k) / ((PetscReal)(mz - 1));
885c4762a1bSJed Brown         PetscReal rad = user->rad + cy * user->height;
886c4762a1bSJed Brown         PetscReal ang = (cx - 0.5) * user->arc;
887c4762a1bSJed Brown         x[k][j][i][0] = rad * PetscSinReal(ang);
888c4762a1bSJed Brown         x[k][j][i][1] = rad * PetscCosReal(ang) - (user->rad + 0.5 * user->height) * PetscCosReal(-0.5 * user->arc);
889c4762a1bSJed Brown         x[k][j][i][2] = user->width * (cz - 0.5);
890c4762a1bSJed Brown       }
891c4762a1bSJed Brown     }
892c4762a1bSJed Brown   }
8939566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, coords, &x));
8949566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinates(da, coords));
8959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coords));
8963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
897c4762a1bSJed Brown }
898c4762a1bSJed Brown 
InitialGuess(DM da,AppCtx * user,Vec X)899d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialGuess(DM da, AppCtx *user, Vec X)
900d71ae5a4SJacob Faibussowitsch {
901c4762a1bSJed Brown   PetscInt i, j, k, xs, ys, zs, xm, ym, zm;
902c4762a1bSJed Brown   PetscInt mx, my, mz;
903c4762a1bSJed Brown   Field ***x;
904c4762a1bSJed Brown 
905c4762a1bSJed Brown   PetscFunctionBegin;
9069566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9079566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9089566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
909c4762a1bSJed Brown 
910c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
911c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
912c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
913c4762a1bSJed Brown         /* reference coordinates */
914f4f49eeaSPierre Jolivet         PetscReal p_x = ((PetscReal)i) / ((PetscReal)(mx - 1));
915f4f49eeaSPierre Jolivet         PetscReal p_y = ((PetscReal)j) / ((PetscReal)(my - 1));
916f4f49eeaSPierre Jolivet         PetscReal p_z = ((PetscReal)k) / ((PetscReal)(mz - 1));
917c4762a1bSJed Brown         PetscReal o_x = p_x;
918c4762a1bSJed Brown         PetscReal o_y = p_y;
919c4762a1bSJed Brown         PetscReal o_z = p_z;
920c4762a1bSJed Brown         x[k][j][i][0] = o_x - p_x;
921c4762a1bSJed Brown         x[k][j][i][1] = o_y - p_y;
922c4762a1bSJed Brown         x[k][j][i][2] = o_z - p_z;
923c4762a1bSJed Brown       }
924c4762a1bSJed Brown     }
925c4762a1bSJed Brown   }
9269566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
9273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
928c4762a1bSJed Brown }
929c4762a1bSJed Brown 
FormRHS(DM da,AppCtx * user,Vec X)930d71ae5a4SJacob Faibussowitsch PetscErrorCode FormRHS(DM da, AppCtx *user, Vec X)
931d71ae5a4SJacob Faibussowitsch {
932c4762a1bSJed Brown   PetscInt i, j, k, xs, ys, zs, xm, ym, zm;
933c4762a1bSJed Brown   PetscInt mx, my, mz;
934c4762a1bSJed Brown   Field ***x;
935c4762a1bSJed Brown 
936c4762a1bSJed Brown   PetscFunctionBegin;
9379566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9389566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9399566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
940c4762a1bSJed Brown 
941c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
942c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
943c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
944c4762a1bSJed Brown         x[k][j][i][0] = 0.;
945c4762a1bSJed Brown         x[k][j][i][1] = 0.;
946c4762a1bSJed Brown         x[k][j][i][2] = 0.;
947c4762a1bSJed Brown         if (i == (mx - 1) / 2 && j == (my - 1)) x[k][j][i][1] = user->ploading / (mz - 1);
948c4762a1bSJed Brown       }
949c4762a1bSJed Brown     }
950c4762a1bSJed Brown   }
9519566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
9523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
953c4762a1bSJed Brown }
954c4762a1bSJed Brown 
DisplayLine(SNES snes,Vec X)955d71ae5a4SJacob Faibussowitsch PetscErrorCode DisplayLine(SNES snes, Vec X)
956d71ae5a4SJacob Faibussowitsch {
957c4762a1bSJed Brown   PetscInt      r, i, j = 0, k = 0, xs, xm, ys, ym, zs, zm, mx, my, mz;
958c4762a1bSJed Brown   Field      ***x;
959c4762a1bSJed Brown   CoordField ***c;
960c4762a1bSJed Brown   DM            da, cda;
961c4762a1bSJed Brown   Vec           C;
962c4762a1bSJed Brown   PetscMPIInt   size, rank;
963c4762a1bSJed Brown 
964c4762a1bSJed Brown   PetscFunctionBegin;
9659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
9669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
9679566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
9689566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9699566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
9709566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &C));
971c4762a1bSJed Brown   j = my / 2;
972c4762a1bSJed Brown   k = mz / 2;
9739566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9749566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
9759566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
976c4762a1bSJed Brown   for (r = 0; r < size; r++) {
977c4762a1bSJed Brown     if (rank == r) {
978c4762a1bSJed Brown       if (j >= ys && j < ys + ym && k >= zs && k < zs + zm) {
979c4762a1bSJed Brown         for (i = xs; i < xs + xm; i++) {
98063a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %d %d: %f %f %f\n", i, 0, 0, (double)PetscRealPart(c[k][j][i][0] + x[k][j][i][0]), (double)PetscRealPart(c[k][j][i][1] + x[k][j][i][1]), (double)PetscRealPart(c[k][j][i][2] + x[k][j][i][2])));
981c4762a1bSJed Brown         }
982c4762a1bSJed Brown       }
983c4762a1bSJed Brown     }
9849566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
985c4762a1bSJed Brown   }
9869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
9879566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
9883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
989c4762a1bSJed Brown }
990c4762a1bSJed Brown 
991c4762a1bSJed Brown /*TEST
992c4762a1bSJed Brown 
993c4762a1bSJed Brown    test:
994c4762a1bSJed Brown       nsize: 2
99597276fddSZach Atkins       args: -da_refine 2 -pc_type mg -rad 10.0 -young 10. -ploading 0.0 -loading -1. -mg_levels_ksp_max_it 2 -snes_monitor_short -ksp_monitor_short
996c4762a1bSJed Brown       requires: !single
997c4762a1bSJed Brown       timeoutfactor: 3
998c4762a1bSJed Brown 
999c4762a1bSJed Brown    test:
1000c4762a1bSJed Brown       suffix: 2
1001a99ef635SJonas Heinzmann       args: -da_refine 2 -pc_type mg -rad 10.0 -young 10. -ploading 0.0 -loading -1. -mg_levels_ksp_max_it 2 -snes_monitor_short -ksp_monitor_short -npc_snes_type fas -npc_fas_levels_snes_type ncg -npc_fas_levels_snes_max_it 3 -npc_snes_monitor_short -npc_fas_levels_snes_linesearch_maxlambda 2.0
1002c4762a1bSJed Brown       requires: !single
1003c4762a1bSJed Brown 
1004c4762a1bSJed Brown    test:
1005c4762a1bSJed Brown       suffix: 3
100697276fddSZach Atkins       args: -da_refine 1 -da_overlap 3 -da_local_subdomains 4 -snes_type aspin -rad 10.0 -young 10. -ploading 0.0 -loading -0.5 -snes_monitor_short -ksp_monitor_short -npc_sub_snes_rtol 1e-2 -ksp_rtol 1e-2 -ksp_max_it 14 -snes_converged_reason -snes_max_linear_solve_fail 100 -npc_sub_ksp_type preonly -npc_sub_pc_type lu
100797276fddSZach Atkins       requires: !single
100897276fddSZach Atkins 
100997276fddSZach Atkins    test:
101097276fddSZach Atkins       suffix: 4
101197276fddSZach Atkins       args: -da_refine 2 -pc_type mg -rad 10.0 -young 10. -ploading -1. -loading -1. -mg_levels_ksp_max_it 2 -snes_monitor_short -ksp_monitor_short -snes_type newtonal -snes_newtonal_step_size 30 -ksp_rtol 1e-4 -info
101297276fddSZach Atkins       requires: !single defined(PETSC_USE_INFO)
101397276fddSZach Atkins       filter: grep -h -e "KSP Residual norm" -e "SNES Function norm" -e "Number of SNES iterations" -e "mu:"
101497276fddSZach Atkins 
101597276fddSZach Atkins    test:
101697276fddSZach Atkins       suffix: 5
101797276fddSZach Atkins       args: -da_refine 2 -pc_type mg -rad 10.0 -young 10. -ploading -1. -loading -1. -mg_levels_ksp_max_it 2 -snes_monitor_short -ksp_monitor_short -snes_type newtonal -snes_newtonal_step_size 30 -snes_newtonal_correction_type normal -ksp_rtol 1e-4
101897276fddSZach Atkins       requires: !single
101997276fddSZach Atkins 
102097276fddSZach Atkins    test:
102197276fddSZach Atkins       suffix: 6
102297276fddSZach Atkins       args: -da_refine 2 -pc_type mg -rad 10.0 -young 10. -ploading -0.5 -loading -1 -mg_levels_ksp_max_it 2 -snes_monitor_short -ksp_monitor_short -snes_type newtonal -snes_newtonal_step_size 0.5 -snes_newtonal_max_continuation_steps 3 -snes_newtonal_scale_rhs false -snes_newtonal_lambda_min -0.1 -ksp_rtol 1e-4
1023c4762a1bSJed Brown       requires: !single
1024c4762a1bSJed Brown 
1025c4762a1bSJed Brown TEST*/
1026