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