Octopus
phonons_lr.F90
Go to the documentation of this file.
1!! Copyright (C) 2007-2012 Xavier Andrade, 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
23 use debug_oct_m
24 use epot_oct_m
25 use forces_oct_m
26 use global_oct_m
27 use grid_oct_m
29 use io_oct_m
31 use ions_oct_m
32 use, intrinsic :: iso_fortran_env
33 use kdotp_oct_m
38 use math_oct_m
39 use mesh_oct_m
42 use mpi_oct_m
45 use parser_oct_m
49 use smear_oct_m
50 use space_oct_m
57 use unit_oct_m
59 use utils_oct_m
61 use v_ks_oct_m
62
63 implicit none
64
65 private
66 public :: &
72
73contains
74
75 ! ---------------------------------------------------------
76 subroutine phonons_lr_run(system, from_scratch)
77 class(*), intent(inout) :: system
78 logical, intent(in) :: from_scratch
79
80 push_sub(phonons_lr_run)
81
82 select type (system)
83 class is (multisystem_basic_t)
84 message(1) = "CalculationMode = vib_modes not implemented for multi-system calculations"
85 call messages_fatal(1, namespace=system%namespace)
86 type is (electrons_t)
87 call phonons_lr_run_legacy(system, from_scratch)
88 end select
89
90 pop_sub(phonons_lr_run)
91 end subroutine phonons_lr_run
92
93 ! ---------------------------------------------------------
94 subroutine phonons_lr_run_legacy(sys, fromscratch)
95 type(electrons_t), target, intent(inout) :: sys
96 logical, intent(in) :: fromscratch
97
98 type(sternheimer_t) :: sh
99 type(lr_t) :: lr(1:1), kdotp_lr(sys%space%dim)
100 type(vibrations_t) :: vib
101 class(perturbation_ionic_t), pointer :: pert
102
103 type(ions_t), pointer :: ions
104 type(states_elec_t), pointer :: st
105 type(grid_t), pointer :: gr
106
107 integer :: natoms, ndim, iatom, idir, jatom, jdir, imat, jmat, iunit_restart, ierr, start_mode
108 complex(real64), allocatable :: force_deriv(:,:)
109 character(len=80) :: str_tmp
110 character(len=300) :: line(1)
111 type(born_charges_t) :: born
112 logical :: normal_mode_wfs, do_infrared, symmetrize
113 type(restart_t) :: restart_load, restart_dump, kdotp_restart, gs_restart
115 push_sub(phonons_lr_run_legacy)
116
117 !some shortcuts
118
119 ions => sys%ions
120 st => sys%st
121 gr => sys%gr
122
123 if (sys%hm%pcm%run_pcm) then
124 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
125 end if
126
127 if (sys%space%is_periodic()) then
128 call messages_not_implemented('linear-response vib_modes for periodic systems', namespace=sys%namespace)
129 end if
130
131 if (sys%hm%ep%nlcc) then
132 call messages_not_implemented('linear-response vib_modes with non-linear core corrections', namespace=sys%namespace)
133 end if
134
135 !%Variable CalcNormalModeWfs
136 !%Type logical
137 !%Default false
138 !%Section Linear Response::Vibrational Modes
139 !%Description
140 !% If set to true, the response wavefunctions for each normal mode will be calculated
141 !% and written in directory <tt>restart/vib_modes/phn_nm_wfs_XXXXX</tt>.
142 !% This part is time-consuming and not parallel, but not needed for most purposes.
143 !%End
144 call parse_variable(sys%namespace, 'CalcNormalModeWfs', .false., normal_mode_wfs)
145
146 !%Variable CalcInfrared
147 !%Type logical
148 !%Default true
149 !%Section Linear Response::Vibrational Modes
150 !%Description
151 !% If set to true, infrared intensities (and born charges) will be calculated
152 !% and written in <tt>vib_modes/infrared</tt>.
153 !%End
154 call parse_variable(sys%namespace, 'CalcInfrared', .true., do_infrared)
155
156 !%Variable SymmetrizeDynamicalMatrix
157 !%Type logical
158 !%Default true
159 !%Section Linear Response::Vibrational Modes
160 !%Description
161 !% If set to true, all entries of the dynamical matrix will be calculated and then
162 !% the matrix will be symmetrized to enforce <math>D_{ij} = D_{ji}</math>. If set to false,
163 !% only the upper half of the matrix will be calculated.
164 !%End
165 call parse_variable(sys%namespace, 'SymmetrizeDynamicalMatrix', .true., symmetrize)
166
167 ! replaced by properly saving and reading the dynamical matrix
168 call messages_obsolete_variable(sys%namespace, 'UseRestartDontSolve')
170 natoms = ions%natoms
171 ndim = sys%space%dim
172
173 call restart_init(gs_restart, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=gr, exact=.true.)
174 if (ierr == 0) then
175 call states_elec_look_and_load(gs_restart, sys%namespace, sys%space, st, sys%gr, sys%kpoints, st%restart_fixed_occ)
176 call restart_end(gs_restart)
177 else
178 message(1) = "Previous gs calculation is required."
179 call messages_fatal(1, namespace=sys%namespace)
180 end if
181
182 ! read kdotp wavefunctions if necessary (for IR intensities)
183 if (sys%space%is_periodic() .and. do_infrared) then
184 message(1) = "Reading kdotp wavefunctions for periodic directions."
185 call messages_info(1, namespace=sys%namespace)
186
187 call restart_init(kdotp_restart, sys%namespace, restart_kdotp, restart_type_load, sys%mc, ierr, mesh=gr)
188 if (ierr /= 0) then
189 message(1) = "Unable to read kdotp wavefunctions."
190 message(2) = "Previous kdotp calculation required."
191 call messages_fatal(2, namespace=sys%namespace)
192 end if
193
194 do idir = 1, sys%space%periodic_dim
195 call lr_init(kdotp_lr(idir))
196 call lr_allocate(kdotp_lr(idir), sys%st, sys%gr)
197
198 ! load wavefunctions
199 str_tmp = trim(kdotp_wfs_tag(idir))
200 call restart_open_dir(kdotp_restart, wfs_tag_sigma(sys%namespace, str_tmp, 1), ierr)
201 if (ierr == 0) then
202 call states_elec_load(kdotp_restart, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, &
203 sys%st%restart_fixed_occ, ierr=ierr, lr=kdotp_lr(idir))
204 end if
205 call restart_close_dir(kdotp_restart)
206
207 if (ierr /= 0) then
208 message(1) = "Unable to read kdotp wavefunctions from '"//trim(wfs_tag_sigma(sys%namespace, str_tmp, 1))//"'."
209 message(2) = "Previous kdotp calculation required."
210 call messages_fatal(2, namespace=sys%namespace)
211 end if
212 end do
213 call restart_end(kdotp_restart)
214 end if
215
216 message(1) = 'Info: Setting up Hamiltonian for linear response.'
217 call messages_info(1, namespace=sys%namespace)
218
219 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
220 call sternheimer_init(sh, sys%namespace, sys%space, sys%gr, sys%st, sys%hm, sys%ks%xc, sys%mc, &
221 wfs_are_cplx = states_are_complex(st))
222
223 call vibrations_init(vib, ions%space, ions%natoms, ions%mass, "lr", sys%namespace)
224
225 call epot_precalc_local_potential(sys%hm%ep, sys%namespace, sys%gr, sys%ions)
226
227 if (do_infrared) then
228 call born_charges_init(born, sys%namespace, ions%natoms, st%val_charge, st%qtot, ndim)
229 end if
230
231 call lr_init(lr(1))
232 call lr_allocate(lr(1), st, gr)
233
234 call restart_init(restart_dump, sys%namespace, restart_vib_modes, restart_type_dump, sys%mc, ierr, mesh=gr)
235 call restart_init(restart_load, sys%namespace, restart_vib_modes, restart_type_load, sys%mc, ierr, mesh=gr)
236
237 !CALCULATE
238
239 ! the ionic contribution, see second term in Eq. 85 in Baroni et al. RMP, 73, 515 (2001)
241
242 ! Compute the -<phi0 | v2 | phi0> term
243 if (states_are_real(st)) then
244 call dionic_pert_matrix_elements_2(sys%gr, sys%namespace, sys%space, sys%ions, sys%hm, 1, st, vib, vib%dyn_matrix)
245 else
246 call zionic_pert_matrix_elements_2(sys%gr, sys%namespace, sys%space, sys%ions, sys%hm, 1, st, vib, vib%dyn_matrix)
247 end if
248
249 if (fromscratch) then
250 start_mode = 1
251 else
252 call phonons_load(restart_load, vib, start_mode)
253 end if
254
255 ! Delete, if fromScratch, or trying to open it failed and there is something wrong with it.
256 if (start_mode == 1) call restart_rm(restart_dump, 'restart')
257
258 ! Output the first start_mode-1 modes
259 do imat = 1, start_mode - 1
260 call vibrations_out_dyn_matrix_row(vib, imat)
261 end do
262
263
264 ! The remaining term is <\psi_1 | v^{(1)} | \psi_0> + <\psi_0 | v^{(1)} | \psi_1>
265 ! We first get \psi_1 from Sternheimer, and then compute the term using X(forces_derivative)
266 ! Note that in order to avoid getting v^{(1)}, we perform an integration by part in this routine
267 pert => perturbation_ionic_t(sys%namespace, ions)
268
269 do imat = start_mode, vib%num_modes
270 iatom = vibrations_get_atom(vib, imat)
271 idir = vibrations_get_dir(vib, imat)
272
273 write(message(1),'(a,i5,a,a1,a)') &
274 "Calculating response to displacement of atom ", iatom, " in ", index2axis(idir), "-direction."
275 call messages_info(1, namespace=sys%namespace)
276
277 ! the converged wfns for the previous mode are probably not a good starting point
278 call lr_zero(lr(1), st)
279
280 if (.not. fromscratch) then
281 message(1) = "Loading restart wavefunctions for linear response."
282 call messages_info(1, namespace=sys%namespace)
283 call restart_open_dir(restart_load, wfs_tag_sigma(sys%namespace, phn_wfs_tag(iatom, idir), 1), ierr)
284 if (ierr == 0) then
285 call states_elec_load(restart_load, sys%namespace, sys%space, st, sys%gr, sys%kpoints, &
286 sys%st%restart_fixed_occ, ierr=ierr, lr = lr(1))
287 end if
288 if (ierr /= 0) then
289 message(1) = "Unable to read response wavefunctions from '"//&
290 trim(wfs_tag_sigma(sys%namespace, phn_wfs_tag(iatom, idir), 1))//"'."
291 call messages_warning(1, namespace=sys%namespace)
292 end if
293 call restart_close_dir(restart_load)
294 end if
295
296 call pert%setup_atom(iatom)
297 call pert%setup_dir(idir)
298
299 ! We now solve the Sternheimber equation in order to get the first-order change in the wavefunction
300 ! This is then used to get the force derivatives
301 safe_allocate(force_deriv(1:ndim, 1:natoms))
302 if (states_are_real(st)) then
303
304 call dsternheimer_solve(sh, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%mc, &
305 lr, 1, m_zero, pert, restart_dump, phn_rho_tag(iatom, idir), phn_wfs_tag(iatom, idir))
306
307 call dforces_derivative(gr, sys%namespace, sys%space, ions, sys%hm%ep, st, sys%kpoints, lr(1), lr(1), force_deriv, &
308 sys%hm%lda_u_level)
309
310 else
311
312 call zsternheimer_solve(sh, sys%namespace, sys%space, sys%gr, sys%kpoints, sys%st, sys%hm, sys%mc, &
313 lr, 1, m_z0, pert, restart_dump, phn_rho_tag(iatom, idir), phn_wfs_tag(iatom, idir))
314
315 call zforces_derivative(gr, sys%namespace, sys%space, ions, sys%hm%ep, st, sys%kpoints, lr(1), lr(1), force_deriv, &
316 sys%hm%lda_u_level)
317
318 end if
319
320 do jmat = 1, vib%num_modes
321 if (.not. symmetrize .and. jmat < imat) then
322 vib%dyn_matrix(jmat, imat) = vib%dyn_matrix(imat, jmat)
323 cycle
324 end if
325
326 jatom = vibrations_get_atom(vib, jmat)
327 jdir = vibrations_get_dir(vib, jmat)
328
329 vib%dyn_matrix(jmat, imat) = vib%dyn_matrix(jmat, imat) + real(force_deriv(jdir, jatom), real64)
330 vib%dyn_matrix(jmat, imat) = vib%dyn_matrix(jmat, imat) * vibrations_norm_factor(vib, iatom, jatom)
331 end do
332 safe_deallocate_a(force_deriv)
333
334 call vibrations_out_dyn_matrix_row(vib, imat)
335
336 if (do_infrared) then
337 if (states_are_real(st)) then
338 call dphonons_lr_infrared(gr, ions, st, lr(1), kdotp_lr, imat, iatom, idir, vib%infrared)
339 else
340 call zphonons_lr_infrared(gr, ions, st, lr(1), kdotp_lr, imat, iatom, idir, vib%infrared)
341 end if
342 end if
343
344 iunit_restart = restart_open(restart_dump, 'restart', position='append')
345 ! open and close makes sure output is not buffered
346 do jmat = 1, vib%num_modes
347 write(line(1), *) jmat, imat, vib%dyn_matrix(jmat, imat)
348 call restart_write(restart_dump, iunit_restart, line, 1, ierr)
349 if (ierr /= 0) then
350 message(1) = "Could not write restart information."
351 call messages_warning(1, namespace=sys%namespace)
352 end if
353 end do
354 write(line(1), *) imat, (vib%infrared(imat, idir), idir = 1, ndim)
355 call restart_write(restart_dump, iunit_restart, line, 1, ierr)
356 if (ierr /= 0) then
357 message(1) = "Could not write restart information."
358 call messages_warning(1, namespace=sys%namespace)
359 end if
360 call restart_close(restart_dump, iunit_restart)
361
362 message(1) = ""
363 call messages_info(1, namespace=sys%namespace)
364 end do
365
366 safe_deallocate_p(pert)
367
368 if (symmetrize) call vibrations_symmetrize_dyn_matrix(vib)
370 call vibrations_output(vib)
371 call axsf_mode_output(vib, ions, gr, sys%namespace)
372
373 if (do_infrared) then
374 if (sys%space%is_periodic() .and. .not. smear_is_semiconducting(st%smear)) then
375 message(1) = "Cannot calculate infrared intensities for periodic system with smearing (i.e. without a gap)."
376 call messages_info(1, namespace=sys%namespace)
377 else
378 call born_from_infrared(vib, born)
379 call born_output_charges(born, ions%atom, ions%charge, ions%natoms, sys%namespace, &
380 ndim, vib_modes_dir, write_real = .true.)
381 call calc_infrared()
382 end if
383
384 call born_charges_end(born)
385 end if
386
387 if (normal_mode_wfs) then
388 message(1) = "Calculating response wavefunctions for normal modes."
389 call messages_info(1, namespace=sys%namespace)
390 if (states_are_real(st)) then
391 call dphonons_lr_wavefunctions(lr(1), sys%namespace, sys%space, st, sys%gr, sys%kpoints, vib, restart_load, &
392 restart_dump)
393 else
394 call zphonons_lr_wavefunctions(lr(1), sys%namespace, sys%space, st, sys%gr, sys%kpoints, vib, restart_load, &
395 restart_dump)
396 end if
397 end if
398
399 !DESTRUCT
400
401 call lr_dealloc(lr(1))
402 call vibrations_end(vib)
403 call sternheimer_end(sh)
405 if (sys%space%is_periodic() .and. do_infrared) then
406 do idir = 1, sys%space%periodic_dim
407 call lr_dealloc(kdotp_lr(idir))
408 end do
409 end if
410 call restart_end(restart_load)
411 call restart_end(restart_dump)
412
413 pop_sub(phonons_lr_run_legacy)
414
415 contains
416
417 ! ---------------------------------------------------------
423 subroutine build_ionic_dyn_matrix()
424 real(real64) :: term, weight, xi(1:ndim), dx(1:ndim), r2
425
427
428 assert(.not. ions%space%is_periodic())
429
430 vib%dyn_matrix(:,:) = m_zero
431
432 do iatom = 1, natoms
433 xi = ions%pos(:, iatom)
434
435 do jatom = 1, natoms
436 if(iatom == jatom) cycle
437
438 dx = xi - ions%pos(:, jatom)
439 r2 = dot_product(dx, dx)
440
441 weight = ions%charge(iatom) * ions%charge(jatom) /(sqrt(r2)**3)
442
443 do idir = 1, ndim
444 do jdir = 1, ndim
445
446 term = weight * (ddelta(idir, jdir) - m_three*dx(idir)*dx(jdir)/r2)
447
448 ! The force is given by F_I = Z_I \sum_K Z_K (R_I-R_K)/|R_I-R_K|^3
449 ! There is therefore a diagonal term from the case J=I, with the sum over K
450 vib%dyn_matrix(vibrations_get_index(vib, iatom, jdir), vibrations_get_index(vib, iatom, idir)) = &
451 vib%dyn_matrix(vibrations_get_index(vib, iatom, jdir), vibrations_get_index(vib, iatom, idir)) + term
452
453 vib%dyn_matrix(vibrations_get_index(vib, jatom, jdir), vibrations_get_index(vib, iatom, idir)) = &
454 vib%dyn_matrix(vibrations_get_index(vib, jatom, jdir), vibrations_get_index(vib, iatom, idir)) - term
455 end do
456 end do
457 end do
458 end do
459
461 end subroutine build_ionic_dyn_matrix
462
463 ! ---------------------------------------------------------
465 subroutine calc_infrared()
466
467 integer :: iunit_ir
468 real(real64) :: lir(1:sys%space%dim+1)
469
471
472 iunit_ir = io_open(vib_modes_dir//'infrared', sys%namespace, action='write')
473
474 write(iunit_ir, '(a)', advance = 'no') '# freq ['//trim(units_abbrev(unit_invcm))//']'
475 do idir = 1, ndim
476 write(iunit_ir, '(a14)', advance = 'no') '<' // index2axis(idir) // '> [' // trim(units_abbrev(units_out%length)) // ']'
477 end do
478 write(iunit_ir, '(a14)') 'average [' // trim(units_abbrev(units_out%length)) // ']'
479
480 do iatom = 1, natoms
481 do idir = 1, ndim
482
483 imat = vibrations_get_index(vib, iatom, idir)
484
485 write(iunit_ir, '(f17.8)', advance = 'no') units_from_atomic(unit_invcm, vib%freq(imat))
486 do jdir = 1, ndim
487 lir(jdir) = dot_product(vib%infrared(:, jdir), vib%normal_mode(:, imat))
488 write(iunit_ir, '(f14.5)', advance = 'no') units_from_atomic(units_out%length, lir(jdir))
489 end do
490
491 lir(ndim+1) = norm2(lir(1:ndim))/sqrt(real(ndim, real64) )
492 write(iunit_ir, '(f17.8)') units_from_atomic(units_out%length, lir(ndim + 1))
493 end do
494 end do
495
496 call io_close(iunit_ir)
498 end subroutine calc_infrared
499
500 end subroutine phonons_lr_run_legacy
501
502
503 ! ---------------------------------------------------------
504 subroutine born_from_infrared(vib, born)
505 type(vibrations_t), intent(in) :: vib
506 type(born_charges_t), intent(inout) :: born
507
508 integer :: imat, idir, iatom
509
510 push_sub(born_from_infrared)
511
512 do imat = 1, vib%num_modes
513 idir = vibrations_get_dir(vib, imat)
514 iatom = vibrations_get_atom(vib, imat)
515 born%charge(1:vib%ndim, idir, iatom) = -vib%infrared(imat, 1:vib%ndim)
516 end do
517
518 pop_sub(born_from_infrared)
519 end subroutine born_from_infrared
520
521
522 ! ---------------------------------------------------------
523 character(len=100) function phn_rho_tag(iatom, dir) result(str)
524 integer, intent(in) :: iatom, dir
525
526 push_sub(phn_rho_tag)
527
528 write(str, '(a,i4.4,a,i1)') 'phn_rho_', iatom, '_', dir
529
530 pop_sub(phn_rho_tag)
531
532 end function phn_rho_tag
533
534
535 ! ---------------------------------------------------------
536 character(len=100) function phn_wfs_tag(iatom, dir) result(str)
537 integer, intent(in) :: iatom, dir
538
539 push_sub(phn_wfs_tag)
540
541 write(str, '(a,i4.4,a,a)') "phn_wfs_", iatom, "_", index2axis(dir)
542
543 pop_sub(phn_wfs_tag)
544
545 end function phn_wfs_tag
546
547
548 ! ---------------------------------------------------------
549 character(len=100) function phn_nm_wfs_tag(inm) result(str)
550 integer, intent(in) :: inm
551
552 push_sub(phn_nm_wfs_tag)
553
554 write(str, '(a,i5.5)') "phn_nm_wfs_", inm
555
556 pop_sub(phn_nm_wfs_tag)
557
558 end function phn_nm_wfs_tag
559
560
561 ! ---------------------------------------------------------
563 subroutine axsf_mode_output(this, ions, mesh, namespace)
564 type(vibrations_t), intent(in) :: this
565 type(ions_t), intent(in) :: ions
566 class(mesh_t), intent(in) :: mesh
567 type(namespace_t), intent(in) :: namespace
568
569 integer :: iunit, iatom, idir, imat, jmat
570 real(real64), allocatable :: forces(:,:)
571 character(len=2) :: suffix
572
573 if (.not. mpi_grp_is_root(mpi_world)) return
574
575 push_sub(axsf_mode_output)
576
577 ! for some reason, direct usage of this%suffix gives an odd result
578 suffix = vibrations_get_suffix(this)
579 iunit = io_open(vib_modes_dir//'normal_modes_'//suffix//'.axsf', namespace, action='write')
580
581 write(iunit, '(a,i6)') 'ANIMSTEPS ', this%num_modes
582 safe_allocate(forces(1:ions%space%dim, 1:ions%natoms))
583 do imat = 1, this%num_modes
584 do jmat = 1, this%num_modes
585 iatom = vibrations_get_atom(this, jmat)
586 idir = vibrations_get_dir(this, jmat)
587 forces(idir, iatom) = this%normal_mode(jmat, imat)
588 end do
589 call write_xsf_geometry(iunit, ions%space, ions%latt, ions%pos, ions%atom, mesh, forces = forces, index = imat)
590 end do
591 safe_deallocate_a(forces)
592 call io_close(iunit)
593
594 pop_sub(axsf_mode_output)
595 end subroutine axsf_mode_output
596
597 ! ---------------------------------------------------------
599 subroutine phonons_load(restart, vib, start_mode)
600 type(restart_t), intent(in) :: restart
601 type(vibrations_t), intent(inout) :: vib
602 integer, intent(out) :: start_mode
603
604 integer :: iunit, ierr, imode, jmode, imode_read, jmode_read
605 character(len=120) :: line(1)
606
607 push_sub(phonons_load)
608
609 iunit = restart_open(restart, 'restart')
610 if (iunit /= -1) then
611 imode_loop: do imode = 1, vib%num_modes
612 do jmode = 1, vib%num_modes
613 call restart_read(restart, iunit, line, 1, ierr)
614 if (ierr /= 0) exit imode_loop
615 read(line(1), fmt=*, iostat=ierr) jmode_read, imode_read, vib%dyn_matrix(jmode, imode)
616 if (imode_read /= imode) then
617 write(message(1),'(a,i9,a,i9)') "Corruption of restart data: row ", imode, " is labeled as ", imode_read
618 call messages_fatal(1)
619 end if
620 if (jmode_read /= jmode) then
621 write(message(1),'(a,i9,a,i9)') "Corruption of restart data: column ", jmode, " is labeled as ", jmode_read
622 call messages_fatal(1)
623 end if
624 end do
625
626 call restart_read(restart, iunit, line, 1, ierr)
627 if (ierr /= 0) exit
628
629 start_mode = imode + 1
630
631 read(line(1), fmt=*, iostat=ierr) imode_read, vib%infrared(imode, 1:vib%ndim)
632 if (imode_read /= imode) then
633 write(message(1),'(a,i9,a,i9)') "Corruption of restart data: infrared row ", imode, " is labeled as ", imode_read
634 call messages_fatal(1)
635 end if
636 end do imode_loop
637
638 write(message(1),'(a,i9,a,i9)') 'Info: Read saved dynamical-matrix rows for ', &
639 start_mode - 1, ' modes out of ', vib%num_modes
640 call messages_info(1)
641
642 call restart_close(restart, iunit)
643 else
644 start_mode = 1
645
646 message(1) = "Could not open restart file 'restart'. Starting from scratch."
647 call messages_warning(1)
648 end if
649
650 pop_sub(phonons_load)
651 end subroutine phonons_load
652
653#include "complex.F90"
654#include "phonons_lr_inc.F90"
655
656#include "undef.F90"
657
658#include "real.F90"
659#include "phonons_lr_inc.F90"
660
661end module phonons_lr_oct_m
662
663!! Local Variables:
664!! mode: f90
665!! coding: utf-8
666!! End:
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)
subroutine, public epot_precalc_local_potential(ep, namespace, gr, ions)
Definition: epot.F90:574
subroutine, public zforces_derivative(gr, namespace, space, ions, ep, st, kpoints, lr, lr2, force_deriv, lda_u_level)
Definition: forces.F90:1892
subroutine, public dforces_derivative(gr, namespace, space, ions, ep, st, kpoints, lr, lr2, force_deriv, lda_u_level)
Definition: forces.F90:1268
real(real64), parameter, public m_zero
Definition: global.F90:188
character(len= *), parameter, public vib_modes_dir
Definition: global.F90:256
complex(real64), parameter, public m_z0
Definition: global.F90:198
real(real64), parameter, public m_three
Definition: global.F90:191
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public write_xsf_geometry(iunit, space, latt, pos, atoms, mesh, forces, index)
for format specification see: http:
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
character(len=100) function, public kdotp_wfs_tag(dir, dir2)
Definition: kdotp_calc.F90:152
subroutine, public lr_zero(lr, st)
subroutine, public lr_allocate(lr, st, mesh, allocate_rho)
subroutine, public lr_init(lr)
subroutine, public lr_dealloc(lr)
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
real(real64) pure function, public ddelta(i, j)
Definition: math.F90:602
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
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.
subroutine, public zionic_pert_matrix_elements_2(gr, namespace, space, ions, hm, ik, st, vib, matrix)
Computes the second order term.
subroutine, public dionic_pert_matrix_elements_2(gr, namespace, space, ions, hm, ik, st, vib, matrix)
Computes the second order term.
subroutine dphonons_lr_infrared(mesh, ions, st, lr, kdotp_lr, imat, iatom, idir, infrared)
Definition: phonons_lr.F90:970
subroutine zphonons_lr_wavefunctions(lr, namespace, space, st, mesh, kpoints, vib, restart_load, restart_dump)
calculate the wavefunction associated with each normal mode
Definition: phonons_lr.F90:830
subroutine, public phonons_lr_run(system, from_scratch)
Definition: phonons_lr.F90:170
subroutine zphonons_lr_infrared(mesh, ions, st, lr, kdotp_lr, imat, iatom, idir, infrared)
Definition: phonons_lr.F90:790
subroutine born_from_infrared(vib, born)
Definition: phonons_lr.F90:598
character(len=100) function, public phn_nm_wfs_tag(inm)
Definition: phonons_lr.F90:643
subroutine phonons_load(restart, vib, start_mode)
Load restart information for a linear-response phonon calculation.
Definition: phonons_lr.F90:693
subroutine dphonons_lr_wavefunctions(lr, namespace, space, st, mesh, kpoints, vib, restart_load, restart_dump)
calculate the wavefunction associated with each normal mode
subroutine phonons_lr_run_legacy(sys, fromscratch)
Definition: phonons_lr.F90:188
subroutine, public axsf_mode_output(this, ions, mesh, namespace)
output eigenvectors as animated XSF file, one per frame, displacements as forces
Definition: phonons_lr.F90:657
character(len=100) function, public phn_rho_tag(iatom, dir)
Definition: phonons_lr.F90:617
character(len=100) function, public phn_wfs_tag(iatom, dir)
Definition: phonons_lr.F90:630
subroutine, public restart_read(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:935
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
Definition: restart.F90:952
integer, parameter, public restart_kdotp
Definition: restart.F90:200
subroutine, public restart_rm(restart, name)
Remove directory or file "name" that is located inside the current restart directory.
Definition: restart.F90:839
integer, parameter, public restart_gs
Definition: restart.F90:200
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:516
integer, parameter, public restart_type_dump
Definition: restart.F90:245
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:908
integer, parameter, public restart_vib_modes
Definition: restart.F90:200
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
Definition: restart.F90:861
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:772
integer, parameter, public restart_type_load
Definition: restart.F90:245
subroutine, public restart_end(restart)
Definition: restart.F90:722
subroutine, public restart_close_dir(restart)
Change back to the base directory. To be called after restart_open_dir.
Definition: restart.F90:806
logical pure function, public smear_is_semiconducting(this)
Definition: smear.F90:845
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, fixed_occ, is_complex, packed)
subroutine, public states_elec_load(restart, namespace, space, st, mesh, kpoints, fixed_occ, 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.
character(len=100) function, public wfs_tag_sigma(namespace, base_name, isigma)
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)
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_invcm
For vibrational frequencies.
type(unit_system_t), public units_out
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
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:691
character(len=2) pure function, public vibrations_get_suffix(this)
Definition: vibrations.F90:222
real(real64) pure function, public vibrations_norm_factor(this, iatom, jatom)
Definition: vibrations.F90:291
subroutine, public vibrations_diag_dyn_matrix(this)
Diagonalize the dynamical matrix.
Definition: vibrations.F90:350
subroutine, public vibrations_out_dyn_matrix_row(this, imat)
Outputs one row of the dynamical matrix.
Definition: vibrations.F90:303
subroutine, public vibrations_init(this, space, natoms, mass, suffix, namespace)
Definition: vibrations.F90:168
integer pure function, public vibrations_get_dir(this, index)
Definition: vibrations.F90:403
subroutine, public vibrations_symmetrize_dyn_matrix(this)
Symmetrize the dynamical matric, which is real symmetric matrix.
Definition: vibrations.F90:231
integer pure function, public vibrations_get_index(this, iatom, idim)
Definition: vibrations.F90:384
subroutine, public vibrations_output(this)
Outputs the eigenvectors and eigenenergies of the dynamical matrix.
Definition: vibrations.F90:413
subroutine, public vibrations_end(this)
Definition: vibrations.F90:207
integer pure function, public vibrations_get_atom(this, index)
Definition: vibrations.F90:394
subroutine calc_infrared()
calculate infrared intensities
Definition: phonons_lr.F90:559
subroutine build_ionic_dyn_matrix()
Computes the ionic contribution to the dynamical matrix.
Definition: phonons_lr.F90:517
Class describing the electron system.
Definition: electrons.F90:218
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Container class for lists of system_oct_m::system_t.
int true(void)