Uses Picard iterations to apply the inverse Jacobian of a system to a vector, to low accuracy. Rather than directly computing the inverse Jacobian, it is more efficient to approximate it. If is the vector being preconditioned, and is the result of applying the preconditioner, then Thus, the preconditioner can be applied by approximately solving this system for . Linearising , this system can be solved efficiently using Picard iteration.
Say that the Jacobian is an system of blocks of the sort implemented in the jacobian_block type (each labeled as ) and that the vector being preconditioned constists of scalar fields (). Then estimate of the solution for the field in the preconditioned vector () is the solution to Depending on the type of fields being used and the direction in which derivatives are being taken, may be tridiaganol, meaning it can be solved efficiently.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=r8), | private | :: | tolerance | = | 1.e-3_r8 | ||
integer, | private | :: | max_iterations | = | 20 |
Create a preconditioner object with the desired tolerance and maximum number of iterations.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=r8), | intent(in), | optional | :: | tolerance | The tolerance within which to apply the inverse Jacobian. Defaults to 0.001. |
|
integer, | intent(in), | optional | :: | max_iterations | The maximum number of iterations to use when applying the preconditioner. Defaults to 20. |
Use Picard iteration to approximately multiply the state vector by the inverse Jacobian. The details for this procedure are in the documentation of the preconditioner type.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(preconditioner), | intent(in) | :: | this | |||
class(jacobian_block), | intent(inout), | dimension(:,:) | :: | jacobian | An matrix approximating the Jacobian for which the preconditioner is used. |
|
class(scalar_field), | intent(in), | dimension(:) | :: | vector | A vector of size which is to be preconditioned. |
|
class(scalar_field), | intent(inout), | dimension(:) | :: | estimate | On entry, an initial guess for the preconditioned vector. On exit, the iteratively determined value of the preconditioned vector. |
type, public :: preconditioner
!* Author: Chris MacMackin
! Date: December 2016
!
! Uses Picard iterations to apply the inverse Jacobian of a system
! to a vector, to low accuracy. Rather than directly computing the
! inverse Jacobian, it is more efficient to approximate it. If
! \(d\) is the vector being preconditioned, and \(z\) is the
! result of applying the preconditioner, then $$ z = J^{-1}d
! \Rightarrow Jz = d.$$ Thus, the preconditioner can be applied by
! approximately solving this system for \(z\). Linearising \(J\),
! this system can be solved efficiently using Picard iteration.
!
! Say that the Jacobian is an \(n\times n\) system of blocks of
! the sort implemented in the [[jacobian_block]] type (each
! labeled as \(J_{j,k}\)) and that the vector being preconditioned
! constists of \(n\) scalar fields (\(d_j\)). Then \(m^{th}\)
! estimate of the solution for the \(j^{th}\) field in the
! preconditioned vector (\(z^m_j\)) is the solution to $$
! J_{j,j}z^m_j = d_j - \sum_{\substack{k=1\ k\ne j}}^n
! J_{j,k}z^{m-1}_k. $$ Depending on the type of fields being used
! and the direction in which derivatives are being taken,
! \(J_{j,j}\) may be tridiaganol, meaning it can be solved
! efficiently. <!--Otherwise, it can be approximated that $$
! J_{j,j}z^m_j = \left(\frac{\partial F}{\partial x_i} +
! F\Delta_i\right)z^{m} \simeq \frac{\partial F}{\partial
! x_i}z^m_j + F\frac{\partial z^{m-1}_j}{\partial x_i}. $$ The
! first term in this approximation corresponds to a diagonal
! matrix (which can be solved trivially), while the second term is
! known and can be subtracted from the right-hand-side of the
! linear system.-->
!
private
real(r8) :: tolerance = 1.e-3_r8
integer :: max_iterations = 20
contains
procedure :: apply => preconditioner_apply
end type preconditioner