FEMSolverClass Module



Interfaces

public interface reverseArray

  • public subroutine reverseArrayReal64(A)

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(inout) :: A(:,:)

public interface reverseVector

  • public subroutine reverseVectorReal64(A)

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(inout) :: A(:)

Derived Types

type, public ::  M_Link_Item_

Components

Type Visibility Attributes Name Initial
integer(kind=int32), public :: rank_and_rowID_1(1:2) = [0, 0]
integer(kind=int32), public :: rank_and_rowID_2(1:2) = [0, 0]

type, public ::  FEMSolver_

Components

Type Visibility Attributes Name Initial
type(FEMDomainp_), public, allocatable :: femdomains(:)
integer(kind=int32), public, allocatable :: DomainIDs(:)
real(kind=real64), public, allocatable :: IfaceElemConnectivity(:,:)
real(kind=real64), public, allocatable :: IfaceElemDomainID(:,:)
logical, public :: initialized = .false.
logical, public :: InterfaceExist = .false.
logical, public :: debug = .false.
real(kind=real64), public, allocatable :: CRS_val(:)
integer(kind=int32), public, allocatable :: CRS_Index_Col(:)
integer(kind=int64), public, allocatable :: CRS_Index_Row(:)
real(kind=real64), public, allocatable :: CRS_RHS(:)
integer(kind=int32), public, allocatable :: number_of_element_per_domain(:)
integer(kind=int32), public, allocatable :: number_of_point_per_domain(:)
real(kind=real64), public, allocatable :: A_CRS_val(:)

General Eigen Value Problem [A]{x} = (lambda)[B]{x}

integer(kind=int32), public, allocatable :: A_CRS_Index_Col(:)
integer(kind=int64), public, allocatable :: A_CRS_Index_Row(:)
logical, public :: A_empty = .true.
real(kind=real64), public, allocatable :: B_CRS_val(:)
integer(kind=int32), public, allocatable :: B_CRS_Index_Col(:)
integer(kind=int64), public, allocatable :: B_CRS_Index_Row(:)
logical, public :: B_empty = .true.
integer(kind=int32), public, allocatable :: fix_eig_IDs(:)
logical, public, allocatable :: fix_lin_exists(:)
real(kind=real64), public, allocatable :: fix_lin_exists_Values(:)
real(kind=real64), public, allocatable :: CRS_x(:)
integer(kind=int32), public, allocatable :: CRS_ID_Starts_From(:)
real(kind=real64), public, allocatable :: A_dense(:,:)
integer(kind=int32), public, allocatable :: Num_nodes_in_Domains(:)
type(MPI_), public, pointer :: MPI_target => null()
type(M_Link_Item_), public, allocatable :: Link_Table(:)
integer(kind=int32), public :: LINK_TABLE_INIT_SIZE = 1000
integer(kind=int32), public :: Link_num = 0
type(CRS_), public :: ILU_MATRIX
integer(kind=int32), public :: itrmax = 100000
real(kind=real64), public :: er0 = dble(1.0e-10)
real(kind=real64), public :: relative_er = dble(1.0e-10)
logical, public :: use_LOBPCG = .false.
integer(kind=int32), public :: LOBPCG_MAX_ITR = 100000
real(kind=real64), public :: LOBPCG_TOL = dble(1.0e-6)
integer(kind=int32), public :: LOBPCG_NUM_MODE = 5

! stack for CRS matrices

type(CRS_), public, allocatable :: CRS_matrix_stack(:)
type(Dictionary_), public :: CRS_matrix_stack_key
integer(kind=int32), public :: LAST_CRS_MATRIX_STACK = 0
integer(kind=int32), public :: MAXSIZE_CRS_MATRIX_STACK = 10

Type-Bound Procedures

