155a74a43SLisandro Dalcin! --------------------------------------------------------------------- 255a74a43SLisandro Dalcin! 355a74a43SLisandro Dalcin! Solid Fuel Ignition (SFI) problem. This problem is modeled by 455a74a43SLisandro Dalcin! the partial differential equation 555a74a43SLisandro Dalcin! 655a74a43SLisandro Dalcin! -Laplacian(u) - lambda * exp(u) = 0, 0 < x,y < 1, 755a74a43SLisandro Dalcin! 855a74a43SLisandro Dalcin! with boundary conditions 955a74a43SLisandro Dalcin! 1055a74a43SLisandro Dalcin! u = 0 for x = 0, x = 1, y = 0, y = 1, 1155a74a43SLisandro Dalcin! 1255a74a43SLisandro Dalcin! A finite difference approximation with the usual 5-point stencil 1355a74a43SLisandro Dalcin! is used to discretize the boundary value problem to obtain a 1455a74a43SLisandro Dalcin! nonlinear system of equations. The problem is solved in a 2D 1555a74a43SLisandro Dalcin! rectangular domain, using distributed arrays (DAs) to partition 1655a74a43SLisandro Dalcin! the parallel grid. 1755a74a43SLisandro Dalcin! 1855a74a43SLisandro Dalcin! -------------------------------------------------------------------- 1955a74a43SLisandro Dalcin 2055a74a43SLisandro Dalcin#include <petsc/finclude/petscdm.h> 2155a74a43SLisandro Dalcin 2255a74a43SLisandro Dalcinmodule Bratu2D 2355a74a43SLisandro Dalcin 2455a74a43SLisandro Dalcin use petsc 2555a74a43SLisandro Dalcin implicit none 2655a74a43SLisandro Dalcin 2755a74a43SLisandro Dalcin type gridinfo 2855a74a43SLisandro Dalcin PetscInt mx, xs, xe, xm, gxs, gxe, gxm 2955a74a43SLisandro Dalcin PetscInt my, ys, ye, ym, gys, gye, gym 3055a74a43SLisandro Dalcin end type gridinfo 3155a74a43SLisandro Dalcin 3255a74a43SLisandro Dalcincontains 3355a74a43SLisandro Dalcin 3455a74a43SLisandro Dalcin subroutine GetGridInfo(da, grd, ierr) 3555a74a43SLisandro Dalcin implicit none 3655a74a43SLisandro Dalcin DM da 3755a74a43SLisandro Dalcin type(gridinfo) grd 3855a74a43SLisandro Dalcin PetscErrorCode ierr 3955a74a43SLisandro Dalcin ! 4055a74a43SLisandro Dalcin PetscCall(DMDAGetInfo(da, PETSC_NULL_INTEGER, grd%mx, grd%my, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, PETSC_NULL_INTEGER, ierr)) 4155a74a43SLisandro Dalcin PetscCall(DMDAGetCorners(da, grd%xs, grd%ys, PETSC_NULL_INTEGER, grd%xm, grd%ym, PETSC_NULL_INTEGER, ierr)) 4255a74a43SLisandro Dalcin PetscCall(DMDAGetGhostCorners(da, grd%gxs, grd%gys, PETSC_NULL_INTEGER, grd%gxm, grd%gym, PETSC_NULL_INTEGER, ierr)) 4355a74a43SLisandro Dalcin 4455a74a43SLisandro Dalcin grd%xs = grd%xs + 1 4555a74a43SLisandro Dalcin grd%ys = grd%ys + 1 4655a74a43SLisandro Dalcin grd%gxs = grd%gxs + 1 4755a74a43SLisandro Dalcin grd%gys = grd%gys + 1 4855a74a43SLisandro Dalcin 4955a74a43SLisandro Dalcin grd%ye = grd%ys + grd%ym - 1 5055a74a43SLisandro Dalcin grd%xe = grd%xs + grd%xm - 1 5155a74a43SLisandro Dalcin grd%gye = grd%gys + grd%gym - 1 5255a74a43SLisandro Dalcin grd%gxe = grd%gxs + grd%gxm - 1 5355a74a43SLisandro Dalcin 5455a74a43SLisandro Dalcin end subroutine GetGridInfo 5555a74a43SLisandro Dalcin 5655a74a43SLisandro Dalcin subroutine InitGuessLocal(grd, x, lambda, ierr) 5755a74a43SLisandro Dalcin implicit none 5855a74a43SLisandro Dalcin type(gridinfo) grd 5955a74a43SLisandro Dalcin PetscScalar x(grd%xs:grd%xe, grd%ys:grd%ye) 6055a74a43SLisandro Dalcin PetscReal lambda 6155a74a43SLisandro Dalcin PetscErrorCode ierr 6255a74a43SLisandro Dalcin ! 6355a74a43SLisandro Dalcin PetscInt i, j 6455a74a43SLisandro Dalcin PetscReal hx, hy, temp, temp1, one 6555a74a43SLisandro Dalcin 6655a74a43SLisandro Dalcin one = 1.0 6755a74a43SLisandro Dalcin hx = one/(dble(grd%mx - 1)) 6855a74a43SLisandro Dalcin hy = one/(dble(grd%my - 1)) 6955a74a43SLisandro Dalcin temp1 = lambda/(lambda + one) 7055a74a43SLisandro Dalcin 7155a74a43SLisandro Dalcin do j = grd%ys, grd%ye 7255a74a43SLisandro Dalcin temp = dble(min(j - 1, grd%my - j))*hy 7355a74a43SLisandro Dalcin do i = grd%xs, grd%xe 7455a74a43SLisandro Dalcin if (i == 1 .or. j == 1 .or. i == grd%mx .or. j == grd%my) then 7555a74a43SLisandro Dalcin ! boundary points 7655a74a43SLisandro Dalcin x(i, j) = 0.0 7755a74a43SLisandro Dalcin else 7855a74a43SLisandro Dalcin ! interior grid points 7955a74a43SLisandro Dalcin x(i, j) = temp1*sqrt(min(dble(min(i - 1, grd%mx - i)*hx), dble(temp))) 8055a74a43SLisandro Dalcin end if 8155a74a43SLisandro Dalcin end do 8255a74a43SLisandro Dalcin end do 8355a74a43SLisandro Dalcin ierr = 0 8455a74a43SLisandro Dalcin 8555a74a43SLisandro Dalcin end subroutine InitGuessLocal 8655a74a43SLisandro Dalcin 8755a74a43SLisandro Dalcin subroutine FunctionLocal(grd, x, f, lambda, ierr) 8855a74a43SLisandro Dalcin implicit none 8955a74a43SLisandro Dalcin type(gridinfo) grd 9055a74a43SLisandro Dalcin PetscScalar x(grd%gxs:grd%gxe, grd%gys:grd%gye) 9155a74a43SLisandro Dalcin PetscScalar f(grd%xs:grd%xe, grd%ys:grd%ye) 9255a74a43SLisandro Dalcin PetscReal lambda 9355a74a43SLisandro Dalcin PetscErrorCode ierr 9455a74a43SLisandro Dalcin ! 9555a74a43SLisandro Dalcin PetscInt i, j 9655a74a43SLisandro Dalcin PetscReal hx, hy, hxdhy, hydhx, sc, one, two 9755a74a43SLisandro Dalcin PetscScalar u, uxx, uyy 9855a74a43SLisandro Dalcin 9955a74a43SLisandro Dalcin one = 1.0 10055a74a43SLisandro Dalcin two = 2.0 10155a74a43SLisandro Dalcin hx = one/dble(grd%mx - 1) 10255a74a43SLisandro Dalcin hy = one/dble(grd%my - 1) 10355a74a43SLisandro Dalcin sc = hx*hy 10455a74a43SLisandro Dalcin hxdhy = hx/hy 10555a74a43SLisandro Dalcin hydhx = hy/hx 10655a74a43SLisandro Dalcin 10755a74a43SLisandro Dalcin do j = grd%ys, grd%ye 10855a74a43SLisandro Dalcin do i = grd%xs, grd%xe 10955a74a43SLisandro Dalcin if (i == 1 .or. j == 1 .or. i == grd%mx .or. j == grd%my) then 11055a74a43SLisandro Dalcin ! boundary points 11155a74a43SLisandro Dalcin f(i, j) = x(i, j) - 0.0 11255a74a43SLisandro Dalcin else 11355a74a43SLisandro Dalcin ! interior grid points 11455a74a43SLisandro Dalcin u = x(i, j) 11555a74a43SLisandro Dalcin uxx = (two*u - x(i - 1, j) - x(i + 1, j))*hydhx 11655a74a43SLisandro Dalcin uyy = (two*u - x(i, j - 1) - x(i, j + 1))*hxdhy 11755a74a43SLisandro Dalcin f(i, j) = uxx + uyy - lambda*exp(u)*sc 11855a74a43SLisandro Dalcin end if 11955a74a43SLisandro Dalcin end do 12055a74a43SLisandro Dalcin end do 12155a74a43SLisandro Dalcin ierr = 0 12255a74a43SLisandro Dalcin 12355a74a43SLisandro Dalcin end subroutine FunctionLocal 12455a74a43SLisandro Dalcin 12555a74a43SLisandro Dalcin subroutine JacobianLocal(grd, x, Jac, lambda, ierr) 12655a74a43SLisandro Dalcin implicit none 12755a74a43SLisandro Dalcin type(gridinfo) grd 12855a74a43SLisandro Dalcin PetscScalar x(grd%gxs:grd%gxe, grd%gys:grd%gye) 12955a74a43SLisandro Dalcin Mat Jac 13055a74a43SLisandro Dalcin PetscReal lambda 13155a74a43SLisandro Dalcin PetscErrorCode ierr 13255a74a43SLisandro Dalcin ! 13355a74a43SLisandro Dalcin PetscInt i, j, row(1), col(5) 13455a74a43SLisandro Dalcin PetscInt ione, ifive 13555a74a43SLisandro Dalcin PetscReal hx, hy, hxdhy, hydhx, sc, v(5), one, two 13655a74a43SLisandro Dalcin 13755a74a43SLisandro Dalcin ione = 1 13855a74a43SLisandro Dalcin ifive = 5 13955a74a43SLisandro Dalcin one = 1.0 14055a74a43SLisandro Dalcin two = 2.0 14155a74a43SLisandro Dalcin hx = one/dble(grd%mx - 1) 14255a74a43SLisandro Dalcin hy = one/dble(grd%my - 1) 14355a74a43SLisandro Dalcin sc = hx*hy 14455a74a43SLisandro Dalcin hxdhy = hx/hy 14555a74a43SLisandro Dalcin hydhx = hy/hx 14655a74a43SLisandro Dalcin 14755a74a43SLisandro Dalcin do j = grd%ys, grd%ye 14855a74a43SLisandro Dalcin row = (j - grd%gys)*grd%gxm + grd%xs - grd%gxs - 1 14955a74a43SLisandro Dalcin do i = grd%xs, grd%xe 15055a74a43SLisandro Dalcin row = row + 1 15155a74a43SLisandro Dalcin if (i == 1 .or. j == 1 .or. i == grd%mx .or. j == grd%my) then 15255a74a43SLisandro Dalcin ! boundary points 15355a74a43SLisandro Dalcin col(1) = row(1) 15455a74a43SLisandro Dalcin v(1) = one 15555a74a43SLisandro Dalcin PetscCall(MatSetValuesLocal(Jac, ione, row, ione, col, v, INSERT_VALUES, ierr)) 15655a74a43SLisandro Dalcin else 15755a74a43SLisandro Dalcin ! interior grid points 15855a74a43SLisandro Dalcin v(1) = -hxdhy 15955a74a43SLisandro Dalcin v(2) = -hydhx 16055a74a43SLisandro Dalcin v(3) = two*(hydhx + hxdhy) - lambda*exp(x(i, j))*sc 16155a74a43SLisandro Dalcin v(4) = -hydhx 16255a74a43SLisandro Dalcin v(5) = -hxdhy 16355a74a43SLisandro Dalcin col(1) = row(1) - grd%gxm 16455a74a43SLisandro Dalcin col(2) = row(1) - 1 16555a74a43SLisandro Dalcin col(3) = row(1) 16655a74a43SLisandro Dalcin col(4) = row(1) + 1 16755a74a43SLisandro Dalcin col(5) = row(1) + grd%gxm 16855a74a43SLisandro Dalcin PetscCall(MatSetValuesLocal(Jac, ione, row, ifive, col, v, INSERT_VALUES, ierr)) 16955a74a43SLisandro Dalcin end if 17055a74a43SLisandro Dalcin end do 17155a74a43SLisandro Dalcin end do 17255a74a43SLisandro Dalcin 17355a74a43SLisandro Dalcin end subroutine JacobianLocal 17455a74a43SLisandro Dalcin 17555a74a43SLisandro Dalcinend module Bratu2D 17655a74a43SLisandro Dalcin 17755a74a43SLisandro Dalcin! -------------------------------------------------------------------- 17855a74a43SLisandro Dalcin 17955a74a43SLisandro Dalcinsubroutine FormInitGuess(da, X, lambda, ierr) 18055a74a43SLisandro Dalcin use Bratu2D 18155a74a43SLisandro Dalcin implicit none 18255a74a43SLisandro Dalcin DM da 18355a74a43SLisandro Dalcin Vec X 18455a74a43SLisandro Dalcin PetscReal lambda 18555a74a43SLisandro Dalcin PetscErrorCode ierr 18655a74a43SLisandro Dalcin ! 18755a74a43SLisandro Dalcin type(gridinfo) :: grd 18855a74a43SLisandro Dalcin PetscScalar, pointer :: xx(:) 18955a74a43SLisandro Dalcin 190*ce78bad3SBarry Smith PetscCall(VecGetArray(X, xx, ierr)) 19155a74a43SLisandro Dalcin PetscCall(GetGridInfo(da, grd, ierr)) 19255a74a43SLisandro Dalcin PetscCall(InitGuessLocal(grd, xx, lambda, ierr)) 193*ce78bad3SBarry Smith PetscCall(VecRestoreArray(X, xx, ierr)) 19455a74a43SLisandro Dalcin 19555a74a43SLisandro Dalcinend subroutine FormInitGuess 19655a74a43SLisandro Dalcin 19755a74a43SLisandro Dalcinsubroutine FormFunction(da, X, F, lambda, ierr) 19855a74a43SLisandro Dalcin use Bratu2D 19955a74a43SLisandro Dalcin implicit none 20055a74a43SLisandro Dalcin DM da 20155a74a43SLisandro Dalcin Vec X 20255a74a43SLisandro Dalcin Vec F 20355a74a43SLisandro Dalcin PetscReal lambda 20455a74a43SLisandro Dalcin PetscErrorCode ierr 20555a74a43SLisandro Dalcin ! 20655a74a43SLisandro Dalcin type(gridinfo) :: grd 20755a74a43SLisandro Dalcin Vec :: localX 20855a74a43SLisandro Dalcin PetscScalar, pointer :: xx(:) 20955a74a43SLisandro Dalcin PetscScalar, pointer :: ff(:) 21055a74a43SLisandro Dalcin 21155a74a43SLisandro Dalcin PetscCall(DMGetLocalVector(da, localX, ierr)) 21255a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX, ierr)) 21355a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX, ierr)) 21455a74a43SLisandro Dalcin 215*ce78bad3SBarry Smith PetscCall(VecGetArray(localX, xx, ierr)) 216*ce78bad3SBarry Smith PetscCall(VecGetArray(F, ff, ierr)) 21755a74a43SLisandro Dalcin 21855a74a43SLisandro Dalcin PetscCall(GetGridInfo(da, grd, ierr)) 21955a74a43SLisandro Dalcin PetscCall(FunctionLocal(grd, xx, ff, lambda, ierr)) 22055a74a43SLisandro Dalcin 221*ce78bad3SBarry Smith PetscCall(VecRestoreArray(F, ff, ierr)) 222*ce78bad3SBarry Smith PetscCall(VecRestoreArray(localX, xx, ierr)) 22355a74a43SLisandro Dalcin PetscCall(DMRestoreLocalVector(da, localX, ierr)) 22455a74a43SLisandro Dalcin 22555a74a43SLisandro Dalcinend subroutine FormFunction 22655a74a43SLisandro Dalcin 22755a74a43SLisandro Dalcinsubroutine FormJacobian(da, X, J, lambda, ierr) 22855a74a43SLisandro Dalcin use Bratu2D 22955a74a43SLisandro Dalcin implicit none 23055a74a43SLisandro Dalcin DM da 23155a74a43SLisandro Dalcin Vec X 23255a74a43SLisandro Dalcin Mat J 23355a74a43SLisandro Dalcin PetscReal lambda 23455a74a43SLisandro Dalcin PetscErrorCode ierr 23555a74a43SLisandro Dalcin ! 23655a74a43SLisandro Dalcin type(gridinfo) :: grd 23755a74a43SLisandro Dalcin Vec :: localX 23855a74a43SLisandro Dalcin PetscScalar, pointer :: xx(:) 23955a74a43SLisandro Dalcin 24055a74a43SLisandro Dalcin PetscCall(DMGetLocalVector(da, localX, ierr)) 24155a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX, ierr)) 24255a74a43SLisandro Dalcin PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX, ierr)) 243*ce78bad3SBarry Smith PetscCall(VecGetArray(localX, xx, ierr)) 24455a74a43SLisandro Dalcin 24555a74a43SLisandro Dalcin PetscCall(GetGridInfo(da, grd, ierr)) 24655a74a43SLisandro Dalcin PetscCall(JacobianLocal(grd, xx, J, lambda, ierr)) 24755a74a43SLisandro Dalcin 248*ce78bad3SBarry Smith PetscCall(VecRestoreArray(localX, xx, ierr)) 24955a74a43SLisandro Dalcin PetscCall(DMRestoreLocalVector(da, localX, ierr)) 25055a74a43SLisandro Dalcin 25155a74a43SLisandro Dalcin PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY, ierr)) 25255a74a43SLisandro Dalcin PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY, ierr)) 25355a74a43SLisandro Dalcin 25455a74a43SLisandro Dalcinend subroutine FormJacobian 25555a74a43SLisandro Dalcin 25655a74a43SLisandro Dalcin! -------------------------------------------------------------------- 25755a74a43SLisandro Dalcin 25855a74a43SLisandro Dalcin! Local Variables: 25955a74a43SLisandro Dalcin! mode: f90 26055a74a43SLisandro Dalcin! End: 261