Actual source code: ex1.c

petsc-3.3-p7 2013-05-11
  2: static char help[] ="Solves the time dependent Bratu problem using pseudo-timestepping.";

  4: /*
  5:    Concepts: TS^pseudo-timestepping
  6:    Concepts: pseudo-timestepping
  7:    Concepts: nonlinear problems
  8:    Processors: 1

 10: */

 12: /* ------------------------------------------------------------------------

 14:     This code demonstrates how one may solve a nonlinear problem 
 15:     with pseudo-timestepping. In this simple example, the pseudo-timestep
 16:     is the same for all grid points, i.e., this is equivalent to using
 17:     the backward Euler method with a variable timestep.

 19:     Note: This example does not require pseudo-timestepping since it
 20:     is an easy nonlinear problem, but it is included to demonstrate how
 21:     the pseudo-timestepping may be done.

 23:     See snes/examples/tutorials/ex4.c[ex4f.F] and 
 24:     snes/examples/tutorials/ex5.c[ex5f.F] where the problem is described
 25:     and solved using Newton's method alone.

 27:   ----------------------------------------------------------------------------- */
 28: /*
 29:     Include "petscts.h" to use the PETSc timestepping routines. Note that
 30:     this file automatically includes "petscsys.h" and other lower-level
 31:     PETSc include files.
 32: */
 33: #include <petscts.h>

 35: /*
 36:   Create an application context to contain data needed by the 
 37:   application-provided call-back routines, FormJacobian() and
 38:   FormFunction().
 39: */
 40: typedef struct {
 41:   PetscReal   param;        /* test problem parameter */
 42:   PetscInt    mx;           /* Discretization in x-direction */
 43:   PetscInt    my;           /* Discretization in y-direction */
 44: } AppCtx;

 46: /* 
 47:    User-defined routines
 48: */
 49: extern PetscErrorCode  FormJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*),
 50:      FormFunction(TS,PetscReal,Vec,Vec,void*),
 51:      FormInitialGuess(Vec,AppCtx*);

 55: int main(int argc,char **argv)
 56: {
 57:   TS             ts;                 /* timestepping context */
 58:   Vec            x,r;               /* solution, residual vectors */
 59:   Mat            J;                  /* Jacobian matrix */
 60:   AppCtx         user;               /* user-defined work context */
 61:   PetscInt       its,N;                /* iterations for convergence */
 63:   PetscReal      param_max = 6.81,param_min = 0.,dt;
 64:   PetscReal      ftime;
 65:   PetscMPIInt    size;

 67:   PetscInitialize(&argc,&argv,PETSC_NULL,help);
 68:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 69:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only");

 71:   user.mx        = 4;
 72:   user.my        = 4;
 73:   user.param     = 6.0;
 74: 
 75:   /*
 76:      Allow user to set the grid dimensions and nonlinearity parameter at run-time
 77:   */
 78:   PetscOptionsGetInt(PETSC_NULL,"-mx",&user.mx,PETSC_NULL);
 79:   PetscOptionsGetInt(PETSC_NULL,"-my",&user.my,PETSC_NULL);
 80:   PetscOptionsGetReal(PETSC_NULL,"-param",&user.param,PETSC_NULL);
 81:   if (user.param >= param_max || user.param <= param_min) SETERRQ(PETSC_COMM_SELF,1,"Parameter is out of range");
 82:   dt = .5/PetscMax(user.mx,user.my);
 83:   PetscOptionsGetReal(PETSC_NULL,"-dt",&dt,PETSC_NULL);
 84:   N          = user.mx*user.my;
 85: 
 86:   /* 
 87:       Create vectors to hold the solution and function value
 88:   */
 89:   VecCreateSeq(PETSC_COMM_SELF,N,&x);
 90:   VecDuplicate(x,&r);

 92:   /*
 93:     Create matrix to hold Jacobian. Preallocate 5 nonzeros per row
 94:     in the sparse matrix. Note that this is not the optimal strategy; see
 95:     the Performance chapter of the users manual for information on 
 96:     preallocating memory in sparse matrices.
 97:   */
 98:   MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,5,0,&J);

100:   /* 
101:      Create timestepper context 
102:   */
103:   TSCreate(PETSC_COMM_WORLD,&ts);
104:   TSSetProblemType(ts,TS_NONLINEAR);

106:   /*
107:      Tell the timestepper context where to compute solutions
108:   */
109:   TSSetSolution(ts,x);

111:   /*
112:      Provide the call-back for the nonlinear function we are 
113:      evaluating. Thus whenever the timestepping routines need the
114:      function they will call this routine. Note the final argument
115:      is the application context used by the call-back functions.
116:   */
117:   TSSetRHSFunction(ts,PETSC_NULL,FormFunction,&user);

119:   /*
120:      Set the Jacobian matrix and the function used to compute 
121:      Jacobians.
122:   */
123:   TSSetRHSJacobian(ts,J,J,FormJacobian,&user);

125:   /*
126:        For the initial guess for the problem
127:   */
128:   FormInitialGuess(x,&user);

130:   /*
131:        This indicates that we are using pseudo timestepping to 
132:      find a steady state solution to the nonlinear problem.
133:   */
134:   TSSetType(ts,TSPSEUDO);

136:   /*
137:        Set the initial time to start at (this is arbitrary for 
138:      steady state problems; and the initial timestep given above
139:   */
140:   TSSetInitialTimeStep(ts,0.0,dt);

142:   /*
143:       Set a large number of timesteps and final duration time
144:      to insure convergence to steady state.
145:   */
146:   TSSetDuration(ts,1000,1.e12);