procedure, public :: init => initFEMSolver
procedure, public :: setDomain => setDomainFEMSolver
procedure, public :: setDomains => setDomainFEMSolver
procedure, public, pass :: setCRS_CRSobjFEMSolver
procedure, public, pass :: setCRSFEMSolver
procedure, public :: setRHS => setRHSFEMSolver
generic, public :: setCRS => setCRSFEMSolver, setCRS_CRSobjFEMSolver
procedure, public :: getCRS => getCRSFEMSolver
procedure, public :: setValue => setValueFEMSolver
procedure, public :: setValues => setValueFEMSolver
procedure, public :: addMatrixValue => addMatrixValueFEMSolver
procedure, public :: setMatrix => setMatrixFEMSolver
procedure, public :: setVector => setVectorFEMSolver
procedure, public :: keepThisMatrixAs => keepThisMatrixAsFEMSolver
procedure, public :: fix => fixFEMSolver
procedure, public :: fix_eig => fix_eigFEMSolver
procedure, public, pass :: importCRSMatrix_FEMSolver
generic, public :: import => importCRSMatrix_FEMSolver
procedure, public :: saveMatrix => saveMatrixFEMSolver
procedure, public :: loadMatrix => loadMatrixFEMSolver
procedure, public :: eig => eigFEMSolver
procedure, public, pass :: solveFEMSolver
procedure, public, pass :: solveFEMSolver_UserDefinedLinearSolver
procedure, public, pass :: solveFEMSolver_UserDefinedLinearSolverAsFunc
generic, public :: solve => solveFEMSolver, solveFEMSolver_UserDefinedLinearSolver, solveFEMSolver_UserDefinedLinearSolverAsFunc
procedure, public :: conditionNumber => conditionNumberFEMSolver
procedure, public :: zeros => zerosFEMSolver
procedure, public :: matmulDiagMatrix => matmulDiagMatrixFEMSolver
procedure, public :: diag => diagFEMSolver
procedure, public :: setEbOM => setEbOMFEMSolver
procedure, public :: get_gap_function => get_gap_function_FEMSolver
procedure, public :: argumented_Lagrangian_RHS => argumented_Lagrangian_RHS_FEMSolver
procedure, public :: num_EbOFEM_projection_point => num_EbOFEM_projection_point_FEMSolver
procedure, public :: get_fix_idx => get_fix_idx_FEMSolver
procedure, public :: get_fix_value => get_fix_value_FEMSolver
procedure, public :: MPI_link => MPI_linkFEMSolver
procedure, public :: MPI_dot_product => MPI_dot_productFEMSolver
procedure, public :: MPI_matmul => MPI_matmulFEMSolver
procedure, public :: MPI_BICGSTAB => MPI_BiCGSTABFEMSolver
procedure, public :: remove => removeFEMSolver

Functions

public function diagFEMSolver(this, cell_centered) result(diag_vector)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this
logical, intent(in), optional :: cell_centered

Return Value real(kind=real64), allocatable, (:)

public function solveFEMSolver_UserDefinedLinearSolver(this, LinearSolver, x0) result(x)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
public subroutine LinearSolver(row_ptr, col_idx, val, rhs, x)
Arguments
Type IntentOptional Attributes Name
integer(kind=int64), intent(in) :: row_ptr(:)
integer(kind=int32), intent(in) :: col_idx(:)
real(kind=real64), intent(in) :: val(:)
real(kind=real64), intent(in) :: rhs(:)
real(kind=real64), intent(inout) :: x(:)
real(kind=real64), intent(in), optional :: x0(:)

Return Value real(kind=real64), allocatable, (:)

public function solveFEMSolver_UserDefinedLinearSolverAsFunc(this, LinearSolver, x0) result(x)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
public function LinearSolver(row_ptr, col_idx, val, rhs, x0) result(x)
Arguments
Type IntentOptional Attributes Name
integer(kind=int64), intent(in) :: row_ptr(:)
integer(kind=int32), intent(in) :: col_idx(:)
real(kind=real64), intent(in) :: val(:)
real(kind=real64), intent(in) :: rhs(:)
real(kind=real64), intent(in) :: x0(:)
Return Value real(kind=real64), allocatable, (:)
real(kind=real64), intent(in) :: x0(:)

Return Value real(kind=real64), allocatable, (:)

public function solveFEMSolver(this, algorithm, preconditioning, x0) result(x)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
character(len=*), intent(in), optional :: algorithm
character(len=*), intent(in), optional :: preconditioning
real(kind=real64), intent(in), optional :: x0(:)

Return Value real(kind=real64), allocatable, (:)

public function num_EbOFEM_projection_point_FEMSolver(this) result(ret)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this

Return Value integer(kind=int32)

public function get_gap_function_FEMSolver(this, DOF, X) result(gap_function)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this
integer(kind=int32), intent(in) :: DOF
real(kind=real64), intent(in) :: X(:)

Return Value real(kind=real64), allocatable, (:)

public function argumented_Lagrangian_RHS_FEMSolver(this, DOF, lambda) result(aL_RHS)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this
integer(kind=int32), intent(in) :: DOF
real(kind=real64), intent(in) :: lambda(:)

