Actual source code: ex1.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

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

  4: /*
  5:    Concepts: TS^pseudo-timestepping
  6:    Concepts: pseudo-timestepping
  7:    Concepts: TS^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,void*), FormFunction(TS,PetscReal,Vec,Vec,void*), FormInitialGuess(Vec,AppCtx*);

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

 63:   PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
 64:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 65:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only");

 67:   user.mx    = 4;
 68:   user.my    = 4;
 69:   user.param = 6.0;

 71:   /*
 72:      Allow user to set the grid dimensions and nonlinearity parameter at run-time
 73:   */
 74:   PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL);
 75:   PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL);
 76:   N  = user.mx*user.my;
 77:   dt = .5/PetscMax(user.mx,user.my);
 78:   PetscOptionsGetReal(NULL,NULL,"-param",&user.param,NULL);
 79:   if (user.param >= param_max || user.param <= param_min) SETERRQ(PETSC_COMM_SELF,1,"Parameter is out of range");

 81:   /*
 82:       Create vectors to hold the solution and function value
 83:   */
 84:   VecCreateSeq(PETSC_COMM_SELF,N,&x);
 85:   VecDuplicate(x,&r);

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

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

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

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

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

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

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

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

137:   /*
138:       Set a large number of timesteps and final duration time
139:      to insure convergence to steady state.
140:   */
141:   TSSetMaxSteps(ts,10000);
142:   TSSetMaxTime(ts,1e12);
143:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

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

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

157:   TSSetUp(ts);

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

165:   /*
166:       Get the number of steps
167:   */
168:   TSGetStepNumber(ts,&its);

170:   PetscPrintf(PETSC_COMM_WORLD,"Number of pseudo timesteps = %D final time %4.2e\n",its,(double)ftime);

172:   /*
173:      Free the data structures constructed above
174:   */
175:   VecDestroy(&x);
176:   VecDestroy(&r);
177:   MatDestroy(&J);
178:   TSDestroy(&ts);
179:   PetscFinalize();
180:   return ierr;
181: }
182: /* ------------------------------------------------------------------ */
183: /*           Bratu (Solid Fuel Ignition) Test Problem                 */
184: /* ------------------------------------------------------------------ */

186: /* --------------------  Form initial approximation ----------------- */

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

196:   mx     = user->mx;
197:   my     = user->my;
198:   lambda = user->param;

200:   hx = one / (PetscReal)(mx-1);
201:   hy = one / (PetscReal)(my-1);

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

221: PetscErrorCode FormFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
222: {
223:   AppCtx            *user = (AppCtx*)ptr;
224:   PetscErrorCode    ierr;
225:   PetscInt          i,j,row,mx,my;
226:   PetscReal         two = 2.0,one = 1.0,lambda;
227:   PetscReal         hx,hy,hxdhy,hydhx;
228:   PetscScalar       ut,ub,ul,ur,u,uxx,uyy,sc,*f;
229:   const PetscScalar *x;

231:   mx     = user->mx;
232:   my     = user->my;
233:   lambda = user->param;

235:   hx    = one / (PetscReal)(mx-1);
236:   hy    = one / (PetscReal)(my-1);
237:   sc    = hx*hy;
238:   hxdhy = hx/hy;
239:   hydhx = hy/hx;

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

266: /*
267:    Calculate the Jacobian matrix J(X,t).

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


284:   mx     = user->mx;
285:   my     = user->my;
286:   lambda = user->param;

288:   hx    = 1.0 / (PetscReal)(mx-1);
289:   hy    = 1.0 / (PetscReal)(my-1);
290:   sc    = hx*hy;
291:   hxdhy = hx/hy;
292:   hydhx = hy/hx;

294:   VecGetArrayRead(X,&x);
295:   for (j=0; j<my; j++) {
296:     for (i=0; i<mx; i++) {
297:       row = i + j*mx;
298:       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
299:         MatSetValues(B,1,&row,1,&row,&one,INSERT_VALUES);
300:         continue;
301:       }
302:       v[0] = hxdhy; col[0] = row - mx;
303:       v[1] = hydhx; col[1] = row - 1;
304:       v[2] = -two*(hydhx + hxdhy) + sc*lambda*PetscExpScalar(x[row]); col[2] = row;
305:       v[3] = hydhx; col[3] = row + 1;
306:       v[4] = hxdhy; col[4] = row + mx;
307:       MatSetValues(B,1,&row,5,col,v,INSERT_VALUES);
308:     }
309:   }
310:   VecRestoreArrayRead(X,&x);
311:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
312:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
313:   if (J != B) {
314:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
315:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
316:   }
317:   return 0;
318: }


321: /*TEST

323:     test:
324:       args: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_monitor_pseudo -snes_atol 1.e-7 -ts_pseudo_frtol 1.e-5 -ts_view draw:tikz:fig.tex

326:     test:
327:       suffix: 2
328:       args: -ts_monitor_pseudo -ts_pseudo_frtol 1.e-5 

330: TEST*/