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