Actual source code: ex1.c

  1: /*$Id: ex1.c,v 1.28 2001/04/10 19:37:12 bsmith Exp $*/

  3: static char help[] ="Solves the time dependent Bratu problem using pseudo-timestepping.";

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

 11: */

 13: /* ------------------------------------------------------------------------

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

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

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

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

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

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

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

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

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

 97:   /* 
 98:      Create timestepper context 
 99:   */
100:   TSCreate(PETSC_COMM_WORLD,TS_NONLINEAR,&ts);

102:   /*
103:      Tell the timestepper context where to compute solutions
104:   */
105:   TSSetSolution(ts,x);

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

115:   /*
116:      Set the Jacobian matrix and the function used to compute 
117:      Jacobians.
118:   */
119:   TSSetRHSJacobian(ts,J,J,FormJacobian,&user);

121:   /*
122:        For the initial guess for the problem
123:   */
124:   FormInitialGuess(x,&user);

126:   /*
127:        This indicates that we are using pseudo timestepping to 
128:      find a steady state solution to the nonlinear problem.
129:   */
130:   TSSetType(ts,TS_PSEUDO);

132:   /*
133:        Set the initial time to start at (this is arbitrary for 
134:      steady state problems; and the initial timestep given above
135:   */
136:   TSSetInitialTimeStep(ts,0.0,dt);

138:   /*
139:       Set a large number of timesteps and final duration time
140:      to insure convergence to steady state.
141:   */
142:   TSSetDuration(ts,1000,1.e12);

144:   /*
145:       Use the default strategy for increasing the timestep
146:   */
147:   TSPseudoSetTimeStep(ts,TSPseudoDefaultTimeStep,0);

149:   /*
150:       Set any additional options from the options database. This
151:      includes all options for the nonlinear and linear solvers used
152:      internally the the timestepping routines.
153:   */
154:   TSSetFromOptions(ts);

156:   TSSetUp(ts);

158:   /*
159:       Perform the solve. This is where the timestepping takes place.
160:   */
161:   TSStep(ts,&its,&ftime);
162: 
163:   printf("Number of pseudo timesteps = %d final time %4.2en",its,ftime);

165:   /* 
166:      Free the data structures constructed above
167:   */
168:   VecDestroy(x);
169:   VecDestroy(r);
170:   MatDestroy(J);
171:   TSDestroy(ts);
172:   PetscFinalize();

174:   return 0;
175: }
176: /* ------------------------------------------------------------------ */
177: /*           Bratu (Solid Fuel Ignition) Test Problem                 */
178: /* ------------------------------------------------------------------ */

180: /* --------------------  Form initial approximation ----------------- */

182: int FormInitialGuess(Vec X,AppCtx *user)
183: {
184:   int     i,j,row,mx,my,ierr;
185:   double  one = 1.0,lambda;
186:   double  temp1,temp,hx,hy;
187:   Scalar  *x;

189:   mx         = user->mx;
190:   my         = user->my;
191:   lambda = user->param;

193:   hx    = one / (double)(mx-1);
194:   hy    = one / (double)(my-1);

196:   VecGetArray(X,&x);
197:   temp1 = lambda/(lambda + one);
198:   for (j=0; j<my; j++) {
199:     temp = (double)(PetscMin(j,my-j-1))*hy;
200:     for (i=0; i<mx; i++) {
201:       row = i + j*mx;
202:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
203:         x[row] = 0.0;
204:         continue;
205:       }
206:       x[row] = temp1*sqrt(PetscMin((double)(PetscMin(i,mx-i-1))*hx,temp));
207:     }
208:   }
209:   VecRestoreArray(X,&x);
210:   return 0;
211: }
212: /* --------------------  Evaluate Function F(x) --------------------- */

214: int FormFunction(TS ts,double t,Vec X,Vec F,void *ptr)
215: {
216:   AppCtx *user = (AppCtx*)ptr;
217:   int     ierr,i,j,row,mx,my;
218:   double  two = 2.0,one = 1.0,lambda;
219:   double  hx,hy,hxdhy,hydhx;
220:   Scalar  ut,ub,ul,ur,u,uxx,uyy,sc,*x,*f;

222:   mx         = user->mx;
223:   my         = user->my;
224:   lambda = user->param;

226:   hx    = one / (double)(mx-1);
227:   hy    = one / (double)(my-1);
228:   sc    = hx*hy;
229:   hxdhy = hx/hy;
230:   hydhx = hy/hx;

232:   VecGetArray(X,&x);
233:   VecGetArray(F,&f);
234:   for (j=0; j<my; j++) {
235:     for (i=0; i<mx; i++) {
236:       row = i + j*mx;
237:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
238:         f[row] = x[row];
239:         continue;
240:       }
241:       u = x[row];
242:       ub = x[row - mx];
243:       ul = x[row - 1];
244:       ut = x[row + mx];
245:       ur = x[row + 1];
246:       uxx = (-ur + two*u - ul)*hydhx;
247:       uyy = (-ut + two*u - ub)*hxdhy;
248:       f[row] = -uxx + -uyy + sc*lambda*PetscExpScalar(u);
249:     }
250:   }
251:   VecRestoreArray(X,&x);
252:   VecRestoreArray(F,&f);
253:   return 0;
254: }
255: /* --------------------  Evaluate Jacobian F'(x) -------------------- */

257: int FormJacobian(TS ts,double t,Vec X,Mat *J,Mat *B,MatStructure *flag,void *ptr)
258: {
259:   AppCtx *user = (AppCtx*)ptr;
260:   Mat     jac = *B;
261:   int     i,j,row,mx,my,col[5],ierr;
262:   Scalar  two = 2.0,one = 1.0,lambda,v[5],sc,*x;
263:   double  hx,hy,hxdhy,hydhx;


266:   mx         = user->mx;
267:   my         = user->my;
268:   lambda = user->param;

270:   hx    = 1.0 / (double)(mx-1);
271:   hy    = 1.0 / (double)(my-1);
272:   sc    = hx*hy;
273:   hxdhy = hx/hy;
274:   hydhx = hy/hx;

276:   VecGetArray(X,&x);
277:   for (j=0; j<my; j++) {
278:     for (i=0; i<mx; i++) {
279:       row = i + j*mx;
280:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
281:         MatSetValues(jac,1,&row,1,&row,&one,INSERT_VALUES);
282:         continue;
283:       }
284:       v[0] = hxdhy; col[0] = row - mx;
285:       v[1] = hydhx; col[1] = row - 1;
286:       v[2] = -two*(hydhx + hxdhy) + sc*lambda*PetscExpScalar(x[row]); col[2] = row;
287:       v[3] = hydhx; col[3] = row + 1;
288:       v[4] = hxdhy; col[4] = row + mx;
289:       MatSetValues(jac,1,&row,5,col,v,INSERT_VALUES);
290:     }
291:   }
292:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
293:   VecRestoreArray(X,&x);
294:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
295:   *flag = SAME_NONZERO_PATTERN;
296:   return 0;
297: }