Octopus
sternheimer.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!! Copyright (C) 2021 Davis Welakuh
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23 use batch_oct_m
25 use comm_oct_m
26 use debug_oct_m
29 use global_oct_m
30 use grid_oct_m
33 use io_oct_m
34 use, intrinsic :: iso_fortran_env
38 use lda_u_oct_m
41 use math_oct_m
42 use mesh_oct_m
45 use mix_oct_m
46 use mpi_oct_m
49 use parser_oct_m
56 use smear_oct_m
57 use space_oct_m
61 use unit_oct_m
63 use types_oct_m
64 use v_ks_oct_m
66 use xc_oct_m
67 use xc_sic_oct_m
68 use xc_f03_lib_m
70
71 implicit none
72
73 private
74 public :: &
95 swap_sigma, &
98 dcalc_hvar, &
99 zcalc_hvar, &
100 dcalc_kvar, &
102
103 character(len=*), public, parameter :: EM_RESP_PHOTONS_DIR = "em_resp_photons/"
104
105 type sternheimer_t
106 private
107 type(linear_solver_t) :: solver
108 type(mix_t) :: mixer
109 type(mixfield_t), pointer :: mixfield
110 type(scf_tol_t) :: scf_tol
111 real(real64), allocatable :: fxc(:,:,:)
112 real(real64), allocatable :: kxc(:,:,:,:)
113 real(real64), pointer, contiguous :: drhs(:, :, :, :) => null()
114 complex(real64), pointer, contiguous :: zrhs(:, :, :, :) => null()
115 real(real64), pointer, contiguous :: dinhomog(:, :, :, :, :) => null()
116 complex(real64), pointer, contiguous :: zinhomog(:, :, :, :, :) => null()
117 logical :: add_fxc_kernel
118 logical :: ok
119 logical :: occ_response
120 logical :: last_occ_response
121 logical :: occ_response_by_sternheimer
122 logical :: preorthogonalization
123 logical, public :: has_photons
124 real(real64) :: domega
125 complex(real64) :: zomega
126 real(real64), allocatable, public :: dphoton_coord_q(:, :)
127 complex(real64), allocatable, public :: zphoton_coord_q(:, :)
128 real(real64) :: pt_eta
129 type(photon_mode_t) :: pt_modes
130
131 real(real64) :: coeff_hartree
132 contains
133 procedure :: add_fxc => sternheimer_add_fxc
134 procedure :: add_hartree => sternheimer_add_hartree
135 end type sternheimer_t
136
137
138contains
139
140 !-----------------------------------------------------------
141 subroutine sternheimer_init(this, namespace, space, gr, st, hm, ks, mc, wfs_are_cplx, set_ham_var, set_occ_response, &
142 set_last_occ_response, occ_response_by_sternheimer)
143 type(sternheimer_t), intent(out) :: this
144 type(namespace_t), intent(in) :: namespace
145 class(space_t), intent(in) :: space
146 type(grid_t), intent(inout) :: gr
147 type(states_elec_t), intent(in) :: st
148 type(hamiltonian_elec_t), intent(in) :: hm
149 type(v_ks_t), intent(in) :: ks
150 type(multicomm_t), intent(in) :: mc
151 logical, intent(in) :: wfs_are_cplx
152 integer, optional, intent(in) :: set_ham_var
153 logical, optional, intent(in) :: set_occ_response
154 logical, optional, intent(in) :: set_last_occ_response
155 logical, optional, intent(in) :: occ_response_by_sternheimer
156
157 logical :: add_hartree
158 integer :: ham_var, iunit
159 logical :: default_preorthog
160
161 push_sub(sternheimer_init)
162
163 if (family_is_mgga(ks%xc%family)) then
164 call messages_not_implemented("MGGA functionals with Sternheimer calculation")
165 end if
166 if (family_is_hybrid(ks%xc)) then
167 call messages_not_implemented("Hybrid functionals with Sternheimer calculation")
168 end if
169
170 if (ks%sic%level /= sic_none) then
171 call messages_not_implemented("SIC with Sternheimer calculation")
172 end if
173
174 if (hm%lda_u_level /= dft_u_none) then
175 call messages_not_implemented("DFT+U with Sterheimer calculation")
176 end if
177
178 if (st%smear%method == smear_fixed_occ) then
179 call messages_experimental("Sternheimer equation for arbitrary occupations", namespace=namespace)
180 end if
181 if (st%smear%method == smear_semiconductor .and. &
182 (abs(st%smear%ef_occ) > m_epsilon) .and. abs(st%smear%ef_occ - m_one) > m_epsilon) then
183 write(message(1),'(a,f12.6)') 'Partial occupation at the Fermi level: ', st%smear%ef_occ
184 message(2) = 'Semiconducting smearing cannot be used for Sternheimer in this situation.'
185 call messages_fatal(2, namespace=namespace)
186 end if
187
188 if (wfs_are_cplx) then
189 call mix_init(this%mixer, namespace, space, gr%der, gr%np, st%d%nspin, func_type_= type_cmplx)
190 else
191 call mix_init(this%mixer, namespace, space, gr%der, gr%np, st%d%nspin, func_type_= type_float)
192 end if
193 call mix_get_field(this%mixer, this%mixfield)
194
195 this%occ_response = optional_default(set_occ_response, .false.)
196 this%occ_response_by_sternheimer = optional_default(occ_response_by_sternheimer, .false.)
197 this%last_occ_response = optional_default(set_last_occ_response, .false.)
199 !%Variable Preorthogonalization
200 !%Type logical
201 !%Section Linear Response::Sternheimer
202 !%Description
203 !% Whether initial linear-response wavefunctions should be orthogonalized
204 !% or not against the occupied states, at the start of each SCF cycle.
205 !% Default is true only if <tt>SmearingFunction = semiconducting</tt>,
206 !% or if the <tt>Occupations</tt> block specifies all full or empty states,
207 !% and we are not solving for linear response in the occupied subspace too.
208 !%End
209 default_preorthog = (st%smear%method == smear_semiconductor .or. &
210 (st%smear%method == smear_fixed_occ .and. st%smear%integral_occs)) &
211 .and. .not. this%occ_response
212 call parse_variable(namespace, 'Preorthogonalization', default_preorthog, this%preorthogonalization)
214 !%Variable HamiltonianVariation
215 !%Type integer
216 !%Default hartree+fxc
217 !%Section Linear Response::Sternheimer
218 !%Description
219 !% The terms to be considered in the variation of the
220 !% Hamiltonian. The external potential (V_ext) is always considered. The default is to include
221 !% also the exchange-correlation and Hartree terms, which fully
222 !% takes into account local fields.
223 !% Just <tt>hartree</tt> gives you the random-phase approximation (RPA).
224 !% If you want to choose the exchange-correlation kernel, use the variable
225 !% <tt>XCKernel</tt>. For <tt>kdotp</tt> and magnetic <tt>em_resp</tt> modes,
226 !% or if <tt>TheoryLevel = independent_particles</tt>,
227 !% the value <tt>V_ext_only</tt> is used and this variable is ignored.
228 !%Option V_ext_only 0
229 !% Neither Hartree nor XC potentials included.
230 !%Option hartree 1
231 !% The variation of the Hartree potential only.
232 !%Option fxc 2
233 !% The exchange-correlation kernel (the variation of the
234 !% exchange-correlation potential) only.
235 !%End
236 if (present(set_ham_var)) then
237 ham_var = set_ham_var
238 else if (hm%theory_level /= independent_particles) then
239 call parse_variable(namespace, 'HamiltonianVariation', 3, ham_var)
240 else
241 ham_var = 0
242 end if
243
244 if (hm%theory_level /= independent_particles) then
245 this%add_fxc_kernel = ((ham_var / 2) == 1)
246 add_hartree = (mod(ham_var, 2) == 1)
247 else
248 this%add_fxc_kernel = .false.
249 add_hartree = .false.
250 end if
251
252
253 message(1) = "Variation of the Hamiltonian in Sternheimer equation: V_ext"
254 if (add_hartree) write(message(1), '(2a)') trim(message(1)), ' + hartree'
255 if (this%add_fxc_kernel) write(message(1), '(2a)') trim(message(1)), ' + fxc'
256
257 message(2) = "Solving Sternheimer equation for"
258 if (this%occ_response) then
259 write(message(2), '(2a)') trim(message(2)), ' full linear response.'
260 else
261 write(message(2), '(2a)') trim(message(2)), ' linear response in unoccupied subspace only.'
262 end if
263
264 message(3) = "Sternheimer preorthogonalization:"
265 if (this%preorthogonalization) then
266 write(message(3), '(2a)') trim(message(3)), ' yes'
267 else
268 write(message(3), '(2a)') trim(message(3)), ' no'
269 end if
270 call messages_info(3, namespace=namespace)
271
272 call linear_solver_init(this%solver, namespace, gr, states_are_real(st), mc, space)
273
274 ! will not converge for non-self-consistent calculation unless LRTolScheme = fixed
275 if (ham_var == 0) then
276 call scf_tol_init(this%scf_tol, namespace, st%qtot, tol_scheme = 0) ! fixed
277 else
278 call scf_tol_init(this%scf_tol, namespace, st%qtot)
279 end if
280
281 ! Build the fxc from the GS density
282 this%coeff_hartree = m_zero
283 if (this%add_fxc_kernel) then
284 this%coeff_hartree = this%coeff_hartree - ks%xc%lrc%alpha / (m_four * m_pi)
285 end if
286 if (add_hartree) this%coeff_hartree = this%coeff_hartree + m_one
287 if (this%add_fxc_kernel) then
288 call sternheimer_build_fxc(this, namespace, gr, st, ks%xc)
289
290 call messages_print_with_emphasis(msg="XC Kernel level", namespace=namespace)
291 call xc_write_fxc_info(ks%xc, namespace=namespace)
292 call xc_sic_write_info(ks%sic, namespace=namespace)
293 call messages_print_with_emphasis(namespace=namespace)
294 end if
295
296
297 ! This variable is documented in xc_oep_init.
298 call parse_variable(namespace, 'EnablePhotons', .false., this%has_photons)
299 call messages_print_var_value('EnablePhotons', this%has_photons, namespace=namespace)
300
301 if (this%has_photons) then
302 call messages_experimental('EnablePhotons = yes', namespace=namespace)
303 call photon_mode_init(this%pt_modes, namespace, space%dim)
304 call photon_mode_set_n_electrons(this%pt_modes, m_zero)
305 call photon_mode_compute_dipoles(this%pt_modes, gr)
306 call io_mkdir(em_resp_photons_dir, namespace)
307 iunit = io_open(em_resp_photons_dir // 'photon_modes', namespace, action='write')
308 call photon_mode_write_info(this%pt_modes, iunit=iunit)
309 safe_allocate(this%zphoton_coord_q(1:this%pt_modes%nmodes, 1:space%dim))
310 end if
311
312 !%Variable PhotonEta
313 !%Type float
314 !%Default 0.0000367
315 !%Section Linear Response::Sternheimer
316 !%Description
317 !% This variable provides the value for the broadening of the photonic spectra
318 !% when the coupling of electrons to photons is enabled in the frequency-dependent Sternheimer equation
319 !%End
320 call parse_variable(namespace, 'PhotonEta', 0.0000367_real64, this%pt_eta, units_inp%energy)
321 call messages_print_var_value('PhotonEta', this%pt_eta, units_inp%energy, namespace=namespace)
322
323
324 pop_sub(sternheimer_init)
325 end subroutine sternheimer_init
326
327 !-----------------------------------------------------------
328 subroutine sternheimer_end(this)
329 type(sternheimer_t), intent(inout) :: this
330
331 push_sub(sternheimer_end)
332
333 safe_deallocate_a(this%zphoton_coord_q)
334
335 call linear_solver_end(this%solver)
336 call scf_tol_end(this%scf_tol)
337 call mix_end(this%mixer)
338
339 nullify(this%mixfield)
340
341 safe_deallocate_a(this%fxc)
342
343 pop_sub(sternheimer_end)
344 end subroutine sternheimer_end
345
346
347 !-----------------------------------------------------------
349 subroutine sternheimer_build_fxc(this, namespace, gr, st, xc)
350 type(sternheimer_t), intent(inout) :: this
351 type(namespace_t), intent(in) :: namespace
352 type(grid_t), intent(in) :: gr
353 type(states_elec_t), intent(in) :: st
354 type(xc_t), intent(in) :: xc
355
356 real(real64), allocatable :: rho(:, :)
357
358 push_sub(sternheimer_build_fxc)
359
360 safe_allocate(this%fxc(1:gr%np, 1:st%d%nspin, 1:st%d%nspin))
361 safe_allocate(rho(1:gr%np_part, 1:st%d%nspin))
362 call states_elec_total_density(st, gr, rho)
363 call xc_get_fxc(xc, gr, namespace, rho, st%d%ispin, this%fxc)
364 safe_deallocate_a(rho)
365
366 pop_sub(sternheimer_build_fxc)
367 end subroutine sternheimer_build_fxc
368
369
370 !-----------------------------------------------------------
371 subroutine sternheimer_build_kxc(this, namespace, mesh, st, xc)
372 type(sternheimer_t), intent(inout) :: this
373 type(namespace_t), intent(in) :: namespace
374 class(mesh_t), intent(in) :: mesh
375 type(states_elec_t), intent(in) :: st
376 type(xc_t), intent(in) :: xc
377
378 real(real64), allocatable :: rho(:, :)
379
380 push_sub(sternheimer_build_kxc)
381
382 if (this%add_fxc()) then
383 safe_allocate(this%kxc(1:mesh%np, 1:st%d%nspin, 1:st%d%nspin, 1:st%d%nspin))
384 this%kxc = m_zero
385
386 safe_allocate(rho(1:mesh%np, 1:st%d%nspin))
387 call states_elec_total_density(st, mesh, rho)
388 call xc_get_kxc(xc, mesh, namespace, rho, st%d%ispin, this%kxc)
389 safe_deallocate_a(rho)
390 end if
391
392 pop_sub(sternheimer_build_kxc)
393
394 end subroutine sternheimer_build_kxc
395
396 !---------------------------------------------
397 subroutine sternheimer_unset_kxc(this)
398 type(sternheimer_t), intent(inout) :: this
399
400 push_sub(sternheimer_unset_kxc)
401
402 safe_deallocate_a(this%kxc)
403
404 pop_sub(sternheimer_unset_kxc)
405 end subroutine sternheimer_unset_kxc
406
407 !-----------------------------------------------------------
408 pure logical function sternheimer_add_fxc(this) result(rr)
409 class(sternheimer_t), intent(in) :: this
410 rr = this%add_fxc_kernel
411 end function sternheimer_add_fxc
412
413
414 !-----------------------------------------------------------
415 pure logical function sternheimer_add_hartree(this) result(rr)
416 class(sternheimer_t), intent(in) :: this
417 rr = abs(this%coeff_hartree) > m_epsilon
418 end function sternheimer_add_hartree
419
420
421 !-----------------------------------------------------------
422 logical function sternheimer_has_converged(this) result(rr)
423 type(sternheimer_t), intent(in) :: this
424 rr = this%ok
425 end function sternheimer_has_converged
426
427 !-----------------------------------------------------------
428 logical pure function sternheimer_have_rhs(this) result(have)
429 type(sternheimer_t), intent(in) :: this
430 have = associated(this%drhs) .or. associated(this%zrhs)
431 end function sternheimer_have_rhs
432
433 !-----------------------------------------------------------
434 subroutine sternheimer_unset_rhs(this)
435 type(sternheimer_t), intent(inout) :: this
436
437 push_sub(sternheimer_unset_rhs)
438
439 nullify(this%drhs)
440 nullify(this%zrhs)
441
442 pop_sub(sternheimer_unset_rhs)
443 end subroutine sternheimer_unset_rhs
445 !-----------------------------------------------------------
446 logical pure function sternheimer_have_inhomog(this) result(have)
447 type(sternheimer_t), intent(in) :: this
448 have = associated(this%dinhomog) .or. associated(this%zinhomog)
449 end function sternheimer_have_inhomog
450
451 !-----------------------------------------------------------
452 subroutine sternheimer_unset_inhomog(this)
453 type(sternheimer_t), intent(inout) :: this
454
456
457 nullify(this%dinhomog)
458 nullify(this%zinhomog)
459
461 end subroutine sternheimer_unset_inhomog
462
463 !-----------------------------------------------------------
464 integer pure function swap_sigma(sigma)
465 integer, intent(in) :: sigma
467 if (sigma == 1) then
468 swap_sigma = 2
469 else
470 swap_sigma = 1
471 end if
472
473 end function swap_sigma
474
475! ---------------------------------------------------------
476 character(len=100) function wfs_tag_sigma(namespace, base_name, isigma) result(str)
477 type(namespace_t), intent(in) :: namespace
478 character(len=*), intent(in) :: base_name
479 integer, intent(in) :: isigma
480
481 character :: sigma_char
482
483 push_sub(wfs_tag_sigma)
484
485 select case (isigma)
486 case (1)
487 sigma_char = '+'
488 case (2)
489 sigma_char = '-'
490 case default
491 write(message(1),'(a,i2)') "Illegal integer isigma passed to wfs_tag_sigma: ", isigma
492 call messages_fatal(1, namespace=namespace)
493 end select
494
495 str = trim(base_name) // sigma_char
496
497 pop_sub(wfs_tag_sigma)
498
499 end function wfs_tag_sigma
500
501 ! --------------------------------------------------------
502
503 subroutine sternheimer_obsolete_variables(namespace, old_prefix, new_prefix)
504 type(namespace_t), intent(in) :: namespace
505 character(len=*), intent(in) :: old_prefix
506 character(len=*), intent(in) :: new_prefix
507
509
510 call messages_obsolete_variable(namespace, trim(old_prefix)//'Preorthogonalization', trim(new_prefix)//'Preorthogonalization')
511 call messages_obsolete_variable(namespace, trim(old_prefix)//'HamiltonianVariation', trim(new_prefix)//'HamiltonianVariation')
512
513 call linear_solver_obsolete_variables(namespace, old_prefix, new_prefix)
514 call scf_tol_obsolete_variables(namespace, old_prefix, new_prefix)
515
518
519 !--------------------------------------------------------------
520 subroutine calc_hvar_photons(this, mesh, nspin, lr_rho, nsigma, hvar, idir)
521 type(sternheimer_t), intent(inout) :: this
522 class(mesh_t), intent(in) :: mesh
523 integer, intent(in) :: nspin
524 integer, intent(in) :: nsigma
525 complex(real64), intent(in) :: lr_rho(:,:)
526 complex(real64), intent(inout) :: hvar(:,:,:)
527 integer, optional, intent(in) :: idir
528
529 real(real64), allocatable :: lambda_dot_r(:)
530 complex(real64), allocatable :: s_lr_rho(:), vp_dip_self_ener(:), vp_bilinear_el_pt(:)
531 integer :: nm, is, ii
532 complex(real64) :: first_moments
533
534 push_sub(calc_hvar_photons)
535 call profiling_in('CALC_HVAR_PHOTONS')
536
537 nm = this%pt_modes%nmodes
538
539 ! photonic terms
540 safe_allocate(s_lr_rho(1:mesh%np))
541 safe_allocate(lambda_dot_r(1:mesh%np))
542 safe_allocate(vp_dip_self_ener(1:mesh%np))
543 safe_allocate(vp_bilinear_el_pt(1:mesh%np))
544
545 ! spin summed density
546 s_lr_rho = m_zero
547 do is = 1, nspin
548 s_lr_rho = s_lr_rho + lr_rho(:, is)
549 end do
550
551 ! Compute electron-photon potentials
552 vp_dip_self_ener = m_zero
553 vp_bilinear_el_pt = m_zero
554 do ii = 1, nm
555 lambda_dot_r(1:mesh%np) = this%pt_modes%lambda(ii)*this%pt_modes%pol_dipole(1:mesh%np, ii)
556 first_moments = zmf_integrate(mesh, lambda_dot_r(1:mesh%np)*s_lr_rho(1:mesh%np))
557
558 ! Compute photon displacement coordinate q_{\alpha}s
559 this%zphoton_coord_q(ii, idir) = (m_one/(m_two*(this%pt_modes%omega(ii))**2)) * &
560 ((m_one/(this%zomega - cmplx(this%pt_modes%omega(ii),-this%pt_eta, real64))) - &
561 (m_one/(this%zomega + cmplx(this%pt_modes%omega(ii), this%pt_eta, real64)))) * &
562 (-(this%pt_modes%omega(ii))**2)*first_moments
563
564 ! Compute potential for bilinear el-pt interaction
565 vp_bilinear_el_pt = vp_bilinear_el_pt - &
566 this%pt_modes%omega(ii)*lambda_dot_r(1:mesh%np)*this%zphoton_coord_q(ii, idir)
567
568 ! Compute potential with dipole-self energy term
569 vp_dip_self_ener = vp_dip_self_ener + first_moments*lambda_dot_r(1:mesh%np)
570 end do
572 hvar(1:mesh%np, 1, 1) = hvar(1:mesh%np, 1, 1) + vp_dip_self_ener(1:mesh%np) + vp_bilinear_el_pt(1:mesh%np)
573
574 safe_deallocate_a(s_lr_rho)
575 safe_deallocate_a(lambda_dot_r)
576 safe_deallocate_a(vp_dip_self_ener)
577 safe_deallocate_a(vp_bilinear_el_pt)
578
579 if (nsigma == 2) hvar(1:mesh%np, 1:nspin, 2) = conjg(hvar(1:mesh%np, 1:nspin, 1))
580
581 call profiling_out('CALC_HVAR_PHOTONS')
582 pop_sub(calc_hvar_photons)
583 end subroutine calc_hvar_photons
584
585
586#include "complex.F90"
587#include "sternheimer_inc.F90"
588
589#include "undef.F90"
590
591#include "real.F90"
592#include "sternheimer_inc.F90"
593
594end module sternheimer_oct_m
595
596!! Local Variables:
597!! mode: f90
598!! coding: utf-8
599!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
subroutine, public states_elec_total_density(st, mesh, total_rho)
This routine calculates the total electronic density.
Definition: density.F90:853
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:237
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_one
Definition: global.F90:192
This module implements the underlying real-space grid.
Definition: grid.F90:119
Definition: io.F90:116
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:360
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:401
A module to handle KS potential, without the external potential.
integer, parameter, public dft_u_none
Definition: lda_u.F90:203
subroutine, public linear_solver_end(this)
subroutine, public linear_solver_init(this, namespace, gr, states_are_real, mc, space)
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1091
character(len=512), private msg
Definition: messages.F90:167
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
subroutine, public mix_get_field(this, mixfield)
Definition: mix.F90:820
subroutine, public mix_init(smix, namespace, space, der, d1, d2, def_, func_type_, prefix_)
Initialise mix_t instance.
Definition: mix.F90:268
subroutine, public mix_end(smix)
Definition: mix.F90:556
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public photon_mode_compute_dipoles(this, mesh)
Computes the polarization dipole.
subroutine, public photon_mode_write_info(this, iunit, namespace)
subroutine, public photon_mode_set_n_electrons(this, qtot)
subroutine, public photon_mode_init(this, namespace, dim, photon_free)
subroutine, public scf_tol_init(this, namespace, qtot, def_maximumiter, tol_scheme)
Definition: scf_tol.F90:162
subroutine, public scf_tol_end(this)
Definition: scf_tol.F90:363
integer, parameter, public smear_semiconductor
Definition: smear.F90:173
integer, parameter, public smear_fixed_occ
Definition: smear.F90:173
pure logical function, public states_are_real(st)
This module handles reading and writing restart information for the states_elec_t.
subroutine, public zsternheimer_set_inhomog(this, inhomog)
subroutine, public zsternheimer_calc_hvar(this, namespace, mesh, hm, lr, nsigma, hvar, idir)
subroutine, public zcalc_hvar(namespace, coeff_hartree, mesh, hm, lr_rho, nsigma, hvar, fxc)
Computes the first-order variation of the Kohn-Sham Hamiltonian.
subroutine, public dsternheimer_solve(this, namespace, space, gr, kpoints, st, hm, mc, lr, nsigma, omega, perturbation, restart, rho_tag, wfs_tag, idir, have_restart_rho, have_exact_freq)
This routine calculates the first-order variations of the wavefunctions for an applied perturbation.
subroutine, public dsternheimer_calc_hvar(this, namespace, mesh, hm, lr, nsigma, hvar, idir)
integer pure function, public swap_sigma(sigma)
subroutine sternheimer_build_fxc(this, namespace, gr, st, xc)
Builds the exchange-correlation kernel for computing the density response.
subroutine, public dcalc_kvar(this, mesh, st, lr_rho1, lr_rho2, nsigma, kvar)
logical function, public sternheimer_has_converged(this)
subroutine, public zsternheimer_set_rhs(this, rhs)
character(len=100) function, public wfs_tag_sigma(namespace, base_name, isigma)
subroutine, public dcalc_hvar(namespace, coeff_hartree, mesh, hm, lr_rho, nsigma, hvar, fxc)
Computes the first-order variation of the Kohn-Sham Hamiltonian.
subroutine calc_hvar_photons(this, mesh, nspin, lr_rho, nsigma, hvar, idir)
logical pure function, public sternheimer_have_rhs(this)
subroutine, public dsternheimer_solve_order2(sh1, sh2, sh_2ndorder, namespace, space, gr, kpoints, st, hm, mc, lr1, lr2, nsigma, omega1, omega2, pert1, pert2, lr_2ndorder, pert_2ndorder, restart, rho_tag, wfs_tag, have_restart_rho, have_exact_freq, give_pert1psi2, give_dl_eig1)
subroutine, public dsternheimer_set_rhs(this, rhs)
pure logical function sternheimer_add_fxc(this)
subroutine, public zsternheimer_solve_order2(sh1, sh2, sh_2ndorder, namespace, space, gr, kpoints, st, hm, mc, lr1, lr2, nsigma, omega1, omega2, pert1, pert2, lr_2ndorder, pert_2ndorder, restart, rho_tag, wfs_tag, have_restart_rho, have_exact_freq, give_pert1psi2, give_dl_eig1)
subroutine, public sternheimer_unset_kxc(this)
subroutine, public zsternheimer_solve(this, namespace, space, gr, kpoints, st, hm, mc, lr, nsigma, omega, perturbation, restart, rho_tag, wfs_tag, idir, have_restart_rho, have_exact_freq)
This routine calculates the first-order variations of the wavefunctions for an applied perturbation.
pure logical function sternheimer_add_hartree(this)
subroutine, public sternheimer_unset_inhomog(this)
logical pure function, public sternheimer_have_inhomog(this)
subroutine, public zcalc_kvar(this, mesh, st, lr_rho1, lr_rho2, nsigma, kvar)
subroutine, public sternheimer_end(this)
subroutine, public sternheimer_obsolete_variables(namespace, old_prefix, new_prefix)
subroutine, public sternheimer_build_kxc(this, namespace, mesh, st, xc)
subroutine, public dsternheimer_set_inhomog(this, inhomog)
subroutine, public sternheimer_init(this, namespace, space, gr, st, hm, ks, mc, wfs_are_cplx, set_ham_var, set_occ_response, set_last_occ_response, occ_response_by_sternheimer)
subroutine, public sternheimer_unset_rhs(this)
type(type_t), public type_float
Definition: types.F90:135
type(type_t), public type_cmplx
Definition: types.F90:136
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
This module provices a simple timer class which can be used to trigger the writing of a restart file ...
Definition: walltimer.F90:123
Definition: xc.F90:116
subroutine, public xc_get_kxc(xcs, mesh, namespace, rho, ispin, kxc)
Definition: xc.F90:2344
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:580
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:608
subroutine, public xc_get_fxc(xcs, gr, namespace, rho, ispin, fxc, fxc_grad, fxc_grad_spin)
Returns the exchange-correlation kernel.
Definition: xc.F90:2002
subroutine, public xc_write_fxc_info(xcs, iunit, namespace)
Definition: xc.F90:255
integer, parameter, public sic_none
no self-interaction correction
Definition: xc_sic.F90:151
subroutine, public xc_sic_write_info(sic, iunit, namespace)
Definition: xc_sic.F90:254
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Describes mesh distribution to nodes.
Definition: mesh.F90:187
The states_elec_t class contains all electronic wave functions.