xref: /petsc/src/snes/tutorials/ex16.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   MPI_Comm       comm;
99c4762a1bSJed Brown   SNES           snes;
100c4762a1bSJed Brown   DM             da;
101c4762a1bSJed Brown   Vec            x,X,b;
102c4762a1bSJed Brown   PetscBool      youngflg,poissonflg,muflg,lambdaflg,view=PETSC_FALSE,viewline=PETSC_FALSE;
103c4762a1bSJed Brown   PetscReal      poisson=0.2,young=4e4;
104c4762a1bSJed Brown   char           filename[PETSC_MAX_PATH_LEN] = "ex16.vts";
105c4762a1bSJed Brown   char           filename_def[PETSC_MAX_PATH_LEN] = "ex16_def.vts";
106c4762a1bSJed Brown 
107*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(FormElements());
109c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(comm,&snes));
1115f80ce2aSJacob 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));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes,(DM)da));
115c4762a1bSJed Brown 
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetNGS(snes,NonlinearGS,&user));
117c4762a1bSJed Brown 
1185f80ce2aSJacob 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));
119c4762a1bSJed Brown   user.loading     = 0.0;
120c4762a1bSJed Brown   user.arc         = PETSC_PI/3.;
121c4762a1bSJed Brown   user.mu          = 4.0;
122c4762a1bSJed Brown   user.lambda      = 1.0;
123c4762a1bSJed Brown   user.rad         = 100.0;
124c4762a1bSJed Brown   user.height      = 3.;
125c4762a1bSJed Brown   user.width       = 1.;
126c4762a1bSJed Brown   user.ploading    = -5e3;
127c4762a1bSJed Brown 
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-arc",&user.arc,NULL));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,&muflg));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-lambda",&user.lambda,&lambdaflg));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-rad",&user.rad,NULL));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-height",&user.height,NULL));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-width",&user.width,NULL));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-loading",&user.loading,NULL));
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ploading",&user.ploading,NULL));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-poisson",&poisson,&poissonflg));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-young",&young,&youngflg));
138c4762a1bSJed Brown   if ((youngflg || poissonflg) || !(muflg || lambdaflg)) {
139c4762a1bSJed Brown     /* set the lame' parameters based upon the poisson ratio and young's modulus */
140c4762a1bSJed Brown     user.lambda = poisson*young / ((1. + poisson)*(1. - 2.*poisson));
141c4762a1bSJed Brown     user.mu     = young/(2.*(1. + poisson));
142c4762a1bSJed Brown   }
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_line",&viewline,NULL));
145c4762a1bSJed Brown 
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,0,"x_disp"));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,1,"y_disp"));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,2,"z_disp"));
149c4762a1bSJed Brown 
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(da,&user));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(FormCoordinates(da,&user));
155c4762a1bSJed Brown 
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&b));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(InitialGuess(da,&user,x));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(FormRHS(da,&user,b));
160c4762a1bSJed Brown 
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"lambda: %f mu: %f\n",(double)user.lambda,(double)user.mu));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* show a cross-section of the initial state */
164c4762a1bSJed Brown   if (viewline) {
1655f80ce2aSJacob Faibussowitsch     CHKERRQ(DisplayLine(snes,x));
166c4762a1bSJed Brown   }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* get the loaded configuration */
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,b,x));
170c4762a1bSJed Brown 
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"Number of SNES iterations = %D\n", its));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetSolution(snes,&X));
174c4762a1bSJed Brown   /* show a cross-section of the final state */
175c4762a1bSJed Brown   if (viewline) {
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(DisplayLine(snes,X));
177c4762a1bSJed Brown   }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   if (view) {
180c4762a1bSJed Brown     PetscViewer viewer;
181c4762a1bSJed Brown     Vec         coords;
1825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer));
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(x,viewer));
1845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
1855f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinates(da,&coords));
1865f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(coords,1.0,x));
1875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename_def,FILE_MODE_WRITE,&viewer));
1885f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(x,viewer));
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerDestroy(&viewer));
190c4762a1bSJed Brown   }
191c4762a1bSJed Brown 
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
196*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
197*b122ec5aSJacob Faibussowitsch   return 0;
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown PetscInt OnBoundary(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   if ((i == 0 || i == mx-1) && j == my-1) return 1;
203c4762a1bSJed Brown   return 0;
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown void BoundaryValue(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz,PetscScalar *val,AppCtx *user)
207c4762a1bSJed Brown {
208c4762a1bSJed Brown   /* reference coordinates */
209c4762a1bSJed Brown   PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx-1)));
210c4762a1bSJed Brown   PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my-1)));
211c4762a1bSJed Brown   PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz-1)));
212c4762a1bSJed Brown   PetscReal o_x = p_x;
213c4762a1bSJed Brown   PetscReal o_y = p_y;
214c4762a1bSJed Brown   PetscReal o_z = p_z;
215c4762a1bSJed Brown   val[0] = o_x - p_x;
216c4762a1bSJed Brown   val[1] = o_y - p_y;
217c4762a1bSJed Brown   val[2] = o_z - p_z;
218c4762a1bSJed Brown }
219c4762a1bSJed Brown 
220c4762a1bSJed Brown void InvertTensor(PetscScalar *t, PetscScalar *ti,PetscReal *dett)
221c4762a1bSJed Brown {
222c4762a1bSJed Brown   const PetscScalar a = t[0];
223c4762a1bSJed Brown   const PetscScalar b = t[1];
224c4762a1bSJed Brown   const PetscScalar c = t[2];
225c4762a1bSJed Brown   const PetscScalar d = t[3];
226c4762a1bSJed Brown   const PetscScalar e = t[4];
227c4762a1bSJed Brown   const PetscScalar f = t[5];
228c4762a1bSJed Brown   const PetscScalar g = t[6];
229c4762a1bSJed Brown   const PetscScalar h = t[7];
230c4762a1bSJed Brown   const PetscScalar i = t[8];
231c4762a1bSJed Brown   const PetscReal   det = PetscRealPart(a*(e*i - f*h) - b*(i*d - f*g) + c*(d*h - e*g));
232c4762a1bSJed Brown   const PetscReal   di = 1. / det;
233c4762a1bSJed Brown   if (ti) {
234c4762a1bSJed Brown     const PetscScalar A = (e*i - f*h);
235c4762a1bSJed Brown     const PetscScalar B = -(d*i - f*g);
236c4762a1bSJed Brown     const PetscScalar C = (d*h - e*g);
237c4762a1bSJed Brown     const PetscScalar D = -(b*i - c*h);
238c4762a1bSJed Brown     const PetscScalar E = (a*i - c*g);
239c4762a1bSJed Brown     const PetscScalar F = -(a*h - b*g);
240c4762a1bSJed Brown     const PetscScalar G = (b*f - c*e);
241c4762a1bSJed Brown     const PetscScalar H = -(a*f - c*d);
242c4762a1bSJed Brown     const PetscScalar II = (a*e - b*d);
243c4762a1bSJed Brown     ti[0] = di*A;
244c4762a1bSJed Brown     ti[1] = di*D;
245c4762a1bSJed Brown     ti[2] = di*G;
246c4762a1bSJed Brown     ti[3] = di*B;
247c4762a1bSJed Brown     ti[4] = di*E;
248c4762a1bSJed Brown     ti[5] = di*H;
249c4762a1bSJed Brown     ti[6] = di*C;
250c4762a1bSJed Brown     ti[7] = di*F;
251c4762a1bSJed Brown     ti[8] = di*II;
252c4762a1bSJed Brown   }
253c4762a1bSJed Brown   if (dett) *dett = det;
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown void TensorTensor(PetscScalar *a,PetscScalar *b,PetscScalar *c)
257c4762a1bSJed Brown {
258c4762a1bSJed Brown   PetscInt i,j,m;
259c4762a1bSJed Brown   for (i=0;i<3;i++) {
260c4762a1bSJed Brown     for (j=0;j<3;j++) {
261c4762a1bSJed Brown       c[i+3*j] = 0;
262c4762a1bSJed Brown       for (m=0;m<3;m++)
263c4762a1bSJed Brown         c[i+3*j] += a[m+3*j]*b[i+3*m];
264c4762a1bSJed Brown     }
265c4762a1bSJed Brown   }
266c4762a1bSJed Brown }
267c4762a1bSJed Brown 
268c4762a1bSJed Brown void TensorTransposeTensor(PetscScalar *a,PetscScalar *b,PetscScalar *c)
269c4762a1bSJed Brown {
270c4762a1bSJed Brown   PetscInt i,j,m;
271c4762a1bSJed Brown   for (i=0;i<3;i++) {
272c4762a1bSJed Brown     for (j=0;j<3;j++) {
273c4762a1bSJed Brown       c[i+3*j] = 0;
274c4762a1bSJed Brown       for (m=0;m<3;m++)
275c4762a1bSJed Brown         c[i+3*j] += a[3*m+j]*b[i+3*m];
276c4762a1bSJed Brown     }
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280c4762a1bSJed Brown void TensorVector(PetscScalar *rot, PetscScalar *vec, PetscScalar *tvec)
281c4762a1bSJed Brown {
282c4762a1bSJed Brown   tvec[0] = rot[0]*vec[0] + rot[1]*vec[1] + rot[2]*vec[2];
283c4762a1bSJed Brown   tvec[1] = rot[3]*vec[0] + rot[4]*vec[1] + rot[5]*vec[2];
284c4762a1bSJed Brown   tvec[2] = rot[6]*vec[0] + rot[7]*vec[1] + rot[8]*vec[2];
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown void DeformationGradient(Field *ex,PetscInt qi,PetscInt qj,PetscInt qk,PetscScalar *invJ,PetscScalar *F)
288c4762a1bSJed Brown {
289c4762a1bSJed Brown   PetscInt ii,jj,kk,l;
290c4762a1bSJed Brown   for (l = 0; l < 9; l++) {
291c4762a1bSJed Brown     F[l] = 0.;
292c4762a1bSJed Brown   }
293c4762a1bSJed Brown   F[0] = 1.;
294c4762a1bSJed Brown   F[4] = 1.;
295c4762a1bSJed Brown   F[8] = 1.;
296c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
297c4762a1bSJed Brown   for (kk=0;kk<NB;kk++) {
298c4762a1bSJed Brown     for (jj=0;jj<NB;jj++) {
299c4762a1bSJed Brown       for (ii=0;ii<NB;ii++) {
300c4762a1bSJed Brown         PetscInt    idx = ii + jj*NB + kk*NB*NB;
301c4762a1bSJed Brown         PetscInt    bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
302c4762a1bSJed Brown         PetscScalar lgrad[3];
303c4762a1bSJed Brown         TensorVector(invJ,&grad[3*bidx],lgrad);
304c4762a1bSJed Brown         F[0] += lgrad[0]*ex[idx][0]; F[1] += lgrad[1]*ex[idx][0]; F[2] += lgrad[2]*ex[idx][0];
305c4762a1bSJed Brown         F[3] += lgrad[0]*ex[idx][1]; F[4] += lgrad[1]*ex[idx][1]; F[5] += lgrad[2]*ex[idx][1];
306c4762a1bSJed Brown         F[6] += lgrad[0]*ex[idx][2]; F[7] += lgrad[1]*ex[idx][2]; F[8] += lgrad[2]*ex[idx][2];
307c4762a1bSJed Brown       }
308c4762a1bSJed Brown     }
309c4762a1bSJed Brown   }
310c4762a1bSJed Brown }
311c4762a1bSJed Brown 
312c4762a1bSJed Brown void DeformationGradientJacobian(PetscInt qi,PetscInt qj,PetscInt qk,PetscInt ii,PetscInt jj,PetscInt kk,PetscInt fld,PetscScalar *invJ,PetscScalar *dF)
313c4762a1bSJed Brown {
314c4762a1bSJed Brown   PetscInt         l;
315c4762a1bSJed Brown   PetscScalar lgrad[3];
316c4762a1bSJed Brown   PetscInt    idx = ii + jj*NB + kk*NB*NB;
317c4762a1bSJed Brown   PetscInt    bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
318c4762a1bSJed Brown   for (l = 0; l < 9; l++) {
319c4762a1bSJed Brown     dF[l] = 0.;
320c4762a1bSJed Brown   }
321c4762a1bSJed Brown   /* form the deformation gradient at this basis function -- loop over element unknowns */
322c4762a1bSJed Brown   TensorVector(invJ,&grad[3*bidx],lgrad);
323c4762a1bSJed Brown   dF[3*fld] = lgrad[0]; dF[3*fld + 1] = lgrad[1]; dF[3*fld + 2] = lgrad[2];
324c4762a1bSJed Brown }
325c4762a1bSJed Brown 
326c4762a1bSJed Brown void LagrangeGreenStrain(PetscScalar *F,PetscScalar *E)
327c4762a1bSJed Brown {
328c4762a1bSJed Brown   PetscInt i,j,m;
329c4762a1bSJed Brown   for (i=0;i<3;i++) {
330c4762a1bSJed Brown     for (j=0;j<3;j++) {
331c4762a1bSJed Brown       E[i+3*j] = 0;
332c4762a1bSJed Brown       for (m=0;m<3;m++)
333c4762a1bSJed Brown         E[i+3*j] += 0.5*F[3*m+j]*F[i+3*m];
334c4762a1bSJed Brown     }
335c4762a1bSJed Brown   }
336c4762a1bSJed Brown   for (i=0;i<3;i++) {
337c4762a1bSJed Brown     E[i+3*i] -= 0.5;
338c4762a1bSJed Brown   }
339c4762a1bSJed Brown }
340c4762a1bSJed Brown 
341c4762a1bSJed Brown void SaintVenantKirchoff(PetscReal lambda,PetscReal mu,PetscScalar *F,PetscScalar *S)
342c4762a1bSJed Brown {
343c4762a1bSJed Brown   PetscInt    i,j;
344c4762a1bSJed Brown   PetscScalar E[9];
345c4762a1bSJed Brown   PetscScalar trE=0;
346c4762a1bSJed Brown   LagrangeGreenStrain(F,E);
347c4762a1bSJed Brown   for (i=0;i<3;i++) {
348c4762a1bSJed Brown     trE += E[i+3*i];
349c4762a1bSJed Brown   }
350c4762a1bSJed Brown   for (i=0;i<3;i++) {
351c4762a1bSJed Brown     for (j=0;j<3;j++) {
352c4762a1bSJed Brown       S[i+3*j] = 2.*mu*E[i+3*j];
353c4762a1bSJed Brown       if (i == j) {
354c4762a1bSJed Brown         S[i+3*j] += trE*lambda;
355c4762a1bSJed Brown       }
356c4762a1bSJed Brown     }
357c4762a1bSJed Brown   }
358c4762a1bSJed Brown }
359c4762a1bSJed Brown 
360c4762a1bSJed Brown void SaintVenantKirchoffJacobian(PetscReal lambda,PetscReal mu,PetscScalar *F,PetscScalar *dF,PetscScalar *dS)
361c4762a1bSJed Brown {
362c4762a1bSJed Brown   PetscScalar FtdF[9],dE[9];
363c4762a1bSJed Brown   PetscInt    i,j;
364c4762a1bSJed Brown   PetscScalar dtrE=0.;
365c4762a1bSJed Brown   TensorTransposeTensor(dF,F,dE);
366c4762a1bSJed Brown   TensorTransposeTensor(F,dF,FtdF);
367c4762a1bSJed Brown   for (i=0;i<9;i++) dE[i] += FtdF[i];
368c4762a1bSJed Brown   for (i=0;i<9;i++) dE[i] *= 0.5;
369c4762a1bSJed Brown   for (i=0;i<3;i++) {
370c4762a1bSJed Brown     dtrE += dE[i+3*i];
371c4762a1bSJed Brown   }
372c4762a1bSJed Brown   for (i=0;i<3;i++) {
373c4762a1bSJed Brown     for (j=0;j<3;j++) {
374c4762a1bSJed Brown       dS[i+3*j] = 2.*mu*dE[i+3*j];
375c4762a1bSJed Brown       if (i == j) {
376c4762a1bSJed Brown         dS[i+3*j] += dtrE*lambda;
377c4762a1bSJed Brown       }
378c4762a1bSJed Brown     }
379c4762a1bSJed Brown   }
380c4762a1bSJed Brown }
381c4762a1bSJed Brown 
382c4762a1bSJed Brown PetscErrorCode FormElements()
383c4762a1bSJed Brown {
384c4762a1bSJed Brown   PetscInt i,j,k,ii,jj,kk;
385c4762a1bSJed Brown   PetscReal bx,by,bz,dbx,dby,dbz;
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   PetscFunctionBegin;
388c4762a1bSJed Brown   /* construct the basis function values and derivatives */
389c4762a1bSJed Brown   for (k = 0; k < NB; k++) {
390c4762a1bSJed Brown     for (j = 0; j < NB; j++) {
391c4762a1bSJed Brown       for (i = 0; i < NB; i++) {
392c4762a1bSJed Brown         /* loop over the quadrature points */
393c4762a1bSJed Brown         for (kk = 0; kk < NQ; kk++) {
394c4762a1bSJed Brown           for (jj = 0; jj < NQ; jj++) {
395c4762a1bSJed Brown             for (ii = 0; ii < NQ; ii++) {
396c4762a1bSJed Brown               PetscInt idx = ii + NQ*jj + NQ*NQ*kk + NEQ*i + NEQ*NB*j + NEQ*NB*NB*k;
397c4762a1bSJed Brown               bx = pts[ii];
398c4762a1bSJed Brown               by = pts[jj];
399c4762a1bSJed Brown               bz = pts[kk];
400c4762a1bSJed Brown               dbx = 1.;
401c4762a1bSJed Brown               dby = 1.;
402c4762a1bSJed Brown               dbz = 1.;
403c4762a1bSJed Brown               if (i == 0) {bx = 1. - bx; dbx = -1;}
404c4762a1bSJed Brown               if (j == 0) {by = 1. - by; dby = -1;}
405c4762a1bSJed Brown               if (k == 0) {bz = 1. - bz; dbz = -1;}
406c4762a1bSJed Brown               vals[idx] = bx*by*bz;
407c4762a1bSJed Brown               grad[3*idx + 0] = dbx*by*bz;
408c4762a1bSJed Brown               grad[3*idx + 1] = dby*bx*bz;
409c4762a1bSJed Brown               grad[3*idx + 2] = dbz*bx*by;
410c4762a1bSJed Brown             }
411c4762a1bSJed Brown           }
412c4762a1bSJed Brown         }
413c4762a1bSJed Brown       }
414c4762a1bSJed Brown     }
415c4762a1bSJed Brown   }
416c4762a1bSJed Brown   PetscFunctionReturn(0);
417c4762a1bSJed Brown }
418c4762a1bSJed Brown 
419c4762a1bSJed 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)
420c4762a1bSJed Brown {
421c4762a1bSJed Brown   PetscInt m;
422c4762a1bSJed Brown   PetscInt ii,jj,kk;
423c4762a1bSJed Brown   /* gather the data -- loop over element unknowns */
424c4762a1bSJed Brown   for (kk=0;kk<NB;kk++) {
425c4762a1bSJed Brown     for (jj=0;jj<NB;jj++) {
426c4762a1bSJed Brown       for (ii=0;ii<NB;ii++) {
427c4762a1bSJed Brown         PetscInt idx = ii + jj*NB + kk*NB*NB;
428c4762a1bSJed Brown         /* decouple the boundary nodes for the displacement variables */
429c4762a1bSJed Brown         if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz)) {
430c4762a1bSJed Brown           BoundaryValue(i+ii,j+jj,k+kk,mx,my,mz,ex[idx],user);
431c4762a1bSJed Brown         } else {
432c4762a1bSJed Brown           for (m=0;m<3;m++) {
433c4762a1bSJed Brown             ex[idx][m] = x[k+kk][j+jj][i+ii][m];
434c4762a1bSJed Brown           }
435c4762a1bSJed Brown         }
436c4762a1bSJed Brown         for (m=0;m<3;m++) {
437c4762a1bSJed Brown           ec[idx][m] = c[k+kk][j+jj][i+ii][m];
438c4762a1bSJed Brown         }
439c4762a1bSJed Brown       }
440c4762a1bSJed Brown     }
441c4762a1bSJed Brown   }
442c4762a1bSJed Brown }
443c4762a1bSJed Brown 
444c4762a1bSJed Brown void QuadraturePointGeometricJacobian(CoordField *ec,PetscInt qi,PetscInt qj,PetscInt qk, PetscScalar *J)
445c4762a1bSJed Brown {
446c4762a1bSJed Brown   PetscInt ii,jj,kk;
447c4762a1bSJed Brown   /* construct the gradient at the given quadrature point named by i,j,k */
448c4762a1bSJed Brown   for (ii = 0; ii < 9; ii++) {
449c4762a1bSJed Brown     J[ii] = 0;
450c4762a1bSJed Brown   }
451c4762a1bSJed Brown   for (kk = 0; kk < NB; kk++) {
452c4762a1bSJed Brown     for (jj = 0; jj < NB; jj++) {
453c4762a1bSJed Brown       for (ii = 0; ii < NB; ii++) {
454c4762a1bSJed Brown         PetscInt idx = ii + jj*NB + kk*NB*NB;
455c4762a1bSJed Brown         PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
456c4762a1bSJed 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];
457c4762a1bSJed 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];
458c4762a1bSJed 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];
459c4762a1bSJed Brown       }
460c4762a1bSJed Brown     }
461c4762a1bSJed Brown   }
462c4762a1bSJed Brown }
463c4762a1bSJed Brown 
464c4762a1bSJed Brown void FormElementJacobian(Field *ex,CoordField *ec,Field *ef,PetscScalar *ej,AppCtx *user)
465c4762a1bSJed Brown {
466c4762a1bSJed Brown   PetscReal   vol;
467c4762a1bSJed Brown   PetscScalar J[9];
468c4762a1bSJed Brown   PetscScalar invJ[9];
469c4762a1bSJed Brown   PetscScalar F[9],S[9],dF[9],dS[9],dFS[9],FdS[9],FS[9];
470c4762a1bSJed Brown   PetscReal   scl;
471c4762a1bSJed Brown   PetscInt    i,j,k,l,ii,jj,kk,ll,qi,qj,qk,m;
472c4762a1bSJed Brown 
473c4762a1bSJed Brown   if (ej) for (i = 0; i < NPB*NPB; i++) ej[i] = 0.;
474c4762a1bSJed Brown   if (ef) for (i = 0; i < NEB; i++) {ef[i][0] = 0.;ef[i][1] = 0.;ef[i][2] = 0.;}
475c4762a1bSJed Brown   /* loop over quadrature */
476c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
477c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
478c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
479c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec,qi,qj,qk,J);
480c4762a1bSJed Brown         InvertTensor(J,invJ,&vol);
481c4762a1bSJed Brown         scl = vol*wts[qi]*wts[qj]*wts[qk];
482c4762a1bSJed Brown         DeformationGradient(ex,qi,qj,qk,invJ,F);
483c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda,user->mu,F,S);
484c4762a1bSJed Brown         /* form the function */
485c4762a1bSJed Brown         if (ef) {
486c4762a1bSJed Brown           TensorTensor(F,S,FS);
487c4762a1bSJed Brown           for (kk=0;kk<NB;kk++) {
488c4762a1bSJed Brown             for (jj=0;jj<NB;jj++) {
489c4762a1bSJed Brown               for (ii=0;ii<NB;ii++) {
490c4762a1bSJed Brown                 PetscInt idx = ii + jj*NB + kk*NB*NB;
491c4762a1bSJed Brown                 PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
492c4762a1bSJed Brown                 PetscScalar lgrad[3];
493c4762a1bSJed Brown                 TensorVector(invJ,&grad[3*bidx],lgrad);
494c4762a1bSJed Brown                 /* mu*F : grad phi_{u,v,w} */
495c4762a1bSJed Brown                 for (m=0;m<3;m++) {
496c4762a1bSJed Brown                   ef[idx][m] += scl*
497c4762a1bSJed Brown                     (lgrad[0]*FS[3*m + 0] + lgrad[1]*FS[3*m + 1] + lgrad[2]*FS[3*m + 2]);
498c4762a1bSJed Brown                 }
499c4762a1bSJed Brown                 ef[idx][1] -= scl*user->loading*vals[bidx];
500c4762a1bSJed Brown               }
501c4762a1bSJed Brown             }
502c4762a1bSJed Brown           }
503c4762a1bSJed Brown         }
504c4762a1bSJed Brown         /* form the jacobian */
505c4762a1bSJed Brown         if (ej) {
506c4762a1bSJed Brown           /* loop over trialfunctions */
507c4762a1bSJed Brown           for (k=0;k<NB;k++) {
508c4762a1bSJed Brown             for (j=0;j<NB;j++) {
509c4762a1bSJed Brown               for (i=0;i<NB;i++) {
510c4762a1bSJed Brown                 for (l=0;l<3;l++) {
511c4762a1bSJed Brown                   PetscInt tridx = l + 3*(i + j*NB + k*NB*NB);
512c4762a1bSJed Brown                   DeformationGradientJacobian(qi,qj,qk,i,j,k,l,invJ,dF);
513c4762a1bSJed Brown                   SaintVenantKirchoffJacobian(user->lambda,user->mu,F,dF,dS);
514c4762a1bSJed Brown                   TensorTensor(dF,S,dFS);
515c4762a1bSJed Brown                   TensorTensor(F,dS,FdS);
516c4762a1bSJed Brown                   for (m=0;m<9;m++) dFS[m] += FdS[m];
517c4762a1bSJed Brown                   /* loop over testfunctions */
518c4762a1bSJed Brown                   for (kk=0;kk<NB;kk++) {
519c4762a1bSJed Brown                     for (jj=0;jj<NB;jj++) {
520c4762a1bSJed Brown                       for (ii=0;ii<NB;ii++) {
521c4762a1bSJed Brown                         PetscInt idx = ii + jj*NB + kk*NB*NB;
522c4762a1bSJed Brown                         PetscInt bidx = 8*idx + qi + NQ*qj + NQ*NQ*qk;
523c4762a1bSJed Brown                         PetscScalar lgrad[3];
524c4762a1bSJed Brown                         TensorVector(invJ,&grad[3*bidx],lgrad);
525c4762a1bSJed Brown                         for (ll=0; ll<3;ll++) {
526c4762a1bSJed Brown                           PetscInt teidx = ll + 3*(ii + jj*NB + kk*NB*NB);
527c4762a1bSJed Brown                           ej[teidx + NPB*tridx] += scl*
528c4762a1bSJed Brown                             (lgrad[0]*dFS[3*ll + 0] + lgrad[1]*dFS[3*ll + 1] + lgrad[2]*dFS[3*ll+2]);
529c4762a1bSJed Brown                         }
530c4762a1bSJed Brown                       }
531c4762a1bSJed Brown                     }
532c4762a1bSJed Brown                   } /* end of testfunctions */
533c4762a1bSJed Brown                 }
534c4762a1bSJed Brown               }
535c4762a1bSJed Brown             }
536c4762a1bSJed Brown           } /* end of trialfunctions */
537c4762a1bSJed Brown         }
538c4762a1bSJed Brown       }
539c4762a1bSJed Brown     }
540c4762a1bSJed Brown   } /* end of quadrature points */
541c4762a1bSJed Brown }
542c4762a1bSJed Brown 
543c4762a1bSJed Brown void FormPBJacobian(PetscInt i,PetscInt j,PetscInt k,Field *ex,CoordField *ec,Field *ef,PetscScalar *ej,AppCtx *user)
544c4762a1bSJed Brown {
545c4762a1bSJed Brown   PetscReal   vol;
546c4762a1bSJed Brown   PetscScalar J[9];
547c4762a1bSJed Brown   PetscScalar invJ[9];
548c4762a1bSJed Brown   PetscScalar F[9],S[9],dF[9],dS[9],dFS[9],FdS[9],FS[9];
549c4762a1bSJed Brown   PetscReal   scl;
550c4762a1bSJed Brown   PetscInt    l,ll,qi,qj,qk,m;
551c4762a1bSJed Brown   PetscInt idx = i + j*NB + k*NB*NB;
552c4762a1bSJed Brown   PetscScalar lgrad[3];
553c4762a1bSJed Brown 
554c4762a1bSJed Brown   if (ej) for (l = 0; l < 9; l++) ej[l] = 0.;
555c4762a1bSJed Brown   if (ef) for (l = 0; l < 1; l++) {ef[l][0] = 0.;ef[l][1] = 0.;ef[l][2] = 0.;}
556c4762a1bSJed Brown   /* loop over quadrature */
557c4762a1bSJed Brown   for (qk = 0; qk < NQ; qk++) {
558c4762a1bSJed Brown     for (qj = 0; qj < NQ; qj++) {
559c4762a1bSJed Brown       for (qi = 0; qi < NQ; qi++) {
560c4762a1bSJed Brown         PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk;
561c4762a1bSJed Brown         QuadraturePointGeometricJacobian(ec,qi,qj,qk,J);
562c4762a1bSJed Brown         InvertTensor(J,invJ,&vol);
563c4762a1bSJed Brown         TensorVector(invJ,&grad[3*bidx],lgrad);
564c4762a1bSJed Brown         scl = vol*wts[qi]*wts[qj]*wts[qk];
565c4762a1bSJed Brown         DeformationGradient(ex,qi,qj,qk,invJ,F);
566c4762a1bSJed Brown         SaintVenantKirchoff(user->lambda,user->mu,F,S);
567c4762a1bSJed Brown         /* form the function */
568c4762a1bSJed Brown         if (ef) {
569c4762a1bSJed Brown           TensorTensor(F,S,FS);
570c4762a1bSJed Brown           for (m=0;m<3;m++) {
571c4762a1bSJed Brown             ef[0][m] += scl*
572c4762a1bSJed Brown               (lgrad[0]*FS[3*m + 0] + lgrad[1]*FS[3*m + 1] + lgrad[2]*FS[3*m + 2]);
573c4762a1bSJed Brown           }
574c4762a1bSJed Brown           ef[0][1] -= scl*user->loading*vals[bidx];
575c4762a1bSJed Brown         }
576c4762a1bSJed Brown         /* form the jacobian */
577c4762a1bSJed Brown         if (ej) {
578c4762a1bSJed Brown           for (l=0;l<3;l++) {
579c4762a1bSJed Brown             DeformationGradientJacobian(qi,qj,qk,i,j,k,l,invJ,dF);
580c4762a1bSJed Brown             SaintVenantKirchoffJacobian(user->lambda,user->mu,F,dF,dS);
581c4762a1bSJed Brown             TensorTensor(dF,S,dFS);
582c4762a1bSJed Brown             TensorTensor(F,dS,FdS);
583c4762a1bSJed Brown             for (m=0;m<9;m++) dFS[m] += FdS[m];
584c4762a1bSJed Brown             for (ll=0; ll<3;ll++) {
585c4762a1bSJed Brown               ej[ll + 3*l] += scl*
586c4762a1bSJed Brown                 (lgrad[0]*dFS[3*ll + 0] + lgrad[1]*dFS[3*ll + 1] + lgrad[2]*dFS[3*ll+2]);
587c4762a1bSJed Brown             }
588c4762a1bSJed Brown           }
589c4762a1bSJed Brown         }
590c4762a1bSJed Brown       }
591c4762a1bSJed Brown     } /* end of quadrature points */
592c4762a1bSJed Brown   }
593c4762a1bSJed Brown }
594c4762a1bSJed Brown 
595c4762a1bSJed Brown void ApplyBCsElement(PetscInt mx,PetscInt my, PetscInt mz, PetscInt i, PetscInt j, PetscInt k,PetscScalar *jacobian)
596c4762a1bSJed Brown {
597c4762a1bSJed Brown   PetscInt ii,jj,kk,ll,ei,ej,ek,el;
598c4762a1bSJed Brown   for (kk=0;kk<NB;kk++) {
599c4762a1bSJed Brown     for (jj=0;jj<NB;jj++) {
600c4762a1bSJed Brown       for (ii=0;ii<NB;ii++) {
601c4762a1bSJed Brown         for (ll = 0;ll<3;ll++) {
602c4762a1bSJed Brown           PetscInt tridx = ll + 3*(ii + jj*NB + kk*NB*NB);
603c4762a1bSJed Brown           for (ek=0;ek<NB;ek++) {
604c4762a1bSJed Brown             for (ej=0;ej<NB;ej++) {
605c4762a1bSJed Brown               for (ei=0;ei<NB;ei++) {
606c4762a1bSJed Brown                 for (el=0;el<3;el++) {
607c4762a1bSJed Brown                   if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz) || OnBoundary(i+ei,j+ej,k+ek,mx,my,mz)) {
608c4762a1bSJed Brown                     PetscInt teidx = el + 3*(ei + ej*NB + ek*NB*NB);
609c4762a1bSJed Brown                     if (teidx == tridx) {
610c4762a1bSJed Brown                       jacobian[tridx + NPB*teidx] = 1.;
611c4762a1bSJed Brown                     } else {
612c4762a1bSJed Brown                       jacobian[tridx + NPB*teidx] = 0.;
613c4762a1bSJed Brown                     }
614c4762a1bSJed Brown                   }
615c4762a1bSJed Brown                 }
616c4762a1bSJed Brown               }
617c4762a1bSJed Brown             }
618c4762a1bSJed Brown           }
619c4762a1bSJed Brown         }
620c4762a1bSJed Brown       }
621c4762a1bSJed Brown     }
622c4762a1bSJed Brown   }
623c4762a1bSJed Brown }
624c4762a1bSJed Brown 
625c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,Field ***x,Mat jacpre,Mat jac,void *ptr)
626c4762a1bSJed Brown {
627c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
628c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
629c4762a1bSJed Brown   PetscInt       i,j,k,m,l;
630c4762a1bSJed Brown   PetscInt       ii,jj,kk;
631c4762a1bSJed Brown   PetscScalar    ej[NPB*NPB];
632c4762a1bSJed Brown   PetscScalar    vals[NPB*NPB];
633c4762a1bSJed Brown   Field          ex[NEB];
634c4762a1bSJed Brown   CoordField     ec[NEB];
635c4762a1bSJed Brown 
636c4762a1bSJed Brown   PetscInt       xs=info->xs,ys=info->ys,zs=info->zs;
637c4762a1bSJed Brown   PetscInt       xm=info->xm,ym=info->ym,zm=info->zm;
638c4762a1bSJed Brown   PetscInt       xes,yes,zes,xee,yee,zee;
639c4762a1bSJed Brown   PetscInt       mx=info->mx,my=info->my,mz=info->mz;
640c4762a1bSJed Brown   DM             cda;
641c4762a1bSJed Brown   CoordField     ***c;
642c4762a1bSJed Brown   Vec            C;
643c4762a1bSJed Brown   PetscInt       nrows;
644c4762a1bSJed Brown   MatStencil     col[NPB],row[NPB];
645c4762a1bSJed Brown   PetscScalar    v[9];
646c4762a1bSJed Brown 
647c4762a1bSJed Brown   PetscFunctionBegin;
6485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(info->da,&cda));
6495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(info->da,&C));
6505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
6515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatScale(jac,0.0));
652c4762a1bSJed Brown 
653c4762a1bSJed Brown   xes = xs;
654c4762a1bSJed Brown   yes = ys;
655c4762a1bSJed Brown   zes = zs;
656c4762a1bSJed Brown   xee = xs+xm;
657c4762a1bSJed Brown   yee = ys+ym;
658c4762a1bSJed Brown   zee = zs+zm;
659c4762a1bSJed Brown   if (xs > 0) xes = xs-1;
660c4762a1bSJed Brown   if (ys > 0) yes = ys-1;
661c4762a1bSJed Brown   if (zs > 0) zes = zs-1;
662c4762a1bSJed Brown   if (xs+xm == mx) xee = xs+xm-1;
663c4762a1bSJed Brown   if (ys+ym == my) yee = ys+ym-1;
664c4762a1bSJed Brown   if (zs+zm == mz) zee = zs+zm-1;
665c4762a1bSJed Brown   for (k=zes; k<zee; k++) {
666c4762a1bSJed Brown     for (j=yes; j<yee; j++) {
667c4762a1bSJed Brown       for (i=xes; i<xee; i++) {
668c4762a1bSJed Brown         GatherElementData(mx,my,mz,x,c,i,j,k,ex,ec,user);
669c4762a1bSJed Brown         FormElementJacobian(ex,ec,NULL,ej,user);
670c4762a1bSJed Brown         ApplyBCsElement(mx,my,mz,i,j,k,ej);
671c4762a1bSJed Brown         nrows = 0.;
672c4762a1bSJed Brown         for (kk=0;kk<NB;kk++) {
673c4762a1bSJed Brown           for (jj=0;jj<NB;jj++) {
674c4762a1bSJed Brown             for (ii=0;ii<NB;ii++) {
675c4762a1bSJed Brown               PetscInt idx = ii + jj*2 + kk*4;
676c4762a1bSJed Brown               for (m=0;m<3;m++) {
677c4762a1bSJed Brown                 col[3*idx+m].i = i+ii;
678c4762a1bSJed Brown                 col[3*idx+m].j = j+jj;
679c4762a1bSJed Brown                 col[3*idx+m].k = k+kk;
680c4762a1bSJed Brown                 col[3*idx+m].c = m;
681c4762a1bSJed Brown                 if (i+ii >= xs && i+ii < xm+xs && j+jj >= ys && j+jj < ys+ym && k+kk >= zs && k+kk < zs+zm) {
682c4762a1bSJed Brown                   row[nrows].i = i+ii;
683c4762a1bSJed Brown                   row[nrows].j = j+jj;
684c4762a1bSJed Brown                   row[nrows].k = k+kk;
685c4762a1bSJed Brown                   row[nrows].c = m;
686c4762a1bSJed Brown                   for (l=0;l<NPB;l++) vals[NPB*nrows + l] = ej[NPB*(3*idx+m) + l];
687c4762a1bSJed Brown                   nrows++;
688c4762a1bSJed Brown                 }
689c4762a1bSJed Brown               }
690c4762a1bSJed Brown             }
691c4762a1bSJed Brown           }
692c4762a1bSJed Brown         }
6935f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(jac,nrows,row,NPB,col,vals,ADD_VALUES));
694c4762a1bSJed Brown       }
695c4762a1bSJed Brown     }
696c4762a1bSJed Brown   }
697c4762a1bSJed Brown 
6985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FLUSH_ASSEMBLY));
6995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FLUSH_ASSEMBLY));
700c4762a1bSJed Brown 
701c4762a1bSJed Brown   /* set the diagonal */
702c4762a1bSJed 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.;
703c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
704c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
705c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
706c4762a1bSJed Brown         if (OnBoundary(i,j,k,mx,my,mz)) {
707c4762a1bSJed Brown           for (m=0; m<3;m++) {
708c4762a1bSJed Brown             col[m].i = i;
709c4762a1bSJed Brown             col[m].j = j;
710c4762a1bSJed Brown             col[m].k = k;
711c4762a1bSJed Brown             col[m].c = m;
712c4762a1bSJed Brown           }
7135f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,3,col,3,col,v,INSERT_VALUES));
714c4762a1bSJed Brown         }
715c4762a1bSJed Brown       }
716c4762a1bSJed Brown     }
717c4762a1bSJed Brown   }
718c4762a1bSJed Brown 
7195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
7205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
721c4762a1bSJed Brown 
7225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
723c4762a1bSJed Brown   PetscFunctionReturn(0);
724c4762a1bSJed Brown }
725c4762a1bSJed Brown 
726c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field ***x,Field ***f,void *ptr)
727c4762a1bSJed Brown {
728c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
729c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
730c4762a1bSJed Brown   PetscInt       i,j,k,l;
731c4762a1bSJed Brown   PetscInt       ii,jj,kk;
732c4762a1bSJed Brown 
733c4762a1bSJed Brown   Field          ef[NEB];
734c4762a1bSJed Brown   Field          ex[NEB];
735c4762a1bSJed Brown   CoordField     ec[NEB];
736c4762a1bSJed Brown 
737c4762a1bSJed Brown   PetscInt       xs=info->xs,ys=info->ys,zs=info->zs;
738c4762a1bSJed Brown   PetscInt       xm=info->xm,ym=info->ym,zm=info->zm;
739c4762a1bSJed Brown   PetscInt       xes,yes,zes,xee,yee,zee;
740c4762a1bSJed Brown   PetscInt       mx=info->mx,my=info->my,mz=info->mz;
741c4762a1bSJed Brown   DM             cda;
742c4762a1bSJed Brown   CoordField     ***c;
743c4762a1bSJed Brown   Vec            C;
744c4762a1bSJed Brown 
745c4762a1bSJed Brown   PetscFunctionBegin;
7465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(info->da,&cda));
7475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(info->da,&C));
7485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
7495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(info->da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
7505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(info->da,&xs,&ys,&zs,&xm,&ym,&zm));
751c4762a1bSJed Brown 
752c4762a1bSJed Brown   /* loop over elements */
753c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
754c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
755c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
756c4762a1bSJed Brown         for (l=0;l<3;l++) {
757c4762a1bSJed Brown           f[k][j][i][l] = 0.;
758c4762a1bSJed Brown         }
759c4762a1bSJed Brown       }
760c4762a1bSJed Brown     }
761c4762a1bSJed Brown   }
762c4762a1bSJed Brown   /* element starts and ends */
763c4762a1bSJed Brown   xes = xs;
764c4762a1bSJed Brown   yes = ys;
765c4762a1bSJed Brown   zes = zs;
766c4762a1bSJed Brown   xee = xs+xm;
767c4762a1bSJed Brown   yee = ys+ym;
768c4762a1bSJed Brown   zee = zs+zm;
769c4762a1bSJed Brown   if (xs > 0) xes = xs - 1;
770c4762a1bSJed Brown   if (ys > 0) yes = ys - 1;
771c4762a1bSJed Brown   if (zs > 0) zes = zs - 1;
772c4762a1bSJed Brown   if (xs+xm == mx) xee = xs+xm-1;
773c4762a1bSJed Brown   if (ys+ym == my) yee = ys+ym-1;
774c4762a1bSJed Brown   if (zs+zm == mz) zee = zs+zm-1;
775c4762a1bSJed Brown   for (k=zes; k<zee; k++) {
776c4762a1bSJed Brown     for (j=yes; j<yee; j++) {
777c4762a1bSJed Brown       for (i=xes; i<xee; i++) {
778c4762a1bSJed Brown         GatherElementData(mx,my,mz,x,c,i,j,k,ex,ec,user);
779c4762a1bSJed Brown         FormElementJacobian(ex,ec,ef,NULL,user);
780c4762a1bSJed Brown         /* put this element's additions into the residuals */
781c4762a1bSJed Brown         for (kk=0;kk<NB;kk++) {
782c4762a1bSJed Brown           for (jj=0;jj<NB;jj++) {
783c4762a1bSJed Brown             for (ii=0;ii<NB;ii++) {
784c4762a1bSJed Brown               PetscInt idx = ii + jj*NB + kk*NB*NB;
785c4762a1bSJed Brown               if (k+kk >= zs && j+jj >= ys && i+ii >= xs && k+kk < zs+zm && j+jj < ys+ym && i+ii < xs+xm) {
786c4762a1bSJed Brown                 if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz)) {
787c4762a1bSJed Brown                   for (l=0;l<3;l++)
788c4762a1bSJed Brown                     f[k+kk][j+jj][i+ii][l] = x[k+kk][j+jj][i+ii][l] - ex[idx][l];
789c4762a1bSJed Brown                 } else {
790c4762a1bSJed Brown                   for (l=0;l<3;l++)
791c4762a1bSJed Brown                     f[k+kk][j+jj][i+ii][l] += ef[idx][l];
792c4762a1bSJed Brown                 }
793c4762a1bSJed Brown               }
794c4762a1bSJed Brown             }
795c4762a1bSJed Brown           }
796c4762a1bSJed Brown         }
797c4762a1bSJed Brown       }
798c4762a1bSJed Brown     }
799c4762a1bSJed Brown   }
8005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
801c4762a1bSJed Brown   PetscFunctionReturn(0);
802c4762a1bSJed Brown }
803c4762a1bSJed Brown 
804c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X,Vec B,void *ptr)
805c4762a1bSJed Brown {
806c4762a1bSJed Brown   /* values for each basis function at each quadrature point */
807c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
808c4762a1bSJed Brown   PetscInt       i,j,k,l,m,n,s;
809c4762a1bSJed Brown   PetscInt       pi,pj,pk;
810c4762a1bSJed Brown   Field          ef[1];
811c4762a1bSJed Brown   Field          ex[8];
812c4762a1bSJed Brown   PetscScalar    ej[9];
813c4762a1bSJed Brown   CoordField     ec[8];
814c4762a1bSJed Brown   PetscScalar    pjac[9],pjinv[9];
815c4762a1bSJed Brown   PetscScalar    pf[3],py[3];
816c4762a1bSJed Brown   PetscInt       xs,ys,zs;
817c4762a1bSJed Brown   PetscInt       xm,ym,zm;
818c4762a1bSJed Brown   PetscInt       mx,my,mz;
819c4762a1bSJed Brown   DM             cda;
820c4762a1bSJed Brown   CoordField     ***c;
821c4762a1bSJed Brown   Vec            C;
822c4762a1bSJed Brown   DM             da;
823c4762a1bSJed Brown   Vec            Xl,Bl;
824c4762a1bSJed Brown   Field          ***x,***b;
825c4762a1bSJed Brown   PetscInt       sweeps,its;
826c4762a1bSJed Brown   PetscReal      atol,rtol,stol;
827c4762a1bSJed Brown   PetscReal      fnorm0 = 0.0,fnorm,ynorm,xnorm = 0.0;
828c4762a1bSJed Brown 
829c4762a1bSJed Brown   PetscFunctionBegin;
8305f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESNGSGetSweeps(snes,&sweeps));
8315f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its));
832c4762a1bSJed Brown 
8335f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
8345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&Xl));
835c4762a1bSJed Brown   if (B) {
8365f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalVector(da,&Bl));
837c4762a1bSJed Brown   }
8385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xl));
8395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xl));
840c4762a1bSJed Brown   if (B) {
8415f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(da,B,INSERT_VALUES,Bl));
8425f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(da,B,INSERT_VALUES,Bl));
843c4762a1bSJed Brown   }
8445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,Xl,&x));
8455f80ce2aSJacob Faibussowitsch   if (B) CHKERRQ(DMDAVecGetArray(da,Bl,&b));
846c4762a1bSJed Brown 
8475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
8485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(da,&C));
8495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
8505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
8515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
852c4762a1bSJed Brown 
853c4762a1bSJed Brown   for (s=0;s<sweeps;s++) {
854c4762a1bSJed Brown     for (k=zs; k<zs+zm; k++) {
855c4762a1bSJed Brown       for (j=ys; j<ys+ym; j++) {
856c4762a1bSJed Brown         for (i=xs; i<xs+xm; i++) {
857c4762a1bSJed Brown           if (OnBoundary(i,j,k,mx,my,mz)) {
858c4762a1bSJed Brown             BoundaryValue(i,j,k,mx,my,mz,x[k][j][i],user);
859c4762a1bSJed Brown           } else {
860c4762a1bSJed Brown             for (n=0;n<its;n++) {
861c4762a1bSJed Brown               for (m=0;m<9;m++) pjac[m] = 0.;
862c4762a1bSJed Brown               for (m=0;m<3;m++) pf[m] = 0.;
863c4762a1bSJed Brown               /* gather the elements for this point */
864c4762a1bSJed Brown               for (pk=-1; pk<1; pk++) {
865c4762a1bSJed Brown                 for (pj=-1; pj<1; pj++) {
866c4762a1bSJed Brown                   for (pi=-1; pi<1; pi++) {
867c4762a1bSJed Brown                     /* check that this element exists */
868c4762a1bSJed Brown                     if (i+pi >= 0 && i+pi < mx-1 && j+pj >= 0 && j+pj < my-1 && k+pk >= 0 && k+pk < mz-1) {
869c4762a1bSJed Brown                       /* create the element function and jacobian */
870c4762a1bSJed Brown                       GatherElementData(mx,my,mz,x,c,i+pi,j+pj,k+pk,ex,ec,user);
871c4762a1bSJed Brown                       FormPBJacobian(-pi,-pj,-pk,ex,ec,ef,ej,user);
872c4762a1bSJed Brown                       /* extract the point named by i,j,k from the whole element jacobian and function */
873c4762a1bSJed Brown                       for (l=0;l<3;l++) {
874c4762a1bSJed Brown                         pf[l] += ef[0][l];
875c4762a1bSJed Brown                         for (m=0;m<3;m++) {
876c4762a1bSJed Brown                           pjac[3*m+l] += ej[3*m+l];
877c4762a1bSJed Brown                         }
878c4762a1bSJed Brown                       }
879c4762a1bSJed Brown                     }
880c4762a1bSJed Brown                   }
881c4762a1bSJed Brown                 }
882c4762a1bSJed Brown               }
883c4762a1bSJed Brown               /* invert */
884c4762a1bSJed Brown               InvertTensor(pjac,pjinv,NULL);
885c4762a1bSJed Brown               /* apply */
886c4762a1bSJed Brown               if (B) for (m=0;m<3;m++) {
887c4762a1bSJed Brown                   pf[m] -= b[k][j][i][m];
888c4762a1bSJed Brown                 }
889c4762a1bSJed Brown               TensorVector(pjinv,pf,py);
890c4762a1bSJed Brown               xnorm=0.;
891c4762a1bSJed Brown               for (m=0;m<3;m++) {
892c4762a1bSJed Brown                 x[k][j][i][m] -= py[m];
893c4762a1bSJed Brown                 xnorm += PetscRealPart(x[k][j][i][m]*x[k][j][i][m]);
894c4762a1bSJed Brown               }
895c4762a1bSJed Brown               fnorm = PetscRealPart(pf[0]*pf[0]+pf[1]*pf[1]+pf[2]*pf[2]);
896c4762a1bSJed Brown               if (n==0) fnorm0 = fnorm;
897c4762a1bSJed Brown               ynorm = PetscRealPart(py[0]*py[0]+py[1]*py[1]+py[2]*py[2]);
898c4762a1bSJed Brown               if (fnorm < atol*atol || fnorm < rtol*rtol*fnorm0 || ynorm < stol*stol*xnorm) break;
899c4762a1bSJed Brown             }
900c4762a1bSJed Brown           }
901c4762a1bSJed Brown         }
902c4762a1bSJed Brown       }
903c4762a1bSJed Brown     }
904c4762a1bSJed Brown   }
9055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,Xl,&x));
9065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(da,Xl,INSERT_VALUES,X));
9075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(da,Xl,INSERT_VALUES,X));
9085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&Xl));
909c4762a1bSJed Brown   if (B) {
9105f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,Bl,&b));
9115f80ce2aSJacob Faibussowitsch     CHKERRQ(DMRestoreLocalVector(da,&Bl));
912c4762a1bSJed Brown   }
9135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
914c4762a1bSJed Brown   PetscFunctionReturn(0);
915c4762a1bSJed Brown }
916c4762a1bSJed Brown 
917c4762a1bSJed Brown PetscErrorCode FormCoordinates(DM da,AppCtx *user)
918c4762a1bSJed Brown {
919c4762a1bSJed Brown   Vec            coords;
920c4762a1bSJed Brown   DM             cda;
921c4762a1bSJed Brown   PetscInt       mx,my,mz;
922c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,zs,xm,ym,zm;
923c4762a1bSJed Brown   CoordField     ***x;
924c4762a1bSJed Brown 
925c4762a1bSJed Brown   PetscFunctionBegin;
9265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
9275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(cda,&coords));
9285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
9295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
9305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,coords,&x));
931c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
932c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
933c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
934c4762a1bSJed Brown         PetscReal cx = ((PetscReal)i) / (((PetscReal)(mx-1)));
935c4762a1bSJed Brown         PetscReal cy = ((PetscReal)j) / (((PetscReal)(my-1)));
936c4762a1bSJed Brown         PetscReal cz = ((PetscReal)k) / (((PetscReal)(mz-1)));
937c4762a1bSJed Brown         PetscReal rad = user->rad + cy*user->height;
938c4762a1bSJed Brown         PetscReal ang = (cx - 0.5)*user->arc;
939c4762a1bSJed Brown         x[k][j][i][0] = rad*PetscSinReal(ang);
940c4762a1bSJed Brown         x[k][j][i][1] = rad*PetscCosReal(ang) - (user->rad + 0.5*user->height)*PetscCosReal(-0.5*user->arc);
941c4762a1bSJed Brown         x[k][j][i][2] = user->width*(cz - 0.5);
942c4762a1bSJed Brown       }
943c4762a1bSJed Brown     }
944c4762a1bSJed Brown   }
9455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,coords,&x));
9465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetCoordinates(da,coords));
9475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&coords));
948c4762a1bSJed Brown   PetscFunctionReturn(0);
949c4762a1bSJed Brown }
950c4762a1bSJed Brown 
951c4762a1bSJed Brown PetscErrorCode InitialGuess(DM da,AppCtx *user,Vec X)
952c4762a1bSJed Brown {
953c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,zs,xm,ym,zm;
954c4762a1bSJed Brown   PetscInt       mx,my,mz;
955c4762a1bSJed Brown   Field          ***x;
956c4762a1bSJed Brown 
957c4762a1bSJed Brown   PetscFunctionBegin;
9585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
9595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
9605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
961c4762a1bSJed Brown 
962c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
963c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
964c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
965c4762a1bSJed Brown         /* reference coordinates */
966c4762a1bSJed Brown         PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx-1)));
967c4762a1bSJed Brown         PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my-1)));
968c4762a1bSJed Brown         PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz-1)));
969c4762a1bSJed Brown         PetscReal o_x = p_x;
970c4762a1bSJed Brown         PetscReal o_y = p_y;
971c4762a1bSJed Brown         PetscReal o_z = p_z;
972c4762a1bSJed Brown         x[k][j][i][0] = o_x - p_x;
973c4762a1bSJed Brown         x[k][j][i][1] = o_y - p_y;
974c4762a1bSJed Brown         x[k][j][i][2] = o_z - p_z;
975c4762a1bSJed Brown       }
976c4762a1bSJed Brown     }
977c4762a1bSJed Brown   }
9785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
979c4762a1bSJed Brown   PetscFunctionReturn(0);
980c4762a1bSJed Brown }
981c4762a1bSJed Brown 
982c4762a1bSJed Brown PetscErrorCode FormRHS(DM da,AppCtx *user,Vec X)
983c4762a1bSJed Brown {
984c4762a1bSJed Brown   PetscInt       i,j,k,xs,ys,zs,xm,ym,zm;
985c4762a1bSJed Brown   PetscInt       mx,my,mz;
986c4762a1bSJed Brown   Field          ***x;
987c4762a1bSJed Brown 
988c4762a1bSJed Brown   PetscFunctionBegin;
9895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
9905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
9915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
992c4762a1bSJed Brown 
993c4762a1bSJed Brown   for (k=zs; k<zs+zm; k++) {
994c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
995c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
996c4762a1bSJed Brown         x[k][j][i][0] = 0.;
997c4762a1bSJed Brown         x[k][j][i][1] = 0.;
998c4762a1bSJed Brown         x[k][j][i][2] = 0.;
999c4762a1bSJed Brown         if (i == (mx-1)/2 && j == (my-1)) x[k][j][i][1] = user->ploading/(mz-1);
1000c4762a1bSJed Brown       }
1001c4762a1bSJed Brown     }
1002c4762a1bSJed Brown   }
10035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
1004c4762a1bSJed Brown   PetscFunctionReturn(0);
1005c4762a1bSJed Brown }
1006c4762a1bSJed Brown 
1007c4762a1bSJed Brown PetscErrorCode DisplayLine(SNES snes,Vec X)
1008c4762a1bSJed Brown {
1009c4762a1bSJed Brown   PetscInt       r,i,j=0,k=0,xs,xm,ys,ym,zs,zm,mx,my,mz;
1010c4762a1bSJed Brown   Field          ***x;
1011c4762a1bSJed Brown   CoordField     ***c;
1012c4762a1bSJed Brown   DM             da,cda;
1013c4762a1bSJed Brown   Vec            C;
1014c4762a1bSJed Brown   PetscMPIInt    size,rank;
1015c4762a1bSJed Brown 
1016c4762a1bSJed Brown   PetscFunctionBegin;
10175f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
10185f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
10195f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
10205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
10215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
10225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(da,&C));
1023c4762a1bSJed Brown   j = my / 2;
1024c4762a1bSJed Brown   k = mz / 2;
10255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
10265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
10275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,C,&c));
1028c4762a1bSJed Brown   for (r=0;r<size;r++) {
1029c4762a1bSJed Brown     if (rank == r) {
1030c4762a1bSJed Brown       if (j >= ys && j < ys+ym && k >= zs && k < zs+zm) {
1031c4762a1bSJed Brown         for (i=xs; i<xs+xm; i++) {
10325f80ce2aSJacob 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])));
1033c4762a1bSJed Brown         }
1034c4762a1bSJed Brown       }
1035c4762a1bSJed Brown     }
10365f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD));
1037c4762a1bSJed Brown   }
10385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
10395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,C,&c));
1040c4762a1bSJed Brown   PetscFunctionReturn(0);
1041c4762a1bSJed Brown }
1042c4762a1bSJed Brown 
1043c4762a1bSJed Brown /*TEST
1044c4762a1bSJed Brown 
1045c4762a1bSJed Brown    test:
1046c4762a1bSJed Brown       nsize: 2
1047c4762a1bSJed 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
1048c4762a1bSJed Brown       requires: !single
1049c4762a1bSJed Brown       timeoutfactor: 3
1050c4762a1bSJed Brown 
1051c4762a1bSJed Brown    test:
1052c4762a1bSJed Brown       suffix: 2
1053c4762a1bSJed 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
1054c4762a1bSJed Brown       requires: !single
1055c4762a1bSJed Brown 
1056c4762a1bSJed Brown    test:
1057c4762a1bSJed Brown       suffix: 3
1058c4762a1bSJed 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
1059c4762a1bSJed Brown       requires: !single
1060c4762a1bSJed Brown 
1061c4762a1bSJed Brown TEST*/
1062