1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] ="Solves the time independent Bratu problem using pseudo-timestepping."; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown /* 5*c4762a1bSJed Brown Concepts: TS^pseudo-timestepping 6*c4762a1bSJed Brown Concepts: pseudo-timestepping 7*c4762a1bSJed Brown Concepts: TS^nonlinear problems 8*c4762a1bSJed Brown Processors: 1 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown */ 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown This code demonstrates how one may solve a nonlinear problem 15*c4762a1bSJed Brown with pseudo-timestepping. In this simple example, the pseudo-timestep 16*c4762a1bSJed Brown is the same for all grid points, i.e., this is equivalent to using 17*c4762a1bSJed Brown the backward Euler method with a variable timestep. 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown Note: This example does not require pseudo-timestepping since it 20*c4762a1bSJed Brown is an easy nonlinear problem, but it is included to demonstrate how 21*c4762a1bSJed Brown the pseudo-timestepping may be done. 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown See snes/tutorials/ex4.c[ex4f.F] and 24*c4762a1bSJed Brown snes/tutorials/ex5.c[ex5f.F] where the problem is described 25*c4762a1bSJed Brown and solved using Newton's method alone. 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown ----------------------------------------------------------------------------- */ 28*c4762a1bSJed Brown /* 29*c4762a1bSJed Brown Include "petscts.h" to use the PETSc timestepping routines. Note that 30*c4762a1bSJed Brown this file automatically includes "petscsys.h" and other lower-level 31*c4762a1bSJed Brown PETSc include files. 32*c4762a1bSJed Brown */ 33*c4762a1bSJed Brown #include <petscts.h> 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown /* 36*c4762a1bSJed Brown Create an application context to contain data needed by the 37*c4762a1bSJed Brown application-provided call-back routines, FormJacobian() and 38*c4762a1bSJed Brown FormFunction(). 39*c4762a1bSJed Brown */ 40*c4762a1bSJed Brown typedef struct { 41*c4762a1bSJed Brown PetscReal param; /* test problem parameter */ 42*c4762a1bSJed Brown PetscInt mx; /* Discretization in x-direction */ 43*c4762a1bSJed Brown PetscInt my; /* Discretization in y-direction */ 44*c4762a1bSJed Brown } AppCtx; 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown /* 47*c4762a1bSJed Brown User-defined routines 48*c4762a1bSJed Brown */ 49*c4762a1bSJed Brown extern PetscErrorCode FormJacobian(TS,PetscReal,Vec,Mat,Mat,void*), FormFunction(TS,PetscReal,Vec,Vec,void*), FormInitialGuess(Vec,AppCtx*); 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown int main(int argc,char **argv) 52*c4762a1bSJed Brown { 53*c4762a1bSJed Brown TS ts; /* timestepping context */ 54*c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 55*c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 56*c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 57*c4762a1bSJed Brown PetscInt its,N; /* iterations for convergence */ 58*c4762a1bSJed Brown PetscErrorCode ierr; 59*c4762a1bSJed Brown PetscReal param_max = 6.81,param_min = 0.,dt; 60*c4762a1bSJed Brown PetscReal ftime; 61*c4762a1bSJed Brown PetscMPIInt size; 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 64*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 65*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only"); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown user.mx = 4; 68*c4762a1bSJed Brown user.my = 4; 69*c4762a1bSJed Brown user.param = 6.0; 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown /* 72*c4762a1bSJed Brown Allow user to set the grid dimensions and nonlinearity parameter at run-time 73*c4762a1bSJed Brown */ 74*c4762a1bSJed Brown PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL); 75*c4762a1bSJed Brown PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL); 76*c4762a1bSJed Brown N = user.mx*user.my; 77*c4762a1bSJed Brown dt = .5/PetscMax(user.mx,user.my); 78*c4762a1bSJed Brown PetscOptionsGetReal(NULL,NULL,"-param",&user.param,NULL); 79*c4762a1bSJed Brown if (user.param >= param_max || user.param <= param_min) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Parameter is out of range"); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown /* 82*c4762a1bSJed Brown Create vectors to hold the solution and function value 83*c4762a1bSJed Brown */ 84*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown /* 88*c4762a1bSJed Brown Create matrix to hold Jacobian. Preallocate 5 nonzeros per row 89*c4762a1bSJed Brown in the sparse matrix. Note that this is not the optimal strategy; see 90*c4762a1bSJed Brown the Performance chapter of the users manual for information on 91*c4762a1bSJed Brown preallocating memory in sparse matrices. 92*c4762a1bSJed Brown */ 93*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,5,0,&J);CHKERRQ(ierr); 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown /* 96*c4762a1bSJed Brown Create timestepper context 97*c4762a1bSJed Brown */ 98*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown /* 102*c4762a1bSJed Brown Tell the timestepper context where to compute solutions 103*c4762a1bSJed Brown */ 104*c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown /* 107*c4762a1bSJed Brown Provide the call-back for the nonlinear function we are 108*c4762a1bSJed Brown evaluating. Thus whenever the timestepping routines need the 109*c4762a1bSJed Brown function they will call this routine. Note the final argument 110*c4762a1bSJed Brown is the application context used by the call-back functions. 111*c4762a1bSJed Brown */ 112*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,FormFunction,&user);CHKERRQ(ierr); 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown /* 115*c4762a1bSJed Brown Set the Jacobian matrix and the function used to compute 116*c4762a1bSJed Brown Jacobians. 117*c4762a1bSJed Brown */ 118*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,J,J,FormJacobian,&user);CHKERRQ(ierr); 119*c4762a1bSJed Brown 120*c4762a1bSJed Brown /* 121*c4762a1bSJed Brown Form the initial guess for the problem 122*c4762a1bSJed Brown */ 123*c4762a1bSJed Brown ierr = FormInitialGuess(x,&user);CHKERRQ(ierr); 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown /* 126*c4762a1bSJed Brown This indicates that we are using pseudo timestepping to 127*c4762a1bSJed Brown find a steady state solution to the nonlinear problem. 128*c4762a1bSJed Brown */ 129*c4762a1bSJed Brown ierr = TSSetType(ts,TSPSEUDO);CHKERRQ(ierr); 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown /* 132*c4762a1bSJed Brown Set the initial time to start at (this is arbitrary for 133*c4762a1bSJed Brown steady state problems); and the initial timestep given above 134*c4762a1bSJed Brown */ 135*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown /* 138*c4762a1bSJed Brown Set a large number of timesteps and final duration time 139*c4762a1bSJed Brown to insure convergence to steady state. 140*c4762a1bSJed Brown */ 141*c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,10000);CHKERRQ(ierr); 142*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,1e12);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown /* 146*c4762a1bSJed Brown Use the default strategy for increasing the timestep 147*c4762a1bSJed Brown */ 148*c4762a1bSJed Brown ierr = TSPseudoSetTimeStep(ts,TSPseudoTimeStepDefault,0);CHKERRQ(ierr); 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown /* 151*c4762a1bSJed Brown Set any additional options from the options database. This 152*c4762a1bSJed Brown includes all options for the nonlinear and linear solvers used 153*c4762a1bSJed Brown internally the timestepping routines. 154*c4762a1bSJed Brown */ 155*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 156*c4762a1bSJed Brown 157*c4762a1bSJed Brown ierr = TSSetUp(ts);CHKERRQ(ierr); 158*c4762a1bSJed Brown 159*c4762a1bSJed Brown /* 160*c4762a1bSJed Brown Perform the solve. This is where the timestepping takes place. 161*c4762a1bSJed Brown */ 162*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 164*c4762a1bSJed Brown 165*c4762a1bSJed Brown /* 166*c4762a1bSJed Brown Get the number of steps 167*c4762a1bSJed Brown */ 168*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&its);CHKERRQ(ierr); 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of pseudo timesteps = %D final time %4.2e\n",its,(double)ftime);CHKERRQ(ierr); 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown /* 173*c4762a1bSJed Brown Free the data structures constructed above 174*c4762a1bSJed Brown */ 175*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 179*c4762a1bSJed Brown ierr = PetscFinalize(); 180*c4762a1bSJed Brown return ierr; 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 183*c4762a1bSJed Brown /* Bratu (Solid Fuel Ignition) Test Problem */ 184*c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 185*c4762a1bSJed Brown 186*c4762a1bSJed Brown /* -------------------- Form initial approximation ----------------- */ 187*c4762a1bSJed Brown 188*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec X,AppCtx *user) 189*c4762a1bSJed Brown { 190*c4762a1bSJed Brown PetscInt i,j,row,mx,my; 191*c4762a1bSJed Brown PetscErrorCode ierr; 192*c4762a1bSJed Brown PetscReal one = 1.0,lambda; 193*c4762a1bSJed Brown PetscReal temp1,temp,hx,hy; 194*c4762a1bSJed Brown PetscScalar *x; 195*c4762a1bSJed Brown 196*c4762a1bSJed Brown mx = user->mx; 197*c4762a1bSJed Brown my = user->my; 198*c4762a1bSJed Brown lambda = user->param; 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown hx = one / (PetscReal)(mx-1); 201*c4762a1bSJed Brown hy = one / (PetscReal)(my-1); 202*c4762a1bSJed Brown 203*c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 204*c4762a1bSJed Brown temp1 = lambda/(lambda + one); 205*c4762a1bSJed Brown for (j=0; j<my; j++) { 206*c4762a1bSJed Brown temp = (PetscReal)(PetscMin(j,my-j-1))*hy; 207*c4762a1bSJed Brown for (i=0; i<mx; i++) { 208*c4762a1bSJed Brown row = i + j*mx; 209*c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 210*c4762a1bSJed Brown x[row] = 0.0; 211*c4762a1bSJed Brown continue; 212*c4762a1bSJed Brown } 213*c4762a1bSJed Brown x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp)); 214*c4762a1bSJed Brown } 215*c4762a1bSJed Brown } 216*c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 217*c4762a1bSJed Brown return 0; 218*c4762a1bSJed Brown } 219*c4762a1bSJed Brown /* -------------------- Evaluate Function F(x) --------------------- */ 220*c4762a1bSJed Brown 221*c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 222*c4762a1bSJed Brown { 223*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 224*c4762a1bSJed Brown PetscErrorCode ierr; 225*c4762a1bSJed Brown PetscInt i,j,row,mx,my; 226*c4762a1bSJed Brown PetscReal two = 2.0,one = 1.0,lambda; 227*c4762a1bSJed Brown PetscReal hx,hy,hxdhy,hydhx; 228*c4762a1bSJed Brown PetscScalar ut,ub,ul,ur,u,uxx,uyy,sc,*f; 229*c4762a1bSJed Brown const PetscScalar *x; 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown mx = user->mx; 232*c4762a1bSJed Brown my = user->my; 233*c4762a1bSJed Brown lambda = user->param; 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown hx = one / (PetscReal)(mx-1); 236*c4762a1bSJed Brown hy = one / (PetscReal)(my-1); 237*c4762a1bSJed Brown sc = hx*hy; 238*c4762a1bSJed Brown hxdhy = hx/hy; 239*c4762a1bSJed Brown hydhx = hy/hx; 240*c4762a1bSJed Brown 241*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 242*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 243*c4762a1bSJed Brown for (j=0; j<my; j++) { 244*c4762a1bSJed Brown for (i=0; i<mx; i++) { 245*c4762a1bSJed Brown row = i + j*mx; 246*c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 247*c4762a1bSJed Brown f[row] = x[row]; 248*c4762a1bSJed Brown continue; 249*c4762a1bSJed Brown } 250*c4762a1bSJed Brown u = x[row]; 251*c4762a1bSJed Brown ub = x[row - mx]; 252*c4762a1bSJed Brown ul = x[row - 1]; 253*c4762a1bSJed Brown ut = x[row + mx]; 254*c4762a1bSJed Brown ur = x[row + 1]; 255*c4762a1bSJed Brown uxx = (-ur + two*u - ul)*hydhx; 256*c4762a1bSJed Brown uyy = (-ut + two*u - ub)*hxdhy; 257*c4762a1bSJed Brown f[row] = -uxx + -uyy + sc*lambda*PetscExpScalar(u); 258*c4762a1bSJed Brown } 259*c4762a1bSJed Brown } 260*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 262*c4762a1bSJed Brown return 0; 263*c4762a1bSJed Brown } 264*c4762a1bSJed Brown /* -------------------- Evaluate Jacobian F'(x) -------------------- */ 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown /* 267*c4762a1bSJed Brown Calculate the Jacobian matrix J(X,t). 268*c4762a1bSJed Brown 269*c4762a1bSJed Brown Note: We put the Jacobian in the preconditioner storage B instead of J. This 270*c4762a1bSJed Brown way we can give the -snes_mf_operator flag to check our work. This replaces 271*c4762a1bSJed Brown J with a finite difference approximation, using our analytic Jacobian B for 272*c4762a1bSJed Brown the preconditioner. 273*c4762a1bSJed Brown */ 274*c4762a1bSJed Brown PetscErrorCode FormJacobian(TS ts,PetscReal t,Vec X,Mat J,Mat B,void *ptr) 275*c4762a1bSJed Brown { 276*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 277*c4762a1bSJed Brown PetscInt i,j,row,mx,my,col[5]; 278*c4762a1bSJed Brown PetscErrorCode ierr; 279*c4762a1bSJed Brown PetscScalar two = 2.0,one = 1.0,lambda,v[5],sc; 280*c4762a1bSJed Brown const PetscScalar *x; 281*c4762a1bSJed Brown PetscReal hx,hy,hxdhy,hydhx; 282*c4762a1bSJed Brown 283*c4762a1bSJed Brown 284*c4762a1bSJed Brown mx = user->mx; 285*c4762a1bSJed Brown my = user->my; 286*c4762a1bSJed Brown lambda = user->param; 287*c4762a1bSJed Brown 288*c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx-1); 289*c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my-1); 290*c4762a1bSJed Brown sc = hx*hy; 291*c4762a1bSJed Brown hxdhy = hx/hy; 292*c4762a1bSJed Brown hydhx = hy/hx; 293*c4762a1bSJed Brown 294*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 295*c4762a1bSJed Brown for (j=0; j<my; j++) { 296*c4762a1bSJed Brown for (i=0; i<mx; i++) { 297*c4762a1bSJed Brown row = i + j*mx; 298*c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx-1 || j == my-1) { 299*c4762a1bSJed Brown ierr = MatSetValues(B,1,&row,1,&row,&one,INSERT_VALUES);CHKERRQ(ierr); 300*c4762a1bSJed Brown continue; 301*c4762a1bSJed Brown } 302*c4762a1bSJed Brown v[0] = hxdhy; col[0] = row - mx; 303*c4762a1bSJed Brown v[1] = hydhx; col[1] = row - 1; 304*c4762a1bSJed Brown v[2] = -two*(hydhx + hxdhy) + sc*lambda*PetscExpScalar(x[row]); col[2] = row; 305*c4762a1bSJed Brown v[3] = hydhx; col[3] = row + 1; 306*c4762a1bSJed Brown v[4] = hxdhy; col[4] = row + mx; 307*c4762a1bSJed Brown ierr = MatSetValues(B,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr); 308*c4762a1bSJed Brown } 309*c4762a1bSJed Brown } 310*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 311*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 312*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 313*c4762a1bSJed Brown if (J != B) { 314*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 315*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 316*c4762a1bSJed Brown } 317*c4762a1bSJed Brown return 0; 318*c4762a1bSJed Brown } 319*c4762a1bSJed Brown 320*c4762a1bSJed Brown 321*c4762a1bSJed Brown /*TEST 322*c4762a1bSJed Brown 323*c4762a1bSJed Brown test: 324*c4762a1bSJed Brown 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 325*c4762a1bSJed Brown 326*c4762a1bSJed Brown test: 327*c4762a1bSJed Brown suffix: 2 328*c4762a1bSJed Brown args: -ts_monitor_pseudo -ts_pseudo_frtol 1.e-5 329*c4762a1bSJed Brown 330*c4762a1bSJed Brown TEST*/ 331*c4762a1bSJed Brown 332