Octopus
em_resp_calc.F90
Go to the documentation of this file.
1!! Copyright (C) 2004-2012 Xavier Andrade, Eugene S. Kadantsev (ekadants@mjs1.phy.queensu.ca), David Strubbe
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use batch_oct_m
23 use comm_oct_m
24 use debug_oct_m
27 use elf_oct_m
28 use grid_oct_m
29 use global_oct_m
31 use ions_oct_m
32 use, intrinsic :: iso_fortran_env
38 use mesh_oct_m
48 use space_oct_m
55 use utils_oct_m
56 use xc_oct_m
58
59 implicit none
60
61 private
62 public :: &
70 dinhomog_b, &
71 zinhomog_b, &
90 freq2str, &
91 magn_dir, &
92 em_wfs_tag, &
93 em_rho_tag, &
98
99
100 type matrix_t
101 private
102 real(real64), allocatable :: dmatrix(:, :)
103 complex(real64), allocatable :: zmatrix(:, :)
104 end type matrix_t
105
106contains
107
108 ! ---------------------------------------------------------
109 subroutine lr_calc_current(st, space, gr, lr, lr_m)
110 type(states_elec_t), intent(inout) :: st
111 class(space_t), intent(in) :: space
112 type(grid_t), intent(in) :: gr
113 type(lr_t), intent(inout) :: lr
114 type(lr_t), optional, intent(inout) :: lr_m
115
116 integer :: idir, ist, ispin, idim, ndim, np
117
118 complex(real64), allocatable :: psi(:, :), gpsi(:,:), gdl_psi(:,:), gdl_psi_m(:,:)
119
120 push_sub(lr_calc_current)
121
122 if (.not. allocated(lr%dl_j)) then
123 safe_allocate(lr%dl_j(1:gr%np, 1:space%dim, 1:st%d%nspin))
124 end if
125
126 np = gr%np
127 ndim = space%dim
128
129 safe_allocate(psi(1:gr%np_part, 1:ndim))
130 safe_allocate(gpsi(1:np, 1:ndim))
131 safe_allocate(gdl_psi(1:np, 1:ndim))
132 if (present(lr_m)) then
133 safe_allocate(gdl_psi_m(1:np, 1:ndim))
134 end if
135
136 lr%dl_j = m_zero
137
138 do ispin = 1, st%d%nspin
139 do ist = 1, st%nst
140
141 call states_elec_set_state(st, gr, ist, ispin, psi)
142
143 do idim = 1, st%d%dim
144
145 call zderivatives_grad(gr%der, lr%zdl_psi(:, idim, ist, ispin), gdl_psi)
146 call zderivatives_grad(gr%der, psi(:, idim), gpsi)
147
148 if (present(lr_m)) then
149
150 call zderivatives_grad(gr%der, lr_m%zdl_psi(:, idim, ist, ispin), gdl_psi_m)
151
152 do idir = 1, space%dim
153
154 lr%dl_j(1:np, idir, ispin) = lr%dl_j(1:np, idir, ispin) + ( &
155 + conjg(psi(1:np, idim))*gdl_psi(1:np, idir) &
156 - psi(1:np, idim)*conjg(gdl_psi_m(1:np, idir)) &
157 + conjg(lr_m%zdl_psi(1:np, idim, ist, ispin)) * gpsi(1:np, idir) &
158 - lr%zdl_psi (1:np, idim, ist, ispin) * conjg(gpsi(1:np, idir)) &
159 ) / (m_two * m_zi)
160 end do
161
162 else
163
164 do idir = 1, space%dim
165
166 lr%dl_j(1:np, idir, ispin) = lr%dl_j(1:np, idir, ispin) + ( &
167 + conjg(psi(1:np, idim))*gdl_psi(1:np, idir) &
168 - psi(1:np, idim)*conjg(gdl_psi(1:np, idir)) &
169 + conjg(lr%zdl_psi(1:np, idim, ist, ispin)) * gpsi(1:np, idir) &
170 - lr%zdl_psi(1:np, idim, ist, ispin) * conjg(gpsi(1:np, idir)) &
171 ) / (m_two * m_zi)
172
173 end do
174
175 end if
176
177 end do
178 end do
179 end do
180
181 safe_deallocate_a(psi)
182 safe_deallocate_a(gpsi)
183 safe_deallocate_a(gdl_psi)
184 if (present(lr_m)) then
185 safe_deallocate_a(gdl_psi_m)
186 end if
187
188 pop_sub(lr_calc_current)
189
190 end subroutine lr_calc_current
191
192
193! ---------------------------------------------------------
194 character(len=12) function freq2str(freq) result(str)
195 real(real64), intent(in) :: freq
196
197 push_sub(freq2str)
199 ! some compilers (xlf) do not put a leading zero when the number
200 ! is smaller than 1. We have to check and correct that behavior.
201
202 write(str, '(f11.4)') freq
203 str = adjustl(str)
204 if (abs(freq) < m_one) then
205 if (freq >= m_zero .and. str(1:1) /= '0') str = "0"//trim(str)
206 if (freq < m_zero .and. str(2:2) /= '0') str = "-0"//trim(str(2:))
207 end if
208 str = trim(adjustl(str))
209
210 pop_sub(freq2str)
211
212 end function freq2str
213
214
215! ---------------------------------------------------------
216 character(len=100) function em_rho_tag(freq, dir, dir2, ipert) result(str)
217 real(real64), intent(in) :: freq
218 integer, intent(in) :: dir
219 integer, optional, intent(in) :: dir2
220 integer, optional, intent(in) :: ipert
221
222 character(len=12) :: str_tmp
223
224 !this function has to be consistent with oct_search_file_lr in liboct/oct_f.c
225
226 push_sub(em_rho_tag)
227
228 str_tmp = freq2str(freq)
229 write(str, '(3a,i1)') 'rho_', trim(str_tmp), '_', dir
230 if (present(dir2)) write(str, '(2a,i1)') trim(str), "_", dir2
231 if (present(ipert)) write(str, '(3a)') trim(str), "_", index2pert(ipert)
232
233 pop_sub(em_rho_tag)
234
235 end function em_rho_tag
236
237
238! ---------------------------------------------------------
239 character(len=100) function em_wfs_tag(idir, ifactor, idir2, ipert) result(str)
240 integer, intent(in) :: idir
241 integer, intent(in) :: ifactor
242 integer, optional, intent(in) :: idir2
243 integer, optional, intent(in) :: ipert
244
245 push_sub(em_wfs_tag)
246
247 write(str, '(3a,i1)') "wfs_", index2axis(idir), "_f", ifactor
248 if (present(idir2)) write(str, '(3a)') trim(str), "_", index2axis(idir2)
249 if (present(ipert)) write(str, '(3a)') trim(str), "_", index2pert(ipert)
250
251 pop_sub(em_wfs_tag)
252
253 end function em_wfs_tag
254
255! ---------------------------------------------------------
256! Provides indices of axes forming the right-hand system
257! to the given axis (to choose components of the position
258! and velocity, r_\alpha and V_\beta, for calculation of
259! the M_\gamma component of the magnetic dipole moment
260! M_\gamma = e_{\alpha \beta \gamma} r_\alpha V_\beta /2)
261 integer pure function magn_dir(dir, ind) result(dir_out)
262 integer, intent(in) :: dir, ind
263
264 select case (dir)
265 case (1)
266 if (ind == 1) then
267 dir_out = 2
268 else
269 dir_out = 3
270 end if
271 case (2)
272 if (ind == 1) then
273 dir_out = 3
274 else
275 dir_out = 1
276 end if
277 case (3)
278 if (ind == 1) then
279 dir_out = 1
280 else
281 dir_out = 2
282 end if
283 case default
284 dir_out = 0
285 end select
286
287 end function magn_dir
288
290! ---------------------------------------------------------
291
292 character(len=2) pure function index2pert(ipert) result(ch)
293 integer, intent(in) :: ipert
294
295 select case (ipert)
296 case (1)
297 ch = 'B'
298 case (2)
299 ch = 'K2'
300 case (3)
301 ch = 'KB'
302 case (4)
303 ch = 'KE'
304 case (5)
305 ch = 'E'
306 end select
307
308 end function index2pert
309
310#include "undef.F90"
311#include "real.F90"
312#include "em_resp_calc_inc.F90"
313
314#include "undef.F90"
315#include "complex.F90"
316#include "em_resp_calc_inc.F90"
317
318end module em_resp_calc_oct_m
319
320!! Local Variables:
321!! mode: f90
322!! coding: utf-8
323!! End:
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public dinhomog_kb_tot(sh, namespace, space, gr, st, hm, ions, idir, idir1, idir2, lr_k, lr_b, lr_k1, lr_k2, lr_kk1, lr_kk2, psi_out)
subroutine, public zlr_calc_magneto_optics_finite(sh, sh_mo, namespace, space, gr, st, hm, ions, nsigma, nfactor, lr_e, lr_b, chi)
subroutine, public dlr_calc_magnetization_periodic(namespace, space, mesh, st, hm, lr_k, magn)
integer pure function, public magn_dir(dir, ind)
subroutine, public dpost_orthogonalize(space, mesh, st, nfactor, nsigma, freq_factor, omega, eta, em_lr, kdotp_em_lr)
subroutine, public zlr_calc_magneto_optics_periodic(sh, sh2, namespace, space, gr, st, hm, ions, nsigma, nfactor, nfactor_ke, freq_factor, lr_e, lr_b, lr_k, lr_ke, lr_kb, frequency, zpol, zpol_kout)
subroutine, public zinhomog_k2_tot(namespace, space, gr, st, hm, ions, idir1, idir2, lr_k1, lr_k2, psi_out)
subroutine, public dlr_calc_beta(sh, namespace, space, gr, st, hm, xc, em_lr, dipole, beta, kdotp_lr, kdotp_em_lr, occ_response, dl_eig)
See (16) in X Andrade et al., J. Chem. Phys. 126, 184106 (2006) for finite systems and (10) in A Dal ...
subroutine, public dlr_calc_susceptibility(namespace, space, gr, st, hm, lr, nsigma, pert, chi_para, chi_dia)
subroutine, public lr_calc_current(st, space, gr, lr, lr_m)
subroutine, public dinhomog_ke_tot(sh, namespace, space, gr, st, hm, ions, idir, nsigma, lr_k, lr_e, lr_kk, psi_out)
subroutine, public dlr_calc_elf(st, space, gr, kpoints, lr, lr_m)
subroutine, public dem_resp_calc_eigenvalues(space, mesh, latt, st, dl_eig)
subroutine, public zinhomog_b(sh, namespace, space, gr, st, hm, ions, idir1, idir2, lr_k1, lr_k2, psi_out)
subroutine, public zcalc_polarizability_periodic(space, mesh, symm, st, em_lr, kdotp_lr, nsigma, zpol, ndir, zpol_k)
alpha_ij(w) = -e sum(m occ, k) [(<u_mk(0)|-id/dk_i)|u_mkj(1)(w)> + <u_mkj(1)(-w)|(-id/dk_i|u_mk(0)>)]
subroutine, public dinhomog_b(sh, namespace, space, gr, st, hm, ions, idir1, idir2, lr_k1, lr_k2, psi_out)
subroutine, public zlr_calc_elf(st, space, gr, kpoints, lr, lr_m)
character(len=2) pure function index2pert(ipert)
subroutine, public dcalc_polarizability_periodic(space, mesh, symm, st, em_lr, kdotp_lr, nsigma, zpol, ndir, zpol_k)
alpha_ij(w) = -e sum(m occ, k) [(<u_mk(0)|-id/dk_i)|u_mkj(1)(w)> + <u_mkj(1)(-w)|(-id/dk_i|u_mk(0)>)]
subroutine, public zpost_orthogonalize(space, mesh, st, nfactor, nsigma, freq_factor, omega, eta, em_lr, kdotp_em_lr)
character(len=100) function, public em_wfs_tag(idir, ifactor, idir2, ipert)
subroutine, public zlr_calc_magnetization_periodic(namespace, space, mesh, st, hm, lr_k, magn)
subroutine, public zinhomog_ke_tot(sh, namespace, space, gr, st, hm, ions, idir, nsigma, lr_k, lr_e, lr_kk, psi_out)
character(len=12) function, public freq2str(freq)
subroutine, public dlr_calc_susceptibility_periodic(namespace, space, symm, mesh, st, hm, lr_k, lr_b, lr_kk, lr_kb, magn)
subroutine, public zlr_calc_susceptibility_periodic(namespace, space, symm, mesh, st, hm, lr_k, lr_b, lr_kk, lr_kb, magn)
subroutine, public zlr_calc_susceptibility(namespace, space, gr, st, hm, lr, nsigma, pert, chi_para, chi_dia)
subroutine, public dlr_calc_magneto_optics_finite(sh, sh_mo, namespace, space, gr, st, hm, ions, nsigma, nfactor, lr_e, lr_b, chi)
character(len=100) function, public em_rho_tag(freq, dir, dir2, ipert)
subroutine, public dinhomog_k2_tot(namespace, space, gr, st, hm, ions, idir1, idir2, lr_k1, lr_k2, psi_out)
subroutine, public zlr_calc_beta(sh, namespace, space, gr, st, hm, xc, em_lr, dipole, beta, kdotp_lr, kdotp_em_lr, occ_response, dl_eig)
See (16) in X Andrade et al., J. Chem. Phys. 126, 184106 (2006) for finite systems and (10) in A Dal ...
subroutine, public zcalc_polarizability_finite(namespace, space, gr, st, hm, lr, nsigma, pert, zpol, doalldirs, ndir)
alpha_ij(w) = - sum(m occ) [<psi_m(0)|r_i|psi_mj(1)(w)> + <psi_mj(1)(-w)|r_i|psi_m(0)>] minus sign is...
subroutine, public dcalc_polarizability_finite(namespace, space, gr, st, hm, lr, nsigma, pert, zpol, doalldirs, ndir)
alpha_ij(w) = - sum(m occ) [<psi_m(0)|r_i|psi_mj(1)(w)> + <psi_mj(1)(-w)|r_i|psi_m(0)>] minus sign is...
subroutine, public dlr_calc_magneto_optics_periodic(sh, sh2, namespace, space, gr, st, hm, ions, nsigma, nfactor, nfactor_ke, freq_factor, lr_e, lr_b, lr_k, lr_ke, lr_kb, frequency, zpol, zpol_kout)
subroutine, public zinhomog_kb_tot(sh, namespace, space, gr, st, hm, ions, idir, idir1, idir2, lr_k, lr_b, lr_k1, lr_k2, lr_kk1, lr_kk2, psi_out)
subroutine, public zem_resp_calc_eigenvalues(space, mesh, latt, st, dl_eig)
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
complex(real64), parameter, public m_zi
Definition: global.F90:205
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:118
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
This module handles spin dimensions of the states and the k-point distribution.
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
character pure function, public index2axis(idir)
Definition: utils.F90:205
Definition: xc.F90:116