Octopus
em_resp.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
21module em_resp_oct_m
23 use box_oct_m
24 use debug_oct_m
27 use forces_oct_m
28 use global_oct_m
29 use grid_oct_m
32 use io_oct_m
33 use ions_oct_m
34 use, intrinsic :: iso_fortran_env
35 use kdotp_oct_m
39 use loct_oct_m
40 use mesh_oct_m
43 use mpi_oct_m
47 use output_oct_m
48 use parser_oct_m
49 use pcm_oct_m
57 use smear_oct_m
58 use sort_oct_m
59 use space_oct_m
65 use string_oct_m
66 use unit_oct_m
68 use utils_oct_m
69 use v_ks_oct_m
70 use xc_oct_m
71
72 implicit none
73
74 private
75
76 public :: &
79
80 integer, parameter :: &
81 PERTURBATION_ELECTRIC = 1, &
84
85
86 type em_resp_t
87 private
88 class(perturbation_t), pointer :: perturbation
89
90 integer :: nsigma
92 integer :: nfactor
95 integer :: nomega
96
97 real(real64) :: eta
98 real(real64) :: freq_factor(3)
99 real(real64), allocatable :: omega(:)
100 type(lr_t), allocatable :: lr(:,:,:)
101 complex(real64), allocatable :: alpha_k(:, :, :, :)
103 complex(real64), allocatable :: alpha_be_k(:, :, :, :)
105 logical :: calc_hyperpol
106 complex(real64) :: alpha(3, 3, 3)
107 complex(real64) :: alpha_be(3, 3, 3)
108 complex(real64) :: alpha0(3, 3, 3)
110 complex(real64) :: alpha_be0(3, 3, 3)
112 complex(real64) :: beta (3, 3, 3)
113
114 complex(real64) :: chi_para(3, 3)
115 complex(real64) :: chi_dia (3, 3)
116 complex(real64) :: magn(3)
117
118 logical :: ok(1:3)
119 logical :: force_no_kdotp
120
121 logical :: calc_rotatory
122 logical :: calc_Born
123 type(Born_charges_t) :: Born_charges(3)
124 logical :: occ_response
125 logical :: wfns_from_scratch
126 logical :: calc_magnetooptics
127 logical :: magnetooptics_nohvar
129 logical :: kpt_output
131 logical :: lrc_kernel
132
133 end type em_resp_t
134
135contains
136
137 ! ---------------------------------------------------------
138 subroutine em_resp_run(system, from_scratch)
139 class(*), intent(inout) :: system
140 logical, intent(in) :: from_scratch
141
142 push_sub(em_resp_run)
143
144 select type (system)
145 class is (multisystem_basic_t)
146 message(1) = "CalculationMode = em_resp not implemented for multi-system calculations"
147 call messages_fatal(1, namespace=system%namespace)
148 type is (electrons_t)
149 call em_resp_run_legacy(system, from_scratch)
150 end select
151
152 pop_sub(em_resp_run)
153 end subroutine em_resp_run
154
155 ! ---------------------------------------------------------
156 subroutine em_resp_run_legacy(sys, fromScratch)
157 type(electrons_t), intent(inout) :: sys
158 logical, intent(in) :: fromScratch
159
160 type(em_resp_t) :: em_vars
161 type(sternheimer_t) :: sh, sh_kdotp, sh2, sh_kmo, sh_mo
162 type(lr_t) :: kdotp_lr(sys%space%dim, 1)
163 type(lr_t), allocatable :: kdotp_em_lr2(:, :, :, :)
164 type(lr_t), allocatable :: b_lr(:, :)
165 type(lr_t), allocatable :: kb_lr(:, :, :), k2_lr(:, :, :)
166 type(lr_t), allocatable :: ke_lr(:, :, :, :)
167 class(perturbation_t), pointer :: pert_kdotp, pert2_none, pert_b
168
169 integer :: sigma, idir, idir2, ierr, iomega, ifactor
170 integer :: ierr_e(3), ierr_e2(3), nfactor_ke
171 character(len=100) :: str_tmp
172 logical :: complex_response, have_to_calculate, use_kdotp, opp_freq, &
173 exact_freq(3), complex_wfs, allocate_rho_em, allocate_rho_mo
174 logical :: magnetic_pert
175
176 real(real64) :: last_omega, frequency, dfrequency_eta
177 real(real64), allocatable :: dl_eig(:,:,:)
178 complex(real64) :: zfrequency_eta, lrc_coef(sys%space%dim, sys%space%dim)
179 type(restart_t) :: gs_restart, kdotp_restart
180
182
183 if (sys%hm%pcm%run_pcm) then
184 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
185 end if
186
187 if (sys%kpoints%use_symmetries) then
188 call messages_experimental("em_resp with k-points symmetries", namespace=sys%namespace)
189 end if
191 if (sys%kpoints%reduced%npoints /= sys%kpoints%full%npoints) then
192 call messages_experimental('em_resp with reduced k-grid', namespace=sys%namespace)
193 end if
195 call parse_input()
197 select type(ptr=>em_vars%perturbation)
199 if (any(abs(em_vars%omega(1:em_vars%nomega)) > m_epsilon)) then
200 call messages_not_implemented('Dynamical magnetic response', namespace=sys%namespace)
201 end if
202 end select
204 em_vars%lrc_kernel = .false.
205 if (abs(sys%ks%xc%kernel_lrc_alpha) > m_epsilon) em_vars%lrc_kernel = .true.
206
207 if (em_vars%lrc_kernel .and. sys%space%periodic_dim < sys%space%dim) then
208 message(1) = 'The use of the LRC kernel for non-periodic dimensions makes no sense.'
209 call messages_warning(1, namespace=sys%namespace)
210 end if
212 complex_wfs = states_are_complex(sys%st)
213 complex_response = (em_vars%eta > m_epsilon) .or. states_are_complex(sys%st)
214 call restart_init(gs_restart, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
215 if (ierr == 0) then
216 call states_elec_look_and_load(gs_restart, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
217 is_complex = complex_response)
218 call restart_end(gs_restart)
219 else
220 message(1) = "Previous gs calculation is required."
221 call messages_fatal(1, namespace=sys%namespace)
222 end if
223
224 ! Use of ForceComplex will make this true after states_elec_look_and_load even if it was not before.
225 ! Otherwise, this line is a tautology.
226 complex_response = states_are_complex(sys%st)
227
228 if (states_are_real(sys%st)) then
229 message(1) = 'Info: Using real wavefunctions.'
230 else
231 message(1) = 'Info: Using complex wavefunctions.'
232 end if
233 call messages_info(1, namespace=sys%namespace)
234
235 ! setup Hamiltonian
236 message(1) = 'Info: Setting up Hamiltonian for linear response'
237 call messages_info(1, namespace=sys%namespace)
238 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
239
240 use_kdotp = sys%space%is_periodic() .and. .not. em_vars%force_no_kdotp
241
242 if (use_kdotp .and. .not. smear_is_semiconducting(sys%st%smear)) then
243 ! there needs to be a gap.
244 message(1) = "em_resp with kdotp can only be used with semiconducting smearing"
245 call messages_fatal(1, namespace=sys%namespace)
246 end if
247
248 ! read kdotp wavefunctions if necessary
249 if (use_kdotp) then
250 message(1) = "Reading kdotp wavefunctions for periodic directions."
251 call messages_info(1, namespace=sys%namespace)
252
253 call restart_init(kdotp_restart, sys%namespace, restart_kdotp, restart_type_load, sys%mc, ierr, mesh=sys%gr)
254 if (ierr /= 0) then
255 message(1) = "Unable to read kdotp wavefunctions."
256 message(2) = "Previous kdotp calculation required."
257 call messages_fatal(2, namespace=sys%namespace)
258 end if
259
260 do idir = 1, sys%space%periodic_dim
261 call lr_init(kdotp_lr(idir, 1))
262 call lr_allocate(kdotp_lr(idir, 1), sys%st, sys%gr, allocate_rho = .false.)
263
264 ! load wavefunctions
265 str_tmp = kdotp_wfs_tag(idir)
266 ! 1 is the sigma index which is used in em_resp
267 call restart_open_dir(kdotp_restart, wfs_tag_sigma(sys%namespace, str_tmp, 1), ierr)
268 if (ierr == 0) then
269 call states_elec_load(kdotp_restart, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, ierr, &
270 lr=kdotp_lr(idir, 1))
271 end if
272 call restart_close_dir(kdotp_restart)
273
274 if (ierr /= 0) then
275 message(1) = "Could not load kdotp wavefunctions from '"//trim(wfs_tag_sigma(sys%namespace, str_tmp, 1))//"'"
276 message(2) = "Previous kdotp calculation required."
277 call messages_fatal(2, namespace=sys%namespace)
278 end if
279 end do
280
281 call restart_end(kdotp_restart)
282 end if
283
284 em_vars%nfactor = 1
285 if (em_vars%calc_hyperpol) em_vars%nfactor = 3
286
287 ! in effect, nsigma = 1 only if hyperpol not being calculated, and the only frequency is zero
288 if (em_vars%calc_hyperpol .or. any(abs(em_vars%omega(1:em_vars%nomega)) > m_epsilon)) then
289 em_vars%nsigma = 2
290 ! positive and negative values of the frequency must be considered
291 else
292 em_vars%nsigma = 1
293 ! only considering positive values
294 end if
295
296 if (em_vars%calc_hyperpol .and. use_kdotp) then
297 pert_kdotp => perturbation_kdotp_t(sys%namespace, sys%ions)
298 pert2_none => perturbation_none_t(sys%namespace)
299 call messages_experimental("Second-order Sternheimer equation", namespace=sys%namespace)
300 call pert2_none%setup_dir(1) ! direction is irrelevant
301 safe_allocate(kdotp_em_lr2(1:sys%space%periodic_dim, 1:sys%space%dim, 1:em_vars%nsigma, 1:em_vars%nfactor))
302 do ifactor = 1, em_vars%nfactor
303 do sigma = 1, em_vars%nsigma
304 do idir = 1, sys%space%periodic_dim
305 do idir2 = 1, sys%space%dim
306 call lr_init(kdotp_em_lr2(idir, idir2, sigma, ifactor))
307 call lr_allocate(kdotp_em_lr2(idir, idir2, sigma, ifactor), sys%st, sys%gr, allocate_rho = .false.)
308 end do
309 end do
310 end do
311 end do
312 call sternheimer_init(sh2, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
313 complex_response, set_ham_var = 0, set_last_occ_response = .false.)
314 call sternheimer_init(sh_kdotp, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
315 complex_response, set_ham_var = 0, set_last_occ_response = .true.)
316 em_vars%occ_response = .true.
317 safe_allocate(dl_eig(1:sys%st%nst, 1:sys%st%nik, 1:sys%space%periodic_dim))
318 end if
319
320 ! Hyperpolarizability requires full corrections to wavefunctions (with projections on occupied states).
321 ! Magnetooptics is implemented only for projections of corrections to wavefunctions on unoccupied states.
322 if (em_vars%calc_magnetooptics) then
323 if (em_vars%calc_hyperpol .and. use_kdotp) then
324 message(1) = "Hyperpolarizability and magnetooptics with kdotp are not compatible."
325 message(2) = "Only calculation of hyperpolarizability will be performed."
326 call messages_warning(2, namespace=sys%namespace)
327 em_vars%calc_magnetooptics = .false.
328 else
329 em_vars%nfactor = 2
330 em_vars%freq_factor(1) = m_one
331 em_vars%freq_factor(2) = -m_one
332 end if
333 end if
334
335 magnetic_pert = .false.
336 select type(ptr => em_vars%perturbation)
338 em_vars%nsigma = 1
339 if (use_kdotp) call messages_experimental("Magnetic perturbation for periodic systems", namespace=sys%namespace)
340 magnetic_pert = .true.
341 end select
342
343 if (em_vars%calc_magnetooptics .or. magnetic_pert) then
344 em_vars%occ_response = .false.
345
346 if (use_kdotp) then
347 pert2_none => perturbation_none_t(sys%namespace)
348 call pert2_none%setup_dir(1)
349
350 safe_allocate(k2_lr(1:sys%space%dim, 1:sys%space%dim, 1))
351 safe_allocate(kb_lr(1:sys%space%dim, 1:sys%space%dim, 1))
352 do idir = 1, sys%space%dim
353 do idir2 = 1, sys%space%dim
354 call lr_init(kb_lr(idir, idir2, 1))
355 call lr_allocate(kb_lr(idir, idir2, 1), sys%st, sys%gr, allocate_rho = .false.)
356 if (idir2 <= idir) then
357 call lr_init(k2_lr(idir, idir2, 1))
358 call lr_allocate(k2_lr(idir, idir2, 1), sys%st, sys%gr, allocate_rho = .false.)
359 end if
360 end do
361 end do
362
363 if (sys%space%periodic_dim < sys%space%dim) then
364 if (magnetic_pert) then
365 message(1) = "All directions should be periodic for magnetic perturbations with kdotp."
366 else
367 message(1) = "All directions should be periodic for magnetooptics with kdotp."
368 end if
369 call messages_fatal(1, namespace=sys%namespace)
370 end if
371 if (.not. complex_response) then
372 do idir = 1, sys%space%dim
373 call dlr_orth_response(sys%gr, sys%st, kdotp_lr(idir, 1), m_zero)
374 end do
375 else
376 do idir = 1, sys%space%dim
377 call zlr_orth_response(sys%gr, sys%st, kdotp_lr(idir, 1), m_z0)
378 end do
379 end if
380 call sternheimer_init(sh_kmo, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
381 complex_response, set_ham_var = 0, set_last_occ_response = em_vars%occ_response)
382 end if
383 end if
384
385 safe_allocate(em_vars%lr(1:sys%space%dim, 1:em_vars%nsigma, 1:em_vars%nfactor))
386 do ifactor = 1, em_vars%nfactor
387 call born_charges_init(em_vars%Born_charges(ifactor), sys%namespace, sys%ions%natoms, &
388 sys%st%val_charge, sys%st%qtot, sys%space%dim)
389 end do
390
391 if (magnetic_pert .and. sys%st%d%nspin == 1 .and. states_are_real(sys%st)) then
392 ! first-order response is zero if there is time-reversal symmetry. F Mauri and SG Louie, PRL 76, 4246 (1996)
393 call sternheimer_init(sh, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
394 complex_response, set_ham_var = 0, set_last_occ_response = em_vars%occ_response)
395 ! set HamiltonianVariation to V_ext_only, in magnetic case
396 else
397 call sternheimer_init(sh, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
398 complex_response, set_last_occ_response = em_vars%occ_response)
399 ! otherwise, use default, which is hartree + fxc
400 end if
401
402 if (em_vars%lrc_kernel .and. (.not. sternheimer_add_hartree(sh)) &
403 .and. (.not. sternheimer_add_fxc(sh))) then
404 message(1) = "Only the G = G'= 0 term of the LRC kernel is taken into account."
405 call messages_warning(1, namespace=sys%namespace)
406 end if
407
408
409 if (mpi_grp_is_root(mpi_world)) then
410 call info()
411 call io_mkdir(em_resp_dir, sys%namespace) ! output
412 end if
413
414 allocate_rho_em = sternheimer_add_fxc(sh) .or. sternheimer_add_hartree(sh)
415 do ifactor = 1, em_vars%nfactor
416 do idir = 1, sys%space%dim
417 do sigma = 1, em_vars%nsigma
418 call lr_init(em_vars%lr(idir, sigma, ifactor))
419 call lr_allocate(em_vars%lr(idir, sigma, ifactor), sys%st, sys%gr, allocate_rho = allocate_rho_em)
420 end do
421 end do
422 end do
423
424 select type(ptr=> em_vars%perturbation)
426 !Do nothing
427 class default
428 em_vars%kpt_output = .false.
429 end select
430 if (.not. use_kdotp .or. sys%st%nik == 1) em_vars%kpt_output = .false.
431
432 if (em_vars%kpt_output) then
433 safe_allocate(em_vars%alpha_k(1:sys%space%dim, 1:sys%space%dim, 1:em_vars%nfactor, 1:sys%st%nik))
434 end if
435
436 if (em_vars%calc_magnetooptics) then
437 if (em_vars%magnetooptics_nohvar) then
438 call sternheimer_init(sh_mo, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
439 complex_response, set_ham_var = 0, set_last_occ_response = em_vars%occ_response)
440 else
441 call sternheimer_init(sh_mo, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
442 complex_response, set_last_occ_response = em_vars%occ_response)
443 call sternheimer_build_kxc(sh_mo, sys%namespace, sys%gr, sys%st, sys%ks%xc)
444 end if
445 call messages_experimental("Magneto-optical response", namespace=sys%namespace)
446 allocate_rho_mo = sternheimer_add_fxc(sh_mo) .or. sternheimer_add_hartree(sh_mo)
447 safe_allocate(b_lr(1:sys%space%dim, 1))
448 do idir = 1, sys%space%dim
449 call lr_init(b_lr(idir, 1))
450 call lr_allocate(b_lr(idir, 1), sys%st, sys%gr, allocate_rho = allocate_rho_mo)
451 end do
452
453 if (use_kdotp) then
454 if (em_vars%kpt_output) then
455 safe_allocate(em_vars%alpha_be_k(1:sys%space%dim, 1:sys%space%dim, 1:sys%space%dim, 1:sys%st%nik))
456 end if
457 nfactor_ke = 1
458 if (sys%kpoints%use_time_reversal .and. sys%kpoints%full%npoints > 1) nfactor_ke = em_vars%nfactor
459 safe_allocate(ke_lr(1:sys%space%dim, 1:sys%space%dim, 1:em_vars%nsigma, 1:nfactor_ke))
460 do idir = 1, sys%space%dim
461 do idir2 = 1, sys%space%dim
462 do sigma = 1, em_vars%nsigma
463 do ifactor = 1, nfactor_ke
464 call lr_init(ke_lr(idir, idir2, sigma, ifactor))
465 call lr_allocate(ke_lr(idir, idir2, sigma, ifactor), sys%st, sys%gr, allocate_rho = .false.)
466 end do
467 end do
468 end do
469 end do
470 else
471 pert_b => perturbation_magnetic_t(sys%namespace, sys%ions)
472 end if
473 end if
474
475
476 last_omega = m_huge
477 do iomega = 1, em_vars%nomega
478
479 em_vars%ok(1:3) = .true.
480
481 do ifactor = 1, em_vars%nfactor
482 frequency = em_vars%freq_factor(ifactor)*em_vars%omega(iomega)
483 zfrequency_eta = frequency + m_zi * em_vars%eta
484 if (em_vars%calc_magnetooptics .and. ifactor == 2) zfrequency_eta = frequency - m_zi * em_vars%eta
485 dfrequency_eta = real(zfrequency_eta, real64)
486
487 if (abs(frequency) < m_epsilon .and. em_vars%calc_magnetooptics .and. use_kdotp) then
488 message(1) = "Magnetooptical response with kdotp requires non-zero frequency."
489 call messages_warning(1, namespace=sys%namespace)
490 end if
491
492 ierr = 0
493 ierr_e(:) = 0
494 ierr_e2(:) = 0
495
496 have_to_calculate = .true.
497 opp_freq = .false.
498
499 ! if this frequency is zero and this is not the first
500 ! iteration we do not have to do anything
501 if (iomega > 1 .and. abs(em_vars%freq_factor(ifactor)) <= m_epsilon) have_to_calculate = .false.
502
503 if (ifactor > 1 .and. (.not. em_vars%calc_magnetooptics)) then
504
505 ! if this frequency is the same as the previous one, just copy it
506 if (have_to_calculate .and. abs(em_vars%freq_factor(ifactor - 1) * em_vars%omega(iomega) &
507 - frequency) < m_epsilon) then
508
509 do idir = 1, sys%space%dim
510 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 1, ifactor - 1), em_vars%lr(idir, 1, ifactor))
511 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 2, ifactor - 1), em_vars%lr(idir, 2, ifactor))
512
513 if (em_vars%calc_hyperpol .and. use_kdotp) then
514 do idir2 = 1, sys%space%periodic_dim
515 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 1, ifactor - 1), &
516 kdotp_em_lr2(idir, idir2, 1, ifactor))
517 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 2, ifactor - 1), &
518 kdotp_em_lr2(idir, idir2, 2, ifactor))
519 end do
520 end if
521 end do
522
523 have_to_calculate = .false.
524
525 end if
526
527 ! if this frequency is minus the previous one, copy it inverted
528 if (have_to_calculate .and. abs(em_vars%freq_factor(ifactor - 1) * em_vars%omega(iomega) &
529 + frequency) < m_epsilon) then
530
531 do idir = 1, sys%space%dim
532 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 1, ifactor - 1), em_vars%lr(idir, 2, ifactor))
533 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 2, ifactor - 1), em_vars%lr(idir, 1, ifactor))
534
535 if (em_vars%calc_hyperpol .and. use_kdotp) then
536 do idir2 = 1, sys%space%periodic_dim
537 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 1, ifactor - 1), &
538 kdotp_em_lr2(idir, idir2, 2, ifactor))
539 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 2, ifactor - 1), &
540 kdotp_em_lr2(idir, idir2, 1, ifactor))
541 end do
542 end if
543 end do
544
545 have_to_calculate = .false.
546
547 end if
548
549 end if
550
551 if (iomega > 1 .and. ifactor == 1 .and. (.not. em_vars%calc_magnetooptics)) then
552
553 ! if this frequency is the same as the previous one, just copy it
554 if (have_to_calculate .and. abs(frequency - last_omega) < m_epsilon) then
555
556 do idir = 1, sys%space%dim
557 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 1, em_vars%nfactor), em_vars%lr(idir, 1, 1))
558 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 2, em_vars%nfactor), em_vars%lr(idir, 2, 1))
559
560 if (em_vars%calc_hyperpol .and. use_kdotp) then
561 do idir2 = 1, sys%space%periodic_dim
562 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 1, em_vars%nfactor), &
563 kdotp_em_lr2(idir, idir2, 1, 1))
564 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 2, em_vars%nfactor), &
565 kdotp_em_lr2(idir, idir2, 2, 1))
566 end do
567 end if
568 end do
569
570 have_to_calculate = .false.
571
572 end if
573
574 ! if this frequency is minus the previous one, copy it inverted
575 if (have_to_calculate .and. abs(frequency + last_omega) < m_epsilon) then
576
577 do idir = 1, sys%space%dim
578 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 1, em_vars%nfactor), em_vars%lr(idir, 2, 1))
579 call lr_copy(sys%st, sys%gr, em_vars%lr(idir, 2, em_vars%nfactor), em_vars%lr(idir, 1, 1))
580
581 if (em_vars%calc_hyperpol .and. use_kdotp) then
582 do idir2 = 1, sys%space%periodic_dim
583 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 1, em_vars%nfactor), &
584 kdotp_em_lr2(idir, idir2, 2, 1))
585 call lr_copy(sys%st, sys%gr, kdotp_em_lr2(idir, idir2, 2, em_vars%nfactor), &
586 kdotp_em_lr2(idir, idir2, 1, 1))
587 end do
588 end if
589 end do
590
591 have_to_calculate = .false.
592
593 end if
594
595 end if
596
597 if (have_to_calculate) then
598
599 exact_freq(:) = .false.
600
601 if (states_are_real(sys%st)) then
602 call drun_sternheimer(em_vars, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%mc, &
603 sys%ions)
604 else
605 call zrun_sternheimer(em_vars, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%mc, &
606 sys%ions)
607 end if
608
609 end if ! have_to_calculate
610
611 if (.not. have_to_calculate) cycle
612
613 if (states_are_real(sys%st)) then
614 call dcalc_properties_linear(em_vars, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%ks%xc, &
615 sys%ions, sys%outp)
616 else
617 call zcalc_properties_linear(em_vars, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%ks%xc, &
618 sys%ions, sys%outp)
619 end if
620
621 end do ! ifactor
622
623 if (states_are_real(sys%st)) then
624 call dcalc_properties_nonlinear(em_vars, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc)
625 else
626 call zcalc_properties_nonlinear(em_vars, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc)
627 end if
628
629 last_omega = em_vars%freq_factor(em_vars%nfactor) * em_vars%omega(iomega)
630
631 end do ! iomega
632
633 do idir = 1, sys%space%dim
634 do sigma = 1, em_vars%nsigma
635 do ifactor = 1, em_vars%nfactor
636 call lr_dealloc(em_vars%lr(idir, sigma, ifactor))
637 end do
638 end do
639 end do
640
641 call sternheimer_end(sh)
642 deallocate(em_vars%perturbation)
643
644 if (use_kdotp) then
645 do idir = 1, sys%space%periodic_dim
646 call lr_dealloc(kdotp_lr(idir, 1))
647 end do
648 end if
649
650 if (em_vars%calc_hyperpol .and. use_kdotp) then
651 call sternheimer_end(sh_kdotp)
652 call sternheimer_end(sh2)
653 safe_deallocate_p(pert_kdotp)
654 safe_deallocate_p(pert2_none)
655 do idir = 1, sys%space%periodic_dim
656 do idir2 = 1, sys%space%periodic_dim
657 do sigma = 1, em_vars%nsigma
658 do ifactor = 1, em_vars%nfactor
659 call lr_dealloc(kdotp_em_lr2(idir, idir2, sigma, ifactor))
660 end do
661 end do
662 end do
663 end do
664 safe_deallocate_a(kdotp_em_lr2)
665 safe_deallocate_a(dl_eig)
666 end if
667
668 if (em_vars%kpt_output) then
669 safe_deallocate_a(em_vars%alpha_k)
670 end if
671
672 if (em_vars%calc_magnetooptics .or. magnetic_pert) then
673 if (use_kdotp) then
674 safe_deallocate_p(pert2_none)
675 call sternheimer_end(sh_kmo)
676 do idir = 1, sys%space%dim
677 do idir2 = 1, sys%space%dim
678 call lr_dealloc(kb_lr(idir, idir2, 1))
679 if (idir2 <= idir) call lr_dealloc(k2_lr(idir, idir2, 1))
680 end do
681 end do
682 safe_deallocate_a(k2_lr)
683 safe_deallocate_a(kb_lr)
684 end if
685 end if
686
687 if (em_vars%calc_magnetooptics) then
688 if (.not. em_vars%magnetooptics_nohvar) call sternheimer_unset_kxc(sh_mo)
689 call sternheimer_end(sh_mo)
690 do idir = 1, sys%space%dim
691 call lr_dealloc(b_lr(idir, 1))
692 end do
693 safe_deallocate_a(b_lr)
694
695 if (use_kdotp) then
696 do idir = 1, sys%space%dim
697 do idir2 = 1, sys%space%dim
698 do sigma = 1, em_vars%nsigma
699 do ifactor = 1, nfactor_ke
700 call lr_dealloc(ke_lr(idir, idir2, sigma, ifactor))
701 end do
702 end do
703 end do
704 end do
705 safe_deallocate_a(ke_lr)
706 if (em_vars%kpt_output) then
707 safe_deallocate_a(em_vars%alpha_be_k)
708 end if
709 else
710 safe_deallocate_p(pert_b)
711 end if
712 end if
713
714 safe_deallocate_a(em_vars%omega)
715 safe_deallocate_a(em_vars%lr)
716 do ifactor = 1, em_vars%nfactor
717 call born_charges_end(em_vars%Born_charges(ifactor))
718 end do
719 call states_elec_deallocate_wfns(sys%st)
720
721 pop_sub(em_resp_run_legacy)
722 contains
723
724 ! ---------------------------------------------------------
725 subroutine parse_input()
726
727 type(block_t) :: blk
728 integer :: nrow, irow, nfreqs_in_row, ifreq, istep, perturb_type
729 real(real64) :: omega_ini, omega_fin, domega
730 logical :: freq_sort
731
733
734 call messages_obsolete_variable(sys%namespace, 'PolFreqs ', 'EMFreqs ')
735 call messages_obsolete_variable(sys%namespace, 'PolHyper ', 'EMHyperpol ')
736 call messages_obsolete_variable(sys%namespace, 'PolEta ', 'EMEta ')
737 call messages_obsolete_variable(sys%namespace, 'PolHamiltonianVariation', 'HamiltonianVariation')
738
739 !%Variable EMFreqs
740 !%Type block
741 !%Section Linear Response::Polarizabilities
742 !%Description
743 !% This block defines for which frequencies the polarizabilities
744 !% will be calculated. If it is not present, the static (<math>\omega = 0</math>) response
745 !% is calculated.
746 !%
747 !% Each row of the block indicates a sequence of frequency values, the
748 !% first column is an integer that indicates the number of steps, the
749 !% second number is the initial frequency, and the third number the final
750 !% frequency. If the first number is one, then only the initial value is
751 !% considered. The block can have any number of rows. Consider the next example:
752 !%
753 !% <tt>%EMFreqs
754 !% <br>31 | 0.0 | 1.0
755 !% <br> 1 | 0.32
756 !% <br>%</tt>
757 !%
758 !%End
759
760 if (parse_block(sys%namespace, 'EMFreqs', blk) == 0) then
761
762 nrow = parse_block_n(blk)
763 em_vars%nomega = 0
764
765 !count the number of frequencies
766 do irow = 0, nrow-1
767 call parse_block_integer(blk, irow, 0, nfreqs_in_row)
768 if (nfreqs_in_row < 1) then
769 message(1) = "EMFreqs: invalid number of frequencies."
770 call messages_fatal(1, namespace=sys%namespace)
771 end if
772 em_vars%nomega = em_vars%nomega + nfreqs_in_row
773 end do
774
775 safe_allocate(em_vars%omega(1:em_vars%nomega))
776
777 !read frequencies
778 ifreq = 1
779 do irow = 0, nrow-1
780 call parse_block_integer(blk, irow, 0, nfreqs_in_row)
781 call parse_block_float(blk, irow, 1, omega_ini)
782 if (nfreqs_in_row > 1) then
783 call parse_block_float(blk, irow, 2, omega_fin)
784 domega = (omega_fin - omega_ini)/(nfreqs_in_row - m_one)
785 do istep = 0, nfreqs_in_row-1
786 em_vars%omega(ifreq + istep) = units_to_atomic(units_inp%energy, omega_ini + domega*istep)
787 end do
788 ifreq = ifreq + nfreqs_in_row
789 else
790 em_vars%omega(ifreq) = units_to_atomic(units_inp%energy, omega_ini)
791 ifreq = ifreq + 1
792 end if
793 end do
794
795 call parse_block_end(blk)
796
797 !%Variable EMFreqsSort
798 !%Type logical
799 !%Default true
800 !%Section Linear Response::Polarizabilities
801 !%Description
802 !% If true, the frequencies specified by the <tt>EMFreqs</tt> block are sorted, so that
803 !% they are calculated in increasing order. Can be set to false to use the order as stated,
804 !% in case this makes better use of available restart information.
805 !%End
806 call parse_variable(sys%namespace, 'EMFreqsSort', .true., freq_sort)
807
808 if (freq_sort) call sort(em_vars%omega)
809
810 else
811 !there is no frequency block, we calculate response for w = 0.0
812 em_vars%nomega = 1
813 safe_allocate(em_vars%omega(1:em_vars%nomega))
814 em_vars%omega(1) = m_zero
815 end if
816
817 !%Variable EMEta
818 !%Type float
819 !%Default 0.0
820 !%Section Linear Response::Polarizabilities
821 !%Description
822 !% The imaginary part of the frequency, effectively a Lorentzian broadening
823 !% for peaks in the spectrum. It can help convergence of the SCF cycle for the
824 !% Sternheimer equation when on a resonance, and it can be used as a positive
825 !% infinitesimal to get the imaginary parts of response functions at poles.
826 !% In units of energy. Cannot be negative.
827 !%End
828
829 call parse_variable(sys%namespace, 'EMEta', m_zero, em_vars%eta, units_inp%energy)
830 if (em_vars%eta < -m_epsilon) then
831 message(1) = "EMEta cannot be negative."
832 call messages_fatal(1, namespace=sys%namespace)
833 end if
834
835 ! reset the values of these variables
836 em_vars%calc_hyperpol = .false.
837 em_vars%freq_factor(1:3) = m_one
838 em_vars%calc_magnetooptics = .false.
839 em_vars%magnetooptics_nohvar = .true.
840 em_vars%kpt_output = .false.
841
842 !%Variable EMPerturbationType
843 !%Type integer
844 !%Default electric
845 !%Section Linear Response::Polarizabilities
846 !%Description
847 !% Which perturbation to consider for electromagnetic linear response.
848 !%Option electric 1
849 !% Electric perturbation used to calculate electric polarizabilities
850 !% and hyperpolarizabilities.
851 !%Option magnetic 2
852 !% Magnetic perturbation used to calculate magnetic susceptibilities.
853 !%Option none 0
854 !% Zero perturbation, for use in testing.
855 !%End
856 call parse_variable(sys%namespace, 'EMPerturbationType', perturbation_electric, perturb_type)
857 call messages_print_var_option('EMPerturbationType', perturb_type, namespace=sys%namespace)
858
859 select case(perturb_type)
860 case(perturbation_electric)
861 em_vars%perturbation => perturbation_electric_t(sys%namespace)
863 em_vars%perturbation => perturbation_magnetic_t(sys%namespace, sys%ions)
865 em_vars%perturbation => perturbation_none_t(sys%namespace)
866 case default
867 assert(.false.)
868 end select
869
870 select type(ptr=>em_vars%perturbation)
872 !%Variable EMCalcRotatoryResponse
873 !%Type logical
874 !%Default false
875 !%Section Linear Response::Polarizabilities
876 !%Description
877 !% Calculate circular-dichroism spectrum from electric perturbation,
878 !% and write to file <tt>rotatory_strength</tt>.
879 !%End
880
881 call parse_variable(sys%namespace, 'EMCalcRotatoryResponse', .false., em_vars%calc_rotatory)
882
883 !%Variable EMHyperpol
884 !%Type block
885 !%Section Linear Response::Polarizabilities
886 !%Description
887 !% This block describes the multiples of the frequency used for
888 !% the dynamic hyperpolarizability. The results are written to the
889 !% file <tt>beta</tt> in the directory for the first multiple.
890 !% There must be three factors, summing to zero: <math>\omega_1 + \omega_2 + \omega_3 = 0</math>.
891 !% For example, for second-harmonic generation, you could use
892 !% <tt>1 | 1 | -2</tt>.
893 !%End
894
895 if (parse_block(sys%namespace, 'EMHyperpol', blk) == 0) then
896 call parse_block_float(blk, 0, 0, em_vars%freq_factor(1))
897 call parse_block_float(blk, 0, 1, em_vars%freq_factor(2))
898 call parse_block_float(blk, 0, 2, em_vars%freq_factor(3))
899
900 call parse_block_end(blk)
901
902 if (abs(sum(em_vars%freq_factor(1:3))) > m_epsilon) then
903 message(1) = "Frequency factors specified by EMHyperpol must sum to zero."
904 call messages_fatal(1, namespace=sys%namespace)
905 end if
906
907 em_vars%calc_hyperpol = .true.
908 end if
909
910 !%Variable EMCalcMagnetooptics
911 !%Type logical
912 !%Default false
913 !%Section Linear Response::Polarizabilities
914 !%Description
915 !% Calculate magneto-optical response.
916 !%End
917 call parse_variable(sys%namespace, 'EMCalcMagnetooptics', .false., em_vars%calc_magnetooptics)
918
919 !%Variable EMMagnetoopticsNoHVar
920 !%Type logical
921 !%Default true
922 !%Section Linear Response::Polarizabilities
923 !%Description
924 !% Exclude corrections to the exchange-correlation and Hartree terms
925 !% from consideration of perturbations induced by a magnetic field
926 !%End
927 call parse_variable(sys%namespace, 'EMMagnetoopticsNoHVar', .true., em_vars%magnetooptics_nohvar)
928
929 !%Variable EMKPointOutput
930 !%Type logical
931 !%Default false
932 !%Section Linear Response::Polarizabilities
933 !%Description
934 !% Give in the output contributions of different k-points to the dielectric constant.
935 !% Can be also used for magneto-optical effects.
936 !%End
937
938 call parse_variable(sys%namespace, 'EMKPointOutput', .false., em_vars%kpt_output)
939
940 end select
941
942 !%Variable EMForceNoKdotP
943 !%Type logical
944 !%Default false
945 !%Section Linear Response::Polarizabilities
946 !%Description
947 !% If the system is periodic, by default wavefunctions from a previous <tt>kdotp</tt> run will
948 !% be read, to be used in the formulas for the polarizability and
949 !% hyperpolarizability in the quantum theory of polarization. For testing purposes,
950 !% you can set this variable to true to disregard the <tt>kdotp</tt> run, and use the formulas
951 !% for the finite system. This variable has no effect for a finite system.
952 !%End
953
954 call parse_variable(sys%namespace, 'EMForceNoKdotP', .false., em_vars%force_no_kdotp)
955
956 !%Variable EMCalcBornCharges
957 !%Type logical
958 !%Default false
959 !%Section Linear Response::Polarizabilities
960 !%Description
961 !% Calculate linear-response Born effective charges from electric perturbation (experimental).
962 !%End
963
964 call parse_variable(sys%namespace, 'EMCalcBornCharges', .false., em_vars%calc_Born)
965 if (em_vars%calc_Born) call messages_experimental("Calculation of Born effective charges", namespace=sys%namespace)
966
967 !%Variable EMOccupiedResponse
968 !%Type logical
969 !%Default false
970 !%Section Linear Response::Polarizabilities
971 !%Description
972 !% Solve for full response without projector into unoccupied subspace.
973 !% Not possible if there are partial occupations.
974 !% When <tt>EMHyperpol</tt> is set for a periodic system, this variable is ignored and
975 !% the full response is always calculated.
976 !%End
977
978 call parse_variable(sys%namespace, 'EMOccupiedResponse', .false., em_vars%occ_response)
979 if (em_vars%occ_response .and. .not. (smear_is_semiconducting(sys%st%smear) .or. sys%st%smear%method == smear_fixed_occ)) then
980 message(1) = "EMOccupiedResponse cannot be used if there are partial occupations."
981 call messages_fatal(1, namespace=sys%namespace)
982 end if
983
984 !%Variable EMWavefunctionsFromScratch
985 !%Type logical
986 !%Default false
987 !%Section Linear Response::Polarizabilities
988 !%Description
989 !% Do not use saved linear-response wavefunctions from a previous run as starting guess.
990 !% Instead initialize to zero as in <tt>FromScratch</tt>, but restart densities will still
991 !% be used. Restart wavefunctions from a very different frequency can hinder convergence.
992 !%End
993
994 call parse_variable(sys%namespace, 'EMWavefunctionsFromScratch', .false., em_vars%wfns_from_scratch)
995
997
998 end subroutine parse_input
999
1000
1001 ! ---------------------------------------------------------
1002 subroutine info()
1003
1004 push_sub(em_resp_run_legacy.info)
1005
1006 call em_vars%perturbation%info()
1007 select type(ptr=>em_vars%perturbation)
1009 if (em_vars%calc_hyperpol) then
1010 call messages_print_with_emphasis(msg='Linear-Response First-Order Hyperpolarizabilities', namespace=sys%namespace)
1011 else
1012 call messages_print_with_emphasis(msg='Linear-Response Polarizabilities', namespace=sys%namespace)
1013 end if
1014 class default
1015 call messages_print_with_emphasis(msg='Magnetic Susceptibilities', namespace=sys%namespace)
1016 end select
1017
1018 if (states_are_real(sys%st)) then
1019 message(1) = 'Wavefunctions type: Real'
1020 else
1021 message(1) = 'Wavefunctions type: Complex'
1022 end if
1023 call messages_info(1, namespace=sys%namespace)
1024
1025 write(message(1),'(a,i3,a)') 'Calculating response for ', em_vars%nomega, ' frequencies.'
1026 call messages_info(1, namespace=sys%namespace)
1027
1028 call messages_print_with_emphasis(namespace=sys%namespace)
1029
1030 pop_sub(em_resp_run_legacy.info)
1031
1032 end subroutine info
1033
1034! Note: unlike the typical usage, here the templates make 'internal procedures'
1035#include "undef.F90"
1036#include "real.F90"
1037#include "em_resp_inc.F90"
1038
1039#include "undef.F90"
1040#include "complex.F90"
1041#include "em_resp_inc.F90"
1042
1043#include "undef.F90"
1044 end subroutine em_resp_run_legacy
1045
1046
1047 ! ---------------------------------------------------------
1048 subroutine em_resp_output(st, namespace, space, gr, hm, ions, outp, sh, em_vars, iomega, ifactor)
1049 type(states_elec_t), intent(inout) :: st
1050 type(namespace_t), intent(in) :: namespace
1051 class(space_t), intent(in) :: space
1052 type(grid_t), intent(in) :: gr
1053 type(hamiltonian_elec_t), intent(inout) :: hm
1054 type(ions_t), intent(in) :: ions
1055 type(output_t), intent(in) :: outp
1056 type(sternheimer_t), intent(in) :: sh
1057 type(em_resp_t), intent(inout) :: em_vars
1058 integer, intent(in) :: iomega
1059 integer, intent(in) :: ifactor
1060
1061 integer :: iunit, idir
1062 character(len=80) :: dirname, str_tmp
1063 logical :: use_kdotp
1064 complex(real64) :: epsilon(space%dim, space%dim)
1065
1066 push_sub(em_resp_output)
1067
1068 use_kdotp = space%is_periodic() .and. .not. em_vars%force_no_kdotp
1069
1070 str_tmp = freq2str(units_from_atomic(units_out%energy, em_vars%freq_factor(ifactor)*em_vars%omega(iomega)))
1071 if (em_vars%calc_magnetooptics) str_tmp = freq2str(units_from_atomic(units_out%energy, em_vars%omega(iomega)))
1072 write(dirname, '(a, a)') em_resp_dir//'freq_', trim(str_tmp)
1073 call io_mkdir(trim(dirname), namespace)
1074
1075 if (sh%has_photons .and. mpi_grp_is_root(mpi_world)) then
1076 iunit = io_open(trim(dirname)//'/photon_coord_q', namespace, action='write')
1077 write(iunit, '(a)') 'Photon coordinate Q [', trim(units_abbrev(units_out%energy)), ']'
1078 write(iunit, '(a)') ' Re Im'
1079 do idir = 1, space%dim
1080 write(iunit, '(f20.6,f20.6)') units_from_atomic(units_out%energy, sum(real(sh%zphoton_coord_q(:, idir)))), &
1081 units_from_atomic(units_out%energy, sum(aimag(sh%zphoton_coord_q(:, idir))))
1082 end do
1083 call io_close(iunit)
1084 end if
1085
1086 call write_eta()
1087
1088 select type(ptr=>em_vars%perturbation)
1090 if ((.not. em_vars%calc_magnetooptics) .or. ifactor == 1) then
1091 call out_polarizability()
1092 if (em_vars%calc_Born) then
1093 call born_output_charges(em_vars%born_charges(ifactor), ions%atom, ions%charge, ions%natoms, &
1094 namespace, space%dim, dirname, write_real = em_vars%eta < m_epsilon)
1095 end if
1096
1097 if (space%periodic_dim == space%dim) then
1099 end if
1100
1101 if ((.not. space%is_periodic() .or. em_vars%force_no_kdotp) .and. em_vars%calc_rotatory) then
1103 end if
1104
1105 else
1106 call out_magnetooptics()
1107 if (iomega == 1) call out_susceptibility()
1108 end if
1109
1111 call out_susceptibility()
1112 end select
1113
1115
1116 pop_sub(em_resp_output)
1117
1118 contains
1119
1120
1121 subroutine write_eta()
1122 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
1123
1124 push_sub(em_resp_output.write_eta)
1125
1126 iunit = io_open(trim(dirname)//'/eta', namespace, action='write')
1127
1128 write(iunit, '(3a)') 'Imaginary part of frequency [', trim(units_abbrev(units_out%energy)), ']'
1129 write(iunit, '(f20.6)') units_from_atomic(units_out%energy, em_vars%eta)
1130
1131 call io_close(iunit)
1132
1133 pop_sub(em_resp_output.write_eta)
1134 end subroutine write_eta
1135
1136
1137 ! ---------------------------------------------------------
1139 subroutine cross_section_header(out_file)
1140 integer, intent(in) :: out_file
1141
1142 character(len=80) :: header_string
1143 integer :: ii, idir, kdir
1144
1145 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
1146
1148
1149 !this header is the same as spectrum.F90
1150 write(out_file, '(a1, a20)', advance = 'no') '#', str_center("Energy", 20)
1151 write(out_file, '(a20)', advance = 'no') str_center("(1/3)*Tr[sigma]", 20)
1152 write(out_file, '(a20)', advance = 'no') str_center("Anisotropy[sigma]", 20)
1153
1154 do idir = 1, space%dim
1155 do kdir = 1, space%dim
1156 write(header_string,'(a6,i1,a1,i1,a1)') 'sigma(', idir, ',', kdir, ')'
1157 write(out_file, '(a20)', advance = 'no') str_center(trim(header_string), 20)
1158 end do
1159 end do
1160
1161 write(out_file, *)
1162 write(out_file, '(a1,a20)', advance = 'no') '#', str_center('['//trim(units_abbrev(units_out%energy)) // ']', 20)
1163 do ii = 1, 2 + space%dim**2
1164 write(out_file, '(a20)', advance = 'no') str_center('['//trim(units_abbrev(units_out%length**2)) // ']', 20)
1165 end do
1166 write(out_file,*)
1167
1169 end subroutine cross_section_header
1170
1171
1172 ! ---------------------------------------------------------
1173 subroutine out_polarizability()
1174 real(real64) :: cross(space%dim, space%dim), crossp(space%dim, space%dim)
1175 real(real64) :: cross_sum, crossp_sum, anisotropy
1176 integer :: idir, idir2
1177
1178 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
1179
1181
1182 iunit = io_open(trim(dirname)//'/alpha', namespace, action='write')
1183
1184 if (.not. em_vars%ok(ifactor)) write(iunit, '(a)') "# WARNING: not converged"
1185
1186 write(iunit, '(3a)') '# Polarizability tensor [', trim(units_abbrev(units_out%polarizability)), ']'
1187 call output_tensor(real(em_vars%alpha(:, :, ifactor), real64), space%dim, units_out%polarizability, iunit=iunit)
1188
1189 call io_close(iunit)
1190
1191 ! CROSS SECTION (THE IMAGINARY PART OF POLARIZABILITY)
1192 if (em_vars%eta > m_epsilon) then
1193 cross(1:space%dim, 1:space%dim) = aimag(em_vars%alpha(1:space%dim, 1:space%dim, ifactor)) * &
1194 em_vars%freq_factor(ifactor) * em_vars%omega(iomega) * (m_four * m_pi / p_c)
1196 do idir = 1, space%dim
1197 do idir2 = 1, space%dim
1198 cross(idir, idir2) = units_from_atomic(units_out%length**2, cross(idir, idir2))
1199 end do
1200 end do
1201
1202 iunit = io_open(trim(dirname)//'/cross_section', namespace, action='write')
1203 if (.not. em_vars%ok(ifactor)) write(iunit, '(a)') "# WARNING: not converged"
1204
1205 crossp(1:space%dim, 1:space%dim) = matmul(cross(1:space%dim, 1:space%dim), cross(1:space%dim, 1:space%dim))
1206
1207 cross_sum = m_zero
1208 crossp_sum = m_zero
1209 do idir = 1, space%dim
1210 cross_sum = cross_sum + cross(idir, idir)
1211 crossp_sum = crossp_sum + crossp(idir, idir)
1212 end do
1213
1214 anisotropy = crossp_sum - m_third * cross_sum**2
1215
1216 call cross_section_header(iunit)
1217 write(iunit,'(3e20.8)', advance = 'no') &
1218 units_from_atomic(units_out%energy, em_vars%freq_factor(ifactor)*em_vars%omega(iomega)), &
1219 cross_sum * m_third, sqrt(max(anisotropy, m_zero))
1220 do idir = 1, space%dim
1221 do idir2 = 1, space%dim
1222 write(iunit,'(e20.8)', advance = 'no') cross(idir, idir2)
1223 end do
1224 end do
1225 write(iunit,'(a)', advance = 'yes')
1226
1227 call io_close(iunit)
1228 end if
1229
1231 end subroutine out_polarizability
1232
1233
1234 ! ---------------------------------------------------------
1236 subroutine out_dielectric_constant()
1237 integer :: idir, idir1, ik
1238 character(len=80) :: header_string
1239 complex(real64), allocatable :: epsilon_k(:, :, :)
1240
1241 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
1242
1244
1245 iunit = io_open(trim(dirname)//'/epsilon', namespace, action='write')
1246 if (.not. em_vars%ok(ifactor)) write(iunit, '(a)') "# WARNING: not converged"
1247
1248 epsilon(1:space%dim, 1:space%dim) = &
1249 4 * m_pi * em_vars%alpha(1:space%dim, 1:space%dim, ifactor) / ions%latt%rcell_volume
1250 do idir = 1, space%dim
1251 epsilon(idir, idir) = epsilon(idir, idir) + m_one
1252 end do
1253
1254 write(iunit, '(a)') '# Real part of dielectric constant'
1255 call output_tensor(real(epsilon(1:space%dim, 1:space%dim), real64), space%dim, unit_one, iunit=iunit)
1256 write(iunit, '(a)')
1257 write(iunit, '(a)') '# Imaginary part of dielectric constant'
1258 call output_tensor(aimag(epsilon(1:space%dim, 1:space%dim)), space%dim, unit_one, iunit=iunit)
1259
1260 if (em_vars%lrc_kernel) then
1261 write(iunit, '(a)')
1262 write(iunit, '(a)') '# Without G = G'' = 0 term of the LRC kernel'
1263
1264 epsilon(1:space%dim, 1:space%dim) = &
1265 4 * m_pi * em_vars%alpha0(1:space%dim, 1:space%dim, ifactor) / ions%latt%rcell_volume
1266 do idir = 1, space%dim
1267 epsilon(idir, idir) = epsilon(idir, idir) + m_one
1268 end do
1269
1270 write(iunit, '(a)') '# Real part of dielectric constant'
1271 call output_tensor(real(epsilon(1:space%dim, 1:space%dim), real64), space%dim, unit_one, iunit=iunit)
1272 write(iunit, '(a)')
1273 write(iunit, '(a)') '# Imaginary part of dielectric constant'
1274 call output_tensor(aimag(epsilon(1:space%dim, 1:space%dim)), space%dim, unit_one, iunit=iunit)
1275 end if
1276
1277 call io_close(iunit)
1278
1279 if (em_vars%kpt_output) then
1280 safe_allocate(epsilon_k(1:space%dim, 1:space%dim, 1:hm%kpoints%reduced%npoints))
1281 do ik = 1, hm%kpoints%reduced%npoints
1282 do idir = 1, space%dim
1283 do idir1 = 1, space%dim
1284 epsilon_k(idir, idir1, ik) = m_four * m_pi * em_vars%alpha_k(idir, idir1, ifactor, ik) / ions%latt%rcell_volume
1285 end do
1286 end do
1287 end do
1288 iunit = io_open(trim(dirname)//'/epsilon_k_re', namespace, action='write')
1289
1290 write(iunit, '(a)') '# Real part of dielectric constant'
1291 write(iunit, '(a10)', advance = 'no') '# index '
1292 write(iunit, '(a20)', advance = 'no') str_center("weight", 20)
1293 write(iunit, '(a20)', advance = 'no') str_center("kx", 20)
1294 write(iunit, '(a20)', advance = 'no') str_center("ky", 20)
1295 write(iunit, '(a20)', advance = 'no') str_center("kz", 20)
1296
1297 do idir = 1, space%dim
1298 do idir1 = 1, space%dim
1299 write(header_string,'(a7,i1,a1,i1,a1)') 'Re eps(', idir, ',', idir1, ')'
1300 write(iunit, '(a20)', advance = 'no') str_center(trim(header_string), 20)
1301 end do
1302 end do
1303 write(iunit, *)
1304
1305 do ik = 1, hm%kpoints%reduced%npoints
1306 write(iunit, '(i8)', advance = 'no') ik
1307 write(iunit, '(e20.8)', advance = 'no') hm%kpoints%reduced%weight(ik)
1308 do idir = 1, space%dim
1309 write(iunit, '(e20.8)', advance = 'no') hm%kpoints%reduced%red_point(idir, ik)
1310 end do
1311 do idir = 1, space%dim
1312 do idir1 = 1, space%dim
1313 write(iunit, '(e20.8)', advance = 'no') real(epsilon_k(idir, idir1, ik), real64)
1314 end do
1315 end do
1316 write(iunit, *)
1317 end do
1318 call io_close(iunit)
1319
1320 iunit = io_open(trim(dirname)//'/epsilon_k_im', namespace, action='write')
1321
1322 write(iunit, '(a)') '# Imaginary part of dielectric constant'
1323 write(iunit, '(a10)', advance = 'no') '# index '
1324 write(iunit, '(a20)', advance = 'no') str_center("weight", 20)
1325 write(iunit, '(a20)', advance = 'no') str_center("kx", 20)
1326 write(iunit, '(a20)', advance = 'no') str_center("ky", 20)
1327 write(iunit, '(a20)', advance = 'no') str_center("kz", 20)
1328
1329 do idir = 1, space%dim
1330 do idir1 = 1, space%dim
1331 write(header_string,'(a7,i1,a1,i1,a1)') 'Im eps(', idir, ',', idir1,')'
1332 write(iunit, '(a20)', advance = 'no') str_center(trim(header_string), 20)
1333 end do
1334 end do
1335 write(iunit, *)
1336
1337 do ik = 1, hm%kpoints%reduced%npoints
1338 write(iunit, '(i8)', advance = 'no') ik
1339 write(iunit, '(e20.8)', advance = 'no') hm%kpoints%reduced%weight(ik)
1340 do idir = 1, space%dim
1341 write(iunit, '(e20.8)', advance = 'no') hm%kpoints%reduced%red_point(idir, ik)
1342 end do
1343 do idir = 1, space%dim
1344 do idir1 = 1, space%dim
1345 write(iunit, '(e20.8)', advance = 'no') aimag(epsilon_k(idir, idir1, ik))
1346 end do
1347 end do
1348 write(iunit, *)
1349 end do
1350 call io_close(iunit)
1351 safe_deallocate_a(epsilon_k)
1352 end if
1353
1355 end subroutine out_dielectric_constant
1356
1357
1358 ! ---------------------------------------------------------
1359 subroutine out_susceptibility()
1360
1361 character(len=80) :: dirname1
1362
1363 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
1364
1366
1367 select type(ptr=>em_vars%perturbation)
1369 write(dirname1, '(a)') em_resp_dir//'freq_0.0000'
1370 call io_mkdir(trim(dirname1), namespace)
1371 iunit = io_open(trim(dirname1)//'/susceptibility', namespace, action='write')
1372 class default
1373 iunit = io_open(trim(dirname)//'/susceptibility', namespace, action='write')
1374 end select
1375
1376 if (.not. em_vars%ok(ifactor)) write(iunit, '(a)') "# WARNING: not converged"
1377
1378 ! There is no separation into the diamagnetic and paramagnetic terms in the expression
1379 ! for periodic systems
1380 if (.not. use_kdotp) then
1381 write(iunit, '(2a)') '# Paramagnetic contribution to the susceptibility tensor [ppm a.u.]'
1382 call output_tensor(real(em_vars%chi_para(:, :), real64), space%dim, unit_ppm, iunit=iunit)
1383 write(iunit, '(1x)')
1384
1385 write(iunit, '(2a)') '# Diamagnetic contribution to the susceptibility tensor [ppm a.u.]'
1386 call output_tensor(real(em_vars%chi_dia(:, :), real64), space%dim, unit_ppm, iunit=iunit)
1387 write(iunit, '(1x)')
1388 end if
1389
1390 write(iunit, '(2a)') '# Total susceptibility tensor [ppm a.u.]'
1391 call output_tensor(real(em_vars%chi_para(:, :) + em_vars%chi_dia(:,:), real64) , &
1392 space%dim, unit_ppm, iunit=iunit)
1393 write(iunit, '(1x)')
1394
1395 write(iunit, '(a)') hyphens
1396
1397 if (.not. use_kdotp) then
1398 write(iunit, '(2a)') '# Paramagnetic contribution to the susceptibility tensor [ppm cgs / mol]'
1399 call output_tensor(real(em_vars%chi_para(:, :), real64), space%dim, unit_susc_ppm_cgs, iunit=iunit)
1400 write(iunit, '(1x)')
1401
1402 write(iunit, '(2a)') '# Diamagnetic contribution to the susceptibility tensor [ppm cgs / mol]'
1403 call output_tensor(real(em_vars%chi_dia(:, :), real64), space%dim, unit_susc_ppm_cgs, iunit=iunit)
1404 write(iunit, '(1x)')
1405 end if
1406
1407 write(iunit, '(2a)') '# Total susceptibility tensor [ppm cgs / mol]'
1408 call output_tensor(real(em_vars%chi_para(:, :) + em_vars%chi_dia(:,:), real64), &
1409 space%dim, unit_susc_ppm_cgs, iunit=iunit)
1410 write(iunit, '(1x)')
1411
1412 if (use_kdotp) then
1413 write(iunit, '(a)') hyphens
1414 write(iunit, '(1a)') '# Magnetization [ppm a.u.]'
1415 write(iunit, '(3f20.8)') units_from_atomic(unit_ppm, real(em_vars%magn(1), real64)), &
1416 units_from_atomic(unit_ppm, real(em_vars%magn(2), real64)), &
1417 units_from_atomic(unit_ppm, real(em_vars%magn(3), real64))
1418 end if
1419
1420 call io_close(iunit)
1422 end subroutine out_susceptibility
1423
1424
1425 ! ---------------------------------------------------------
1426 subroutine out_wfn_and_densities()
1427 integer :: idir, isigma
1428
1430
1431 do idir = 1, space%dim
1432 if (states_are_complex(st)) then
1433
1434 if (space%dim == 3 .and. outp%what(option__output__elf)) then
1435 if (em_vars%nsigma == 1) then
1436 call zlr_calc_elf(st, space, gr, hm%kpoints, em_vars%lr(idir, 1, ifactor))
1437 else
1438 call zlr_calc_elf(st, space, gr, hm%kpoints, em_vars%lr(idir, 1, ifactor), em_vars%lr(idir, 2, ifactor))
1439 end if
1440 end if
1441 do isigma = 1, em_vars%nsigma
1442 call zoutput_lr(outp, namespace, space, dirname, st, gr, em_vars%lr(idir, isigma, ifactor), idir, isigma, ions, &
1443 units_out%force)
1444 end do
1445 else
1446
1447 if (space%dim == 3 .and. outp%what(option__output__elf)) then
1448 if (em_vars%nsigma == 1) then
1449 call dlr_calc_elf(st, space, gr, hm%kpoints, em_vars%lr(idir, 1, ifactor))
1450 else
1451 call dlr_calc_elf(st, space, gr, hm%kpoints, em_vars%lr(idir, 1, ifactor), em_vars%lr(idir, 2, ifactor))
1452 end if
1453 end if
1454
1455 do isigma = 1, em_vars%nsigma
1456 call doutput_lr(outp, namespace, space, dirname, st, gr, em_vars%lr(idir, isigma, ifactor), idir, isigma, ions, &
1457 units_out%force)
1458 end do
1459
1460 end if
1461 end do
1462
1464
1465 end subroutine out_wfn_and_densities
1466
1467
1468 ! ---------------------------------------------------------
1471 subroutine out_circular_dichroism()
1472
1473 type(perturbation_magnetic_t), pointer :: angular_momentum
1474 integer :: idir
1475 real(real64) :: ff
1476 complex(real64) :: dic
1477 complex(real64), allocatable :: psi(:, :, :, :)
1478
1480
1481 if (states_are_complex(st) .and. em_vars%nsigma == 2) then
1482
1483 message(1) = "Info: Calculating rotatory response."
1484 call messages_info(1, namespace=namespace)
1485
1486 angular_momentum => perturbation_magnetic_t(namespace, ions)
1487
1488 safe_allocate(psi(1:gr%np_part, 1:st%d%dim, st%st_start:st%st_end, st%d%kpt%start:st%d%kpt%end))
1489
1490 call states_elec_get_state(st, gr, psi)
1491
1492 dic = m_zero
1493 do idir = 1, space%dim
1494 call angular_momentum%setup_dir(idir)
1495 dic = dic &
1496 + angular_momentum%zexpectation_value(namespace, space, gr, hm, st, psi, &
1497 em_vars%lr(idir, 1, ifactor)%zdl_psi) &
1498 + angular_momentum%zexpectation_value(namespace, space, gr, hm, st, &
1499 em_vars%lr(idir, 2, ifactor)%zdl_psi, psi)
1500 end do
1501
1502 safe_deallocate_a(psi)
1503
1504 safe_deallocate_p(angular_momentum)
1505
1506 dic = dic*m_zi*m_half
1507
1508 if (mpi_grp_is_root(mpi_world)) then ! only first node outputs
1509
1510 iunit = io_open(trim(dirname)//'/rotatory_strength', namespace, action='write')
1511
1512 ! print header
1513 write(iunit, '(a1,a20,a20,a20)') '#', str_center("Energy", 20), str_center("R", 20), str_center("Re[beta]", 20)
1514 write(iunit, '(a1,a20,a20,a20)') '#', str_center('['//trim(units_abbrev(units_out%energy)) // ']', 20), &
1515 str_center('['//trim(units_abbrev(units_out%length**3)) //']', 20), &
1516 str_center('['//trim(units_abbrev(units_out%length**4)) //']', 20)
1517
1518 ff = m_zero
1519 if (abs(em_vars%omega(iomega)) > m_epsilon) ff = real(dic, real64) /(m_three*em_vars%omega(iomega))
1520
1521 write(iunit, '(3e20.8)') units_from_atomic(units_out%energy, em_vars%omega(iomega)), &
1522 units_from_atomic(units_out%length**3, aimag(dic)/(p_c*m_pi)), units_from_atomic(units_out%length**4, ff)
1523
1524 call io_close(iunit)
1525 end if
1526 end if
1527
1529
1530 end subroutine out_circular_dichroism
1531
1532 ! ---------------------------------------------------------
1533 subroutine out_magnetooptics
1534 integer :: idir, ik
1535 complex(real64) :: epsilon_m(4), diff(4), eps_mk(space%dim)
1536
1537 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
1538
1540
1541 ! This code assumes 3D
1542 assert(space%dim == 3)
1543
1544 diff(:) = m_zero
1545 do idir = 1, space%dim
1546 diff(idir) = m_half * (em_vars%alpha_be(magn_dir(idir, 1), magn_dir(idir, 2), idir) - &
1547 em_vars%alpha_be(magn_dir(idir, 2), magn_dir(idir, 1), idir))
1548 end do
1549 diff(4) = (diff(1) + diff(2) + diff(3)) / m_three
1550
1551 iunit = io_open(trim(dirname)//'/alpha_mo', namespace, action='write')
1552
1553 if (.not. em_vars%ok(ifactor)) write(iunit, '(a)') "# WARNING: not converged"
1554
1555 write(iunit, '(a1, a25)', advance = 'no') '#', str_center(" ", 25)
1556 write(iunit, '(a20)', advance = 'no') str_center(" yz,x = -zy,x", 20)
1557 write(iunit, '(a20)', advance = 'no') str_center(" zx,y = -xz,y", 20)
1558 write(iunit, '(a20)', advance = 'no') str_center(" xy,z = -yx,z", 20)
1559 write(iunit, '(a20)', advance = 'no') str_center(" Average", 20)
1560 write(iunit, *)
1561
1562 write(iunit, '(a25)', advance = 'no') str_center("Re alpha [a.u.]", 25)
1563 do idir = 1, space%dim + 1
1564 write(iunit, '(e20.8)', advance = 'no') real(diff(idir), real64)
1565 end do
1566 write(iunit, *)
1567
1568 write(iunit, '(a25)', advance = 'no') str_center("Im alpha [a.u.]", 25)
1569 do idir = 1, space%dim + 1
1570 write(iunit, '(e20.8)', advance = 'no') aimag(diff(idir))
1571 end do
1572 write(iunit, *)
1573
1574 if (space%is_periodic()) then
1575 ! This code assumes 3D periodic
1576 assert(space%periodic_dim == 3)
1577
1578 do idir = 1, space%dim
1579 epsilon_m(idir) = 4 * m_pi * diff(idir) / ions%latt%rcell_volume
1580 end do
1581 epsilon_m(4) = 4 * m_pi * diff(4) / ions%latt%rcell_volume
1582
1583 write(iunit, '(a25)', advance = 'no') str_center("Re epsilon (B = 1 a.u.)", 25)
1584 do idir = 1, space%dim + 1
1585 write(iunit, '(e20.8)', advance = 'no') real(epsilon_m(idir), real64)
1586 end do
1587 write(iunit, *)
1588
1589 write(iunit, '(a25)', advance = 'no') str_center("Im epsilon (B = 1 a.u.)", 25)
1590 do idir = 1, space%dim + 1
1591 write(iunit, '(e20.8)', advance = 'no') aimag(epsilon_m(idir))
1592 end do
1593 write(iunit, *)
1594
1595 if (em_vars%lrc_kernel) then
1596 write(iunit, '(a)')
1597 write(iunit, '(a)') '# Without the G = G'' = 0 term of the LRC kernel'
1598
1599 diff(:) = m_zero
1600 epsilon_m(:) = m_zero
1601 do idir = 1, space%dim
1602 diff(idir) = m_half * (em_vars%alpha_be0(magn_dir(idir, 1), magn_dir(idir, 2), idir) - &
1603 em_vars%alpha_be0(magn_dir(idir, 2), magn_dir(idir, 1), idir))
1604
1605 epsilon_m(idir) = 4 * m_pi * diff(idir) / ions%latt%rcell_volume
1606 end do
1607 diff(4) = (diff(1) + diff(2) + diff(3)) / m_three
1608 epsilon_m(4) = 4 * m_pi * diff(4) / ions%latt%rcell_volume
1609
1610 write(iunit, '(a1, a25)', advance = 'no') '#', str_center(" ", 25)
1611 write(iunit, '(a20)', advance = 'no') str_center(" yz,x = -zy,x", 20)
1612 write(iunit, '(a20)', advance = 'no') str_center(" zx,y = -xz,y", 20)
1613 write(iunit, '(a20)', advance = 'no') str_center(" xy,z = -yx,z", 20)
1614 write(iunit, '(a20)', advance = 'no') str_center(" Average", 20)
1615 write(iunit, *)
1616
1617 write(iunit, '(a25)', advance = 'no') str_center("Re alpha [a.u.]", 25)
1618 do idir = 1, space%dim + 1
1619 write(iunit, '(e20.8)', advance = 'no') real(diff(idir), real64)
1620 end do
1621 write(iunit, *)
1622
1623 write(iunit, '(a25)', advance = 'no') str_center("Im alpha [a.u.]", 25)
1624 do idir = 1, space%dim + 1
1625 write(iunit, '(e20.8)', advance = 'no') aimag(diff(idir))
1626 end do
1627 write(iunit, *)
1628
1629 write(iunit, '(a25)', advance = 'no') str_center("Re epsilon (B = 1 a.u.)", 25)
1630 do idir = 1, space%dim + 1
1631 write(iunit, '(e20.8)', advance = 'no') real(epsilon_m(idir), real64)
1632 end do
1633 write(iunit, *)
1634
1635 write(iunit, '(a25)', advance = 'no') str_center("Im epsilon (B = 1 a.u.)", 25)
1636 do idir = 1, space%dim + 1
1637 write(iunit, '(e20.8)', advance = 'no') aimag(epsilon_m(idir))
1638 end do
1639 write(iunit, *)
1640 end if
1641 end if
1642 call io_close(iunit)
1643
1644 if (space%is_periodic() .and. em_vars%kpt_output) then
1645 iunit = io_open(trim(dirname)//'/epsilon_mo_k', namespace, action='write')
1646
1647 write(iunit, '(a)') '# Contribution to dielectric tensor for B = 1 a.u.'
1648 write(iunit, '(a10)', advance = 'no') '# index '
1649 write(iunit, '(a20)', advance = 'no') str_center("weight", 20)
1650 write(iunit, '(a20)', advance = 'no') str_center("kx", 20)
1651 write(iunit, '(a20)', advance = 'no') str_center("ky", 20)
1652 write(iunit, '(a20)', advance = 'no') str_center("kz", 20)
1653 write(iunit, '(a20)', advance = 'no') str_center("Re eps_yz,x", 20)
1654 write(iunit, '(a20)', advance = 'no') str_center("Re eps_zx,y", 20)
1655 write(iunit, '(a20)', advance = 'no') str_center("Re eps_xy,z", 20)
1656 write(iunit, '(a20)', advance = 'no') str_center("Im eps_yz,x", 20)
1657 write(iunit, '(a20)', advance = 'no') str_center("Im eps_zx,y", 20)
1658 write(iunit, '(a20)', advance = 'no') str_center("Im eps_xy,z", 20)
1659 write(iunit, *)
1660
1661 do ik = 1, hm%kpoints%reduced%npoints
1662 write(iunit, '(i8)', advance = 'no') ik
1663 write(iunit, '(e20.8)', advance = 'no') hm%kpoints%reduced%weight(ik)
1664 do idir = 1, space%dim
1665 eps_mk(idir) = m_two * m_pi * (em_vars%alpha_be_k(magn_dir(idir, 1), magn_dir(idir, 2), idir, ik) - &
1666 em_vars%alpha_be_k(magn_dir(idir, 2), magn_dir(idir, 1), idir, ik)) / ions%latt%rcell_volume
1667 end do
1668
1669 do idir = 1, space%dim
1670 write(iunit, '(e20.8)', advance = 'no') hm%kpoints%reduced%red_point(idir, ik)
1671 end do
1672 do idir = 1, space%dim
1673 write(iunit, '(e20.8)', advance = 'no') real(eps_mk(idir), real64)
1674 end do
1675 do idir = 1, space%dim
1676 write(iunit, '(e20.8)', advance = 'no') aimag(eps_mk(idir))
1677 end do
1678 write(iunit, *)
1679 end do
1680 call io_close(iunit)
1681 end if
1682
1684 end subroutine out_magnetooptics
1685
1686 end subroutine em_resp_output
1687
1688 ! ---------------------------------------------------------
1692 subroutine out_hyperpolarizability(box, beta, freq_factor, converged, dirname, namespace)
1693 class(box_t), intent(in) :: box
1694 complex(real64), intent(in) :: beta(:, :, :)
1695 real(real64), intent(in) :: freq_factor(:)
1696 logical, intent(in) :: converged
1697 character(len=*), intent(in) :: dirname
1698 type(namespace_t), intent(in) :: namespace
1699
1700 complex(real64) :: bpar(1:box%dim), bper(1:box%dim), bk(1:box%dim)
1701 complex(real64) :: HRS_VV, HRS_HV
1702 integer :: ii, jj, kk, iunit
1703
1704 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
1705
1706 push_sub(out_hyperpolarizability)
1707
1708 ! Output first hyperpolarizability (beta)
1709 iunit = io_open(trim(dirname)//'/beta', namespace, action='write')
1710
1711 if (.not. converged) write(iunit, '(a)') "# WARNING: not converged"
1712
1713 write(iunit, '(a,3(f4.1,a),2a)', advance='no') 'First hyperpolarizability tensor: beta(', &
1714 freq_factor(1), ', ', freq_factor(2), ', ', freq_factor(3), ') [', &
1715 trim(units_abbrev(units_out%hyperpolarizability)), ']'
1716
1717 write(iunit, '()')
1718
1719 do ii = 1, box%dim
1720 do jj = 1, box%dim
1721 do kk = 1, box%dim
1722 write(iunit,'(a,e20.8,e20.8)') 'beta '// &
1723 index2axis(ii)//index2axis(jj)//index2axis(kk)//' ', &
1724 units_from_atomic(units_out%hyperpolarizability, real( beta(ii, jj, kk), real64)) , &
1725 units_from_atomic(units_out%hyperpolarizability, aimag(beta(ii, jj, kk)))
1726 end do
1727 end do
1728 end do
1729
1730 if (box%dim == 3) then
1731 bpar = m_zero
1732 bper = m_zero
1733
1734 do ii = 1, box%dim
1735 do jj = 1, box%dim
1736 bpar(ii) = bpar(ii) + beta(ii, jj, jj) + beta(jj, ii, jj) + beta(jj, jj, ii)
1737 bper(ii) = bper(ii) + m_two*beta(ii, jj, jj) - m_three*beta(jj, ii, jj) + m_two*beta(jj, jj, ii)
1738 end do
1739 end do
1740
1741 write(iunit, '()')
1742
1743 bpar = bpar / m_five
1744 bper = bper / m_five
1745 bk(1:box%dim) = m_three*m_half*(bpar(1:box%dim) - bper(1:box%dim))
1746
1747 do ii = 1, box%dim
1748 write(iunit, '(a, 2e20.8)') 'beta // '//index2axis(ii), &
1749 units_from_atomic(units_out%hyperpolarizability, real(bpar(ii), real64)), &
1750 units_from_atomic(units_out%hyperpolarizability, aimag(bpar(ii)))
1751 end do
1752
1753 write(iunit, '()')
1754
1755 do ii = 1, box%dim
1756 write(iunit, '(a, 2e20.8)') 'beta _L '//index2axis(ii), &
1757 units_from_atomic(units_out%hyperpolarizability, real(bper(ii), real64)), &
1758 units_from_atomic(units_out%hyperpolarizability, aimag(bper(ii)))
1759 end do
1760
1761 write(iunit, '()')
1762
1763 do ii = 1, box%dim
1764 write(iunit, '(a, 2e20.8)') 'beta k '//index2axis(ii), &
1765 units_from_atomic(units_out%hyperpolarizability, real(bk(ii), real64)), &
1766 units_from_atomic(units_out%hyperpolarizability, aimag(bk(ii)))
1767 end do
1768
1769 call calc_beta_hrs(box, beta, hrs_vv, hrs_hv)
1770
1771 write(iunit, '()')
1772 write(iunit, '(a)') 'beta for liquid- or gas-phase hyper-Rayleigh scattering:'
1773 write(iunit, '(a, 2e20.8)') 'VV polarization ', &
1774 units_from_atomic(units_out%hyperpolarizability, real(sqrt(hrs_vv), real64)), &
1775 units_from_atomic(units_out%hyperpolarizability, aimag(sqrt(hrs_vv)))
1776 write(iunit, '(a, 2e20.8)') 'HV polarization ', &
1777 units_from_atomic(units_out%hyperpolarizability, real(sqrt(hrs_hv), real64)), &
1778 units_from_atomic(units_out%hyperpolarizability, aimag(sqrt(hrs_hv)))
1779 end if
1780
1781 call io_close(iunit)
1783
1784 contains
1785
1786 ! ---------------------------------------------------------
1791 subroutine calc_beta_hrs(box, beta, HRS_VV, HRS_HV)
1792 class(box_t), intent(in) :: box
1793 complex(real64), intent(in) :: beta(:, :, :)
1794 complex(real64), intent(out) :: HRS_VV, HRS_HV
1795
1796 complex(real64) :: HRS_A, HRS_B, HRS_C, HRS_D, HRS_E
1797 complex(real64) :: HRS_B1, HRS_B2, HRS_C1, HRS_C2, HRS_C3, HRS_D1, HRS_D2, HRS_D3, HRS_E1, HRS_E2
1798 integer :: ii, jj
1799
1801
1802 ! first calculate VV (vertical-vertical) polarization, FFF in Decius et al.
1803 hrs_a = m_zero
1804 do ii = 1, box%dim
1805 hrs_a = hrs_a + beta(ii,ii,ii)**2
1806 end do
1807
1808 hrs_b = m_zero
1809 hrs_c = m_zero
1810 do ii = 1, box%dim
1811 do jj = 1, box%dim
1812 if (ii /= jj) then
1813 hrs_b = hrs_b + beta(ii,ii,ii) * (beta(ii,jj,jj) + beta(jj,ii,jj) + beta(jj,jj,ii))
1814 hrs_c = hrs_c + (beta(ii,ii,jj) + beta(ii,jj,ii) + beta(jj,ii,ii))**2
1815 end if
1816 end do
1817 end do
1818
1819 hrs_d = (beta(1,1,2) + beta(1,2,1) + beta(2,1,1)) * (beta(2,3,3) + beta(3,2,3) + beta(3,3,2)) &
1820 + (beta(2,2,3) + beta(2,3,2) + beta(3,2,2)) * (beta(3,1,1) + beta(1,3,1) + beta(1,1,3)) &
1821 + (beta(3,3,1) + beta(3,1,3) + beta(1,3,3)) * (beta(1,2,2) + beta(2,1,2) + beta(2,2,1))
1822
1823 hrs_e = (beta(1,2,3) + beta(1,3,2) + beta(2,1,3) + beta(2,3,1) + beta(3,1,2) + beta(3,2,1))**2
1824
1825 hrs_vv = (m_one / 7.0_real64) * hrs_a &
1826 + (m_two / 35.0_real64) * hrs_b &
1827 + (m_one / 35.0_real64) * hrs_c &
1828 + (m_two / 105.0_real64) * hrs_d &
1829 + (m_one / 105.0_real64) * hrs_e
1830
1831 ! now calculate HV (horizontal-vertical) polarization, FGG in Decius et al.
1832 hrs_b1 = m_zero
1833 hrs_b2 = m_zero
1834 hrs_c1 = m_zero
1835 hrs_c2 = m_zero
1836 hrs_c3 = m_zero
1837 do ii = 1, box%dim
1838 do jj = 1, box%dim
1839 if (ii /= jj) then
1840 hrs_b1 = hrs_b1 + beta(ii,ii,ii) * beta(ii,jj,jj)
1841 hrs_b2 = hrs_b2 + beta(ii,ii,ii) * (beta(jj,ii,jj) + beta(jj,jj,ii))
1842 hrs_c1 = hrs_c1 + (beta(ii,ii,jj) + beta(ii,jj,ii))**2
1843 hrs_c2 = hrs_c2 + beta(jj,ii,ii) * (beta(ii,ii,jj) + beta(ii,jj,ii))
1844 hrs_c3 = hrs_c3 + beta(jj,ii,ii)**2
1845 end if
1846 end do
1847 end do
1848
1849 hrs_d1 = (beta(1,1,2) + beta(1,2,1) + beta(2,1,1)) * (beta(3,2,3) + beta(3,3,2)) &
1850 + (beta(2,2,3) + beta(2,3,2) + beta(3,2,2)) * (beta(1,3,1) + beta(1,1,3)) &
1851 + (beta(3,3,1) + beta(3,1,3) + beta(1,3,3)) * (beta(2,1,2) + beta(2,2,1))
1852 hrs_d2 = (beta(1,1,2) + beta(1,2,1)) * beta(2,3,3) &
1853 + (beta(2,2,3) + beta(2,3,2)) * beta(3,1,1) &
1854 + (beta(3,3,1) + beta(3,1,3)) * beta(1,2,2)
1855 hrs_d3 = beta(2,1,1) * beta(2,3,3) &
1856 + beta(3,2,2) * beta(3,1,1) &
1857 + beta(1,3,3) * beta(1,2,2)
1858
1859 hrs_e1 = (beta(1,2,3) + beta(1,3,2))**2 &
1860 + (beta(2,1,3) + beta(2,3,1))**2 &
1861 + (beta(3,1,2) + beta(3,2,1))**2
1862
1863 hrs_e2 = (beta(1,2,3) + beta(1,3,2)) * (beta(2,1,3) + beta(2,3,1)) &
1864 + (beta(2,1,3) + beta(2,3,1)) * (beta(3,1,2) + beta(3,2,1)) &
1865 + (beta(3,1,2) + beta(3,2,1)) * (beta(1,2,3) + beta(1,3,2))
1866
1867 hrs_hv = (m_one / 35.0_real64) * hrs_a &
1868 + (m_four / 105.0_real64) * hrs_b1 &
1869 - (m_one / 35.0_real64) * hrs_b2 &
1870 + (m_two / 105.0_real64) * hrs_c1 &
1871 - (m_one / 35.0_real64) * hrs_c2 &
1872 + (m_three / 35.0_real64) * hrs_c3 &
1873 - (m_one / 105.0_real64) * hrs_d1 &
1874 - (m_one / 105.0_real64) * hrs_d2 &
1875 + (m_two / 35.0_real64) * hrs_d3 &
1876 + (m_one / 35.0_real64) * hrs_e1 &
1877 - (m_one / 105.0_real64) * hrs_e2
1878
1880 end subroutine calc_beta_hrs
1881
1882 end subroutine out_hyperpolarizability
1883
1884end module em_resp_oct_m
1885
1886!! Local Variables:
1887!! mode: f90
1888!! coding: utf-8
1889!! End:
subroutine out_dielectric_constant()
epsilon = 1 + 4 * pi * alpha/volume
Definition: em_resp.F90:2803
subroutine write_eta()
Definition: em_resp.F90:2688
subroutine calc_beta_hrs(box, beta, HRS_VV, HRS_HV)
calculate hyper-Rayleigh scattering hyperpolarizabilities SJ Cyvin, JE Rauch, and JC Decius,...
Definition: em_resp.F90:3358
subroutine out_circular_dichroism()
See D Varsano, LA Espinosa Leal, Xavier Andrade, MAL Marques, Rosa di Felice, Angel Rubio,...
Definition: em_resp.F90:3038
subroutine out_magnetooptics
Definition: em_resp.F90:3100
subroutine out_susceptibility()
Definition: em_resp.F90:2926
subroutine dcalc_properties_linear(em_vars, namespace, space, gr, kpoints, st, hm, xc, ions, outp)
Definition: em_resp.F90:1691
subroutine out_wfn_and_densities()
Definition: em_resp.F90:2993
subroutine drun_sternheimer(em_vars, namespace, space, gr, kpoints, st, hm, mc, ions)
Definition: em_resp.F90:1196
subroutine out_polarizability()
Definition: em_resp.F90:2740
subroutine dcalc_properties_nonlinear(em_vars, namespace, space, gr, st, hm, xc)
Definition: em_resp.F90:1820
subroutine zrun_sternheimer(em_vars, namespace, space, gr, kpoints, st, hm, mc, ions)
Definition: em_resp.F90:1925
subroutine parse_input()
Definition: em_resp.F90:819
subroutine cross_section_header(out_file)
Note: this should be in spectrum.F90.
Definition: em_resp.F90:2706
subroutine zcalc_properties_nonlinear(em_vars, namespace, space, gr, st, hm, xc)
Definition: em_resp.F90:2549
subroutine zcalc_properties_linear(em_vars, namespace, space, gr, kpoints, st, hm, xc, ions, outp)
Definition: em_resp.F90:2420
This is the common interface to a sorting routine. It performs the shell algorithm,...
Definition: sort.F90:149
subroutine, public born_charges_end(this)
subroutine, public born_output_charges(this, atom, charge, natoms, namespace, dim, dirname, write_real)
subroutine, public born_charges_init(this, namespace, natoms, val_charge, qtot, dim)
integer pure function, public magn_dir(dir, ind)
subroutine, public dlr_calc_elf(st, space, gr, kpoints, lr, lr_m)
subroutine, public zlr_calc_elf(st, space, gr, kpoints, lr, lr_m)
character(len=12) function, public freq2str(freq)
integer, parameter perturbation_magnetic
Definition: em_resp.F90:173
subroutine em_resp_run_legacy(sys, fromScratch)
Definition: em_resp.F90:250
subroutine, public out_hyperpolarizability(box, beta, freq_factor, converged, dirname, namespace)
Ref: David M Bishop, Rev Mod Phys 62, 343 (1990) beta generalized to lack of Kleinman symmetry.
Definition: em_resp.F90:3259
subroutine, public em_resp_run(system, from_scratch)
Definition: em_resp.F90:232
subroutine em_resp_output(st, namespace, space, gr, hm, ions, outp, sh, em_vars, iomega, ifactor)
Definition: em_resp.F90:2615
integer, parameter perturbation_none
Definition: em_resp.F90:173
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_huge
Definition: global.F90:205
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_third
Definition: global.F90:194
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
complex(real64), parameter, public m_z0
Definition: global.F90:197
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
character(len= *), parameter, public em_resp_dir
Definition: global.F90:252
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:223
real(real64), parameter, public m_one
Definition: global.F90:188
real(real64), parameter, public m_three
Definition: global.F90:190
real(real64), parameter, public m_five
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
character(len=100) function, public kdotp_wfs_tag(dir, dir2)
Definition: kdotp_calc.F90:152
subroutine, public zlr_orth_response(mesh, st, lr, omega)
subroutine, public lr_copy(st, mesh, src, dest)
subroutine, public lr_allocate(lr, st, mesh, allocate_rho)
subroutine, public lr_init(lr)
subroutine, public lr_dealloc(lr)
subroutine, public dlr_orth_response(mesh, st, lr, omega)
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_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=68), parameter, public hyphens
Definition: messages.F90:161
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
This module implements the basic mulsisystem class, a container system for other systems.
this module contains the output system
Definition: output_low.F90:115
this module contains the output system
Definition: output.F90:115
subroutine, public zoutput_lr(outp, namespace, space, dir, st, mesh, lr, idir, isigma, ions, pert_unit)
Definition: output.F90:2676
subroutine, public doutput_lr(outp, namespace, space, dir, st, mesh, lr, idir, isigma, ions, pert_unit)
Definition: output.F90:2945
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
integer, parameter, public restart_kdotp
Definition: restart.F90:229
integer, parameter, public restart_gs
Definition: restart.F90:229
type(restart_data_t), dimension(restart_n_data_types) info
Definition: restart.F90:254
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:514
subroutine, public restart_open_dir(restart, dirname, ierr)
Change the restart directory to dirname, where "dirname" is a subdirectory of the base restart direct...
Definition: restart.F90:770
integer, parameter, public restart_type_load
Definition: restart.F90:225
subroutine, public restart_end(restart)
Definition: restart.F90:720
subroutine, public restart_close_dir(restart)
Change back to the base directory. To be called after restart_open_dir.
Definition: restart.F90:804
integer, parameter, public smear_fixed_occ
Definition: smear.F90:171
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:829
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
This module handles reading and writing restart information for the states_elec_t.
subroutine, public states_elec_look_and_load(restart, namespace, space, st, mesh, kpoints, is_complex, packed)
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, ierr, iter, lr, lowest_missing, label, verbose, skip)
returns in ierr: <0 => Fatal error, or nothing read =0 => read all wavefunctions >0 => could only rea...
character(len=100) function, public wfs_tag_sigma(namespace, base_name, isigma)
logical function, public sternheimer_add_fxc(this)
subroutine, public sternheimer_unset_kxc(this)
subroutine, public sternheimer_init(this, namespace, space, gr, st, hm, xc, mc, wfs_are_cplx, set_ham_var, set_occ_response, set_last_occ_response, occ_response_by_sternheimer)
subroutine, public sternheimer_end(this)
subroutine, public sternheimer_build_kxc(this, namespace, mesh, st, xc)
logical function, public sternheimer_add_hartree(this)
character(len=80) function, public str_center(s_in, l_in)
puts space around string, so that it is centered
Definition: string.F90:197
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_t), public unit_ppm
Parts per million.
type(unit_system_t), public units_out
type(unit_t), public unit_susc_ppm_cgs
Some magnetic stuff.
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
subroutine, public output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
Definition: utils.F90:242
character pure function, public index2axis(idir)
Definition: utils.F90:202
subroutine, public v_ks_h_setup(namespace, space, gr, ions, ext_partners, st, ks, hm, calc_eigenval, calc_current)
Definition: v_ks.F90:679
Definition: xc.F90:114
class to tell whether a point is inside or outside
Definition: box.F90:141
Class describing the electron system.
Definition: electrons.F90:214
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Container class for lists of system_oct_m::system_t.
output handler class
Definition: output_low.F90:163
The states_elec_t class contains all electronic wave functions.
int true(void)