1*55a74a43SLisandro Dalcin /* ------------------------------------------------------------------------ 2*55a74a43SLisandro Dalcin 3*55a74a43SLisandro Dalcin Solid Fuel Ignition (SFI) problem. This problem is modeled by the 4*55a74a43SLisandro Dalcin partial differential equation 5*55a74a43SLisandro Dalcin 6*55a74a43SLisandro Dalcin -Laplacian(u) - lambda * exp(u) = 0, 0 < x,y,z < 1, 7*55a74a43SLisandro Dalcin 8*55a74a43SLisandro Dalcin with boundary conditions 9*55a74a43SLisandro Dalcin 10*55a74a43SLisandro Dalcin u = 0 for x = 0, x = 1, y = 0, y = 1, z = 0, z = 1 11*55a74a43SLisandro Dalcin 12*55a74a43SLisandro Dalcin A finite difference approximation with the usual 7-point stencil 13*55a74a43SLisandro Dalcin is used to discretize the boundary value problem to obtain a 14*55a74a43SLisandro Dalcin nonlinear system of equations. The problem is solved in a 3D 15*55a74a43SLisandro Dalcin rectangular domain, using distributed arrays (DAs) to partition 16*55a74a43SLisandro Dalcin the parallel grid. 17*55a74a43SLisandro Dalcin 18*55a74a43SLisandro Dalcin ------------------------------------------------------------------------- */ 19*55a74a43SLisandro Dalcin 20*55a74a43SLisandro Dalcin #include "Bratu3D.h" 21*55a74a43SLisandro Dalcin 22*55a74a43SLisandro Dalcin PetscErrorCode FormInitGuess(DM da, Vec X, Params *p) 23*55a74a43SLisandro Dalcin { 24*55a74a43SLisandro Dalcin PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm; 25*55a74a43SLisandro Dalcin PetscReal lambda,temp1,hx,hy,hz,tempk,tempj; 26*55a74a43SLisandro Dalcin PetscScalar ***x; 27*55a74a43SLisandro Dalcin 28*55a74a43SLisandro Dalcin PetscFunctionBegin; 29*55a74a43SLisandro Dalcin PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 30*55a74a43SLisandro Dalcin lambda = p->lambda_; 31*55a74a43SLisandro Dalcin hx = 1.0/(PetscReal)(Mx-1); 32*55a74a43SLisandro Dalcin hy = 1.0/(PetscReal)(My-1); 33*55a74a43SLisandro Dalcin hz = 1.0/(PetscReal)(Mz-1); 34*55a74a43SLisandro Dalcin temp1 = lambda/(lambda + 1.0); 35*55a74a43SLisandro Dalcin 36*55a74a43SLisandro Dalcin /* 37*55a74a43SLisandro Dalcin Get a pointer to vector data. 38*55a74a43SLisandro Dalcin 39*55a74a43SLisandro Dalcin - For default PETSc vectors, VecGetArray() returns a pointer to 40*55a74a43SLisandro Dalcin the data array. Otherwise, the routine is implementation 41*55a74a43SLisandro Dalcin dependent. 42*55a74a43SLisandro Dalcin 43*55a74a43SLisandro Dalcin - You MUST call VecRestoreArray() when you no longer need access 44*55a74a43SLisandro Dalcin to the array. 45*55a74a43SLisandro Dalcin */ 46*55a74a43SLisandro Dalcin PetscCall(DMDAVecGetArray(da,X,&x)); 47*55a74a43SLisandro Dalcin 48*55a74a43SLisandro Dalcin /* 49*55a74a43SLisandro Dalcin Get local grid boundaries (for 3-dimensional DMDA): 50*55a74a43SLisandro Dalcin 51*55a74a43SLisandro Dalcin - xs, ys, zs: starting grid indices (no ghost points) 52*55a74a43SLisandro Dalcin 53*55a74a43SLisandro Dalcin - xm, ym, zm: widths of local grid (no ghost points) 54*55a74a43SLisandro Dalcin */ 55*55a74a43SLisandro Dalcin PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 56*55a74a43SLisandro Dalcin 57*55a74a43SLisandro Dalcin /* 58*55a74a43SLisandro Dalcin Compute initial guess over the locally owned part of the grid 59*55a74a43SLisandro Dalcin */ 60*55a74a43SLisandro Dalcin for (k=zs; k<zs+zm; k++) { 61*55a74a43SLisandro Dalcin tempk = (PetscReal)(PetscMin(k,Mz-k-1))*hz; 62*55a74a43SLisandro Dalcin for (j=ys; j<ys+ym; j++) { 63*55a74a43SLisandro Dalcin tempj = PetscMin((PetscReal)(PetscMin(j,My-j-1))*hy,tempk); 64*55a74a43SLisandro Dalcin for (i=xs; i<xs+xm; i++) { 65*55a74a43SLisandro Dalcin if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) { 66*55a74a43SLisandro Dalcin /* boundary conditions are all zero Dirichlet */ 67*55a74a43SLisandro Dalcin x[k][j][i] = 0.0; 68*55a74a43SLisandro Dalcin } else { 69*55a74a43SLisandro Dalcin x[k][j][i] = temp1*sqrt(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,tempj)); 70*55a74a43SLisandro Dalcin } 71*55a74a43SLisandro Dalcin } 72*55a74a43SLisandro Dalcin } 73*55a74a43SLisandro Dalcin } 74*55a74a43SLisandro Dalcin 75*55a74a43SLisandro Dalcin /* 76*55a74a43SLisandro Dalcin Restore vector 77*55a74a43SLisandro Dalcin */ 78*55a74a43SLisandro Dalcin PetscCall(DMDAVecRestoreArray(da,X,&x)); 79*55a74a43SLisandro Dalcin PetscFunctionReturn(PETSC_SUCCESS); 80*55a74a43SLisandro Dalcin } 81*55a74a43SLisandro Dalcin 82*55a74a43SLisandro Dalcin PetscErrorCode FormFunction(DM da, Vec X, Vec F, Params *p) 83*55a74a43SLisandro Dalcin { 84*55a74a43SLisandro Dalcin PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm; 85*55a74a43SLisandro Dalcin PetscReal two = 2.0,lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc; 86*55a74a43SLisandro Dalcin PetscScalar u_north,u_south,u_east,u_west,u_up,u_down,u,u_xx,u_yy,u_zz,***x,***f; 87*55a74a43SLisandro Dalcin Vec localX; 88*55a74a43SLisandro Dalcin 89*55a74a43SLisandro Dalcin PetscFunctionBegin; 90*55a74a43SLisandro Dalcin PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 91*55a74a43SLisandro Dalcin lambda = p->lambda_; 92*55a74a43SLisandro Dalcin hx = 1.0/(PetscReal)(Mx-1); 93*55a74a43SLisandro Dalcin hy = 1.0/(PetscReal)(My-1); 94*55a74a43SLisandro Dalcin hz = 1.0/(PetscReal)(Mz-1); 95*55a74a43SLisandro Dalcin sc = hx*hy*hz*lambda; 96*55a74a43SLisandro Dalcin hxhzdhy = hx*hz/hy; 97*55a74a43SLisandro Dalcin hyhzdhx = hy*hz/hx; 98*55a74a43SLisandro Dalcin hxhydhz = hx*hy/hz; 99*55a74a43SLisandro Dalcin 100*55a74a43SLisandro Dalcin PetscCall(DMGetLocalVector(da,&localX)); 101*55a74a43SLisandro Dalcin 102*55a74a43SLisandro Dalcin /* 103*55a74a43SLisandro Dalcin Scatter ghost points to local vector,using the 2-step process 104*55a74a43SLisandro Dalcin DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). By placing code 105*55a74a43SLisandro Dalcin between these two statements, computations can be done while 106*55a74a43SLisandro Dalcin messages are in transition. 107*55a74a43SLisandro Dalcin */ 108*55a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 109*55a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 110*55a74a43SLisandro Dalcin 111*55a74a43SLisandro Dalcin /* 112*55a74a43SLisandro Dalcin Get pointers to vector data. 113*55a74a43SLisandro Dalcin */ 114*55a74a43SLisandro Dalcin PetscCall(DMDAVecGetArray(da,localX,&x)); 115*55a74a43SLisandro Dalcin PetscCall(DMDAVecGetArray(da,F,&f)); 116*55a74a43SLisandro Dalcin 117*55a74a43SLisandro Dalcin /* 118*55a74a43SLisandro Dalcin Get local grid boundaries. 119*55a74a43SLisandro Dalcin */ 120*55a74a43SLisandro Dalcin PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 121*55a74a43SLisandro Dalcin 122*55a74a43SLisandro Dalcin /* 123*55a74a43SLisandro Dalcin Compute function over the locally owned part of the grid. 124*55a74a43SLisandro Dalcin */ 125*55a74a43SLisandro Dalcin for (k=zs; k<zs+zm; k++) { 126*55a74a43SLisandro Dalcin for (j=ys; j<ys+ym; j++) { 127*55a74a43SLisandro Dalcin for (i=xs; i<xs+xm; i++) { 128*55a74a43SLisandro Dalcin if (i == 0 || j == 0 || k == 0 || i == Mx-1 || j == My-1 || k == Mz-1) { 129*55a74a43SLisandro Dalcin /* boundary points */ 130*55a74a43SLisandro Dalcin f[k][j][i] = x[k][j][i] - 0.0; 131*55a74a43SLisandro Dalcin } else { 132*55a74a43SLisandro Dalcin /* interior grid points */ 133*55a74a43SLisandro Dalcin u = x[k][j][i]; 134*55a74a43SLisandro Dalcin u_east = x[k][j][i+1]; 135*55a74a43SLisandro Dalcin u_west = x[k][j][i-1]; 136*55a74a43SLisandro Dalcin u_north = x[k][j+1][i]; 137*55a74a43SLisandro Dalcin u_south = x[k][j-1][i]; 138*55a74a43SLisandro Dalcin u_up = x[k+1][j][i]; 139*55a74a43SLisandro Dalcin u_down = x[k-1][j][i]; 140*55a74a43SLisandro Dalcin u_xx = (-u_east + two*u - u_west)*hyhzdhx; 141*55a74a43SLisandro Dalcin u_yy = (-u_north + two*u - u_south)*hxhzdhy; 142*55a74a43SLisandro Dalcin u_zz = (-u_up + two*u - u_down)*hxhydhz; 143*55a74a43SLisandro Dalcin f[k][j][i] = u_xx + u_yy + u_zz - sc*PetscExpScalar(u); 144*55a74a43SLisandro Dalcin } 145*55a74a43SLisandro Dalcin } 146*55a74a43SLisandro Dalcin } 147*55a74a43SLisandro Dalcin } 148*55a74a43SLisandro Dalcin 149*55a74a43SLisandro Dalcin /* 150*55a74a43SLisandro Dalcin Restore vectors. 151*55a74a43SLisandro Dalcin */ 152*55a74a43SLisandro Dalcin PetscCall(DMDAVecRestoreArray(da,F,&f)); 153*55a74a43SLisandro Dalcin PetscCall(DMDAVecRestoreArray(da,localX,&x)); 154*55a74a43SLisandro Dalcin PetscCall(DMRestoreLocalVector(da,&localX)); 155*55a74a43SLisandro Dalcin PetscCall(PetscLogFlops(11.0*ym*xm)); 156*55a74a43SLisandro Dalcin PetscFunctionReturn(PETSC_SUCCESS); 157*55a74a43SLisandro Dalcin } 158*55a74a43SLisandro Dalcin 159*55a74a43SLisandro Dalcin PetscErrorCode FormJacobian(DM da, Vec X, Mat J, Params *p) 160*55a74a43SLisandro Dalcin { 161*55a74a43SLisandro Dalcin PetscInt i,j,k,Mx,My,Mz,xs,ys,zs,xm,ym,zm; 162*55a74a43SLisandro Dalcin PetscReal lambda,hx,hy,hz,hxhzdhy,hyhzdhx,hxhydhz,sc; 163*55a74a43SLisandro Dalcin PetscScalar v[7],***x; 164*55a74a43SLisandro Dalcin MatStencil col[7],row; 165*55a74a43SLisandro Dalcin Vec localX; 166*55a74a43SLisandro Dalcin 167*55a74a43SLisandro Dalcin PetscFunctionBegin; 168*55a74a43SLisandro Dalcin PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,&Mz,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 169*55a74a43SLisandro Dalcin lambda = p->lambda_; 170*55a74a43SLisandro Dalcin hx = 1.0/(PetscReal)(Mx-1); 171*55a74a43SLisandro Dalcin hy = 1.0/(PetscReal)(My-1); 172*55a74a43SLisandro Dalcin hz = 1.0/(PetscReal)(Mz-1); 173*55a74a43SLisandro Dalcin sc = hx*hy*hz*lambda; 174*55a74a43SLisandro Dalcin hxhzdhy = hx*hz/hy; 175*55a74a43SLisandro Dalcin hyhzdhx = hy*hz/hx; 176*55a74a43SLisandro Dalcin hxhydhz = hx*hy/hz; 177*55a74a43SLisandro Dalcin 178*55a74a43SLisandro Dalcin PetscCall(DMGetLocalVector(da,&localX)); 179*55a74a43SLisandro Dalcin 180*55a74a43SLisandro Dalcin /* 181*55a74a43SLisandro Dalcin Scatter ghost points to local vector, using the 2-step process 182*55a74a43SLisandro Dalcin DMGlobalToLocalBegin(), DMGlobalToLocalEnd(). By placing code 183*55a74a43SLisandro Dalcin between these two statements, computations can be done while 184*55a74a43SLisandro Dalcin messages are in transition. 185*55a74a43SLisandro Dalcin */ 186*55a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 187*55a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 188*55a74a43SLisandro Dalcin 189*55a74a43SLisandro Dalcin /* 190*55a74a43SLisandro Dalcin Get pointer to vector data. 191*55a74a43SLisandro Dalcin */ 192*55a74a43SLisandro Dalcin PetscCall(DMDAVecGetArray(da,localX,&x)); 193*55a74a43SLisandro Dalcin 194*55a74a43SLisandro Dalcin /* 195*55a74a43SLisandro Dalcin Get local grid boundaries. 196*55a74a43SLisandro Dalcin */ 197*55a74a43SLisandro Dalcin PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 198*55a74a43SLisandro Dalcin 199*55a74a43SLisandro Dalcin /* 200*55a74a43SLisandro Dalcin Compute entries for the locally owned part of the Jacobian. 201*55a74a43SLisandro Dalcin 202*55a74a43SLisandro Dalcin - Currently, all PETSc parallel matrix formats are partitioned by 203*55a74a43SLisandro Dalcin contiguous chunks of rows across the processors. 204*55a74a43SLisandro Dalcin 205*55a74a43SLisandro Dalcin - Each processor needs to insert only elements that it owns 206*55a74a43SLisandro Dalcin locally (but any non-local elements will be sent to the 207*55a74a43SLisandro Dalcin appropriate processor during matrix assembly). 208*55a74a43SLisandro Dalcin 209*55a74a43SLisandro Dalcin - Here, we set all entries for a particular row at once. 210*55a74a43SLisandro Dalcin 211*55a74a43SLisandro Dalcin - We can set matrix entries either using either 212*55a74a43SLisandro Dalcin MatSetValuesLocal() or MatSetValues(), as discussed above. 213*55a74a43SLisandro Dalcin */ 214*55a74a43SLisandro Dalcin for (k=zs; k<zs+zm; k++) { 215*55a74a43SLisandro Dalcin for (j=ys; j<ys+ym; j++) { 216*55a74a43SLisandro Dalcin for (i=xs; i<xs+xm; i++) { 217*55a74a43SLisandro Dalcin row.k = k; row.j = j; row.i = i; 218*55a74a43SLisandro Dalcin /* boundary points */ 219*55a74a43SLisandro Dalcin if (i == 0 || j == 0 || k == 0|| i == Mx-1 || j == My-1 || k == Mz-1) { 220*55a74a43SLisandro Dalcin v[0] = 1.0; 221*55a74a43SLisandro Dalcin PetscCall(MatSetValuesStencil(J,1,&row,1,&row,v,INSERT_VALUES)); 222*55a74a43SLisandro Dalcin } else { 223*55a74a43SLisandro Dalcin /* interior grid points */ 224*55a74a43SLisandro Dalcin v[0] = -hxhydhz; col[0].k=k-1;col[0].j=j; col[0].i = i; 225*55a74a43SLisandro Dalcin v[1] = -hxhzdhy; col[1].k=k; col[1].j=j-1;col[1].i = i; 226*55a74a43SLisandro Dalcin v[2] = -hyhzdhx; col[2].k=k; col[2].j=j; col[2].i = i-1; 227*55a74a43SLisandro Dalcin v[3] = 2.0*(hyhzdhx+hxhzdhy+hxhydhz)-sc*PetscExpScalar(x[k][j][i]);col[3].k=row.k;col[3].j=row.j;col[3].i = row.i; 228*55a74a43SLisandro Dalcin v[4] = -hyhzdhx; col[4].k=k; col[4].j=j; col[4].i = i+1; 229*55a74a43SLisandro Dalcin v[5] = -hxhzdhy; col[5].k=k; col[5].j=j+1;col[5].i = i; 230*55a74a43SLisandro Dalcin v[6] = -hxhydhz; col[6].k=k+1;col[6].j=j; col[6].i = i; 231*55a74a43SLisandro Dalcin PetscCall(MatSetValuesStencil(J,1,&row,7,col,v,INSERT_VALUES)); 232*55a74a43SLisandro Dalcin } 233*55a74a43SLisandro Dalcin } 234*55a74a43SLisandro Dalcin } 235*55a74a43SLisandro Dalcin } 236*55a74a43SLisandro Dalcin PetscCall(DMDAVecRestoreArray(da,localX,&x)); 237*55a74a43SLisandro Dalcin PetscCall(DMRestoreLocalVector(da,&localX)); 238*55a74a43SLisandro Dalcin 239*55a74a43SLisandro Dalcin /* 240*55a74a43SLisandro Dalcin Assemble matrix, using the 2-step process: MatAssemblyBegin(), 241*55a74a43SLisandro Dalcin MatAssemblyEnd(). 242*55a74a43SLisandro Dalcin */ 243*55a74a43SLisandro Dalcin PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 244*55a74a43SLisandro Dalcin PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 245*55a74a43SLisandro Dalcin 246*55a74a43SLisandro Dalcin /* 247*55a74a43SLisandro Dalcin Tell the matrix we will never add a new nonzero location to the 248*55a74a43SLisandro Dalcin matrix. If we do, it will generate an error. 249*55a74a43SLisandro Dalcin */ 250*55a74a43SLisandro Dalcin PetscCall(MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 251*55a74a43SLisandro Dalcin PetscFunctionReturn(PETSC_SUCCESS); 252*55a74a43SLisandro Dalcin } 253