xref: /petsc/src/snes/tutorials/ex16.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Large-deformation Elasticity Buckling Example";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*F-----------------------------------------------------------------------
5c4762a1bSJed Brown 
6c4762a1bSJed Brown     This example solves the 3D large deformation elasticity problem
7c4762a1bSJed Brown 
8c4762a1bSJed Brown \begin{equation}
9c4762a1bSJed Brown  \int_{\Omega}F \cdot S : \nabla v d\Omega + \int_{\Omega} (loading)\mathbf{e}_y \cdot v d\Omega = 0
10c4762a1bSJed Brown \end{equation}
11c4762a1bSJed Brown 
12c4762a1bSJed Brown     F is the deformation gradient, and S is the second Piola-Kirchhoff tensor from the Saint Venant-Kirchhoff model of
13c4762a1bSJed Brown     hyperelasticity.  \Omega is a (arc) angle subsection of a cylindrical shell of thickness (height), inner radius
14c4762a1bSJed Brown     (rad) and width (width).  The problem is discretized using Q1 finite elements on a logically structured grid.
15c4762a1bSJed Brown     Homogenous Dirichlet boundary conditions are applied at the centers of the ends of the sphere.
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     This example is tunable with the following options:
18c4762a1bSJed Brown     -rad : the radius of the circle
19c4762a1bSJed Brown     -arc : set the angle of the arch represented
20c4762a1bSJed Brown     -loading : set the bulk loading (the mass)
21c4762a1bSJed Brown     -ploading : set the point loading at the top
22c4762a1bSJed Brown     -height : set the height of the arch
23c4762a1bSJed Brown     -width : set the width of the arch
24c4762a1bSJed Brown     -view_line : print initial and final offsets of the centerline of the
25c4762a1bSJed Brown                  beam along the x direction
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     The material properties may be modified using either:
28c4762a1bSJed Brown     -mu      : the first lame' parameter
29c4762a1bSJed Brown     -lambda  : the second lame' parameter
30c4762a1bSJed Brown 
31c4762a1bSJed Brown     Or:
32c4762a1bSJed Brown     -young   : the Young's modulus
33c4762a1bSJed Brown     -poisson : the poisson ratio
34c4762a1bSJed Brown 
35c4762a1bSJed Brown     This example is meant to show the strain placed upon the nonlinear solvers when trying to "snap through" the arch
36c4762a1bSJed Brown     using the loading.  Under certain parameter regimes, the arch will invert under the load, and the number of Newton
37c4762a1bSJed Brown     steps will jump considerably.  Composed nonlinear solvers may be used to mitigate this difficulty.
38c4762a1bSJed Brown 
39c4762a1bSJed Brown     The initial setup follows the example in pg. 268 of "Nonlinear Finite Element Methods" by Peter Wriggers, but is a
40c4762a1bSJed Brown     3D extension.
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   ------------------------------------------------------------------------F*/
43c4762a1bSJed Brown 
44c4762a1bSJed Brown #include <petscsnes.h>
45c4762a1bSJed Brown #include <petscdm.h>
46c4762a1bSJed Brown #include <petscdmda.h>
47c4762a1bSJed Brown 
48c4762a1bSJed Brown #define QP0 0.2113248654051871
49c4762a1bSJed Brown #define QP1 0.7886751345948129
50c4762a1bSJed Brown #define NQ  2
51c4762a1bSJed Brown #define NB  2
52c4762a1bSJed Brown #define NEB 8
53c4762a1bSJed Brown #define NEQ 8
54c4762a1bSJed Brown #define NPB 24
55c4762a1bSJed Brown 
56c4762a1bSJed Brown #define NVALS NEB *NEQ
57c4762a1bSJed Brown const PetscReal pts[NQ] = {QP0, QP1};
58c4762a1bSJed Brown const PetscReal wts[NQ] = {0.5, 0.5};
59c4762a1bSJed Brown 
60c4762a1bSJed Brown PetscScalar vals[NVALS];
61c4762a1bSJed Brown PetscScalar grad[3 * NVALS];
62c4762a1bSJed Brown 
63c4762a1bSJed Brown typedef PetscScalar Field[3];
64c4762a1bSJed Brown typedef PetscScalar CoordField[3];
65c4762a1bSJed Brown 
66c4762a1bSJed Brown typedef PetscScalar JacField[9];
67c4762a1bSJed Brown 
68c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *, Field ***, Field ***, void *);
69c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *, Field ***, Mat, Mat, void *);
70c4762a1bSJed Brown PetscErrorCode DisplayLine(SNES, Vec);
71c4762a1bSJed Brown PetscErrorCode FormElements();
72c4762a1bSJed Brown 
73c4762a1bSJed Brown typedef struct {
74c4762a1bSJed Brown   PetscReal loading;
75c4762a1bSJed Brown   PetscReal mu;
76c4762a1bSJed Brown   PetscReal lambda;
77c4762a1bSJed Brown   PetscReal rad;
78c4762a1bSJed Brown   PetscReal height;
79c4762a1bSJed Brown   PetscReal width;
80c4762a1bSJed Brown   PetscReal arc;
81c4762a1bSJed Brown   PetscReal ploading;
82c4762a1bSJed Brown } AppCtx;
83c4762a1bSJed Brown 
84c4762a1bSJed Brown PetscErrorCode        InitialGuess(DM, AppCtx *, Vec);
85c4762a1bSJed Brown PetscErrorCode        FormRHS(DM, AppCtx *, Vec);
86c4762a1bSJed Brown PetscErrorCode        FormCoordinates(DM, AppCtx *);
87c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES, Vec, Vec, void *);
88c4762a1bSJed Brown 
899371c9d4SSatish Balay int main(int argc, char **argv) {
90c4762a1bSJed Brown   AppCtx    user; /* user-defined work context */
91c4762a1bSJed Brown   PetscInt  mx, my, its;
92c4762a1bSJed Brown   MPI_Comm  comm;
93c4762a1bSJed Brown   SNES      snes;
94c4762a1bSJed Brown   DM        da;
95c4762a1bSJed Brown   Vec       x, X, b;
96c4762a1bSJed Brown   PetscBool youngflg, poissonflg, muflg, lambdaflg, view = PETSC_FALSE, viewline = PETSC_FALSE;
97c4762a1bSJed Brown   PetscReal poisson = 0.2, young = 4e4;
98c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN]     = "ex16.vts";
99c4762a1bSJed Brown   char      filename_def[PETSC_MAX_PATH_LEN] = "ex16_def.vts";
100c4762a1bSJed Brown 
101327415f7SBarry Smith   PetscFunctionBeginUser;
1029566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1039566063dSJacob Faibussowitsch   PetscCall(FormElements());
104c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1059566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &snes));
1069566063dSJacob 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));
1079566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1089566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1099566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, (DM)da));
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch   PetscCall(SNESSetNGS(snes, NonlinearGS, &user));
112c4762a1bSJed Brown 
1139566063dSJacob 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));
114c4762a1bSJed Brown   user.loading  = 0.0;
115c4762a1bSJed Brown   user.arc      = PETSC_PI / 3.;
116c4762a1bSJed Brown   user.mu       = 4.0;
117c4762a1bSJed Brown   user.lambda   = 1.0;
118c4762a1bSJed Brown   user.rad      = 100.0;
119c4762a1bSJed Brown   user.height   = 3.;
120c4762a1bSJed Brown   user.width    = 1.;
121c4762a1bSJed Brown   user.ploading = -5e3;
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-arc", &user.arc, NULL));
1249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, &muflg));
1259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lambda", &user.lambda, &lambdaflg));
1269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-rad", &user.rad, NULL));
1279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-height", &user.height, NULL));
1289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-width", &user.width, NULL));
1299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-loading", &user.loading, NULL));
1309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-ploading", &user.ploading, NULL));
1319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-poisson", &poisson, &poissonflg));
1329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-young", &young, &youngflg));
133c4762a1bSJed Brown   if ((youngflg || poissonflg) || !(muflg || lambdaflg)) {
134c4762a1bSJed Brown     /* set the lame' parameters based upon the poisson ratio and young's modulus */
135c4762a1bSJed Brown     user.lambda = poisson * young / ((1. + poisson) * (1. - 2. * poisson));
136c4762a1bSJed Brown     user.mu     = young / (2. * (1. + poisson));
137c4762a1bSJed Brown   }
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view", &view, NULL));
1399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_line", &viewline, NULL));
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "x_disp"));
1429566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "y_disp"));
1439566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 2, "z_disp"));
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
1469566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode(*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user));
1479566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobian)FormJacobianLocal, &user));
1489566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
1499566063dSJacob Faibussowitsch   PetscCall(FormCoordinates(da, &user));
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
1529566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &b));
1539566063dSJacob Faibussowitsch   PetscCall(InitialGuess(da, &user, x));
1549566063dSJacob Faibussowitsch   PetscCall(FormRHS(da, &user, b));
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "lambda: %f mu: %f\n", (double)user.lambda, (double)user.mu));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* show a cross-section of the initial state */
1591baa6e33SBarry Smith   if (viewline) PetscCall(DisplayLine(snes, x));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* get the loaded configuration */
1629566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, b, x));
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
16563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
1669566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &X));
167c4762a1bSJed Brown   /* show a cross-section of the final state */
1681baa6e33SBarry Smith   if (viewline) PetscCall(DisplayLine(snes, X));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   if (view) {
171c4762a1bSJed Brown     PetscViewer viewer;
172c4762a1bSJed Brown     Vec         coords;
1739566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD, filename, FILE_MODE_WRITE, &viewer));
1749566063dSJacob Faibussowitsch     PetscCall(VecView(x, viewer));
1759566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1769566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(da, &coords));
1779566063dSJacob Faibussowitsch     PetscCall(VecAXPY(coords, 1.0, x));
1789566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD, filename_def, FILE_MODE_WRITE, &viewer));
1799566063dSJacob Faibussowitsch     PetscCall(VecView(x, viewer));
1809566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown 
1839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1859566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1869566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1879566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
188b122ec5aSJacob Faibussowitsch   return 0;
189c4762a1bSJed Brown }
190c4762a1bSJed Brown 
1919371c9d4SSatish Balay PetscInt OnBoundary(PetscInt i, PetscInt j, PetscInt k, PetscInt mx, PetscInt my, PetscInt mz) {
192c4762a1bSJed Brown   if ((i == 0 || i == mx - 1) && j == my - 1) return 1;
193c4762a1bSJed Brown   return 0;
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
1969371c9d4SSatish Balay void BoundaryValue(PetscInt i, PetscInt j, PetscInt k, PetscInt mx, PetscInt my, PetscInt mz, PetscScalar *val, AppCtx *user) {
197c4762a1bSJed Brown   /* reference coordinates */
198c4762a1bSJed Brown   PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx - 1)));
199c4762a1bSJed Brown   PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my - 1)));
200c4762a1bSJed Brown   PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz - 1)));
201c4762a1bSJed Brown   PetscReal o_x = p_x;
202c4762a1bSJed Brown   PetscReal o_y = p_y;
203c4762a1bSJed Brown   PetscReal o_z = p_z;
204c4762a1bSJed Brown   val[0]        = o_x - p_x;
205c4762a1bSJed Brown   val[1]        = o_y - p_y;
206c4762a1bSJed Brown   val[2]        = o_z - p_z;
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
2099371c9d4SSatish Balay void InvertTensor(PetscScalar *t, PetscScalar *ti, PetscReal *dett) {
210c4762a1bSJed Brown   const PetscScalar a   = t[0];
211c4762a1bSJed Brown   const PetscScalar b   = t[1];
212c4762a1bSJed Brown   const PetscScalar c   = t[2];
213c4762a1bSJed Brown   const PetscScalar d   = t[3];
214c4762a1bSJed Brown   const PetscScalar e   = t[4];
215c4762a1bSJed Brown   const PetscScalar f   = t[5];
216c4762a1bSJed Brown   const PetscScalar g   = t[6];
217c4762a1bSJed Brown   const PetscScalar h   = t[7];
218c4762a1bSJed Brown   const PetscScalar i   = t[8];
219c4762a1bSJed Brown   const PetscReal   det = PetscRealPart(a * (e * i - f * h) - b * (i * d - f * g) + c * (d * h - e * g));
220c4762a1bSJed Brown   const PetscReal   di  = 1. / det;
221c4762a1bSJed Brown   if (ti) {
222c4762a1bSJed Brown     const PetscScalar A  = (e * i - f * h);
223c4762a1bSJed Brown     const PetscScalar B  = -(d * i - f * g);
224c4762a1bSJed Brown     const PetscScalar C  = (d * h - e * g);
225c4762a1bSJed Brown     const PetscScalar D  = -(b * i - c * h);
226c4762a1bSJed Brown     const PetscScalar E  = (a * i - c * g);
227c4762a1bSJed Brown     const PetscScalar F  = -(a * h - b * g);
228c4762a1bSJed Brown     const PetscScalar G  = (b * f - c * e);
229c4762a1bSJed Brown     const PetscScalar H  = -(a * f - c * d);
230c4762a1bSJed Brown     const PetscScalar II = (a * e - b * d);
231c4762a1bSJed Brown     ti[0]                = di * A;
232c4762a1bSJed Brown     ti[1]                = di * D;
233c4762a1bSJed Brown     ti[2]                = di * G;
234c4762a1bSJed Brown     ti[3]                = di * B;
235c4762a1bSJed Brown     ti[4]                = di * E;
236c4762a1bSJed Brown     ti[5]                = di * H;
237c4762a1bSJed Brown     ti[6]                = di * C;
238c4762a1bSJed Brown     ti[7]                = di * F;
239c4762a1bSJed Brown     ti[8]                = di * II;
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown   if (dett) *dett = det;
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
2449371c9d4SSatish Balay void TensorTensor(PetscScalar *a, PetscScalar *b, PetscScalar *c) {
245c4762a1bSJed Brown   PetscInt i, j, m;
246c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
247c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
248c4762a1bSJed Brown       c[i + 3 * j] = 0;
2499371c9d4SSatish Balay       for (m = 0; m < 3; m++) c[i + 3 * j] += a[m + 3 * j] * b[i + 3 * m];
250c4762a1bSJed Brown     }
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown }
253c4762a1bSJed Brown 
2549371c9d4SSatish Balay void TensorTransposeTensor(PetscScalar *a, PetscScalar *b, PetscScalar *c) {
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[3 * m + j] * b[i + 3 * m];
260c4762a1bSJed Brown     }
261c4762a1bSJed Brown   }
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
2649371c9d4SSatish Balay void TensorVector(PetscScalar *rot, PetscScalar *vec, PetscScalar *tvec) {
265c4762a1bSJed Brown   tvec[0] = rot[0] * vec[0] + rot[1] * vec[1] + rot[2] * vec[2];
266c4762a1bSJed Brown   tvec[1] = rot[3] * vec[0] + rot[4] * vec[1] + rot[5] * vec[2];
267c4762a1bSJed Brown   tvec[2] = rot[6] * vec[0] + rot[7] * vec[1] + rot[8] * vec[2];
268c4762a1bSJed Brown }
269c4762a1bSJed Brown 
2709371c9d4SSatish Balay void DeformationGradient(Field *ex, PetscInt qi, PetscInt qj, PetscInt qk, PetscScalar *invJ, PetscScalar *F) {
271c4762a1bSJed Brown   PetscInt ii, jj, kk, l;
2729371c9d4SSatish Balay   for (l = 0; l < 9; l++) { F[l] = 0.; }
273c4762a1bSJed Brown   F[0] = 1.;
274c4762a1bSJed Brown   F[4] = 1.;
275c4762a1bSJed Brown   F[8] = 1.;
276c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
277c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
278c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
279c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
280c4762a1bSJed Brown         PetscInt    idx  = ii + jj * NB + kk * NB * NB;
281c4762a1bSJed Brown         PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
282c4762a1bSJed Brown         PetscScalar lgrad[3];
283c4762a1bSJed Brown         TensorVector(invJ, &grad[3 * bidx], lgrad);
2849371c9d4SSatish Balay         F[0] += lgrad[0] * ex[idx][0];
2859371c9d4SSatish Balay         F[1] += lgrad[1] * ex[idx][0];
2869371c9d4SSatish Balay         F[2] += lgrad[2] * ex[idx][0];
2879371c9d4SSatish Balay         F[3] += lgrad[0] * ex[idx][1];
2889371c9d4SSatish Balay         F[4] += lgrad[1] * ex[idx][1];
2899371c9d4SSatish Balay         F[5] += lgrad[2] * ex[idx][1];
2909371c9d4SSatish Balay         F[6] += lgrad[0] * ex[idx][2];
2919371c9d4SSatish Balay         F[7] += lgrad[1] * ex[idx][2];
2929371c9d4SSatish Balay         F[8] += lgrad[2] * ex[idx][2];
293c4762a1bSJed Brown       }
294c4762a1bSJed Brown     }
295c4762a1bSJed Brown   }
296c4762a1bSJed Brown }
297c4762a1bSJed Brown 
2989371c9d4SSatish Balay void DeformationGradientJacobian(PetscInt qi, PetscInt qj, PetscInt qk, PetscInt ii, PetscInt jj, PetscInt kk, PetscInt fld, PetscScalar *invJ, PetscScalar *dF) {
299c4762a1bSJed Brown   PetscInt    l;
300c4762a1bSJed Brown   PetscScalar lgrad[3];
301c4762a1bSJed Brown   PetscInt    idx  = ii + jj * NB + kk * NB * NB;
302c4762a1bSJed Brown   PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
3039371c9d4SSatish Balay   for (l = 0; l < 9; l++) { dF[l] = 0.; }
304c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
305c4762a1bSJed Brown   TensorVector(invJ, &grad[3 * bidx], lgrad);
3069371c9d4SSatish Balay   dF[3 * fld]     = lgrad[0];
3079371c9d4SSatish Balay   dF[3 * fld + 1] = lgrad[1];
3089371c9d4SSatish Balay   dF[3 * fld + 2] = lgrad[2];
309c4762a1bSJed Brown }
310c4762a1bSJed Brown 
3119371c9d4SSatish Balay void LagrangeGreenStrain(PetscScalar *F, PetscScalar *E) {
312c4762a1bSJed Brown   PetscInt i, j, m;
313c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
314c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
315c4762a1bSJed Brown       E[i + 3 * j] = 0;
3169371c9d4SSatish Balay       for (m = 0; m < 3; m++) E[i + 3 * j] += 0.5 * F[3 * m + j] * F[i + 3 * m];
317c4762a1bSJed Brown     }
318c4762a1bSJed Brown   }
3199371c9d4SSatish Balay   for (i = 0; i < 3; i++) { E[i + 3 * i] -= 0.5; }
320c4762a1bSJed Brown }
321c4762a1bSJed Brown 
3229371c9d4SSatish Balay void SaintVenantKirchoff(PetscReal lambda, PetscReal mu, PetscScalar *F, PetscScalar *S) {
323c4762a1bSJed Brown   PetscInt    i, j;
324c4762a1bSJed Brown   PetscScalar E[9];
325c4762a1bSJed Brown   PetscScalar trE = 0;
326c4762a1bSJed Brown   LagrangeGreenStrain(F, E);
3279371c9d4SSatish Balay   for (i = 0; i < 3; i++) { trE += E[i + 3 * i]; }
328c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
329c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
330c4762a1bSJed Brown       S[i + 3 * j] = 2. * mu * E[i + 3 * j];
3319371c9d4SSatish Balay       if (i == j) { S[i + 3 * j] += trE * lambda; }
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown   }
334c4762a1bSJed Brown }
335c4762a1bSJed Brown 
3369371c9d4SSatish Balay void SaintVenantKirchoffJacobian(PetscReal lambda, PetscReal mu, PetscScalar *F, PetscScalar *dF, PetscScalar *dS) {
337c4762a1bSJed Brown   PetscScalar FtdF[9], dE[9];
338c4762a1bSJed Brown   PetscInt    i, j;
339c4762a1bSJed Brown   PetscScalar dtrE = 0.;
340c4762a1bSJed Brown   TensorTransposeTensor(dF, F, dE);
341c4762a1bSJed Brown   TensorTransposeTensor(F, dF, FtdF);
342c4762a1bSJed Brown   for (i = 0; i < 9; i++) dE[i] += FtdF[i];
343c4762a1bSJed Brown   for (i = 0; i < 9; i++) dE[i] *= 0.5;
3449371c9d4SSatish Balay   for (i = 0; i < 3; i++) { dtrE += dE[i + 3 * i]; }
345c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
346c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
347c4762a1bSJed Brown       dS[i + 3 * j] = 2. * mu * dE[i + 3 * j];
3489371c9d4SSatish Balay       if (i == j) { dS[i + 3 * j] += dtrE * lambda; }
349c4762a1bSJed Brown     }
350c4762a1bSJed Brown   }
351c4762a1bSJed Brown }
352c4762a1bSJed Brown 
3539371c9d4SSatish Balay PetscErrorCode FormElements() {
354c4762a1bSJed Brown   PetscInt  i, j, k, ii, jj, kk;
355c4762a1bSJed Brown   PetscReal bx, by, bz, dbx, dby, dbz;
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   PetscFunctionBegin;
358c4762a1bSJed Brown   /* construct the basis function values and derivatives */
359c4762a1bSJed Brown   for (k = 0; k < NB; k++) {
360c4762a1bSJed Brown     for (j = 0; j < NB; j++) {
361c4762a1bSJed Brown       for (i = 0; i < NB; i++) {
362c4762a1bSJed Brown         /* loop over the quadrature points */
363c4762a1bSJed Brown         for (kk = 0; kk < NQ; kk++) {
364c4762a1bSJed Brown           for (jj = 0; jj < NQ; jj++) {
365c4762a1bSJed Brown             for (ii = 0; ii < NQ; ii++) {
366c4762a1bSJed Brown               PetscInt idx = ii + NQ * jj + NQ * NQ * kk + NEQ * i + NEQ * NB * j + NEQ * NB * NB * k;
367c4762a1bSJed Brown               bx           = pts[ii];
368c4762a1bSJed Brown               by           = pts[jj];
369c4762a1bSJed Brown               bz           = pts[kk];
370c4762a1bSJed Brown               dbx          = 1.;
371c4762a1bSJed Brown               dby          = 1.;
372c4762a1bSJed Brown               dbz          = 1.;
3739371c9d4SSatish Balay               if (i == 0) {
3749371c9d4SSatish Balay                 bx  = 1. - bx;
3759371c9d4SSatish Balay                 dbx = -1;
3769371c9d4SSatish Balay               }
3779371c9d4SSatish Balay               if (j == 0) {
3789371c9d4SSatish Balay                 by  = 1. - by;
3799371c9d4SSatish Balay                 dby = -1;
3809371c9d4SSatish Balay               }
3819371c9d4SSatish Balay               if (k == 0) {
3829371c9d4SSatish Balay                 bz  = 1. - bz;
3839371c9d4SSatish Balay                 dbz = -1;
3849371c9d4SSatish Balay               }
385c4762a1bSJed Brown               vals[idx]         = bx * by * bz;
386c4762a1bSJed Brown               grad[3 * idx + 0] = dbx * by * bz;
387c4762a1bSJed Brown               grad[3 * idx + 1] = dby * bx * bz;
388c4762a1bSJed Brown               grad[3 * idx + 2] = dbz * bx * by;
389c4762a1bSJed Brown             }
390c4762a1bSJed Brown           }
391c4762a1bSJed Brown         }
392c4762a1bSJed Brown       }
393c4762a1bSJed Brown     }
394c4762a1bSJed Brown   }
395c4762a1bSJed Brown   PetscFunctionReturn(0);
396c4762a1bSJed Brown }
397c4762a1bSJed Brown 
3989371c9d4SSatish Balay void GatherElementData(PetscInt mx, PetscInt my, PetscInt mz, Field ***x, CoordField ***c, PetscInt i, PetscInt j, PetscInt k, Field *ex, CoordField *ec, AppCtx *user) {
399c4762a1bSJed Brown   PetscInt m;
400c4762a1bSJed Brown   PetscInt ii, jj, kk;
401c4762a1bSJed Brown   /* gather the data -- loop over element unknowns */
402c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
403c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
404c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
405c4762a1bSJed Brown         PetscInt idx = ii + jj * NB + kk * NB * NB;
406c4762a1bSJed Brown         /* decouple the boundary nodes for the displacement variables */
407c4762a1bSJed Brown         if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz)) {
408c4762a1bSJed Brown           BoundaryValue(i + ii, j + jj, k + kk, mx, my, mz, ex[idx], user);
409c4762a1bSJed Brown         } else {
4109371c9d4SSatish Balay           for (m = 0; m < 3; m++) { ex[idx][m] = x[k + kk][j + jj][i + ii][m]; }
411c4762a1bSJed Brown         }
4129371c9d4SSatish Balay         for (m = 0; m < 3; m++) { ec[idx][m] = c[k + kk][j + jj][i + ii][m]; }
413c4762a1bSJed Brown       }
414c4762a1bSJed Brown     }
415c4762a1bSJed Brown   }
416c4762a1bSJed Brown }
417c4762a1bSJed Brown 
4189371c9d4SSatish Balay void QuadraturePointGeometricJacobian(CoordField *ec, PetscInt qi, PetscInt qj, PetscInt qk, PetscScalar *J) {
419c4762a1bSJed Brown   PetscInt ii, jj, kk;
420c4762a1bSJed Brown   /* construct the gradient at the given quadrature point named by i,j,k */
4219371c9d4SSatish Balay   for (ii = 0; ii < 9; ii++) { J[ii] = 0; }
422c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
423c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
424c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
425c4762a1bSJed Brown         PetscInt idx  = ii + jj * NB + kk * NB * NB;
426c4762a1bSJed Brown         PetscInt bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
4279371c9d4SSatish Balay         J[0] += grad[3 * bidx + 0] * ec[idx][0];
4289371c9d4SSatish Balay         J[1] += grad[3 * bidx + 1] * ec[idx][0];
4299371c9d4SSatish Balay         J[2] += grad[3 * bidx + 2] * ec[idx][0];
4309371c9d4SSatish Balay         J[3] += grad[3 * bidx + 0] * ec[idx][1];
4319371c9d4SSatish Balay         J[4] += grad[3 * bidx + 1] * ec[idx][1];
4329371c9d4SSatish Balay         J[5] += grad[3 * bidx + 2] * ec[idx][1];
4339371c9d4SSatish Balay         J[6] += grad[3 * bidx + 0] * ec[idx][2];
4349371c9d4SSatish Balay         J[7] += grad[3 * bidx + 1] * ec[idx][2];
4359371c9d4SSatish Balay         J[8] += grad[3 * bidx + 2] * ec[idx][2];
436c4762a1bSJed Brown       }
437c4762a1bSJed Brown     }
438c4762a1bSJed Brown   }
439c4762a1bSJed Brown }
440c4762a1bSJed Brown 
4419371c9d4SSatish Balay void FormElementJacobian(Field *ex, CoordField *ec, Field *ef, PetscScalar *ej, AppCtx *user) {
442c4762a1bSJed Brown   PetscReal   vol;
443c4762a1bSJed Brown   PetscScalar J[9];
444c4762a1bSJed Brown   PetscScalar invJ[9];
445c4762a1bSJed Brown   PetscScalar F[9], S[9], dF[9], dS[9], dFS[9], FdS[9], FS[9];
446c4762a1bSJed Brown   PetscReal   scl;
447c4762a1bSJed Brown   PetscInt    i, j, k, l, ii, jj, kk, ll, qi, qj, qk, m;
448c4762a1bSJed Brown 
4499371c9d4SSatish Balay   if (ej)
4509371c9d4SSatish Balay     for (i = 0; i < NPB * NPB; i++) ej[i] = 0.;
4519371c9d4SSatish Balay   if (ef)
4529371c9d4SSatish Balay     for (i = 0; i < NEB; i++) {
4539371c9d4SSatish Balay       ef[i][0] = 0.;
4549371c9d4SSatish Balay       ef[i][1] = 0.;
4559371c9d4SSatish Balay       ef[i][2] = 0.;
4569371c9d4SSatish Balay     }
457c4762a1bSJed Brown   /* loop over quadrature */
458c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
459c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
460c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
461c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec, qi, qj, qk, J);
462c4762a1bSJed Brown         InvertTensor(J, invJ, &vol);
463c4762a1bSJed Brown         scl = vol * wts[qi] * wts[qj] * wts[qk];
464c4762a1bSJed Brown         DeformationGradient(ex, qi, qj, qk, invJ, F);
465c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda, user->mu, F, S);
466c4762a1bSJed Brown         /* form the function */
467c4762a1bSJed Brown         if (ef) {
468c4762a1bSJed Brown           TensorTensor(F, S, FS);
469c4762a1bSJed Brown           for (kk = 0; kk < NB; kk++) {
470c4762a1bSJed Brown             for (jj = 0; jj < NB; jj++) {
471c4762a1bSJed Brown               for (ii = 0; ii < NB; ii++) {
472c4762a1bSJed Brown                 PetscInt    idx  = ii + jj * NB + kk * NB * NB;
473c4762a1bSJed Brown                 PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
474c4762a1bSJed Brown                 PetscScalar lgrad[3];
475c4762a1bSJed Brown                 TensorVector(invJ, &grad[3 * bidx], lgrad);
476c4762a1bSJed Brown                 /* mu*F : grad phi_{u,v,w} */
4779371c9d4SSatish Balay                 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]); }
478c4762a1bSJed Brown                 ef[idx][1] -= scl * user->loading * vals[bidx];
479c4762a1bSJed Brown               }
480c4762a1bSJed Brown             }
481c4762a1bSJed Brown           }
482c4762a1bSJed Brown         }
483c4762a1bSJed Brown         /* form the jacobian */
484c4762a1bSJed Brown         if (ej) {
485c4762a1bSJed Brown           /* loop over trialfunctions */
486c4762a1bSJed Brown           for (k = 0; k < NB; k++) {
487c4762a1bSJed Brown             for (j = 0; j < NB; j++) {
488c4762a1bSJed Brown               for (i = 0; i < NB; i++) {
489c4762a1bSJed Brown                 for (l = 0; l < 3; l++) {
490c4762a1bSJed Brown                   PetscInt tridx = l + 3 * (i + j * NB + k * NB * NB);
491c4762a1bSJed Brown                   DeformationGradientJacobian(qi, qj, qk, i, j, k, l, invJ, dF);
492c4762a1bSJed Brown                   SaintVenantKirchoffJacobian(user->lambda, user->mu, F, dF, dS);
493c4762a1bSJed Brown                   TensorTensor(dF, S, dFS);
494c4762a1bSJed Brown                   TensorTensor(F, dS, FdS);
495c4762a1bSJed Brown                   for (m = 0; m < 9; m++) dFS[m] += FdS[m];
496c4762a1bSJed Brown                   /* loop over testfunctions */
497c4762a1bSJed Brown                   for (kk = 0; kk < NB; kk++) {
498c4762a1bSJed Brown                     for (jj = 0; jj < NB; jj++) {
499c4762a1bSJed Brown                       for (ii = 0; ii < NB; ii++) {
500c4762a1bSJed Brown                         PetscInt    idx  = ii + jj * NB + kk * NB * NB;
501c4762a1bSJed Brown                         PetscInt    bidx = 8 * idx + qi + NQ * qj + NQ * NQ * qk;
502c4762a1bSJed Brown                         PetscScalar lgrad[3];
503c4762a1bSJed Brown                         TensorVector(invJ, &grad[3 * bidx], lgrad);
504c4762a1bSJed Brown                         for (ll = 0; ll < 3; ll++) {
505c4762a1bSJed Brown                           PetscInt teidx = ll + 3 * (ii + jj * NB + kk * NB * NB);
5069371c9d4SSatish Balay                           ej[teidx + NPB * tridx] += scl * (lgrad[0] * dFS[3 * ll + 0] + lgrad[1] * dFS[3 * ll + 1] + lgrad[2] * dFS[3 * ll + 2]);
507c4762a1bSJed Brown                         }
508c4762a1bSJed Brown                       }
509c4762a1bSJed Brown                     }
510c4762a1bSJed Brown                   } /* end of testfunctions */
511c4762a1bSJed Brown                 }
512c4762a1bSJed Brown               }
513c4762a1bSJed Brown             }
514c4762a1bSJed Brown           } /* end of trialfunctions */
515c4762a1bSJed Brown         }
516c4762a1bSJed Brown       }
517c4762a1bSJed Brown     }
518c4762a1bSJed Brown   } /* end of quadrature points */
519c4762a1bSJed Brown }
520c4762a1bSJed Brown 
5219371c9d4SSatish Balay void FormPBJacobian(PetscInt i, PetscInt j, PetscInt k, Field *ex, CoordField *ec, Field *ef, PetscScalar *ej, AppCtx *user) {
522c4762a1bSJed Brown   PetscReal   vol;
523c4762a1bSJed Brown   PetscScalar J[9];
524c4762a1bSJed Brown   PetscScalar invJ[9];
525c4762a1bSJed Brown   PetscScalar F[9], S[9], dF[9], dS[9], dFS[9], FdS[9], FS[9];
526c4762a1bSJed Brown   PetscReal   scl;
527c4762a1bSJed Brown   PetscInt    l, ll, qi, qj, qk, m;
528c4762a1bSJed Brown   PetscInt    idx = i + j * NB + k * NB * NB;
529c4762a1bSJed Brown   PetscScalar lgrad[3];
530c4762a1bSJed Brown 
5319371c9d4SSatish Balay   if (ej)
5329371c9d4SSatish Balay     for (l = 0; l < 9; l++) ej[l] = 0.;
5339371c9d4SSatish Balay   if (ef)
5349371c9d4SSatish Balay     for (l = 0; l < 1; l++) {
5359371c9d4SSatish Balay       ef[l][0] = 0.;
5369371c9d4SSatish Balay       ef[l][1] = 0.;
5379371c9d4SSatish Balay       ef[l][2] = 0.;
5389371c9d4SSatish Balay     }
539c4762a1bSJed Brown   /* loop over quadrature */
540c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
541c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
542c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
543c4762a1bSJed Brown         PetscInt bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
544c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec, qi, qj, qk, J);
545c4762a1bSJed Brown         InvertTensor(J, invJ, &vol);
546c4762a1bSJed Brown         TensorVector(invJ, &grad[3 * bidx], lgrad);
547c4762a1bSJed Brown         scl = vol * wts[qi] * wts[qj] * wts[qk];
548c4762a1bSJed Brown         DeformationGradient(ex, qi, qj, qk, invJ, F);
549c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda, user->mu, F, S);
550c4762a1bSJed Brown         /* form the function */
551c4762a1bSJed Brown         if (ef) {
552c4762a1bSJed Brown           TensorTensor(F, S, FS);
5539371c9d4SSatish Balay           for (m = 0; m < 3; m++) { ef[0][m] += scl * (lgrad[0] * FS[3 * m + 0] + lgrad[1] * FS[3 * m + 1] + lgrad[2] * FS[3 * m + 2]); }
554c4762a1bSJed Brown           ef[0][1] -= scl * user->loading * vals[bidx];
555c4762a1bSJed Brown         }
556c4762a1bSJed Brown         /* form the jacobian */
557c4762a1bSJed Brown         if (ej) {
558c4762a1bSJed Brown           for (l = 0; l < 3; l++) {
559c4762a1bSJed Brown             DeformationGradientJacobian(qi, qj, qk, i, j, k, l, invJ, dF);
560c4762a1bSJed Brown             SaintVenantKirchoffJacobian(user->lambda, user->mu, F, dF, dS);
561c4762a1bSJed Brown             TensorTensor(dF, S, dFS);
562c4762a1bSJed Brown             TensorTensor(F, dS, FdS);
563c4762a1bSJed Brown             for (m = 0; m < 9; m++) dFS[m] += FdS[m];
5649371c9d4SSatish Balay             for (ll = 0; ll < 3; ll++) { ej[ll + 3 * l] += scl * (lgrad[0] * dFS[3 * ll + 0] + lgrad[1] * dFS[3 * ll + 1] + lgrad[2] * dFS[3 * ll + 2]); }
565c4762a1bSJed Brown           }
566c4762a1bSJed Brown         }
567c4762a1bSJed Brown       }
568c4762a1bSJed Brown     } /* end of quadrature points */
569c4762a1bSJed Brown   }
570c4762a1bSJed Brown }
571c4762a1bSJed Brown 
5729371c9d4SSatish Balay void ApplyBCsElement(PetscInt mx, PetscInt my, PetscInt mz, PetscInt i, PetscInt j, PetscInt k, PetscScalar *jacobian) {
573c4762a1bSJed Brown   PetscInt ii, jj, kk, ll, ei, ej, ek, el;
574c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
575c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
576c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
577c4762a1bSJed Brown         for (ll = 0; ll < 3; ll++) {
578c4762a1bSJed Brown           PetscInt tridx = ll + 3 * (ii + jj * NB + kk * NB * NB);
579c4762a1bSJed Brown           for (ek = 0; ek < NB; ek++) {
580c4762a1bSJed Brown             for (ej = 0; ej < NB; ej++) {
581c4762a1bSJed Brown               for (ei = 0; ei < NB; ei++) {
582c4762a1bSJed Brown                 for (el = 0; el < 3; el++) {
583c4762a1bSJed Brown                   if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz) || OnBoundary(i + ei, j + ej, k + ek, mx, my, mz)) {
584c4762a1bSJed Brown                     PetscInt teidx = el + 3 * (ei + ej * NB + ek * NB * NB);
585c4762a1bSJed Brown                     if (teidx == tridx) {
586c4762a1bSJed Brown                       jacobian[tridx + NPB * teidx] = 1.;
587c4762a1bSJed Brown                     } else {
588c4762a1bSJed Brown                       jacobian[tridx + NPB * teidx] = 0.;
589c4762a1bSJed Brown                     }
590c4762a1bSJed Brown                   }
591c4762a1bSJed Brown                 }
592c4762a1bSJed Brown               }
593c4762a1bSJed Brown             }
594c4762a1bSJed Brown           }
595c4762a1bSJed Brown         }
596c4762a1bSJed Brown       }
597c4762a1bSJed Brown     }
598c4762a1bSJed Brown   }
599c4762a1bSJed Brown }
600c4762a1bSJed Brown 
6019371c9d4SSatish Balay PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, Field ***x, Mat jacpre, Mat jac, void *ptr) {
602c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
603c4762a1bSJed Brown   AppCtx     *user = (AppCtx *)ptr;
604c4762a1bSJed Brown   PetscInt    i, j, k, m, l;
605c4762a1bSJed Brown   PetscInt    ii, jj, kk;
606c4762a1bSJed Brown   PetscScalar ej[NPB * NPB];
607c4762a1bSJed Brown   PetscScalar vals[NPB * NPB];
608c4762a1bSJed Brown   Field       ex[NEB];
609c4762a1bSJed Brown   CoordField  ec[NEB];
610c4762a1bSJed Brown 
611c4762a1bSJed Brown   PetscInt      xs = info->xs, ys = info->ys, zs = info->zs;
612c4762a1bSJed Brown   PetscInt      xm = info->xm, ym = info->ym, zm = info->zm;
613c4762a1bSJed Brown   PetscInt      xes, yes, zes, xee, yee, zee;
614c4762a1bSJed Brown   PetscInt      mx = info->mx, my = info->my, mz = info->mz;
615c4762a1bSJed Brown   DM            cda;
616c4762a1bSJed Brown   CoordField ***c;
617c4762a1bSJed Brown   Vec           C;
618c4762a1bSJed Brown   PetscInt      nrows;
619c4762a1bSJed Brown   MatStencil    col[NPB], row[NPB];
620c4762a1bSJed Brown   PetscScalar   v[9];
621c4762a1bSJed Brown 
622c4762a1bSJed Brown   PetscFunctionBegin;
6239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(info->da, &cda));
6249566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(info->da, &C));
6259566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
6269566063dSJacob Faibussowitsch   PetscCall(MatScale(jac, 0.0));
627c4762a1bSJed Brown 
628c4762a1bSJed Brown   xes = xs;
629c4762a1bSJed Brown   yes = ys;
630c4762a1bSJed Brown   zes = zs;
631c4762a1bSJed Brown   xee = xs + xm;
632c4762a1bSJed Brown   yee = ys + ym;
633c4762a1bSJed Brown   zee = zs + zm;
634c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
635c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
636c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
637c4762a1bSJed Brown   if (xs + xm == mx) xee = xs + xm - 1;
638c4762a1bSJed Brown   if (ys + ym == my) yee = ys + ym - 1;
639c4762a1bSJed Brown   if (zs + zm == mz) zee = zs + zm - 1;
640c4762a1bSJed Brown   for (k = zes; k < zee; k++) {
641c4762a1bSJed Brown     for (j = yes; j < yee; j++) {
642c4762a1bSJed Brown       for (i = xes; i < xee; i++) {
643c4762a1bSJed Brown         GatherElementData(mx, my, mz, x, c, i, j, k, ex, ec, user);
644c4762a1bSJed Brown         FormElementJacobian(ex, ec, NULL, ej, user);
645c4762a1bSJed Brown         ApplyBCsElement(mx, my, mz, i, j, k, ej);
646c4762a1bSJed Brown         nrows = 0.;
647c4762a1bSJed Brown         for (kk = 0; kk < NB; kk++) {
648c4762a1bSJed Brown           for (jj = 0; jj < NB; jj++) {
649c4762a1bSJed Brown             for (ii = 0; ii < NB; ii++) {
650c4762a1bSJed Brown               PetscInt idx = ii + jj * 2 + kk * 4;
651c4762a1bSJed Brown               for (m = 0; m < 3; m++) {
652c4762a1bSJed Brown                 col[3 * idx + m].i = i + ii;
653c4762a1bSJed Brown                 col[3 * idx + m].j = j + jj;
654c4762a1bSJed Brown                 col[3 * idx + m].k = k + kk;
655c4762a1bSJed Brown                 col[3 * idx + m].c = m;
656c4762a1bSJed Brown                 if (i + ii >= xs && i + ii < xm + xs && j + jj >= ys && j + jj < ys + ym && k + kk >= zs && k + kk < zs + zm) {
657c4762a1bSJed Brown                   row[nrows].i = i + ii;
658c4762a1bSJed Brown                   row[nrows].j = j + jj;
659c4762a1bSJed Brown                   row[nrows].k = k + kk;
660c4762a1bSJed Brown                   row[nrows].c = m;
661c4762a1bSJed Brown                   for (l = 0; l < NPB; l++) vals[NPB * nrows + l] = ej[NPB * (3 * idx + m) + l];
662c4762a1bSJed Brown                   nrows++;
663c4762a1bSJed Brown                 }
664c4762a1bSJed Brown               }
665c4762a1bSJed Brown             }
666c4762a1bSJed Brown           }
667c4762a1bSJed Brown         }
6689566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jac, nrows, row, NPB, col, vals, ADD_VALUES));
669c4762a1bSJed Brown       }
670c4762a1bSJed Brown     }
671c4762a1bSJed Brown   }
672c4762a1bSJed Brown 
6739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FLUSH_ASSEMBLY));
6749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FLUSH_ASSEMBLY));
675c4762a1bSJed Brown 
676c4762a1bSJed Brown   /* set the diagonal */
6779371c9d4SSatish Balay   v[0] = 1.;
6789371c9d4SSatish Balay   v[1] = 0.;
6799371c9d4SSatish Balay   v[2] = 0.;
6809371c9d4SSatish Balay   v[3] = 0.;
6819371c9d4SSatish Balay   v[4] = 1.;
6829371c9d4SSatish Balay   v[5] = 0.;
6839371c9d4SSatish Balay   v[6] = 0.;
6849371c9d4SSatish Balay   v[7] = 0.;
6859371c9d4SSatish Balay   v[8] = 1.;
686c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
687c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
688c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
689c4762a1bSJed Brown         if (OnBoundary(i, j, k, mx, my, mz)) {
690c4762a1bSJed Brown           for (m = 0; m < 3; m++) {
691c4762a1bSJed Brown             col[m].i = i;
692c4762a1bSJed Brown             col[m].j = j;
693c4762a1bSJed Brown             col[m].k = k;
694c4762a1bSJed Brown             col[m].c = m;
695c4762a1bSJed Brown           }
6969566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(jac, 3, col, 3, col, v, INSERT_VALUES));
697c4762a1bSJed Brown         }
698c4762a1bSJed Brown       }
699c4762a1bSJed Brown     }
700c4762a1bSJed Brown   }
701c4762a1bSJed Brown 
7029566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
7039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
704c4762a1bSJed Brown 
7059566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
706c4762a1bSJed Brown   PetscFunctionReturn(0);
707c4762a1bSJed Brown }
708c4762a1bSJed Brown 
7099371c9d4SSatish Balay PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field ***x, Field ***f, void *ptr) {
710c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
711c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
712c4762a1bSJed Brown   PetscInt i, j, k, l;
713c4762a1bSJed Brown   PetscInt ii, jj, kk;
714c4762a1bSJed Brown 
715c4762a1bSJed Brown   Field      ef[NEB];
716c4762a1bSJed Brown   Field      ex[NEB];
717c4762a1bSJed Brown   CoordField ec[NEB];
718c4762a1bSJed Brown 
719c4762a1bSJed Brown   PetscInt      xs = info->xs, ys = info->ys, zs = info->zs;
720c4762a1bSJed Brown   PetscInt      xm = info->xm, ym = info->ym, zm = info->zm;
721c4762a1bSJed Brown   PetscInt      xes, yes, zes, xee, yee, zee;
722c4762a1bSJed Brown   PetscInt      mx = info->mx, my = info->my, mz = info->mz;
723c4762a1bSJed Brown   DM            cda;
724c4762a1bSJed Brown   CoordField ***c;
725c4762a1bSJed Brown   Vec           C;
726c4762a1bSJed Brown 
727c4762a1bSJed Brown   PetscFunctionBegin;
7289566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(info->da, &cda));
7299566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(info->da, &C));
7309566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
7319566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(info->da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
7329566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(info->da, &xs, &ys, &zs, &xm, &ym, &zm));
733c4762a1bSJed Brown 
734c4762a1bSJed Brown   /* loop over elements */
735c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
736c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
737c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
7389371c9d4SSatish Balay         for (l = 0; l < 3; l++) { f[k][j][i][l] = 0.; }
739c4762a1bSJed Brown       }
740c4762a1bSJed Brown     }
741c4762a1bSJed Brown   }
742c4762a1bSJed Brown   /* element starts and ends */
743c4762a1bSJed Brown   xes = xs;
744c4762a1bSJed Brown   yes = ys;
745c4762a1bSJed Brown   zes = zs;
746c4762a1bSJed Brown   xee = xs + xm;
747c4762a1bSJed Brown   yee = ys + ym;
748c4762a1bSJed Brown   zee = zs + zm;
749c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
750c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
751c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
752c4762a1bSJed Brown   if (xs + xm == mx) xee = xs + xm - 1;
753c4762a1bSJed Brown   if (ys + ym == my) yee = ys + ym - 1;
754c4762a1bSJed Brown   if (zs + zm == mz) zee = zs + zm - 1;
755c4762a1bSJed Brown   for (k = zes; k < zee; k++) {
756c4762a1bSJed Brown     for (j = yes; j < yee; j++) {
757c4762a1bSJed Brown       for (i = xes; i < xee; i++) {
758c4762a1bSJed Brown         GatherElementData(mx, my, mz, x, c, i, j, k, ex, ec, user);
759c4762a1bSJed Brown         FormElementJacobian(ex, ec, ef, NULL, user);
760c4762a1bSJed Brown         /* put this element's additions into the residuals */
761c4762a1bSJed Brown         for (kk = 0; kk < NB; kk++) {
762c4762a1bSJed Brown           for (jj = 0; jj < NB; jj++) {
763c4762a1bSJed Brown             for (ii = 0; ii < NB; ii++) {
764c4762a1bSJed Brown               PetscInt idx = ii + jj * NB + kk * NB * NB;
765c4762a1bSJed Brown               if (k + kk >= zs && j + jj >= ys && i + ii >= xs && k + kk < zs + zm && j + jj < ys + ym && i + ii < xs + xm) {
766c4762a1bSJed Brown                 if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz)) {
7679371c9d4SSatish 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];
768c4762a1bSJed Brown                 } else {
7699371c9d4SSatish Balay                   for (l = 0; l < 3; l++) f[k + kk][j + jj][i + ii][l] += ef[idx][l];
770c4762a1bSJed Brown                 }
771c4762a1bSJed Brown               }
772c4762a1bSJed Brown             }
773c4762a1bSJed Brown           }
774c4762a1bSJed Brown         }
775c4762a1bSJed Brown       }
776c4762a1bSJed Brown     }
777c4762a1bSJed Brown   }
7789566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
779c4762a1bSJed Brown   PetscFunctionReturn(0);
780c4762a1bSJed Brown }
781c4762a1bSJed Brown 
7829371c9d4SSatish Balay PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ptr) {
783c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
784c4762a1bSJed Brown   AppCtx       *user = (AppCtx *)ptr;
785c4762a1bSJed Brown   PetscInt      i, j, k, l, m, n, s;
786c4762a1bSJed Brown   PetscInt      pi, pj, pk;
787c4762a1bSJed Brown   Field         ef[1];
788c4762a1bSJed Brown   Field         ex[8];
789c4762a1bSJed Brown   PetscScalar   ej[9];
790c4762a1bSJed Brown   CoordField    ec[8];
791c4762a1bSJed Brown   PetscScalar   pjac[9], pjinv[9];
792c4762a1bSJed Brown   PetscScalar   pf[3], py[3];
793c4762a1bSJed Brown   PetscInt      xs, ys, zs;
794c4762a1bSJed Brown   PetscInt      xm, ym, zm;
795c4762a1bSJed Brown   PetscInt      mx, my, mz;
796c4762a1bSJed Brown   DM            cda;
797c4762a1bSJed Brown   CoordField ***c;
798c4762a1bSJed Brown   Vec           C;
799c4762a1bSJed Brown   DM            da;
800c4762a1bSJed Brown   Vec           Xl, Bl;
801c4762a1bSJed Brown   Field      ***x, ***b;
802c4762a1bSJed Brown   PetscInt      sweeps, its;
803c4762a1bSJed Brown   PetscReal     atol, rtol, stol;
804c4762a1bSJed Brown   PetscReal     fnorm0 = 0.0, fnorm, ynorm, xnorm = 0.0;
805c4762a1bSJed Brown 
806c4762a1bSJed Brown   PetscFunctionBegin;
8079566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetSweeps(snes, &sweeps));
8089566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its));
809c4762a1bSJed Brown 
8109566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
8119566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xl));
812*48a46eb9SPierre Jolivet   if (B) PetscCall(DMGetLocalVector(da, &Bl));
8139566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xl));
8149566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xl));
815c4762a1bSJed Brown   if (B) {
8169566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, Bl));
8179566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, Bl));
818c4762a1bSJed Brown   }
8199566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xl, &x));
8209566063dSJacob Faibussowitsch   if (B) PetscCall(DMDAVecGetArray(da, Bl, &b));
821c4762a1bSJed Brown 
8229566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
8239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &C));
8249566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
8259566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
8269566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
827c4762a1bSJed Brown 
828c4762a1bSJed Brown   for (s = 0; s < sweeps; s++) {
829c4762a1bSJed Brown     for (k = zs; k < zs + zm; k++) {
830c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
831c4762a1bSJed Brown         for (i = xs; i < xs + xm; i++) {
832c4762a1bSJed Brown           if (OnBoundary(i, j, k, mx, my, mz)) {
833c4762a1bSJed Brown             BoundaryValue(i, j, k, mx, my, mz, x[k][j][i], user);
834c4762a1bSJed Brown           } else {
835c4762a1bSJed Brown             for (n = 0; n < its; n++) {
836c4762a1bSJed Brown               for (m = 0; m < 9; m++) pjac[m] = 0.;
837c4762a1bSJed Brown               for (m = 0; m < 3; m++) pf[m] = 0.;
838c4762a1bSJed Brown               /* gather the elements for this point */
839c4762a1bSJed Brown               for (pk = -1; pk < 1; pk++) {
840c4762a1bSJed Brown                 for (pj = -1; pj < 1; pj++) {
841c4762a1bSJed Brown                   for (pi = -1; pi < 1; pi++) {
842c4762a1bSJed Brown                     /* check that this element exists */
843c4762a1bSJed Brown                     if (i + pi >= 0 && i + pi < mx - 1 && j + pj >= 0 && j + pj < my - 1 && k + pk >= 0 && k + pk < mz - 1) {
844c4762a1bSJed Brown                       /* create the element function and jacobian */
845c4762a1bSJed Brown                       GatherElementData(mx, my, mz, x, c, i + pi, j + pj, k + pk, ex, ec, user);
846c4762a1bSJed Brown                       FormPBJacobian(-pi, -pj, -pk, ex, ec, ef, ej, user);
847c4762a1bSJed Brown                       /* extract the point named by i,j,k from the whole element jacobian and function */
848c4762a1bSJed Brown                       for (l = 0; l < 3; l++) {
849c4762a1bSJed Brown                         pf[l] += ef[0][l];
8509371c9d4SSatish Balay                         for (m = 0; m < 3; m++) { pjac[3 * m + l] += ej[3 * m + l]; }
851c4762a1bSJed Brown                       }
852c4762a1bSJed Brown                     }
853c4762a1bSJed Brown                   }
854c4762a1bSJed Brown                 }
855c4762a1bSJed Brown               }
856c4762a1bSJed Brown               /* invert */
857c4762a1bSJed Brown               InvertTensor(pjac, pjinv, NULL);
858c4762a1bSJed Brown               /* apply */
8599371c9d4SSatish Balay               if (B)
8609371c9d4SSatish Balay                 for (m = 0; m < 3; m++) { pf[m] -= b[k][j][i][m]; }
861c4762a1bSJed Brown               TensorVector(pjinv, pf, py);
862c4762a1bSJed Brown               xnorm = 0.;
863c4762a1bSJed Brown               for (m = 0; m < 3; m++) {
864c4762a1bSJed Brown                 x[k][j][i][m] -= py[m];
865c4762a1bSJed Brown                 xnorm += PetscRealPart(x[k][j][i][m] * x[k][j][i][m]);
866c4762a1bSJed Brown               }
867c4762a1bSJed Brown               fnorm = PetscRealPart(pf[0] * pf[0] + pf[1] * pf[1] + pf[2] * pf[2]);
868c4762a1bSJed Brown               if (n == 0) fnorm0 = fnorm;
869c4762a1bSJed Brown               ynorm = PetscRealPart(py[0] * py[0] + py[1] * py[1] + py[2] * py[2]);
870c4762a1bSJed Brown               if (fnorm < atol * atol || fnorm < rtol * rtol * fnorm0 || ynorm < stol * stol * xnorm) break;
871c4762a1bSJed Brown             }
872c4762a1bSJed Brown           }
873c4762a1bSJed Brown         }
874c4762a1bSJed Brown       }
875c4762a1bSJed Brown     }
876c4762a1bSJed Brown   }
8779566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xl, &x));
8789566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, Xl, INSERT_VALUES, X));
8799566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, Xl, INSERT_VALUES, X));
8809566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xl));
881c4762a1bSJed Brown   if (B) {
8829566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, Bl, &b));
8839566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da, &Bl));
884c4762a1bSJed Brown   }
8859566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
886c4762a1bSJed Brown   PetscFunctionReturn(0);
887c4762a1bSJed Brown }
888c4762a1bSJed Brown 
8899371c9d4SSatish Balay PetscErrorCode FormCoordinates(DM da, AppCtx *user) {
890c4762a1bSJed Brown   Vec           coords;
891c4762a1bSJed Brown   DM            cda;
892c4762a1bSJed Brown   PetscInt      mx, my, mz;
893c4762a1bSJed Brown   PetscInt      i, j, k, xs, ys, zs, xm, ym, zm;
894c4762a1bSJed Brown   CoordField ***x;
895c4762a1bSJed Brown 
896c4762a1bSJed Brown   PetscFunctionBegin;
8979566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
8989566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(cda, &coords));
8999566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9009566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9019566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, coords, &x));
902c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
903c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
904c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
905c4762a1bSJed Brown         PetscReal cx  = ((PetscReal)i) / (((PetscReal)(mx - 1)));
906c4762a1bSJed Brown         PetscReal cy  = ((PetscReal)j) / (((PetscReal)(my - 1)));
907c4762a1bSJed Brown         PetscReal cz  = ((PetscReal)k) / (((PetscReal)(mz - 1)));
908c4762a1bSJed Brown         PetscReal rad = user->rad + cy * user->height;
909c4762a1bSJed Brown         PetscReal ang = (cx - 0.5) * user->arc;
910c4762a1bSJed Brown         x[k][j][i][0] = rad * PetscSinReal(ang);
911c4762a1bSJed Brown         x[k][j][i][1] = rad * PetscCosReal(ang) - (user->rad + 0.5 * user->height) * PetscCosReal(-0.5 * user->arc);
912c4762a1bSJed Brown         x[k][j][i][2] = user->width * (cz - 0.5);
913c4762a1bSJed Brown       }
914c4762a1bSJed Brown     }
915c4762a1bSJed Brown   }
9169566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, coords, &x));
9179566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinates(da, coords));
9189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coords));
919c4762a1bSJed Brown   PetscFunctionReturn(0);
920c4762a1bSJed Brown }
921c4762a1bSJed Brown 
9229371c9d4SSatish Balay PetscErrorCode InitialGuess(DM da, AppCtx *user, Vec X) {
923c4762a1bSJed Brown   PetscInt i, j, k, xs, ys, zs, xm, ym, zm;
924c4762a1bSJed Brown   PetscInt mx, my, mz;
925c4762a1bSJed Brown   Field ***x;
926c4762a1bSJed Brown 
927c4762a1bSJed Brown   PetscFunctionBegin;
9289566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9299566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9309566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
931c4762a1bSJed Brown 
932c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
933c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
934c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
935c4762a1bSJed Brown         /* reference coordinates */
936c4762a1bSJed Brown         PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx - 1)));
937c4762a1bSJed Brown         PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my - 1)));
938c4762a1bSJed Brown         PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz - 1)));
939c4762a1bSJed Brown         PetscReal o_x = p_x;
940c4762a1bSJed Brown         PetscReal o_y = p_y;
941c4762a1bSJed Brown         PetscReal o_z = p_z;
942c4762a1bSJed Brown         x[k][j][i][0] = o_x - p_x;
943c4762a1bSJed Brown         x[k][j][i][1] = o_y - p_y;
944c4762a1bSJed Brown         x[k][j][i][2] = o_z - p_z;
945c4762a1bSJed Brown       }
946c4762a1bSJed Brown     }
947c4762a1bSJed Brown   }
9489566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
949c4762a1bSJed Brown   PetscFunctionReturn(0);
950c4762a1bSJed Brown }
951c4762a1bSJed Brown 
9529371c9d4SSatish Balay PetscErrorCode FormRHS(DM da, AppCtx *user, Vec X) {
953c4762a1bSJed Brown   PetscInt i, j, k, xs, ys, zs, xm, ym, zm;
954c4762a1bSJed Brown   PetscInt mx, my, mz;
955c4762a1bSJed Brown   Field ***x;
956c4762a1bSJed Brown 
957c4762a1bSJed Brown   PetscFunctionBegin;
9589566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9599566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
961c4762a1bSJed Brown 
962c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
963c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
964c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
965c4762a1bSJed Brown         x[k][j][i][0] = 0.;
966c4762a1bSJed Brown         x[k][j][i][1] = 0.;
967c4762a1bSJed Brown         x[k][j][i][2] = 0.;
968c4762a1bSJed Brown         if (i == (mx - 1) / 2 && j == (my - 1)) x[k][j][i][1] = user->ploading / (mz - 1);
969c4762a1bSJed Brown       }
970c4762a1bSJed Brown     }
971c4762a1bSJed Brown   }
9729566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
973c4762a1bSJed Brown   PetscFunctionReturn(0);
974c4762a1bSJed Brown }
975c4762a1bSJed Brown 
9769371c9d4SSatish Balay PetscErrorCode DisplayLine(SNES snes, Vec X) {
977c4762a1bSJed Brown   PetscInt      r, i, j = 0, k = 0, xs, xm, ys, ym, zs, zm, mx, my, mz;
978c4762a1bSJed Brown   Field      ***x;
979c4762a1bSJed Brown   CoordField ***c;
980c4762a1bSJed Brown   DM            da, cda;
981c4762a1bSJed Brown   Vec           C;
982c4762a1bSJed Brown   PetscMPIInt   size, rank;
983c4762a1bSJed Brown 
984c4762a1bSJed Brown   PetscFunctionBegin;
9859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
9869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
9879566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
9889566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9899566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
9909566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &C));
991c4762a1bSJed Brown   j = my / 2;
992c4762a1bSJed Brown   k = mz / 2;
9939566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9949566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
9959566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
996c4762a1bSJed Brown   for (r = 0; r < size; r++) {
997c4762a1bSJed Brown     if (rank == r) {
998c4762a1bSJed Brown       if (j >= ys && j < ys + ym && k >= zs && k < zs + zm) {
999c4762a1bSJed Brown         for (i = xs; i < xs + xm; i++) {
100063a3b9bcSJacob 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])));
1001c4762a1bSJed Brown         }
1002c4762a1bSJed Brown       }
1003c4762a1bSJed Brown     }
10049566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
1005c4762a1bSJed Brown   }
10069566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
10079566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
1008c4762a1bSJed Brown   PetscFunctionReturn(0);
1009c4762a1bSJed Brown }
1010c4762a1bSJed Brown 
1011c4762a1bSJed Brown /*TEST
1012c4762a1bSJed Brown 
1013c4762a1bSJed Brown    test:
1014c4762a1bSJed Brown       nsize: 2
1015c4762a1bSJed Brown       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 -snes_max_it 7
1016c4762a1bSJed Brown       requires: !single
1017c4762a1bSJed Brown       timeoutfactor: 3
1018c4762a1bSJed Brown 
1019c4762a1bSJed Brown    test:
1020c4762a1bSJed Brown       suffix: 2
1021c4762a1bSJed Brown       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 -snes_max_it 2
1022c4762a1bSJed Brown       requires: !single
1023c4762a1bSJed Brown 
1024c4762a1bSJed Brown    test:
1025c4762a1bSJed Brown       suffix: 3
1026c4762a1bSJed Brown       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 -snes_max_it 4
1027c4762a1bSJed Brown       requires: !single
1028c4762a1bSJed Brown 
1029c4762a1bSJed Brown TEST*/
1030