Octopus
propagation.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
23 use batch_oct_m
26 use debug_oct_m
32 use epot_oct_m
35 use forces_oct_m
37 use global_oct_m
38 use grid_oct_m
42 use ions_oct_m
43 use, intrinsic :: iso_fortran_env
46 use lasers_oct_m
48 use loct_oct_m
49 use mesh_oct_m
52 use mpi_oct_m
57 use parser_oct_m
64 use space_oct_m
68 use target_oct_m
70 use td_oct_m
72 use v_ks_oct_m
73
74 implicit none
75
76 private
77
78 public :: &
82 fwd_step, &
83 bwd_step, &
84 bwd_step_2, &
85 oct_prop_t, &
89
90
91 type oct_prop_t
92 private
93 integer :: number_checkpoints
94 integer, allocatable :: iter(:)
95 integer :: niter
96 character(len=100) :: dirname
97 type(restart_t) :: restart_load
98 type(restart_t) :: restart_dump
99 end type oct_prop_t
100
101
103 integer :: niter_
104 integer :: number_checkpoints_
105 real(real64) :: eta_
106 real(real64) :: delta_
107 logical :: zbr98_
108 logical :: gradients_
109
110contains
111
117 subroutine propagation_mod_init(niter, eta, delta, number_checkpoints, zbr98, gradients)
118 integer, intent(in) :: niter
119 real(real64), intent(in) :: eta
120 real(real64), intent(in) :: delta
121 integer, intent(in) :: number_checkpoints
122 logical, intent(in) :: zbr98
123 logical, intent(in) :: gradients
124
125 assert(.not. (zbr98 .and. gradients))
126
127 push_sub(propagation_mod_init)
128
129 niter_ = niter
130 eta_ = eta
131 delta_ = delta
132 number_checkpoints_ = number_checkpoints
133 zbr98_ = zbr98
134 gradients_ = gradients
135
136 pop_sub(propagation_mod_init)
137 end subroutine propagation_mod_init
138 ! ---------------------------------------------------------
139
140
146 subroutine propagate_forward(sys, td, par, tg, qcpsi, prop, write_iter)
147 type(electrons_t), intent(inout) :: sys
148 type(td_t), intent(inout) :: td
149 type(controlfunction_t), intent(in) :: par
150 type(target_t), intent(inout) :: tg
151 type(opt_control_state_t), intent(inout) :: qcpsi
152 type(oct_prop_t), optional, intent(inout) :: prop
153 logical, optional, intent(in) :: write_iter
154
155 integer :: ii, istep, ierr
156 logical :: write_iter_ = .false.
157 type(td_write_t) :: write_handler
158 real(real64), allocatable :: x_initial(:,:)
159 logical :: vel_target_ = .false.
160 type(states_elec_t), pointer :: psi
161
162 real(real64) :: init_time, final_time
163
164 push_sub(propagate_forward)
165
166 message(1) = "Info: Forward propagation."
167 call messages_info(1, namespace=sys%namespace)
168
169 call controlfunction_to_h(par, sys%ext_partners)
170
171 write_iter_ = .false.
172 if (present(write_iter)) write_iter_ = write_iter
173
174 psi => opt_control_point_qs(qcpsi)
175 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%pack()
176 call opt_control_get_classical(sys%ions, qcpsi)
177
178 if (write_iter_) then
179 call td_write_init(write_handler, sys%namespace, sys%space, sys%outp, sys%gr, sys%st, sys%hm, &
180 sys%ions, sys%ext_partners, sys%ks, ion_dynamics_ions_move(td%ions_dyn), &
181 list_has_gauge_field(sys%ext_partners), sys%hm%kick, td%iter, td%max_iter, &
182 td%dt, sys%mc)
183 call td_write_data(write_handler)
184 end if
185
188 ! setup the Hamiltonian
189 call density_calc(psi, sys%gr, psi%rho)
190 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners, time = m_zero)
191 call sys%hm%ks_pot%run_zero_iter(td%tr%vks_old)
192 if (ion_dynamics_ions_move(td%ions_dyn)) then
193 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, &
194 sys%ext_partners, psi, time = m_zero)
195 call forces_calculate(sys%gr, sys%namespace, sys%ions, sys%hm, sys%ext_partners, psi, sys%ks, t = m_zero, dt = td%dt)
196 end if
199 if (target_type(tg) == oct_tg_hhgnew) then
201 end if
202
203 if (target_type(tg) == oct_tg_velocity .or. target_type(tg) == oct_tg_hhgnew) then
204 safe_allocate_source_a(x_initial, sys%ions%pos)
205 vel_target_ = .true.
206 sys%ions%vel = m_zero
207 sys%ions%tot_force = m_zero
208 end if
209
210 if (.not. target_move_ions(tg)) call epot_precalc_local_potential(sys%hm%ep, sys%namespace, sys%gr, sys%ions)
211
212 call target_tdcalc(tg, sys%namespace, sys%space, sys%hm, sys%gr, sys%ions, sys%ext_partners, psi, 0, td%max_iter)
213
214 if (present(prop)) then
215 call oct_prop_dump_states(prop, sys%space, 0, psi, sys%gr, sys%kpoints, ierr)
216 if (ierr /= 0) then
217 message(1) = "Unable to write OCT states restart."
218 call messages_warning(1, namespace=sys%namespace)
219 end if
220 end if
221
222 init_time = loct_clock()
223 if (mpi_grp_is_root(mpi_world)) call loct_progress_bar(-1, td%max_iter)
224
225 ii = 1
226 do istep = 1, td%max_iter
227 ! time-iterate wavefunctions
228
229 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, psi, td%tr, istep*td%dt, td%dt, istep, &
230 td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
231
232 if (present(prop)) then
233 call oct_prop_dump_states(prop, sys%space, istep, psi, sys%gr, sys%kpoints, ierr)
234 if (ierr /= 0) then
235 message(1) = "Unable to write OCT states restart."
236 call messages_warning(1, namespace=sys%namespace)
237 end if
238 end if
240 ! update
241 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners, time = istep*td%dt)
242 call energy_calc_total(sys%namespace, sys%space, sys%hm, sys%gr, psi, sys%ext_partners)
243
244 if (sys%hm%abs_boundaries%abtype == mask_absorbing) call zvmask(sys%gr, sys%hm, psi)
245
246 ! if td_target
247 call target_tdcalc(tg, sys%namespace, sys%space, sys%hm, sys%gr, sys%ions, sys%ext_partners, psi, istep, td%max_iter)
248
249 ! only write in final run
250 if (write_iter_) then
251 call td_write_iter(write_handler, sys%namespace, sys%space, sys%outp, sys%gr, psi, sys%hm, sys%ions, &
252 sys%ext_partners, sys%hm%kick, sys%ks, td%dt, istep, sys%mc, sys%td%recalculate_gs)
253 ii = ii + 1
254 if (any(ii == sys%outp%output_interval + 1) .or. istep == td%max_iter) then ! output
255 if (istep == td%max_iter) sys%outp%output_interval = ii - 1
256 ii = istep
257 call td_write_data(write_handler)
258 end if
259 end if
260
261 if ((mod(istep, 100) == 0) .and. mpi_grp_is_root(mpi_world)) call loct_progress_bar(istep, td%max_iter)
262 end do
263 if (mpi_grp_is_root(mpi_world)) write(stdout, '(1x)')
264
265 final_time = loct_clock()
266 write(message(1),'(a,f12.2,a)') 'Propagation time: ', final_time - init_time, ' seconds.'
267 call messages_info(1, namespace=sys%namespace)
268
269 if (vel_target_) then
270 sys%ions%pos = x_initial
271 safe_deallocate_a(x_initial)
272 end if
273
274 call opt_control_set_classical(sys%ions, qcpsi)
275
276 if (write_iter_) call td_write_end(write_handler)
277 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%unpack()
278 nullify(psi)
279 pop_sub(propagate_forward)
280 end subroutine propagate_forward
281 ! ---------------------------------------------------------
282
283
288 subroutine propagate_backward(sys, td, qcpsi, prop)
289 type(electrons_t), intent(inout) :: sys
290 type(td_t), intent(inout) :: td
291 type(opt_control_state_t), intent(inout) :: qcpsi
292 type(oct_prop_t), intent(inout) :: prop
293
294 integer :: istep, ierr
295 type(states_elec_t), pointer :: psi
296
297 push_sub(propagate_backward)
298
299 message(1) = "Info: Backward propagation."
300 call messages_info(1, namespace=sys%namespace)
301
302 psi => opt_control_point_qs(qcpsi)
303 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%pack()
304
305 call hamiltonian_elec_adjoint(sys%hm)
306
307 ! setup the Hamiltonian
308 call density_calc(psi, sys%gr, psi%rho)
309 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
310 call sys%hm%ks_pot%run_zero_iter(td%tr%vks_old)
311
312 call oct_prop_dump_states(prop, sys%space, td%max_iter, psi, sys%gr, sys%kpoints, ierr)
313 if (ierr /= 0) then
314 message(1) = "Unable to write OCT states restart."
315 call messages_warning(1, namespace=sys%namespace)
316 end if
317
318 if (mpi_grp_is_root(mpi_world)) call loct_progress_bar(-1, td%max_iter)
319
320 do istep = td%max_iter, 1, -1
321 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, psi, td%tr, &
322 (istep - 1)*td%dt, -td%dt, istep-1, td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
323
324 call oct_prop_dump_states(prop, sys%space, istep - 1, psi, sys%gr, sys%kpoints, ierr)
325 if (ierr /= 0) then
326 message(1) = "Unable to write OCT states restart."
327 call messages_warning(1, namespace=sys%namespace)
328 end if
329
330 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
331 if (mod(istep, 100) == 0 .and. mpi_grp_is_root(mpi_world)) call loct_progress_bar(td%max_iter - istep + 1, td%max_iter)
332 end do
333
334 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%unpack()
335 nullify(psi)
336 pop_sub(propagate_backward)
337 end subroutine propagate_backward
338 ! ---------------------------------------------------------
339
340
354 subroutine fwd_step(sys, td, tg, par, par_chi, qcpsi, prop_chi, prop_psi)
355 type(electrons_t), intent(inout) :: sys
356 type(td_t), intent(inout) :: td
357 type(target_t), intent(inout) :: tg
358 type(controlfunction_t), intent(inout) :: par
359 type(controlfunction_t), intent(in) :: par_chi
360 type(opt_control_state_t), intent(inout) :: qcpsi
361 type(oct_prop_t), intent(inout) :: prop_chi
362 type(oct_prop_t), intent(inout) :: prop_psi
363
364 integer :: i, ierr
365 logical :: aux_fwd_propagation
366 type(states_elec_t) :: psi2
367 type(opt_control_state_t) :: qcchi
368 type(controlfunction_t) :: par_prev
369 type(propagator_base_t) :: tr_chi, tr_psi2
370 type(states_elec_t), pointer :: psi, chi
371
372 push_sub(fwd_step)
373
374 message(1) = "Info: Forward propagation."
375 call messages_info(1, namespace=sys%namespace)
376
378
379 call opt_control_state_null(qcchi)
380 call opt_control_state_copy(qcchi, qcpsi)
382 psi => opt_control_point_qs(qcpsi)
383 chi => opt_control_point_qs(qcchi)
384 call propagator_elec_copy(tr_chi, td%tr)
385 ! The propagation of chi should not be self-consistent, because the Kohn-Sham
386 ! potential used is the one created by psi. Note, however, that it is likely that
387 ! the first two iterations are done self-consistently nonetheless.
389
390 aux_fwd_propagation = (target_mode(tg) == oct_targetmode_td .or. &
391 (sys%hm%theory_level /= independent_particles .and. &
392 .not. sys%ks%frozen_hxc))
393 if (aux_fwd_propagation) then
394 call states_elec_copy(psi2, psi)
395 call controlfunction_copy(par_prev, par)
396 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi2%pack()
397 end if
398
399 ! setup forward propagation
400 call density_calc(psi, sys%gr, psi%rho)
401 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
402 call sys%hm%ks_pot%run_zero_iter(td%tr%vks_old)
403 call sys%hm%ks_pot%run_zero_iter(tr_chi%vks_old)
404 if (aux_fwd_propagation) then
405 call propagator_elec_copy(tr_psi2, td%tr)
406 call sys%hm%ks_pot%run_zero_iter(tr_psi2%vks_old)
407 end if
408
409 call oct_prop_dump_states(prop_psi, sys%space, 0, psi, sys%gr, sys%kpoints, ierr)
410 if (ierr /= 0) then
411 message(1) = "Unable to write OCT states restart."
412 call messages_warning(1, namespace=sys%namespace)
413 end if
414 call oct_prop_load_states(prop_chi, sys%namespace, sys%space, chi, sys%gr, sys%kpoints, 0, ierr)
415 if (ierr /= 0) then
416 message(1) = "Unable to read OCT states restart."
417 call messages_fatal(1, namespace=sys%namespace)
418 end if
419
420 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%pack()
421 if (sys%st%pack_states .and. sys%hm%apply_packed()) call chi%pack()
422
423 do i = 1, td%max_iter
424 call update_field(i, par, sys%space, sys%gr, sys%hm, sys%ext_partners, sys%ions, qcpsi, qcchi, par_chi, dir = 'f')
425 call update_hamiltonian_elec_chi(i, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
426 td, tg, par_chi, sys%ions, psi2)
427 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners, time = (i - 1)*td%dt)
428 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, chi, tr_chi, i*td%dt, td%dt, i, &
429 td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
430 if (aux_fwd_propagation) then
431 call update_hamiltonian_elec_psi(i, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
432 td, tg, par_prev, psi2, sys%ions)
433 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, psi2, tr_psi2, i*td%dt, td%dt, i, &
434 td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
435 end if
436 call update_hamiltonian_elec_psi(i, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, &
437 sys%ext_partners, td, tg, par, psi, sys%ions)
438 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners, time = (i - 1)*td%dt)
439 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, psi, td%tr, i*td%dt, td%dt, i, &
440 td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
441 call target_tdcalc(tg, sys%namespace, sys%space, sys%hm, sys%gr, sys%ions, sys%ext_partners, psi, i, td%max_iter)
442
443 call oct_prop_dump_states(prop_psi, sys%space, i, psi, sys%gr, sys%kpoints, ierr)
444 if (ierr /= 0) then
445 message(1) = "Unable to write OCT states restart."
446 call messages_warning(1, namespace=sys%namespace)
447 end if
448 call oct_prop_check(prop_chi, sys%namespace, sys%space, chi, sys%gr, sys%kpoints, i)
449 end do
450 call update_field(td%max_iter+1, par, sys%space, sys%gr, sys%hm, sys%ext_partners, sys%ions, qcpsi, qcchi, par_chi, dir = 'f')
451
452 call density_calc(psi, sys%gr, psi%rho)
453 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
454
455 if (target_mode(tg) == oct_targetmode_td .or. &
456 (sys%hm%theory_level /= independent_particles .and. (.not. sys%ks%frozen_hxc))) then
457 call states_elec_end(psi2)
458 call controlfunction_end(par_prev)
459 end if
460
462 if (aux_fwd_propagation) call propagator_elec_end(tr_psi2)
463 call states_elec_end(chi)
464 call propagator_elec_end(tr_chi)
465 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%unpack()
466 nullify(psi)
467 nullify(chi)
468 pop_sub(fwd_step)
469 end subroutine fwd_step
470 ! ---------------------------------------------------------
471
472
482 subroutine bwd_step(sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
483 type(electrons_t), intent(inout) :: sys
484 type(td_t), intent(inout) :: td
485 type(target_t), intent(inout) :: tg
486 type(controlfunction_t), intent(in) :: par
487 type(controlfunction_t), intent(inout) :: par_chi
488 type(opt_control_state_t), intent(inout) :: qcchi
489 type(oct_prop_t), intent(inout) :: prop_chi
490 type(oct_prop_t), intent(inout) :: prop_psi
491
492 integer :: i, ierr
493 type(propagator_base_t) :: tr_chi
494 type(opt_control_state_t) :: qcpsi
495 type(states_elec_t), pointer :: chi, psi
496
497 push_sub(bwd_step)
498
499 message(1) = "Info: Backward propagation."
500 call messages_info(1, namespace=sys%namespace)
501
502 call controlfunction_to_realtime(par_chi)
503
504 chi => opt_control_point_qs(qcchi)
505 psi => opt_control_point_qs(qcpsi)
506
507 call propagator_elec_copy(tr_chi, td%tr)
508 ! The propagation of chi should not be self-consistent, because the Kohn-Sham
509 ! potential used is the one created by psi. Note, however, that it is likely that
510 ! the first two iterations are done self-consistently nonetheless.
512
513 call states_elec_copy(psi, chi)
514 call oct_prop_load_states(prop_psi, sys%namespace, sys%space, psi, sys%gr, sys%kpoints, td%max_iter, ierr)
515 if (ierr /= 0) then
516 message(1) = "Unable to read OCT states restart."
517 call messages_fatal(1, namespace=sys%namespace)
518 end if
519 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%pack()
520 if (sys%st%pack_states .and. sys%hm%apply_packed()) call chi%pack()
521
522 call density_calc(psi, sys%gr, psi%rho)
523 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
524 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners)
525 call sys%hm%ks_pot%run_zero_iter(td%tr%vks_old)
526 call sys%hm%ks_pot%run_zero_iter(tr_chi%vks_old)
527
528 td%dt = -td%dt
529 call oct_prop_dump_states(prop_chi, sys%space, td%max_iter, chi, sys%gr, sys%kpoints, ierr)
530 if (ierr /= 0) then
531 message(1) = "Unable to write OCT states restart."
532 call messages_warning(1, namespace=sys%namespace)
533 end if
534
535 do i = td%max_iter, 1, -1
536 call oct_prop_check(prop_psi, sys%namespace, sys%space, psi, sys%gr, sys%kpoints, i)
537 call update_field(i, par_chi, sys%space, sys%gr, sys%hm, sys%ext_partners, sys%ions, qcpsi, qcchi, par, dir = 'b')
538 call update_hamiltonian_elec_chi(i-1, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
539 td, tg, par_chi, sys%ions, psi)
540 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners, time = abs(i*td%dt))
541 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, chi, tr_chi, abs((i-1)*td%dt), td%dt, &
542 i-1, td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
543 call oct_prop_dump_states(prop_chi, sys%space, i-1, chi, sys%gr, sys%kpoints, ierr)
544 if (ierr /= 0) then
545 message(1) = "Unable to write OCT states restart."
546 call messages_warning(1, namespace=sys%namespace)
547 end if
548 call update_hamiltonian_elec_psi(i-1, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
549 td, tg, par, psi, sys%ions)
550 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners, time = abs(i*td%dt))
551 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, psi, td%tr, abs((i-1)*td%dt), td%dt, &
552 i-1, td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
553 end do
554 td%dt = -td%dt
555 call update_field(0, par_chi, sys%space, sys%gr, sys%hm, sys%ext_partners, sys%ions, qcpsi, qcchi, par, dir = 'b')
556
557 call density_calc(psi, sys%gr, psi%rho)
558 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
559 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners)
560
561 call controlfunction_to_basis(par_chi)
562 call states_elec_end(psi)
563 call propagator_elec_end(tr_chi)
564 if (sys%st%pack_states .and. sys%hm%apply_packed()) call chi%unpack()
565 nullify(chi)
566 nullify(psi)
567 pop_sub(bwd_step)
568 end subroutine bwd_step
569 ! ---------------------------------------------------------
570
571
584 subroutine bwd_step_2(sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
585 type(electrons_t), intent(inout) :: sys
586 type(td_t), intent(inout) :: td
587 type(target_t), intent(inout) :: tg
588 type(controlfunction_t), intent(in) :: par
589 type(controlfunction_t), intent(inout) :: par_chi
590 type(opt_control_state_t), intent(inout) :: qcchi
591 type(oct_prop_t), intent(inout) :: prop_chi
592 type(oct_prop_t), intent(inout) :: prop_psi
593
594 integer :: i, ierr, ik, ib
595 logical :: freeze
596 type(propagator_base_t) :: tr_chi
597 type(opt_control_state_t) :: qcpsi
598 type(states_elec_t) :: st_ref
599 type(states_elec_t), pointer :: chi, psi
600 real(real64), pointer :: q(:, :), p(:, :)
601 real(real64), allocatable :: qtildehalf(:, :), qinitial(:, :)
602 real(real64), allocatable :: vhxc(:, :)
603 real(real64), allocatable :: fold(:, :), fnew(:, :)
604 type(ion_state_t) :: ions_state_initial, ions_state_final
605
606 real(real64) :: init_time, final_time
607
608 push_sub(bwd_step_2)
609
610 chi => opt_control_point_qs(qcchi)
611 q => opt_control_point_q(qcchi)
612 p => opt_control_point_p(qcchi)
613 safe_allocate(qtildehalf(1:sys%ions%natoms, 1:sys%ions%space%dim))
614 safe_allocate(qinitial(1:sys%ions%space%dim, 1:sys%ions%natoms))
615
616 call propagator_elec_copy(tr_chi, td%tr)
617 ! The propagation of chi should not be self-consistent, because the Kohn-Sham
618 ! potential used is the one created by psi. Note, however, that it is likely that
619 ! the first two iterations are done self-consistently nonetheless.
621
622 call opt_control_state_null(qcpsi)
623 call opt_control_state_copy(qcpsi, qcchi)
624 psi => opt_control_point_qs(qcpsi)
625 call oct_prop_load_states(prop_psi, sys%namespace, sys%space, psi, sys%gr, sys%kpoints, td%max_iter, ierr)
626 if (ierr /= 0) then
627 message(1) = "Unable to read OCT states restart."
628 call messages_fatal(1, namespace=sys%namespace)
629 end if
630
631 safe_allocate(vhxc(1:sys%gr%np, 1:sys%hm%d%nspin))
632
633 call density_calc(psi, sys%gr, psi%rho)
634 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
635 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners)
636 call sys%hm%ks_pot%run_zero_iter(td%tr%vks_old)
637 call sys%hm%ks_pot%run_zero_iter(tr_chi%vks_old)
638 td%dt = -td%dt
639 call oct_prop_dump_states(prop_chi, sys%space, td%max_iter, chi, sys%gr, sys%kpoints, ierr)
640 if (ierr /= 0) then
641 message(1) = "Unable to write OCT states restart."
642 call messages_warning(1, namespace=sys%namespace)
643 end if
644
645 call states_elec_copy(st_ref, psi)
646 if (sys%st%pack_states .and. sys%hm%apply_packed()) call psi%pack()
647 if (sys%st%pack_states .and. sys%hm%apply_packed()) call st_ref%pack()
648
649 if (ion_dynamics_ions_move(td%ions_dyn)) then
650 call forces_calculate(sys%gr, sys%namespace, sys%ions, sys%hm, sys%ext_partners, &
651 psi, sys%ks, t = td%max_iter*abs(td%dt), dt = td%dt)
652 end if
653
654 message(1) = "Info: Backward propagation."
655 call messages_info(1, namespace=sys%namespace)
656 if (mpi_grp_is_root(mpi_world)) call loct_progress_bar(-1, td%max_iter)
657
658 init_time = loct_clock()
659
660 do i = td%max_iter, 1, -1
661
662 call oct_prop_check(prop_psi, sys%namespace, sys%space, psi, sys%gr, sys%kpoints, i)
663 call update_field(i, par_chi, sys%space, sys%gr, sys%hm, sys%ext_partners, sys%ions, qcpsi, qcchi, par, dir = 'b')
664
665 select case (td%tr%method)
666
668
669 call update_hamiltonian_elec_psi(i-1, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
670 td, tg, par, psi, sys%ions)
671 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, psi, td%tr, abs((i-1)*td%dt), td%dt, &
672 i-1, td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler, qcchi = qcchi)
673
674 case default
675
676 if (ion_dynamics_ions_move(td%ions_dyn)) then
677 qtildehalf = q
678
679 call ion_dynamics_save_state(td%ions_dyn, sys%ions, ions_state_initial)
680 call ion_dynamics_propagate(td%ions_dyn, sys%ions, abs((i-1)*td%dt), m_half * td%dt, sys%namespace)
681 qinitial = sys%ions%pos
682 call ion_dynamics_restore_state(td%ions_dyn, sys%ions, ions_state_initial)
683
684 safe_allocate(fold(1:sys%ions%natoms, 1:sys%ions%space%dim))
685 safe_allocate(fnew(1:sys%ions%natoms, 1:sys%ions%space%dim))
686 call forces_costate_calculate(sys%gr, sys%namespace, sys%ions, sys%hm, psi, chi, fold, q)
687
688 call ion_dynamics_verlet_step1(sys%ions, qtildehalf, p, fold, m_half * td%dt)
689 call ion_dynamics_verlet_step1(sys%ions, q, p, fold, td%dt)
690 end if
691
692 ! Here propagate psi one full step, and then simply interpolate to get the state
693 ! at half the time interval. Perhaps one could gain some accuracy by performing two
694 ! successive propagations of half time step.
695 call update_hamiltonian_elec_psi(i-1, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
696 td, tg, par, psi, sys%ions)
697
698 do ik = psi%d%kpt%start, psi%d%kpt%end
699 do ib = psi%group%block_start, psi%group%block_end
700 call psi%group%psib(ib, ik)%copy_data_to(sys%gr%np, st_ref%group%psib(ib, ik))
701 end do
702 end do
703
704 vhxc(:, :) = sys%hm%ks_pot%vhxc(:, :)
705 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, psi, td%tr, abs((i-1)*td%dt), td%dt, &
706 i-1, td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
707
708 if (ion_dynamics_ions_move(td%ions_dyn)) then
709 call ion_dynamics_save_state(td%ions_dyn, sys%ions, ions_state_final)
710 sys%ions%pos = qinitial
711 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, &
712 sys%ext_partners, psi, time = abs((i-1)*td%dt))
713 end if
714
715 do ik = psi%d%kpt%start, psi%d%kpt%end
716 do ib = psi%group%block_start, psi%group%block_end
717 call batch_scal(sys%gr%np, cmplx(m_half, m_zero, real64) , &
718 st_ref%group%psib(ib, ik))
719 call batch_axpy(sys%gr%np, cmplx(m_half, m_zero, real64) , &
720 psi%group%psib(ib, ik), st_ref%group%psib(ib, ik))
721 end do
722 end do
723
724 sys%hm%ks_pot%vhxc(:, :) = m_half * (sys%hm%ks_pot%vhxc(:, :) + vhxc(:, :))
725 call update_hamiltonian_elec_chi(i-1, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
726 td, tg, par, sys%ions, st_ref, qtildehalf)
727 freeze = ion_dynamics_freeze(td%ions_dyn)
728 call propagator_elec_dt(sys%ks, sys%namespace, sys%space, sys%hm, sys%gr, chi, tr_chi, abs((i-1)*td%dt), td%dt, &
729 i-1, td%ions_dyn, sys%ions, sys%ext_partners, sys%outp, td%write_handler)
730 if (freeze) call ion_dynamics_unfreeze(td%ions_dyn)
731
732 if (ion_dynamics_ions_move(td%ions_dyn)) then
733 call ion_dynamics_restore_state(td%ions_dyn, sys%ions, ions_state_final)
734 call hamiltonian_elec_epot_generate(sys%hm, sys%namespace, sys%space, sys%gr, sys%ions, &
735 sys%ext_partners, psi, time = abs((i-1)*td%dt))
736 call forces_calculate(sys%gr, sys%namespace, sys%ions, sys%hm, sys%ext_partners, &
737 psi, sys%ks, t = abs((i-1)*td%dt), dt = td%dt)
738 call forces_costate_calculate(sys%gr, sys%namespace, sys%ions, sys%hm, psi, chi, fnew, q)
739 call ion_dynamics_verlet_step2(sys%ions, p, fold, fnew, td%dt)
740 safe_deallocate_a(fold)
741 safe_deallocate_a(fnew)
742 end if
743
744 sys%hm%ks_pot%vhxc(:, :) = vhxc(:, :)
745
746 end select
747
748 call oct_prop_dump_states(prop_chi, sys%space, i-1, chi, sys%gr, sys%kpoints, ierr)
749 if (ierr /= 0) then
750 message(1) = "Unable to write OCT states restart."
751 call messages_warning(1, namespace=sys%namespace)
752 end if
753
754 if ((mod(i, 100) == 0).and. mpi_grp_is_root(mpi_world)) call loct_progress_bar(td%max_iter-i, td%max_iter)
755 end do
756 if (mpi_grp_is_root(mpi_world)) then
757 call loct_progress_bar(td%max_iter, td%max_iter)
758 write(stdout, '(1x)')
759 end if
760
761 final_time = loct_clock()
762 write(message(1),'(a,f12.2,a)') 'Propagation time: ', final_time - init_time, ' seconds.'
763 call messages_info(1, namespace=sys%namespace)
764
765 call states_elec_end(st_ref)
766
767 td%dt = -td%dt
768 call update_hamiltonian_elec_psi(0, sys%namespace, sys%space, sys%gr, sys%ks, sys%hm, sys%ext_partners, &
769 td, tg, par, psi, sys%ions)
770 call update_field(0, par_chi, sys%space, sys%gr, sys%hm, sys%ext_partners, sys%ions, qcpsi, qcchi, par, dir = 'b')
771
772 call density_calc(psi, sys%gr, psi%rho)
773 call v_ks_calc(sys%ks, sys%namespace, sys%space, sys%hm, psi, sys%ions, sys%ext_partners)
774 call sys%hm%update(sys%gr, sys%namespace, sys%space, sys%ext_partners)
775
776 call propagator_elec_end(tr_chi)
777
778 safe_deallocate_a(vhxc)
779 call states_elec_end(psi)
780
781 nullify(chi)
782 nullify(psi)
783 nullify(q)
784 nullify(p)
785 safe_deallocate_a(qtildehalf)
786 safe_deallocate_a(qinitial)
787 pop_sub(bwd_step_2)
788 end subroutine bwd_step_2
789 ! ----------------------------------------------------------
790
791
792 ! ----------------------------------------------------------
793 !
794 ! ----------------------------------------------------------
795 subroutine update_hamiltonian_elec_chi(iter, namespace, space, gr, ks, hm, ext_partners, td, tg, par_chi, ions, st, qtildehalf)
796 integer, intent(in) :: iter
797 type(namespace_t), intent(in) :: namespace
798 class(space_t), intent(in) :: space
799 type(grid_t), intent(inout) :: gr
800 type(v_ks_t), intent(inout) :: ks
801 type(hamiltonian_elec_t), intent(inout) :: hm
802 type(partner_list_t), intent(in) :: ext_partners
803 type(td_t), intent(inout) :: td
804 type(target_t), intent(inout) :: tg
805 type(controlfunction_t), intent(in) :: par_chi
806 type(ions_t), intent(in) :: ions
807 type(states_elec_t), intent(inout) :: st
808 real(real64), optional, intent(in) :: qtildehalf(:, :)
809
810 type(states_elec_t) :: inh
811 type(perturbation_ionic_t), pointer :: pert
812 integer :: j, iatom, idim
813 complex(real64), allocatable :: dvpsi(:, :, :), zpsi(:, :), inhzpsi(:, :)
814 integer :: ist, ik, ib
815
817
818 assert(.not. st%parallel_in_states)
819 assert(.not. st%d%kpt%parallel)
820
821 if (target_mode(tg) == oct_targetmode_td) then
822 call states_elec_copy(inh, st)
823 call target_inh(st, gr, hm%kpoints, tg, abs(td%dt)*iter, inh, iter)
824 call hamiltonian_elec_set_inh(hm, inh)
825 call states_elec_end(inh)
826 end if
827
828 if (ion_dynamics_ions_move(td%ions_dyn)) then
829 call states_elec_copy(inh, st)
830 safe_allocate(dvpsi(1:gr%np_part, 1:st%d%dim, 1:space%dim))
831 do ik = inh%d%kpt%start, inh%d%kpt%end
832 do ib = inh%group%block_start, inh%group%block_end
833 call batch_set_zero(inh%group%psib(ib, ik))
834 end do
835 end do
836 safe_allocate(zpsi(1:gr%np_part, 1:st%d%dim))
837 safe_allocate(inhzpsi(1:gr%np_part, 1:st%d%dim))
838 do ist = 1, st%nst
839 do ik = 1, st%nik
840
841 call states_elec_get_state(st, gr, ist, ik, zpsi)
842 call states_elec_get_state(inh, gr, ist, ik, inhzpsi)
843
844 pert => perturbation_ionic_t(namespace, ions)
845 do iatom = 1, ions%natoms
846 call pert%setup_atom(iatom)
847 do idim = 1, space%dim
848 call pert%setup_dir(idim)
849 call pert%zapply(namespace, space, gr, hm, ik, zpsi(:, :), dvpsi(:, :, idim))
850
851 call lalg_axpy(gr%np, st%d%dim, -st%occ(ist, ik)*qtildehalf(iatom, idim), &
852 dvpsi(:, :, idim), inhzpsi(:, :))
853 end do
854 end do
855 safe_deallocate_p(pert)
856 call states_elec_set_state(inh, gr, ist, ik, inhzpsi)
857 end do
858 end do
859
860 safe_deallocate_a(zpsi)
861 safe_deallocate_a(inhzpsi)
862 safe_deallocate_a(dvpsi)
863 call hamiltonian_elec_set_inh(hm, inh)
864 call states_elec_end(inh)
865 end if
866
867 if (hm%theory_level /= independent_particles .and. (.not. ks%frozen_hxc)) then
868 call density_calc(st, gr, st%rho)
869 call oct_exchange_set(hm%oct_exchange, st, gr)
870 end if
871
873
874 do j = iter - 2, iter + 2
875 if (j >= 0 .and. j <= td%max_iter) then
876 call controlfunction_to_h_val(par_chi, ext_partners, j+1)
877 end if
878 end do
879
881 end subroutine update_hamiltonian_elec_chi
882 ! ---------------------------------------------------------
883
884
885 ! ----------------------------------------------------------
886 !
887 ! ----------------------------------------------------------
888 subroutine update_hamiltonian_elec_psi(iter, namespace, space, gr, ks, hm, ext_partners, td, tg, par, st, ions)
889 integer, intent(in) :: iter
890 type(namespace_t), intent(in) :: namespace
891 type(electron_space_t), intent(in) :: space
892 type(grid_t), intent(inout) :: gr
893 type(v_ks_t), intent(inout) :: ks
894 type(hamiltonian_elec_t), intent(inout) :: hm
895 type(partner_list_t), intent(in) :: ext_partners
896 type(td_t), intent(inout) :: td
897 type(target_t), intent(inout) :: tg
898 type(controlfunction_t), intent(in) :: par
899 type(states_elec_t), intent(inout) :: st
900 type(ions_t), intent(in) :: ions
901
902 integer :: j
903
905
906 if (target_mode(tg) == oct_targetmode_td) then
908 end if
909
910 if (ion_dynamics_ions_move(td%ions_dyn)) then
912 end if
913
914 if (hm%theory_level /= independent_particles .and. (.not. ks%frozen_hxc)) then
915 call oct_exchange_remove(hm%oct_exchange)
916 end if
917
919
920 do j = iter - 2, iter + 2
921 if (j >= 0 .and. j <= td%max_iter) then
922 call controlfunction_to_h_val(par, ext_partners, j+1)
923 end if
924 end do
925 if (hm%theory_level /= independent_particles .and. (.not. ks%frozen_hxc)) then
926 call density_calc(st, gr, st%rho)
927 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners)
928 call hm%update(gr, namespace, space, ext_partners)
929 end if
930
932 end subroutine update_hamiltonian_elec_psi
933 ! ---------------------------------------------------------
934
935
936 ! ---------------------------------------------------------
937 subroutine calculate_g(space, gr, hm, lasers, psi, chi, dl, dq)
938 class(space_t), intent(in) :: space
939 type(grid_t), intent(inout) :: gr
940 type(hamiltonian_elec_t), intent(in) :: hm
941 type(lasers_t), intent(in) :: lasers
942 type(states_elec_t), intent(inout) :: psi
943 type(states_elec_t), intent(inout) :: chi
944 complex(real64), intent(inout) :: dl(:), dq(:)
945
946 complex(real64), allocatable :: zpsi(:, :), zoppsi(:, :)
947 integer :: no_parameters, j, ik, p
948
949 push_sub(calculate_g)
950
951 no_parameters = lasers%no_lasers
952
953 safe_allocate(zpsi(1:gr%np_part, 1:chi%d%dim))
954 safe_allocate(zoppsi(1:gr%np_part, 1:chi%d%dim))
955
956 do j = 1, no_parameters
957
958 dl(j) = m_z0
959 do ik = 1, psi%nik
960 do p = 1, psi%nst
961
962 call states_elec_get_state(psi, gr, p, ik, zpsi)
963
964 zoppsi = m_z0
965 if (allocated(hm%ep%a_static)) then
966 call vlaser_operator_linear(lasers%lasers(j), gr%der, hm%d, zpsi, &
967 zoppsi, ik, hm%ep%gyromagnetic_ratio, hm%ep%a_static)
968 else
969 call vlaser_operator_linear(lasers%lasers(j), gr%der, hm%d, zpsi, &
970 zoppsi, ik, hm%ep%gyromagnetic_ratio)
971 end if
972
973 call states_elec_get_state(chi, gr, p, ik, zpsi)
974 dl(j) = dl(j) + zmf_dotp(gr, psi%d%dim, zpsi, zoppsi)
975 end do
976 end do
977
978 ! The quadratic part should only be computed if necessary.
979 if (laser_kind(lasers%lasers(j)) == e_field_magnetic) then
980
981 dq(j) = m_z0
982 do ik = 1, psi%nik
983 do p = 1, psi%nst
984 zoppsi = m_z0
985
986 call states_elec_get_state(psi, gr, p, ik, zpsi)
987 call vlaser_operator_quadratic(lasers%lasers(j), gr, space, zpsi, zoppsi)
988
989 call states_elec_get_state(chi, gr, p, ik, zpsi)
990 dq(j) = dq(j) + zmf_dotp(gr, psi%d%dim, zpsi, zoppsi)
991
992 end do
993 end do
994
995 else
996 dq(j) = m_z0
997 end if
998 end do
999
1000 safe_deallocate_a(zpsi)
1001 safe_deallocate_a(zoppsi)
1002
1003 pop_sub(calculate_g)
1004 end subroutine calculate_g
1005 ! ---------------------------------------------------------
1006
1007
1008
1009
1022 subroutine update_field(iter, cp, space, gr, hm, ext_partners, ions, qcpsi, qcchi, cpp, dir)
1023 class(space_t), intent(in) :: space
1024 integer, intent(in) :: iter
1025 type(controlfunction_t), intent(inout) :: cp
1026 type(grid_t), intent(inout) :: gr
1027 type(hamiltonian_elec_t), intent(in) :: hm
1028 type(partner_list_t), intent(in) :: ext_partners
1029 type(ions_t), intent(in) :: ions
1030 type(opt_control_state_t), intent(inout) :: qcpsi
1031 type(opt_control_state_t), intent(inout) :: qcchi
1032 type(controlfunction_t), intent(in) :: cpp
1033 character(len=1), intent(in) :: dir
1034
1035 complex(real64) :: d1, pol(3)
1036 complex(real64), allocatable :: dl(:), dq(:), zpsi(:, :), zchi(:, :)
1037 real(real64), allocatable :: d(:)
1038 integer :: j, no_parameters, iatom
1039 type(states_elec_t), pointer :: psi, chi
1040 real(real64), pointer :: q(:, :)
1041 type(lasers_t), pointer :: lasers
1042
1043 push_sub(update_field)
1044
1045 psi => opt_control_point_qs(qcpsi)
1046 chi => opt_control_point_qs(qcchi)
1047 q => opt_control_point_q(qcchi)
1048
1049 no_parameters = controlfunction_number(cp)
1050
1051 safe_allocate(dl(1:no_parameters))
1052 safe_allocate(dq(1:no_parameters))
1053 safe_allocate( d(1:no_parameters))
1054
1055
1056 lasers => list_get_lasers(ext_partners)
1057 if(associated(lasers)) then
1058 call calculate_g(space, gr, hm, lasers, psi, chi, dl, dq)
1059 end if
1060
1061 d1 = m_z1
1062 if (zbr98_) then
1063 safe_allocate(zpsi(1:gr%np, 1:psi%d%dim))
1064 safe_allocate(zchi(1:gr%np, 1:chi%d%dim))
1065
1066 call states_elec_get_state(psi, gr, 1, 1, zpsi)
1067 call states_elec_get_state(chi, gr, 1, 1, zchi)
1068
1069 d1 = zmf_dotp(gr, psi%d%dim, zpsi, zchi)
1070 do j = 1, no_parameters
1071 d(j) = aimag(d1*dl(j)) / controlfunction_alpha(cp, j)
1072 end do
1073
1074 safe_deallocate_a(zpsi)
1075 safe_deallocate_a(zchi)
1076
1077 elseif (gradients_) then
1078 do j = 1, no_parameters
1079 d(j) = m_two * aimag(dl(j))
1080 end do
1081 else
1082 do j = 1, no_parameters
1083 d(j) = aimag(dl(j)) / controlfunction_alpha(cp, j)
1084 end do
1085 end if
1086
1087 ! This is for the classical target.
1088 if (dir == 'b' .and. associated(lasers)) then
1089 pol = laser_polarization(lasers%lasers(1))
1090 do iatom = 1, ions%natoms
1091 d(1) = d(1) - ions%charge(iatom) * real(sum(pol(1:ions%space%dim)*q(iatom, 1:ions%space%dim)), real64)
1092 end do
1093 end if
1094
1095
1096 if (dir == 'f') then
1097 call controlfunction_update(cp, cpp, dir, iter, delta_, d, dq)
1098 else
1099 call controlfunction_update(cp, cpp, dir, iter, eta_, d, dq)
1100 end if
1101
1102 nullify(q)
1103 nullify(psi)
1104 nullify(chi)
1105 safe_deallocate_a(d)
1106 safe_deallocate_a(dl)
1107 safe_deallocate_a(dq)
1108 pop_sub(update_field)
1109 end subroutine update_field
1110 ! ---------------------------------------------------------
1111
1112
1113 ! ---------------------------------------------------------
1114 subroutine oct_prop_init(prop, namespace, dirname, mesh, mc)
1115 type(oct_prop_t), intent(inout) :: prop
1116 type(namespace_t), intent(in) :: namespace
1117 character(len=*), intent(in) :: dirname
1118 class(mesh_t), intent(in) :: mesh
1119 type(multicomm_t), intent(in) :: mc
1120
1121 integer :: j, ierr
1122
1123 push_sub(oct_prop_init)
1124
1125 prop%dirname = dirname
1126 prop%niter = niter_
1127 prop%number_checkpoints = number_checkpoints_
1128
1129 ! The OCT_DIR//trim(dirname) will be used to write and read information during the calculation,
1130 ! so they need to use the same path.
1131 call restart_init(prop%restart_dump, namespace, restart_oct, restart_type_dump, mc, ierr, mesh=mesh)
1132 call restart_init(prop%restart_load, namespace, restart_oct, restart_type_load, mc, ierr, mesh=mesh)
1133
1134 safe_allocate(prop%iter(1:prop%number_checkpoints+2))
1135 prop%iter(1) = 0
1136 do j = 1, prop%number_checkpoints
1137 prop%iter(j+1) = nint(real(niter_, real64) /(prop%number_checkpoints+1) * j)
1138 end do
1139 prop%iter(prop%number_checkpoints+2) = niter_
1140
1141 pop_sub(oct_prop_init)
1142 end subroutine oct_prop_init
1143 ! ---------------------------------------------------------
1144
1145
1146 ! ---------------------------------------------------------
1147 subroutine oct_prop_end(prop)
1148 type(oct_prop_t), intent(inout) :: prop
1149
1150 push_sub(oct_prop_end)
1151
1152 call restart_end(prop%restart_load)
1153 call restart_end(prop%restart_dump)
1154
1155 safe_deallocate_a(prop%iter)
1156 ! This routine should maybe delete the files?
1157
1158 pop_sub(oct_prop_end)
1159 end subroutine oct_prop_end
1160 ! ---------------------------------------------------------
1161
1162
1163 ! ---------------------------------------------------------
1164 subroutine oct_prop_check(prop, namespace, space, psi, mesh, kpoints, iter)
1165 type(oct_prop_t), intent(inout) :: prop
1166 type(namespace_t), intent(in) :: namespace
1167 class(space_t), intent(in) :: space
1168 type(states_elec_t), intent(inout) :: psi
1169 class(mesh_t), intent(in) :: mesh
1170 type(kpoints_t), intent(in) :: kpoints
1171 integer, intent(in) :: iter
1172
1173 type(states_elec_t) :: stored_st
1174 character(len=80) :: dirname
1175 integer :: j, ierr
1176 complex(real64) :: overlap, prev_overlap
1177 real(real64), parameter :: WARNING_THRESHOLD = 1.0e-2_real64
1178
1179 push_sub(oct_prop_check)
1180
1181 do j = 1, prop%number_checkpoints + 2
1182 if (prop%iter(j) == iter) then
1183 call states_elec_copy(stored_st, psi)
1184 write(dirname,'(a, i4.4)') trim(prop%dirname), j
1185 call restart_open_dir(prop%restart_load, dirname, ierr)
1186 if (ierr == 0) then
1187 call states_elec_load(prop%restart_load, namespace, space, stored_st, mesh, kpoints, ierr, verbose=.false.)
1188 end if
1189 if (ierr /= 0) then
1190 message(1) = "Unable to read wavefunctions from '"//trim(dirname)//"'."
1191 call messages_fatal(1, namespace=namespace)
1192 end if
1193 call restart_close_dir(prop%restart_load)
1194 prev_overlap = zstates_elec_mpdotp(namespace, mesh, stored_st, stored_st)
1195 overlap = zstates_elec_mpdotp(namespace, mesh, stored_st, psi)
1196 if (abs(overlap - prev_overlap) > warning_threshold) then
1197 write(message(1), '(a,es13.4)') &
1198 "Forward-backward propagation produced an error of", abs(overlap-prev_overlap)
1199 write(message(2), '(a,i8)') "Iter = ", iter
1200 call messages_warning(2, namespace=namespace)
1201 end if
1202 ! Restore state only if the number of checkpoints is larger than zero.
1203 if (prop%number_checkpoints > 0) then
1204 call states_elec_end(psi)
1205 call states_elec_copy(psi, stored_st)
1206 end if
1207 call states_elec_end(stored_st)
1208 end if
1209 end do
1210 pop_sub(oct_prop_check)
1211 end subroutine oct_prop_check
1212 ! ---------------------------------------------------------
1213
1214
1215 ! ---------------------------------------------------------
1216 subroutine oct_prop_dump_states(prop, space, iter, psi, mesh, kpoints, ierr)
1217 type(oct_prop_t), intent(inout) :: prop
1218 class(space_t), intent(in) :: space
1219 integer, intent(in) :: iter
1220 type(states_elec_t), intent(in) :: psi
1221 class(mesh_t), intent(in) :: mesh
1222 type(kpoints_t), intent(in) :: kpoints
1223 integer, intent(out) :: ierr
1224
1225 integer :: j, err
1226 character(len=80) :: dirname
1227
1228 push_sub(oct_prop_dump_states)
1229
1230 ierr = 0
1231
1232 if (restart_skip(prop%restart_dump)) then
1233 pop_sub(oct_prop_dump_states)
1234 return
1235 end if
1236
1237 message(1) = "Debug: Writing OCT propagation states restart."
1238 call messages_info(1, debug_only=.true.)
1239
1240 do j = 1, prop%number_checkpoints + 2
1241 if (prop%iter(j) == iter) then
1242 write(dirname,'(a,i4.4)') trim(prop%dirname), j
1243 call restart_open_dir(prop%restart_dump, dirname, err)
1244 if (err == 0) then
1245 call states_elec_dump(prop%restart_dump, space, psi, mesh, kpoints, err, iter, verbose = .false.)
1246 end if
1247 if (err /= 0) then
1248 message(1) = "Unable to write wavefunctions to '"//trim(dirname)//"'."
1249 call messages_warning(1)
1250 ierr = ierr + 2**j
1251 end if
1252 call restart_close_dir(prop%restart_dump)
1253 end if
1254 end do
1255
1256 message(1) = "Debug: Writing OCT propagation states restart done."
1257 call messages_info(1, debug_only=.true.)
1258
1259 pop_sub(oct_prop_dump_states)
1260 end subroutine oct_prop_dump_states
1261 ! ---------------------------------------------------------
1262
1263
1264 ! ---------------------------------------------------------
1265 subroutine oct_prop_load_states(prop, namespace, space, psi, mesh, kpoints, iter, ierr)
1266 type(oct_prop_t), intent(inout) :: prop
1267 type(namespace_t), intent(in) :: namespace
1268 class(space_t), intent(in) :: space
1269 type(states_elec_t), intent(inout) :: psi
1270 class(mesh_t), intent(in) :: mesh
1271 type(kpoints_t), intent(in) :: kpoints
1272 integer, intent(in) :: iter
1273 integer, intent(out) :: ierr
1274
1275 integer :: j, err
1276 character(len=80) :: dirname
1277
1278 push_sub(oct_prop_load_states)
1279
1280 ierr = 0
1281
1282 if (restart_skip(prop%restart_load)) then
1283 ierr = -1
1284 pop_sub(oct_prop_load_states)
1285 return
1286 end if
1287
1288 message(1) = "Debug: Reading OCT propagation states restart."
1289 call messages_info(1, namespace=namespace, debug_only=.true.)
1290
1291 do j = 1, prop%number_checkpoints + 2
1292 if (prop%iter(j) == iter) then
1293 write(dirname,'(a, i4.4)') trim(prop%dirname), j
1294 call restart_open_dir(prop%restart_load, dirname, err)
1295 if (err == 0) then
1296 call states_elec_load(prop%restart_load, namespace, space, psi, mesh, kpoints, err, verbose=.false.)
1297 end if
1298 if (err /= 0) then
1299 message(1) = "Unable to read wavefunctions from '"//trim(dirname)//"'."
1300 call messages_warning(1, namespace=namespace)
1301 ierr = ierr + 2**j
1302 end if
1303 call restart_close_dir(prop%restart_load)
1304 end if
1305 end do
1306
1307 message(1) = "Debug: Reading OCT propagation states restart done."
1308 call messages_info(1, namespace=namespace, debug_only=.true.)
1310 pop_sub(oct_prop_load_states)
1311 end subroutine oct_prop_load_states
1312 ! ---------------------------------------------------------
1313
1314 ! ---------------------------------------------------------
1315 subroutine vlaser_operator_quadratic(laser, mesh, space, psi, hpsi)
1316 type(laser_t), intent(in) :: laser
1317 class(mesh_t), intent(in) :: mesh
1318 class(space_t), intent(in) :: space
1319 complex(real64), intent(inout) :: psi(:,:)
1320 complex(real64), intent(inout) :: hpsi(:,:)
1321
1322 integer :: ip
1323 logical :: vector_potential, magnetic_field
1324
1325 real(real64) :: a_field(3), a_field_prime(3), bb(3), b_prime(3)
1326 real(real64), allocatable :: aa(:, :), a_prime(:, :)
1327
1329
1330 a_field = m_zero
1331
1332 vector_potential = .false.
1333 magnetic_field = .false.
1334
1335 select case (laser_kind(laser))
1336 case (e_field_electric) ! do nothing
1337 case (e_field_magnetic)
1338 if (.not. allocated(aa)) then
1339 safe_allocate(aa(1:mesh%np_part, 1:space%dim))
1340 aa = m_zero
1341 safe_allocate(a_prime(1:mesh%np_part, 1:space%dim))
1342 end if
1343 a_prime = m_zero
1344 call laser_vector_potential(laser, mesh, a_prime)
1345 aa = aa + a_prime
1346 b_prime = m_zero
1347 call laser_field(laser, b_prime(1:space%dim))
1348 bb = bb + b_prime
1349 magnetic_field = .true.
1351 a_field_prime = m_zero
1352 call laser_field(laser, a_field_prime(1:space%dim))
1353 a_field = a_field + a_field_prime
1354 vector_potential = .true.
1355 end select
1356
1357 if (magnetic_field) then
1358 do ip = 1, mesh%np
1359 hpsi(ip, :) = hpsi(ip, :) + m_half * &
1360 dot_product(aa(ip, 1:space%dim), aa(ip, 1:space%dim)) * psi(ip, :) / p_c**2
1361 end do
1362 safe_deallocate_a(aa)
1363 safe_deallocate_a(a_prime)
1364 end if
1365 if (vector_potential) then
1366 do ip = 1, mesh%np
1367 hpsi(ip, :) = hpsi(ip, :) + m_half * &
1368 dot_product(a_field(1:space%dim), a_field(1:space%dim))*psi(ip, :) / p_c**2
1369 end do
1370 end if
1371
1373 end subroutine vlaser_operator_quadratic
1374
1375 ! ---------------------------------------------------------
1376 subroutine vlaser_operator_linear(laser, der, std, psi, hpsi, ik, gyromagnetic_ratio, a_static)
1377 type(laser_t), intent(in) :: laser
1378 type(derivatives_t), intent(in) :: der
1379 type(states_elec_dim_t), intent(in) :: std
1380 complex(real64), intent(inout) :: psi(:,:)
1381 complex(real64), intent(inout) :: hpsi(:,:)
1382 integer, intent(in) :: ik
1383 real(real64), intent(in) :: gyromagnetic_ratio
1384 real(real64), optional, intent(in) :: a_static(:,:)
1385
1386 integer :: ip, idim
1387 logical :: electric_field, vector_potential, magnetic_field
1388 complex(real64), allocatable :: grad(:, :, :), lhpsi(:, :)
1389
1390 real(real64) :: a_field(3), a_field_prime(3), bb(3), b_prime(3)
1391
1392 real(real64), allocatable :: vv(:), pot(:), aa(:, :), a_prime(:, :)
1393
1394 push_sub(vlaser_operator_linear)
1395
1396 a_field = m_zero
1397
1398 electric_field = .false.
1399 vector_potential = .false.
1400 magnetic_field = .false.
1401
1402 select case (laser_kind(laser))
1404 if (.not. allocated(vv)) then
1405 safe_allocate(vv(1:der%mesh%np))
1406 end if
1407 vv = m_zero
1408 call laser_potential(laser, der%mesh, vv)
1409 electric_field = .true.
1410
1411 case (e_field_electric)
1412 if (.not. allocated(vv)) then
1413 safe_allocate(vv(1:der%mesh%np))
1414 vv = m_zero
1415 safe_allocate(pot(1:der%mesh%np))
1416 end if
1417 pot = m_zero
1418 call laser_potential(laser, der%mesh, pot)
1419 vv = vv + pot
1420 electric_field = .true.
1421 safe_deallocate_a(pot)
1422
1423 case (e_field_magnetic)
1424 if (.not. allocated(aa)) then
1425 safe_allocate(aa(1:der%mesh%np_part, 1:der%dim))
1426 aa = m_zero
1427 safe_allocate(a_prime(1:der%mesh%np_part, 1:der%dim))
1428 end if
1429 a_prime = m_zero
1430 call laser_vector_potential(laser, der%mesh, a_prime)
1431 aa = aa + a_prime
1432 b_prime = m_zero
1433 call laser_field(laser, b_prime(1:der%dim))
1434 bb = bb + b_prime
1435 magnetic_field = .true.
1437 a_field_prime = m_zero
1438 call laser_field(laser, a_field_prime(1:der%dim))
1439 a_field = a_field + a_field_prime
1440 vector_potential = .true.
1441 end select
1442
1443 if (electric_field) then
1444 do idim = 1, std%dim
1445 hpsi(1:der%mesh%np, idim)= hpsi(1:der%mesh%np, idim) + vv(1:der%mesh%np) * psi(1:der%mesh%np, idim)
1446 end do
1447 safe_deallocate_a(vv)
1448 end if
1449
1450 if (magnetic_field) then
1451 safe_allocate(grad(1:der%mesh%np_part, 1:der%dim, 1:std%dim))
1452
1453 do idim = 1, std%dim
1454 call zderivatives_grad(der, psi(:, idim), grad(:, :, idim))
1455 end do
1456
1457 ! If there is a static magnetic field, its associated vector potential is coupled with
1458 ! the time-dependent one defined as a "laser" (ideally one should just add them all and
1459 ! do the calculation only once...). Note that h%ep%a_static already has been divided
1460 ! by P_c, and therefore here we only divide by P_c, and not P_c**2.
1461 !
1462 ! We put a minus sign, since for the moment vector potential for
1463 ! lasers and for the static magnetic field use a different
1464 ! convetion.
1465 if (present(a_static)) then
1466 do ip = 1, der%mesh%np
1467 hpsi(ip, :) = hpsi(ip, :) - dot_product(aa(ip, 1:der%dim), a_static(ip, 1:der%dim)) * psi(ip, :) / p_c
1468 end do
1469 end if
1470
1471 select case (std%ispin)
1473 do ip = 1, der%mesh%np
1474 hpsi(ip, 1) = hpsi(ip, 1) - m_zi * dot_product(aa(ip, 1:der%dim), grad(ip, 1:der%dim, 1)) / p_c
1475 end do
1476 case (spinors)
1477 do ip = 1, der%mesh%np
1478 do idim = 1, std%dim
1479 hpsi(ip, idim) = hpsi(ip, idim) - m_zi * &
1480 dot_product(aa(ip, 1:der%dim), grad(ip, 1:der%dim, idim)) / p_c
1481 end do
1482 end do
1483 end select
1484
1485
1486 select case (std%ispin)
1487 case (spin_polarized)
1488 safe_allocate(lhpsi(1:der%mesh%np, 1:std%dim))
1489 if (modulo(ik+1, 2) == 0) then ! we have a spin down
1490 lhpsi(1:der%mesh%np, 1) = - m_half / p_c * norm2(bb) * psi(1:der%mesh%np, 1)
1491 else
1492 lhpsi(1:der%mesh%np, 1) = + m_half / p_c * norm2(bb) * psi(1:der%mesh%np, 1)
1493 end if
1494 hpsi(1:der%mesh%np, :) = hpsi(1:der%mesh%np, :) + (gyromagnetic_ratio * m_half) * lhpsi(1:der%mesh%np, :)
1495 safe_deallocate_a(lhpsi)
1496
1497 case (spinors)
1498 safe_allocate(lhpsi(1:der%mesh%np, 1:std%dim))
1499 lhpsi(1:der%mesh%np, 1) = m_half / p_c * (bb(3) * psi(1:der%mesh%np, 1) &
1500 + cmplx(bb(1), -bb(2), real64) * psi(1:der%mesh%np, 2))
1501 lhpsi(1:der%mesh%np, 2) = m_half / p_c * (-bb(3) * psi(1:der%mesh%np, 2) &
1502 + cmplx(bb(1), bb(2), real64) * psi(1:der%mesh%np, 1))
1503 hpsi(1:der%mesh%np, :) = hpsi(1:der%mesh%np, :) + (gyromagnetic_ratio * m_half) * lhpsi(1:der%mesh%np, :)
1504 safe_deallocate_a(lhpsi)
1505 end select
1506
1507 safe_deallocate_a(grad)
1508 safe_deallocate_a(aa)
1509 safe_deallocate_a(a_prime)
1510 end if
1511
1512 if (vector_potential) then
1513 safe_allocate(grad(1:der%mesh%np_part, 1:der%dim, 1:std%dim))
1514
1515 do idim = 1, std%dim
1516 call zderivatives_grad(der, psi(:, idim), grad(:, :, idim))
1517 end do
1518
1519 select case (std%ispin)
1521 do ip = 1, der%mesh%np
1522 hpsi(ip, 1) = hpsi(ip, 1) - m_zi * dot_product(a_field(1:der%dim), grad(ip, 1:der%dim, 1)) / p_c
1523 end do
1524 case (spinors)
1525 do ip = 1, der%mesh%np
1526 do idim = 1, std%dim
1527 hpsi(ip, idim) = hpsi(ip, idim) - m_zi * &
1528 dot_product(a_field(1:der%dim), grad(ip, 1:der%dim, idim)) / p_c
1529 end do
1530 end do
1531 end select
1532 safe_deallocate_a(grad)
1533 end if
1534
1535 pop_sub(vlaser_operator_linear)
1536 end subroutine vlaser_operator_linear
1537
1538
1539end module propagation_oct_m
1540
1541!! Local Variables:
1542!! mode: f90
1543!! coding: utf-8
1544!! End:
batchified version of the BLAS axpy routine:
Definition: batch_ops.F90:154
scale a batch by a constant or vector
Definition: batch_ops.F90:162
constant times a vector plus a vector
Definition: lalg_basic.F90:171
integer, parameter, public mask_absorbing
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
Definition: batch_ops.F90:242
This module contains the definition of the data type that holds a "control function" used for OCT run...
subroutine, public controlfunction_to_h_val(cp, ext_partners, val)
subroutine, public controlfunction_to_realtime(par)
subroutine, public controlfunction_update(cp, cpp, dir, iter, mu, dd, dq)
Update the control function(s) given in "cp", according to the formula cp = (1 - mu) * cpp + mu * dd ...
subroutine, public controlfunction_end(cp)
real(real64) pure function, public controlfunction_alpha(par, ipar)
integer pure function, public controlfunction_number(par)
subroutine, public controlfunction_copy(cp_out, cp_in)
subroutine, public controlfunction_to_h(cp, ext_partners)
subroutine, public controlfunction_to_basis(par)
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:609
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
subroutine, public energy_calc_total(namespace, space, hm, gr, st, ext_partners, iunit, full)
This subroutine calculates the total energy of the system. Basically, it adds up the KS eigenvalues,...
subroutine, public epot_precalc_local_potential(ep, namespace, gr, ions)
Definition: epot.F90:574
logical function, public list_has_gauge_field(partners)
type(lasers_t) function, pointer, public list_get_lasers(partners)
subroutine, public forces_costate_calculate(gr, namespace, ions, hm, psi, chi, ff, qq)
Definition: forces.F90:218
subroutine, public forces_calculate(gr, namespace, ions, hm, ext_partners, st, ks, vhxc_old, t, dt)
Definition: forces.F90:339
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
complex(real64), parameter, public m_zi
Definition: global.F90:202
complex(real64), parameter, public m_z1
Definition: global.F90:199
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
Definition: global.F90:224
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public hamiltonian_elec_set_inh(hm, st)
subroutine, public zvmask(mesh, hm, st)
subroutine, public hamiltonian_elec_adjoint(hm)
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
subroutine, public hamiltonian_elec_remove_inh(hm)
subroutine, public hamiltonian_elec_not_adjoint(hm)
This module defines classes and functions for interaction partners.
subroutine, public ion_dynamics_verlet_step2(ions, v, fold, fnew, dt)
A bare verlet integrator.
subroutine, public ion_dynamics_verlet_step1(ions, q, v, fold, dt)
A bare verlet integrator.
subroutine, public ion_dynamics_restore_state(this, ions, state)
subroutine, public ion_dynamics_propagate(this, ions, time, dt, namespace)
subroutine, public ion_dynamics_save_state(this, ions, state)
logical pure function, public ion_dynamics_ions_move(this)
subroutine, public ion_dynamics_unfreeze(this)
Unfreezes the ionic movement.
logical function, public ion_dynamics_freeze(this)
Freezes the ionic movement.
A module to handle KS potential, without the external potential.
integer, parameter, public independent_particles
complex(real64) function, dimension(3), public laser_polarization(laser)
Definition: lasers.F90:726
subroutine, public laser_vector_potential(laser, mesh, aa, time)
Definition: lasers.F90:1078
integer, parameter, public e_field_electric
Definition: lasers.F90:177
integer, parameter, public e_field_vector_potential
Definition: lasers.F90:177
subroutine, public laser_potential(laser, mesh, pot, time)
Definition: lasers.F90:1043
integer, parameter, public e_field_scalar_potential
Definition: lasers.F90:177
integer pure elemental function, public laser_kind(laser)
Definition: lasers.F90:715
subroutine, public laser_field(laser, field, time)
Retrieves the value of either the electric or the magnetic field. If the laser is given by a scalar p...
Definition: lasers.F90:1115
integer, parameter, public e_field_magnetic
Definition: lasers.F90:177
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_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
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:414
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
subroutine, public oct_exchange_remove(this)
subroutine, public oct_exchange_set(this, st, mesh)
This module holds the "opt_control_state_t" datatype, which contains a quantum-classical state.
real(real64) function, dimension(:, :), pointer, public opt_control_point_p(ocs)
subroutine, public opt_control_set_classical(ions, ocs)
real(real64) function, dimension(:, :), pointer, public opt_control_point_q(ocs)
type(states_elec_t) function, pointer, public opt_control_point_qs(ocs)
subroutine, public opt_control_state_null(ocs)
subroutine, public opt_control_state_copy(ocsout, ocsin)
subroutine, public opt_control_get_classical(ions, ocs)
subroutine, public propagation_mod_init(niter, eta, delta, number_checkpoints, zbr98, gradients)
This subroutine must be called before any QOCT propagations are done. It simply stores in the module ...
subroutine, public oct_prop_end(prop)
subroutine update_hamiltonian_elec_chi(iter, namespace, space, gr, ks, hm, ext_partners, td, tg, par_chi, ions, st, qtildehalf)
subroutine update_hamiltonian_elec_psi(iter, namespace, space, gr, ks, hm, ext_partners, td, tg, par, st, ions)
subroutine oct_prop_load_states(prop, namespace, space, psi, mesh, kpoints, iter, ierr)
subroutine, public bwd_step(sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
subroutine update_field(iter, cp, space, gr, hm, ext_partners, ions, qcpsi, qcchi, cpp, dir)
Calculates the value of the control functions at iteration iter, from the state psi and the Lagrange-...
subroutine vlaser_operator_linear(laser, der, std, psi, hpsi, ik, gyromagnetic_ratio, a_static)
subroutine, public propagate_forward(sys, td, par, tg, qcpsi, prop, write_iter)
subroutine calculate_g(space, gr, hm, lasers, psi, chi, dl, dq)
subroutine, public oct_prop_check(prop, namespace, space, psi, mesh, kpoints, iter)
subroutine, public fwd_step(sys, td, tg, par, par_chi, qcpsi, prop_chi, prop_psi)
subroutine, public propagate_backward(sys, td, qcpsi, prop)
subroutine vlaser_operator_quadratic(laser, mesh, space, psi, hpsi)
subroutine, public oct_prop_init(prop, namespace, dirname, mesh, mc)
subroutine oct_prop_dump_states(prop, space, iter, psi, mesh, kpoints, ierr)
subroutine, public bwd_step_2(sys, td, tg, par, par_chi, qcchi, prop_chi, prop_psi)
integer, parameter, public prop_explicit_runge_kutta4
subroutine, public propagator_elec_copy(tro, tri)
subroutine, public propagator_elec_dt(ks, namespace, space, hm, gr, st, tr, time, dt, nt, ions_dyn, ions, ext_partners, outp, write_handler, scsteps, update_energy, qcchi)
Propagates st from time - dt to t. If dt<0, it propagates backwards from t+|dt| to t.
subroutine, public propagator_elec_remove_scf_prop(tr)
subroutine, public propagator_elec_end(tr)
integer, parameter, public restart_oct
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
logical pure function, public restart_skip(restart)
Returns true if the restart information should neither be read nor written. This might happen because...
Definition: restart.F90:969
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
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
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...
integer, parameter, public oct_tg_velocity
Definition: target_low.F90:212
integer, parameter, public oct_tg_hhgnew
Definition: target_low.F90:212
integer pure function, public target_type(tg)
Definition: target_low.F90:265
integer, parameter, public oct_targetmode_td
Definition: target_low.F90:227
integer pure function, public target_mode(tg)
Definition: target_low.F90:243
logical pure function, public target_move_ions(tg)
Definition: target_low.F90:281
subroutine, public target_init_propagation(tg)
This routine performs all the things that must be initialized prior to a forward evolution,...
Definition: target.F90:186
subroutine, public target_inh(psi, gr, kpoints, tg, time, inh, iter)
Calculates the inhomogeneous term that appears in the equation for chi, and places it into inh.
Definition: target.F90:487
subroutine, public target_tdcalc(tg, namespace, space, hm, gr, ions, ext_partners, psi, time, max_time)
Calculates, at a given point in time marked by the integer index, the integrand of the target functio...
Definition: target.F90:444
Definition: td.F90:114
subroutine, public td_write_data(writ)
Definition: td_write.F90:1190
subroutine, public td_write_end(writ)
Definition: td_write.F90:978
subroutine, public td_write_iter(writ, namespace, space, outp, gr, st, hm, ions, ext_partners, kick, ks, dt, iter, mc, recalculate_gs)
Definition: td_write.F90:1025
subroutine, public td_write_init(writ, namespace, space, outp, gr, st, hm, ions, ext_partners, ks, ions_move, with_gauge_field, kick, iter, max_iter, dt, mc)
Initialize files to write when prograting in time.
Definition: td_write.F90:371
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:730
This is the data type used to hold a control function.
class representing derivatives
Extension of space that contains the knowledge of the spin dimension.
Class describing the electron system.
Definition: electrons.F90:218
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Stores all communicators and groups.
Definition: multicomm.F90:206
This is the datatype that contains the objects that are propagated: in principle this could be both t...
class for organizing spins and k-points
The states_elec_t class contains all electronic wave functions.
int true(void)