The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines
oft_native_solvers Module Reference

Detailed Description

Abstract solver interfaces and select native implementations.

Abstract interface definitions

Native solver implementations

Preconditioner implementations

See also
oft_cg, oft_gmres, oft_petsc_solvers
Authors
Chris Hansen
Date
August 2011

Data Types

type  oft_bjprecond
 Block-Jacobi preconditioner. More...
 
type  oft_diag_cscale
 Diagonal preconditioner. More...
 
type  oft_diag_scale
 Diagonal preconditioner. More...
 
type  oft_identity_inv
 Identity matrix inversion. More...
 
type  oft_jblock_precond
 Symmetric Jacobi smoother. More...
 
type  oft_ml_precond
 Multi-level preconditioner level. More...
 
type  oft_ml_trans
 Multi-level transfer level. More...
 
type  oft_native_cg_eigsolver
 CG eigensolver class. More...
 
type  oft_native_cg_solver
 CG solver class. More...
 
type  oft_native_gmres_csolver
 GMRES solver class (complex) More...
 
type  oft_native_gmres_solver
 GMRES solver class. More...
 
type  oft_nksolver
 Native Newton solver. More...
 
interface  oft_update_jacobian
 Needs docs. More...
 

Functions/Subroutines

recursive subroutine bjprecond_apply (self, u, g)
 Precondition a linear system using a Block-Jacobi method.
 
subroutine bjprecond_delete (self)
 Destroy Block-Jacobi preconditioner and deallocate all internal storage.
 
subroutine bjprecond_setup_xml (self, solver_node, level)
 Setup block-Jacobi smoother from XML definition.
 
recursive subroutine bjprecond_update (self, new_pattern)
 Update solver after changing settings/operators.
 
subroutine cdiag_scale_apply (self, u, g)
 Solve the linear system for a diagonal matrix (complex)
 
subroutine cdiag_scale_delete (self)
 Destroy diagonal preconditioner and deallocate all internal storage.
 
subroutine cg_delete (self)
 Destroy diagonal preconditioner and deallocate all internal storage.
 
subroutine cg_eig_delete (self)
 Destroy diagonal preconditioner and deallocate all internal storage.
 
subroutine cg_eigsolver_apply (self, u, alam)
 Solve a general eigenvalue system using the Conjugate-Gradient method.
 
subroutine cg_setup_xml (self, solver_node, level)
 Setup CG solver from XML definition.
 
recursive subroutine cg_solver_apply (self, u, g)
 Solve a linear system using the Conjugate-Gradient method.
 
subroutine cgmres_delete (self)
 Destroy diagonal preconditioner and deallocate all internal storage.
 
subroutine cgmres_setup_xml (self, solver_node, level)
 Setup GMRES solver from XML definition.
 
recursive subroutine cgmres_solver_apply (self, u, g)
 Solve a linear system using the flexible GMRES method (complex)
 
subroutine diag_scale_apply (self, u, g)
 Solve the linear system for a diagonal matrix.
 
logical function, public diag_scale_cast (self, source)
 Cast a solver object to a oft_diag_scale.
 
subroutine diag_scale_delete (self)
 Destroy diagonal preconditioner and deallocate all internal storage.
 
subroutine gmres_delete (self)
 Destroy diagonal preconditioner and deallocate all internal storage.
 
subroutine gmres_setup_xml (self, solver_node, level)
 Setup GMRES solver from XML definition.
 
recursive subroutine gmres_solver_apply (self, u, g)
 Solve a linear system using the flexible GMRES method.
 
subroutine identity_inv_apply (self, u, g)
 Solve the trivial linear system for the identity matrix.
 
subroutine identity_inv_delete (self)
 Destroy diagonal preconditioner and deallocate all internal storage.
 
subroutine jblock_precond_apply (self, u, g)
 Apply 1-step of a symmetric Jacobi smoother with native CRS matrices.
 
