FMS  2024.03
Flexible Modeling System
tridiagonal_mod

Solves a tridiagonal system of equations. More...

Data Types

interface  tri_invert
 Interface to solve tridiagonal systems of equations for either kind value. Module level variables will be deallocated and allocated for every Since this relies on the state of module variables (unless A,B,C are specified) the values stored are distinct for each kind call unless the added optional argument store_both_kinds is true. More...
 
type  tridiag_reals_r8
 

Functions/Subroutines

subroutine close_tridiagonal
 Releases memory used by the solver.
 
subroutine tri_invert_ (x, d, a, b, c, store_both_kinds)
 Sets up and solves the tridiagonal system of equations. More...
 
 tri_invert_r4
 
 tri_invert_r8
 

Variables

real(r4_kind), dimension(:,:), allocatable, private bb
 
real(r8_kind), dimension(:,:), allocatable, private bb
 
real(r4_kind), dimension(:,:,:), allocatable, private cc
 
real(r8_kind), dimension(:,:,:), allocatable, private cc
 
real(r4_kind), dimension(:,:,:), allocatable, private e
 
real(r8_kind), dimension(:,:,:), allocatable, private e
 
real(r4_kind), dimension(:,:,:), allocatable, private g
 
real(r8_kind), dimension(:,:,:), allocatable, private g
 
logical, private init_tridiagonal_r4 = .false.
 true when fields in tridiag_r4 are allocated
 
logical, private init_tridiagonal_r8 = .false.
 true when fields in tridiag_r8 are allocated
 
type(tridiag_reals_r4tridiag_r4
 holds reals stored from r4_kind calls to tri_invert
 
type(tridiag_reals_r8tridiag_r8
 holds reals stored from r8_kind calls to tri_invert
 

Detailed Description

Solves a tridiagonal system of equations.

The following schematic represents the system of equations solved, where X is the solution.

     | B(1)  A(1)   0     0                .......            0    |  |X(1)|   |D(1)|
     | C(2)  B(2)  A(2)   0                .......            0    |  |X(2)|   |D(2)|
     |  0    C(3)  B(3)  A(3)  0           .......            0    |  | .. |   | .. |
     |  ..........................................                 |  | .. | = | .. |
     |  ..........................................                 |  | .. |   | .. |
     |                                  C(N-2) B(N-2) A(N-2)  0    |  | .. |   | .. |
     |                                    0    C(N-1) B(N-1) A(N-1)|  | .. |   | .. |
     |                                    0      0    C(N)   B(N)  |  |X(N)|   |D(N)|

 

To solve this system

   call tri_invert(X,D,A,B,C)

       real, intent(out), dimension(:,:,:) :: X
       real, intent(in),  dimension(:,:,:) :: D
       real, optional,    dimension(:,:,:) :: A,B,C
 

For simplicity (?), A and C are assumed to be dimensioned the same size as B, D, and X, although any input values for A(N) and C(1) are ignored. (some checks are needed here)

If A is not present, it is assumed that the matrix (A,B.C) has not been changed since the last call to tri_invert.

To release memory,

    call close_tridiagonal
 

Arguments A, B, and C are optional, and are saved as module variables if one recalls tri_invert without changing (A,B,C)

Note
Optional arguments A,B,C have no intent declaration, so the default intent is inout. The value of A(N) is modified on output, and B and C are unchanged.

The following private allocatable arrays save the relevant information if one recalls tri_invert without changing (A,B,C):

        allocate ( e  (size(x,1), size(x,2), size(x,3)) )
        allocate ( g  (size(x,1), size(x,2), size(x,3)) )
        allocate ( cc (size(x,1), size(x,2), size(x,3)) )
        allocate ( bb (size(x,1), size(x,2)) )
 

This storage is deallocated when close_tridiagonal is called.


Data Type Documentation

◆ tridiagonal_mod::tri_invert

interface tridiagonal_mod::tri_invert

Interface to solve tridiagonal systems of equations for either kind value. Module level variables will be deallocated and allocated for every Since this relies on the state of module variables (unless A,B,C are specified) the values stored are distinct for each kind call unless the added optional argument store_both_kinds is true.

Definition at line 104 of file tridiagonal.F90.

Public Member Functions

 tri_invert_r4
 
 tri_invert_r8
 

◆ tridiagonal_mod::tridiag_reals_r8

type tridiagonal_mod::tridiag_reals_r8

Definition at line 88 of file tridiagonal.F90.

Collaboration diagram for tridiag_reals_r8:
[legend]

Private Attributes

real(r8_kind), dimension(:,:), allocatable, private bb
 
real(r8_kind), dimension(:,:,:), allocatable, private cc
 
real(r8_kind), dimension(:,:,:), allocatable, private e
 
real(r8_kind), dimension(:,:,:), allocatable, private g
 

Function/Subroutine Documentation

◆ tri_invert_()

subroutine tri_invert_ ( real(fms_trid_kind_), dimension(:,:,:), intent(out)  x,
real(fms_trid_kind_), dimension(:,:,:), intent(in)  d,
real(fms_trid_kind_), dimension(:,:,:), optional  a,
real(fms_trid_kind_), dimension(:,:,:), optional  b,
real(fms_trid_kind_), dimension(:,:,:), optional  c,
logical, optional  store_both_kinds 
)

Sets up and solves the tridiagonal system of equations.

For simplicity, A and C are assumed to be dimensioned the same size as B, D, and X, although any input values for A(N) and C(1) are ignored. There are no checks to make sure the sizes agree.

The value of A(N) is modified on output, and B and C are unchanged.

For mixed precision, this routine uses the kind size macro(FMS_TRID_KIND_) to determine which module variables are used/stored. This means a,b, and c values will only be stored for calls of the same real kind value unless store_both_kinds is present and .true..

Parameters
[out]xSolution to the tridiagonal system of equations
[in]dThe right-hand side term, see the schematic above.
cLeft hand side terms(see schematic on module page). If not provided, values from last call are used
store_both_kindsWill save module state variables for both kind types in order to be used in subsequent calls with either kind.

Definition at line 34 of file tridiagonal.inc.