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