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 102*327415f7SBarry Smith PetscFunctionBeginUser; 1039566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 1049566063dSJacob Faibussowitsch PetscCall(FormElements()); 105c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1069566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm,&snes)); 1079566063dSJacob 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)); 1089566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1099566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1109566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes,(DM)da)); 111c4762a1bSJed Brown 1129566063dSJacob Faibussowitsch PetscCall(SNESSetNGS(snes,NonlinearGS,&user)); 113c4762a1bSJed Brown 1149566063dSJacob 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)); 115c4762a1bSJed Brown user.loading = 0.0; 116c4762a1bSJed Brown user.arc = PETSC_PI/3.; 117c4762a1bSJed Brown user.mu = 4.0; 118c4762a1bSJed Brown user.lambda = 1.0; 119c4762a1bSJed Brown user.rad = 100.0; 120c4762a1bSJed Brown user.height = 3.; 121c4762a1bSJed Brown user.width = 1.; 122c4762a1bSJed Brown user.ploading = -5e3; 123c4762a1bSJed Brown 1249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-arc",&user.arc,NULL)); 1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,&muflg)); 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-lambda",&user.lambda,&lambdaflg)); 1279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-rad",&user.rad,NULL)); 1289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-height",&user.height,NULL)); 1299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-width",&user.width,NULL)); 1309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-loading",&user.loading,NULL)); 1319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-ploading",&user.ploading,NULL)); 1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-poisson",&poisson,&poissonflg)); 1339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-young",&young,&youngflg)); 134c4762a1bSJed Brown if ((youngflg || poissonflg) || !(muflg || lambdaflg)) { 135c4762a1bSJed Brown /* set the lame' parameters based upon the poisson ratio and young's modulus */ 136c4762a1bSJed Brown user.lambda = poisson*young / ((1. + poisson)*(1. - 2.*poisson)); 137c4762a1bSJed Brown user.mu = young/(2.*(1. + poisson)); 138c4762a1bSJed Brown } 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-view",&view,NULL)); 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_line",&viewline,NULL)); 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,0,"x_disp")); 1439566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,1,"y_disp")); 1449566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,2,"z_disp")); 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(da,&user)); 1479566063dSJacob Faibussowitsch PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user)); 1489566063dSJacob Faibussowitsch PetscCall(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user)); 1499566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 1509566063dSJacob Faibussowitsch PetscCall(FormCoordinates(da,&user)); 151c4762a1bSJed Brown 1529566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&x)); 1539566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&b)); 1549566063dSJacob Faibussowitsch PetscCall(InitialGuess(da,&user,x)); 1559566063dSJacob Faibussowitsch PetscCall(FormRHS(da,&user,b)); 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"lambda: %f mu: %f\n",(double)user.lambda,(double)user.mu)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* show a cross-section of the initial state */ 1601baa6e33SBarry Smith if (viewline) PetscCall(DisplayLine(snes,x)); 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* get the loaded configuration */ 1639566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,b,x)); 164c4762a1bSJed Brown 1659566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 16663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Number of SNES iterations = %" PetscInt_FMT "\n", its)); 1679566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes,&X)); 168c4762a1bSJed Brown /* show a cross-section of the final state */ 1691baa6e33SBarry Smith if (viewline) PetscCall(DisplayLine(snes,X)); 170c4762a1bSJed Brown 171c4762a1bSJed Brown if (view) { 172c4762a1bSJed Brown PetscViewer viewer; 173c4762a1bSJed Brown Vec coords; 1749566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&viewer)); 1759566063dSJacob Faibussowitsch PetscCall(VecView(x,viewer)); 1769566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&coords)); 1789566063dSJacob Faibussowitsch PetscCall(VecAXPY(coords,1.0,x)); 1799566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(PETSC_COMM_WORLD,filename_def,FILE_MODE_WRITE,&viewer)); 1809566063dSJacob Faibussowitsch PetscCall(VecView(x,viewer)); 1819566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1879566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1889566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 189b122ec5aSJacob Faibussowitsch return 0; 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192c4762a1bSJed Brown PetscInt OnBoundary(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz) 193c4762a1bSJed Brown { 194c4762a1bSJed Brown if ((i == 0 || i == mx-1) && j == my-1) return 1; 195c4762a1bSJed Brown return 0; 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown void BoundaryValue(PetscInt i,PetscInt j,PetscInt k,PetscInt mx,PetscInt my,PetscInt mz,PetscScalar *val,AppCtx *user) 199c4762a1bSJed Brown { 200c4762a1bSJed Brown /* reference coordinates */ 201c4762a1bSJed Brown PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx-1))); 202c4762a1bSJed Brown PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my-1))); 203c4762a1bSJed Brown PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz-1))); 204c4762a1bSJed Brown PetscReal o_x = p_x; 205c4762a1bSJed Brown PetscReal o_y = p_y; 206c4762a1bSJed Brown PetscReal o_z = p_z; 207c4762a1bSJed Brown val[0] = o_x - p_x; 208c4762a1bSJed Brown val[1] = o_y - p_y; 209c4762a1bSJed Brown val[2] = o_z - p_z; 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown void InvertTensor(PetscScalar *t, PetscScalar *ti,PetscReal *dett) 213c4762a1bSJed Brown { 214c4762a1bSJed Brown const PetscScalar a = t[0]; 215c4762a1bSJed Brown const PetscScalar b = t[1]; 216c4762a1bSJed Brown const PetscScalar c = t[2]; 217c4762a1bSJed Brown const PetscScalar d = t[3]; 218c4762a1bSJed Brown const PetscScalar e = t[4]; 219c4762a1bSJed Brown const PetscScalar f = t[5]; 220c4762a1bSJed Brown const PetscScalar g = t[6]; 221c4762a1bSJed Brown const PetscScalar h = t[7]; 222c4762a1bSJed Brown const PetscScalar i = t[8]; 223c4762a1bSJed Brown const PetscReal det = PetscRealPart(a*(e*i - f*h) - b*(i*d - f*g) + c*(d*h - e*g)); 224c4762a1bSJed Brown const PetscReal di = 1. / det; 225c4762a1bSJed Brown if (ti) { 226c4762a1bSJed Brown const PetscScalar A = (e*i - f*h); 227c4762a1bSJed Brown const PetscScalar B = -(d*i - f*g); 228c4762a1bSJed Brown const PetscScalar C = (d*h - e*g); 229c4762a1bSJed Brown const PetscScalar D = -(b*i - c*h); 230c4762a1bSJed Brown const PetscScalar E = (a*i - c*g); 231c4762a1bSJed Brown const PetscScalar F = -(a*h - b*g); 232c4762a1bSJed Brown const PetscScalar G = (b*f - c*e); 233c4762a1bSJed Brown const PetscScalar H = -(a*f - c*d); 234c4762a1bSJed Brown const PetscScalar II = (a*e - b*d); 235c4762a1bSJed Brown ti[0] = di*A; 236c4762a1bSJed Brown ti[1] = di*D; 237c4762a1bSJed Brown ti[2] = di*G; 238c4762a1bSJed Brown ti[3] = di*B; 239c4762a1bSJed Brown ti[4] = di*E; 240c4762a1bSJed Brown ti[5] = di*H; 241c4762a1bSJed Brown ti[6] = di*C; 242c4762a1bSJed Brown ti[7] = di*F; 243c4762a1bSJed Brown ti[8] = di*II; 244c4762a1bSJed Brown } 245c4762a1bSJed Brown if (dett) *dett = det; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown void TensorTensor(PetscScalar *a,PetscScalar *b,PetscScalar *c) 249c4762a1bSJed Brown { 250c4762a1bSJed Brown PetscInt i,j,m; 251c4762a1bSJed Brown for (i=0;i<3;i++) { 252c4762a1bSJed Brown for (j=0;j<3;j++) { 253c4762a1bSJed Brown c[i+3*j] = 0; 254c4762a1bSJed Brown for (m=0;m<3;m++) 255c4762a1bSJed Brown c[i+3*j] += a[m+3*j]*b[i+3*m]; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown void TensorTransposeTensor(PetscScalar *a,PetscScalar *b,PetscScalar *c) 261c4762a1bSJed Brown { 262c4762a1bSJed Brown PetscInt i,j,m; 263c4762a1bSJed Brown for (i=0;i<3;i++) { 264c4762a1bSJed Brown for (j=0;j<3;j++) { 265c4762a1bSJed Brown c[i+3*j] = 0; 266c4762a1bSJed Brown for (m=0;m<3;m++) 267c4762a1bSJed Brown c[i+3*j] += a[3*m+j]*b[i+3*m]; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown } 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown void TensorVector(PetscScalar *rot, PetscScalar *vec, PetscScalar *tvec) 273c4762a1bSJed Brown { 274c4762a1bSJed Brown tvec[0] = rot[0]*vec[0] + rot[1]*vec[1] + rot[2]*vec[2]; 275c4762a1bSJed Brown tvec[1] = rot[3]*vec[0] + rot[4]*vec[1] + rot[5]*vec[2]; 276c4762a1bSJed Brown tvec[2] = rot[6]*vec[0] + rot[7]*vec[1] + rot[8]*vec[2]; 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown void DeformationGradient(Field *ex,PetscInt qi,PetscInt qj,PetscInt qk,PetscScalar *invJ,PetscScalar *F) 280c4762a1bSJed Brown { 281c4762a1bSJed Brown PetscInt ii,jj,kk,l; 282c4762a1bSJed Brown for (l = 0; l < 9; l++) { 283c4762a1bSJed Brown F[l] = 0.; 284c4762a1bSJed Brown } 285c4762a1bSJed Brown F[0] = 1.; 286c4762a1bSJed Brown F[4] = 1.; 287c4762a1bSJed Brown F[8] = 1.; 288c4762a1bSJed Brown /* form the deformation gradient at this basis function -- loop over element unknowns */ 289c4762a1bSJed Brown for (kk=0;kk<NB;kk++) { 290c4762a1bSJed Brown for (jj=0;jj<NB;jj++) { 291c4762a1bSJed Brown for (ii=0;ii<NB;ii++) { 292c4762a1bSJed Brown PetscInt idx = ii + jj*NB + kk*NB*NB; 293c4762a1bSJed Brown PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk; 294c4762a1bSJed Brown PetscScalar lgrad[3]; 295c4762a1bSJed Brown TensorVector(invJ,&grad[3*bidx],lgrad); 296c4762a1bSJed Brown F[0] += lgrad[0]*ex[idx][0]; F[1] += lgrad[1]*ex[idx][0]; F[2] += lgrad[2]*ex[idx][0]; 297c4762a1bSJed Brown F[3] += lgrad[0]*ex[idx][1]; F[4] += lgrad[1]*ex[idx][1]; F[5] += lgrad[2]*ex[idx][1]; 298c4762a1bSJed Brown F[6] += lgrad[0]*ex[idx][2]; F[7] += lgrad[1]*ex[idx][2]; F[8] += lgrad[2]*ex[idx][2]; 299c4762a1bSJed Brown } 300c4762a1bSJed Brown } 301c4762a1bSJed Brown } 302c4762a1bSJed Brown } 303c4762a1bSJed Brown 304c4762a1bSJed Brown void DeformationGradientJacobian(PetscInt qi,PetscInt qj,PetscInt qk,PetscInt ii,PetscInt jj,PetscInt kk,PetscInt fld,PetscScalar *invJ,PetscScalar *dF) 305c4762a1bSJed Brown { 306c4762a1bSJed Brown PetscInt l; 307c4762a1bSJed Brown PetscScalar lgrad[3]; 308c4762a1bSJed Brown PetscInt idx = ii + jj*NB + kk*NB*NB; 309c4762a1bSJed Brown PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk; 310c4762a1bSJed Brown for (l = 0; l < 9; l++) { 311c4762a1bSJed Brown dF[l] = 0.; 312c4762a1bSJed Brown } 313c4762a1bSJed Brown /* form the deformation gradient at this basis function -- loop over element unknowns */ 314c4762a1bSJed Brown TensorVector(invJ,&grad[3*bidx],lgrad); 315c4762a1bSJed Brown dF[3*fld] = lgrad[0]; dF[3*fld + 1] = lgrad[1]; dF[3*fld + 2] = lgrad[2]; 316c4762a1bSJed Brown } 317c4762a1bSJed Brown 318c4762a1bSJed Brown void LagrangeGreenStrain(PetscScalar *F,PetscScalar *E) 319c4762a1bSJed Brown { 320c4762a1bSJed Brown PetscInt i,j,m; 321c4762a1bSJed Brown for (i=0;i<3;i++) { 322c4762a1bSJed Brown for (j=0;j<3;j++) { 323c4762a1bSJed Brown E[i+3*j] = 0; 324c4762a1bSJed Brown for (m=0;m<3;m++) 325c4762a1bSJed Brown E[i+3*j] += 0.5*F[3*m+j]*F[i+3*m]; 326c4762a1bSJed Brown } 327c4762a1bSJed Brown } 328c4762a1bSJed Brown for (i=0;i<3;i++) { 329c4762a1bSJed Brown E[i+3*i] -= 0.5; 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown void SaintVenantKirchoff(PetscReal lambda,PetscReal mu,PetscScalar *F,PetscScalar *S) 334c4762a1bSJed Brown { 335c4762a1bSJed Brown PetscInt i,j; 336c4762a1bSJed Brown PetscScalar E[9]; 337c4762a1bSJed Brown PetscScalar trE=0; 338c4762a1bSJed Brown LagrangeGreenStrain(F,E); 339c4762a1bSJed Brown for (i=0;i<3;i++) { 340c4762a1bSJed Brown trE += E[i+3*i]; 341c4762a1bSJed Brown } 342c4762a1bSJed Brown for (i=0;i<3;i++) { 343c4762a1bSJed Brown for (j=0;j<3;j++) { 344c4762a1bSJed Brown S[i+3*j] = 2.*mu*E[i+3*j]; 345c4762a1bSJed Brown if (i == j) { 346c4762a1bSJed Brown S[i+3*j] += trE*lambda; 347c4762a1bSJed Brown } 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 350c4762a1bSJed Brown } 351c4762a1bSJed Brown 352c4762a1bSJed Brown void SaintVenantKirchoffJacobian(PetscReal lambda,PetscReal mu,PetscScalar *F,PetscScalar *dF,PetscScalar *dS) 353c4762a1bSJed Brown { 354c4762a1bSJed Brown PetscScalar FtdF[9],dE[9]; 355c4762a1bSJed Brown PetscInt i,j; 356c4762a1bSJed Brown PetscScalar dtrE=0.; 357c4762a1bSJed Brown TensorTransposeTensor(dF,F,dE); 358c4762a1bSJed Brown TensorTransposeTensor(F,dF,FtdF); 359c4762a1bSJed Brown for (i=0;i<9;i++) dE[i] += FtdF[i]; 360c4762a1bSJed Brown for (i=0;i<9;i++) dE[i] *= 0.5; 361c4762a1bSJed Brown for (i=0;i<3;i++) { 362c4762a1bSJed Brown dtrE += dE[i+3*i]; 363c4762a1bSJed Brown } 364c4762a1bSJed Brown for (i=0;i<3;i++) { 365c4762a1bSJed Brown for (j=0;j<3;j++) { 366c4762a1bSJed Brown dS[i+3*j] = 2.*mu*dE[i+3*j]; 367c4762a1bSJed Brown if (i == j) { 368c4762a1bSJed Brown dS[i+3*j] += dtrE*lambda; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown } 371c4762a1bSJed Brown } 372c4762a1bSJed Brown } 373c4762a1bSJed Brown 374c4762a1bSJed Brown PetscErrorCode FormElements() 375c4762a1bSJed Brown { 376c4762a1bSJed Brown PetscInt i,j,k,ii,jj,kk; 377c4762a1bSJed Brown PetscReal bx,by,bz,dbx,dby,dbz; 378c4762a1bSJed Brown 379c4762a1bSJed Brown PetscFunctionBegin; 380c4762a1bSJed Brown /* construct the basis function values and derivatives */ 381c4762a1bSJed Brown for (k = 0; k < NB; k++) { 382c4762a1bSJed Brown for (j = 0; j < NB; j++) { 383c4762a1bSJed Brown for (i = 0; i < NB; i++) { 384c4762a1bSJed Brown /* loop over the quadrature points */ 385c4762a1bSJed Brown for (kk = 0; kk < NQ; kk++) { 386c4762a1bSJed Brown for (jj = 0; jj < NQ; jj++) { 387c4762a1bSJed Brown for (ii = 0; ii < NQ; ii++) { 388c4762a1bSJed Brown PetscInt idx = ii + NQ*jj + NQ*NQ*kk + NEQ*i + NEQ*NB*j + NEQ*NB*NB*k; 389c4762a1bSJed Brown bx = pts[ii]; 390c4762a1bSJed Brown by = pts[jj]; 391c4762a1bSJed Brown bz = pts[kk]; 392c4762a1bSJed Brown dbx = 1.; 393c4762a1bSJed Brown dby = 1.; 394c4762a1bSJed Brown dbz = 1.; 395c4762a1bSJed Brown if (i == 0) {bx = 1. - bx; dbx = -1;} 396c4762a1bSJed Brown if (j == 0) {by = 1. - by; dby = -1;} 397c4762a1bSJed Brown if (k == 0) {bz = 1. - bz; dbz = -1;} 398c4762a1bSJed Brown vals[idx] = bx*by*bz; 399c4762a1bSJed Brown grad[3*idx + 0] = dbx*by*bz; 400c4762a1bSJed Brown grad[3*idx + 1] = dby*bx*bz; 401c4762a1bSJed Brown grad[3*idx + 2] = dbz*bx*by; 402c4762a1bSJed Brown } 403c4762a1bSJed Brown } 404c4762a1bSJed Brown } 405c4762a1bSJed Brown } 406c4762a1bSJed Brown } 407c4762a1bSJed Brown } 408c4762a1bSJed Brown PetscFunctionReturn(0); 409c4762a1bSJed Brown } 410c4762a1bSJed Brown 411c4762a1bSJed 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) 412c4762a1bSJed Brown { 413c4762a1bSJed Brown PetscInt m; 414c4762a1bSJed Brown PetscInt ii,jj,kk; 415c4762a1bSJed Brown /* gather the data -- loop over element unknowns */ 416c4762a1bSJed Brown for (kk=0;kk<NB;kk++) { 417c4762a1bSJed Brown for (jj=0;jj<NB;jj++) { 418c4762a1bSJed Brown for (ii=0;ii<NB;ii++) { 419c4762a1bSJed Brown PetscInt idx = ii + jj*NB + kk*NB*NB; 420c4762a1bSJed Brown /* decouple the boundary nodes for the displacement variables */ 421c4762a1bSJed Brown if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz)) { 422c4762a1bSJed Brown BoundaryValue(i+ii,j+jj,k+kk,mx,my,mz,ex[idx],user); 423c4762a1bSJed Brown } else { 424c4762a1bSJed Brown for (m=0;m<3;m++) { 425c4762a1bSJed Brown ex[idx][m] = x[k+kk][j+jj][i+ii][m]; 426c4762a1bSJed Brown } 427c4762a1bSJed Brown } 428c4762a1bSJed Brown for (m=0;m<3;m++) { 429c4762a1bSJed Brown ec[idx][m] = c[k+kk][j+jj][i+ii][m]; 430c4762a1bSJed Brown } 431c4762a1bSJed Brown } 432c4762a1bSJed Brown } 433c4762a1bSJed Brown } 434c4762a1bSJed Brown } 435c4762a1bSJed Brown 436c4762a1bSJed Brown void QuadraturePointGeometricJacobian(CoordField *ec,PetscInt qi,PetscInt qj,PetscInt qk, PetscScalar *J) 437c4762a1bSJed Brown { 438c4762a1bSJed Brown PetscInt ii,jj,kk; 439c4762a1bSJed Brown /* construct the gradient at the given quadrature point named by i,j,k */ 440c4762a1bSJed Brown for (ii = 0; ii < 9; ii++) { 441c4762a1bSJed Brown J[ii] = 0; 442c4762a1bSJed Brown } 443c4762a1bSJed Brown for (kk = 0; kk < NB; kk++) { 444c4762a1bSJed Brown for (jj = 0; jj < NB; jj++) { 445c4762a1bSJed Brown for (ii = 0; ii < NB; ii++) { 446c4762a1bSJed Brown PetscInt idx = ii + jj*NB + kk*NB*NB; 447c4762a1bSJed Brown PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk; 448c4762a1bSJed 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]; 449c4762a1bSJed 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]; 450c4762a1bSJed 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]; 451c4762a1bSJed Brown } 452c4762a1bSJed Brown } 453c4762a1bSJed Brown } 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 456c4762a1bSJed Brown void FormElementJacobian(Field *ex,CoordField *ec,Field *ef,PetscScalar *ej,AppCtx *user) 457c4762a1bSJed Brown { 458c4762a1bSJed Brown PetscReal vol; 459c4762a1bSJed Brown PetscScalar J[9]; 460c4762a1bSJed Brown PetscScalar invJ[9]; 461c4762a1bSJed Brown PetscScalar F[9],S[9],dF[9],dS[9],dFS[9],FdS[9],FS[9]; 462c4762a1bSJed Brown PetscReal scl; 463c4762a1bSJed Brown PetscInt i,j,k,l,ii,jj,kk,ll,qi,qj,qk,m; 464c4762a1bSJed Brown 465c4762a1bSJed Brown if (ej) for (i = 0; i < NPB*NPB; i++) ej[i] = 0.; 466c4762a1bSJed Brown if (ef) for (i = 0; i < NEB; i++) {ef[i][0] = 0.;ef[i][1] = 0.;ef[i][2] = 0.;} 467c4762a1bSJed Brown /* loop over quadrature */ 468c4762a1bSJed Brown for (qk = 0; qk < NQ; qk++) { 469c4762a1bSJed Brown for (qj = 0; qj < NQ; qj++) { 470c4762a1bSJed Brown for (qi = 0; qi < NQ; qi++) { 471c4762a1bSJed Brown QuadraturePointGeometricJacobian(ec,qi,qj,qk,J); 472c4762a1bSJed Brown InvertTensor(J,invJ,&vol); 473c4762a1bSJed Brown scl = vol*wts[qi]*wts[qj]*wts[qk]; 474c4762a1bSJed Brown DeformationGradient(ex,qi,qj,qk,invJ,F); 475c4762a1bSJed Brown SaintVenantKirchoff(user->lambda,user->mu,F,S); 476c4762a1bSJed Brown /* form the function */ 477c4762a1bSJed Brown if (ef) { 478c4762a1bSJed Brown TensorTensor(F,S,FS); 479c4762a1bSJed Brown for (kk=0;kk<NB;kk++) { 480c4762a1bSJed Brown for (jj=0;jj<NB;jj++) { 481c4762a1bSJed Brown for (ii=0;ii<NB;ii++) { 482c4762a1bSJed Brown PetscInt idx = ii + jj*NB + kk*NB*NB; 483c4762a1bSJed Brown PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk; 484c4762a1bSJed Brown PetscScalar lgrad[3]; 485c4762a1bSJed Brown TensorVector(invJ,&grad[3*bidx],lgrad); 486c4762a1bSJed Brown /* mu*F : grad phi_{u,v,w} */ 487c4762a1bSJed Brown for (m=0;m<3;m++) { 488c4762a1bSJed Brown ef[idx][m] += scl* 489c4762a1bSJed Brown (lgrad[0]*FS[3*m + 0] + lgrad[1]*FS[3*m + 1] + lgrad[2]*FS[3*m + 2]); 490c4762a1bSJed Brown } 491c4762a1bSJed Brown ef[idx][1] -= scl*user->loading*vals[bidx]; 492c4762a1bSJed Brown } 493c4762a1bSJed Brown } 494c4762a1bSJed Brown } 495c4762a1bSJed Brown } 496c4762a1bSJed Brown /* form the jacobian */ 497c4762a1bSJed Brown if (ej) { 498c4762a1bSJed Brown /* loop over trialfunctions */ 499c4762a1bSJed Brown for (k=0;k<NB;k++) { 500c4762a1bSJed Brown for (j=0;j<NB;j++) { 501c4762a1bSJed Brown for (i=0;i<NB;i++) { 502c4762a1bSJed Brown for (l=0;l<3;l++) { 503c4762a1bSJed Brown PetscInt tridx = l + 3*(i + j*NB + k*NB*NB); 504c4762a1bSJed Brown DeformationGradientJacobian(qi,qj,qk,i,j,k,l,invJ,dF); 505c4762a1bSJed Brown SaintVenantKirchoffJacobian(user->lambda,user->mu,F,dF,dS); 506c4762a1bSJed Brown TensorTensor(dF,S,dFS); 507c4762a1bSJed Brown TensorTensor(F,dS,FdS); 508c4762a1bSJed Brown for (m=0;m<9;m++) dFS[m] += FdS[m]; 509c4762a1bSJed Brown /* loop over testfunctions */ 510c4762a1bSJed Brown for (kk=0;kk<NB;kk++) { 511c4762a1bSJed Brown for (jj=0;jj<NB;jj++) { 512c4762a1bSJed Brown for (ii=0;ii<NB;ii++) { 513c4762a1bSJed Brown PetscInt idx = ii + jj*NB + kk*NB*NB; 514c4762a1bSJed Brown PetscInt bidx = 8*idx + qi + NQ*qj + NQ*NQ*qk; 515c4762a1bSJed Brown PetscScalar lgrad[3]; 516c4762a1bSJed Brown TensorVector(invJ,&grad[3*bidx],lgrad); 517c4762a1bSJed Brown for (ll=0; ll<3;ll++) { 518c4762a1bSJed Brown PetscInt teidx = ll + 3*(ii + jj*NB + kk*NB*NB); 519c4762a1bSJed Brown ej[teidx + NPB*tridx] += scl* 520c4762a1bSJed Brown (lgrad[0]*dFS[3*ll + 0] + lgrad[1]*dFS[3*ll + 1] + lgrad[2]*dFS[3*ll+2]); 521c4762a1bSJed Brown } 522c4762a1bSJed Brown } 523c4762a1bSJed Brown } 524c4762a1bSJed Brown } /* end of testfunctions */ 525c4762a1bSJed Brown } 526c4762a1bSJed Brown } 527c4762a1bSJed Brown } 528c4762a1bSJed Brown } /* end of trialfunctions */ 529c4762a1bSJed Brown } 530c4762a1bSJed Brown } 531c4762a1bSJed Brown } 532c4762a1bSJed Brown } /* end of quadrature points */ 533c4762a1bSJed Brown } 534c4762a1bSJed Brown 535c4762a1bSJed Brown void FormPBJacobian(PetscInt i,PetscInt j,PetscInt k,Field *ex,CoordField *ec,Field *ef,PetscScalar *ej,AppCtx *user) 536c4762a1bSJed Brown { 537c4762a1bSJed Brown PetscReal vol; 538c4762a1bSJed Brown PetscScalar J[9]; 539c4762a1bSJed Brown PetscScalar invJ[9]; 540c4762a1bSJed Brown PetscScalar F[9],S[9],dF[9],dS[9],dFS[9],FdS[9],FS[9]; 541c4762a1bSJed Brown PetscReal scl; 542c4762a1bSJed Brown PetscInt l,ll,qi,qj,qk,m; 543c4762a1bSJed Brown PetscInt idx = i + j*NB + k*NB*NB; 544c4762a1bSJed Brown PetscScalar lgrad[3]; 545c4762a1bSJed Brown 546c4762a1bSJed Brown if (ej) for (l = 0; l < 9; l++) ej[l] = 0.; 547c4762a1bSJed Brown if (ef) for (l = 0; l < 1; l++) {ef[l][0] = 0.;ef[l][1] = 0.;ef[l][2] = 0.;} 548c4762a1bSJed Brown /* loop over quadrature */ 549c4762a1bSJed Brown for (qk = 0; qk < NQ; qk++) { 550c4762a1bSJed Brown for (qj = 0; qj < NQ; qj++) { 551c4762a1bSJed Brown for (qi = 0; qi < NQ; qi++) { 552c4762a1bSJed Brown PetscInt bidx = NEB*idx + qi + NQ*qj + NQ*NQ*qk; 553c4762a1bSJed Brown QuadraturePointGeometricJacobian(ec,qi,qj,qk,J); 554c4762a1bSJed Brown InvertTensor(J,invJ,&vol); 555c4762a1bSJed Brown TensorVector(invJ,&grad[3*bidx],lgrad); 556c4762a1bSJed Brown scl = vol*wts[qi]*wts[qj]*wts[qk]; 557c4762a1bSJed Brown DeformationGradient(ex,qi,qj,qk,invJ,F); 558c4762a1bSJed Brown SaintVenantKirchoff(user->lambda,user->mu,F,S); 559c4762a1bSJed Brown /* form the function */ 560c4762a1bSJed Brown if (ef) { 561c4762a1bSJed Brown TensorTensor(F,S,FS); 562c4762a1bSJed Brown for (m=0;m<3;m++) { 563c4762a1bSJed Brown ef[0][m] += scl* 564c4762a1bSJed Brown (lgrad[0]*FS[3*m + 0] + lgrad[1]*FS[3*m + 1] + lgrad[2]*FS[3*m + 2]); 565c4762a1bSJed Brown } 566c4762a1bSJed Brown ef[0][1] -= scl*user->loading*vals[bidx]; 567c4762a1bSJed Brown } 568c4762a1bSJed Brown /* form the jacobian */ 569c4762a1bSJed Brown if (ej) { 570c4762a1bSJed Brown for (l=0;l<3;l++) { 571c4762a1bSJed Brown DeformationGradientJacobian(qi,qj,qk,i,j,k,l,invJ,dF); 572c4762a1bSJed Brown SaintVenantKirchoffJacobian(user->lambda,user->mu,F,dF,dS); 573c4762a1bSJed Brown TensorTensor(dF,S,dFS); 574c4762a1bSJed Brown TensorTensor(F,dS,FdS); 575c4762a1bSJed Brown for (m=0;m<9;m++) dFS[m] += FdS[m]; 576c4762a1bSJed Brown for (ll=0; ll<3;ll++) { 577c4762a1bSJed Brown ej[ll + 3*l] += scl* 578c4762a1bSJed Brown (lgrad[0]*dFS[3*ll + 0] + lgrad[1]*dFS[3*ll + 1] + lgrad[2]*dFS[3*ll+2]); 579c4762a1bSJed Brown } 580c4762a1bSJed Brown } 581c4762a1bSJed Brown } 582c4762a1bSJed Brown } 583c4762a1bSJed Brown } /* end of quadrature points */ 584c4762a1bSJed Brown } 585c4762a1bSJed Brown } 586c4762a1bSJed Brown 587c4762a1bSJed Brown void ApplyBCsElement(PetscInt mx,PetscInt my, PetscInt mz, PetscInt i, PetscInt j, PetscInt k,PetscScalar *jacobian) 588c4762a1bSJed Brown { 589c4762a1bSJed Brown PetscInt ii,jj,kk,ll,ei,ej,ek,el; 590c4762a1bSJed Brown for (kk=0;kk<NB;kk++) { 591c4762a1bSJed Brown for (jj=0;jj<NB;jj++) { 592c4762a1bSJed Brown for (ii=0;ii<NB;ii++) { 593c4762a1bSJed Brown for (ll = 0;ll<3;ll++) { 594c4762a1bSJed Brown PetscInt tridx = ll + 3*(ii + jj*NB + kk*NB*NB); 595c4762a1bSJed Brown for (ek=0;ek<NB;ek++) { 596c4762a1bSJed Brown for (ej=0;ej<NB;ej++) { 597c4762a1bSJed Brown for (ei=0;ei<NB;ei++) { 598c4762a1bSJed Brown for (el=0;el<3;el++) { 599c4762a1bSJed Brown if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz) || OnBoundary(i+ei,j+ej,k+ek,mx,my,mz)) { 600c4762a1bSJed Brown PetscInt teidx = el + 3*(ei + ej*NB + ek*NB*NB); 601c4762a1bSJed Brown if (teidx == tridx) { 602c4762a1bSJed Brown jacobian[tridx + NPB*teidx] = 1.; 603c4762a1bSJed Brown } else { 604c4762a1bSJed Brown jacobian[tridx + NPB*teidx] = 0.; 605c4762a1bSJed Brown } 606c4762a1bSJed Brown } 607c4762a1bSJed Brown } 608c4762a1bSJed Brown } 609c4762a1bSJed Brown } 610c4762a1bSJed Brown } 611c4762a1bSJed Brown } 612c4762a1bSJed Brown } 613c4762a1bSJed Brown } 614c4762a1bSJed Brown } 615c4762a1bSJed Brown } 616c4762a1bSJed Brown 617c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,Field ***x,Mat jacpre,Mat jac,void *ptr) 618c4762a1bSJed Brown { 619c4762a1bSJed Brown /* values for each basis function at each quadrature point */ 620c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 621c4762a1bSJed Brown PetscInt i,j,k,m,l; 622c4762a1bSJed Brown PetscInt ii,jj,kk; 623c4762a1bSJed Brown PetscScalar ej[NPB*NPB]; 624c4762a1bSJed Brown PetscScalar vals[NPB*NPB]; 625c4762a1bSJed Brown Field ex[NEB]; 626c4762a1bSJed Brown CoordField ec[NEB]; 627c4762a1bSJed Brown 628c4762a1bSJed Brown PetscInt xs=info->xs,ys=info->ys,zs=info->zs; 629c4762a1bSJed Brown PetscInt xm=info->xm,ym=info->ym,zm=info->zm; 630c4762a1bSJed Brown PetscInt xes,yes,zes,xee,yee,zee; 631c4762a1bSJed Brown PetscInt mx=info->mx,my=info->my,mz=info->mz; 632c4762a1bSJed Brown DM cda; 633c4762a1bSJed Brown CoordField ***c; 634c4762a1bSJed Brown Vec C; 635c4762a1bSJed Brown PetscInt nrows; 636c4762a1bSJed Brown MatStencil col[NPB],row[NPB]; 637c4762a1bSJed Brown PetscScalar v[9]; 638c4762a1bSJed Brown 639c4762a1bSJed Brown PetscFunctionBegin; 6409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(info->da,&cda)); 6419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(info->da,&C)); 6429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda,C,&c)); 6439566063dSJacob Faibussowitsch PetscCall(MatScale(jac,0.0)); 644c4762a1bSJed Brown 645c4762a1bSJed Brown xes = xs; 646c4762a1bSJed Brown yes = ys; 647c4762a1bSJed Brown zes = zs; 648c4762a1bSJed Brown xee = xs+xm; 649c4762a1bSJed Brown yee = ys+ym; 650c4762a1bSJed Brown zee = zs+zm; 651c4762a1bSJed Brown if (xs > 0) xes = xs-1; 652c4762a1bSJed Brown if (ys > 0) yes = ys-1; 653c4762a1bSJed Brown if (zs > 0) zes = zs-1; 654c4762a1bSJed Brown if (xs+xm == mx) xee = xs+xm-1; 655c4762a1bSJed Brown if (ys+ym == my) yee = ys+ym-1; 656c4762a1bSJed Brown if (zs+zm == mz) zee = zs+zm-1; 657c4762a1bSJed Brown for (k=zes; k<zee; k++) { 658c4762a1bSJed Brown for (j=yes; j<yee; j++) { 659c4762a1bSJed Brown for (i=xes; i<xee; i++) { 660c4762a1bSJed Brown GatherElementData(mx,my,mz,x,c,i,j,k,ex,ec,user); 661c4762a1bSJed Brown FormElementJacobian(ex,ec,NULL,ej,user); 662c4762a1bSJed Brown ApplyBCsElement(mx,my,mz,i,j,k,ej); 663c4762a1bSJed Brown nrows = 0.; 664c4762a1bSJed Brown for (kk=0;kk<NB;kk++) { 665c4762a1bSJed Brown for (jj=0;jj<NB;jj++) { 666c4762a1bSJed Brown for (ii=0;ii<NB;ii++) { 667c4762a1bSJed Brown PetscInt idx = ii + jj*2 + kk*4; 668c4762a1bSJed Brown for (m=0;m<3;m++) { 669c4762a1bSJed Brown col[3*idx+m].i = i+ii; 670c4762a1bSJed Brown col[3*idx+m].j = j+jj; 671c4762a1bSJed Brown col[3*idx+m].k = k+kk; 672c4762a1bSJed Brown col[3*idx+m].c = m; 673c4762a1bSJed Brown if (i+ii >= xs && i+ii < xm+xs && j+jj >= ys && j+jj < ys+ym && k+kk >= zs && k+kk < zs+zm) { 674c4762a1bSJed Brown row[nrows].i = i+ii; 675c4762a1bSJed Brown row[nrows].j = j+jj; 676c4762a1bSJed Brown row[nrows].k = k+kk; 677c4762a1bSJed Brown row[nrows].c = m; 678c4762a1bSJed Brown for (l=0;l<NPB;l++) vals[NPB*nrows + l] = ej[NPB*(3*idx+m) + l]; 679c4762a1bSJed Brown nrows++; 680c4762a1bSJed Brown } 681c4762a1bSJed Brown } 682c4762a1bSJed Brown } 683c4762a1bSJed Brown } 684c4762a1bSJed Brown } 6859566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,nrows,row,NPB,col,vals,ADD_VALUES)); 686c4762a1bSJed Brown } 687c4762a1bSJed Brown } 688c4762a1bSJed Brown } 689c4762a1bSJed Brown 6909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FLUSH_ASSEMBLY)); 6919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FLUSH_ASSEMBLY)); 692c4762a1bSJed Brown 693c4762a1bSJed Brown /* set the diagonal */ 694c4762a1bSJed 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.; 695c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 696c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 697c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 698c4762a1bSJed Brown if (OnBoundary(i,j,k,mx,my,mz)) { 699c4762a1bSJed Brown for (m=0; m<3;m++) { 700c4762a1bSJed Brown col[m].i = i; 701c4762a1bSJed Brown col[m].j = j; 702c4762a1bSJed Brown col[m].k = k; 703c4762a1bSJed Brown col[m].c = m; 704c4762a1bSJed Brown } 7059566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(jac,3,col,3,col,v,INSERT_VALUES)); 706c4762a1bSJed Brown } 707c4762a1bSJed Brown } 708c4762a1bSJed Brown } 709c4762a1bSJed Brown } 710c4762a1bSJed Brown 7119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 7129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 713c4762a1bSJed Brown 7149566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda,C,&c)); 715c4762a1bSJed Brown PetscFunctionReturn(0); 716c4762a1bSJed Brown } 717c4762a1bSJed Brown 718c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field ***x,Field ***f,void *ptr) 719c4762a1bSJed Brown { 720c4762a1bSJed Brown /* values for each basis function at each quadrature point */ 721c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 722c4762a1bSJed Brown PetscInt i,j,k,l; 723c4762a1bSJed Brown PetscInt ii,jj,kk; 724c4762a1bSJed Brown 725c4762a1bSJed Brown Field ef[NEB]; 726c4762a1bSJed Brown Field ex[NEB]; 727c4762a1bSJed Brown CoordField ec[NEB]; 728c4762a1bSJed Brown 729c4762a1bSJed Brown PetscInt xs=info->xs,ys=info->ys,zs=info->zs; 730c4762a1bSJed Brown PetscInt xm=info->xm,ym=info->ym,zm=info->zm; 731c4762a1bSJed Brown PetscInt xes,yes,zes,xee,yee,zee; 732c4762a1bSJed Brown PetscInt mx=info->mx,my=info->my,mz=info->mz; 733c4762a1bSJed Brown DM cda; 734c4762a1bSJed Brown CoordField ***c; 735c4762a1bSJed Brown Vec C; 736c4762a1bSJed Brown 737c4762a1bSJed Brown PetscFunctionBegin; 7389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(info->da,&cda)); 7399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(info->da,&C)); 7409566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda,C,&c)); 7419566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(info->da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 7429566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(info->da,&xs,&ys,&zs,&xm,&ym,&zm)); 743c4762a1bSJed Brown 744c4762a1bSJed Brown /* loop over elements */ 745c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 746c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 747c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 748c4762a1bSJed Brown for (l=0;l<3;l++) { 749c4762a1bSJed Brown f[k][j][i][l] = 0.; 750c4762a1bSJed Brown } 751c4762a1bSJed Brown } 752c4762a1bSJed Brown } 753c4762a1bSJed Brown } 754c4762a1bSJed Brown /* element starts and ends */ 755c4762a1bSJed Brown xes = xs; 756c4762a1bSJed Brown yes = ys; 757c4762a1bSJed Brown zes = zs; 758c4762a1bSJed Brown xee = xs+xm; 759c4762a1bSJed Brown yee = ys+ym; 760c4762a1bSJed Brown zee = zs+zm; 761c4762a1bSJed Brown if (xs > 0) xes = xs - 1; 762c4762a1bSJed Brown if (ys > 0) yes = ys - 1; 763c4762a1bSJed Brown if (zs > 0) zes = zs - 1; 764c4762a1bSJed Brown if (xs+xm == mx) xee = xs+xm-1; 765c4762a1bSJed Brown if (ys+ym == my) yee = ys+ym-1; 766c4762a1bSJed Brown if (zs+zm == mz) zee = zs+zm-1; 767c4762a1bSJed Brown for (k=zes; k<zee; k++) { 768c4762a1bSJed Brown for (j=yes; j<yee; j++) { 769c4762a1bSJed Brown for (i=xes; i<xee; i++) { 770c4762a1bSJed Brown GatherElementData(mx,my,mz,x,c,i,j,k,ex,ec,user); 771c4762a1bSJed Brown FormElementJacobian(ex,ec,ef,NULL,user); 772c4762a1bSJed Brown /* put this element's additions into the residuals */ 773c4762a1bSJed Brown for (kk=0;kk<NB;kk++) { 774c4762a1bSJed Brown for (jj=0;jj<NB;jj++) { 775c4762a1bSJed Brown for (ii=0;ii<NB;ii++) { 776c4762a1bSJed Brown PetscInt idx = ii + jj*NB + kk*NB*NB; 777c4762a1bSJed Brown if (k+kk >= zs && j+jj >= ys && i+ii >= xs && k+kk < zs+zm && j+jj < ys+ym && i+ii < xs+xm) { 778c4762a1bSJed Brown if (OnBoundary(i+ii,j+jj,k+kk,mx,my,mz)) { 779c4762a1bSJed Brown for (l=0;l<3;l++) 780c4762a1bSJed Brown f[k+kk][j+jj][i+ii][l] = x[k+kk][j+jj][i+ii][l] - ex[idx][l]; 781c4762a1bSJed Brown } else { 782c4762a1bSJed Brown for (l=0;l<3;l++) 783c4762a1bSJed Brown f[k+kk][j+jj][i+ii][l] += ef[idx][l]; 784c4762a1bSJed Brown } 785c4762a1bSJed Brown } 786c4762a1bSJed Brown } 787c4762a1bSJed Brown } 788c4762a1bSJed Brown } 789c4762a1bSJed Brown } 790c4762a1bSJed Brown } 791c4762a1bSJed Brown } 7929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda,C,&c)); 793c4762a1bSJed Brown PetscFunctionReturn(0); 794c4762a1bSJed Brown } 795c4762a1bSJed Brown 796c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X,Vec B,void *ptr) 797c4762a1bSJed Brown { 798c4762a1bSJed Brown /* values for each basis function at each quadrature point */ 799c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 800c4762a1bSJed Brown PetscInt i,j,k,l,m,n,s; 801c4762a1bSJed Brown PetscInt pi,pj,pk; 802c4762a1bSJed Brown Field ef[1]; 803c4762a1bSJed Brown Field ex[8]; 804c4762a1bSJed Brown PetscScalar ej[9]; 805c4762a1bSJed Brown CoordField ec[8]; 806c4762a1bSJed Brown PetscScalar pjac[9],pjinv[9]; 807c4762a1bSJed Brown PetscScalar pf[3],py[3]; 808c4762a1bSJed Brown PetscInt xs,ys,zs; 809c4762a1bSJed Brown PetscInt xm,ym,zm; 810c4762a1bSJed Brown PetscInt mx,my,mz; 811c4762a1bSJed Brown DM cda; 812c4762a1bSJed Brown CoordField ***c; 813c4762a1bSJed Brown Vec C; 814c4762a1bSJed Brown DM da; 815c4762a1bSJed Brown Vec Xl,Bl; 816c4762a1bSJed Brown Field ***x,***b; 817c4762a1bSJed Brown PetscInt sweeps,its; 818c4762a1bSJed Brown PetscReal atol,rtol,stol; 819c4762a1bSJed Brown PetscReal fnorm0 = 0.0,fnorm,ynorm,xnorm = 0.0; 820c4762a1bSJed Brown 821c4762a1bSJed Brown PetscFunctionBegin; 8229566063dSJacob Faibussowitsch PetscCall(SNESNGSGetSweeps(snes,&sweeps)); 8239566063dSJacob Faibussowitsch PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its)); 824c4762a1bSJed Brown 8259566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&da)); 8269566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Xl)); 827c4762a1bSJed Brown if (B) { 8289566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&Bl)); 829c4762a1bSJed Brown } 8309566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xl)); 8319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xl)); 832c4762a1bSJed Brown if (B) { 8339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,Bl)); 8349566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,Bl)); 835c4762a1bSJed Brown } 8369566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,Xl,&x)); 8379566063dSJacob Faibussowitsch if (B) PetscCall(DMDAVecGetArray(da,Bl,&b)); 838c4762a1bSJed Brown 8399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da,&cda)); 8409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(da,&C)); 8419566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda,C,&c)); 8429566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 8439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 844c4762a1bSJed Brown 845c4762a1bSJed Brown for (s=0;s<sweeps;s++) { 846c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 847c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 848c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 849c4762a1bSJed Brown if (OnBoundary(i,j,k,mx,my,mz)) { 850c4762a1bSJed Brown BoundaryValue(i,j,k,mx,my,mz,x[k][j][i],user); 851c4762a1bSJed Brown } else { 852c4762a1bSJed Brown for (n=0;n<its;n++) { 853c4762a1bSJed Brown for (m=0;m<9;m++) pjac[m] = 0.; 854c4762a1bSJed Brown for (m=0;m<3;m++) pf[m] = 0.; 855c4762a1bSJed Brown /* gather the elements for this point */ 856c4762a1bSJed Brown for (pk=-1; pk<1; pk++) { 857c4762a1bSJed Brown for (pj=-1; pj<1; pj++) { 858c4762a1bSJed Brown for (pi=-1; pi<1; pi++) { 859c4762a1bSJed Brown /* check that this element exists */ 860c4762a1bSJed Brown if (i+pi >= 0 && i+pi < mx-1 && j+pj >= 0 && j+pj < my-1 && k+pk >= 0 && k+pk < mz-1) { 861c4762a1bSJed Brown /* create the element function and jacobian */ 862c4762a1bSJed Brown GatherElementData(mx,my,mz,x,c,i+pi,j+pj,k+pk,ex,ec,user); 863c4762a1bSJed Brown FormPBJacobian(-pi,-pj,-pk,ex,ec,ef,ej,user); 864c4762a1bSJed Brown /* extract the point named by i,j,k from the whole element jacobian and function */ 865c4762a1bSJed Brown for (l=0;l<3;l++) { 866c4762a1bSJed Brown pf[l] += ef[0][l]; 867c4762a1bSJed Brown for (m=0;m<3;m++) { 868c4762a1bSJed Brown pjac[3*m+l] += ej[3*m+l]; 869c4762a1bSJed Brown } 870c4762a1bSJed Brown } 871c4762a1bSJed Brown } 872c4762a1bSJed Brown } 873c4762a1bSJed Brown } 874c4762a1bSJed Brown } 875c4762a1bSJed Brown /* invert */ 876c4762a1bSJed Brown InvertTensor(pjac,pjinv,NULL); 877c4762a1bSJed Brown /* apply */ 878c4762a1bSJed Brown if (B) for (m=0;m<3;m++) { 879c4762a1bSJed Brown pf[m] -= b[k][j][i][m]; 880c4762a1bSJed Brown } 881c4762a1bSJed Brown TensorVector(pjinv,pf,py); 882c4762a1bSJed Brown xnorm=0.; 883c4762a1bSJed Brown for (m=0;m<3;m++) { 884c4762a1bSJed Brown x[k][j][i][m] -= py[m]; 885c4762a1bSJed Brown xnorm += PetscRealPart(x[k][j][i][m]*x[k][j][i][m]); 886c4762a1bSJed Brown } 887c4762a1bSJed Brown fnorm = PetscRealPart(pf[0]*pf[0]+pf[1]*pf[1]+pf[2]*pf[2]); 888c4762a1bSJed Brown if (n==0) fnorm0 = fnorm; 889c4762a1bSJed Brown ynorm = PetscRealPart(py[0]*py[0]+py[1]*py[1]+py[2]*py[2]); 890c4762a1bSJed Brown if (fnorm < atol*atol || fnorm < rtol*rtol*fnorm0 || ynorm < stol*stol*xnorm) break; 891c4762a1bSJed Brown } 892c4762a1bSJed Brown } 893c4762a1bSJed Brown } 894c4762a1bSJed Brown } 895c4762a1bSJed Brown } 896c4762a1bSJed Brown } 8979566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Xl,&x)); 8989566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da,Xl,INSERT_VALUES,X)); 8999566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da,Xl,INSERT_VALUES,X)); 9009566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Xl)); 901c4762a1bSJed Brown if (B) { 9029566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,Bl,&b)); 9039566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&Bl)); 904c4762a1bSJed Brown } 9059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda,C,&c)); 906c4762a1bSJed Brown PetscFunctionReturn(0); 907c4762a1bSJed Brown } 908c4762a1bSJed Brown 909c4762a1bSJed Brown PetscErrorCode FormCoordinates(DM da,AppCtx *user) 910c4762a1bSJed Brown { 911c4762a1bSJed Brown Vec coords; 912c4762a1bSJed Brown DM cda; 913c4762a1bSJed Brown PetscInt mx,my,mz; 914c4762a1bSJed Brown PetscInt i,j,k,xs,ys,zs,xm,ym,zm; 915c4762a1bSJed Brown CoordField ***x; 916c4762a1bSJed Brown 917c4762a1bSJed Brown PetscFunctionBegin; 9189566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da,&cda)); 9199566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(cda,&coords)); 9209566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 9219566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 9229566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,coords,&x)); 923c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 924c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 925c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 926c4762a1bSJed Brown PetscReal cx = ((PetscReal)i) / (((PetscReal)(mx-1))); 927c4762a1bSJed Brown PetscReal cy = ((PetscReal)j) / (((PetscReal)(my-1))); 928c4762a1bSJed Brown PetscReal cz = ((PetscReal)k) / (((PetscReal)(mz-1))); 929c4762a1bSJed Brown PetscReal rad = user->rad + cy*user->height; 930c4762a1bSJed Brown PetscReal ang = (cx - 0.5)*user->arc; 931c4762a1bSJed Brown x[k][j][i][0] = rad*PetscSinReal(ang); 932c4762a1bSJed Brown x[k][j][i][1] = rad*PetscCosReal(ang) - (user->rad + 0.5*user->height)*PetscCosReal(-0.5*user->arc); 933c4762a1bSJed Brown x[k][j][i][2] = user->width*(cz - 0.5); 934c4762a1bSJed Brown } 935c4762a1bSJed Brown } 936c4762a1bSJed Brown } 9379566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,coords,&x)); 9389566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da,coords)); 9399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 940c4762a1bSJed Brown PetscFunctionReturn(0); 941c4762a1bSJed Brown } 942c4762a1bSJed Brown 943c4762a1bSJed Brown PetscErrorCode InitialGuess(DM da,AppCtx *user,Vec X) 944c4762a1bSJed Brown { 945c4762a1bSJed Brown PetscInt i,j,k,xs,ys,zs,xm,ym,zm; 946c4762a1bSJed Brown PetscInt mx,my,mz; 947c4762a1bSJed Brown Field ***x; 948c4762a1bSJed Brown 949c4762a1bSJed Brown PetscFunctionBegin; 9509566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 9519566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 9529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,X,&x)); 953c4762a1bSJed Brown 954c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 955c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 956c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 957c4762a1bSJed Brown /* reference coordinates */ 958c4762a1bSJed Brown PetscReal p_x = ((PetscReal)i) / (((PetscReal)(mx-1))); 959c4762a1bSJed Brown PetscReal p_y = ((PetscReal)j) / (((PetscReal)(my-1))); 960c4762a1bSJed Brown PetscReal p_z = ((PetscReal)k) / (((PetscReal)(mz-1))); 961c4762a1bSJed Brown PetscReal o_x = p_x; 962c4762a1bSJed Brown PetscReal o_y = p_y; 963c4762a1bSJed Brown PetscReal o_z = p_z; 964c4762a1bSJed Brown x[k][j][i][0] = o_x - p_x; 965c4762a1bSJed Brown x[k][j][i][1] = o_y - p_y; 966c4762a1bSJed Brown x[k][j][i][2] = o_z - p_z; 967c4762a1bSJed Brown } 968c4762a1bSJed Brown } 969c4762a1bSJed Brown } 9709566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,X,&x)); 971c4762a1bSJed Brown PetscFunctionReturn(0); 972c4762a1bSJed Brown } 973c4762a1bSJed Brown 974c4762a1bSJed Brown PetscErrorCode FormRHS(DM da,AppCtx *user,Vec X) 975c4762a1bSJed Brown { 976c4762a1bSJed Brown PetscInt i,j,k,xs,ys,zs,xm,ym,zm; 977c4762a1bSJed Brown PetscInt mx,my,mz; 978c4762a1bSJed Brown Field ***x; 979c4762a1bSJed Brown 980c4762a1bSJed Brown PetscFunctionBegin; 9819566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 9829566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 9839566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,X,&x)); 984c4762a1bSJed Brown 985c4762a1bSJed Brown for (k=zs; k<zs+zm; k++) { 986c4762a1bSJed Brown for (j=ys; j<ys+ym; j++) { 987c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 988c4762a1bSJed Brown x[k][j][i][0] = 0.; 989c4762a1bSJed Brown x[k][j][i][1] = 0.; 990c4762a1bSJed Brown x[k][j][i][2] = 0.; 991c4762a1bSJed Brown if (i == (mx-1)/2 && j == (my-1)) x[k][j][i][1] = user->ploading/(mz-1); 992c4762a1bSJed Brown } 993c4762a1bSJed Brown } 994c4762a1bSJed Brown } 9959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,X,&x)); 996c4762a1bSJed Brown PetscFunctionReturn(0); 997c4762a1bSJed Brown } 998c4762a1bSJed Brown 999c4762a1bSJed Brown PetscErrorCode DisplayLine(SNES snes,Vec X) 1000c4762a1bSJed Brown { 1001c4762a1bSJed Brown PetscInt r,i,j=0,k=0,xs,xm,ys,ym,zs,zm,mx,my,mz; 1002c4762a1bSJed Brown Field ***x; 1003c4762a1bSJed Brown CoordField ***c; 1004c4762a1bSJed Brown DM da,cda; 1005c4762a1bSJed Brown Vec C; 1006c4762a1bSJed Brown PetscMPIInt size,rank; 1007c4762a1bSJed Brown 1008c4762a1bSJed Brown PetscFunctionBegin; 10099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 10109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 10119566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&da)); 10129566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&mx,&my,&mz,0,0,0,0,0,0,0,0,0)); 10139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da,&cda)); 10149566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&C)); 1015c4762a1bSJed Brown j = my / 2; 1016c4762a1bSJed Brown k = mz / 2; 10179566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 10189566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,X,&x)); 10199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda,C,&c)); 1020c4762a1bSJed Brown for (r=0;r<size;r++) { 1021c4762a1bSJed Brown if (rank == r) { 1022c4762a1bSJed Brown if (j >= ys && j < ys+ym && k >= zs && k < zs+zm) { 1023c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 102463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %d %d: %f %f %f\n",i,0,0,(double)PetscRealPart(c[k][j][i][0] + x[k][j][i][0]),(double)PetscRealPart(c[k][j][i][1] + x[k][j][i][1]),(double)PetscRealPart(c[k][j][i][2] + x[k][j][i][2]))); 1025c4762a1bSJed Brown } 1026c4762a1bSJed Brown } 1027c4762a1bSJed Brown } 10289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 1029c4762a1bSJed Brown } 10309566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,X,&x)); 10319566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda,C,&c)); 1032c4762a1bSJed Brown PetscFunctionReturn(0); 1033c4762a1bSJed Brown } 1034c4762a1bSJed Brown 1035c4762a1bSJed Brown /*TEST 1036c4762a1bSJed Brown 1037c4762a1bSJed Brown test: 1038c4762a1bSJed Brown nsize: 2 1039c4762a1bSJed Brown args: -da_refine 2 -pc_type mg -rad 10.0 -young 10. -ploading 0.0 -loading -1. -mg_levels_ksp_max_it 2 -snes_monitor_short -ksp_monitor_short -snes_max_it 7 1040c4762a1bSJed Brown requires: !single 1041c4762a1bSJed Brown timeoutfactor: 3 1042c4762a1bSJed Brown 1043c4762a1bSJed Brown test: 1044c4762a1bSJed Brown suffix: 2 1045c4762a1bSJed Brown args: -da_refine 2 -pc_type mg -rad 10.0 -young 10. -ploading 0.0 -loading -1. -mg_levels_ksp_max_it 2 -snes_monitor_short -ksp_monitor_short -npc_snes_type fas -npc_fas_levels_snes_type ncg -npc_fas_levels_snes_max_it 3 -npc_snes_monitor_short -snes_max_it 2 1046c4762a1bSJed Brown requires: !single 1047c4762a1bSJed Brown 1048c4762a1bSJed Brown test: 1049c4762a1bSJed Brown suffix: 3 1050c4762a1bSJed 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 1051c4762a1bSJed Brown requires: !single 1052c4762a1bSJed Brown 1053c4762a1bSJed Brown TEST*/ 1054