Octopus
td_calc.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
21module td_calc_oct_m
22 use comm_oct_m
23 use debug_oct_m
27 use forces_oct_m
28 use global_oct_m
29 use grid_oct_m
33 use ions_oct_m
34 use iso_c_binding
35 use, intrinsic :: iso_fortran_env
37 use lasers_oct_m
38 use lda_u_oct_m
41 use mesh_oct_m
45 use space_oct_m
49 use xc_oct_m
50
51 implicit none
52
53 private
54 public :: &
58
59contains
60
61 ! ---------------------------------------------------------
71 subroutine td_calc_tacc(namespace, space, gr, ions, ext_partners, st, hm, acc, time)
72 type(namespace_t), intent(in) :: namespace
73 class(space_t), intent(in) :: space
74 type(grid_t), intent(in) :: gr
75 type(ions_t), intent(in) :: ions
76 type(partner_list_t), intent(in) :: ext_partners
77 type(states_elec_t), intent(in) :: st
78 type(hamiltonian_elec_t), intent(in) :: hm
79 real(real64), intent(in) :: time
80 real(real64), intent(out) :: acc(:)
81
82 real(real64) :: field(space%dim), x(space%dim)
83 complex(real64), allocatable :: zpsi(:, :), hzpsi(:,:), hhzpsi(:,:), xzpsi(:,:,:), vnl_xzpsi(:,:)
84 integer :: j, k, ik, ist, idim
85 type(lasers_t), pointer :: lasers
86
87 push_sub(td_calc_tacc)
88
89 ! See #683
90 if (family_is_mgga(hm%xc%family) .or. family_is_hybrid(hm%xc)) then
91 call messages_not_implemented("Dipole acceleration with mGGAs and hybrid functionals")
92 end if
93
94 if (hm%lda_u_level /= dft_u_none) then
95 call messages_not_implemented("Dipole acceleration with DFT+U")
96 end if
97
98 ! This is due to an ASSERT in X(total_force_from_potential)
99 if (st%symmetrize_density) then
100 call messages_not_implemented("Dipole acceleration with SymmetrizeDensity = yes")
101 end if
102
103 ! The term i<[V_l,p]> + i<[V_nl,p]> may be considered as equal but opposite to the
104 ! force exerted by the electrons on the ions. COMMENT: This has to be thought about.
105 ! Maybe we are forgetting something....
106 call total_force_calculate(space, gr, ions, hm%ep, st, hm%phase, acc, hm%lda_u_level)
107
108 ! Adds the laser contribution : i<[V_laser, p]>
109 ! WARNING: this ignores the possibility of non-electric td external fields.
110 field = m_zero
111 lasers => list_get_lasers(ext_partners)
112 if(associated(lasers)) then
113 do j = 1, lasers%no_lasers
114 call laser_electric_field(lasers%lasers(j), field(1:space%dim), time, 0.001_real64)
115 acc(1:space%dim) = acc(1:space%dim) - st%qtot*field(1:space%dim)
116 end do
117 end if
118
119 if (.not. hm%ep%non_local) then
120 pop_sub(td_calc_tacc)
121 return
122 end if
123
124 ! And now, i<[H,[V_nl,x]]>
125 x = m_zero
126 safe_allocate(zpsi(1:gr%np_part, 1:st%d%dim))
127 safe_allocate(hzpsi(1:gr%np_part, 1:st%d%dim))
128 safe_allocate(hhzpsi(1:3, 1:gr%np_part))
129
130 do ik = st%d%kpt%start, st%d%kpt%end
131 do ist = st%st_start, st%st_end
132
133 call states_elec_get_state(st, gr, ist, ik, zpsi)
134
135 call zhamiltonian_elec_apply_single(hm, namespace, gr, zpsi, hzpsi, ist, ik)
136
137 safe_allocate(xzpsi(1:gr%np_part, 1:st%d%dim, 1:3))
138 safe_allocate(vnl_xzpsi(1:gr%np_part, 1:st%d%dim))
139 xzpsi = m_z0
140 do k = 1, gr%np
141 do j = 1, space%dim
142 xzpsi(k, 1:st%d%dim, j) = gr%x(k, j)*zpsi(k, 1:st%d%dim)
143 end do
144 end do
145
146 do j = 1, space%dim
147 call zhamiltonian_elec_apply_single(hm, namespace, gr, xzpsi(:, :, j), vnl_xzpsi, ist, ik, &
149
150 do idim = 1, st%d%dim
151 x(j) = x(j) - 2*st%occ(ist, ik)*real(zmf_dotp(gr, hzpsi(1:gr%np, idim), vnl_xzpsi(:, idim)), real64)
152 end do
153 end do
154
155 xzpsi = m_z0
156 do k = 1, gr%np
157 do j = 1, space%dim
158 xzpsi(k, 1:st%d%dim, j) = gr%x(k, j)*hzpsi(k, 1:st%d%dim)
159 end do
160 end do
161
162 do j = 1, space%dim
163 call zhamiltonian_elec_apply_single(hm, namespace, gr, xzpsi(:, :, j), vnl_xzpsi, ist, ik, &
165
166 do idim = 1, st%d%dim
167 x(j) = x(j) + 2*st%occ(ist, ik)*real(zmf_dotp(gr, zpsi(:, idim), vnl_xzpsi(:, idim)), real64)
168 end do
169 end do
170 safe_deallocate_a(xzpsi)
171 safe_deallocate_a(vnl_xzpsi)
172
173 end do
174 end do
175 safe_deallocate_a(hzpsi)
176 safe_deallocate_a(hhzpsi)
177
178 if (st%parallel_in_states) then
179 call comm_allreduce(st%mpi_grp, x)
180 end if
181 acc = acc + x
182
183 pop_sub(td_calc_tacc)
184 end subroutine td_calc_tacc
185
186! ---------------------------------------------------------
192! ---------------------------------------------------------
193 subroutine td_calc_tvel(gr, st, space, kpoints, vel)
194 type(grid_t), intent(in) :: gr
195 type(states_elec_t), intent(in) :: st
196 type(space_t), intent(in) :: space
197 type(kpoints_t), intent(in) :: kpoints
198 real(real64), intent(out) :: vel(:)
199
200 real(real64), allocatable :: momentum(:,:,:)
201 integer :: st_start, st_end
202
203 push_sub(td_calc_tvel)
204
205 st_start = st%st_start
206 st_end = st%st_end
207
208 safe_allocate(momentum(1:space%dim, st_start:st_end, st%d%kpt%start:st%d%kpt%end))
209 call elec_momentum_me(gr, st, space, kpoints, momentum)
210
211 momentum(1:space%dim, st_start:st_end, 1) = &
212 sum(momentum(1:space%dim, st_start:st_end, st%d%kpt%start:st%d%kpt%end), 3)
213 momentum(1:space%dim, 1, 1) = sum(momentum(1:space%dim, st_start:st_end, 1), 2)
214 vel = momentum(:, 1, 1)
215
216 safe_deallocate_a(momentum)
217
218 pop_sub(td_calc_tvel)
219 end subroutine td_calc_tvel
220
221! ---------------------------------------------------------
224! ---------------------------------------------------------
225 subroutine td_calc_ionch(mesh, st, ch, Nch)
226 type(mesh_t), intent(in) :: mesh
227 type(states_elec_t), intent(in) :: st
228 integer, intent(in) :: nch
229 real(real64), intent(out) :: ch(0:nch)
230
231 integer :: ik, ist, ii, jj, idim, nid
232 real(real64) :: prod, prod0
233 real(real64), allocatable :: n(:), nnot(:)
234 complex(real64), allocatable :: zpsi(:)
235 !combinations
236 integer :: next
237 type(c_ptr) :: c
238 integer, allocatable :: idx0(:), idx(:), idxref(:)
239
240 push_sub(td_calc_ionch)
241
242 safe_allocate( n(1: nch))
243 safe_allocate(nnot(1: nch))
244 safe_allocate(zpsi(1:mesh%np))
245
246 n(:) = m_zero
247 nnot(:)= m_zero
248
249
250 ii = 1
251 do ik = 1, st%nik
252 do ist = 1, st%nst
253 do idim = 1, st%d%dim
254
255 if (st%st_start <= ist .and. ist <= st%st_end .and. &
256 st%d%kpt%start <= ik .and. ik <= st%d%kpt%end) then
257 call states_elec_get_state(st, mesh, idim, ist, ik, zpsi)
258 n(ii) = zmf_nrm2(mesh, zpsi)
259 nnot(ii) = m_one - n(ii)**2
260 end if
261
262 ii = ii + 1
263
264 end do
265 end do
266 end do
267
268 if (st%parallel_in_states) then
269 call comm_allreduce(st%mpi_grp, n)
270 call comm_allreduce(st%mpi_grp, nnot)
271 end if
272
273! print * ,mpi_world%rank, "N =", N(:)
274! print * ,mpi_world%rank, "Nnot =", Nnot(:)
275
276 ch(:) = m_zero
277
278 nid = nch - 1
279 safe_allocate(idx(0:nid))
280 safe_allocate(idxref(0:nid))
281 idxref = (/ (ii, ii = 0, nid) /)
282 do ii = 0, nch
283! print *, "Nch= ", Nch, "ii", ii
284 call loct_combination_init(c, nch, ii)
285 if (ii == 0) then
286 safe_allocate(idx0(0:0))
287 else
288 safe_allocate(idx0(0:ii-1))
289 end if
290! print *,"size(idx0,1)=",size(idx0,1)
291 if (debug%info) then
292 call messages_write("P(")
293 call messages_write(ii)
294 call messages_write(") = 0")
295 call messages_new_line()
296 end if
297 ! Loop over combinations
298 do
299 prod = m_one
300 prod0 = m_one
301 if (debug%info) then
302 call messages_write("P(")
303 call messages_write(ii)
304 call messages_write(") += ")
305 end if
306 call loct_get_combination(c, idx0)
307 idx(:) = idxref(:)
308 do jj = 0, ii-1
309 idx(idx0(jj))= -1
310 if (debug%info) then
311 call messages_write(" No(")
312 call messages_write(idx0(jj)+1)
313 call messages_write(") ")
314 end if
315 prod0 = prod0 * nnot(idx0(jj)+1)
316 end do
317
318 do jj = 0, nid
319 if (idx(jj) >= 0) then
320 if (debug%info) then
321 call messages_write(" N(")
322 call messages_write(idx(jj)+1)
323 call messages_write(") ")
324 end if
325 prod = prod * n(idx(jj)+1)
326 end if
327 end do
328
329 if (debug%info) call messages_new_line()
330
331 ch(ii) = ch(ii) + prod*prod0
332
333 call loct_combination_next(c, next)
334 if (next /= 0) exit
335 end do
336 safe_deallocate_a(idx0)
338
339 if (debug%info) call messages_info()
340 end do
341
342 safe_deallocate_a(idx)
343 safe_deallocate_a(idxref)
344
345
346 safe_deallocate_a(n)
347 safe_deallocate_a(nnot)
348 safe_deallocate_a(zpsi)
349
350 pop_sub(td_calc_ionch)
351 end subroutine td_calc_ionch
352
353
354end module td_calc_oct_m
355
356!! Local Variables:
357!! mode: f90
358!! coding: utf-8
359!! End:
Functions to generate combinations.
Definition: loct_math.F90:293
type(debug_t), save, public debug
Definition: debug.F90:156
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public elec_momentum_me(gr, st, space, kpoints, momentum)
type(lasers_t) function, pointer, public list_get_lasers(partners)
subroutine, public total_force_calculate(space, gr, ions, ep, st, phase, x, lda_u)
This computes the total forces on the ions created by the electrons (it excludes the force due to pos...
Definition: forces.F90:191
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the underlying real-space grid.
Definition: grid.F90:117
integer, parameter, public term_non_local_potential
subroutine, public zhamiltonian_elec_apply_single(hm, namespace, mesh, psi, hpsi, ist, ik, terms, set_bc, set_phase)
This module defines classes and functions for interaction partners.
subroutine, public laser_electric_field(laser, field, time, dt)
Returns a vector with the electric field, no matter whether the laser is described directly as an ele...
Definition: lasers.F90:1191
integer, parameter, public dft_u_none
Definition: lda_u.F90:201
subroutine, public loct_get_combination(c, comb)
Definition: loct_math.F90:377
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1113
subroutine, public messages_new_line()
Definition: messages.F90:1134
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public td_calc_tacc(namespace, space, gr, ions, ext_partners, st, hm, acc, time)
Electronic acceleration (to calculate harmonic spectrum...) It is calculated as:
Definition: td_calc.F90:165
subroutine, public td_calc_ionch(mesh, st, ch, Nch)
Multiple ionization probabilities calculated form the KS orbital densities C. Ullrich,...
Definition: td_calc.F90:319
subroutine, public td_calc_tvel(gr, st, space, kpoints, vel)
Electronic velocity (to calculate harmonic spectrum...) It is calculated as:
Definition: td_calc.F90:287
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.