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