PROGRESS  master
prg_quantumdynamics_mod.F90
Go to the documentation of this file.
1 
10  use bml
11  implicit none
12 
13  private
14  integer, parameter :: dp = kind(1.0d0)
15 
16 
21 
22 contains
23 
41  subroutine prg_kick_density(kick_direc,kick_mag,dens,norbs,mdim,S,SINV,&
42  which_atom,r,bmltype,thresh)
43  integer :: i
44  integer, intent(in) :: kick_direc, norbs, mdim
45  integer, allocatable, intent(in) :: which_atom(:)
46  real(dp), allocatable :: r(:,:)
47  real(dp) :: kick_mag, thresh
48  complex(dp), allocatable, intent(inout) :: dens(:,:)
49  complex(dp), allocatable :: tmat1(:), tmat2(:), s(:,:), sinv(:,:)
50  complex(dp) :: telem
51  type(bml_matrix_t) :: t1, t2, rho_bml, s_bml, sinv_bml
52  character(len=*), intent(in) :: bmltype
53 
54  allocate (tmat1(norbs))
55  allocate (tmat2(norbs))
56 
57  do i=1,norbs
58  telem = cmplx(0.0_dp,-kick_mag*r(which_atom(i),kick_direc))
59  tmat1(i) = exp(telem)
60  tmat2(i) = exp(-telem)
61  enddo
62 
63  call bml_zero_matrix(bmltype, bml_element_complex,dp,norbs,mdim, t1)
64  call bml_set_diagonal(t1,tmat1,thresh)
65  deallocate (tmat1)
66  call bml_import_from_dense(bmltype,dens,rho_bml,thresh,mdim)
67  call bml_zero_matrix(bmltype, bml_element_complex,dp,norbs,mdim, t2)
68  call bml_multiply(t1,rho_bml,t2)
69  call bml_import_from_dense(bmltype,s,s_bml,thresh,mdim)
70  call bml_multiply(t2,s_bml,rho_bml)
71  call bml_deallocate(s_bml)
72 
73  call bml_set_diagonal(t1,tmat2,thresh)
74  deallocate (tmat2)
75  call bml_multiply(rho_bml,t1,t2)
76  call bml_deallocate(t1)
77  call bml_import_from_dense(bmltype,sinv,sinv_bml,thresh,mdim)
78  call bml_multiply(t2,sinv_bml,rho_bml)
79  call bml_deallocate(sinv_bml)
80  call bml_deallocate(t2)
81 
82  call bml_export_to_dense(rho_bml,dens)
83 
84  end subroutine prg_kick_density
85 
86 
97  subroutine prg_get_sparsity_cplxmat(matrix_type,element_type,thresh,a_dense)
98  character(len=*), intent(in) :: matrix_type ,element_type
99  complex(dp),intent(in) :: a_dense(:,:)
100  type(bml_matrix_t) :: a
101  real(dp),intent(in) :: thresh
102  integer :: asize
103  character(len=20) :: test_type
104 
105  asize=size(a_dense,1)
106  call bml_import_from_dense(matrix_type, a_dense, a, thresh,asize)
107  call bml_deallocate(a)
108 
109  end subroutine prg_get_sparsity_cplxmat
110 
121  subroutine prg_get_sparsity_realmat(matrix_type,element_type,thresh,a_dense)
122  character(len=*), intent(in) :: matrix_type ,element_type
123  real(dp),intent(in) :: a_dense(:,:)
124  type(bml_matrix_t) :: a
125  real(dp),intent(in) :: thresh
126  integer :: asize
127 
128  asize=size(a_dense,1)
129  call bml_import_from_dense(matrix_type,a_dense,a,thresh,asize)
130  call bml_deallocate(a)
131 
132  end subroutine prg_get_sparsity_realmat
133 
134 
152  subroutine prg_kick_density_bml(kick_direc,kick_mag,rho_bml,s_bml,sinv_bml,mdim,which_atom,&
153  r,matrix_type,thresh)
154  integer :: i, norbs, mdim
155  integer, intent(in) :: kick_direc
156  integer, allocatable, intent(in) :: which_atom(:)
157  real(dp), allocatable :: r(:,:)
158  real(dp) :: kick_mag, thresh
159  complex(dp), allocatable :: tmat1(:), tmat2(:)
160  complex(dp) :: telem
161  type(bml_matrix_t) :: t1, t2, rho_bml,s_bml,sinv_bml
162  character(len=*), intent(in) :: matrix_type
163 
164  norbs=bml_get_n(rho_bml)
165  allocate (tmat1(norbs))
166  allocate (tmat2(norbs))
167 
168  do i=1,norbs
169  telem = cmplx(0.0_dp,-kick_mag*r(which_atom(i),kick_direc))
170  tmat1(i) = exp(telem)
171  tmat2(i) = exp(-telem)
172  enddo
173 
174  call bml_zero_matrix(matrix_type, bml_element_complex,dp,norbs,mdim, t1)
175  call bml_zero_matrix(matrix_type, bml_element_complex,dp,norbs,mdim, t2)
176  call bml_set_diagonal(t1,tmat1,thresh)
177  deallocate (tmat1)
178 
179  call bml_multiply(t1,rho_bml,t2)
180  call bml_multiply(t2,s_bml,rho_bml)
181  call bml_set_diagonal(t1,tmat2,thresh)
182  deallocate (tmat2)
183  call bml_multiply(rho_bml,t1,t2)
184  call bml_deallocate(t1)
185 
186  call bml_multiply(t2,sinv_bml,rho_bml)
187 
188  call bml_deallocate(t2)
189 
190  end subroutine prg_kick_density_bml
191 
192 
210  subroutine prg_lvni_bml(h1_bml, sinv_bml, dt, hbar, rhoold_bml, rho_bml,aux_bml,matrix_type,mdim,thresh)
211  integer :: norbs, mdim
212  real(dp) :: hbar, dt
213  type(bml_matrix_t) :: h1_bml, sinv_bml
214  type(bml_matrix_t) :: rho_bml, rhoold_bml, aux_bml,aux_1
215  complex(dp) :: dr
216  real(dp), intent(in) :: thresh
217  character(len=*), intent(in) :: matrix_type
218 
219  norbs=bml_get_n(rho_bml)
220  dr = 2.0_dp*dt*cmplx(0.0_dp,-1.0_dp/hbar)
221  call bml_multiply(h1_bml,rho_bml,aux_bml,1.0_dp,0.0_dp,thresh)
222  call bml_zero_matrix(matrix_type,bml_element_complex,dp,norbs,mdim,aux_1)
223  call bml_multiply(sinv_bml,aux_bml,aux_1,1.0_dp,0.0_dp,thresh)
224  call bml_multiply(rho_bml,h1_bml,aux_bml,1.0_dp,0.0_dp,thresh)
225  call bml_multiply(aux_bml,sinv_bml,aux_1,-1.0_dp,1.0d0)
226  call bml_scale(dr,aux_1,aux_bml)
227  call bml_deallocate(aux_1)
228  call bml_add(aux_bml,rhoold_bml,1.0_dp,1.0_dp,thresh)
229  call bml_copy(rho_bml,rhoold_bml)
230  call bml_copy(aux_bml,rho_bml)
231 
232  end subroutine prg_lvni_bml
233 
234 
246  subroutine prg_getcharge(rho_bml,s_bml,charges,aux_bml,z,spindex,N,nats,thresh)
247  integer :: i, j, k, nats,norbs
248  integer, allocatable,intent(in) :: n(:), spindex(:)
249  complex(dp),allocatable :: auxd(:,:)
250  real(dp), allocatable :: charges(:)
251  real(dp), intent(in) :: thresh,z(:)
252  type(bml_matrix_t), intent(in) :: s_bml, rho_bml
253  type(bml_matrix_t) :: aux_bml
254 
255  norbs=bml_get_n(rho_bml)
256  allocate(auxd(norbs,norbs))
257  auxd=cmplx(0.0_dp,0.0_dp)
258  if(.not.allocated(charges)) allocate(charges(nats))
259  call bml_multiply(rho_bml,s_bml,aux_bml,1.0_dp,0.0_dp,thresh)
260  call bml_export_to_dense(aux_bml,auxd)
261  k=0
262  do i = 1,nats
263  charges(i)=0.0_dp
264  do j = 1, n(spindex(i))
265  k = k+1
266  charges(i) = charges(i) + 2.0_dp*auxd(k,k)
267  enddo
268  charges(i) = -charges(i) + z(spindex(i))
269  enddo
270  deallocate(auxd)
271 
272  end subroutine prg_getcharge
273 
274 
281  subroutine prg_getdipole(charges,r,mu)
282  integer :: i, norbs
283  real(dp), intent(in) :: charges(:), r(:,:)
284  real(dp), intent(inout) :: mu(3)
285 
286  norbs = size(charges,1)
287  mu = 0.0_dp
288 
289  do i=1,norbs
290  mu=mu+r(:,i)*charges(i)
291  enddo
292 
293  end subroutine prg_getdipole
294 
295 
297  !simulate
298  !! photo-excitation.
299  !! This routine does:
300  !! Manually induces electronic excitation based on the initial filling matrix
301  !! by promoting an electron.
302  !! \param fill_mat the initial filling matrix
303  !! \param orbit_orig the origin orbital with an electron to be promoted
304  !! \param orbit_exci the destination orbital of the promoted electron
305  !!
306  subroutine prg_excitation(fill_mat, orbit_orig, orbit_exci)
307  integer, intent(inout) :: fill_mat(:)
308  integer, intent(in) :: orbit_orig, orbit_exci
309 
310  fill_mat(orbit_orig) = fill_mat(orbit_orig) - 1.0d0
311  fill_mat(orbit_exci) = fill_mat(orbit_exci) + 1.0d0
312 
313  end subroutine prg_excitation
314 
315 end module prg_quantumdynamics_mod
A module to add in common quantum dynamical operations.
subroutine, public prg_kick_density(kick_direc, kick_mag, dens, norbs, mdim, S, SINV, which_atom, r, bmltype, thresh)
Provides perturbation to initial density matrix in the form of an electric field kick....
subroutine, public prg_kick_density_bml(kick_direc, kick_mag, rho_bml, s_bml, sinv_bml, mdim, which_atom, r, matrix_type, thresh)
Provides perturbation to initial density matrix in the form of an electric field kick given input mat...
subroutine, public prg_get_sparsity_cplxmat(matrix_type, element_type, thresh, a_dense)
This computes the sparsity of a complex matrix given a threshold value This routine does: where is ...
subroutine, public prg_getcharge(rho_bml, s_bml, charges, aux_bml, z, spindex, N, nats, thresh)
Constructs the charges from the density matrix.
subroutine, public prg_excitation(fill_mat, orbit_orig, orbit_exci)
Produce an excitation in the initially calculated density matrix to.
subroutine, public prg_lvni_bml(h1_bml, sinv_bml, dt, hbar, rhoold_bml, rho_bml, aux_bml, matrix_type, mdim, thresh)
Performs Liouville-von Neumann integration using leap-frog method. This routine does: where the time...
subroutine, public prg_getdipole(charges, r, mu)
This routine computes the dipole moment of the system with units determined by the units of the coord...
subroutine, public prg_get_sparsity_realmat(matrix_type, element_type, thresh, a_dense)
This computes the sparsity of a real matrix given a threshold value This routine does: where is the...