Stellarator-Tools
blocktridiagonalsolver_s Module Reference

Module contains subroutines for solver for block tri-diagonal matrices. More...

Functions/Subroutines

subroutine, public initialize (inN, inM)
 To be invoked, before solve, by user. More...
 
subroutine, public getranks (irank, inranks)
 
subroutine, public setmatrixrowcoll (globrow, Lj, j)
 To be invoked, before solve, after initialize, by user, to set up the matrix. More...
 
subroutine, public setmatrixrowcold (globrow, Dj, j)
 To be invoked, before solve, after initialize, by user, to set up the matrix. More...
 
subroutine, public setmatrixrowcolu (globrow, Uj, j)
 To be invoked, before solve, after initialize, by user, to set up the matrix. More...
 
subroutine, public setmatrixrhs (globrow, b)
 To be invoked, before solve, after initialize, by user, to set up the matrix. More...
 
subroutine, public getmatrixrowcoll (globrow, Lj, j)
 Can be invoked after initialize, by user, to get a copy of the matrix. More...
 
subroutine, public getmatrixrowcold (globrow, Dj, j)
 Can be invoked after initialize, by user, to get a copy of the matrix. More...
 
subroutine, public getmatrixrowcolu (globrow, Uj, j)
 Can be invoked after initialize, by user, to get a copy of the matrix. More...
 
subroutine, public getmatrixrhs (globrow, b)
 Can be invoked after initialize, by user, to get a copy of the RHS. More...
 
subroutine, public getsolutionvector (globrow, x)
 To be invoked, after solve, by user to get the solution vector of a local row. More...
 
subroutine, public finalize (do_mpifinalize)
 To be invoked, after solve, by user. More...
 
subroutine, public forwardsolve
 BCYCLIC forward phase; to be called after Initialize.
 
subroutine, public backwardsolve
 BCYCLIC backward phase; to be called after ForwardSolve.
 
subroutine, public checkconditionnumber (nblkrow, bsize, anorm, rcond, info)
 
subroutine, public checksymmetry (asymIndx)
 
subroutine, public storediagonal (blkrow, colnum, buf)
 
subroutine, public parmatvec (invec, outvec, outveclength)
 
subroutine, public applyparallelscaling (lmp, colscale, lverbose)
 
subroutine, public findminmax_tri (lmp)
 
subroutine, public refactorhessian (lmp)
 
subroutine, public padsides (arrin, blksize, top, bot)
 

Variables

real(dp), public maxeigen_tri
 maximum eigenvalue from Gerschgorin theorem
 
real(dp), public dmin_tri
 minimum non-zero diagonal element, used to scale Levenberg parameter (SPH 08-31-16)
 

Detailed Description

Module contains subroutines for solver for block tri-diagonal matrices.

Author
Kalyan S. Perumalla and Sudip K. Seal, ORNL
Date
2009-2017

Function/Subroutine Documentation

◆ finalize()

subroutine, public blocktridiagonalsolver_s::finalize ( logical, intent(in)  do_mpifinalize)

To be invoked, after solve, by user.

Parameters
[in]do_mpifinalizeInvoke MPI_Finalize?

◆ getmatrixrhs()

subroutine, public blocktridiagonalsolver_s::getmatrixrhs ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(out)  b 
)

Can be invoked after initialize, by user, to get a copy of the RHS.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[out]bRHS column corr. to globrow

◆ getmatrixrowcold()

subroutine, public blocktridiagonalsolver_s::getmatrixrowcold ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(out)  Dj,
integer, intent(in)  j 
)

Can be invoked after initialize, by user, to get a copy of the matrix.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[out]djj'th colum of D at globrow
[in]jcolumn number of L at globrow that is being set

◆ getmatrixrowcoll()

subroutine, public blocktridiagonalsolver_s::getmatrixrowcoll ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(out)  Lj,
integer, intent(in)  j 
)

Can be invoked after initialize, by user, to get a copy of the matrix.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[out]ljj'th colum of L at globrow; 1st L is always 0
[in]jcolumn number of L at globrow that is being set

◆ getmatrixrowcolu()

subroutine, public blocktridiagonalsolver_s::getmatrixrowcolu ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(out)  Uj,
integer, intent(in)  j 
)

Can be invoked after initialize, by user, to get a copy of the matrix.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[out]ujj'th colum of L at globrow; Nth U is always 0
[in]jcolumn number of U at globrow that is being set

◆ getsolutionvector()

subroutine, public blocktridiagonalsolver_s::getsolutionvector ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(out)  x 
)

To be invoked, after solve, by user to get the solution vector of a local row.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[out]xSolution column corr. to globrow

◆ initialize()

subroutine, public blocktridiagonalsolver_s::initialize ( integer, intent(in)  inN,
integer, intent(in)  inM 
)

To be invoked, before solve, by user.

Parameters
[in]innNum of row blocks in input block tri-diag matrix
[in]inmSize of each square sub-matrix block

◆ setmatrixrhs()

subroutine, public blocktridiagonalsolver_s::setmatrixrhs ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(in)  b 
)

To be invoked, before solve, after initialize, by user, to set up the matrix.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[in]bRHS column corr. to globrow

◆ setmatrixrowcold()

subroutine, public blocktridiagonalsolver_s::setmatrixrowcold ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(in)  Dj,
integer, intent(in)  j 
)

To be invoked, before solve, after initialize, by user, to set up the matrix.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[in]djj'th colum of D at globrow
[in]jcolumn number of L at globrow that is being set

◆ setmatrixrowcoll()

subroutine, public blocktridiagonalsolver_s::setmatrixrowcoll ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(in)  Lj,
integer, intent(in)  j 
)

To be invoked, before solve, after initialize, by user, to set up the matrix.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[in]ljj'th colum of L at globrow; 1st L is always 0
[in]jcolumn number of L at globrow that is being set

◆ setmatrixrowcolu()

subroutine, public blocktridiagonalsolver_s::setmatrixrowcolu ( integer, intent(in)  globrow,
real(dp), dimension(:), intent(in)  Uj,
integer, intent(in)  j 
)

To be invoked, before solve, after initialize, by user, to set up the matrix.

Parameters
[in]globrowOriginal/global block-row num in [1..N]
[in]ujj'th colum of U at globrow; Nth U is always 0
[in]jcolumn number of U at globrow that is being set