1d8606c27SBarry Smith! 2d8606c27SBarry Smith! Solves the time dependent Bratu problem using pseudo-timestepping 3d8606c27SBarry Smith! 4d8606c27SBarry Smith! This code demonstrates how one may solve a nonlinear problem 5d8606c27SBarry Smith! with pseudo-timestepping. In this simple example, the pseudo-timestep 6d8606c27SBarry Smith! is the same for all grid points, i.e., this is equivalent to using 7d8606c27SBarry Smith! the backward Euler method with a variable timestep. 8d8606c27SBarry Smith! 9d8606c27SBarry Smith! Note: This example does not require pseudo-timestepping since it 10d8606c27SBarry Smith! is an easy nonlinear problem, but it is included to demonstrate how 11d8606c27SBarry Smith! the pseudo-timestepping may be done. 12d8606c27SBarry Smith! 13d8606c27SBarry Smith! See snes/tutorials/ex4.c[ex4f.F] and 14d8606c27SBarry Smith! snes/tutorials/ex5.c[ex5f.F] where the problem is described 15d8606c27SBarry Smith! and solved using the method of Newton alone. 16d8606c27SBarry Smith! 17d8606c27SBarry Smith! 18d8606c27SBarry Smith!23456789012345678901234567890123456789012345678901234567890123456789012 19d8606c27SBarry Smith#include <petsc/finclude/petscts.h> 20*01fa2b5aSMartin Diehlmodule ex1fmodule 21d8606c27SBarry Smith use petscts 22d8606c27SBarry Smith implicit none 23e7a95102SMartin Diehlcontains 24e7a95102SMartin Diehl! 25e7a95102SMartin Diehl! -------------------- Evaluate Function F(x) --------------------- 26e7a95102SMartin Diehl! 27e7a95102SMartin Diehl subroutine FormFunction(ts, t, X, F, user, ierr) 28d8606c27SBarry Smith 29e7a95102SMartin Diehl TS ts 30e7a95102SMartin Diehl PetscReal t 31e7a95102SMartin Diehl Vec X, F 32e7a95102SMartin Diehl PetscReal user(3) 33e7a95102SMartin Diehl PetscErrorCode ierr 34e7a95102SMartin Diehl PetscInt i, j, row, mx, my 35e7a95102SMartin Diehl PetscReal two, lambda 36e7a95102SMartin Diehl PetscReal hx, hy, hxdhy, hydhx 37e7a95102SMartin Diehl PetscScalar ut, ub, ul, ur, u 38e7a95102SMartin Diehl PetscScalar uxx, uyy, sc 39e7a95102SMartin Diehl PetscScalar, pointer :: xx(:), ff(:) 40e7a95102SMartin Diehl PetscInt, parameter :: param = 1, lmx = 2, lmy = 3 41e7a95102SMartin Diehl 42e7a95102SMartin Diehl two = 2.0 43e7a95102SMartin Diehl 44e7a95102SMartin Diehl mx = int(user(lmx)) 45e7a95102SMartin Diehl my = int(user(lmy)) 46e7a95102SMartin Diehl lambda = user(param) 47e7a95102SMartin Diehl 48e7a95102SMartin Diehl hx = 1.0/real(mx - 1) 49e7a95102SMartin Diehl hy = 1.0/real(my - 1) 50e7a95102SMartin Diehl sc = hx*hy 51e7a95102SMartin Diehl hxdhy = hx/hy 52e7a95102SMartin Diehl hydhx = hy/hx 53e7a95102SMartin Diehl 54e7a95102SMartin Diehl PetscCall(VecGetArrayRead(X, xx, ierr)) 55e7a95102SMartin Diehl PetscCall(VecGetArray(F, ff, ierr)) 56e7a95102SMartin Diehl do j = 1, my 57e7a95102SMartin Diehl do i = 1, mx 58e7a95102SMartin Diehl row = i + (j - 1)*mx 59e7a95102SMartin Diehl if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then 60e7a95102SMartin Diehl ff(row) = xx(row) 61e7a95102SMartin Diehl else 62e7a95102SMartin Diehl u = xx(row) 63e7a95102SMartin Diehl ub = xx(row - mx) 64e7a95102SMartin Diehl ul = xx(row - 1) 65e7a95102SMartin Diehl ut = xx(row + mx) 66e7a95102SMartin Diehl ur = xx(row + 1) 67e7a95102SMartin Diehl uxx = (-ur + two*u - ul)*hydhx 68e7a95102SMartin Diehl uyy = (-ut + two*u - ub)*hxdhy 69e7a95102SMartin Diehl ff(row) = -uxx - uyy + sc*lambda*exp(u) 70e7a95102SMartin Diehl end if 71e7a95102SMartin Diehl end do 72e7a95102SMartin Diehl end do 73e7a95102SMartin Diehl 74e7a95102SMartin Diehl PetscCall(VecRestoreArrayRead(X, xx, ierr)) 75e7a95102SMartin Diehl PetscCall(VecRestoreArray(F, ff, ierr)) 76e7a95102SMartin Diehl end 77e7a95102SMartin Diehl! 78e7a95102SMartin Diehl! -------------------- Evaluate Jacobian of F(x) -------------------- 79e7a95102SMartin Diehl! 80e7a95102SMartin Diehl subroutine FormJacobian(ts, ctime, X, JJ, B, user, ierr) 81e7a95102SMartin Diehl 82e7a95102SMartin Diehl TS ts 83e7a95102SMartin Diehl Vec X 84e7a95102SMartin Diehl Mat JJ, B 85e7a95102SMartin Diehl PetscReal user(3), ctime 86e7a95102SMartin Diehl PetscErrorCode ierr 87e7a95102SMartin Diehl Mat jac 88e7a95102SMartin Diehl PetscInt i, j, row(1), mx, my 89e7a95102SMartin Diehl PetscInt col(5), i1, i5 90e7a95102SMartin Diehl PetscScalar two, one, lambda 91e7a95102SMartin Diehl PetscScalar v(5), sc 92e7a95102SMartin Diehl PetscScalar, pointer :: xx(:) 93e7a95102SMartin Diehl PetscReal hx, hy, hxdhy, hydhx 94e7a95102SMartin Diehl 95e7a95102SMartin Diehl PetscInt, parameter :: param = 1, lmx = 2, lmy = 3 96e7a95102SMartin Diehl 97e7a95102SMartin Diehl i1 = 1 98e7a95102SMartin Diehl i5 = 5 99e7a95102SMartin Diehl jac = B 100e7a95102SMartin Diehl two = 2.0 101e7a95102SMartin Diehl one = 1.0 102e7a95102SMartin Diehl 103e7a95102SMartin Diehl mx = int(user(lmx)) 104e7a95102SMartin Diehl my = int(user(lmy)) 105e7a95102SMartin Diehl lambda = user(param) 106e7a95102SMartin Diehl 107e7a95102SMartin Diehl hx = 1.0/real(mx - 1) 108e7a95102SMartin Diehl hy = 1.0/real(my - 1) 109e7a95102SMartin Diehl sc = hx*hy 110e7a95102SMartin Diehl hxdhy = hx/hy 111e7a95102SMartin Diehl hydhx = hy/hx 112e7a95102SMartin Diehl 113e7a95102SMartin Diehl PetscCall(VecGetArrayRead(X, xx, ierr)) 114e7a95102SMartin Diehl do j = 1, my 115e7a95102SMartin Diehl do i = 1, mx 116e7a95102SMartin Diehl! 117e7a95102SMartin Diehl! When inserting into PETSc matrices, indices start at 0 118e7a95102SMartin Diehl! 119e7a95102SMartin Diehl row(1) = i - 1 + (j - 1)*mx 120e7a95102SMartin Diehl if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then 121e7a95102SMartin Diehl PetscCall(MatSetValues(jac, i1, [row], i1, [row], [one], INSERT_VALUES, ierr)) 122e7a95102SMartin Diehl else 123e7a95102SMartin Diehl v(1) = hxdhy 124e7a95102SMartin Diehl col(1) = row(1) - mx 125e7a95102SMartin Diehl v(2) = hydhx 126e7a95102SMartin Diehl col(2) = row(1) - 1 127e7a95102SMartin Diehl v(3) = -two*(hydhx + hxdhy) + sc*lambda*exp(xx(row(1))) 128e7a95102SMartin Diehl col(3) = row(1) 129e7a95102SMartin Diehl v(4) = hydhx 130e7a95102SMartin Diehl col(4) = row(1) + 1 131e7a95102SMartin Diehl v(5) = hxdhy 132e7a95102SMartin Diehl col(5) = row(1) + mx 133e7a95102SMartin Diehl PetscCall(MatSetValues(jac, i1, [row], i5, col, v, INSERT_VALUES, ierr)) 134e7a95102SMartin Diehl end if 135e7a95102SMartin Diehl end do 136e7a95102SMartin Diehl end do 137e7a95102SMartin Diehl PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY, ierr)) 138e7a95102SMartin Diehl PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY, ierr)) 139e7a95102SMartin Diehl PetscCall(VecRestoreArray(X, xx, ierr)) 140e7a95102SMartin Diehl end 141e7a95102SMartin Diehl! 142e7a95102SMartin Diehl! -------------------- Form initial approximation ----------------- 143e7a95102SMartin Diehl! 144e7a95102SMartin Diehl subroutine FormInitialGuess(X, user, ierr) 145e7a95102SMartin Diehl 146e7a95102SMartin Diehl Vec X 147e7a95102SMartin Diehl PetscReal user(3) 148e7a95102SMartin Diehl PetscInt i, j, row, mx, my 149e7a95102SMartin Diehl PetscErrorCode ierr 150e7a95102SMartin Diehl PetscReal one, lambda 151e7a95102SMartin Diehl PetscReal temp1, temp, hx, hy 152e7a95102SMartin Diehl PetscScalar, pointer :: xx(:) 153e7a95102SMartin Diehl PetscInt, parameter :: param = 1, lmx = 2, lmy = 3 154e7a95102SMartin Diehl 155e7a95102SMartin Diehl one = 1.0 156e7a95102SMartin Diehl 157e7a95102SMartin Diehl mx = int(user(lmx)) 158e7a95102SMartin Diehl my = int(user(lmy)) 159e7a95102SMartin Diehl lambda = user(param) 160e7a95102SMartin Diehl 161e7a95102SMartin Diehl hy = one/(my - 1) 162e7a95102SMartin Diehl hx = one/(mx - 1) 163e7a95102SMartin Diehl 164e7a95102SMartin Diehl PetscCall(VecGetArray(X, xx, ierr)) 165e7a95102SMartin Diehl temp1 = lambda/(lambda + one) 166e7a95102SMartin Diehl do j = 1, my 167e7a95102SMartin Diehl temp = min(j - 1, my - j)*hy 168e7a95102SMartin Diehl do i = 1, mx 169e7a95102SMartin Diehl row = i + (j - 1)*mx 170e7a95102SMartin Diehl if (i == 1 .or. j == 1 .or. i == mx .or. j == my) then 171e7a95102SMartin Diehl xx(row) = 0.0 172e7a95102SMartin Diehl else 173e7a95102SMartin Diehl xx(row) = temp1*sqrt(min(min(i - 1, mx - i)*hx, temp)) 174e7a95102SMartin Diehl end if 175e7a95102SMartin Diehl end do 176e7a95102SMartin Diehl end do 177e7a95102SMartin Diehl PetscCall(VecRestoreArray(X, xx, ierr)) 178e7a95102SMartin Diehl end 179*01fa2b5aSMartin Diehlend module ex1fmodule 180e7a95102SMartin Diehlprogram main 181e7a95102SMartin Diehl use petscts 182*01fa2b5aSMartin Diehl use ex1fmodule 183e7a95102SMartin Diehl implicit none 184d8606c27SBarry Smith! 185d8606c27SBarry Smith! Create an application context to contain data needed by the 186d8606c27SBarry Smith! application-provided call-back routines, FormJacobian() and 187d8606c27SBarry Smith! FormFunction(). We use a double precision array with three 188d8606c27SBarry Smith! entries indexed by param, lmx, lmy. 189d8606c27SBarry Smith! 190d8606c27SBarry Smith PetscReal user(3) 191e7a95102SMartin Diehl PetscInt i5 192e7a95102SMartin Diehl PetscInt, parameter :: param = 1, lmx = 2, lmy = 3 193d8606c27SBarry Smith! 194d8606c27SBarry Smith! Data for problem 195d8606c27SBarry Smith! 196d8606c27SBarry Smith TS ts 197d8606c27SBarry Smith Vec x, r 198d8606c27SBarry Smith Mat J 199d8606c27SBarry Smith PetscInt its, N, i1000, itmp 200d8606c27SBarry Smith PetscBool flg 201d8606c27SBarry Smith PetscErrorCode ierr 202d8606c27SBarry Smith PetscReal param_max, param_min, dt 203d8606c27SBarry Smith PetscReal tmax 204d8606c27SBarry Smith PetscReal ftime 205d8606c27SBarry Smith TSConvergedReason reason 206d8606c27SBarry Smith 207d8606c27SBarry Smith i5 = 5 208d8606c27SBarry Smith param_max = 6.81 209d8606c27SBarry Smith param_min = 0 210d8606c27SBarry Smith 211d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 212d8606c27SBarry Smith user(lmx) = 4 213d8606c27SBarry Smith user(lmy) = 4 214d8606c27SBarry Smith user(param) = 6.0 215d8606c27SBarry Smith 216d8606c27SBarry Smith! 217d8606c27SBarry Smith! Allow user to set the grid dimensions and nonlinearity parameter at run-time 218d8606c27SBarry Smith! 219d8606c27SBarry Smith PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-mx', user(lmx), flg, ierr)) 220d8606c27SBarry Smith itmp = 4 221d8606c27SBarry Smith PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-my', itmp, flg, ierr)) 222d8606c27SBarry Smith user(lmy) = itmp 223d8606c27SBarry Smith PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-param', user(param), flg, ierr)) 2244820e4eaSBarry Smith if (user(param) >= param_max .or. user(param) <= param_min) then 225d8606c27SBarry Smith print *, 'Parameter is out of range' 226d8606c27SBarry Smith end if 2274820e4eaSBarry Smith if (user(lmx) > user(lmy)) then 228d8606c27SBarry Smith dt = .5/user(lmx) 229d8606c27SBarry Smith else 230d8606c27SBarry Smith dt = .5/user(lmy) 231d8606c27SBarry Smith end if 232d8606c27SBarry Smith PetscCallA(PetscOptionsGetReal(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, '-dt', dt, flg, ierr)) 233d8606c27SBarry Smith N = int(user(lmx)*user(lmy)) 234d8606c27SBarry Smith 235d8606c27SBarry Smith! 236d8606c27SBarry Smith! Create vectors to hold the solution and function value 237d8606c27SBarry Smith! 238d8606c27SBarry Smith PetscCallA(VecCreateSeq(PETSC_COMM_SELF, N, x, ierr)) 239d8606c27SBarry Smith PetscCallA(VecDuplicate(x, r, ierr)) 240d8606c27SBarry Smith 241d8606c27SBarry Smith! 242d8606c27SBarry Smith! Create matrix to hold Jacobian. Preallocate 5 nonzeros per row 243d8606c27SBarry Smith! in the sparse matrix. Note that this is not the optimal strategy see 244d8606c27SBarry Smith! the Performance chapter of the users manual for information on 245d8606c27SBarry Smith! preallocating memory in sparse matrices. 246d8606c27SBarry Smith! 247d8606c27SBarry Smith i5 = 5 2485d83a8b1SBarry Smith PetscCallA(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, i5, PETSC_NULL_INTEGER_ARRAY, J, ierr)) 249d8606c27SBarry Smith 250d8606c27SBarry Smith! 251d8606c27SBarry Smith! Create timestepper context 252d8606c27SBarry Smith! 253d8606c27SBarry Smith 254d8606c27SBarry Smith PetscCallA(TSCreate(PETSC_COMM_WORLD, ts, ierr)) 255d8606c27SBarry Smith PetscCallA(TSSetProblemType(ts, TS_NONLINEAR, ierr)) 256d8606c27SBarry Smith 257d8606c27SBarry Smith! 258d8606c27SBarry Smith! Tell the timestepper context where to compute solutions 259d8606c27SBarry Smith! 260d8606c27SBarry Smith 261d8606c27SBarry Smith PetscCallA(TSSetSolution(ts, x, ierr)) 262d8606c27SBarry Smith 263d8606c27SBarry Smith! 264d8606c27SBarry Smith! Provide the call-back for the nonlinear function we are 265d8606c27SBarry Smith! evaluating. Thus whenever the timestepping routines need the 266d8606c27SBarry Smith! function they will call this routine. Note the final argument 267d8606c27SBarry Smith! is the application context used by the call-back functions. 268d8606c27SBarry Smith! 269d8606c27SBarry Smith 270d8606c27SBarry Smith PetscCallA(TSSetRHSFunction(ts, PETSC_NULL_VEC, FormFunction, user, ierr)) 271d8606c27SBarry Smith 272d8606c27SBarry Smith! 273d8606c27SBarry Smith! Set the Jacobian matrix and the function used to compute 274d8606c27SBarry Smith! Jacobians. 275d8606c27SBarry Smith! 276d8606c27SBarry Smith 277d8606c27SBarry Smith PetscCallA(TSSetRHSJacobian(ts, J, J, FormJacobian, user, ierr)) 278d8606c27SBarry Smith 279d8606c27SBarry Smith! 280d8606c27SBarry Smith! For the initial guess for the problem 281d8606c27SBarry Smith! 282d8606c27SBarry Smith 283d8606c27SBarry Smith PetscCallA(FormInitialGuess(x, user, ierr)) 284d8606c27SBarry Smith 285d8606c27SBarry Smith! 286d8606c27SBarry Smith! This indicates that we are using pseudo timestepping to 287d8606c27SBarry Smith! find a steady state solution to the nonlinear problem. 288d8606c27SBarry Smith! 289d8606c27SBarry Smith 290d8606c27SBarry Smith PetscCallA(TSSetType(ts, TSPSEUDO, ierr)) 291d8606c27SBarry Smith 292d8606c27SBarry Smith! 293d8606c27SBarry Smith! Set the initial time to start at (this is arbitrary for 294d8606c27SBarry Smith! steady state problems and the initial timestep given above 295d8606c27SBarry Smith! 296d8606c27SBarry Smith 297d8606c27SBarry Smith PetscCallA(TSSetTimeStep(ts, dt, ierr)) 298d8606c27SBarry Smith 299d8606c27SBarry Smith! 300d8606c27SBarry Smith! Set a large number of timesteps and final duration time 301d8606c27SBarry Smith! to insure convergence to steady state. 302d8606c27SBarry Smith! 303d8606c27SBarry Smith i1000 = 1000 304d8606c27SBarry Smith tmax = 1.e12 305d8606c27SBarry Smith PetscCallA(TSSetMaxSteps(ts, i1000, ierr)) 306d8606c27SBarry Smith PetscCallA(TSSetMaxTime(ts, tmax, ierr)) 307d8606c27SBarry Smith PetscCallA(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER, ierr)) 308d8606c27SBarry Smith 309d8606c27SBarry Smith! 310d8606c27SBarry Smith! Set any additional options from the options database. This 311d8606c27SBarry Smith! includes all options for the nonlinear and linear solvers used 312d8606c27SBarry Smith! internally the timestepping routines. 313d8606c27SBarry Smith! 314d8606c27SBarry Smith 315d8606c27SBarry Smith PetscCallA(TSSetFromOptions(ts, ierr)) 316d8606c27SBarry Smith 317d8606c27SBarry Smith PetscCallA(TSSetUp(ts, ierr)) 318d8606c27SBarry Smith 319d8606c27SBarry Smith! 320d8606c27SBarry Smith! Perform the solve. This is where the timestepping takes place. 321d8606c27SBarry Smith! 322d8606c27SBarry Smith PetscCallA(TSSolve(ts, x, ierr)) 323d8606c27SBarry Smith PetscCallA(TSGetSolveTime(ts, ftime, ierr)) 324d8606c27SBarry Smith PetscCallA(TSGetStepNumber(ts, its, ierr)) 325d8606c27SBarry Smith PetscCallA(TSGetConvergedReason(ts, reason, ierr)) 326d8606c27SBarry Smith 327d8606c27SBarry Smith write (6, 100) its, ftime, reason 328d8606c27SBarry Smith100 format('Number of pseudo time-steps ', i5, ' final time ', 1pe9.2, ' reason ', i3) 329d8606c27SBarry Smith 330d8606c27SBarry Smith! 331d8606c27SBarry Smith! Free the data structures constructed above 332d8606c27SBarry Smith! 333d8606c27SBarry Smith 334d8606c27SBarry Smith PetscCallA(VecDestroy(x, ierr)) 335d8606c27SBarry Smith PetscCallA(VecDestroy(r, ierr)) 336d8606c27SBarry Smith PetscCallA(MatDestroy(J, ierr)) 337d8606c27SBarry Smith PetscCallA(TSDestroy(ts, ierr)) 338d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 339d8606c27SBarry Smithend 340d8606c27SBarry Smith!/*TEST 341d8606c27SBarry Smith! 342d8606c27SBarry Smith! test: 343d8606c27SBarry Smith! TODO: broken 344d8606c27SBarry Smith! args: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_monitor_pseudo -ts_max_snes_failures 3 -ts_pseudo_frtol 1.e-5 -snes_stol 1e-5 345d8606c27SBarry Smith! 346d8606c27SBarry Smith!TEST*/ 347