xref: /petsc/src/snes/tutorials/ex16.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Large-deformation Elasticity Buckling Example";
3c4762a1bSJed Brown /*T
4c4762a1bSJed Brown    Concepts: SNES^solving a system of nonlinear equations;
5c4762a1bSJed Brown    Concepts: DMDA^using distributed arrays;
6c4762a1bSJed Brown    Processors: n
7c4762a1bSJed Brown T*/
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /*F-----------------------------------------------------------------------
10c4762a1bSJed Brown 
11c4762a1bSJed Brown     This example solves the 3D large deformation elasticity problem
12c4762a1bSJed Brown 
13c4762a1bSJed Brown \begin{equation}
14c4762a1bSJed Brown  \int_{\Omega}F \cdot S : \nabla v d\Omega + \int_{\Omega} (loading)\mathbf{e}_y \cdot v d\Omega = 0
15c4762a1bSJed Brown \end{equation}
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     F is the deformation gradient, and S is the second Piola-Kirchhoff tensor from the Saint Venant-Kirchhoff model of
18c4762a1bSJed Brown     hyperelasticity.  \Omega is a (arc) angle subsection of a cylindrical shell of thickness (height), inner radius
19c4762a1bSJed Brown     (rad) and width (width).  The problem is discretized using Q1 finite elements on a logically structured grid.
20c4762a1bSJed Brown     Homogenous Dirichlet boundary conditions are applied at the centers of the ends of the sphere.
21c4762a1bSJed Brown 
22c4762a1bSJed Brown     This example is tunable with the following options:
23c4762a1bSJed Brown     -rad : the radius of the circle
24c4762a1bSJed Brown     -arc : set the angle of the arch represented
25c4762a1bSJed Brown     -loading : set the bulk loading (the mass)
26c4762a1bSJed Brown     -ploading : set the point loading at the top
27c4762a1bSJed Brown     -height : set the height of the arch
28c4762a1bSJed Brown     -width : set the width of the arch
29c4762a1bSJed Brown     -view_line : print initial and final offsets of the centerline of the
30c4762a1bSJed Brown                  beam along the x direction
31c4762a1bSJed Brown 
32c4762a1bSJed Brown     The material properties may be modified using either:
33c4762a1bSJed Brown     -mu      : the first lame' parameter
34c4762a1bSJed Brown     -lambda  : the second lame' parameter
35c4762a1bSJed Brown 
36c4762a1bSJed Brown     Or:
37c4762a1bSJed Brown     -young   : the Young's modulus
38c4762a1bSJed Brown     -poisson : the poisson ratio
39c4762a1bSJed Brown 
40c4762a1bSJed Brown     This example is meant to show the strain placed upon the nonlinear solvers when trying to "snap through" the arch
41c4762a1bSJed Brown     using the loading.  Under certain parameter regimes, the arch will invert under the load, and the number of Newton
42c4762a1bSJed Brown     steps will jump considerably.  Composed nonlinear solvers may be used to mitigate this difficulty.
43c4762a1bSJed Brown 
44c4762a1bSJed Brown     The initial setup follows the example in pg. 268 of "Nonlinear Finite Element Methods" by Peter Wriggers, but is a
45c4762a1bSJed Brown     3D extension.
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   ------------------------------------------------------------------------F*/
48c4762a1bSJed Brown 
49c4762a1bSJed Brown #include <petscsnes.h>
50c4762a1bSJed Brown #include <petscdm.h>
51c4762a1bSJed Brown #include <petscdmda.h>
52c4762a1bSJed Brown 
53c4762a1bSJed Brown #define QP0 0.2113248654051871
54c4762a1bSJed Brown #define QP1 0.7886751345948129
55c4762a1bSJed Brown #define NQ  2
56c4762a1bSJed Brown #define NB  2
57c4762a1bSJed Brown #define NEB 8
58c4762a1bSJed Brown #define NEQ 8
59c4762a1bSJed Brown #define NPB 24
60c4762a1bSJed Brown 
61c4762a1bSJed Brown #define NVALS NEB*NEQ
62c4762a1bSJed Brown const PetscReal pts[NQ] = {QP0,QP1};
63c4762a1bSJed Brown const PetscReal wts[NQ] = {0.5,0.5};
64c4762a1bSJed Brown 
65c4762a1bSJed Brown PetscScalar vals[NVALS];
66c4762a1bSJed Brown PetscScalar grad[3*NVALS];
67c4762a1bSJed Brown 
68c4762a1bSJed Brown typedef PetscScalar Field[3];
69c4762a1bSJed Brown typedef PetscScalar CoordField[3];
70c4762a1bSJed Brown 
71c4762a1bSJed Brown typedef PetscScalar JacField[9];
72c4762a1bSJed Brown 
73c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field***,Field***,void*);
74c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *,Field ***,Mat,Mat,void *);
75c4762a1bSJed Brown PetscErrorCode DisplayLine(SNES,Vec);
76c4762a1bSJed Brown PetscErrorCode FormElements();
77c4762a1bSJed Brown 
78c4762a1bSJed Brown typedef struct {
79c4762a1bSJed Brown   PetscReal loading;
80c4762a1bSJed Brown   PetscReal mu;
81c4762a1bSJed Brown   PetscReal lambda;
82c4762a1bSJed Brown   PetscReal rad;
83c4762a1bSJed Brown   PetscReal height;
84c4762a1bSJed Brown   PetscReal width;
85c4762a1bSJed Brown   PetscReal arc;
86c4762a1bSJed Brown   PetscReal ploading;
87c4762a1bSJed Brown } AppCtx;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown PetscErrorCode InitialGuess(DM,AppCtx *,Vec);
90c4762a1bSJed Brown PetscErrorCode FormRHS(DM,AppCtx *,Vec);
91c4762a1bSJed Brown PetscErrorCode FormCoordinates(DM,AppCtx *);
92c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown int main(int argc,char **argv)
95c4762a1bSJed Brown {
96c4762a1bSJed Brown   AppCtx         user;                /* user-defined work context */
97c4762a1bSJed Brown   PetscInt       mx,my,its;
98c4762a1bSJed Brown   PetscErrorCode ierr;
99c4762a1bSJed Brown   MPI_Comm       comm;
100c4762a1bSJed Brown   SNES           snes;
101c4762a1bSJed Brown   DM             da;
102c4762a1bSJed Brown   Vec            x,X,b;
103c4762a1bSJed Brown   PetscBool      youngflg,poissonflg,muflg,lambdaflg,view=PETSC_FALSE,viewline=PETSC_FALSE;
104c4762a1bSJed Brown   PetscReal      poisson=0.2,young=4e4;
105c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN] = "ex16.vts";
106c4762a1bSJed Brown   char           filename_def[PETSC_MAX_PATH_LEN] = "ex16_def.vts";
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormElements());
110c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(comm,&snes));
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes,(DM)da));
116c4762a1bSJed Brown 
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetNGS(snes,NonlinearGS,&user));
118c4762a1bSJed Brown 
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
120c4762a1bSJed Brown   user.loading     = 0.0;
121c4762a1bSJed Brown   user.arc         = PETSC_PI/3.;
122c4762a1bSJed Brown   user.mu          = 4.0;
123c4762a1bSJed Brown   user.lambda      = 1.0;
124c4762a1bSJed Brown   user.rad         = 100.0;
125c4762a1bSJed Brown   user.height      = 3.;
126c4762a1bSJed Brown   user.width       = 1.;
127c4762a1bSJed Brown   user.ploading    = -5e3;
128c4762a1bSJed Brown 
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-arc",&user.arc,NULL));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,&muflg));
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-lambda",&user.lambda,&lambdaflg));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-rad",&user.rad,NULL));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-height",&user.height,NULL));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-width",&user.width,NULL));
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-loading",&user.loading,NULL));
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ploading",&user.ploading,NULL));
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-poisson",&poisson,&poissonflg));
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-young",&young,&youngflg));
139c4762a1bSJed Brown   if ((youngflg || poissonflg) || !(muflg || lambdaflg)) {
140c4762a1bSJed Brown     /* set the lame' parameters based upon the poisson ratio and young's modulus */
141c4762a1bSJed Brown     user.lambda = poisson*young / ((1. + poisson)*(1. - 2.*poisson));
142c4762a1bSJed Brown     user.mu     = young/(2.*(1. + poisson));
143c4762a1bSJed Brown   }
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL));
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_line",&viewline,NULL));
146c4762a1bSJed Brown 
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,0,"x_disp"));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,1,"y_disp"));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,2,"z_disp"));
150c4762a1bSJed Brown 
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(da,&user));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user));
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormCoordinates(da,&user));
156c4762a1bSJed Brown 
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&b));
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(InitialGuess(da,&user,x));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormRHS(da,&user,b));
161c4762a1bSJed Brown 
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"lambda: %f mu: %f\n",(double)user.lambda,(double)user.mu));
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* show a cross-section of the initial state */
165c4762a1bSJed Brown   if (viewline) {
166*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DisplayLine(snes,x));
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* get the loaded configuration */
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,b,x));
171c4762a1bSJed Brown 
172*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"Number of SNES iterations = %D\n", its));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes,&X));
175c4762a1bSJed Brown   /* show a cross-section of the final state */
176c4762a1bSJed Brown   if (viewline) {
177*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DisplayLine(snes,X));
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   if (view) {
181c4762a1bSJed Brown     PetscViewer viewer;
182c4762a1bSJed Brown     Vec         coords;
183*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer));
184*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(x,viewer));
185*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinates(da,&coords));
187*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(coords,1.0,x));
188*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename_def,FILE_MODE_WRITE,&viewer));
189*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(x,viewer));
190*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
191c4762a1bSJed Brown   }
192c4762a1bSJed Brown 
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
197c4762a1bSJed Brown   ierr = PetscFinalize();
198c4762a1bSJed Brown   return ierr;
199c4762a1bSJed Brown }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown PetscInt OnBoundary(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz)
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   if ((i == 0 || i == mx-1) && j == my-1) return 1;
204c4762a1bSJed Brown   return 0;
205c4762a1bSJed Brown }
206c4762a1bSJed Brown 
207c4762a1bSJed Brown void BoundaryValue(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz,PetscScalar *val,AppCtx *user)
208c4762a1bSJed Brown {
209c4762a1bSJed Brown   /* reference coordinates */
210c4762a1bSJed Brown   PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx-1)));
211c4762a1bSJed Brown   PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my-1)));
212c4762a1bSJed Brown   PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz-1)));
213c4762a1bSJed Brown   PetscReal o_x = p_x;
214c4762a1bSJed Brown   PetscReal o_y = p_y;
215c4762a1bSJed Brown   PetscReal o_z = p_z;
216c4762a1bSJed Brown   val[0] = o_x - p_x;
217c4762a1bSJed Brown   val[1] = o_y - p_y;
218c4762a1bSJed Brown   val[2] = o_z - p_z;
219c4762a1bSJed Brown }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown void InvertTensor(PetscScalar *t, PetscScalar *ti,PetscReal *dett)
222c4762a1bSJed Brown {
223c4762a1bSJed Brown   const PetscScalar a = t[0];
224c4762a1bSJed Brown   const PetscScalar b = t[1];
225c4762a1bSJed Brown   const PetscScalar c = t[2];
226c4762a1bSJed Brown   const PetscScalar d = t[3];
227c4762a1bSJed Brown   const PetscScalar e = t[4];
228c4762a1bSJed Brown   const PetscScalar f = t[5];
229c4762a1bSJed Brown   const PetscScalar g = t[6];
230c4762a1bSJed Brown   const PetscScalar h = t[7];
231c4762a1bSJed Brown   const PetscScalar i = t[8];
232c4762a1bSJed Brown   const PetscReal   det = PetscRealPart(a*(e*i - f*h) - b*(i*d - f*g) + c*(d*h - e*g));
233c4762a1bSJed Brown   const PetscReal   di = 1. / det;
234c4762a1bSJed Brown   if (ti) {
235c4762a1bSJed Brown     const PetscScalar A = (e*i - f*h);
236c4762a1bSJed Brown     const PetscScalar B = -(d*i - f*g);
237c4762a1bSJed Brown     const PetscScalar C = (d*h - e*g);
238c4762a1bSJed Brown     const PetscScalar D = -(b*i - c*h);
239c4762a1bSJed Brown     const PetscScalar E = (a*i - c*g);
240c4762a1bSJed Brown     const PetscScalar F = -(a*h - b*g);
241c4762a1bSJed Brown     const PetscScalar G = (b*f - c*e);
242c4762a1bSJed Brown     const PetscScalar H = -(a*f - c*d);
243c4762a1bSJed Brown     const PetscScalar II = (a*e - b*d);
244c4762a1bSJed Brown     ti[0] = di*A;
245c4762a1bSJed Brown     ti[1] = di*D;
246c4762a1bSJed Brown     ti[2] = di*G;
247c4762a1bSJed Brown     ti[3] = di*B;
248c4762a1bSJed Brown     ti[4] = di*E;
249c4762a1bSJed Brown     ti[5] = di*H;
250c4762a1bSJed Brown     ti[6] = di*C;
251c4762a1bSJed Brown     ti[7] = di*F;
252c4762a1bSJed Brown     ti[8] = di*II;
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown   if (dett) *dett = det;
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown void TensorTensor(PetscScalar *a,PetscScalar *b,PetscScalar *c)
258c4762a1bSJed Brown {
259c4762a1bSJed Brown   PetscInt i,j,m;
260c4762a1bSJed Brown   for (i=0;i<3;i++) {
261c4762a1bSJed Brown     for (j=0;j<3;j++) {
262c4762a1bSJed Brown       c[i+3*j] = 0;
263c4762a1bSJed Brown       for (m=0;m<3;m++)
264c4762a1bSJed Brown         c[i+3*j] += a[m+3*j]*b[i+3*m];
265c4762a1bSJed Brown     }
266c4762a1bSJed Brown   }
267c4762a1bSJed Brown }
268c4762a1bSJed Brown 
269c4762a1bSJed Brown void TensorTransposeTensor(PetscScalar *a,PetscScalar *b,PetscScalar *c)
270c4762a1bSJed Brown {
271c4762a1bSJed Brown   PetscInt i,j,m;
272c4762a1bSJed Brown   for (i=0;i<3;i++) {
273c4762a1bSJed Brown     for (j=0;j<3;j++) {
274c4762a1bSJed Brown       c[i+3*j] = 0;
275c4762a1bSJed Brown       for (m=0;m<3;m++)
276c4762a1bSJed Brown         c[i+3*j] += a[3*m+j]*b[i+3*m];
277c4762a1bSJed Brown     }
278c4762a1bSJed Brown   }
279c4762a1bSJed Brown }
280c4762a1bSJed Brown 
281c4762a1bSJed Brown void TensorVector(PetscScalar *rot, PetscScalar *vec, PetscScalar *tvec)
282c4762a1bSJed Brown {
283c4762a1bSJed Brown   tvec[0] = rot[0]*vec[0] + rot[1]*vec[1] + rot[2]*vec[2];
284c4762a1bSJed Brown   tvec[1] = rot[3]*vec[0] + rot[4]*vec[1] + rot[5]*vec[2];
285c4762a1bSJed Brown   tvec[2] = rot[6]*vec[0] + rot[7]*vec[1] + rot[8]*vec[2];
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown void DeformationGradient(Field *ex,PetscInt qi,PetscInt qj,PetscInt qk,PetscScalar *invJ,PetscScalar *F)
289c4762a1bSJed Brown {
290c4762a1bSJed Brown   PetscInt ii,jj,kk,l;
291c4762a1bSJed Brown   for (l = 0; l < 9; l++) {
292c4762a1bSJed Brown     F[l] = 0.;
293c4762a1bSJed Brown   }
294c4762a1bSJed Brown   F[0] = 1.;
295c4762a1bSJed Brown   F[4] = 1.;
296c4762a1bSJed Brown   F[8] = 1.;
297c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
298c4762a1bSJed Brown   for (kk=0;kk<NB;kk++) {
299c4762a1bSJed Brown     for (jj=0;jj<NB;jj++) {
300c4762a1bSJed Brown       for (ii=0;ii<NB;ii++) {
301c4762a1bSJed Brown         PetscInt    idx = ii + jj*NB + kk*NB*NB;
302c4762a1bSJed Brown         PetscInt    bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
303c4762a1bSJed Brown         PetscScalar lgrad[3];
304c4762a1bSJed Brown         TensorVector(invJ,&grad[3*bidx],lgrad);
305c4762a1bSJed Brown         F[0] += lgrad[0]*ex[idx][0]; F[1] += lgrad[1]*ex[idx][0]; F[2] += lgrad[2]*ex[idx][0];
306c4762a1bSJed Brown         F[3] += lgrad[0]*ex[idx][1]; F[4] += lgrad[1]*ex[idx][1]; F[5] += lgrad[2]*ex[idx][1];
307c4762a1bSJed Brown         F[6] += lgrad[0]*ex[idx][2]; F[7] += lgrad[1]*ex[idx][2]; F[8] += lgrad[2]*ex[idx][2];
308c4762a1bSJed Brown       }
309c4762a1bSJed Brown     }
310c4762a1bSJed Brown   }
311c4762a1bSJed Brown }
312c4762a1bSJed Brown 
313c4762a1bSJed Brown void DeformationGradientJacobian(PetscInt qi,PetscInt qj,PetscInt qk,PetscInt ii,PetscInt jj,PetscInt kk,PetscInt fld,PetscScalar *invJ,PetscScalar *dF)
314c4762a1bSJed Brown {
315c4762a1bSJed Brown   PetscInt         l;
316c4762a1bSJed Brown   PetscScalar lgrad[3];
317c4762a1bSJed Brown   PetscInt    idx = ii + jj*NB + kk*NB*NB;
318c4762a1bSJed Brown   PetscInt    bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
319c4762a1bSJed Brown   for (l = 0; l < 9; l++) {
320c4762a1bSJed Brown     dF[l] = 0.;
321c4762a1bSJed Brown   }
322c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
323c4762a1bSJed Brown   TensorVector(invJ,&grad[3*bidx],lgrad);
324c4762a1bSJed Brown   dF[3*fld] = lgrad[0]; dF[3*fld + 1] = lgrad[1]; dF[3*fld + 2] = lgrad[2];
325c4762a1bSJed Brown }
326c4762a1bSJed Brown 
327c4762a1bSJed Brown void LagrangeGreenStrain(PetscScalar *F,PetscScalar *E)
328c4762a1bSJed Brown {
329c4762a1bSJed Brown   PetscInt i,j,m;
330c4762a1bSJed Brown   for (i=0;i<3;i++) {
331c4762a1bSJed Brown     for (j=0;j<3;j++) {
332c4762a1bSJed Brown       E[i+3*j] = 0;
333c4762a1bSJed Brown       for (m=0;m<3;m++)
334c4762a1bSJed Brown         E[i+3*j] += 0.5*F[3*m+j]*F[i+3*m];
335c4762a1bSJed Brown     }
336c4762a1bSJed Brown   }
337c4762a1bSJed Brown   for (i=0;i<3;i++) {
338c4762a1bSJed Brown     E[i+3*i] -= 0.5;
339c4762a1bSJed Brown   }
340c4762a1bSJed Brown }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown void SaintVenantKirchoff(PetscReal lambda,PetscReal mu,PetscScalar *F,PetscScalar *S)
343c4762a1bSJed Brown {
344c4762a1bSJed Brown   PetscInt    i,j;
345c4762a1bSJed Brown   PetscScalar E[9];
346c4762a1bSJed Brown   PetscScalar trE=0;
347c4762a1bSJed Brown   LagrangeGreenStrain(F,E);
348c4762a1bSJed Brown   for (i=0;i<3;i++) {
349c4762a1bSJed Brown     trE += E[i+3*i];
350c4762a1bSJed Brown   }
351c4762a1bSJed Brown   for (i=0;i<3;i++) {
352c4762a1bSJed Brown     for (j=0;j<3;j++) {
353c4762a1bSJed Brown       S[i+3*j] = 2.*mu*E[i+3*j];
354c4762a1bSJed Brown       if (i == j) {
355c4762a1bSJed Brown         S[i+3*j] += trE*lambda;
356c4762a1bSJed Brown       }
357c4762a1bSJed Brown     }
358c4762a1bSJed Brown   }
359c4762a1bSJed Brown }
360c4762a1bSJed Brown 
361c4762a1bSJed Brown void SaintVenantKirchoffJacobian(PetscReal lambda,PetscReal mu,PetscScalar *F,PetscScalar *dF,PetscScalar *dS)
362c4762a1bSJed Brown {
363c4762a1bSJed Brown   PetscScalar FtdF[9],dE[9];
364c4762a1bSJed Brown   PetscInt    i,j;
365c4762a1bSJed Brown   PetscScalar dtrE=0.;
366c4762a1bSJed Brown   TensorTransposeTensor(dF,F,dE);
367c4762a1bSJed Brown   TensorTransposeTensor(F,dF,FtdF);
368c4762a1bSJed Brown   for (i=0;i<9;i++) dE[i] += FtdF[i];
369c4762a1bSJed Brown   for (i=0;i<9;i++) dE[i] *= 0.5;
370c4762a1bSJed Brown   for (i=0;i<3;i++) {
371c4762a1bSJed Brown     dtrE += dE[i+3*i];
372c4762a1bSJed Brown   }
373c4762a1bSJed Brown   for (i=0;i<3;i++) {
374c4762a1bSJed Brown     for (j=0;j<3;j++) {
375c4762a1bSJed Brown       dS[i+3*j] = 2.*mu*dE[i+3*j];
376c4762a1bSJed Brown       if (i == j) {
377c4762a1bSJed Brown         dS[i+3*j] += dtrE*lambda;
378c4762a1bSJed Brown       }
379c4762a1bSJed Brown     }
380c4762a1bSJed Brown   }
381c4762a1bSJed Brown }
382c4762a1bSJed Brown 
383c4762a1bSJed Brown PetscErrorCode FormElements()
384c4762a1bSJed Brown {
385c4762a1bSJed Brown   PetscInt i,j,k,ii,jj,kk;
386c4762a1bSJed Brown   PetscReal bx,by,bz,dbx,dby,dbz;
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   PetscFunctionBegin;
389c4762a1bSJed Brown   /* construct the basis function values and derivatives */
390c4762a1bSJed Brown   for (k = 0; k < NB; k++) {
391c4762a1bSJed Brown     for (j = 0; j < NB; j++) {
392c4762a1bSJed Brown       for (i = 0; i < NB; i++) {
393c4762a1bSJed Brown         /* loop over the quadrature points */
394c4762a1bSJed Brown         for (kk = 0; kk < NQ; kk++) {
395c4762a1bSJed Brown           for (jj = 0; jj < NQ; jj++) {
396c4762a1bSJed Brown             for (ii = 0; ii < NQ; ii++) {
397c4762a1bSJed Brown               PetscInt idx = ii + NQ*jj + NQ*NQ*kk + NEQ*i + NEQ*NB*j + NEQ*NB*NB*k;
398c4762a1bSJed Brown               bx = pts[ii];
399c4762a1bSJed Brown               by = pts[jj];
400c4762a1bSJed Brown               bz = pts[kk];
401c4762a1bSJed Brown               dbx = 1.;
402c4762a1bSJed Brown               dby = 1.;
403c4762a1bSJed Brown               dbz = 1.;
404c4762a1bSJed Brown               if (i == 0) {bx = 1. - bx; dbx = -1;}
405c4762a1bSJed Brown               if (j == 0) {by = 1. - by; dby = -1;}
406c4762a1bSJed Brown               if (k == 0) {bz = 1. - bz; dbz = -1;}
407c4762a1bSJed Brown               vals[idx] = bx*by*bz;
408c4762a1bSJed Brown               grad[3*idx + 0] = dbx*by*bz;
409c4762a1bSJed Brown               grad[3*idx + 1] = dby*bx*bz;
410c4762a1bSJed Brown               grad[3*idx + 2] = dbz*bx*by;
411c4762a1bSJed Brown             }
412c4762a1bSJed Brown           }
413c4762a1bSJed Brown         }
414c4762a1bSJed Brown       }
415c4762a1bSJed Brown     }
416c4762a1bSJed Brown   }
417c4762a1bSJed Brown   PetscFunctionReturn(0);
418c4762a1bSJed Brown }
419c4762a1bSJed Brown 
420c4762a1bSJed Brown void GatherElementData(PetscInt mx,PetscInt my,PetscInt mz,Field ***x,CoordField ***c,PetscInt i,PetscInt j,PetscInt k,Field *ex,CoordField *ec,AppCtx *user)
421c4762a1bSJed Brown {
422c4762a1bSJed Brown   PetscInt m;
423c4762a1bSJed Brown   PetscInt ii,jj,kk;
424c4762a1bSJed Brown   /* gather the data -- loop over element unknowns */
425c4762a1bSJed Brown   for (kk=0;kk<NB;kk++) {
426c4762a1bSJed Brown     for (jj=0;jj<NB;jj++) {
427c4762a1bSJed Brown       for (ii=0;ii<NB;ii++) {
428c4762a1bSJed Brown         PetscInt idx = ii + jj*NB + kk*NB*NB;
429c4762a1bSJed Brown         /* decouple the boundary nodes for the displacement variables */
430c4762a1bSJed Brown         if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz)) {
431c4762a1bSJed Brown           BoundaryValue(i+ii,j+jj,k+kk,mx,my,mz,ex[idx],user);
432c4762a1bSJed Brown         } else {
433c4762a1bSJed Brown           for (m=0;m<3;m++) {
434c4762a1bSJed Brown             ex[idx][m] = x[k+kk][j+jj][i+ii][m];
435c4762a1bSJed Brown           }
436c4762a1bSJed Brown         }
437c4762a1bSJed Brown         for (m=0;m<3;m++) {
438c4762a1bSJed Brown           ec[idx][m] = c[k+kk][j+jj][i+ii][m];
439c4762a1bSJed Brown         }
440c4762a1bSJed Brown       }
441c4762a1bSJed Brown     }
442c4762a1bSJed Brown   }
443c4762a1bSJed Brown }
444c4762a1bSJed Brown 
445c4762a1bSJed Brown void QuadraturePointGeometricJacobian(CoordField *ec,PetscInt qi,PetscInt qj,PetscInt qk, PetscScalar *J)
446c4762a1bSJed Brown {
447c4762a1bSJed Brown   PetscInt ii,jj,kk;
448c4762a1bSJed Brown   /* construct the gradient at the given quadrature point named by i,j,k */
449c4762a1bSJed Brown   for (ii = 0; ii < 9; ii++) {
450c4762a1bSJed Brown     J[ii] = 0;
451c4762a1bSJed Brown   }
452c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
453c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
454c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
455c4762a1bSJed Brown         PetscInt idx = ii + jj*NB + kk*NB*NB;
456c4762a1bSJed Brown         PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
457c4762a1bSJed Brown         J[0] += grad[3*bidx + 0]*ec[idx][0]; J[1] += grad[3*bidx + 1]*ec[idx][0]; J[2] += grad[3*bidx + 2]*ec[idx][0];
458c4762a1bSJed Brown         J[3] += grad[3*bidx + 0]*ec[idx][1]; J[4] += grad[3*bidx + 1]*ec[idx][1]; J[5] += grad[3*bidx + 2]*ec[idx][1];
459c4762a1bSJed Brown         J[6] += grad[3*bidx + 0]*ec[idx][2]; J[7] += grad[3*bidx + 1]*ec[idx][2]; J[8] += grad[3*bidx + 2]*ec[idx][2];
460c4762a1bSJed Brown       }
461c4762a1bSJed Brown     }
462c4762a1bSJed Brown   }
463c4762a1bSJed Brown }
464c4762a1bSJed Brown 
465c4762a1bSJed Brown void FormElementJacobian(Field *ex,CoordField *ec,Field *ef,PetscScalar *ej,AppCtx *user)
466c4762a1bSJed Brown {
467c4762a1bSJed Brown   PetscReal   vol;
468c4762a1bSJed Brown   PetscScalar J[9];
469c4762a1bSJed Brown   PetscScalar invJ[9];
470c4762a1bSJed Brown   PetscScalar F[9],S[9],dF[9],dS[9],dFS[9],FdS[9],FS[9];
471c4762a1bSJed Brown   PetscReal   scl;
472c4762a1bSJed Brown   PetscInt    i,j,k,l,ii,jj,kk,ll,qi,qj,qk,m;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   if (ej) for (i = 0; i < NPB*NPB; i++) ej[i] = 0.;
475c4762a1bSJed Brown   if (ef) for (i = 0; i < NEB; i++) {ef[i][0] = 0.;ef[i][1] = 0.;ef[i][2] = 0.;}
476c4762a1bSJed Brown   /* loop over quadrature */
477c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
478c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
479c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
480c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec,qi,qj,qk,J);
481c4762a1bSJed Brown         InvertTensor(J,invJ,&vol);
482c4762a1bSJed Brown         scl = vol*wts[qi]*wts[qj]*wts[qk];
483c4762a1bSJed Brown         DeformationGradient(ex,qi,qj,qk,invJ,F);
484c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda,user->mu,F,S);
485c4762a1bSJed Brown         /* form the function */
486c4762a1bSJed Brown         if (ef) {
487c4762a1bSJed Brown           TensorTensor(F,S,FS);
488c4762a1bSJed Brown           for (kk=0;kk<NB;kk++) {
489c4762a1bSJed Brown             for (jj=0;jj<NB;jj++) {
490c4762a1bSJed Brown               for (ii=0;ii<NB;ii++) {
491c4762a1bSJed Brown                 PetscInt idx = ii + jj*NB + kk*NB*NB;
492c4762a1bSJed Brown                 PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
493c4762a1bSJed Brown                 PetscScalar lgrad[3];
494c4762a1bSJed Brown                 TensorVector(invJ,&grad[3*bidx],lgrad);
495c4762a1bSJed Brown                 /* mu*F : grad phi_{u,v,w} */
496c4762a1bSJed Brown                 for (m=0;m<3;m++) {
497c4762a1bSJed Brown                   ef[idx][m] += scl*
498c4762a1bSJed Brown                     (lgrad[0]*FS[3*m + 0] + lgrad[1]*FS[3*m + 1] + lgrad[2]*FS[3*m + 2]);
499c4762a1bSJed Brown                 }
500c4762a1bSJed Brown                 ef[idx][1] -= scl*user->loading*vals[bidx];
501c4762a1bSJed Brown               }
502c4762a1bSJed Brown             }
503c4762a1bSJed Brown           }
504c4762a1bSJed Brown         }
505c4762a1bSJed Brown         /* form the jacobian */
506c4762a1bSJed Brown         if (ej) {
507c4762a1bSJed Brown           /* loop over trialfunctions */
508c4762a1bSJed Brown           for (k=0;k<NB;k++) {
509c4762a1bSJed Brown             for (j=0;j<NB;j++) {
510c4762a1bSJed Brown               for (i=0;i<NB;i++) {
511c4762a1bSJed Brown                 for (l=0;l<3;l++) {
512c4762a1bSJed Brown                   PetscInt tridx = l + 3*(i + j*NB + k*NB*NB);
513c4762a1bSJed Brown                   DeformationGradientJacobian(qi,qj,qk,i,j,k,l,invJ,dF);
514c4762a1bSJed Brown                   SaintVenantKirchoffJacobian(user->lambda,user->mu,F,dF,dS);
515c4762a1bSJed Brown                   TensorTensor(dF,S,dFS);
516c4762a1bSJed Brown                   TensorTensor(F,dS,FdS);
517c4762a1bSJed Brown                   for (m=0;m<9;m++) dFS[m] += FdS[m];
518c4762a1bSJed Brown                   /* loop over testfunctions */
519c4762a1bSJed Brown                   for (kk=0;kk<NB;kk++) {
520c4762a1bSJed Brown                     for (jj=0;jj<NB;jj++) {
521c4762a1bSJed Brown                       for (ii=0;ii<NB;ii++) {
522c4762a1bSJed Brown                         PetscInt idx = ii + jj*NB + kk*NB*NB;
523c4762a1bSJed Brown                         PetscInt bidx = 8*idx + qi + NQ*qj + NQ*NQ*qk;
524c4762a1bSJed Brown                         PetscScalar lgrad[3];
525c4762a1bSJed Brown                         TensorVector(invJ,&grad[3*bidx],lgrad);
526c4762a1bSJed Brown                         for (ll=0; ll<3;ll++) {
527c4762a1bSJed Brown                           PetscInt teidx = ll + 3*(ii + jj*NB + kk*NB*NB);
528c4762a1bSJed Brown                           ej[teidx + NPB*tridx] += scl*
529c4762a1bSJed Brown                             (lgrad[0]*dFS[3*ll + 0] + lgrad[1]*dFS[3*ll + 1] + lgrad[2]*dFS[3*ll+2]);
530c4762a1bSJed Brown                         }
531c4762a1bSJed Brown                       }
532c4762a1bSJed Brown                     }
533c4762a1bSJed Brown                   } /* end of testfunctions */
534c4762a1bSJed Brown                 }
535c4762a1bSJed Brown               }
536c4762a1bSJed Brown             }
537c4762a1bSJed Brown           } /* end of trialfunctions */
538c4762a1bSJed Brown         }
539c4762a1bSJed Brown       }
540c4762a1bSJed Brown     }
541c4762a1bSJed Brown   } /* end of quadrature points */
542c4762a1bSJed Brown }
543c4762a1bSJed Brown 
544c4762a1bSJed Brown void FormPBJacobian(PetscInt i,PetscInt j,PetscInt k,Field *ex,CoordField *ec,Field *ef,PetscScalar *ej,AppCtx *user)
545c4762a1bSJed Brown {
546c4762a1bSJed Brown   PetscReal   vol;
547c4762a1bSJed Brown   PetscScalar J[9];
548c4762a1bSJed Brown   PetscScalar invJ[9];
549c4762a1bSJed Brown   PetscScalar F[9],S[9],dF[9],dS[9],dFS[9],FdS[9],FS[9];
550c4762a1bSJed Brown   PetscReal   scl;
551c4762a1bSJed Brown   PetscInt    l,ll,qi,qj,qk,m;
552c4762a1bSJed Brown   PetscInt idx = i + j*NB + k*NB*NB;
553c4762a1bSJed Brown   PetscScalar lgrad[3];
554c4762a1bSJed Brown 
555c4762a1bSJed Brown   if (ej) for (l = 0; l < 9; l++) ej[l] = 0.;
556c4762a1bSJed Brown   if (ef) for (l = 0; l < 1; l++) {ef[l][0] = 0.;ef[l][1] = 0.;ef[l][2] = 0.;}
557c4762a1bSJed Brown   /* loop over quadrature */
558c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
559c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
560c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
561c4762a1bSJed Brown         PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
562c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec,qi,qj,qk,J);
563c4762a1bSJed Brown         InvertTensor(J,invJ,&vol);
564c4762a1bSJed Brown         TensorVector(invJ,&grad[3*bidx],lgrad);
565c4762a1bSJed Brown         scl = vol*wts[qi]*wts[qj]*wts[qk];
566c4762a1bSJed Brown         DeformationGradient(ex,qi,qj,qk,invJ,F);
567c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda,user->mu,F,S);
568c4762a1bSJed Brown         /* form the function */
569c4762a1bSJed Brown         if (ef) {
570c4762a1bSJed Brown           TensorTensor(F,S,FS);
571c4762a1bSJed Brown           for (m=0;m<3;m++) {
572c4762a1bSJed Brown             ef[0][m] += scl*
573c4762a1bSJed Brown               (lgrad[0]*FS[3*m + 0] + lgrad[1]*FS[3*m + 1] + lgrad[2]*FS[3*m + 2]);
574c4762a1bSJed Brown           }
575c4762a1bSJed Brown           ef[0][1] -= scl*user->loading*vals[bidx];
576c4762a1bSJed Brown         }
577c4762a1bSJed Brown         /* form the jacobian */
578c4762a1bSJed Brown         if (ej) {
579c4762a1bSJed Brown           for (l=0;l<3;l++) {
580c4762a1bSJed Brown             DeformationGradientJacobian(qi,qj,qk,i,j,k,l,invJ,dF);
581c4762a1bSJed Brown             SaintVenantKirchoffJacobian(user->lambda,user->mu,F,dF,dS);
582c4762a1bSJed Brown             TensorTensor(dF,S,dFS);
583c4762a1bSJed Brown             TensorTensor(F,dS,FdS);
584c4762a1bSJed Brown             for (m=0;m<9;m++) dFS[m] += FdS[m];
585c4762a1bSJed Brown             for (ll=0; ll<3;ll++) {
586c4762a1bSJed Brown               ej[ll + 3*l] += scl*
587c4762a1bSJed Brown                 (lgrad[0]*dFS[3*ll + 0] + lgrad[1]*dFS[3*ll + 1] + lgrad[2]*dFS[3*ll+2]);
588c4762a1bSJed Brown             }
589c4762a1bSJed Brown           }
590c4762a1bSJed Brown         }
591c4762a1bSJed Brown       }
592c4762a1bSJed Brown     } /* end of quadrature points */
593c4762a1bSJed Brown   }
594c4762a1bSJed Brown }
595c4762a1bSJed Brown 
596c4762a1bSJed Brown void ApplyBCsElement(PetscInt mx,PetscInt my, PetscInt mz, PetscInt i, PetscInt j, PetscInt k,PetscScalar *jacobian)
597c4762a1bSJed Brown {
598c4762a1bSJed Brown   PetscInt ii,jj,kk,ll,ei,ej,ek,el;
599c4762a1bSJed Brown   for (kk=0;kk<NB;kk++) {
600c4762a1bSJed Brown     for (jj=0;jj<NB;jj++) {
601c4762a1bSJed Brown       for (ii=0;ii<NB;ii++) {
602c4762a1bSJed Brown         for (ll = 0;ll<3;ll++) {
603c4762a1bSJed Brown           PetscInt tridx = ll + 3*(ii + jj*NB + kk*NB*NB);
604c4762a1bSJed Brown           for (ek=0;ek<NB;ek++) {
605c4762a1bSJed Brown             for (ej=0;ej<NB;ej++) {
606c4762a1bSJed Brown               for (ei=0;ei<NB;ei++) {
607c4762a1bSJed Brown                 for (el=0;el<3;el++) {
608c4762a1bSJed Brown                   if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz) || OnBoundary(i+ei,j+ej,k+ek,mx,my,mz)) {
609c4762a1bSJed Brown                     PetscInt teidx = el + 3*(ei + ej*NB + ek*NB*NB);
610c4762a1bSJed Brown                     if (teidx == tridx) {
611c4762a1bSJed Brown                       jacobian[tridx + NPB*teidx] = 1.;
612c4762a1bSJed Brown                     } else {
613c4762a1bSJed Brown                       jacobian[tridx + NPB*teidx] = 0.;
614c4762a1bSJed Brown                     }
615c4762a1bSJed Brown                   }
616c4762a1bSJed Brown                 }
617c4762a1bSJed Brown               }
618c4762a1bSJed Brown             }
619c4762a1bSJed Brown           }
620c4762a1bSJed Brown         }
621c4762a1bSJed Brown       }
622c4762a1bSJed Brown     }
623c4762a1bSJed Brown   }
624c4762a1bSJed Brown }
625c4762a1bSJed Brown 
626c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,Field ***x,Mat jacpre,Mat jac,void *ptr)
627c4762a1bSJed Brown {
628c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
629c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
630c4762a1bSJed Brown   PetscInt       i,j,k,m,l;
631c4762a1bSJed Brown   PetscInt       ii,jj,kk;
632c4762a1bSJed Brown   PetscScalar    ej[NPB*NPB];
633c4762a1bSJed Brown   PetscScalar    vals[NPB*NPB];
634c4762a1bSJed Brown   Field          ex[NEB];
635c4762a1bSJed Brown   CoordField     ec[NEB];
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   PetscInt       xs=info->xs,ys=info->ys,zs=info->zs;
638c4762a1bSJed Brown   PetscInt       xm=info->xm,ym=info->ym,zm=info->zm;
639c4762a1bSJed Brown   PetscInt       xes,yes,zes,xee,yee,zee;
640c4762a1bSJed Brown   PetscInt       mx=info->mx,my=info->my,mz=info->mz;
641c4762a1bSJed Brown   DM             cda;
642c4762a1bSJed Brown   CoordField     ***c;
643c4762a1bSJed Brown   Vec            C;
644c4762a1bSJed Brown   PetscInt       nrows;
645c4762a1bSJed Brown   MatStencil     col[NPB],row[NPB];
646c4762a1bSJed Brown   PetscScalar    v[9];
647c4762a1bSJed Brown 
648c4762a1bSJed Brown   PetscFunctionBegin;
649*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(info->da,&cda));
650*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(info->da,&C));
651*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
652*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(jac,0.0));
653c4762a1bSJed Brown 
654c4762a1bSJed Brown   xes = xs;
655c4762a1bSJed Brown   yes = ys;
656c4762a1bSJed Brown   zes = zs;
657c4762a1bSJed Brown   xee = xs+xm;
658c4762a1bSJed Brown   yee = ys+ym;
659c4762a1bSJed Brown   zee = zs+zm;
660c4762a1bSJed Brown   if (xs > 0) xes = xs-1;
661c4762a1bSJed Brown   if (ys > 0) yes = ys-1;
662c4762a1bSJed Brown   if (zs > 0) zes = zs-1;
663c4762a1bSJed Brown   if (xs+xm == mx) xee = xs+xm-1;
664c4762a1bSJed Brown   if (ys+ym == my) yee = ys+ym-1;
665c4762a1bSJed Brown   if (zs+zm == mz) zee = zs+zm-1;
666c4762a1bSJed Brown   for (k=zes; k<zee; k++) {
667c4762a1bSJed Brown     for (j=yes; j<yee; j++) {
668c4762a1bSJed Brown       for (i=xes; i<xee; i++) {
669c4762a1bSJed Brown         GatherElementData(mx,my,mz,x,c,i,j,k,ex,ec,user);
670c4762a1bSJed Brown         FormElementJacobian(ex,ec,NULL,ej,user);
671c4762a1bSJed Brown         ApplyBCsElement(mx,my,mz,i,j,k,ej);
672c4762a1bSJed Brown         nrows = 0.;
673c4762a1bSJed Brown         for (kk=0;kk<NB;kk++) {
674c4762a1bSJed Brown           for (jj=0;jj<NB;jj++) {
675c4762a1bSJed Brown             for (ii=0;ii<NB;ii++) {
676c4762a1bSJed Brown               PetscInt idx = ii + jj*2 + kk*4;
677c4762a1bSJed Brown               for (m=0;m<3;m++) {
678c4762a1bSJed Brown                 col[3*idx+m].i = i+ii;
679c4762a1bSJed Brown                 col[3*idx+m].j = j+jj;
680c4762a1bSJed Brown                 col[3*idx+m].k = k+kk;
681c4762a1bSJed Brown                 col[3*idx+m].c = m;
682c4762a1bSJed Brown                 if (i+ii >= xs && i+ii < xm+xs && j+jj >= ys && j+jj < ys+ym && k+kk >= zs && k+kk < zs+zm) {
683c4762a1bSJed Brown                   row[nrows].i = i+ii;
684c4762a1bSJed Brown                   row[nrows].j = j+jj;
685c4762a1bSJed Brown                   row[nrows].k = k+kk;
686c4762a1bSJed Brown                   row[nrows].c = m;
687c4762a1bSJed Brown                   for (l=0;l<NPB;l++) vals[NPB*nrows + l] = ej[NPB*(3*idx+m) + l];
688c4762a1bSJed Brown                   nrows++;
689c4762a1bSJed Brown                 }
690c4762a1bSJed Brown               }
691c4762a1bSJed Brown             }
692c4762a1bSJed Brown           }
693c4762a1bSJed Brown         }
694*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(jac,nrows,row,NPB,col,vals,ADD_VALUES));
695c4762a1bSJed Brown       }
696c4762a1bSJed Brown     }
697c4762a1bSJed Brown   }
698c4762a1bSJed Brown 
699*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FLUSH_ASSEMBLY));
700*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FLUSH_ASSEMBLY));
701c4762a1bSJed Brown 
702c4762a1bSJed Brown   /* set the diagonal */
703c4762a1bSJed Brown   v[0] = 1.;v[1] = 0.;v[2] = 0.;v[3] = 0.;v[4] = 1.;v[5] = 0.;v[6] = 0.;v[7] = 0.;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           }
714*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,3,col,3,col,v,INSERT_VALUES));
715c4762a1bSJed Brown         }
716c4762a1bSJed Brown       }
717c4762a1bSJed Brown     }
718c4762a1bSJed Brown   }
719c4762a1bSJed Brown 
720*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
721*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
722c4762a1bSJed Brown 
723*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
724c4762a1bSJed Brown   PetscFunctionReturn(0);
725c4762a1bSJed Brown }
726c4762a1bSJed Brown 
727c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field ***x,Field ***f,void *ptr)
728c4762a1bSJed Brown {
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;
747*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(info->da,&cda));
748*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(info->da,&C));
749*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
750*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(info->da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
751*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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++) {
757c4762a1bSJed Brown         for (l=0;l<3;l++) {
758c4762a1bSJed Brown           f[k][j][i][l] = 0.;
759c4762a1bSJed Brown         }
760c4762a1bSJed Brown       }
761c4762a1bSJed Brown     }
762c4762a1bSJed Brown   }
763c4762a1bSJed Brown   /* element starts and ends */
764c4762a1bSJed Brown   xes = xs;
765c4762a1bSJed Brown   yes = ys;
766c4762a1bSJed Brown   zes = zs;
767c4762a1bSJed Brown   xee = xs+xm;
768c4762a1bSJed Brown   yee = ys+ym;
769c4762a1bSJed Brown   zee = zs+zm;
770c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
771c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
772c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
773c4762a1bSJed Brown   if (xs+xm == mx) xee = xs+xm-1;
774c4762a1bSJed Brown   if (ys+ym == my) yee = ys+ym-1;
775c4762a1bSJed Brown   if (zs+zm == mz) zee = zs+zm-1;
776c4762a1bSJed Brown   for (k=zes; k<zee; k++) {
777c4762a1bSJed Brown     for (j=yes; j<yee; j++) {
778c4762a1bSJed Brown       for (i=xes; i<xee; i++) {
779c4762a1bSJed Brown         GatherElementData(mx,my,mz,x,c,i,j,k,ex,ec,user);
780c4762a1bSJed Brown         FormElementJacobian(ex,ec,ef,NULL,user);
781c4762a1bSJed Brown         /* put this element's additions into the residuals */
782c4762a1bSJed Brown         for (kk=0;kk<NB;kk++) {
783c4762a1bSJed Brown           for (jj=0;jj<NB;jj++) {
784c4762a1bSJed Brown             for (ii=0;ii<NB;ii++) {
785c4762a1bSJed Brown               PetscInt idx = ii + jj*NB + kk*NB*NB;
786c4762a1bSJed Brown               if (k+kk >= zs && j+jj >= ys && i+ii >= xs && k+kk < zs+zm && j+jj < ys+ym && i+ii < xs+xm) {
787c4762a1bSJed Brown                 if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz)) {
788c4762a1bSJed Brown                   for (l=0;l<3;l++)
789c4762a1bSJed Brown                     f[k+kk][j+jj][i+ii][l] = x[k+kk][j+jj][i+ii][l] - ex[idx][l];
790c4762a1bSJed Brown                 } else {
791c4762a1bSJed Brown                   for (l=0;l<3;l++)
792c4762a1bSJed Brown                     f[k+kk][j+jj][i+ii][l] += ef[idx][l];
793c4762a1bSJed Brown                 }
794c4762a1bSJed Brown               }
795c4762a1bSJed Brown             }
796c4762a1bSJed Brown           }
797c4762a1bSJed Brown         }
798c4762a1bSJed Brown       }
799c4762a1bSJed Brown     }
800c4762a1bSJed Brown   }
801*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
802c4762a1bSJed Brown   PetscFunctionReturn(0);
803c4762a1bSJed Brown }
804c4762a1bSJed Brown 
805c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X,Vec B,void *ptr)
806c4762a1bSJed Brown {
807c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
808c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
809c4762a1bSJed Brown   PetscInt       i,j,k,l,m,n,s;
810c4762a1bSJed Brown   PetscInt       pi,pj,pk;
811c4762a1bSJed Brown   Field          ef[1];
812c4762a1bSJed Brown   Field          ex[8];
813c4762a1bSJed Brown   PetscScalar    ej[9];
814c4762a1bSJed Brown   CoordField     ec[8];
815c4762a1bSJed Brown   PetscScalar    pjac[9],pjinv[9];
816c4762a1bSJed Brown   PetscScalar    pf[3],py[3];
817c4762a1bSJed Brown   PetscInt       xs,ys,zs;
818c4762a1bSJed Brown   PetscInt       xm,ym,zm;
819c4762a1bSJed Brown   PetscInt       mx,my,mz;
820c4762a1bSJed Brown   DM             cda;
821c4762a1bSJed Brown   CoordField     ***c;
822c4762a1bSJed Brown   Vec            C;
823c4762a1bSJed Brown   DM             da;
824c4762a1bSJed Brown   Vec            Xl,Bl;
825c4762a1bSJed Brown   Field          ***x,***b;
826c4762a1bSJed Brown   PetscInt       sweeps,its;
827c4762a1bSJed Brown   PetscReal      atol,rtol,stol;
828c4762a1bSJed Brown   PetscReal      fnorm0 = 0.0,fnorm,ynorm,xnorm = 0.0;
829c4762a1bSJed Brown 
830c4762a1bSJed Brown   PetscFunctionBegin;
831*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESNGSGetSweeps(snes,&sweeps));
832*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its));
833c4762a1bSJed Brown 
834*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
835*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xl));
836c4762a1bSJed Brown   if (B) {
837*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(da,&Bl));
838c4762a1bSJed Brown   }
839*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xl));
840*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xl));
841c4762a1bSJed Brown   if (B) {
842*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(da,B,INSERT_VALUES,Bl));
843*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(da,B,INSERT_VALUES,Bl));
844c4762a1bSJed Brown   }
845*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xl,&x));
846*5f80ce2aSJacob Faibussowitsch   if (B) CHKERRQ(DMDAVecGetArray(da,Bl,&b));
847c4762a1bSJed Brown 
848*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
849*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(da,&C));
850*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
851*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
852*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
853c4762a1bSJed Brown 
854c4762a1bSJed Brown   for (s=0;s<sweeps;s++) {
855c4762a1bSJed Brown     for (k=zs; k<zs+zm; k++) {
856c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
857c4762a1bSJed Brown         for (i=xs; i<xs+xm; i++) {
858c4762a1bSJed Brown           if (OnBoundary(i,j,k,mx,my,mz)) {
859c4762a1bSJed Brown             BoundaryValue(i,j,k,mx,my,mz,x[k][j][i],user);
860c4762a1bSJed Brown           } else {
861c4762a1bSJed Brown             for (n=0;n<its;n++) {
862c4762a1bSJed Brown               for (m=0;m<9;m++) pjac[m] = 0.;
863c4762a1bSJed Brown               for (m=0;m<3;m++) pf[m] = 0.;
864c4762a1bSJed Brown               /* gather the elements for this point */
865c4762a1bSJed Brown               for (pk=-1; pk<1; pk++) {
866c4762a1bSJed Brown                 for (pj=-1; pj<1; pj++) {
867c4762a1bSJed Brown                   for (pi=-1; pi<1; pi++) {
868c4762a1bSJed Brown                     /* check that this element exists */
869c4762a1bSJed Brown                     if (i+pi >= 0 && i+pi < mx-1 && j+pj >= 0 && j+pj < my-1 && k+pk >= 0 && k+pk < mz-1) {
870c4762a1bSJed Brown                       /* create the element function and jacobian */
871c4762a1bSJed Brown                       GatherElementData(mx,my,mz,x,c,i+pi,j+pj,k+pk,ex,ec,user);
872c4762a1bSJed Brown                       FormPBJacobian(-pi,-pj,-pk,ex,ec,ef,ej,user);
873c4762a1bSJed Brown                       /* extract the point named by i,j,k from the whole element jacobian and function */
874c4762a1bSJed Brown                       for (l=0;l<3;l++) {
875c4762a1bSJed Brown                         pf[l] += ef[0][l];
876c4762a1bSJed Brown                         for (m=0;m<3;m++) {
877c4762a1bSJed Brown                           pjac[3*m+l] += ej[3*m+l];
878c4762a1bSJed Brown                         }
879c4762a1bSJed Brown                       }
880c4762a1bSJed Brown                     }
881c4762a1bSJed Brown                   }
882c4762a1bSJed Brown                 }
883c4762a1bSJed Brown               }
884c4762a1bSJed Brown               /* invert */
885c4762a1bSJed Brown               InvertTensor(pjac,pjinv,NULL);
886c4762a1bSJed Brown               /* apply */
887c4762a1bSJed Brown               if (B) for (m=0;m<3;m++) {
888c4762a1bSJed Brown                   pf[m] -= b[k][j][i][m];
889c4762a1bSJed Brown                 }
890c4762a1bSJed Brown               TensorVector(pjinv,pf,py);
891c4762a1bSJed Brown               xnorm=0.;
892c4762a1bSJed Brown               for (m=0;m<3;m++) {
893c4762a1bSJed Brown                 x[k][j][i][m] -= py[m];
894c4762a1bSJed Brown                 xnorm += PetscRealPart(x[k][j][i][m]*x[k][j][i][m]);
895c4762a1bSJed Brown               }
896c4762a1bSJed Brown               fnorm = PetscRealPart(pf[0]*pf[0]+pf[1]*pf[1]+pf[2]*pf[2]);
897c4762a1bSJed Brown               if (n==0) fnorm0 = fnorm;
898c4762a1bSJed Brown               ynorm = PetscRealPart(py[0]*py[0]+py[1]*py[1]+py[2]*py[2]);
899c4762a1bSJed Brown               if (fnorm < atol*atol || fnorm < rtol*rtol*fnorm0 || ynorm < stol*stol*xnorm) break;
900c4762a1bSJed Brown             }
901c4762a1bSJed Brown           }
902c4762a1bSJed Brown         }
903c4762a1bSJed Brown       }
904c4762a1bSJed Brown     }
905c4762a1bSJed Brown   }
906*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xl,&x));
907*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(da,Xl,INSERT_VALUES,X));
908*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(da,Xl,INSERT_VALUES,X));
909*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xl));
910c4762a1bSJed Brown   if (B) {
911*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,Bl,&b));
912*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(da,&Bl));
913c4762a1bSJed Brown   }
914*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
915c4762a1bSJed Brown   PetscFunctionReturn(0);
916c4762a1bSJed Brown }
917c4762a1bSJed Brown 
918c4762a1bSJed Brown PetscErrorCode FormCoordinates(DM da,AppCtx *user)
919c4762a1bSJed Brown {
920c4762a1bSJed Brown   Vec            coords;
921c4762a1bSJed Brown   DM             cda;
922c4762a1bSJed Brown   PetscInt       mx,my,mz;
923c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,zs,xm,ym,zm;
924c4762a1bSJed Brown   CoordField     ***x;
925c4762a1bSJed Brown 
926c4762a1bSJed Brown   PetscFunctionBegin;
927*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
928*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(cda,&coords));
929*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
930*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
931*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,coords,&x));
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         PetscReal cx = ((PetscReal)i) / (((PetscReal)(mx-1)));
936c4762a1bSJed Brown         PetscReal cy = ((PetscReal)j) / (((PetscReal)(my-1)));
937c4762a1bSJed Brown         PetscReal cz = ((PetscReal)k) / (((PetscReal)(mz-1)));
938c4762a1bSJed Brown         PetscReal rad = user->rad + cy*user->height;
939c4762a1bSJed Brown         PetscReal ang = (cx - 0.5)*user->arc;
940c4762a1bSJed Brown         x[k][j][i][0] = rad*PetscSinReal(ang);
941c4762a1bSJed Brown         x[k][j][i][1] = rad*PetscCosReal(ang) - (user->rad + 0.5*user->height)*PetscCosReal(-0.5*user->arc);
942c4762a1bSJed Brown         x[k][j][i][2] = user->width*(cz - 0.5);
943c4762a1bSJed Brown       }
944c4762a1bSJed Brown     }
945c4762a1bSJed Brown   }
946*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,coords,&x));
947*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetCoordinates(da,coords));
948*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&coords));
949c4762a1bSJed Brown   PetscFunctionReturn(0);
950c4762a1bSJed Brown }
951c4762a1bSJed Brown 
952c4762a1bSJed Brown PetscErrorCode InitialGuess(DM da,AppCtx *user,Vec X)
953c4762a1bSJed Brown {
954c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,zs,xm,ym,zm;
955c4762a1bSJed Brown   PetscInt       mx,my,mz;
956c4762a1bSJed Brown   Field          ***x;
957c4762a1bSJed Brown 
958c4762a1bSJed Brown   PetscFunctionBegin;
959*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
960*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
961*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
962c4762a1bSJed Brown 
963c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
964c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
965c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
966c4762a1bSJed Brown         /* reference coordinates */
967c4762a1bSJed Brown         PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx-1)));
968c4762a1bSJed Brown         PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my-1)));
969c4762a1bSJed Brown         PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz-1)));
970c4762a1bSJed Brown         PetscReal o_x = p_x;
971c4762a1bSJed Brown         PetscReal o_y = p_y;
972c4762a1bSJed Brown         PetscReal o_z = p_z;
973c4762a1bSJed Brown         x[k][j][i][0] = o_x - p_x;
974c4762a1bSJed Brown         x[k][j][i][1] = o_y - p_y;
975c4762a1bSJed Brown         x[k][j][i][2] = o_z - p_z;
976c4762a1bSJed Brown       }
977c4762a1bSJed Brown     }
978c4762a1bSJed Brown   }
979*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
980c4762a1bSJed Brown   PetscFunctionReturn(0);
981c4762a1bSJed Brown }
982c4762a1bSJed Brown 
983c4762a1bSJed Brown PetscErrorCode FormRHS(DM da,AppCtx *user,Vec X)
984c4762a1bSJed Brown {
985c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,zs,xm,ym,zm;
986c4762a1bSJed Brown   PetscInt       mx,my,mz;
987c4762a1bSJed Brown   Field          ***x;
988c4762a1bSJed Brown 
989c4762a1bSJed Brown   PetscFunctionBegin;
990*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
991*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
992*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
993c4762a1bSJed Brown 
994c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
995c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
996c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
997c4762a1bSJed Brown         x[k][j][i][0] = 0.;
998c4762a1bSJed Brown         x[k][j][i][1] = 0.;
999c4762a1bSJed Brown         x[k][j][i][2] = 0.;
1000c4762a1bSJed Brown         if (i == (mx-1)/2 && j == (my-1)) x[k][j][i][1] = user->ploading/(mz-1);
1001c4762a1bSJed Brown       }
1002c4762a1bSJed Brown     }
1003c4762a1bSJed Brown   }
1004*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
1005c4762a1bSJed Brown   PetscFunctionReturn(0);
1006c4762a1bSJed Brown }
1007c4762a1bSJed Brown 
1008c4762a1bSJed Brown PetscErrorCode DisplayLine(SNES snes,Vec X)
1009c4762a1bSJed Brown {
1010c4762a1bSJed Brown   PetscInt       r,i,j=0,k=0,xs,xm,ys,ym,zs,zm,mx,my,mz;
1011c4762a1bSJed Brown   Field          ***x;
1012c4762a1bSJed Brown   CoordField     ***c;
1013c4762a1bSJed Brown   DM             da,cda;
1014c4762a1bSJed Brown   Vec            C;
1015c4762a1bSJed Brown   PetscMPIInt    size,rank;
1016c4762a1bSJed Brown 
1017c4762a1bSJed Brown   PetscFunctionBegin;
1018*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1019*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1020*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
1021*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
1022*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
1023*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(da,&C));
1024c4762a1bSJed Brown   j = my / 2;
1025c4762a1bSJed Brown   k = mz / 2;
1026*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
1027*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
1028*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
1029c4762a1bSJed Brown   for (r=0;r<size;r++) {
1030c4762a1bSJed Brown     if (rank == r) {
1031c4762a1bSJed Brown       if (j >= ys && j < ys+ym && k >= zs && k < zs+zm) {
1032c4762a1bSJed Brown         for (i=xs; i<xs+xm; i++) {
1033*5f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%D %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])));
1034c4762a1bSJed Brown         }
1035c4762a1bSJed Brown       }
1036c4762a1bSJed Brown     }
1037*5f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
1038c4762a1bSJed Brown   }
1039*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
1040*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
1041c4762a1bSJed Brown   PetscFunctionReturn(0);
1042c4762a1bSJed Brown }
1043c4762a1bSJed Brown 
1044c4762a1bSJed Brown /*TEST
1045c4762a1bSJed Brown 
1046c4762a1bSJed Brown    test:
1047c4762a1bSJed Brown       nsize: 2
1048c4762a1bSJed 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
1049c4762a1bSJed Brown       requires: !single
1050c4762a1bSJed Brown       timeoutfactor: 3
1051c4762a1bSJed Brown 
1052c4762a1bSJed Brown    test:
1053c4762a1bSJed Brown       suffix: 2
1054c4762a1bSJed 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
1055c4762a1bSJed Brown       requires: !single
1056c4762a1bSJed Brown 
1057c4762a1bSJed Brown    test:
1058c4762a1bSJed Brown       suffix: 3
1059c4762a1bSJed 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
1060c4762a1bSJed Brown       requires: !single
1061c4762a1bSJed Brown 
1062c4762a1bSJed Brown TEST*/
1063