Octopus
static_pol.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20#define RESTART_FILE 'dipoles'
21
24 use debug_oct_m
25 use elf_oct_m
27 use global_oct_m
29 use io_oct_m
31 use ions_oct_m
32 use, intrinsic :: iso_fortran_env
33 use lcao_oct_m
36 use mpi_oct_m
39 use parser_oct_m
42 use scf_oct_m
43 use space_oct_m
48 use unit_oct_m
50 use utils_oct_m
51 use v_ks_oct_m
52
53 implicit none
54
55 private
56 public :: &
58
59
60contains
61
62 ! ---------------------------------------------------------
63 subroutine static_pol_run(system, from_scratch)
64 class(*), intent(inout) :: system
65 logical, intent(in) :: from_scratch
66
67 push_sub(static_pol_run)
68
69 select type (system)
70 class is (multisystem_basic_t)
71 message(1) = "CalculationMode = static_pol not implemented for multi-system calculations"
72 call messages_fatal(1, namespace=system%namespace)
73 type is (electrons_t)
74 call static_pol_run_legacy(system, from_scratch)
75 end select
76
77 pop_sub(static_pol_run)
78 end subroutine static_pol_run
79
80 ! ---------------------------------------------------------
81 subroutine static_pol_run_legacy(sys, fromScratch)
82 type(electrons_t), intent(inout) :: sys
83 logical, intent(in) :: fromScratch
84
85 type(scf_t) :: scfv
86 integer :: iunit, ios, i_start, ii, jj, is, isign, ierr, read_count, verbosity
87 real(real64) :: e_field, e_field_saved
88 real(real64), allocatable :: Vpsl_save(:), trrho(:), dipole(:, :, :)
89 real(real64), allocatable :: elf(:,:), lr_elf(:,:), elfd(:,:), lr_elfd(:,:)
90 real(real64), allocatable :: lr_rho(:,:), lr_rho2(:,:), gs_rho(:,:), tmp_rho(:,:)
91 real(real64) :: center_dipole(1:sys%space%dim), diag_dipole(1:sys%space%dim), ionic_dipole(1:sys%space%dim), &
92 print_dipole(1:sys%space%dim)
93 type(born_charges_t) :: born_charges
94 logical :: calc_born, start_density_is_zero_field, write_restart_densities, calc_diagonal, verbose
95 logical :: diagonal_done, center_written, fromScratch_local, field_written
96 character(len=MAX_PATH_LEN) :: fname, dir_name
97 character(len=120) :: line(1)
98 character :: sign_char
99 type(restart_t) :: gs_restart, restart_load, restart_dump
100
101 push_sub(static_pol_run_legacy)
102
103 if (sys%hm%pcm%run_pcm) then
104 call messages_not_implemented("PCM for CalculationMode /= gs or td", namespace=sys%namespace)
105 end if
106
107 if (sys%kpoints%use_symmetries) then
108 call messages_experimental("KPoints symmetries with CalculationMode = em_resp", namespace=sys%namespace)
109 end if
110
111 call init_()
112
113 ! load wavefunctions
114 call restart_init(gs_restart, sys%namespace, restart_gs, restart_type_load, sys%mc, ierr, mesh=sys%gr, exact=.true.)
115 if (ierr == 0) then
116 call states_elec_load(gs_restart, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, ierr)
117 end if
118 if (ierr /= 0) then
119 message(1) = "Unable to read wavefunctions."
120 call messages_fatal(1, namespace=sys%namespace)
121 end if
122 call restart_end(gs_restart)
123
124 if (sys%space%is_periodic()) then
125 message(1) = "Electric field cannot be applied to a periodic system (currently)."
126 call messages_fatal(1, namespace=sys%namespace)
127 end if
128
129 ! set up Hamiltonian
130 message(1) = 'Info: Setting up Hamiltonian.'
131 call messages_info(1, namespace=sys%namespace)
132 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
133 sys%hm, calc_eigenval = .false.) ! we read them from restart
134
135 ! Allocate the dipole
136 safe_allocate(dipole(1:sys%space%dim, 1:sys%space%dim, 1:2))
137 dipole = m_zero
138
139 i_start = 1
140 center_written = .false.
141 diagonal_done = .false.
142 field_written = .false.
143
144 call restart_init(restart_dump, sys%namespace, restart_em_resp_fd, restart_type_dump, sys%mc, ierr, mesh=sys%gr)
145
146 if (.not. fromscratch) then
147 call restart_init(restart_load, sys%namespace, restart_em_resp_fd, restart_type_load, sys%mc, ierr, mesh=sys%gr)
148 if (ierr == 0) then
149 iunit = restart_open(restart_load, restart_file)
150 else
151 iunit = 0
152 end if
153
154 if (iunit > 0) then
155 ! Finds out how many dipoles have already been written.
156 rewind(iunit)
157
158 read(iunit, fmt=*, iostat = ios) e_field_saved
159 field_written = (ios == 0)
160
161 read(iunit, fmt=*, iostat = ios) (center_dipole(jj), jj = 1, sys%space%dim)
162 center_written = (ios == 0)
163
164 do ii = 1, 3
165 read(iunit, fmt=*, iostat = ios) ((dipole(ii, jj, isign), jj = 1, sys%space%dim), isign = 1, 2)
166 if (ios /= 0) exit
167 i_start = i_start + 1
168 end do
169
170 read(iunit, fmt=*, iostat = ios) (diag_dipole(jj), jj = 1, sys%space%dim)
171 diagonal_done = (ios == 0)
172
173 call restart_close(restart_load, iunit)
174 end if
175
176 ! if saved dipoles used a different e_field, we cannot use them
177 if (field_written .and. abs(e_field_saved - e_field) > 1e-15_real64) then
178 message(1) = "Saved dipoles are from a different electric field, cannot use them."
179 call messages_warning(1, namespace=sys%namespace)
180 center_written = .false.
181 diagonal_done = .false.
182 i_start = 1
183 end if
184
185 read_count = (i_start - 1) * 2
186 if (center_written) read_count = read_count + 1
187 if (diagonal_done) read_count = read_count + 1
188 write(message(1),'(a,i1,a)') "Using ", read_count, " dipole(s) from file."
189 call messages_info(1, namespace=sys%namespace)
190 end if
191
192 if (sys%outp%what(option__output__density) .or. &
193 sys%outp%what(option__output__pol_density)) then
194 if (i_start > 2 .and. calc_diagonal) then
195 i_start = 2
196 diagonal_done = .false.
197 !FIXME: take derivatives between yz and z (not y) so can restart from only last (z) calc
198 end if
199 end if
200
201 if (i_start == 1) then
202 ! open new file, erase old data, write e_field
203 iunit = restart_open(restart_dump, restart_file, status='replace')
204 write(line(1), fmt='(e20.12)') e_field
205 call restart_write(restart_dump, iunit, line, 1, ierr)
206 if (ierr /= 0) then
207 message(1) = "Unsuccessful write of electric field."
208 call messages_warning(1, namespace=sys%namespace)
209 end if
210 call restart_close(restart_dump, iunit)
211 center_written = .false.
212 end if
213
214 ! Save local potential
215 safe_allocate(vpsl_save(1:sys%gr%np))
216 vpsl_save = sys%hm%ep%Vpsl
217
218 ! Allocate the trrho to contain the trace of the density.
219 safe_allocate(trrho(1:sys%gr%np))
220 safe_allocate(gs_rho(1:sys%gr%np, 1:sys%st%d%nspin))
221 safe_allocate(tmp_rho(1:sys%gr%np, 1:sys%st%d%nspin))
222 trrho = m_zero
223 gs_rho = m_zero
224
225 call output_init_()
226 call scf_init(scfv, sys%namespace, sys%gr, sys%ions, sys%st, sys%mc, sys%hm, sys%space)
227 call born_charges_init(born_charges, sys%namespace, sys%ions%natoms, sys%st%val_charge, &
228 sys%st%qtot, sys%space%dim)
229
230 ! now calculate the dipole without field
231
232 sys%hm%ep%vpsl(1:sys%gr%np) = vpsl_save(1:sys%gr%np)
233 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners)
234
235 write(message(1), '(a)')
236 write(message(2), '(a)') 'Info: Calculating dipole moment for zero field.'
237 call messages_info(2, namespace=sys%namespace)
238 call scf_run(scfv, sys%namespace, sys%space, sys%mc, sys%gr, sys%ions, sys%ext_partners, sys%st, &
239 sys%ks, sys%hm, verbosity = verbosity)
240
241 gs_rho(1:sys%gr%np, 1:sys%st%d%nspin) = sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin)
242 trrho = m_zero
243 do is = 1, sys%st%d%spin_channels
244 trrho(1:sys%gr%np) = trrho(1:sys%gr%np) + gs_rho(1:sys%gr%np, is)
245 end do
246
247 ! calculate dipole
248 do jj = 1, sys%space%dim
249 center_dipole(jj) = dmf_moment(sys%gr, trrho, jj, 1)
250 end do
251
252 ! Writes the dipole to file
253 if (.not. center_written) then
254 iunit = restart_open(restart_dump, restart_file, position='append')
255 write(line(1), fmt='(6e20.12)') (center_dipole(jj), jj = 1, sys%space%dim)
256 call restart_write(restart_dump, iunit, line, 1, ierr)
257 if (ierr /= 0) then
258 message(1) = "Unsuccessful write of center dipole."
259 call messages_warning(1, namespace=sys%namespace)
260 end if
261 call restart_close(restart_dump, iunit)
262 end if
263
264 if (mpi_grp_is_root(mpi_world)) then
265 ionic_dipole(1:sys%space%dim) = sys%ions%dipole()
266 print_dipole(1:sys%space%dim) = center_dipole(1:sys%space%dim) + ionic_dipole(1:sys%space%dim)
267 call output_dipole(print_dipole, sys%space%dim, namespace=sys%namespace)
268 end if
269
270 do ii = i_start, sys%space%dim
271 do isign = 1, 2
272 write(message(1), '(a)')
273 write(message(2), '(a,f6.4,5a)') 'Info: Calculating dipole moment for field ', &
274 units_from_atomic(units_out%force, (-1)**isign * e_field), ' ', &
275 trim(units_abbrev(units_out%force)), ' in the ', index2axis(ii), '-direction.'
276 call messages_info(2, namespace=sys%namespace)
277 ! there would be an extra factor of -1 in here that is for the electronic charge
278 ! except that we treat electrons as positive
279
280 sys%hm%ep%vpsl(1:sys%gr%np) = vpsl_save(1:sys%gr%np) + (-1)**isign * sys%gr%x(1:sys%gr%np, ii) * e_field
281 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners)
282
283 if (isign == 1) then
284 sign_char = '+'
285 else
286 sign_char = '-'
287 end if
288
289 write(dir_name,'(a)') "field_"//index2axis(ii)//sign_char
290 fromscratch_local = fromscratch
291
292 if (.not. fromscratch) then
293 call restart_open_dir(restart_load, trim(dir_name), ierr)
294 if (ierr == 0) then
295 call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, ierr)
296 end if
297 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
298 if (ierr /= 0) fromscratch_local = .true.
299 call restart_close_dir(restart_load)
300 end if
301
302 if (fromscratch_local) then
303 if (start_density_is_zero_field) then
304 sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin) = gs_rho(1:sys%gr%np, 1:sys%st%d%nspin)
305 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
306 else
307 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
308 sys%hm, lmm_r = scfv%lmm_r)
309 end if
310 end if
311
312 call scf_mix_clear(scfv)
313 call scf_run(scfv, sys%namespace, sys%space, sys%mc, sys%gr, sys%ions, sys%ext_partners, sys%st, &
314 sys%ks, sys%hm, verbosity = verbosity)
315
316 trrho = m_zero
317 do is = 1, sys%st%d%spin_channels
318 trrho(1:sys%gr%np) = trrho(1:sys%gr%np) + sys%st%rho(1:sys%gr%np, is)
319 end do
320
321 ! calculate dipole
322 do jj = 1, sys%space%dim
323 dipole(ii, jj, isign) = dmf_moment(sys%gr, trrho, jj, 1)
324 end do
325
326 if (mpi_grp_is_root(mpi_world)) then
327 print_dipole(1:sys%space%dim) = dipole(ii, 1:sys%space%dim, isign) + ionic_dipole(1:sys%space%dim)
328 call output_dipole(print_dipole, sys%space%dim, namespace=sys%namespace)
329 end if
330
331 call output_cycle_()
332
333 if (write_restart_densities) then
334 call restart_open_dir(restart_dump, trim(dir_name), ierr)
335 if (ierr == 0) then
336 call states_elec_dump(restart_dump, sys%space, sys%st, sys%gr, sys%kpoints, ierr)
337 end if
338 call restart_close_dir(restart_dump)
339 if (ierr /= 0) then
340 message(1) = 'Unable to write states wavefunctions.'
341 call messages_warning(1, namespace=sys%namespace)
342 end if
343 end if
344 end do
345
346 ! Writes the dipole to file
347 iunit = restart_open(restart_dump, restart_file, position='append')
348 write(line(1), '(6e20.12)') ((dipole(ii, jj, isign), jj = 1, sys%space%dim), isign = 1, 2)
349 call restart_write(restart_dump, iunit, line, 1, ierr)
350 if (ierr /= 0) then
351 message(1) = "Unsuccessful write of dipole."
352 call messages_warning(1, namespace=sys%namespace)
353 end if
354 call restart_close(restart_dump, iunit)
355 end do
356
357 if (.not. diagonal_done .and. calc_diagonal) then
358 write(message(1), '(a)')
359 write(message(2), '(a,f6.4,3a, f6.4, 3a)') 'Info: Calculating dipole moment for field ', &
360 units_from_atomic(units_out%force, e_field), ' ', &
361 trim(units_abbrev(units_out%force)), ' in the '//index2axis(2)//'-direction plus ', &
362 units_from_atomic(units_out%force, e_field), ' ', &
363 trim(units_abbrev(units_out%force)), ' in the '//index2axis(3)//'-direction.'
364 call messages_info(2, namespace=sys%namespace)
365
366 sys%hm%ep%vpsl(1:sys%gr%np) = vpsl_save(1:sys%gr%np) &
367 - (sys%gr%x(1:sys%gr%np, 2) + sys%gr%x(1:sys%gr%np, 3)) * e_field
368 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners)
369
370 if (isign == 1) then
371 sign_char = '+'
372 else
373 sign_char = '-'
374 end if
375
376 fromscratch_local = fromscratch
377
378 if (.not. fromscratch) then
379 call restart_open_dir(restart_load, "field_yz+", ierr)
380 if (ierr == 0) then
381 call states_elec_load(restart_load, sys%namespace, sys%space, sys%st, sys%gr, sys%kpoints, ierr)
382 end if
383 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
384 if (ierr /= 0) fromscratch_local = .true.
385 call restart_close_dir(restart_load)
386 end if
387
388 if (fromscratch_local) then
389 if (start_density_is_zero_field) then
390 sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin) = gs_rho(1:sys%gr%np, 1:sys%st%d%nspin)
391 call v_ks_h_setup(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, sys%hm)
392 else
393 call lcao_run(sys%namespace, sys%space, sys%gr, sys%ions, sys%ext_partners, sys%st, sys%ks, &
394 sys%hm, lmm_r = scfv%lmm_r)
395 end if
396 end if
397
398 call scf_mix_clear(scfv)
399 call scf_run(scfv, sys%namespace, sys%space, sys%mc, sys%gr, sys%ions, sys%ext_partners, sys%st, &
400 sys%ks, sys%hm, verbosity = verbosity)
401
402 trrho = m_zero
403 do is = 1, sys%st%d%spin_channels
404 trrho(1:sys%gr%np) = trrho(1:sys%gr%np) + sys%st%rho(1:sys%gr%np, is)
405 end do
406
407 ! calculate dipole
408 do jj = 1, sys%space%dim
409 diag_dipole(jj) = dmf_moment(sys%gr, trrho, jj, 1)
410 end do
411
412 if (mpi_grp_is_root(mpi_world)) then
413 print_dipole(1:sys%space%dim) = diag_dipole(1:sys%space%dim) + ionic_dipole(1:sys%space%dim)
414 call output_dipole(print_dipole, sys%space%dim, namespace=sys%namespace)
415 end if
416
417 ! Writes the dipole to file
418 iunit = restart_open(restart_dump, restart_file, position='append')
419 write(line(1), fmt='(3e20.12)') (diag_dipole(jj), jj = 1, sys%space%dim)
420 call restart_write(restart_dump, iunit, line, 1, ierr)
421 if (ierr /= 0) then
422 message(1) = "Unsuccessful write of dipole."
423 call messages_warning(1, namespace=sys%namespace)
424 end if
425 call restart_close(restart_dump, iunit)
426
427 if (write_restart_densities) then
428 call restart_open_dir(restart_dump, "field_yz+", ierr)
429 if (ierr == 0) then
430 call states_elec_dump(restart_dump, sys%space, sys%st, sys%gr, sys%kpoints, ierr)
431 end if
432 call restart_close_dir(restart_dump)
433 if (ierr /= 0) then
434 message(1) = 'Unable to write states wavefunctions.'
435 call messages_warning(1, namespace=sys%namespace)
436 end if
437 end if
438
439 end if
440
441 if (.not. fromscratch) call restart_end(restart_load)
442 call scf_end(scfv)
443 call output_end_()
444 call born_charges_end(born_charges)
445
446 safe_deallocate_a(vpsl_save)
447 safe_deallocate_a(trrho)
448 safe_deallocate_a(gs_rho)
449 safe_deallocate_a(tmp_rho)
450 safe_deallocate_a(dipole)
451 call end_()
452 pop_sub(static_pol_run_legacy)
453
454 contains
455
456 ! ---------------------------------------------------------
457 subroutine init_()
459
460 call states_elec_allocate_wfns(sys%st, sys%gr)
461
462 call messages_obsolete_variable(sys%namespace, "EMStaticField", "EMStaticElectricField")
463 !%Variable EMStaticElectricField
464 !%Type float
465 !%Default 0.01 a.u.
466 !%Section Linear Response::Static Polarization
467 !%Description
468 !% Magnitude of the static electric field used to calculate the static polarizability,
469 !% if <tt>ResponseMethod = finite_differences</tt>.
470 !%End
471 call parse_variable(sys%namespace, 'EMStaticElectricField', 0.01_real64, e_field, units_inp%force)
472 if (e_field <= m_zero) then
473 write(message(1), '(a,e14.6,a)') "Input: '", e_field, "' is not a valid EMStaticElectricField."
474 message(2) = '(Must have EMStaticElectricField > 0)'
475 call messages_fatal(2, namespace=sys%namespace)
476 end if
477
478 ! variable defined in em_resp
479 call parse_variable(sys%namespace, 'EMCalcbornCharges', .false., calc_born)
480 if (calc_born) call messages_experimental("Calculation of born effective charges", namespace=sys%namespace)
481
482 !%Variable EMStartDensityIsZeroField
483 !%Type logical
484 !%Default true
485 !%Section Linear Response::Static Polarization
486 !%Description
487 !% Use the charge density from the zero-field calculation as the starting density for
488 !% SCF calculations with applied fields. For small fields, this will be fastest.
489 !% If there is trouble converging with larger fields, set to false,
490 !% to initialize the calculation for each field from scratch, as specified by the LCAO variables.
491 !% Only applies if <tt>ResponseMethod = finite_differences</tt>.
492 !%End
493 call parse_variable(sys%namespace, 'EMStartDensityIsZeroField', .true., start_density_is_zero_field)
494
495 !%Variable EMCalcDiagonalField
496 !%Type logical
497 !%Default true
498 !%Section Linear Response::Static Polarization
499 !%Description
500 !% Calculate <i>yz</i>-field for <math>\beta_{xyz}</math> hyperpolarizability, which is sometimes harder to converge.
501 !% Only applies if <tt>ResponseMethod = finite_differences</tt>.
502 !%End
503 call parse_variable(sys%namespace, 'EMCalcDiagonalField', .true., calc_diagonal)
504
505 !%Variable EMWriteRestartDensities
506 !%Type logical
507 !%Default true
508 !%Section Linear Response::Static Polarization
509 !%Description
510 !% Write density after each calculation for restart, rather than just the resulting electronic dipole moment.
511 !% Only applies if <tt>ResponseMethod = finite_differences</tt>. Restarting from calculations at smaller
512 !% fields can be helpful if there are convergence problems.
513 !%End
514 call parse_variable(sys%namespace, 'EMWriteRestartDensities', .true., write_restart_densities)
515
516 !%Variable EMVerbose
517 !%Type logical
518 !%Default false
519 !%Section Linear Response::Static Polarization
520 !%Description
521 !% Write full SCF output.
522 !% Only applies if <tt>ResponseMethod = finite_differences</tt>.
523 !%End
524 call parse_variable(sys%namespace, 'EMVerbose', .false., verbose)
525
526 if (verbose) then
527 verbosity = verb_full
528 else
529 verbosity = verb_compact
530 end if
531
533 end subroutine init_
534
535 ! ---------------------------------------------------------
536 subroutine end_()
537
538 push_sub(end_)
539
540 call states_elec_deallocate_wfns(sys%st)
541
542 pop_sub(end_)
543 end subroutine end_
544
545
546 !-------------------------------------------------------------
547 subroutine output_init_()
548 push_sub(output_init_)
549
550 !allocate memory for what we want to output
551 if (sys%outp%what(option__output__density) .or. &
552 sys%outp%what(option__output__pol_density)) then
553 safe_allocate(lr_rho(1:sys%gr%np, 1:sys%st%d%nspin))
554 safe_allocate(lr_rho2(1:sys%gr%np, 1:sys%st%d%nspin))
555 end if
556
557 if (sys%outp%what(option__output__elf)) then
558 safe_allocate( elf(1:sys%gr%np, 1:sys%st%d%nspin))
559 safe_allocate( lr_elf(1:sys%gr%np, 1:sys%st%d%nspin))
560 safe_allocate( elfd(1:sys%gr%np, 1:sys%st%d%nspin))
561 safe_allocate(lr_elfd(1:sys%gr%np, 1:sys%st%d%nspin))
562 end if
563
564 pop_sub(output_init_)
565 end subroutine output_init_
566
567
568 !-------------------------------------------------------------
569 subroutine output_cycle_()
570 integer :: iatom
571 type(unit_t) :: fn_unit
572
573 push_sub(output_cycle_)
574
575 ! BORN CHARGES
576 if (calc_born) then
577 do iatom = 1, sys%ions%natoms
578 if (isign == 1) then
579 ! temporary assignment for use in next cycle when isign == 2
580 born_charges%charge(ii, 1:sys%space%dim, iatom) = sys%ions%tot_force(:, iatom)
581 else
582 born_charges%charge(ii, 1:sys%space%dim, iatom) = &
583 (sys%ions%tot_force(:, iatom) - born_charges%charge(ii, 1:sys%space%dim, iatom)) &
584 / (m_two*e_field)
585 born_charges%charge(ii, ii, iatom) = born_charges%charge(ii, ii, iatom) + sys%ions%charge(iatom)
586 ! since the efield is applied in the SCF calculation by just altering the external potential felt by the electrons,
587 ! the ionic force due to the efield is not included in the forces returned by the SCF run, and so the ionic
588 ! contribution to the born charge must be added by hand here
589 end if
590 end do
591 end if
592
593 !DENSITY AND POLARIZABILITY DENSITY
594 if (sys%outp%what(option__output__density) .or. &
595 sys%outp%what(option__output__pol_density)) then
596
597 if (isign == 1 .and. ii == 2) then
598 tmp_rho(1:sys%gr%np, 1:sys%st%d%nspin) = sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin)
599 ! for use in off-diagonal non-linear densities
600 end if
601
602 if (isign == 1) then
603 ! temporary assignment for use in next cycle when isign == 2
604 lr_rho(1:sys%gr%np, 1:sys%st%d%nspin) = sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin)
605 else
606 lr_rho2(1:sys%gr%np, 1:sys%st%d%nspin) = &
607 -(sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin) + lr_rho(1:sys%gr%np, 1:sys%st%d%nspin) - &
608 2 * gs_rho(1:sys%gr%np, 1:sys%st%d%nspin)) / e_field**2
609
610 lr_rho(1:sys%gr%np, 1:sys%st%d%nspin) = &
611 (sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin) - lr_rho(1:sys%gr%np, 1:sys%st%d%nspin)) / (m_two*e_field)
612
613 !write
614 do is = 1, sys%st%d%nspin
615 if (sys%outp%what(option__output__density)) then
616 fn_unit = units_out%length**(1-sys%space%dim) / units_out%energy
617 write(fname, '(a,i1,2a)') 'fd_density-sp', is, '-', index2axis(ii)
618 call dio_function_output(sys%outp%how(option__output__density), em_resp_fd_dir, trim(fname),&
619 sys%namespace, sys%space, sys%gr, lr_rho(:, is), fn_unit, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
620
621 ! save the trouble of writing many copies of each density, since ii,jj = jj,ii
622 fn_unit = units_out%length**(2-sys%space%dim) / units_out%energy**2
623 do jj = ii, sys%space%dim
624 write(fname, '(a,i1,4a)') 'fd2_density-sp', is, '-', index2axis(ii), '-', index2axis(jj)
625 call dio_function_output(sys%outp%how(option__output__density), em_resp_fd_dir, trim(fname),&
626 sys%namespace, sys%space, sys%gr, lr_rho2(:, is), fn_unit, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
627 end do
628 end if
630 if (sys%outp%what(option__output__pol_density)) then
631 do jj = ii, sys%space%dim
632 fn_unit = units_out%length**(2-sys%space%dim) / units_out%energy
633 write(fname, '(a,i1,4a)') 'alpha_density-sp', is, '-', index2axis(ii), '-', index2axis(jj)
634 call dio_function_output(sys%outp%how(option__output__pol_density), em_resp_fd_dir, trim(fname), &
635 sys%namespace, sys%space, sys%gr, -sys%gr%x(:, jj) * lr_rho(:, is), &
636 fn_unit, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
637
638 fn_unit = units_out%length**(3-sys%space%dim) / units_out%energy**2
639 write(fname, '(a,i1,6a)') 'beta_density-sp', is, '-', index2axis(ii), &
640 '-', index2axis(ii), '-', index2axis(jj)
641 call dio_function_output(sys%outp%how(option__output__pol_density), em_resp_fd_dir, trim(fname), &
642 sys%namespace, sys%space, sys%gr, -sys%gr%x(:, jj) * lr_rho2(:, is), &
643 fn_unit, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
644 end do
645 end if
646 end do
647
648 end if
649 end if
650
651 !ELF
652 if (sys%outp%what(option__output__elf)) then
653
654 if (isign == 1) then
655 call elf_calc(sys%space, sys%st, sys%gr, sys%kpoints, elf, elfd)
656 else
657 call elf_calc(sys%space, sys%st, sys%gr, sys%kpoints, lr_elf, lr_elfd)
658
659 !numerical derivative
660 lr_elf(1:sys%gr%np, 1:sys%st%d%nspin) = &
661 ( lr_elf(1:sys%gr%np, 1:sys%st%d%nspin) - elf(1:sys%gr%np, 1:sys%st%d%nspin)) / (m_two * e_field)
662 lr_elfd(1:sys%gr%np, 1:sys%st%d%nspin) = &
663 (lr_elfd(1:sys%gr%np, 1:sys%st%d%nspin) - elfd(1:sys%gr%np, 1:sys%st%d%nspin)) / (m_two * e_field)
664
665 !write
666 do is = 1, sys%st%d%nspin
667 write(fname, '(a,i1,2a)') 'lr_elf-sp', is, '-', index2axis(ii)
668 call dio_function_output(sys%outp%how(option__output__elf), em_resp_fd_dir, trim(fname),&
669 sys%namespace, sys%space, sys%gr, lr_elf(:, is), unit_one, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
670 write(fname, '(a,i1,2a)') 'lr_elf_D-sp', is, '-', index2axis(ii)
671 call dio_function_output(sys%outp%how(option__output__elf), em_resp_fd_dir, trim(fname),&
672 sys%namespace, sys%space, sys%gr, lr_elfd(:, is), unit_one, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
673 end do
674 end if
675
676 end if
677
678 pop_sub(output_cycle_)
679 end subroutine output_cycle_
680
681
682 !-------------------------------------------------------------
683 subroutine output_end_()
684 real(real64) :: alpha(sys%space%dim, sys%space%dim)
685 complex(real64) :: beta(sys%space%dim, sys%space%dim, sys%space%dim)
686 integer :: iunit, idir
687 type(unit_t) :: fn_unit
688 real(real64) :: freq_factor(3)
689
690 push_sub(output_end_)
691
692 call io_mkdir(em_resp_fd_dir, sys%namespace)
693
694 if ((sys%outp%what(option__output__density) .or. &
695 sys%outp%what(option__output__pol_density)) .and. calc_diagonal) then
696 lr_rho2(1:sys%gr%np, 1:sys%st%d%nspin) = &
697 -(sys%st%rho(1:sys%gr%np, 1:sys%st%d%nspin) - lr_rho(1:sys%gr%np, 1:sys%st%d%nspin) &
698 - tmp_rho(1:sys%gr%np, 1:sys%st%d%nspin) + gs_rho(1:sys%gr%np, 1:sys%st%d%nspin)) / e_field**2
699
700 do is = 1, sys%st%d%nspin
701 if (sys%outp%what(option__output__density)) then
702 fn_unit = units_out%length**(2-sys%space%dim) / units_out%energy**2
703 write(fname, '(a,i1,a)') 'fd2_density-sp', is, '-y-z'
704 call dio_function_output(sys%outp%how(option__output__density), em_resp_fd_dir, trim(fname),&
705 sys%namespace, sys%space, sys%gr, lr_rho2(:, is), fn_unit, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
706 end if
707
708 if (sys%outp%what(option__output__pol_density)) then
709 fn_unit = units_out%length**(3-sys%space%dim) / units_out%energy**2
710 write(fname, '(a,i1,a)') 'beta_density-sp', is, '-x-y-z'
711 call dio_function_output(sys%outp%how(option__output__pol_density), em_resp_fd_dir, trim(fname),&
712 sys%namespace, sys%space, sys%gr, -sys%gr%x(:, 1) * lr_rho2(:, is), &
713 fn_unit, ierr, pos=sys%ions%pos, atoms=sys%ions%atom)
714 end if
715 end do
716 end if
717
718 if (mpi_grp_is_root(mpi_world)) then ! output pol file
719 iunit = io_open(em_resp_fd_dir//'alpha', sys%namespace, action='write')
720 write(iunit, '(3a)') '# Polarizability tensor [', trim(units_abbrev(units_out%polarizability)), ']'
721
722 alpha(1:sys%space%dim, 1:sys%space%dim) = (dipole(1:sys%space%dim, 1:sys%space%dim, 1) - &
723 dipole(1:sys%space%dim, 1:sys%space%dim, 2)) / (m_two * e_field)
724
725 beta = m_zero
726
727 do idir = 1, sys%space%dim
728 beta(1:sys%space%dim, idir, idir) = &
729 -(dipole(idir, 1:sys%space%dim, 1) + dipole(idir, 1:sys%space%dim, 2) - &
730 m_two * center_dipole(1:sys%space%dim)) / e_field**2
731 beta(idir, 1:sys%space%dim, idir) = beta(1:sys%space%dim, idir, idir)
732 beta(idir, idir, 1:sys%space%dim) = beta(1:sys%space%dim, idir, idir)
733 end do
734
735 if (calc_diagonal) then
736 beta(1, 2, 3) = -(diag_dipole(1) - dipole(2, 1, 1) - dipole(3, 1, 1) + center_dipole(1)) / e_field**2
737 else
738 beta(1, 2, 3) = m_zero
739 end if
740
741 beta(2, 3, 1) = beta(1, 2, 3)
742 beta(3, 1, 2) = beta(1, 2, 3)
743 beta(3, 2, 1) = beta(1, 2, 3)
744 beta(1, 3, 2) = beta(1, 2, 3)
745 beta(2, 1, 3) = beta(1, 2, 3)
746
747 call output_tensor(alpha, sys%space%dim, units_out%polarizability, iunit=iunit)
748 call io_close(iunit)
749
750 freq_factor(1:3) = m_zero ! for compatibility with em_resp version
751 call out_hyperpolarizability(sys%gr%box, beta, freq_factor(1:3), .true., em_resp_fd_dir, sys%namespace)
752
753 if (calc_born) then
754 call born_output_charges(born_charges, sys%ions%atom, sys%ions%charge, sys%ions%natoms, &
755 sys%namespace, sys%space%dim, em_resp_fd_dir, states_are_real(sys%st))
756 end if
757 end if
758
759 if (sys%outp%what(option__output__density) .or. &
760 sys%outp%what(option__output__pol_density)) then
761 safe_deallocate_a(lr_rho)
762 safe_deallocate_a(lr_rho2)
763 end if
764
765 if (sys%outp%what(option__output__elf)) then
766 safe_deallocate_a(lr_elf)
767 safe_deallocate_a(elf)
768 safe_deallocate_a(lr_elfd)
769 safe_deallocate_a(elfd)
770 end if
771
772 pop_sub(output_end_)
773 end subroutine output_end_
774
775 end subroutine static_pol_run_legacy
777end module static_pol_oct_m
778
779!! Local Variables:
780!! mode: f90
781!! coding: utf-8
782!! End:
subroutine init_(fromscratch)
Definition: geom_opt.F90:354
subroutine end_()
Definition: geom_opt.F90:804
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 elf_calc(space, st, gr, kpoints, elf, de)
(time-dependent) electron localization function, (TD)ELF.
Definition: elf.F90:167
subroutine, public out_hyperpolarizability(box, beta, freq_factor, converged, dirname, namespace)
Ref: David M Bishop, Rev Mod Phys 62, 343 (1990) beta generalized to lack of Kleinman symmetry.
Definition: em_resp.F90:3259
character(len= *), parameter, public em_resp_fd_dir
Definition: global.F90:253
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
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 lcao_run(namespace, space, gr, ions, ext_partners, st, ks, hm, st_start, lmm_r)
Definition: lcao.F90:799
This module defines various routines, operating on mesh functions.
real(real64) function, public dmf_moment(mesh, ff, idir, order)
This function calculates the "order" moment of the function ff.
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
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.
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
Definition: restart.F90:950
integer, parameter, public restart_gs
Definition: restart.F90:229
integer, parameter, public restart_em_resp_fd
Definition: restart.F90:229
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_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:906
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:859
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
subroutine, public scf_mix_clear(scf)
Definition: scf.F90:556
subroutine, public scf_run(scf, namespace, space, mc, gr, ions, ext_partners, st, ks, hm, outp, verbosity, iters_done, restart_load, restart_dump)
Definition: scf.F90:574
integer, parameter, public verb_full
Definition: scf.F90:197
integer, parameter, public verb_compact
Definition: scf.F90:197
subroutine, public scf_init(scf, namespace, gr, ions, st, mc, hm, space)
Definition: scf.F90:237
subroutine, public scf_end(scf)
Definition: scf.F90:526
pure logical function, public states_are_real(st)
subroutine, public states_elec_deallocate_wfns(st)
Deallocates the KS wavefunctions defined within a states_elec_t structure.
subroutine, public states_elec_allocate_wfns(st, mesh, wfs_type, skip, packed)
Allocates 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_dump(restart, space, st, mesh, kpoints, ierr, iter, lr, st_start_writing, verbose)
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 static_pol_run_legacy(sys, fromScratch)
Definition: static_pol.F90:175
subroutine, public static_pol_run(system, from_scratch)
Definition: static_pol.F90:157
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
subroutine, public output_dipole(dipole, ndim, iunit, namespace)
Definition: utils.F90:278
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
subroutine output_init_()
Definition: static_pol.F90:641
subroutine output_end_()
Definition: static_pol.F90:777
subroutine output_cycle_()
Definition: static_pol.F90:663
Class describing the electron system.
Definition: electrons.F90:214
Container class for lists of system_oct_m::system_t.
int true(void)