logical function, public jblock_precond_cast (self, source)
 Cast a solver object to a jblock_precond_cast.
 
subroutine jblock_precond_delete (self)
 Destroy symmetric Jacobi preconditioner and deallocate all internal storage.
 
subroutine jblock_setup_xml (self, solver_node, level)
 Setup symmetric Jacobi smoother from XML definition.
 
recursive subroutine ml_precond_apply (self, u, g)
 Apply 1-step of a multi-level preconditioner.
 
logical function, public ml_precond_cast (self, source)
 Cast a solver object to a oft_ml_precond.
 
recursive subroutine ml_precond_delete (self)
 Destroy Multi-Level preconditioner and deallocate all internal storage.
 
recursive subroutine ml_precond_update (self, new_pattern)
 Update solver after changing settings/operators.
 
recursive subroutine ml_precond_view (self)
 Print solver configuration.
 
recursive subroutine ml_trans_apply (self, u, g)
 Transfer solution between distributed and shared levels as part of a ML preconditioner.
 
logical function, public ml_trans_cast (self, source)
 Cast a solver object to a oft_ml_trans.
 
subroutine ml_trans_delete (self)
 Destroy Multi-Level preconditioner and deallocate all internal storage.
 
logical function native_cg_eigsolver_cast (self, source)
 Cast an eigensolver object to a oft_native_cg_eigsolver.
 
logical function, public native_cg_solver_cast (self, source)
 Cast a solver object to a oft_native_cg_solver.
 
logical function, public native_gmres_solver_cast (self, source)
 Cast a solver object to a oft_native_gmres_solver.
 
subroutine nksolver_apply (self, u, g)
 Solve a nonlinear system Newton's method.
 
subroutine nksolver_delete (self)
 Destroy Newton solver and deallocate all internal storage.
 

Function/Subroutine Documentation

◆ bjprecond_apply()

recursive subroutine bjprecond_apply ( class(oft_bjprecond), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Precondition a linear system using a Block-Jacobi method.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ bjprecond_delete()

subroutine bjprecond_delete ( class(oft_bjprecond), intent(inout)  self)
private

Destroy Block-Jacobi preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ bjprecond_setup_xml()

subroutine bjprecond_setup_xml ( class(oft_bjprecond), intent(inout)  self,
type(xml_node), intent(in), pointer  solver_node,
integer(i4), intent(in), optional  level 
)
private

Setup block-Jacobi smoother from XML definition.

Parameters
[in,out]selfSolver object
[in]solver_nodeXML element containing solver definition
[in]levelLevel in MG hierarchy (optional)

◆ bjprecond_update()

recursive subroutine bjprecond_update ( class(oft_bjprecond), intent(inout)  self,
logical, intent(in), optional  new_pattern 
)
private

Update solver after changing settings/operators.

Parameters
[in,out]selfSolver object
[in]new_patternUpdate matrix non-zero pattern (optional)

◆ cdiag_scale_apply()

subroutine cdiag_scale_apply ( class(oft_diag_cscale), intent(inout)  self,
class(oft_cvector), intent(inout)  u,
class(oft_cvector), intent(inout)  g 
)
private

Solve the linear system for a diagonal matrix (complex)

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ cdiag_scale_delete()

subroutine cdiag_scale_delete ( class(oft_diag_cscale), intent(inout)  self)
private

Destroy diagonal preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ cg_delete()

subroutine cg_delete ( class(oft_native_cg_solver), intent(inout)  self)
private

Destroy diagonal preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ cg_eig_delete()

subroutine cg_eig_delete ( class(oft_native_cg_eigsolver), intent(inout)  self)
private

Destroy diagonal preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ cg_eigsolver_apply()

subroutine cg_eigsolver_apply ( class(oft_native_cg_eigsolver), intent(inout)  self,
class(oft_vector), intent(inout)  u,
real(r8), intent(inout)  alam 
)
private

Solve a general eigenvalue system using the Conjugate-Gradient method.

This solver uses a Non-Linear Conjugate-Gradient method to minimize the Rayleigh Quotient ( \( R = \frac{x*A*x}{x*M*x} \)).

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]alamEigenvalue

