xref: /petsc/src/binding/petsc4py/demo/legacy/bratu2d/bratu2df90.f90 (revision 5a48edb989d3ea10d6aff6c0e26d581c18691deb)
1*55a74a43SLisandro Dalcin! file: bratu2df90.f90
2*55a74a43SLisandro Dalcin! to build a Python module, use this:
3*55a74a43SLisandro Dalcin! $$ f2py -m bratu2df90 -c bratu2df90.f90
4*55a74a43SLisandro Dalcin
5*55a74a43SLisandro Dalcinsubroutine bratu2d (m, n, alpha, x, f)
6*55a74a43SLisandro Dalcin  !f2py intent(hide) :: m = shape(x,0)
7*55a74a43SLisandro Dalcin  !f2py intent(hide) :: n = shape(x,1)
8*55a74a43SLisandro Dalcin  integer                          :: m, n
9*55a74a43SLisandro Dalcin  real(kind=8)                     :: alpha
10*55a74a43SLisandro Dalcin  real(kind=8), intent(in), target :: x(m,n)
11*55a74a43SLisandro Dalcin  real(kind=8), intent(inout)      :: f(m,n)
12*55a74a43SLisandro Dalcin  real(kind=8) :: hx, hy
13*55a74a43SLisandro Dalcin  real(kind=8), pointer, &
14*55a74a43SLisandro Dalcin       dimension(:,:) :: u, uN, uS, uE, uW
15*55a74a43SLisandro Dalcin  ! setup 5-points stencil
16*55a74a43SLisandro Dalcin  u  => x(2:m-1, 2:n-1) ! center
17*55a74a43SLisandro Dalcin  uN => x(2:m-1, 1:n-2) ! north
18*55a74a43SLisandro Dalcin  uS => x(2:m-1, 3:n  ) ! south
19*55a74a43SLisandro Dalcin  uW => x(1:m-2, 2:n-1) ! west
20*55a74a43SLisandro Dalcin  uE => x(3:m,   2:n-1) ! east
21*55a74a43SLisandro Dalcin  ! compute nonlinear function
22*55a74a43SLisandro Dalcin  hx = 1.0/(m-1) ! x grid spacing
23*55a74a43SLisandro Dalcin  hy = 1.0/(n-1) ! y grid spacing
24*55a74a43SLisandro Dalcin  f(:,:) = x
25*55a74a43SLisandro Dalcin  f(2:m-1, 2:n-1) =  &
26*55a74a43SLisandro Dalcin         (2*u - uE - uW) * (hy/hx) &
27*55a74a43SLisandro Dalcin       + (2*u - uN - uS) * (hx/hy) &
28*55a74a43SLisandro Dalcin       - alpha * exp(u)  * (hx*hy)
29*55a74a43SLisandro Dalcinend subroutine bratu2d
30