[table of contents] [master index] [comments] [modules] [programs] [variables] [types] [procedures]
PURPOSE
Integrator for the first-order perturbation theory for q(r,s).
Propagates both q(r,s) and delq(r,s)
SOURCE
module response_step_mod use const_mod use chemistry_mod use fft_mod private public :: response_step_startup public :: fft_step public :: pspropagate public :: alloc_exparrays public :: calc_exparrays public :: calc_exp_grid public :: ps_propagate
SUBROUTINE
response_step_startup(ngrid,order)
PURPOSE
Allocate the module variable exponential arrays
ARGUMENTS
integer ngrid - number of grid points
integer order - order of Richardson extrapolation
SOURCE
subroutine response_step_startup(Ngrid, order) integer,intent(IN) :: ngrid(3) integer,intent(IN) :: order
SUBROUTINE
ps_propagate(plan,planc,ds,delqin,q0,exp_gpdotr,order,delta_omega, &
monomer,p_mon,delqout,dsloc,qout)
PURPOSE
Propagate q and delq by solving the inhomog pde.
call fft_step and fftc_step for unperturbed and perturbation
partition functions.
ARGUMENTS
fft_plan planc - see fft_mod for details
complex delqin - delta q at beginning of the block
complex delqout - Output delta q array
complex exp_gpdotr - exponential array of G.r for perturbation at G+k
complex f_inhomo - = d_w * q_0(1:2**order)
integer p_mon - index of the perturbed monomer
integer monomer - index of the monomer of the concerned chain site
over which the propagation is being implemented
integer order - order of extrapolation
real dsloc - contour step size for this block
SOURCE
subroutine ps_propagate(delqin, delqout, p_mon, exp_gpdotr, &
f_inhomo, monomer, planc, order, dsloc)
complex(long), intent(IN) :: delqin(:,:,:)
complex(long), intent(OUT) :: delqout(:,:,:)
complex(long), intent(IN) :: exp_gpdotr(:,:,:)
real(long), intent(IN) :: f_inhomo(:,:,:,:)
integer, intent(IN) :: p_mon
integer, intent(IN) :: monomer
type(fft_plan) :: planc
real(long), intent(IN) :: dsloc
integer, intent(IN) :: order
SUBROUTINE
fftc_step(plan,order,q_in,q_out)
PURPOSE
solve for perturbed partition function of one step
call FFT and inverse FFT
ARGUMENTS
fft_plan plan - see fft_mod for details
integer order - extrapolation order
real q_in - input partition
real q_out - output partition
SOURCE
subroutine fftc_step(plan,order,q_in,q_out) implicit none type(fft_plan), intent(IN) :: plan integer, intent(IN) :: order complex(long), intent(IN) :: q_in(0:,0:,0:) complex(long), intent(OUT) :: q_out(0:,0:,0:)
SUBROUTINE
calc_exp_grid(ksq_grid, ksq_grid_pert, omega_grid,
monomer, ds, order, dsarray, pertb)
PURPOSE
Calculate the module variable exponential arrays
ARGUMENTS
real ksq_grid - |G|^2 for reciprocal lattice vectors
real omega_grid - omega fields at grid points
real ksq_grid_pert - |k+G|^2 on grids
integer monomer - monomer type for block
real ds - chain step
integer order - order of extrapolation
real dsarray - output array that has stepsizes as elements
logical pertb - calculate the perturbed exponential if true
SOURCE
subroutine calc_exp_grid(ksq_grid, ksq_grid_pert, omega_grid, &
monomer, ds, order, dsarray, pertb)
implicit none
real(long), intent(IN) :: ksq_grid(:,:,:)
real(long), intent(IN) :: ksq_grid_pert(:,:,:)
real(long), intent(IN) :: omega_grid(:,:,:,:)
integer, intent(IN) :: monomer
real(long), intent(IN) :: ds
integer, intent(IN) :: order
real(long), intent(out) :: dsarray(:)
logical, intent(IN) :: pertb