◆ cg_setup_xml()

subroutine cg_setup_xml ( class(oft_native_cg_solver), intent(inout)  self,
type(xml_node), intent(in), pointer  solver_node,
integer(i4), intent(in), optional  level 
)
private

Setup CG solver from XML definition.

Parameters
[in,out]selfSolver object
[in]solver_nodeXML element containing solver definition
[in]levelLevel in MG hierarchy (optional)

◆ cg_solver_apply()

recursive subroutine cg_solver_apply ( class(oft_native_cg_solver), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Solve a linear system using the Conjugate-Gradient method.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ cgmres_delete()

subroutine cgmres_delete ( class(oft_native_gmres_csolver), intent(inout)  self)
private

Destroy diagonal preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ cgmres_setup_xml()

subroutine cgmres_setup_xml ( class(oft_native_gmres_csolver), intent(inout)  self,
type(xml_node), intent(in), pointer  solver_node,
integer(i4), intent(in), optional  level 
)
private

Setup GMRES solver from XML definition.

Parameters
[in,out]selfSolver object
[in]solver_nodeXML element containing solver definition
[in]levelLevel in MG hierarchy (optional)

◆ cgmres_solver_apply()

recursive subroutine cgmres_solver_apply ( class(oft_native_gmres_csolver), intent(inout)  self,
class(oft_cvector), intent(inout)  u,
class(oft_cvector), intent(inout)  g 
)
private

Solve a linear system using the flexible GMRES method (complex)

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ diag_scale_apply()

subroutine diag_scale_apply ( class(oft_diag_scale), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Solve the linear system for a diagonal matrix.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ diag_scale_cast()

logical function, public diag_scale_cast ( type(oft_diag_scale), intent(out), pointer  self,
class(oft_solver), intent(in), target  source 
)

Cast a solver object to a oft_diag_scale.

The source matrix must be oft_diag_scale or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ diag_scale_delete()

subroutine diag_scale_delete ( class(oft_diag_scale), intent(inout)  self)
private

Destroy diagonal preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ gmres_delete()

subroutine gmres_delete ( class(oft_native_gmres_solver), intent(inout)  self)
private

Destroy diagonal preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ gmres_setup_xml()

subroutine gmres_setup_xml ( class(oft_native_gmres_solver), intent(inout)  self,
type(xml_node), intent(in), pointer  solver_node,
integer(i4), intent(in), optional  level 
)
private

Setup GMRES solver from XML definition.

Parameters
[in,out]selfSolver object
[in]solver_nodeXML element containing solver definition
[in]levelLevel in MG hierarchy (optional)

◆ gmres_solver_apply()

recursive subroutine gmres_solver_apply ( class(oft_native_gmres_solver), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Solve a linear system using the flexible GMRES method.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ identity_inv_apply()

subroutine identity_inv_apply ( class(oft_identity_inv), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Solve the trivial linear system for the identity matrix.

Note
Typically used for simple cases when a solver wrapper is required but the trivial case of \( u = g \) is the desired result
Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ identity_inv_delete()

subroutine identity_inv_delete ( class(oft_identity_inv), intent(inout)  self)
private

Destroy diagonal preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ jblock_precond_apply()

subroutine jblock_precond_apply ( class(oft_jblock_precond), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Apply 1-step of a symmetric Jacobi smoother with native CRS matrices.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ jblock_precond_cast()

logical function, public jblock_precond_cast ( type(oft_jblock_precond), intent(out), pointer  self,
class(oft_solver), intent(in), target  source 
)

Cast a solver object to a jblock_precond_cast.

The source matrix must be jblock_precond_cast or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ jblock_precond_delete()

subroutine jblock_precond_delete ( class(oft_jblock_precond), intent(inout)  self)
private

Destroy symmetric Jacobi preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ jblock_setup_xml()

subroutine jblock_setup_xml ( class(oft_jblock_precond), intent(inout)  self,
type(xml_node), intent(in), pointer  solver_node,
integer(i4), intent(in), optional  level 
)
private

Setup symmetric Jacobi smoother from XML definition.

Parameters
[in,out]selfSolver object
[in]solver_nodeXML element containing solver definition
[in]levelLevel in MG hierarchy (optional)

◆ ml_precond_apply()

recursive subroutine ml_precond_apply ( class(oft_ml_precond), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Apply 1-step of a multi-level preconditioner.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ ml_precond_cast()

logical function, public ml_precond_cast ( class(oft_ml_precond), intent(out), pointer  self,
class(oft_solver), intent(in), target  source 
)

Cast a solver object to a oft_ml_precond.

The source matrix must be oft_ml_precond or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ ml_precond_delete()

recursive subroutine ml_precond_delete ( class(oft_ml_precond), intent(inout)  self)
private

Destroy Multi-Level preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ ml_precond_update()

recursive subroutine ml_precond_update ( class(oft_ml_precond), intent(inout)  self,
logical, intent(in), optional  new_pattern 
)
private

Update solver after changing settings/operators.

Parameters
[in,out]selfSolver object

◆ ml_precond_view()

recursive subroutine ml_precond_view ( class(oft_ml_precond), intent(inout)  self)
private

Print solver configuration.

Parameters
[in,out]selfSolver object

◆ ml_trans_apply()

recursive subroutine ml_trans_apply ( class(oft_ml_trans), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Transfer solution between distributed and shared levels as part of a ML preconditioner.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ ml_trans_cast()

logical function, public ml_trans_cast ( class(oft_ml_trans), intent(out), pointer  self,
class(oft_solver), intent(in), target  source 
)

Cast a solver object to a oft_ml_trans.

The source matrix must be oft_ml_trans or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ ml_trans_delete()

subroutine ml_trans_delete ( class(oft_ml_trans), intent(inout)  self)
private

Destroy Multi-Level preconditioner and deallocate all internal storage.

Parameters
[in,out]selfSolver object

◆ native_cg_eigsolver_cast()

logical function native_cg_eigsolver_cast ( type(oft_native_cg_eigsolver), intent(out), pointer  self,
class(oft_eigsolver), intent(in), target  source 
)
private

Cast an eigensolver object to a oft_native_cg_eigsolver.

The source matrix must be oft_native_cg_eigsolver or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ native_cg_solver_cast()

logical function, public native_cg_solver_cast ( type(oft_native_cg_solver), intent(out), pointer  self,
class(oft_solver), intent(in), target  source 
)

Cast a solver object to a oft_native_cg_solver.

The source matrix must be oft_native_cg_solver or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ native_gmres_solver_cast()

logical function, public native_gmres_solver_cast ( type(oft_native_gmres_solver), intent(out), pointer  self,
class(oft_solver), intent(in), target  source 
)

Cast a solver object to a oft_native_gmres_solver.

The source matrix must be oft_native_gmres_solver or a child class, otherwise pointer will be returned as null and success == .FALSE.

Parameters
[out]selfReference to source object with desired class
[in]sourceSource solver to cast
Returns
Cast success flag

◆ nksolver_apply()

subroutine nksolver_apply ( class(oft_nksolver), intent(inout)  self,
class(oft_vector), intent(inout)  u,
class(oft_vector), intent(inout)  g 
)
private

Solve a nonlinear system Newton's method.

Parameters
[in,out]selfSolver object
[in,out]uGuess (input), Solution (output)
[in,out]gRHS (input), Residual (output)

◆ nksolver_delete()

subroutine nksolver_delete ( class(oft_nksolver), intent(inout)  self)
private

Destroy Newton solver and deallocate all internal storage.

Parameters
[in,out]selfSolver object