double_gram_schmidt_step Interface

public interface double_gram_schmidt_step

Description

Given an array of abstract_vector and an abstract_vector (or array of abstract_vectors) , this subroutine returns a modified vector orthogonal to all columns of , i.e.

using a double Gram-Schmidt process. On exit, is orthogonal to albeit does not have unit-norm. Note moreover that is assumed to be an orthonormal set of vectors. The function can also return the projection coefficients .

Syntax

    call double_gram_schmidt_step(y, X, info [, if_chk_orthonormal] [, beta])

Arguments

y : abstract_vector (or array of abstract_vector) that needs to be orthogonalize in-place against .

X : Array of abstract_vector against which needs to be orthogonalized. Note the function assumes that is an orthonormal set of vectors, i.e. . If it this is not the case, the result are meaningless.

info : integer Information flag.

if_chk_orthonormal (optional) : logical flag (default .true.) to check whether is an orthonormal set of vectors or not. If the orthonormality returns .false., the function throws an error. Note that this check is however computationally expensive and can be disable for the sake of performances.

beta (optional) : real or complex array containing the coefficients .


Module Procedures

private subroutine DGS_vector_against_basis_rsp(y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_rsp), intent(inout) :: y

Input abstract_vector to orthogonalize

class(abstract_vector_rsp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
real(kind=sp), intent(out), optional :: beta(:)

Projection coefficients if requested

private subroutine DGS_basis_against_basis_rsp(Y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_rsp), intent(inout) :: Y(:)

Input abstract_vector basis to orthogonalize

class(abstract_vector_rsp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
real(kind=sp), intent(out), optional :: beta(:,:)

Projection coefficients if requested

private subroutine DGS_vector_against_basis_rdp(y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_rdp), intent(inout) :: y

Input abstract_vector to orthogonalize

class(abstract_vector_rdp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
real(kind=dp), intent(out), optional :: beta(:)

Projection coefficients if requested

private subroutine DGS_basis_against_basis_rdp(Y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_rdp), intent(inout) :: Y(:)

Input abstract_vector basis to orthogonalize

class(abstract_vector_rdp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
real(kind=dp), intent(out), optional :: beta(:,:)

Projection coefficients if requested

private subroutine DGS_vector_against_basis_csp(y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_csp), intent(inout) :: y

Input abstract_vector to orthogonalize

class(abstract_vector_csp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
complex(kind=sp), intent(out), optional :: beta(:)

Projection coefficients if requested

private subroutine DGS_basis_against_basis_csp(Y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_csp), intent(inout) :: Y(:)

Input abstract_vector basis to orthogonalize

class(abstract_vector_csp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
complex(kind=sp), intent(out), optional :: beta(:,:)

Projection coefficients if requested

private subroutine DGS_vector_against_basis_cdp(y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_cdp), intent(inout) :: y

Input abstract_vector to orthogonalize

class(abstract_vector_cdp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
complex(kind=dp), intent(out), optional :: beta(:)

Projection coefficients if requested

private subroutine DGS_basis_against_basis_cdp(Y, X, info, if_chk_orthonormal, beta)

Computes one step of the double Gram-Schmidt orthogonalization process of the abstract_vector y against the abstract_vector basis X

Arguments

Type IntentOptional Attributes Name
class(abstract_vector_cdp), intent(inout) :: Y(:)

Input abstract_vector basis to orthogonalize

class(abstract_vector_cdp), intent(in) :: X(:)

Input abstract_vector basis to orthogonalize against

integer, intent(out) :: info

Information flag.

logical, intent(in), optional :: if_chk_orthonormal
complex(kind=dp), intent(out), optional :: beta(:,:)

Projection coefficients if requested