Return Value real(kind=real64), allocatable, (:)

public function getCRSval(val, row_ptr, col_idx, row, col) result(ret)

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in) :: val(:)
integer(kind=int64), intent(in) :: row_ptr(:)
integer(kind=int32), intent(in) :: col_idx(:)
integer(kind=int32), intent(in) :: row
integer(kind=int32), intent(in) :: col

Return Value real(kind=real64)

public function sortByIDreal64ColisVector(vectors, ID) result(new_vectors)

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in) :: vectors(:,:)
integer(kind=int32), intent(in) :: ID(:)

Return Value real(kind=real64), allocatable, (:,:)

public pure function getCRSFEMSolver(this, name) result(ret)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this
character(len=*), intent(in), optional :: name

Return Value type(CRS_)

public function MPI_dot_productFEMSolver(this, a, b) result(dp)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this
real(kind=real64), intent(in) :: a(:)
real(kind=real64), intent(in) :: b(:)

Return Value real(kind=real64)

public function MPI_matmulFEMSolver(this, A, b) result(my_c)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this
type(CRS_), intent(in), optional :: A
real(kind=real64), intent(in) :: b(:)

Return Value real(kind=real64), allocatable, (:)

public function conditionNumberFEMSolver(this) result(RCOND)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this

Return Value real(kind=real64)

public function get_fix_idx_FEMSolver(this) result(fix_idx)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this

Return Value integer(kind=int32), allocatable, (:)

public function get_fix_value_FEMSolver(this) result(fix_value)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this

Return Value integer(kind=int32), allocatable, (:)


Subroutines

public recursive subroutine initFEMSolver(this, numDomain, FEMDomains, DomainIDs, DOF, MPI_target, NUM_CRS_MATRIX_STACK)

General Eigen Value Problem [A]{x} = (lambda)[B]{x}

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in), optional :: numDomain
type(FEMDomain_), intent(in), optional :: FEMDomains(:)
integer(kind=int32), intent(in), optional :: DomainIDs(:)
integer(kind=int32), intent(in), optional :: DOF
type(MPI_), intent(in), optional, target :: MPI_target
integer(kind=int32), intent(in), optional :: NUM_CRS_MATRIX_STACK

public recursive subroutine setDomainFEMSolver(this, FEMDomain, FEMDomains, FEMDomainPointers, DomainID, DomainIDs)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
type(FEMDomain_), intent(in), optional, target :: FEMDomain
type(FEMDomain_), intent(in), optional, target :: FEMDomains(:)
type(FEMDomainp_), intent(in), optional, target :: FEMDomainPointers(:)
integer(kind=int32), intent(in), optional :: DomainID
integer(kind=int32), intent(in), optional :: DomainIDs(:)

public subroutine setCRS_CRSobjFEMSolver(this, CRS)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
type(CRS_), intent(in) :: CRS

public subroutine setRHSFEMSolver(this, RHS)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
real(kind=real64), intent(in) :: RHS(:)

public recursive subroutine setCRSFEMSolver(this, DOF, debug)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in) :: DOF
logical, intent(in), optional :: debug

public subroutine setMatrixFEMSolver(this, DomainID, ElementID, DOF, Matrix, as_Dense)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in) :: DomainID
integer(kind=int32), intent(in) :: ElementID
integer(kind=int32), intent(in) :: DOF
real(kind=real64), intent(in) :: Matrix(:,:)
logical, intent(in), optional :: as_Dense

public subroutine setVectorFEMSolver(this, DomainID, ElementID, DOF, Vector)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in) :: DomainID
integer(kind=int32), intent(in) :: ElementID
integer(kind=int32), intent(in) :: DOF
real(kind=real64), intent(in) :: Vector(:)

public subroutine addMatrixValueFEMSolver(this, row_id, col_id, SingleValue, as_Dense)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in) :: row_id
integer(kind=int32), intent(in) :: col_id
real(kind=real64), intent(in) :: SingleValue
logical, intent(in), optional :: as_Dense

public subroutine setValueFEMSolver(this, DomainID, ElementID, DOF, Matrix, Vector, as_Dense)

$OMP atomic $OMP atomic

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in) :: DomainID
integer(kind=int32), intent(in) :: ElementID
integer(kind=int32), intent(in) :: DOF
real(kind=real64), intent(in), optional :: Matrix(:,:)
real(kind=real64), intent(in), optional :: Vector(:)
logical, intent(in), optional :: as_Dense

