Octopus
kdotp.F90
Go to the documentation of this file.
1!!! Copyright (C) 2008-2010 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 kdotp_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use grid_oct_m
25 use output_oct_m
27 use io_oct_m
28 use, intrinsic :: iso_fortran_env
33 use mesh_oct_m
36 use mpi_oct_m
39 use parser_oct_m
45 use smear_oct_m
52 use unit_oct_m
54 use utils_oct_m
55 use v_ks_oct_m
56
57 implicit none
58
59 private
60
61 public :: &
64
65 type kdotp_t
66 private
67 class(perturbation_t), pointer :: pert => null()
68 class(perturbation_t), pointer :: pert2 => null()
69
70 real(real64), allocatable :: eff_mass_inv(:,:,:,:)
72 real(real64), allocatable :: velocity(:,:,:)
73
74 type(lr_t), allocatable :: lr(:,:)
77
78 type(lr_t), allocatable :: lr2(:,:,:)
80
81 logical :: converged
82 integer :: occ_solution_method
83 real(real64) :: degen_thres
84 real(real64) :: eta
85 end type kdotp_t
86
87contains
88
89 ! ---------------------------------------------------------
90 subroutine kdotp_lr_run(system, from_scratch)
91 class(*), intent(inout) :: system
92 logical, intent(in) :: from_scratch
93
94 push_sub(kdotp_lr_run)
95
96 select type (system)
97 class is (multisystem_basic_t)
98 message(1) = "CalculationMode = kdotp not implemented for multi-system calculations"
99 call messages_fatal(1, namespace=system%namespace)
100 type is (electrons_t)
101 call kdotp_lr_run_legacy(system, from_scratch)
102 end select
103
104 pop_sub(kdotp_lr_run)
105 end subroutine kdotp_lr_run
106
107 ! ---------------------------------------------------------
108 subroutine kdotp_lr_run_legacy(sys, fromScratch)
109 type(electrons_t), intent(inout) :: sys
110 logical, intent(in) :: fromScratch
111
112 type(kdotp_t) :: kdotp_vars
113 type(sternheimer_t) :: sh, sh2
114 logical :: calc_eff_mass, calc_2nd_order, complex_response
115
116 integer :: idir, idir2, ierr, pdim, ispin
117 character(len=100) :: str_tmp
118 real(real64) :: errornorm
119 type(restart_t) :: restart_load, restart_dump
120
121 class(perturbation_t), pointer :: pert2 ! for the second direction in second-order kdotp
122
123 push_sub(kdotp_lr_run_legacy)
124
125 call messages_experimental("k.p perturbation and calculation of effective masses", namespace=sys%namespace)
126
127 if (sys%hm%pcm%run_pcm) then
128 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
129 end if
130
131 !TODO: This test belongs to the pert.F90 file, in the case of the velocity operator not been defined
132 ! from the Hamiltonian
133 ! In this case, there are other terms missing (MGGA, DFT+U for instance).
134 if (sys%hm%theory_level == hartree_fock) then
135 call messages_not_implemented('Commutator of Fock operator', namespace=sys%namespace)
136 end if
137 if (sys%hm%theory_level == generalized_kohn_sham_dft) then
138 call messages_not_implemented('k.p with generalized Kohn-Sham DFT', namespace=sys%namespace)
139 end if
140
141
142 if (sys%kpoints%use_symmetries) then
143 call messages_experimental("KPoints symmetries with CalculationMode = kdotp", namespace=sys%namespace)
144 end if
145
146 pdim = sys%space%periodic_dim
147
148 if (.not. sys%space%is_periodic()) then
149 message(1) = "k.p perturbation cannot be used for a finite system."
150 call messages_fatal(1, namespace=sys%namespace)
151 end if
152
153 kdotp_vars%pert => perturbation_kdotp_t(sys%namespace, sys%ions)
154 safe_allocate(kdotp_vars%lr(1, 1:pdim))
155
156 call parse_input()
157
158 if (calc_2nd_order) then
159 kdotp_vars%pert2 => perturbation_none_t(sys%namespace)
160 pert2 => perturbation_kdotp_t(sys%namespace, sys%ions)
161 call kdotp_vars%pert2%setup_dir(1) ! direction is irrelevant
162 safe_allocate(kdotp_vars%lr2(1, 1:pdim, 1:pdim))
163 end if
164
165 !Read ground-state wavefunctions
166 complex_response = (abs(kdotp_vars%eta) > m_epsilon) .or. states_are_complex(sys%st)
167 call restart_init(restart_load, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
168 if (ierr == 0) then
169 call states_elec_look_and_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
170 is_complex = complex_response)
171 call restart_end(restart_load)
172 else
173 message(1) = "A previous gs calculation is required."
174 call messages_fatal(1, namespace=sys%namespace)
175 end if
177 ! Use of ForceComplex will make this true after states_elec_look_and_load even if it was not before.
178 ! Otherwise, this line is a tautology.
179 complex_response = states_are_complex(sys%st)
180
181 ! Start restart. Note: we are going to use the same directory to read and write.
182 ! Therefore, restart_dump must be initialized first to make sure the directory
183 ! exists when we initialize restart_load.
184 call restart_init(restart_dump, sys%namespace, restart_kdotp, restart_type_dump, sys%mc, ierr, mesh=sys%gr)
185 ! no problem if this fails
186 call restart_init(restart_load, sys%namespace, restart_kdotp, restart_type_load, sys%mc, ierr, mesh=sys%gr)
187
188 ! setup Hamiltonian
189 message(1) = 'Info: Setting up Hamiltonian for linear response.'
190 call messages_info(1, namespace=sys%namespace)
191 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
192
193 if (states_are_real(sys%st)) then
194 message(1) = 'Info: Using real wavefunctions.'
195 else
196 message(1) = 'Info: Using complex wavefunctions.'
197 end if
198 call messages_info(1, namespace=sys%namespace)
199
200 ! Band velocities
201 message(1) = 'Calculating band velocities.'
202 call messages_info(1, namespace=sys%namespace)
203
204 ! Compute band velocities
205 safe_allocate(kdotp_vars%velocity(1:pdim, 1:sys%st%nst, 1:sys%st%nik))
206 call calc_band_velocity(sys%namespace, sys%space, sys%gr, sys%st, sys%hm, kdotp_vars%pert, &
207 kdotp_vars%velocity(:,:,:))
208
209 ! Output of the band velocities
210 if (mpi_grp_is_root(mpi_world)) then
211 call io_mkdir(kdotp_dir, sys%namespace) ! data output
212 call kdotp_write_band_velocity(sys%st, pdim, kdotp_vars%velocity(:,:,:), sys%namespace)
213 end if
214 safe_deallocate_a(kdotp_vars%velocity)
215
216
217 call sternheimer_obsolete_variables(sys%namespace, 'KdotP_', 'KdotP')
218 call sternheimer_init(sh, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, complex_response, &
219 set_ham_var = 0, set_occ_response = (kdotp_vars%occ_solution_method == 0), &
220 set_last_occ_response = (kdotp_vars%occ_solution_method == 0), occ_response_by_sternheimer = .true.)
221 ! ham_var_set = 0 results in HamiltonianVariation = V_ext_only
222 if (calc_2nd_order) then
223 call sternheimer_init(sh2, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, complex_response, &
224 set_ham_var = 0, set_occ_response = .false., set_last_occ_response = .false.)
225 end if
226
227 do idir = 1, pdim
228 call lr_init(kdotp_vars%lr(1, idir))
229 call lr_allocate(kdotp_vars%lr(1, idir), sys%st, sys%gr)
230
231 if (calc_2nd_order) then
232 do idir2 = idir, pdim
233 call lr_init(kdotp_vars%lr2(1, idir, idir2))
234 call lr_allocate(kdotp_vars%lr2(1, idir, idir2), sys%st, sys%gr)
235 end do
236 end if
237
238 ! load wavefunctions
239 if (.not. fromscratch) then
240 str_tmp = kdotp_wfs_tag(idir)
241 call restart_open_dir(restart_load, wfs_tag_sigma(sys%namespace, str_tmp, 1), ierr)
242 if (ierr == 0) call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
243 ierr, lr=kdotp_vars%lr(1, idir))
244 call restart_close_dir(restart_load)
245
246 if (ierr /= 0) then
247 message(1) = "Unable to read response wavefunctions from '"//trim(wfs_tag_sigma(sys%namespace, str_tmp, 1))//"'."
248 call messages_warning(1, namespace=sys%namespace)
249 end if
250
251 if (calc_2nd_order) then
252 do idir2 = idir, pdim
253 str_tmp = kdotp_wfs_tag(idir, idir2)
254 call restart_open_dir(restart_load, wfs_tag_sigma(sys%namespace, str_tmp, 1), ierr)
255 if (ierr == 0) then
256 call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
257 ierr, lr=kdotp_vars%lr2(1, idir, idir2))
258 end if
259 call restart_close_dir(restart_load)
260
261 if (ierr /= 0) then
262 message(1) = "Unable to read response wavefunctions from '"//trim(wfs_tag_sigma(sys%namespace, str_tmp, 1))//"'."
263 call messages_warning(1, namespace=sys%namespace)
264 end if
265 end do
266 end if
267 end if
268 end do
269
270 call info()
271 message(1) = "Info: Calculating k.p linear response of ground-state wavefunctions."
272 call messages_info(1, namespace=sys%namespace)
273 kdotp_vars%converged = .true.
274
275 ! solve the Sternheimer equation
276 do idir = 1, pdim
277 write(message(1), '(3a)') 'Info: Calculating response for the ', index2axis(idir), &
278 '-direction.'
279 call messages_info(1, namespace=sys%namespace)
280 call kdotp_vars%pert%setup_dir(idir)
281
282 if (states_are_real(sys%st)) then
283 call dsternheimer_solve(sh, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%mc, &
284 kdotp_vars%lr(1:1, idir), 1, m_zero, kdotp_vars%pert, restart_dump, "", kdotp_wfs_tag(idir), &
285 have_restart_rho = .false.)
286 if (kdotp_vars%occ_solution_method == 1) then
287 call dkdotp_add_occ(sys%namespace, sys%space, sys%gr, sys%st, sys%hm, kdotp_vars%pert, &
288 kdotp_vars%lr(1, idir), kdotp_vars%degen_thres)
289 end if
290 else
291 call zsternheimer_solve(sh, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%mc, &
292 kdotp_vars%lr(1:1, idir), 1, m_zi * kdotp_vars%eta, kdotp_vars%pert, restart_dump, "", &
293 kdotp_wfs_tag(idir), have_restart_rho = .false.)
294 if (kdotp_vars%occ_solution_method == 1) then
295 call zkdotp_add_occ(sys%namespace, sys%space, sys%gr, sys%st, sys%hm, kdotp_vars%pert, &
296 kdotp_vars%lr(1, idir), kdotp_vars%degen_thres)
297 end if
298 end if
299
300 kdotp_vars%converged = kdotp_vars%converged .and. sternheimer_has_converged(sh)
301
302 errornorm = m_zero
303 if (states_are_real(sys%st)) then
304 call doutput_lr(sys%outp, sys%namespace, sys%space, kdotp_dir, sys%st, sys%gr, kdotp_vars%lr(1, idir), idir, 1, &
305 sys%ions, units_out%force)
306
307 do ispin = 1, sys%st%d%nspin
308 errornorm = hypot(errornorm, dmf_nrm2(sys%gr, kdotp_vars%lr(1, idir)%ddl_rho(:, ispin)))
309 end do
310 else
311 call zoutput_lr(sys%outp, sys%namespace, sys%space, kdotp_dir, sys%st, sys%gr, kdotp_vars%lr(1, idir), idir, 1, &
312 sys%ions, units_out%force)
313
314 do ispin = 1, sys%st%d%nspin
315 errornorm = hypot(errornorm, zmf_nrm2(sys%gr, kdotp_vars%lr(1, idir)%zdl_rho(:, ispin)))
316 end do
317 end if
318
319 write(message(1),'(a,g13.6)') "Norm of relative density variation = ", errornorm / sys%st%qtot
320 call messages_info(1, namespace=sys%namespace)
321
322 if (calc_2nd_order) then
323 ! by equality of mixed partial derivatives, kdotp_vars%lr2(idir, idir2) = kdotp_vars%lr2(idir2, idir)
324 do idir2 = idir, pdim
325 write(message(1), '(3a)') 'Info: Calculating second-order response in the ', index2axis(idir2), &
326 '-direction.'
327 call messages_info(1, namespace=sys%namespace)
328
329 call pert2%setup_dir(idir2)
330
331 if (states_are_real(sys%st)) then
332 call dsternheimer_solve_order2(sh, sh, sh2, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, &
333 sys%mc, kdotp_vars%lr(1:1, idir), kdotp_vars%lr(1:1, idir2), &
334 1, m_zero, m_zero, kdotp_vars%pert, pert2, &
335 kdotp_vars%lr2(1:1, idir, idir2), kdotp_vars%pert2, restart_dump, "", kdotp_wfs_tag(idir, idir2), &
336 have_restart_rho = .false., have_exact_freq = .true.)
337 else
338 call zsternheimer_solve_order2(sh, sh, sh2, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, &
339 sys%mc, kdotp_vars%lr(1:1, idir), kdotp_vars%lr(1:1, idir2), &
340 1, m_zi * kdotp_vars%eta, m_zi * kdotp_vars%eta, kdotp_vars%pert, pert2, &
341 kdotp_vars%lr2(1:1, idir, idir2), kdotp_vars%pert2, restart_dump, "", kdotp_wfs_tag(idir, idir2), &
342 have_restart_rho = .false., have_exact_freq = .true.)
343 end if
344
345 end do
346 message(1) = ""
347 call messages_info(1, namespace=sys%namespace)
348 end if
349 end do ! idir
350
351 ! calculate effective masses
352 if (calc_eff_mass) then
353 message(1) = "Info: Calculating effective masses."
354 call messages_info(1, namespace=sys%namespace)
355
356 safe_allocate(kdotp_vars%eff_mass_inv(1:pdim, 1:pdim, 1:sys%st%nst, 1:sys%st%nik))
357 kdotp_vars%eff_mass_inv(:,:,:,:) = m_zero
358
359
360 if (states_are_real(sys%st)) then
361 call dcalc_eff_mass_inv(sys%namespace, sys%space, sys%gr, sys%st, sys%hm, kdotp_vars%lr, &
362 kdotp_vars%pert, kdotp_vars%eff_mass_inv, kdotp_vars%degen_thres)
363 else
364 call zcalc_eff_mass_inv(sys%namespace, sys%space, sys%gr, sys%st, sys%hm, kdotp_vars%lr, &
365 kdotp_vars%pert, kdotp_vars%eff_mass_inv, kdotp_vars%degen_thres)
366 end if
367
368 call kdotp_write_degeneracies(sys%st, kdotp_vars%degen_thres, sys%namespace)
369 call kdotp_write_eff_mass(sys%st, sys%kpoints, kdotp_vars, sys%namespace, sys%space%periodic_dim)
370
371 safe_deallocate_a(kdotp_vars%eff_mass_inv)
372 end if
373
374 ! clean up some things
375 do idir = 1, pdim
376 call lr_dealloc(kdotp_vars%lr(1, idir))
377
378 if (calc_2nd_order) then
379 do idir2 = idir, pdim
380 call lr_dealloc(kdotp_vars%lr2(1, idir, idir2))
381 end do
382 end if
383 end do
384
385 call restart_end(restart_load)
386 call restart_end(restart_dump)
387 call sternheimer_end(sh)
388 deallocate(kdotp_vars%pert)
389 safe_deallocate_a(kdotp_vars%lr)
390
391 if (calc_2nd_order) then
392 call sternheimer_end(sh2)
393 safe_deallocate_p(pert2)
394 deallocate(kdotp_vars%pert2)
395 safe_deallocate_a(kdotp_vars%lr2)
396 end if
397
398 call states_elec_deallocate_wfns(sys%st)
399
400 pop_sub(kdotp_lr_run_legacy)
401
402 contains
403
404 ! ---------------------------------------------------------
405
406 subroutine parse_input()
407
409
410 !%Variable KdotPOccupiedSolutionMethod
411 !%Type integer
412 !%Default sternheimer_eqn
413 !%Section Linear Response::KdotP
414 !%Description
415 !% Method of calculating the contribution of the projection of the
416 !% linear-response wavefunctions in the occupied subspace.
417 !%Option sternheimer_eqn 0
418 !% The Sternheimer equation is solved including the occupied subspace,
419 !% to get the full linear-response wavefunctions.
420 !%Option sum_over_states 1
421 !% The Sternheimer equation is solved only in the unoccupied subspace,
422 !% and a sum-over-states perturbation-theory expression is used to
423 !% evaluate the contributions in the occupied subspace.
424 !%End
425
426 call messages_obsolete_variable(sys%namespace, 'KdotP_OccupiedSolutionMethod', 'KdotPOccupiedSolutionMethod')
427
428 call parse_variable(sys%namespace, 'KdotPOccupiedSolutionMethod', 0, kdotp_vars%occ_solution_method)
429 if (kdotp_vars%occ_solution_method == 1 .and. .not. smear_is_semiconducting(sys%st%smear)) then
430 call messages_not_implemented('KdotPOccupiedSolutionMethod = sum_over_states for non-semiconducting smearing', &
431 namespace=sys%namespace)
432 end if
433
434 !%Variable DegeneracyThreshold
435 !%Type float
436 !%Default 1e-5
437 !%Section States
438 !%Description
439 !% States with energy <math>E_i</math> and <math>E_j</math> will be considered degenerate
440 !% if <math> \left| E_i - E_j \right| < </math><tt>DegeneracyThreshold</tt>.
441 !%End
442 call parse_variable(sys%namespace, 'DegeneracyThreshold', units_from_atomic(units_inp%energy, 1e-5_real64), &
443 kdotp_vars%degen_thres)
444 kdotp_vars%degen_thres = units_to_atomic(units_inp%energy, kdotp_vars%degen_thres)
445
446 !%Variable KdotPEta
447 !%Type float
448 !%Default 0.0
449 !%Section Linear Response::KdotP
450 !%Description
451 !% Imaginary frequency added to Sternheimer equation which may improve convergence.
452 !% Not recommended.
453 !%End
454 call messages_obsolete_variable(sys%namespace, 'KdotP_Eta', 'KdotPEta')
455 call parse_variable(sys%namespace, 'KdotPEta', m_zero, kdotp_vars%eta)
456 kdotp_vars%eta = units_to_atomic(units_inp%energy, kdotp_vars%eta)
457
458 !%Variable KdotPCalculateEffectiveMasses
459 !%Type logical
460 !%Default true
461 !%Section Linear Response::KdotP
462 !%Description
463 !% If true, uses <tt>kdotp</tt> perturbations of ground-state wavefunctions
464 !% to calculate effective masses. It is not correct for degenerate states.
465 !%End
466 call messages_obsolete_variable(sys%namespace, 'KdotP_CalculateEffectiveMasses', 'KdotPCalculateEffectiveMasses')
467 call parse_variable(sys%namespace, 'KdotPCalculateEffectiveMasses', .true., calc_eff_mass)
468
469 !%Variable KdotPCalcSecondOrder
470 !%Type logical
471 !%Default false
472 !%Section Linear Response::KdotP
473 !%Description
474 !% If true, calculates second-order response of wavefunctions as well as first-order response.
475 !% Note that the second derivative of the Hamiltonian is NOT included in this calculation.
476 !% This is needed for a subsequent run in <tt>CalculationMode = em_resp</tt> with <tt>EMHyperpol</tt>.
477 !%End
478 call parse_variable(sys%namespace, 'KdotPCalcSecondOrder', .false., calc_2nd_order)
479
481
482 end subroutine parse_input
483
484 ! ---------------------------------------------------------
485 subroutine info()
486
487 push_sub(kdotp_lr_run_legacy.info)
488
489 call kdotp_vars%pert%info()
490
491 call messages_print_with_emphasis(msg='k.p perturbation theory', namespace=sys%namespace)
492
493 if (kdotp_vars%occ_solution_method == 0) then
494 message(1) = 'Occupied solution method: Sternheimer equation.'
495 else
496 message(1) = 'Occupied solution method: sum over states.'
497 end if
498
499 call messages_info(1, namespace=sys%namespace)
500
501 call messages_print_with_emphasis(namespace=sys%namespace)
502
504
505 end subroutine info
506
507 end subroutine kdotp_lr_run_legacy
508
509 ! ---------------------------------------------------------
510 subroutine kdotp_write_band_velocity(st, periodic_dim, velocity, namespace)
511 type(states_elec_t), intent(inout) :: st
512 integer, intent(in) :: periodic_dim
513 real(real64), intent(in) :: velocity(:,:,:)
514 type(namespace_t), intent(inout) :: namespace
515
516 character(len=80) :: filename, tmp
517 integer :: iunit, ik, ist, ik2, ispin, idir
518
519 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
520
522
523 write(filename, '(a)') kdotp_dir//'velocity'
524 iunit = io_open(trim(filename), namespace, action='write')
525 write(iunit,'(a)') '# Band velocities'
526
527 do ik = 1, st%nik
528 ispin = st%d%get_spin_index(ik)
529 ik2 = st%d%get_kpoint_index(ik)
530 tmp = int2str(ik2)
531
532 write(iunit,'(a,i1,a,a)') '# spin = ', ispin, ', k-point = ', trim(tmp)
533
534 write(iunit,'(a)',advance='no') '# state energy '
535 do idir = 1, periodic_dim
536 write(iunit,'(3a)',advance='no') 'vg(', trim(index2axis(idir)), ') '
537 end do
538 write(iunit,'(a)')
539
540 write(iunit,'(3a)',advance='no') '# [', trim(units_abbrev(units_out%energy)), '] '
541 do idir = 1, periodic_dim
542 write(iunit,'(3a)',advance='no') '[', trim(units_abbrev(units_out%velocity)), '] '
543 end do
544 write(iunit,'(a)')
545
546 do ist = 1, st%nst
547 write(iunit,'(i5,f12.5,3f12.5)') ist, units_from_atomic(units_out%energy, st%eigenval(ist, ik)), &
548 velocity(1:periodic_dim, ist, ik)
549 end do
550 end do
551
552 call io_close(iunit)
554 end subroutine kdotp_write_band_velocity
555
556 ! ---------------------------------------------------------
557 subroutine kdotp_write_eff_mass(st, kpoints, kdotp_vars, namespace, periodic_dim)
558 type(states_elec_t), intent(inout) :: st
559 type(kpoints_t), intent(in) :: kpoints
560 type(kdotp_t), intent(inout) :: kdotp_vars
561 type(namespace_t), intent(in) :: namespace
562 integer, intent(in) :: periodic_dim
563
564 character(len=80) :: filename, tmp
565 integer :: iunit, ik, ist, ik2, ispin
566
567 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
568
569 push_sub(kdotp_write_eff_mass)
570
571 do ik = 1, st%nik
572 ispin = st%d%get_spin_index(ik)
573 ik2 = st%d%get_kpoint_index(ik)
574
575 tmp = int2str(ik2)
576 write(filename, '(3a, i1)') kdotp_dir//'kpoint_', trim(tmp), '_', ispin
577 iunit = io_open(trim(filename), namespace, action='write')
579 write(iunit,'(a, i10)') '# spin index = ', ispin
580 write(iunit,'(a, i10)') '# k-point index = ', ik2
581 write(iunit,'(a, 99f12.8)') '# k-point coordinates = ', kpoints%get_point(ik2)
582 if (.not. kdotp_vars%converged) write(iunit, '(a)') "# WARNING: not converged"
583
584 write(iunit,'(a)')
585 write(iunit,'(a)') '# Inverse effective-mass tensors'
586 do ist = 1, st%nst
587 write(iunit,'(a)')
588 tmp = int2str(ist)
589 write(iunit,'(a, a, a, f12.8, a, a)') 'State #', trim(tmp), ', Energy = ', &
590 units_from_atomic(units_out%energy, st%eigenval(ist, ik)), ' ', units_abbrev(units_out%energy)
591 call output_tensor(kdotp_vars%eff_mass_inv(:, :, ist, ik), periodic_dim, unit_one, iunit=iunit)
592 end do
593
594 write(iunit,'(a)')
595 write(iunit,'(a)') '# Effective-mass tensors'
596 do ist = 1, st%nst
597 write(iunit,'(a)')
598 tmp = int2str(ist)
599 write(iunit,'(a, a, a, f12.8, a, a)') 'State #', trim(tmp), ', Energy = ', &
600 units_from_atomic(units_out%energy, st%eigenval(ist, ik)), ' ', units_abbrev(units_out%energy)
601 call lalg_inverse(periodic_dim, kdotp_vars%eff_mass_inv(:, :, ist, ik), 'dir')
602 call output_tensor(kdotp_vars%eff_mass_inv(:, :, ist, ik), periodic_dim, unit_one, iunit=iunit)
603 end do
604
605 call io_close(iunit)
606 end do
607
608 pop_sub(kdotp_write_eff_mass)
609 end subroutine kdotp_write_eff_mass
610
611 ! ---------------------------------------------------------
612 subroutine kdotp_write_degeneracies(st, threshold, namespace)
613 type(states_elec_t), intent(inout) :: st
614 real(real64), intent(in) :: threshold
615 type(namespace_t), intent(in) :: namespace
616
617 character(len=80) :: tmp
618 integer :: ik, ist, ist2, ik2, ispin
619
620 if (.not. mpi_grp_is_root(mpi_world)) return ! only first node outputs
621
623
624 call messages_print_with_emphasis(msg='Degenerate subspaces', namespace=namespace)
625
626 do ik = 1, st%nik
627 ispin = st%d%get_spin_index(ik)
628 ik2 = st%d%get_kpoint_index(ik)
629
630 tmp = int2str(ik2)
631 write(message(1), '(3a, i1)') 'k-point ', trim(tmp), ', spin ', ispin
632 call messages_info(1, namespace=namespace)
633
634 ist = 1
635 do while (ist <= st%nst)
636 ! test for degeneracies
637 write(message(1),'(a)') '===='
638 tmp = int2str(ist)
639 write(message(2),'(a, a, a, f12.8, a, a)') 'State #', trim(tmp), ', Energy = ', &
640 units_from_atomic(units_out%energy, st%eigenval(ist, ik)), ' ', units_abbrev(units_out%energy)
641 call messages_info(2, namespace=namespace)
642
643 ist2 = ist + 1
644 do while (ist2 <= st%nst .and. &
645 abs(st%eigenval(min(ist2, st%nst), ik) - st%eigenval(ist, ik)) < threshold)
646 tmp = int2str(ist2)
647 write(message(1),'(a, a, a, f12.8, a, a)') 'State #', trim(tmp), ', Energy = ', &
648 units_from_atomic(units_out%energy, st%eigenval(ist2, ik)), ' ', units_abbrev(units_out%energy)
649 call messages_info(1, namespace=namespace)
650 ist2 = ist2 + 1
651 end do
652
653 ist = ist2
654 end do
655
656 write(message(1),'()')
657 call messages_info(1, namespace=namespace)
658 end do
659
660 message(1) = "Velocities and effective masses are not correct within degenerate subspaces."
661 call messages_warning(1, namespace=namespace)
662
664 end subroutine kdotp_write_degeneracies
665
666 ! ---------------------------------------------------------
667 character(len=12) pure function int2str(ii) result(str)
668 integer, intent(in) :: ii
669
670 write(str, '(i11)') ii
671 str = trim(adjustl(str))
672
673 end function int2str
674
675end module kdotp_oct_m
676
677!! Local Variables:
678!! mode: f90
679!! coding: utf-8
680!! End:
subroutine parse_input()
Definition: em_resp.F90:819
double hypot(double __x, double __y) __attribute__((__nothrow__
real(real64), parameter, public m_zero
Definition: global.F90:187
character(len= *), parameter, public kdotp_dir
Definition: global.F90:254
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
This module implements the underlying real-space grid.
Definition: grid.F90:117
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public hartree_fock
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
subroutine, public calc_band_velocity(namespace, space, gr, st, hm, pert, velocity)
Computes the k-point and band-resolved velocity.
Definition: kdotp_calc.F90:174
subroutine, public dkdotp_add_occ(namespace, space, gr, st, hm, pert, kdotp_lr, degen_thres)
add projection onto occupied states, by sum over states
Definition: kdotp_calc.F90:400
subroutine, public zkdotp_add_occ(namespace, space, gr, st, hm, pert, kdotp_lr, degen_thres)
add projection onto occupied states, by sum over states
Definition: kdotp_calc.F90:639
subroutine, public dcalc_eff_mass_inv(namespace, space, gr, st, hm, lr, pert, eff_mass_inv, degen_thres)
Computes the effective mass tensor.
Definition: kdotp_calc.F90:302
character(len=100) function, public kdotp_wfs_tag(dir, dir2)
Definition: kdotp_calc.F90:152
subroutine, public zcalc_eff_mass_inv(namespace, space, gr, st, hm, lr, pert, eff_mass_inv, degen_thres)
Computes the effective mass tensor.
Definition: kdotp_calc.F90:541
subroutine kdotp_write_band_velocity(st, periodic_dim, velocity, namespace)
Definition: kdotp.F90:604
subroutine kdotp_lr_run_legacy(sys, fromScratch)
Definition: kdotp.F90:202
subroutine kdotp_write_eff_mass(st, kpoints, kdotp_vars, namespace, periodic_dim)
Definition: kdotp.F90:651
subroutine kdotp_write_degeneracies(st, threshold, namespace)
Definition: kdotp.F90:706
subroutine, public kdotp_lr_run(system, from_scratch)
Definition: kdotp.F90:184
character(len=12) pure function, public int2str(ii)
Definition: kdotp.F90:761
subroutine, public lr_allocate(lr, st, mesh, allocate_rho)
subroutine, public lr_init(lr)
subroutine, public lr_dealloc(lr)
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=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 implements the basic mulsisystem class, a container system for other systems.
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, 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
integer, parameter, public restart_type_dump
Definition: restart.F90:225
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
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:829
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...
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.
logical function, public sternheimer_has_converged(this)
character(len=100) function, public wfs_tag_sigma(namespace, base_name, isigma)
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 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_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 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.
subroutine, public sternheimer_end(this)
subroutine, public sternheimer_obsolete_variables(namespace, old_prefix, new_prefix)
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_system_t), public units_out
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
Class describing the electron system.
Definition: electrons.F90:214
Container class for lists of system_oct_m::system_t.
The states_elec_t class contains all electronic wave functions.
int true(void)