ddx 0.6.0
Libary for domain-decomposition methods for polarizable continuum models
ddx_lpb_core Module Reference

Core routines and parameters specific to ddLPB. More...

Functions/Subroutines

subroutine wghpot_f (params, constants, workspace, gradphi, f, ddx_error)
 Find intermediate F0 in the RHS of the ddLPB model given in Eq.(82) More...
 
subroutine calcv2_lpb (params, constants, isph, pot, x)
 Intermediate computation of Bx. More...
 
subroutine convert_ddcosmo (params, constants, direction, vector)
 Scale the ddCOSMO solution vector. More...
 
subroutine inthsp (params, constants, rijn, isph, basloc, fac_hsp, work)
 Intermediate calculation in calcv2_lpb subroutine. More...
 
subroutine adjrhs_lpb (params, constants, isph, xi, vlm, basloc, vplm, vcos, vsin, tmp_bessel)
 Intermediate for computing Bstarx Compute the Adjoint matix B*x Taken from ddx_core routine adjrhs. More...
 
subroutine inthsp_adj (params, constants, rjin, jsph, basloc, fac_hsp, work)
 Intermediate calculation in adjrhs_lpb subroutine. More...
 

Detailed Description

Core routines and parameters specific to ddLPB.

ddX software

Version
1.0.0
Author
Abhinav Jha and Michele Nottoli
Date
2022-03-31

Function/Subroutine Documentation

◆ wghpot_f()

subroutine ddx_lpb_core::wghpot_f ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
type(ddx_workspace_type), intent(inout)  workspace,
real(dp), dimension(3, constants % ncav), intent(in)  gradphi,
real(dp), dimension(params % ngrid, params % nsph), intent(out)  f,
type(ddx_error_type), intent(inout)  ddx_error 
)

Find intermediate F0 in the RHS of the ddLPB model given in Eq.(82)

Parameters
[in]params: Input parameter file
[in]constants: Input constants file
[in]workspace: Input workspace
[in]gradphi: Gradient of psi_0
[out]f: Intermediate calculation of F0
[in,out]ddx_errorddX error

Definition at line 27 of file ddx_lpb_core.f90.

◆ calcv2_lpb()

subroutine ddx_lpb_core::calcv2_lpb ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  isph,
real(dp), dimension(params % ngrid), intent(inout)  pot,
real(dp), dimension(constants % nbasis, params % nsph), intent(in)  x 
)

Intermediate computation of Bx.

Parameters
[in]params: Input parameter file
[in]constants: Input constants file
[in]isph: Number of the sphere
[in,out]pot: Array of size ngrid
[in]x: Input vector (Usually X_e)

Definition at line 154 of file ddx_lpb_core.f90.

◆ convert_ddcosmo()

subroutine ddx_lpb_core::convert_ddcosmo ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  direction,
real(dp), dimension(constants % nbasis, params % nsph), intent(inout)  vector 
)

Scale the ddCOSMO solution vector.

Parameters
[in]params: Input parameter file
[in]constants: Input constants file
[in]direction: Direction of the scaling
[in,out]vector: ddCOSMO solution vector

Definition at line 204 of file ddx_lpb_core.f90.

◆ inthsp()

subroutine ddx_lpb_core::inthsp ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  rijn,
integer, intent(in)  isph,
real(dp), dimension(constants % nbasis), intent(in)  basloc,
real(dp), dimension(constants % nbasis), intent(inout)  fac_hsp,
complex(dp), dimension(max(2, params % lmax+1)), intent(out)  work 
)

Intermediate calculation in calcv2_lpb subroutine.

Parameters
[in]params: Input parameter file
[in]constants: Input constants file
[in]rijn: Radius of sphers x_ijn
[in]ri: Radius of sphers x_i
[in]isph: Index of sphere
[in]basloc: Spherical Harmonic
[out]fac_hsp: Return bessel function ratio multiplied by the spherical harmonic Y_l'm'. Array of size nylm

Definition at line 239 of file ddx_lpb_core.f90.

◆ adjrhs_lpb()

subroutine ddx_lpb_core::adjrhs_lpb ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
integer, intent(in)  isph,
real(dp), dimension(params % ngrid, params % nsph), intent(in)  xi,
real(dp), dimension((params % lmax+1)**2), intent(out)  vlm,
real(dp), dimension((params % lmax+1)**2), intent(out)  basloc,
real(dp), dimension((params % lmax+1)**2), intent(out)  vplm,
real(dp), dimension(params % lmax+1), intent(out)  vcos,
real(dp), dimension(params % lmax+1), intent(out)  vsin,
complex(dp), dimension(max(2, params % lmax+1)), intent(out)  tmp_bessel 
)

Intermediate for computing Bstarx Compute the Adjoint matix B*x Taken from ddx_core routine adjrhs.

Parameters
[in]params: Input parameter file
[in]constants: Input constants file
[in]isph: Number of the sphere
[in,out]basloc: Used to compute spherical harmonic
[in,out]vplm: Used to compute spherical harmonic
[in,out]vcos: Used to compute spherical harmonic
[in,out]vsin: Used to compute spherical harmonic

Definition at line 276 of file ddx_lpb_core.f90.

◆ inthsp_adj()

subroutine ddx_lpb_core::inthsp_adj ( type(ddx_params_type), intent(in)  params,
type(ddx_constants_type), intent(in)  constants,
real(dp), intent(in)  rjin,
integer, intent(in)  jsph,
real(dp), dimension(constants % nbasis), intent(in)  basloc,
real(dp), dimension(constants % nbasis), intent(inout)  fac_hsp,
complex(dp), dimension(max(2, params % lmax+1)), intent(out)  work 
)

Intermediate calculation in adjrhs_lpb subroutine.

Parameters
[in]params: Input parameter file
[in]constants: Input constants file
[in]rjin: Radius of sphers x_jin
[in]jsph: Index of sphere
[in]basloc: Spherical Harmonic
[out]fac_hsp: Return bessel function ratio multiplied by the spherical harmonic Y_l'm'. Array of size nylm

Definition at line 352 of file ddx_lpb_core.f90.