public subroutine fixFEMSolver(this, DomainID, IDs, FixValue, FixValues)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in), optional :: DomainID
integer(kind=int32), intent(in) :: IDs(:)
real(kind=real64), intent(in), optional :: FixValue
real(kind=real64), intent(in), optional :: FixValues(:)

public subroutine matmulDiagMatrixFEMSolver(this, diagMat)

diag is diagonal component of n x n matrix

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
real(kind=real64), intent(in) :: diagMat(:)

public subroutine loadMatrixFEMSolver(this, from)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
character(len=1), intent(in) :: from

public subroutine saveMatrixFEMSolver(this, name, CRS_as_dense, if_dense_exists, zero_or_nonzero)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(in) :: this
character(len=*), intent(in) :: name
logical, intent(in), optional :: CRS_as_dense
logical, intent(in), optional :: if_dense_exists
logical, intent(in), optional :: zero_or_nonzero

public subroutine zerosFEMSolver(this)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this

public subroutine LanczosMethod(this, eigen_value, Eigen_vectors, max_itr)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
real(kind=real64), allocatable :: eigen_value(:)
real(kind=real64), allocatable :: Eigen_vectors(:,:)
integer(kind=int32), intent(in) :: max_itr

public recursive subroutine eigFEMSolver(this, num_eigen, eigen_value, eigen_vectors)

INPUT INPUT/OUTPUT INPUT OUTPUT

Read more…

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in), optional :: num_eigen
real(kind=real64), intent(inout), allocatable :: eigen_value(:)
real(kind=real64), intent(inout), allocatable :: eigen_vectors(:,:)

public subroutine keepThisMatrixAsFEMSolver(this, As)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
character(len=*), intent(in) :: As

public subroutine fix_eigFEMSolver(this, IDs)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in) :: IDs(:)

public subroutine removeFEMSolver(this, only_matrices)

General Eigen Value Problem [A]{x} = (lambda)[B]{x}

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
logical, intent(in), optional :: only_matrices

public subroutine setEbOMFEMSolver(this, penalty, DOF)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
real(kind=real64), intent(in) :: penalty
integer(kind=int32), intent(in) :: DOF

public subroutine bicgstab_CRS_2(a, ptr_i, index_j, x, b, itrmax, er, relative_er, debug)

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(inout) :: a(:)
integer(kind=int64), intent(inout) :: ptr_i(:)
integer(kind=int32), intent(inout) :: index_j(:)
real(kind=real64), intent(inout) :: x(:)
real(kind=real64), intent(inout) :: b(:)
integer(kind=int32), intent(inout) :: itrmax
real(kind=real64), intent(inout) :: er
real(kind=real64), intent(in), optional :: relative_er
logical, intent(in), optional :: debug

public subroutine MPI_BiCGSTABFEMSolver(this, x)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
real(kind=real64), intent(inout), allocatable :: x(:)

public subroutine gauss_jordan_crs(a0_val, a0_row_ptr, a0_col_idx, x, b, n)

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in) :: a0_val(n)
integer(kind=int64), intent(in) :: a0_row_ptr(:)
integer(kind=int32), intent(in) :: a0_col_idx(:)
real(kind=real64), intent(inout) :: x(n)
real(kind=real64), intent(in) :: b(n)
integer(kind=int32), intent(in) :: n

public subroutine JacobiPreconditionerCRS(val, row_ptr, col_idx, rhs)

$OMP atomic $OMP end do $OMP end parallel

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(inout) :: val(:)
integer(kind=int64), intent(inout) :: row_ptr(:)
integer(kind=int32), intent(inout) :: col_idx(:)
real(kind=real64), intent(inout) :: rhs(:)

public subroutine reverseVectorReal64(A)

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(inout) :: A(:)

public subroutine reverseArrayReal64(A)

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(inout) :: A(:,:)

public subroutine MPI_linkFEMSolver(this, rank_and_rowID_1, rank_and_rowID_2)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
integer(kind=int32), intent(in) :: rank_and_rowID_1(1:2)
integer(kind=int32), intent(in) :: rank_and_rowID_2(1:2)

public subroutine importCRSMatrix_FEMSolver(this, CRS)

Arguments

Type IntentOptional Attributes Name
class(FEMSolver_), intent(inout) :: this
type(CRS_), intent(in) :: CRS