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