148:   /*
149:       Use the default strategy for increasing the timestep
150:   */
151:   TSPseudoSetTimeStep(ts,TSPseudoDefaultTimeStep,0);

153:   /*
154:       Set any additional options from the options database. This
155:      includes all options for the nonlinear and linear solvers used
156:      internally the the timestepping routines.
157:   */
158:   TSSetFromOptions(ts);

160:   TSSetUp(ts);

162:   /*
163:       Perform the solve. This is where the timestepping takes place.
164:   */
165:   TSSolve(ts,x,&ftime);

167:   /*
168:       Get the number of steps
169:   */
170:   TSGetTimeStepNumber(ts,&its);

172:   printf("Number of pseudo timesteps = %d final time %4.2e\n",(int)its,ftime);

174:   /* 
175:      Free the data structures constructed above
176:   */
177:   VecDestroy(&x);
178:   VecDestroy(&r);
179:   MatDestroy(&J);
180:   TSDestroy(&ts);
181:   PetscFinalize();

183:   return 0;
184: }
185: /* ------------------------------------------------------------------ */
186: /*           Bratu (Solid Fuel Ignition) Test Problem                 */
187: /* ------------------------------------------------------------------ */

189: /* --------------------  Form initial approximation ----------------- */

193: PetscErrorCode FormInitialGuess(Vec X,AppCtx *user)
194: {
195:   PetscInt       i,j,row,mx,my;
197:   PetscReal      one = 1.0,lambda;
198:   PetscReal      temp1,temp,hx,hy;
199:   PetscScalar    *x;

201:   mx         = user->mx;
202:   my         = user->my;
203:   lambda = user->param;

205:   hx    = one / (PetscReal)(mx-1);
206:   hy    = one / (PetscReal)(my-1);

208:   VecGetArray(X,&x);
209:   temp1 = lambda/(lambda + one);
210:   for (j=0; j<my; j++) {
211:     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
212:     for (i=0; i<mx; i++) {
213:       row = i + j*mx;
214:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
215:         x[row] = 0.0;
216:         continue;
217:       }
218:       x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
219:     }
220:   }
221:   VecRestoreArray(X,&x);
222:   return 0;
223: }
224: /* --------------------  Evaluate Function F(x) --------------------- */

228: PetscErrorCode FormFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
229: {
230:   AppCtx         *user = (AppCtx*)ptr;
232:   PetscInt       i,j,row,mx,my;
233:   PetscReal      two = 2.0,one = 1.0,lambda;
234:   PetscReal      hx,hy,hxdhy,hydhx;
235:   PetscScalar    ut,ub,ul,ur,u,uxx,uyy,sc,*x,*f;

237:   mx         = user->mx;
238:   my         = user->my;
239:   lambda = user->param;

241:   hx    = one / (PetscReal)(mx-1);
242:   hy    = one / (PetscReal)(my-1);
243:   sc    = hx*hy;
244:   hxdhy = hx/hy;
245:   hydhx = hy/hx;

247:   VecGetArray(X,&x);
248:   VecGetArray(F,&f);
249:   for (j=0; j<my; j++) {
250:     for (i=0; i<mx; i++) {
251:       row = i + j*mx;
252:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
253:         f[row] = x[row];
254:         continue;
255:       }
256:       u = x[row];
257:       ub = x[row - mx];
258:       ul = x[row - 1];
259:       ut = x[row + mx];
260:       ur = x[row + 1];
261:       uxx = (-ur + two*u - ul)*hydhx;
262:       uyy = (-ut + two*u - ub)*hxdhy;
263:       f[row] = -uxx + -uyy + sc*lambda*PetscExpScalar(u);
264:     }
265:   }
266:   VecRestoreArray(X,&x);
267:   VecRestoreArray(F,&f);
268:   return 0;
269: }
270: /* --------------------  Evaluate Jacobian F'(x) -------------------- */

274: /*
275:    Calculate the Jacobian matrix J(X,t).

277:    Note: We put the Jacobian in the preconditioner storage B instead of J. This
278:    way we can give the -snes_mf_operator flag to check our work. This replaces
279:    J with a finite difference approximation, using our analytic Jacobian B for
280:    the preconditioner.
281: */
282: PetscErrorCode FormJacobian(TS ts,PetscReal t,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
283: {
284:   AppCtx         *user = (AppCtx*)ptr;
285:   Mat            jac = *B;
286:   PetscInt       i,j,row,mx,my,col[5];
288:   PetscScalar    two = 2.0,one = 1.0,lambda,v[5],sc,*x;
289:   PetscReal      hx,hy,hxdhy,hydhx;


292:   mx         = user->mx;
293:   my         = user->my;
294:   lambda = user->param;

296:   hx    = 1.0 / (PetscReal)(mx-1);
297:   hy    = 1.0 / (PetscReal)(my-1);
298:   sc    = hx*hy;
299:   hxdhy = hx/hy;
300:   hydhx = hy/hx;

302:   VecGetArray(X,&x);
303:   for (j=0; j<my; j++) {
304:     for (i=0; i<mx; i++) {
305:       row = i + j*mx;
306:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
307:         MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);
308:         continue;
309:       }
310:       v[0] = hxdhy; col[0] = row - mx;
311:       v[1] = hydhx; col[1] = row - 1;
312:       v[2] = -two*(hydhx + hxdhy) + sc*lambda*PetscExpScalar(x[row]); col[2] = row;
313:       v[3] = hydhx; col[3] = row + 1;
314:       v[4] = hxdhy; col[4] = row + mx;
315:       MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);
316:     }
317:   }
318:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
319:   VecRestoreArray(X,&x);
320:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
321:   *flag = SAME_NONZERO_PATTERN;
322:   return 0;
323: }