xref: /petsc/src/snes/tutorials/ex16.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
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
36c4762a1bSJed Brown     steps will jump considerably.  Composed nonlinear solvers may be used to mitigate this difficulty.
37c4762a1bSJed Brown 
38c4762a1bSJed Brown     The initial setup follows the example in pg. 268 of "Nonlinear Finite Element Methods" by Peter Wriggers, but is a
39c4762a1bSJed Brown     3D extension.
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   ------------------------------------------------------------------------F*/
42c4762a1bSJed Brown 
43c4762a1bSJed Brown #include <petscsnes.h>
44c4762a1bSJed Brown #include <petscdm.h>
45c4762a1bSJed Brown #include <petscdmda.h>
46c4762a1bSJed Brown 
47c4762a1bSJed Brown #define QP0 0.2113248654051871
48c4762a1bSJed Brown #define QP1 0.7886751345948129
49c4762a1bSJed Brown #define NQ  2
50c4762a1bSJed Brown #define NB  2
51c4762a1bSJed Brown #define NEB 8
52c4762a1bSJed Brown #define NEQ 8
53c4762a1bSJed Brown #define NPB 24
54c4762a1bSJed Brown 
55c4762a1bSJed Brown #define NVALS NEB *NEQ
56c4762a1bSJed Brown const PetscReal pts[NQ] = {QP0, QP1};
57c4762a1bSJed Brown const PetscReal wts[NQ] = {0.5, 0.5};
58c4762a1bSJed Brown 
59c4762a1bSJed Brown PetscScalar vals[NVALS];
60c4762a1bSJed Brown PetscScalar grad[3 * NVALS];
61c4762a1bSJed Brown 
62c4762a1bSJed Brown typedef PetscScalar Field[3];
63c4762a1bSJed Brown typedef PetscScalar CoordField[3];
64c4762a1bSJed Brown 
65c4762a1bSJed Brown typedef PetscScalar JacField[9];
66c4762a1bSJed Brown 
67c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *, Field ***, Field ***, void *);
68c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *, Field ***, Mat, Mat, void *);
69c4762a1bSJed Brown PetscErrorCode DisplayLine(SNES, Vec);
703274405dSPierre Jolivet PetscErrorCode FormElements(void);
71c4762a1bSJed Brown 
72c4762a1bSJed Brown typedef struct {
73c4762a1bSJed Brown   PetscReal loading;
74c4762a1bSJed Brown   PetscReal mu;
75c4762a1bSJed Brown   PetscReal lambda;
76c4762a1bSJed Brown   PetscReal rad;
77c4762a1bSJed Brown   PetscReal height;
78c4762a1bSJed Brown   PetscReal width;
79c4762a1bSJed Brown   PetscReal arc;
80c4762a1bSJed Brown   PetscReal ploading;
81c4762a1bSJed Brown } AppCtx;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown PetscErrorCode        InitialGuess(DM, AppCtx *, Vec);
84c4762a1bSJed Brown PetscErrorCode        FormRHS(DM, AppCtx *, Vec);
85c4762a1bSJed Brown PetscErrorCode        FormCoordinates(DM, AppCtx *);
86c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES, Vec, Vec, void *);
87c4762a1bSJed Brown 
88d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
89d71ae5a4SJacob Faibussowitsch {
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));
1478434afd1SBarry Smith   PetscCall(DMDASNESSetJacobianLocal(da, (DMDASNESJacobianFn *)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 
191d71ae5a4SJacob Faibussowitsch PetscInt OnBoundary(PetscInt i, PetscInt j, PetscInt k, PetscInt mx, PetscInt my, PetscInt mz)
192d71ae5a4SJacob Faibussowitsch {
193c4762a1bSJed Brown   if ((i == 0 || i == mx - 1) && j == my - 1) return 1;
194c4762a1bSJed Brown   return 0;
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
197d71ae5a4SJacob Faibussowitsch void BoundaryValue(PetscInt i, PetscInt j, PetscInt k, PetscInt mx, PetscInt my, PetscInt mz, PetscScalar *val, AppCtx *user)
198d71ae5a4SJacob Faibussowitsch {
199c4762a1bSJed Brown   /* reference coordinates */
200*f4f49eeaSPierre Jolivet   PetscReal p_x = ((PetscReal)i) / ((PetscReal)(mx - 1));
201*f4f49eeaSPierre Jolivet   PetscReal p_y = ((PetscReal)j) / ((PetscReal)(my - 1));
202*f4f49eeaSPierre Jolivet   PetscReal p_z = ((PetscReal)k) / ((PetscReal)(mz - 1));
203c4762a1bSJed Brown   PetscReal o_x = p_x;
204c4762a1bSJed Brown   PetscReal o_y = p_y;
205c4762a1bSJed Brown   PetscReal o_z = p_z;
206c4762a1bSJed Brown   val[0]        = o_x - p_x;
207c4762a1bSJed Brown   val[1]        = o_y - p_y;
208c4762a1bSJed Brown   val[2]        = o_z - p_z;
209c4762a1bSJed Brown }
210c4762a1bSJed Brown 
211d71ae5a4SJacob Faibussowitsch void InvertTensor(PetscScalar *t, PetscScalar *ti, PetscReal *dett)
212d71ae5a4SJacob Faibussowitsch {
213c4762a1bSJed Brown   const PetscScalar a   = t[0];
214c4762a1bSJed Brown   const PetscScalar b   = t[1];
215c4762a1bSJed Brown   const PetscScalar c   = t[2];
216c4762a1bSJed Brown   const PetscScalar d   = t[3];
217c4762a1bSJed Brown   const PetscScalar e   = t[4];
218c4762a1bSJed Brown   const PetscScalar f   = t[5];
219c4762a1bSJed Brown   const PetscScalar g   = t[6];
220c4762a1bSJed Brown   const PetscScalar h   = t[7];
221c4762a1bSJed Brown   const PetscScalar i   = t[8];
222c4762a1bSJed Brown   const PetscReal   det = PetscRealPart(a * (e * i - f * h) - b * (i * d - f * g) + c * (d * h - e * g));
223c4762a1bSJed Brown   const PetscReal   di  = 1. / det;
224c4762a1bSJed Brown   if (ti) {
225c4762a1bSJed Brown     const PetscScalar A  = (e * i - f * h);
226c4762a1bSJed Brown     const PetscScalar B  = -(d * i - f * g);
227c4762a1bSJed Brown     const PetscScalar C  = (d * h - e * g);
228c4762a1bSJed Brown     const PetscScalar D  = -(b * i - c * h);
229c4762a1bSJed Brown     const PetscScalar E  = (a * i - c * g);
230c4762a1bSJed Brown     const PetscScalar F  = -(a * h - b * g);
231c4762a1bSJed Brown     const PetscScalar G  = (b * f - c * e);
232c4762a1bSJed Brown     const PetscScalar H  = -(a * f - c * d);
233c4762a1bSJed Brown     const PetscScalar II = (a * e - b * d);
234c4762a1bSJed Brown     ti[0]                = di * A;
235c4762a1bSJed Brown     ti[1]                = di * D;
236c4762a1bSJed Brown     ti[2]                = di * G;
237c4762a1bSJed Brown     ti[3]                = di * B;
238c4762a1bSJed Brown     ti[4]                = di * E;
239c4762a1bSJed Brown     ti[5]                = di * H;
240c4762a1bSJed Brown     ti[6]                = di * C;
241c4762a1bSJed Brown     ti[7]                = di * F;
242c4762a1bSJed Brown     ti[8]                = di * II;
243c4762a1bSJed Brown   }
244c4762a1bSJed Brown   if (dett) *dett = det;
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
247d71ae5a4SJacob Faibussowitsch void TensorTensor(PetscScalar *a, PetscScalar *b, PetscScalar *c)
248d71ae5a4SJacob Faibussowitsch {
249c4762a1bSJed Brown   PetscInt i, j, m;
250c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
251c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
252c4762a1bSJed Brown       c[i + 3 * j] = 0;
2539371c9d4SSatish Balay       for (m = 0; m < 3; m++) c[i + 3 * j] += a[m + 3 * j] * b[i + 3 * m];
254c4762a1bSJed Brown     }
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258d71ae5a4SJacob Faibussowitsch void TensorTransposeTensor(PetscScalar *a, PetscScalar *b, PetscScalar *c)
259d71ae5a4SJacob Faibussowitsch {
260c4762a1bSJed Brown   PetscInt i, j, m;
261c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
262c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
263c4762a1bSJed Brown       c[i + 3 * j] = 0;
2649371c9d4SSatish Balay       for (m = 0; m < 3; m++) c[i + 3 * j] += a[3 * m + j] * b[i + 3 * m];
265c4762a1bSJed Brown     }
266c4762a1bSJed Brown   }
267c4762a1bSJed Brown }
268c4762a1bSJed Brown 
269d71ae5a4SJacob Faibussowitsch void TensorVector(PetscScalar *rot, PetscScalar *vec, PetscScalar *tvec)
270d71ae5a4SJacob Faibussowitsch {
271c4762a1bSJed Brown   tvec[0] = rot[0] * vec[0] + rot[1] * vec[1] + rot[2] * vec[2];
272c4762a1bSJed Brown   tvec[1] = rot[3] * vec[0] + rot[4] * vec[1] + rot[5] * vec[2];
273c4762a1bSJed Brown   tvec[2] = rot[6] * vec[0] + rot[7] * vec[1] + rot[8] * vec[2];
274c4762a1bSJed Brown }
275c4762a1bSJed Brown 
276d71ae5a4SJacob Faibussowitsch void DeformationGradient(Field *ex, PetscInt qi, PetscInt qj, PetscInt qk, PetscScalar *invJ, PetscScalar *F)
277d71ae5a4SJacob Faibussowitsch {
278c4762a1bSJed Brown   PetscInt ii, jj, kk, l;
279ad540459SPierre Jolivet   for (l = 0; l < 9; l++) F[l] = 0.;
280c4762a1bSJed Brown   F[0] = 1.;
281c4762a1bSJed Brown   F[4] = 1.;
282c4762a1bSJed Brown   F[8] = 1.;
283c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
284c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
285c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
286c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
287c4762a1bSJed Brown         PetscInt    idx  = ii + jj * NB + kk * NB * NB;
288c4762a1bSJed Brown         PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
289c4762a1bSJed Brown         PetscScalar lgrad[3];
290c4762a1bSJed Brown         TensorVector(invJ, &grad[3 * bidx], lgrad);
2919371c9d4SSatish Balay         F[0] += lgrad[0] * ex[idx][0];
2929371c9d4SSatish Balay         F[1] += lgrad[1] * ex[idx][0];
2939371c9d4SSatish Balay         F[2] += lgrad[2] * ex[idx][0];
2949371c9d4SSatish Balay         F[3] += lgrad[0] * ex[idx][1];
2959371c9d4SSatish Balay         F[4] += lgrad[1] * ex[idx][1];
2969371c9d4SSatish Balay         F[5] += lgrad[2] * ex[idx][1];
2979371c9d4SSatish Balay         F[6] += lgrad[0] * ex[idx][2];
2989371c9d4SSatish Balay         F[7] += lgrad[1] * ex[idx][2];
2999371c9d4SSatish Balay         F[8] += lgrad[2] * ex[idx][2];
300c4762a1bSJed Brown       }
301c4762a1bSJed Brown     }
302c4762a1bSJed Brown   }
303c4762a1bSJed Brown }
304c4762a1bSJed Brown 
305d71ae5a4SJacob Faibussowitsch void DeformationGradientJacobian(PetscInt qi, PetscInt qj, PetscInt qk, PetscInt ii, PetscInt jj, PetscInt kk, PetscInt fld, PetscScalar *invJ, PetscScalar *dF)
306d71ae5a4SJacob Faibussowitsch {
307c4762a1bSJed Brown   PetscInt    l;
308c4762a1bSJed Brown   PetscScalar lgrad[3];
309c4762a1bSJed Brown   PetscInt    idx  = ii + jj * NB + kk * NB * NB;
310c4762a1bSJed Brown   PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
311ad540459SPierre Jolivet   for (l = 0; l < 9; l++) dF[l] = 0.;
312c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
313c4762a1bSJed Brown   TensorVector(invJ, &grad[3 * bidx], lgrad);
3149371c9d4SSatish Balay   dF[3 * fld]     = lgrad[0];
3159371c9d4SSatish Balay   dF[3 * fld + 1] = lgrad[1];
3169371c9d4SSatish Balay   dF[3 * fld + 2] = lgrad[2];
317c4762a1bSJed Brown }
318c4762a1bSJed Brown 
319d71ae5a4SJacob Faibussowitsch void LagrangeGreenStrain(PetscScalar *F, PetscScalar *E)
320d71ae5a4SJacob Faibussowitsch {
321c4762a1bSJed Brown   PetscInt i, j, m;
322c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
323c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
324c4762a1bSJed Brown       E[i + 3 * j] = 0;
3259371c9d4SSatish Balay       for (m = 0; m < 3; m++) E[i + 3 * j] += 0.5 * F[3 * m + j] * F[i + 3 * m];
326c4762a1bSJed Brown     }
327c4762a1bSJed Brown   }
328ad540459SPierre Jolivet   for (i = 0; i < 3; i++) E[i + 3 * i] -= 0.5;
329c4762a1bSJed Brown }
330c4762a1bSJed Brown 
331d71ae5a4SJacob Faibussowitsch void SaintVenantKirchoff(PetscReal lambda, PetscReal mu, PetscScalar *F, PetscScalar *S)
332d71ae5a4SJacob Faibussowitsch {
333c4762a1bSJed Brown   PetscInt    i, j;
334c4762a1bSJed Brown   PetscScalar E[9];
335c4762a1bSJed Brown   PetscScalar trE = 0;
336c4762a1bSJed Brown   LagrangeGreenStrain(F, E);
337ad540459SPierre Jolivet   for (i = 0; i < 3; i++) trE += E[i + 3 * i];
338c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
339c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
340c4762a1bSJed Brown       S[i + 3 * j] = 2. * mu * E[i + 3 * j];
341ad540459SPierre Jolivet       if (i == j) S[i + 3 * j] += trE * lambda;
342c4762a1bSJed Brown     }
343c4762a1bSJed Brown   }
344c4762a1bSJed Brown }
345c4762a1bSJed Brown 
346d71ae5a4SJacob Faibussowitsch void SaintVenantKirchoffJacobian(PetscReal lambda, PetscReal mu, PetscScalar *F, PetscScalar *dF, PetscScalar *dS)
347d71ae5a4SJacob Faibussowitsch {
348c4762a1bSJed Brown   PetscScalar FtdF[9], dE[9];
349c4762a1bSJed Brown   PetscInt    i, j;
350c4762a1bSJed Brown   PetscScalar dtrE = 0.;
351c4762a1bSJed Brown   TensorTransposeTensor(dF, F, dE);
352c4762a1bSJed Brown   TensorTransposeTensor(F, dF, FtdF);
353c4762a1bSJed Brown   for (i = 0; i < 9; i++) dE[i] += FtdF[i];
354c4762a1bSJed Brown   for (i = 0; i < 9; i++) dE[i] *= 0.5;
355ad540459SPierre Jolivet   for (i = 0; i < 3; i++) dtrE += dE[i + 3 * i];
356c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
357c4762a1bSJed Brown     for (j = 0; j < 3; j++) {
358c4762a1bSJed Brown       dS[i + 3 * j] = 2. * mu * dE[i + 3 * j];
359ad540459SPierre Jolivet       if (i == j) dS[i + 3 * j] += dtrE * lambda;
360c4762a1bSJed Brown     }
361c4762a1bSJed Brown   }
362c4762a1bSJed Brown }
363c4762a1bSJed Brown 
364d71ae5a4SJacob Faibussowitsch PetscErrorCode FormElements()
365d71ae5a4SJacob Faibussowitsch {
366c4762a1bSJed Brown   PetscInt  i, j, k, ii, jj, kk;
367c4762a1bSJed Brown   PetscReal bx, by, bz, dbx, dby, dbz;
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   PetscFunctionBegin;
370c4762a1bSJed Brown   /* construct the basis function values and derivatives */
371c4762a1bSJed Brown   for (k = 0; k < NB; k++) {
372c4762a1bSJed Brown     for (j = 0; j < NB; j++) {
373c4762a1bSJed Brown       for (i = 0; i < NB; i++) {
374c4762a1bSJed Brown         /* loop over the quadrature points */
375c4762a1bSJed Brown         for (kk = 0; kk < NQ; kk++) {
376c4762a1bSJed Brown           for (jj = 0; jj < NQ; jj++) {
377c4762a1bSJed Brown             for (ii = 0; ii < NQ; ii++) {
378c4762a1bSJed Brown               PetscInt idx = ii + NQ * jj + NQ * NQ * kk + NEQ * i + NEQ * NB * j + NEQ * NB * NB * k;
379c4762a1bSJed Brown               bx           = pts[ii];
380c4762a1bSJed Brown               by           = pts[jj];
381c4762a1bSJed Brown               bz           = pts[kk];
382c4762a1bSJed Brown               dbx          = 1.;
383c4762a1bSJed Brown               dby          = 1.;
384c4762a1bSJed Brown               dbz          = 1.;
3859371c9d4SSatish Balay               if (i == 0) {
3869371c9d4SSatish Balay                 bx  = 1. - bx;
3879371c9d4SSatish Balay                 dbx = -1;
3889371c9d4SSatish Balay               }
3899371c9d4SSatish Balay               if (j == 0) {
3909371c9d4SSatish Balay                 by  = 1. - by;
3919371c9d4SSatish Balay                 dby = -1;
3929371c9d4SSatish Balay               }
3939371c9d4SSatish Balay               if (k == 0) {
3949371c9d4SSatish Balay                 bz  = 1. - bz;
3959371c9d4SSatish Balay                 dbz = -1;
3969371c9d4SSatish Balay               }
397c4762a1bSJed Brown               vals[idx]         = bx * by * bz;
398c4762a1bSJed Brown               grad[3 * idx + 0] = dbx * by * bz;
399c4762a1bSJed Brown               grad[3 * idx + 1] = dby * bx * bz;
400c4762a1bSJed Brown               grad[3 * idx + 2] = dbz * bx * by;
401c4762a1bSJed Brown             }
402c4762a1bSJed Brown           }
403c4762a1bSJed Brown         }
404c4762a1bSJed Brown       }
405c4762a1bSJed Brown     }
406c4762a1bSJed Brown   }
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408c4762a1bSJed Brown }
409c4762a1bSJed Brown 
410d71ae5a4SJacob 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)
411d71ae5a4SJacob Faibussowitsch {
412c4762a1bSJed Brown   PetscInt m;
413c4762a1bSJed Brown   PetscInt ii, jj, kk;
414c4762a1bSJed Brown   /* gather the data -- loop over element unknowns */
415c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
416c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
417c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
418c4762a1bSJed Brown         PetscInt idx = ii + jj * NB + kk * NB * NB;
419c4762a1bSJed Brown         /* decouple the boundary nodes for the displacement variables */
420c4762a1bSJed Brown         if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz)) {
421c4762a1bSJed Brown           BoundaryValue(i + ii, j + jj, k + kk, mx, my, mz, ex[idx], user);
422c4762a1bSJed Brown         } else {
423ad540459SPierre Jolivet           for (m = 0; m < 3; m++) ex[idx][m] = x[k + kk][j + jj][i + ii][m];
424c4762a1bSJed Brown         }
425ad540459SPierre Jolivet         for (m = 0; m < 3; m++) ec[idx][m] = c[k + kk][j + jj][i + ii][m];
426c4762a1bSJed Brown       }
427c4762a1bSJed Brown     }
428c4762a1bSJed Brown   }
429c4762a1bSJed Brown }
430c4762a1bSJed Brown 
431d71ae5a4SJacob Faibussowitsch void QuadraturePointGeometricJacobian(CoordField *ec, PetscInt qi, PetscInt qj, PetscInt qk, PetscScalar *J)
432d71ae5a4SJacob Faibussowitsch {
433c4762a1bSJed Brown   PetscInt ii, jj, kk;
434c4762a1bSJed Brown   /* construct the gradient at the given quadrature point named by i,j,k */
435ad540459SPierre Jolivet   for (ii = 0; ii < 9; ii++) J[ii] = 0;
436c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
437c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
438c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
439c4762a1bSJed Brown         PetscInt idx  = ii + jj * NB + kk * NB * NB;
440c4762a1bSJed Brown         PetscInt bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
4419371c9d4SSatish Balay         J[0] += grad[3 * bidx + 0] * ec[idx][0];
4429371c9d4SSatish Balay         J[1] += grad[3 * bidx + 1] * ec[idx][0];
4439371c9d4SSatish Balay         J[2] += grad[3 * bidx + 2] * ec[idx][0];
4449371c9d4SSatish Balay         J[3] += grad[3 * bidx + 0] * ec[idx][1];
4459371c9d4SSatish Balay         J[4] += grad[3 * bidx + 1] * ec[idx][1];
4469371c9d4SSatish Balay         J[5] += grad[3 * bidx + 2] * ec[idx][1];
4479371c9d4SSatish Balay         J[6] += grad[3 * bidx + 0] * ec[idx][2];
4489371c9d4SSatish Balay         J[7] += grad[3 * bidx + 1] * ec[idx][2];
4499371c9d4SSatish Balay         J[8] += grad[3 * bidx + 2] * ec[idx][2];
450c4762a1bSJed Brown       }
451c4762a1bSJed Brown     }
452c4762a1bSJed Brown   }
453c4762a1bSJed Brown }
454c4762a1bSJed Brown 
455d71ae5a4SJacob Faibussowitsch void FormElementJacobian(Field *ex, CoordField *ec, Field *ef, PetscScalar *ej, AppCtx *user)
456d71ae5a4SJacob Faibussowitsch {
457c4762a1bSJed Brown   PetscReal   vol;
458c4762a1bSJed Brown   PetscScalar J[9];
459c4762a1bSJed Brown   PetscScalar invJ[9];
460c4762a1bSJed Brown   PetscScalar F[9], S[9], dF[9], dS[9], dFS[9], FdS[9], FS[9];
461c4762a1bSJed Brown   PetscReal   scl;
462c4762a1bSJed Brown   PetscInt    i, j, k, l, ii, jj, kk, ll, qi, qj, qk, m;
463c4762a1bSJed Brown 
4649371c9d4SSatish Balay   if (ej)
4659371c9d4SSatish Balay     for (i = 0; i < NPB * NPB; i++) ej[i] = 0.;
4669371c9d4SSatish Balay   if (ef)
4679371c9d4SSatish Balay     for (i = 0; i < NEB; i++) {
4689371c9d4SSatish Balay       ef[i][0] = 0.;
4699371c9d4SSatish Balay       ef[i][1] = 0.;
4709371c9d4SSatish Balay       ef[i][2] = 0.;
4719371c9d4SSatish Balay     }
472c4762a1bSJed Brown   /* loop over quadrature */
473c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
474c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
475c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
476c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec, qi, qj, qk, J);
477c4762a1bSJed Brown         InvertTensor(J, invJ, &vol);
478c4762a1bSJed Brown         scl = vol * wts[qi] * wts[qj] * wts[qk];
479c4762a1bSJed Brown         DeformationGradient(ex, qi, qj, qk, invJ, F);
480c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda, user->mu, F, S);
481c4762a1bSJed Brown         /* form the function */
482c4762a1bSJed Brown         if (ef) {
483c4762a1bSJed Brown           TensorTensor(F, S, FS);
484c4762a1bSJed Brown           for (kk = 0; kk < NB; kk++) {
485c4762a1bSJed Brown             for (jj = 0; jj < NB; jj++) {
486c4762a1bSJed Brown               for (ii = 0; ii < NB; ii++) {
487c4762a1bSJed Brown                 PetscInt    idx  = ii + jj * NB + kk * NB * NB;
488c4762a1bSJed Brown                 PetscInt    bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
489c4762a1bSJed Brown                 PetscScalar lgrad[3];
490c4762a1bSJed Brown                 TensorVector(invJ, &grad[3 * bidx], lgrad);
491c4762a1bSJed Brown                 /* mu*F : grad phi_{u,v,w} */
492ad540459SPierre 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]);
493c4762a1bSJed Brown                 ef[idx][1] -= scl * user->loading * vals[bidx];
494c4762a1bSJed Brown               }
495c4762a1bSJed Brown             }
496c4762a1bSJed Brown           }
497c4762a1bSJed Brown         }
498c4762a1bSJed Brown         /* form the jacobian */
499c4762a1bSJed Brown         if (ej) {
500c4762a1bSJed Brown           /* loop over trialfunctions */
501c4762a1bSJed Brown           for (k = 0; k < NB; k++) {
502c4762a1bSJed Brown             for (j = 0; j < NB; j++) {
503c4762a1bSJed Brown               for (i = 0; i < NB; i++) {
504c4762a1bSJed Brown                 for (l = 0; l < 3; l++) {
505c4762a1bSJed Brown                   PetscInt tridx = l + 3 * (i + j * NB + k * NB * NB);
506c4762a1bSJed Brown                   DeformationGradientJacobian(qi, qj, qk, i, j, k, l, invJ, dF);
507c4762a1bSJed Brown                   SaintVenantKirchoffJacobian(user->lambda, user->mu, F, dF, dS);
508c4762a1bSJed Brown                   TensorTensor(dF, S, dFS);
509c4762a1bSJed Brown                   TensorTensor(F, dS, FdS);
510c4762a1bSJed Brown                   for (m = 0; m < 9; m++) dFS[m] += FdS[m];
511c4762a1bSJed Brown                   /* loop over testfunctions */
512c4762a1bSJed Brown                   for (kk = 0; kk < NB; kk++) {
513c4762a1bSJed Brown                     for (jj = 0; jj < NB; jj++) {
514c4762a1bSJed Brown                       for (ii = 0; ii < NB; ii++) {
515c4762a1bSJed Brown                         PetscInt    idx  = ii + jj * NB + kk * NB * NB;
516c4762a1bSJed Brown                         PetscInt    bidx = 8 * idx + qi + NQ * qj + NQ * NQ * qk;
517c4762a1bSJed Brown                         PetscScalar lgrad[3];
518c4762a1bSJed Brown                         TensorVector(invJ, &grad[3 * bidx], lgrad);
519c4762a1bSJed Brown                         for (ll = 0; ll < 3; ll++) {
520c4762a1bSJed Brown                           PetscInt teidx = ll + 3 * (ii + jj * NB + kk * NB * NB);
5219371c9d4SSatish Balay                           ej[teidx + NPB * tridx] += scl * (lgrad[0] * dFS[3 * ll + 0] + lgrad[1] * dFS[3 * ll + 1] + lgrad[2] * dFS[3 * ll + 2]);
522c4762a1bSJed Brown                         }
523c4762a1bSJed Brown                       }
524c4762a1bSJed Brown                     }
525c4762a1bSJed Brown                   } /* end of testfunctions */
526c4762a1bSJed Brown                 }
527c4762a1bSJed Brown               }
528c4762a1bSJed Brown             }
529c4762a1bSJed Brown           } /* end of trialfunctions */
530c4762a1bSJed Brown         }
531c4762a1bSJed Brown       }
532c4762a1bSJed Brown     }
533c4762a1bSJed Brown   } /* end of quadrature points */
534c4762a1bSJed Brown }
535c4762a1bSJed Brown 
536d71ae5a4SJacob Faibussowitsch void FormPBJacobian(PetscInt i, PetscInt j, PetscInt k, Field *ex, CoordField *ec, Field *ef, PetscScalar *ej, AppCtx *user)
537d71ae5a4SJacob Faibussowitsch {
538c4762a1bSJed Brown   PetscReal   vol;
539c4762a1bSJed Brown   PetscScalar J[9];
540c4762a1bSJed Brown   PetscScalar invJ[9];
541c4762a1bSJed Brown   PetscScalar F[9], S[9], dF[9], dS[9], dFS[9], FdS[9], FS[9];
542c4762a1bSJed Brown   PetscReal   scl;
543c4762a1bSJed Brown   PetscInt    l, ll, qi, qj, qk, m;
544c4762a1bSJed Brown   PetscInt    idx = i + j * NB + k * NB * NB;
545c4762a1bSJed Brown   PetscScalar lgrad[3];
546c4762a1bSJed Brown 
5479371c9d4SSatish Balay   if (ej)
5489371c9d4SSatish Balay     for (l = 0; l < 9; l++) ej[l] = 0.;
5499371c9d4SSatish Balay   if (ef)
5509371c9d4SSatish Balay     for (l = 0; l < 1; l++) {
5519371c9d4SSatish Balay       ef[l][0] = 0.;
5529371c9d4SSatish Balay       ef[l][1] = 0.;
5539371c9d4SSatish Balay       ef[l][2] = 0.;
5549371c9d4SSatish Balay     }
555c4762a1bSJed Brown   /* loop over quadrature */
556c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
557c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
558c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
559c4762a1bSJed Brown         PetscInt bidx = NEB * idx + qi + NQ * qj + NQ * NQ * qk;
560c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec, qi, qj, qk, J);
561c4762a1bSJed Brown         InvertTensor(J, invJ, &vol);
562c4762a1bSJed Brown         TensorVector(invJ, &grad[3 * bidx], lgrad);
563c4762a1bSJed Brown         scl = vol * wts[qi] * wts[qj] * wts[qk];
564c4762a1bSJed Brown         DeformationGradient(ex, qi, qj, qk, invJ, F);
565c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda, user->mu, F, S);
566c4762a1bSJed Brown         /* form the function */
567c4762a1bSJed Brown         if (ef) {
568c4762a1bSJed Brown           TensorTensor(F, S, FS);
569ad540459SPierre Jolivet           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]);
570c4762a1bSJed Brown           ef[0][1] -= scl * user->loading * vals[bidx];
571c4762a1bSJed Brown         }
572c4762a1bSJed Brown         /* form the jacobian */
573c4762a1bSJed Brown         if (ej) {
574c4762a1bSJed Brown           for (l = 0; l < 3; l++) {
575c4762a1bSJed Brown             DeformationGradientJacobian(qi, qj, qk, i, j, k, l, invJ, dF);
576c4762a1bSJed Brown             SaintVenantKirchoffJacobian(user->lambda, user->mu, F, dF, dS);
577c4762a1bSJed Brown             TensorTensor(dF, S, dFS);
578c4762a1bSJed Brown             TensorTensor(F, dS, FdS);
579c4762a1bSJed Brown             for (m = 0; m < 9; m++) dFS[m] += FdS[m];
580ad540459SPierre Jolivet             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]);
581c4762a1bSJed Brown           }
582c4762a1bSJed Brown         }
583c4762a1bSJed Brown       }
584c4762a1bSJed Brown     } /* end of quadrature points */
585c4762a1bSJed Brown   }
586c4762a1bSJed Brown }
587c4762a1bSJed Brown 
588d71ae5a4SJacob Faibussowitsch void ApplyBCsElement(PetscInt mx, PetscInt my, PetscInt mz, PetscInt i, PetscInt j, PetscInt k, PetscScalar *jacobian)
589d71ae5a4SJacob Faibussowitsch {
590c4762a1bSJed Brown   PetscInt ii, jj, kk, ll, ei, ej, ek, el;
591c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
592c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
593c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
594c4762a1bSJed Brown         for (ll = 0; ll < 3; ll++) {
595c4762a1bSJed Brown           PetscInt tridx = ll + 3 * (ii + jj * NB + kk * NB * NB);
596c4762a1bSJed Brown           for (ek = 0; ek < NB; ek++) {
597c4762a1bSJed Brown             for (ej = 0; ej < NB; ej++) {
598c4762a1bSJed Brown               for (ei = 0; ei < NB; ei++) {
599c4762a1bSJed Brown                 for (el = 0; el < 3; el++) {
600c4762a1bSJed Brown                   if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz) || OnBoundary(i + ei, j + ej, k + ek, mx, my, mz)) {
601c4762a1bSJed Brown                     PetscInt teidx = el + 3 * (ei + ej * NB + ek * NB * NB);
602c4762a1bSJed Brown                     if (teidx == tridx) {
603c4762a1bSJed Brown                       jacobian[tridx + NPB * teidx] = 1.;
604c4762a1bSJed Brown                     } else {
605c4762a1bSJed Brown                       jacobian[tridx + NPB * teidx] = 0.;
606c4762a1bSJed Brown                     }
607c4762a1bSJed Brown                   }
608c4762a1bSJed Brown                 }
609c4762a1bSJed Brown               }
610c4762a1bSJed Brown             }
611c4762a1bSJed Brown           }
612c4762a1bSJed Brown         }
613c4762a1bSJed Brown       }
614c4762a1bSJed Brown     }
615c4762a1bSJed Brown   }
616c4762a1bSJed Brown }
617c4762a1bSJed Brown 
618d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobianLocal(DMDALocalInfo *info, Field ***x, Mat jacpre, Mat jac, void *ptr)
619d71ae5a4SJacob Faibussowitsch {
620c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
621c4762a1bSJed Brown   AppCtx     *user = (AppCtx *)ptr;
622c4762a1bSJed Brown   PetscInt    i, j, k, m, l;
623c4762a1bSJed Brown   PetscInt    ii, jj, kk;
624c4762a1bSJed Brown   PetscScalar ej[NPB * NPB];
625c4762a1bSJed Brown   PetscScalar vals[NPB * NPB];
626c4762a1bSJed Brown   Field       ex[NEB];
627c4762a1bSJed Brown   CoordField  ec[NEB];
628c4762a1bSJed Brown 
629c4762a1bSJed Brown   PetscInt      xs = info->xs, ys = info->ys, zs = info->zs;
630c4762a1bSJed Brown   PetscInt      xm = info->xm, ym = info->ym, zm = info->zm;
631c4762a1bSJed Brown   PetscInt      xes, yes, zes, xee, yee, zee;
632c4762a1bSJed Brown   PetscInt      mx = info->mx, my = info->my, mz = info->mz;
633c4762a1bSJed Brown   DM            cda;
634c4762a1bSJed Brown   CoordField ***c;
635c4762a1bSJed Brown   Vec           C;
636c4762a1bSJed Brown   PetscInt      nrows;
637c4762a1bSJed Brown   MatStencil    col[NPB], row[NPB];
638c4762a1bSJed Brown   PetscScalar   v[9];
639c4762a1bSJed Brown 
640c4762a1bSJed Brown   PetscFunctionBegin;
6419566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(info->da, &cda));
6429566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(info->da, &C));
6439566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
6449566063dSJacob Faibussowitsch   PetscCall(MatScale(jac, 0.0));
645c4762a1bSJed Brown 
646c4762a1bSJed Brown   xes = xs;
647c4762a1bSJed Brown   yes = ys;
648c4762a1bSJed Brown   zes = zs;
649c4762a1bSJed Brown   xee = xs + xm;
650c4762a1bSJed Brown   yee = ys + ym;
651c4762a1bSJed Brown   zee = zs + zm;
652c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
653c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
654c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
655c4762a1bSJed Brown   if (xs + xm == mx) xee = xs + xm - 1;
656c4762a1bSJed Brown   if (ys + ym == my) yee = ys + ym - 1;
657c4762a1bSJed Brown   if (zs + zm == mz) zee = zs + zm - 1;
658c4762a1bSJed Brown   for (k = zes; k < zee; k++) {
659c4762a1bSJed Brown     for (j = yes; j < yee; j++) {
660c4762a1bSJed Brown       for (i = xes; i < xee; i++) {
661c4762a1bSJed Brown         GatherElementData(mx, my, mz, x, c, i, j, k, ex, ec, user);
662c4762a1bSJed Brown         FormElementJacobian(ex, ec, NULL, ej, user);
663c4762a1bSJed Brown         ApplyBCsElement(mx, my, mz, i, j, k, ej);
664c4762a1bSJed Brown         nrows = 0.;
665c4762a1bSJed Brown         for (kk = 0; kk < NB; kk++) {
666c4762a1bSJed Brown           for (jj = 0; jj < NB; jj++) {
667c4762a1bSJed Brown             for (ii = 0; ii < NB; ii++) {
668c4762a1bSJed Brown               PetscInt idx = ii + jj * 2 + kk * 4;
669c4762a1bSJed Brown               for (m = 0; m < 3; m++) {
670c4762a1bSJed Brown                 col[3 * idx + m].i = i + ii;
671c4762a1bSJed Brown                 col[3 * idx + m].j = j + jj;
672c4762a1bSJed Brown                 col[3 * idx + m].k = k + kk;
673c4762a1bSJed Brown                 col[3 * idx + m].c = m;
674c4762a1bSJed Brown                 if (i + ii >= xs && i + ii < xm + xs && j + jj >= ys && j + jj < ys + ym && k + kk >= zs && k + kk < zs + zm) {
675c4762a1bSJed Brown                   row[nrows].i = i + ii;
676c4762a1bSJed Brown                   row[nrows].j = j + jj;
677c4762a1bSJed Brown                   row[nrows].k = k + kk;
678c4762a1bSJed Brown                   row[nrows].c = m;
679c4762a1bSJed Brown                   for (l = 0; l < NPB; l++) vals[NPB * nrows + l] = ej[NPB * (3 * idx + m) + l];
680c4762a1bSJed Brown                   nrows++;
681c4762a1bSJed Brown                 }
682c4762a1bSJed Brown               }
683c4762a1bSJed Brown             }
684c4762a1bSJed Brown           }
685c4762a1bSJed Brown         }
6869566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(jac, nrows, row, NPB, col, vals, ADD_VALUES));
687c4762a1bSJed Brown       }
688c4762a1bSJed Brown     }
689c4762a1bSJed Brown   }
690c4762a1bSJed Brown 
6919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FLUSH_ASSEMBLY));
6929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FLUSH_ASSEMBLY));
693c4762a1bSJed Brown 
694c4762a1bSJed Brown   /* set the diagonal */
6959371c9d4SSatish Balay   v[0] = 1.;
6969371c9d4SSatish Balay   v[1] = 0.;
6979371c9d4SSatish Balay   v[2] = 0.;
6989371c9d4SSatish Balay   v[3] = 0.;
6999371c9d4SSatish Balay   v[4] = 1.;
7009371c9d4SSatish Balay   v[5] = 0.;
7019371c9d4SSatish Balay   v[6] = 0.;
7029371c9d4SSatish Balay   v[7] = 0.;
7039371c9d4SSatish Balay   v[8] = 1.;
704c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
705c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
706c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
707c4762a1bSJed Brown         if (OnBoundary(i, j, k, mx, my, mz)) {
708c4762a1bSJed Brown           for (m = 0; m < 3; m++) {
709c4762a1bSJed Brown             col[m].i = i;
710c4762a1bSJed Brown             col[m].j = j;
711c4762a1bSJed Brown             col[m].k = k;
712c4762a1bSJed Brown             col[m].c = m;
713c4762a1bSJed Brown           }
7149566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(jac, 3, col, 3, col, v, INSERT_VALUES));
715c4762a1bSJed Brown         }
716c4762a1bSJed Brown       }
717c4762a1bSJed Brown     }
718c4762a1bSJed Brown   }
719c4762a1bSJed Brown 
7209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
7219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
722c4762a1bSJed Brown 
7239566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
7243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
725c4762a1bSJed Brown }
726c4762a1bSJed Brown 
727d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field ***x, Field ***f, void *ptr)
728d71ae5a4SJacob Faibussowitsch {
729c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
730c4762a1bSJed Brown   AppCtx  *user = (AppCtx *)ptr;
731c4762a1bSJed Brown   PetscInt i, j, k, l;
732c4762a1bSJed Brown   PetscInt ii, jj, kk;
733c4762a1bSJed Brown 
734c4762a1bSJed Brown   Field      ef[NEB];
735c4762a1bSJed Brown   Field      ex[NEB];
736c4762a1bSJed Brown   CoordField ec[NEB];
737c4762a1bSJed Brown 
738c4762a1bSJed Brown   PetscInt      xs = info->xs, ys = info->ys, zs = info->zs;
739c4762a1bSJed Brown   PetscInt      xm = info->xm, ym = info->ym, zm = info->zm;
740c4762a1bSJed Brown   PetscInt      xes, yes, zes, xee, yee, zee;
741c4762a1bSJed Brown   PetscInt      mx = info->mx, my = info->my, mz = info->mz;
742c4762a1bSJed Brown   DM            cda;
743c4762a1bSJed Brown   CoordField ***c;
744c4762a1bSJed Brown   Vec           C;
745c4762a1bSJed Brown 
746c4762a1bSJed Brown   PetscFunctionBegin;
7479566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(info->da, &cda));
7489566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(info->da, &C));
7499566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
7509566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(info->da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
7519566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(info->da, &xs, &ys, &zs, &xm, &ym, &zm));
752c4762a1bSJed Brown 
753c4762a1bSJed Brown   /* loop over elements */
754c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
755c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
756c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
757ad540459SPierre Jolivet         for (l = 0; l < 3; l++) f[k][j][i][l] = 0.;
758c4762a1bSJed Brown       }
759c4762a1bSJed Brown     }
760c4762a1bSJed Brown   }
761c4762a1bSJed Brown   /* element starts and ends */
762c4762a1bSJed Brown   xes = xs;
763c4762a1bSJed Brown   yes = ys;
764c4762a1bSJed Brown   zes = zs;
765c4762a1bSJed Brown   xee = xs + xm;
766c4762a1bSJed Brown   yee = ys + ym;
767c4762a1bSJed Brown   zee = zs + zm;
768c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
769c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
770c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
771c4762a1bSJed Brown   if (xs + xm == mx) xee = xs + xm - 1;
772c4762a1bSJed Brown   if (ys + ym == my) yee = ys + ym - 1;
773c4762a1bSJed Brown   if (zs + zm == mz) zee = zs + zm - 1;
774c4762a1bSJed Brown   for (k = zes; k < zee; k++) {
775c4762a1bSJed Brown     for (j = yes; j < yee; j++) {
776c4762a1bSJed Brown       for (i = xes; i < xee; i++) {
777c4762a1bSJed Brown         GatherElementData(mx, my, mz, x, c, i, j, k, ex, ec, user);
778c4762a1bSJed Brown         FormElementJacobian(ex, ec, ef, NULL, user);
779c4762a1bSJed Brown         /* put this element's additions into the residuals */
780c4762a1bSJed Brown         for (kk = 0; kk < NB; kk++) {
781c4762a1bSJed Brown           for (jj = 0; jj < NB; jj++) {
782c4762a1bSJed Brown             for (ii = 0; ii < NB; ii++) {
783c4762a1bSJed Brown               PetscInt idx = ii + jj * NB + kk * NB * NB;
784c4762a1bSJed Brown               if (k + kk >= zs && j + jj >= ys && i + ii >= xs && k + kk < zs + zm && j + jj < ys + ym && i + ii < xs + xm) {
785c4762a1bSJed Brown                 if (OnBoundary(i + ii, j + jj, k + kk, mx, my, mz)) {
7869371c9d4SSatish 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];
787c4762a1bSJed Brown                 } else {
7889371c9d4SSatish Balay                   for (l = 0; l < 3; l++) f[k + kk][j + jj][i + ii][l] += ef[idx][l];
789c4762a1bSJed Brown                 }
790c4762a1bSJed Brown               }
791c4762a1bSJed Brown             }
792c4762a1bSJed Brown           }
793c4762a1bSJed Brown         }
794c4762a1bSJed Brown       }
795c4762a1bSJed Brown     }
796c4762a1bSJed Brown   }
7979566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
7983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
799c4762a1bSJed Brown }
800c4762a1bSJed Brown 
801d71ae5a4SJacob Faibussowitsch PetscErrorCode NonlinearGS(SNES snes, Vec X, Vec B, void *ptr)
802d71ae5a4SJacob Faibussowitsch {
803c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
804c4762a1bSJed Brown   AppCtx       *user = (AppCtx *)ptr;
805c4762a1bSJed Brown   PetscInt      i, j, k, l, m, n, s;
806c4762a1bSJed Brown   PetscInt      pi, pj, pk;
807c4762a1bSJed Brown   Field         ef[1];
808c4762a1bSJed Brown   Field         ex[8];
809c4762a1bSJed Brown   PetscScalar   ej[9];
810c4762a1bSJed Brown   CoordField    ec[8];
811c4762a1bSJed Brown   PetscScalar   pjac[9], pjinv[9];
812c4762a1bSJed Brown   PetscScalar   pf[3], py[3];
813c4762a1bSJed Brown   PetscInt      xs, ys, zs;
814c4762a1bSJed Brown   PetscInt      xm, ym, zm;
815c4762a1bSJed Brown   PetscInt      mx, my, mz;
816c4762a1bSJed Brown   DM            cda;
817c4762a1bSJed Brown   CoordField ***c;
818c4762a1bSJed Brown   Vec           C;
819c4762a1bSJed Brown   DM            da;
820c4762a1bSJed Brown   Vec           Xl, Bl;
821c4762a1bSJed Brown   Field      ***x, ***b;
822c4762a1bSJed Brown   PetscInt      sweeps, its;
823c4762a1bSJed Brown   PetscReal     atol, rtol, stol;
824c4762a1bSJed Brown   PetscReal     fnorm0 = 0.0, fnorm, ynorm, xnorm = 0.0;
825c4762a1bSJed Brown 
826c4762a1bSJed Brown   PetscFunctionBegin;
8279566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetSweeps(snes, &sweeps));
8289566063dSJacob Faibussowitsch   PetscCall(SNESNGSGetTolerances(snes, &atol, &rtol, &stol, &its));
829c4762a1bSJed Brown 
8309566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
8319566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xl));
83248a46eb9SPierre Jolivet   if (B) PetscCall(DMGetLocalVector(da, &Bl));
8339566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xl));
8349566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xl));
835c4762a1bSJed Brown   if (B) {
8369566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(da, B, INSERT_VALUES, Bl));
8379566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(da, B, INSERT_VALUES, Bl));
838c4762a1bSJed Brown   }
8399566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Xl, &x));
8409566063dSJacob Faibussowitsch   if (B) PetscCall(DMDAVecGetArray(da, Bl, &b));
841c4762a1bSJed Brown 
8429566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
8439566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &C));
8449566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
8459566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
8469566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
847c4762a1bSJed Brown 
848c4762a1bSJed Brown   for (s = 0; s < sweeps; s++) {
849c4762a1bSJed Brown     for (k = zs; k < zs + zm; k++) {
850c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
851c4762a1bSJed Brown         for (i = xs; i < xs + xm; i++) {
852c4762a1bSJed Brown           if (OnBoundary(i, j, k, mx, my, mz)) {
853c4762a1bSJed Brown             BoundaryValue(i, j, k, mx, my, mz, x[k][j][i], user);
854c4762a1bSJed Brown           } else {
855c4762a1bSJed Brown             for (n = 0; n < its; n++) {
856c4762a1bSJed Brown               for (m = 0; m < 9; m++) pjac[m] = 0.;
857c4762a1bSJed Brown               for (m = 0; m < 3; m++) pf[m] = 0.;
858c4762a1bSJed Brown               /* gather the elements for this point */
859c4762a1bSJed Brown               for (pk = -1; pk < 1; pk++) {
860c4762a1bSJed Brown                 for (pj = -1; pj < 1; pj++) {
861c4762a1bSJed Brown                   for (pi = -1; pi < 1; pi++) {
862c4762a1bSJed Brown                     /* check that this element exists */
863c4762a1bSJed Brown                     if (i + pi >= 0 && i + pi < mx - 1 && j + pj >= 0 && j + pj < my - 1 && k + pk >= 0 && k + pk < mz - 1) {
864c4762a1bSJed Brown                       /* create the element function and jacobian */
865c4762a1bSJed Brown                       GatherElementData(mx, my, mz, x, c, i + pi, j + pj, k + pk, ex, ec, user);
866c4762a1bSJed Brown                       FormPBJacobian(-pi, -pj, -pk, ex, ec, ef, ej, user);
867c4762a1bSJed Brown                       /* extract the point named by i,j,k from the whole element jacobian and function */
868c4762a1bSJed Brown                       for (l = 0; l < 3; l++) {
869c4762a1bSJed Brown                         pf[l] += ef[0][l];
870ad540459SPierre Jolivet                         for (m = 0; m < 3; m++) pjac[3 * m + l] += ej[3 * m + l];
871c4762a1bSJed Brown                       }
872c4762a1bSJed Brown                     }
873c4762a1bSJed Brown                   }
874c4762a1bSJed Brown                 }
875c4762a1bSJed Brown               }
876c4762a1bSJed Brown               /* invert */
877c4762a1bSJed Brown               InvertTensor(pjac, pjinv, NULL);
878c4762a1bSJed Brown               /* apply */
8799371c9d4SSatish Balay               if (B)
880ad540459SPierre Jolivet                 for (m = 0; m < 3; m++) pf[m] -= b[k][j][i][m];
881c4762a1bSJed Brown               TensorVector(pjinv, pf, py);
882c4762a1bSJed Brown               xnorm = 0.;
883c4762a1bSJed Brown               for (m = 0; m < 3; m++) {
884c4762a1bSJed Brown                 x[k][j][i][m] -= py[m];
885c4762a1bSJed Brown                 xnorm += PetscRealPart(x[k][j][i][m] * x[k][j][i][m]);
886c4762a1bSJed Brown               }
887c4762a1bSJed Brown               fnorm = PetscRealPart(pf[0] * pf[0] + pf[1] * pf[1] + pf[2] * pf[2]);
888c4762a1bSJed Brown               if (n == 0) fnorm0 = fnorm;
889c4762a1bSJed Brown               ynorm = PetscRealPart(py[0] * py[0] + py[1] * py[1] + py[2] * py[2]);
890c4762a1bSJed Brown               if (fnorm < atol * atol || fnorm < rtol * rtol * fnorm0 || ynorm < stol * stol * xnorm) break;
891c4762a1bSJed Brown             }
892c4762a1bSJed Brown           }
893c4762a1bSJed Brown         }
894c4762a1bSJed Brown       }
895c4762a1bSJed Brown     }
896c4762a1bSJed Brown   }
8979566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Xl, &x));
8989566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da, Xl, INSERT_VALUES, X));
8999566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da, Xl, INSERT_VALUES, X));
9009566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xl));
901c4762a1bSJed Brown   if (B) {
9029566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, Bl, &b));
9039566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(da, &Bl));
904c4762a1bSJed Brown   }
9059566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
9063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
907c4762a1bSJed Brown }
908c4762a1bSJed Brown 
909d71ae5a4SJacob Faibussowitsch PetscErrorCode FormCoordinates(DM da, AppCtx *user)
910d71ae5a4SJacob Faibussowitsch {
911c4762a1bSJed Brown   Vec           coords;
912c4762a1bSJed Brown   DM            cda;
913c4762a1bSJed Brown   PetscInt      mx, my, mz;
914c4762a1bSJed Brown   PetscInt      i, j, k, xs, ys, zs, xm, ym, zm;
915c4762a1bSJed Brown   CoordField ***x;
916c4762a1bSJed Brown 
917c4762a1bSJed Brown   PetscFunctionBegin;
9189566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
9199566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(cda, &coords));
9209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9219566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9229566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, coords, &x));
923c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
924c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
925c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
926*f4f49eeaSPierre Jolivet         PetscReal cx  = ((PetscReal)i) / ((PetscReal)(mx - 1));
927*f4f49eeaSPierre Jolivet         PetscReal cy  = ((PetscReal)j) / ((PetscReal)(my - 1));
928*f4f49eeaSPierre Jolivet         PetscReal cz  = ((PetscReal)k) / ((PetscReal)(mz - 1));
929c4762a1bSJed Brown         PetscReal rad = user->rad + cy * user->height;
930c4762a1bSJed Brown         PetscReal ang = (cx - 0.5) * user->arc;
931c4762a1bSJed Brown         x[k][j][i][0] = rad * PetscSinReal(ang);
932c4762a1bSJed Brown         x[k][j][i][1] = rad * PetscCosReal(ang) - (user->rad + 0.5 * user->height) * PetscCosReal(-0.5 * user->arc);
933c4762a1bSJed Brown         x[k][j][i][2] = user->width * (cz - 0.5);
934c4762a1bSJed Brown       }
935c4762a1bSJed Brown     }
936c4762a1bSJed Brown   }
9379566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, coords, &x));
9389566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinates(da, coords));
9399566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coords));
9403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
941c4762a1bSJed Brown }
942c4762a1bSJed Brown 
943d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialGuess(DM da, AppCtx *user, Vec X)
944d71ae5a4SJacob Faibussowitsch {
945c4762a1bSJed Brown   PetscInt i, j, k, xs, ys, zs, xm, ym, zm;
946c4762a1bSJed Brown   PetscInt mx, my, mz;
947c4762a1bSJed Brown   Field ***x;
948c4762a1bSJed Brown 
949c4762a1bSJed Brown   PetscFunctionBegin;
9509566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9519566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9529566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
953c4762a1bSJed Brown 
954c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
955c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
956c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
957c4762a1bSJed Brown         /* reference coordinates */
958*f4f49eeaSPierre Jolivet         PetscReal p_x = ((PetscReal)i) / ((PetscReal)(mx - 1));
959*f4f49eeaSPierre Jolivet         PetscReal p_y = ((PetscReal)j) / ((PetscReal)(my - 1));
960*f4f49eeaSPierre Jolivet         PetscReal p_z = ((PetscReal)k) / ((PetscReal)(mz - 1));
961c4762a1bSJed Brown         PetscReal o_x = p_x;
962c4762a1bSJed Brown         PetscReal o_y = p_y;
963c4762a1bSJed Brown         PetscReal o_z = p_z;
964c4762a1bSJed Brown         x[k][j][i][0] = o_x - p_x;
965c4762a1bSJed Brown         x[k][j][i][1] = o_y - p_y;
966c4762a1bSJed Brown         x[k][j][i][2] = o_z - p_z;
967c4762a1bSJed Brown       }
968c4762a1bSJed Brown     }
969c4762a1bSJed Brown   }
9709566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
9713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
972c4762a1bSJed Brown }
973c4762a1bSJed Brown 
974d71ae5a4SJacob Faibussowitsch PetscErrorCode FormRHS(DM da, AppCtx *user, Vec X)
975d71ae5a4SJacob Faibussowitsch {
976c4762a1bSJed Brown   PetscInt i, j, k, xs, ys, zs, xm, ym, zm;
977c4762a1bSJed Brown   PetscInt mx, my, mz;
978c4762a1bSJed Brown   Field ***x;
979c4762a1bSJed Brown 
980c4762a1bSJed Brown   PetscFunctionBegin;
9819566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
9829566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
9839566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
984c4762a1bSJed Brown 
985c4762a1bSJed Brown   for (k = zs; k < zs + zm; k++) {
986c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
987c4762a1bSJed Brown       for (i = xs; i < xs + xm; i++) {
988c4762a1bSJed Brown         x[k][j][i][0] = 0.;
989c4762a1bSJed Brown         x[k][j][i][1] = 0.;
990c4762a1bSJed Brown         x[k][j][i][2] = 0.;
991c4762a1bSJed Brown         if (i == (mx - 1) / 2 && j == (my - 1)) x[k][j][i][1] = user->ploading / (mz - 1);
992c4762a1bSJed Brown       }
993c4762a1bSJed Brown     }
994c4762a1bSJed Brown   }
9959566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
9963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
997c4762a1bSJed Brown }
998c4762a1bSJed Brown 
999d71ae5a4SJacob Faibussowitsch PetscErrorCode DisplayLine(SNES snes, Vec X)
1000d71ae5a4SJacob Faibussowitsch {
1001c4762a1bSJed Brown   PetscInt      r, i, j = 0, k = 0, xs, xm, ys, ym, zs, zm, mx, my, mz;
1002c4762a1bSJed Brown   Field      ***x;
1003c4762a1bSJed Brown   CoordField ***c;
1004c4762a1bSJed Brown   DM            da, cda;
1005c4762a1bSJed Brown   Vec           C;
1006c4762a1bSJed Brown   PetscMPIInt   size, rank;
1007c4762a1bSJed Brown 
1008c4762a1bSJed Brown   PetscFunctionBegin;
10099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
10119566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
10129566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
10139566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
10149566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &C));
1015c4762a1bSJed Brown   j = my / 2;
1016c4762a1bSJed Brown   k = mz / 2;
10179566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
10189566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
10199566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, C, &c));
1020c4762a1bSJed Brown   for (r = 0; r < size; r++) {
1021c4762a1bSJed Brown     if (rank == r) {
1022c4762a1bSJed Brown       if (j >= ys && j < ys + ym && k >= zs && k < zs + zm) {
1023c4762a1bSJed Brown         for (i = xs; i < xs + xm; i++) {
102463a3b9bcSJacob 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])));
1025c4762a1bSJed Brown         }
1026c4762a1bSJed Brown       }
1027c4762a1bSJed Brown     }
10289566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
1029c4762a1bSJed Brown   }
10309566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
10319566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, C, &c));
10323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1033c4762a1bSJed Brown }
1034c4762a1bSJed Brown 
1035c4762a1bSJed Brown /*TEST
1036c4762a1bSJed Brown 
1037c4762a1bSJed Brown    test:
1038c4762a1bSJed Brown       nsize: 2
1039c4762a1bSJed 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
1040c4762a1bSJed Brown       requires: !single
1041c4762a1bSJed Brown       timeoutfactor: 3
1042c4762a1bSJed Brown 
1043c4762a1bSJed Brown    test:
1044c4762a1bSJed Brown       suffix: 2
1045c4762a1bSJed 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
1046c4762a1bSJed Brown       requires: !single
1047c4762a1bSJed Brown 
1048c4762a1bSJed Brown    test:
1049c4762a1bSJed Brown       suffix: 3
1050de54d9edSStefano Zampini       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 -npc_sub_ksp_type preonly -npc_sub_pc_type lu
1051c4762a1bSJed Brown       requires: !single
1052c4762a1bSJed Brown 
1053c4762a1bSJed Brown TEST*/
1054