Octopus
propagator_etrs.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
22 use accel_oct_m
23 use batch_oct_m
24 use debug_oct_m
29 use grid_oct_m
31 use global_oct_m
36 use ions_oct_m
39 use lda_u_oct_m
41 use math_oct_m
46 use parser_oct_m
49 use space_oct_m
53 use types_oct_m
54 use v_ks_oct_m
57 use xc_oct_m
58
59 implicit none
60
61 private
62
63 public :: &
64 td_etrs, &
65 td_etrs_sc, &
66 td_aetrs, &
69
70contains
71
72 ! ---------------------------------------------------------
74 subroutine td_etrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, &
75 ions_dyn, ions, mc)
76 type(v_ks_t), intent(inout) :: ks
77 type(namespace_t), intent(in) :: namespace
78 type(electron_space_t), intent(in) :: space
79 type(hamiltonian_elec_t), intent(inout) :: hm
80 type(partner_list_t), intent(in) :: ext_partners
81 type(grid_t), intent(inout) :: gr
82 type(states_elec_t), intent(inout) :: st
83 type(propagator_base_t), intent(inout) :: tr
84 real(real64), intent(in) :: time
85 real(real64), intent(in) :: dt
86 type(ion_dynamics_t), intent(inout) :: ions_dyn
87 type(ions_t), intent(inout) :: ions
88 type(multicomm_t), intent(inout) :: mc
89
90 type(xc_copied_potentials_t) :: vhxc_t1, vhxc_t2
91 type(gauge_field_t), pointer :: gfield
92
93 push_sub(td_etrs)
94
95 if (hm%theory_level /= independent_particles) then
96
97 call hm%ks_pot%store_potentials(vhxc_t1)
98
99 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt, dt)
100
101 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
102 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
103
104 call hm%ks_pot%store_potentials(vhxc_t2)
105 call hm%ks_pot%restore_potentials(vhxc_t1)
106 call hm%update(gr, namespace, space, ext_partners, time = time - dt)
107
108 else
109
110 call propagation_ops_elec_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
111
112 end if
113
114 ! propagate dt/2 with H(t)
115
116 ! first move the ions to time t
117 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
118 ext_partners, mc, time, dt)
119
120 gfield => list_get_gauge_field(ext_partners)
121 if(associated(gfield)) then
122 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, dt, time)
123 end if
124
125 call hm%ks_pot%restore_potentials(vhxc_t2)
126
127 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time)
128
129 ! propagate dt/2 with H(time - dt)
130 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
131
132 pop_sub(td_etrs)
133 end subroutine td_etrs
134
135 ! ---------------------------------------------------------
137 subroutine td_etrs_sc(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, &
138 ions_dyn, ions, mc, sctol, scsteps)
139 type(v_ks_t), intent(inout) :: ks
140 type(namespace_t), intent(in) :: namespace
141 type(electron_space_t), intent(in) :: space
142 type(hamiltonian_elec_t), intent(inout) :: hm
143 type(partner_list_t), intent(in) :: ext_partners
144 type(grid_t), intent(inout) :: gr
145 type(states_elec_t), intent(inout) :: st
146 type(propagator_base_t), intent(inout) :: tr
147 real(real64), intent(in) :: time
148 real(real64), intent(in) :: dt
149 type(ion_dynamics_t), intent(inout) :: ions_dyn
150 type(ions_t), intent(inout) :: ions
151 type(multicomm_t), intent(inout) :: mc
152 real(real64), intent(in) :: sctol
153 integer, optional, intent(out) :: scsteps
154
155 real(real64) :: diff
156 integer :: ik, ib, iter
157 class(wfs_elec_t), allocatable :: psi2(:, :)
158 ! these are hardcoded for the moment
159 integer, parameter :: niter = 10
160 type(gauge_field_t), pointer :: gfield
161 type(xc_copied_potentials_t) :: vhxc_t1, vhxc_t2
162
163 push_sub(td_etrs_sc)
164
165 assert(hm%theory_level /= independent_particles)
166
167 call hm%ks_pot%store_potentials(vhxc_t1)
168
170 call messages_write(' Self-consistency iteration:')
171 call messages_info(namespace=namespace)
172
173 !Propagate the states to t+dt/2 and compute the density at t+dt
174 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt, dt)
175
176 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
177 calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
178
179
180 call hm%ks_pot%store_potentials(vhxc_t2)
181 call hm%ks_pot%restore_potentials(vhxc_t1)
182
183 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time - dt)
184
185 ! propagate dt/2 with H(t)
186
187 ! first move the ions to time t
188 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
189 ext_partners, mc, time, dt)
190
191 gfield => list_get_gauge_field(ext_partners)
192 if(associated(gfield)) then
193 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, dt, time)
194 end if
195
196 call hm%ks_pot%restore_potentials(vhxc_t2)
197
198 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time)
199
200 safe_allocate_type_array(wfs_elec_t, psi2, (st%group%block_start:st%group%block_end, st%d%kpt%start:st%d%kpt%end))
201
202 ! store the state at half iteration
203 do ik = st%d%kpt%start, st%d%kpt%end
204 do ib = st%group%block_start, st%group%block_end
205 call st%group%psib(ib, ik)%copy_to(psi2(ib, ik), copy_data=.true.)
206 end do
207 end do
208
209 do iter = 1, niter
210
211 call hm%ks_pot%store_potentials(vhxc_t2)
212
213 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, m_half * dt)
214
215 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
216 time = time, calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
217 call lda_u_update_occ_matrices(hm%lda_u, namespace, gr, st, hm%phase, hm%energy)
218
219 ! now check how much the potential changed
220 diff = hm%ks_pot%check_convergence(vhxc_t2, gr, st%rho, st%qtot)
221
222 call messages_write(' step ')
223 call messages_write(iter)
224 call messages_write(', residue = ')
225 call messages_write(abs(diff), fmt = '(1x,es9.2)')
226 call messages_info(namespace=namespace)
227
228 if (diff <= sctol) exit
229
230 if (iter /= niter) then
231 ! we are not converged, restore the states
232 do ik = st%d%kpt%start, st%d%kpt%end
233 do ib = st%group%block_start, st%group%block_end
234 call psi2(ib, ik)%copy_data_to(gr%np, st%group%psib(ib, ik))
235 end do
236 end do
237 end if
238
239 end do
240
241 if (hm%lda_u_level /= dft_u_none) then
242 call lda_u_write_u(hm%lda_u, namespace=namespace)
243 call lda_u_write_v(hm%lda_u, namespace=namespace)
244 end if
245
246 ! print an empty line
247 call messages_info(namespace=namespace)
248
249 if (present(scsteps)) scsteps = iter
250
251 do ik = st%d%kpt%start, st%d%kpt%end
252 do ib = st%group%block_start, st%group%block_end
253 call psi2(ib, ik)%end()
254 end do
255 end do
256
257 safe_deallocate_a(psi2)
258
259 pop_sub(td_etrs_sc)
260 end subroutine td_etrs_sc
261
262 ! ---------------------------------------------------------
264 subroutine td_aetrs(namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
265 type(namespace_t), intent(in) :: namespace
266 type(electron_space_t), intent(in) :: space
267 type(hamiltonian_elec_t), intent(inout) :: hm
268 type(grid_t), intent(inout) :: gr
269 type(states_elec_t), intent(inout) :: st
270 type(propagator_base_t), intent(inout) :: tr
271 real(real64), intent(in) :: time
272 real(real64), intent(in) :: dt
273 type(ion_dynamics_t), intent(inout) :: ions_dyn
274 type(ions_t), intent(inout) :: ions
275 type(partner_list_t), intent(in) :: ext_partners
276 type(multicomm_t), intent(inout) :: mc
277
278 type(gauge_field_t), pointer :: gfield
279
280 push_sub(td_aetrs)
281
282 ! propagate half of the time step with H(time - dt)
283 call propagation_ops_elec_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
284
285 !Get the potentials from the interpolation
286 call propagation_ops_elec_interpolate_get(hm, tr%vks_old)
287
288 ! move the ions to time t
289 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, &
290 ions, ext_partners, mc, time, dt)
291
292 !Propagate gauge field
293 gfield => list_get_gauge_field(ext_partners)
294 if(associated(gfield)) then
295 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, dt, time)
296 end if
297
298 !Update Hamiltonian
299 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time)
300
301 !Do the time propagation for the second half of the time step
302 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
303
304 pop_sub(td_aetrs)
305 end subroutine td_aetrs
306
307 ! ---------------------------------------------------------
309 subroutine td_aetrs_sc(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, &
310 ions_dyn, ions, mc, sctol, scsteps)
311 type(v_ks_t), intent(inout) :: ks
312 type(namespace_t), intent(in) :: namespace
313 type(electron_space_t), intent(in) :: space
314 type(hamiltonian_elec_t), intent(inout) :: hm
315 type(partner_list_t), intent(in) :: ext_partners
316 type(grid_t), intent(inout) :: gr
317 type(states_elec_t), intent(inout) :: st
318 type(propagator_base_t), intent(inout) :: tr
319 real(real64), intent(in) :: time
320 real(real64), intent(in) :: dt
321 type(ion_dynamics_t), intent(inout) :: ions_dyn
322 type(ions_t), intent(inout) :: ions
323 type(multicomm_t), intent(inout) :: mc
324 real(real64), intent(in) :: sctol
325 integer, optional, intent(out) :: scsteps
326
327 real(real64) :: diff
328 integer :: iter
329 type(states_elec_t) :: st0
330 integer, parameter :: niter = 10
331 type(gauge_field_t), pointer :: gfield
332 type(xc_copied_potentials_t) :: vhxc_t2
333
334 push_sub(td_aetrs_sc)
335
336 assert(hm%theory_level /= independent_particles)
337
338 call messages_new_line()
339 call messages_write(' Self-consistency iteration:')
340 call messages_info(namespace=namespace)
341
342 ! First half-step with H(t - dt)
343 call propagation_ops_elec_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
344
345 ! Initial guess for H(t) from extrapolated potentials
346 call propagation_ops_elec_interpolate_get(hm, tr%vks_old)
347
348 ! Move ions and gauge to t only once before SC loop
349 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, ions, &
350 ext_partners, mc, time, dt)
351
352 gfield => list_get_gauge_field(ext_partners)
353 if (associated(gfield)) then
354 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, dt, time)
355 end if
356
357 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time)
358
359 ! Store the half-step state to restart each SC iteration
360 call states_elec_copy(st0, st)
361
362 do iter = 1, niter
363 call hm%ks_pot%store_potentials(vhxc_t2)
364
365 call propagation_ops_elec_fuse_density_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
366
367 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
368 time = time, calc_current = .false., calc_energy = .false., calc_eigenval = .false.)
369 call lda_u_update_occ_matrices(hm%lda_u, namespace, gr, st, hm%phase, hm%energy)
370
371 diff = hm%ks_pot%check_convergence(vhxc_t2, gr, st%rho, st%qtot)
372
373 call messages_write(' step ')
374 call messages_write(iter)
375 call messages_write(', residue = ')
376 call messages_write(abs(diff), fmt = '(1x,es9.2)')
377 call messages_info(namespace=namespace)
378
379 if (diff <= sctol) exit
380
381 if (iter /= niter) then
382 ! restore initial states for next iteration
383 call states_elec_group_copy(st0%d, st0%group, st%group)
384 end if
385 end do
386
387 if (hm%lda_u_level /= dft_u_none) then
388 call lda_u_write_u(hm%lda_u, namespace=namespace)
389 call lda_u_write_v(hm%lda_u, namespace=namespace)
390 end if
391
392 call messages_info(namespace=namespace)
393
394 if (present(scsteps)) scsteps = iter
395
396 call states_elec_end(st0)
397
398 pop_sub(td_aetrs_sc)
399 end subroutine td_aetrs_sc
400
401 ! ---------------------------------------------------------
403 subroutine td_caetrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, &
404 ions_dyn, ions, mc)
405 type(v_ks_t), intent(inout) :: ks
406 type(namespace_t), intent(in) :: namespace
407 type(electron_space_t), intent(in) :: space
408 type(hamiltonian_elec_t), intent(inout) :: hm
409 type(partner_list_t), intent(in) :: ext_partners
410 type(grid_t), intent(inout) :: gr
411 type(states_elec_t), intent(inout) :: st
412 type(propagator_base_t), intent(inout) :: tr
413 real(real64), intent(in) :: time
414 real(real64), intent(in) :: dt
415 type(ion_dynamics_t), intent(inout) :: ions_dyn
416 type(ions_t), intent(inout) :: ions
417 type(multicomm_t), intent(inout) :: mc
418
419 integer :: ik, ispin, ip, ist, ib
420 real(real64) :: vv
421 complex(real64) :: phase
422 type(density_calc_t) :: dens_calc
423 integer(int64) :: pnp
424 integer(int64), dimension(3) :: gsizes, bsizes
425 type(accel_mem_t) :: phase_buff
426 type(gauge_field_t), pointer :: gfield
427 type(xc_copied_potentials_t) :: vold
428
429 push_sub(td_caetrs)
430
431 ! Get the interpolated KS potentials into vold
432 call hm%ks_pot%get_interpolated_potentials(tr%vks_old, 2, storage=vold)
433 ! And set it to the Hamiltonian
434 call hm%ks_pot%restore_potentials(vold)
435
436 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time - dt)
437
438 call v_ks_calc_start(ks, namespace, space, hm, st, ions, ions%latt, ext_partners, &
439 time = time - dt, calc_energy = .false.)
440
441 ! propagate half of the time step with H(time - dt)
442 call propagation_ops_elec_exp_apply(tr%te, namespace, st, gr, hm, m_half*dt)
443
444 call v_ks_calc_finish(ks, hm, namespace, space, ions%latt, st, ext_partners)
445
446 call hm%ks_pot%set_interpolated_potentials(tr%vks_old, 1)
447
448 call hm%ks_pot%perform_interpolation(tr%vks_old, (/time - dt, time - m_two*dt, time - m_three*dt/), time)
449
450 ! Replace vold by 0.5(vhxc+vold)
451 call hm%ks_pot%mix_potentials(vold, dt)
452
453 ! copy vold to a cl buffer
454 if (accel_is_enabled() .and. hm%apply_packed()) then
455 if (family_is_mgga_with_exc(hm%xc)) then
456 call messages_not_implemented('CAETRS propagator with accel and MGGA with energy functionals', namespace=namespace)
457 end if
458 pnp = accel_padded_size(gr%np)
459 call accel_create_buffer(phase_buff, accel_mem_read_only, type_float, pnp*st%d%nspin)
460 call vold%copy_vhxc_to_buffer(int(gr%np, int64), st%d%nspin, pnp, phase_buff)
461 end if
462
463 !Get the potentials from the interpolator
464 call propagation_ops_elec_interpolate_get(hm, tr%vks_old)
465
466 ! move the ions to time t
467 call propagation_ops_elec_move_ions(tr%propagation_ops_elec, gr, hm, st, namespace, space, ions_dyn, &
468 ions, ext_partners, mc, time, dt)
469
470 gfield => list_get_gauge_field(ext_partners)
471 if(associated(gfield)) then
472 call propagation_ops_elec_propagate_gauge_field(tr%propagation_ops_elec, gfield, dt, time)
473 end if
474
475 call propagation_ops_elec_update_hamiltonian(namespace, space, st, gr, hm, ext_partners, time)
476
477 call density_calc_init(dens_calc, st, gr, st%rho)
478
479 ! propagate the other half with H(t)
480 do ik = st%d%kpt%start, st%d%kpt%end
481 ispin = st%d%get_spin_index(ik)
482
483 do ib = st%group%block_start, st%group%block_end
484 if (hm%apply_packed()) then
485 call st%group%psib(ib, ik)%do_pack()
486 if (hamiltonian_elec_inh_term(hm)) call hm%inh_st%group%psib(ib, ik)%do_pack()
487 end if
488
489 call profiling_in("CAETRS_PHASE")
490 select case (st%group%psib(ib, ik)%status())
491 case (batch_not_packed)
492 do ip = 1, gr%np
493 vv = vold%vhxc(ip, ispin)
494 phase = cmplx(cos(vv), -sin(vv), real64)
495 do ist = 1, st%group%psib(ib, ik)%nst_linear
496 st%group%psib(ib, ik)%zff_linear(ip, ist) = st%group%psib(ib, ik)%zff_linear(ip, ist)*phase
497 end do
498 end do
499 case (batch_packed)
500 do ip = 1, gr%np
501 vv = vold%vhxc(ip, ispin)
502 phase = cmplx(cos(vv), -sin(vv), real64)
503 do ist = 1, st%group%psib(ib, ik)%nst_linear
504 st%group%psib(ib, ik)%zff_pack(ist, ip) = st%group%psib(ib, ik)%zff_pack(ist, ip)*phase
505 end do
506 end do
508 call accel_set_kernel_arg(kernel_phase, 0, pnp*(ispin - 1))
509 call accel_set_kernel_arg(kernel_phase, 1, phase_buff)
510 call accel_set_kernel_arg(kernel_phase, 2, st%group%psib(ib, ik)%ff_device)
511 call accel_set_kernel_arg(kernel_phase, 3, log2(st%group%psib(ib, ik)%pack_size(1)))
512
513 ! Compute the grid (extend to another dimensions if the size of the problem is too big)
514 call accel_grid_size_extend_dim(pnp, st%group%psib(ib, ik)%pack_size(1), gsizes, bsizes)
515
516 call accel_kernel_run(kernel_phase, gsizes, bsizes)
517
518 end select
519 call profiling_out("CAETRS_PHASE")
520
521 call hm%phase%set_phase_corr(gr, st%group%psib(ib, ik))
522 if (hamiltonian_elec_inh_term(hm)) then
523 call tr%te%apply_batch(namespace, gr, hm, st%group%psib(ib, ik), m_half*dt, &
524 inh_psib = hm%inh_st%group%psib(ib, ik))
525 else
526 call tr%te%apply_batch(namespace, gr, hm, st%group%psib(ib, ik), m_half*dt)
527 end if
528 call hm%phase%unset_phase_corr(gr, st%group%psib(ib, ik))
529
530 call density_calc_accumulate(dens_calc, st%group%psib(ib, ik))
531
532 if (hm%apply_packed()) then
533 call st%group%psib(ib, ik)%do_unpack()
534 if (hamiltonian_elec_inh_term(hm)) call hm%inh_st%group%psib(ib, ik)%do_unpack()
535 end if
536 end do
537 end do
538
539 call density_calc_end(dens_calc)
540
541 if (accel_is_enabled() .and. hm%apply_packed()) then
542 call accel_free_buffer(phase_buff)
543 end if
544
545 pop_sub(td_caetrs)
546 end subroutine td_caetrs
547
548end module propagator_etrs_oct_m
549
550!! Local Variables:
551!! mode: f90
552!! coding: utf-8
553!! End:
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
subroutine, public accel_free_buffer(this, async)
Definition: accel.F90:1005
type(accel_kernel_t), target, save, public kernel_phase
Definition: accel.F90:274
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
integer, parameter, public accel_mem_read_only
Definition: accel.F90:185
This module implements batches of mesh functions.
Definition: batch.F90:135
integer, parameter, public batch_not_packed
functions are stored in CPU memory, unpacked order
Definition: batch.F90:286
integer, parameter, public batch_device_packed
functions are stored in device memory in packed order
Definition: batch.F90:286
integer, parameter, public batch_packed
functions are stored in CPU memory, in transposed (packed) order
Definition: batch.F90:286
This module implements a calculator for the density and defines related functions.
Definition: density.F90:122
subroutine, public density_calc_end(this, allreduce, symmetrize)
Finalize the density calculation.
Definition: density.F90:565
subroutine, public density_calc_accumulate(this, psib, async)
Accumulate weighted orbital densities for a batch psib.
Definition: density.F90:396
subroutine, public density_calc_init(this, st, gr, density)
initialize the density calculator
Definition: density.F90:190
type(gauge_field_t) function, pointer, public list_get_gauge_field(partners)
real(real64), parameter, public m_two
Definition: global.F90:193
integer, parameter, public independent_particles
Theory level.
Definition: global.F90:241
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_three
Definition: global.F90:194
This module implements the underlying real-space grid.
Definition: grid.F90:119
pure logical function, public hamiltonian_elec_inh_term(hm)
This module defines classes and functions for interaction partners.
A module to handle KS potential, without the external potential.
subroutine, public lda_u_write_u(this, iunit, namespace)
Definition: lda_u_io.F90:533
subroutine, public lda_u_write_v(this, iunit, namespace)
Definition: lda_u_io.F90:582
integer, parameter, public dft_u_none
Definition: lda_u.F90:205
subroutine, public lda_u_update_occ_matrices(this, namespace, mesh, st, phase, energy)
Definition: lda_u.F90:892
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
subroutine, public messages_new_line()
Definition: messages.F90:1088
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
subroutine, public propagation_ops_elec_propagate_gauge_field(wo, gfield, dt, time, save_gf)
subroutine, public propagation_ops_elec_update_hamiltonian(namespace, space, st, mesh, hm, ext_partners, time)
subroutine, public propagation_ops_elec_interpolate_get(hm, vks_old)
subroutine, public propagation_ops_elec_fuse_density_exp_apply(te, namespace, st, gr, hm, dt, dt2, op)
subroutine, public propagation_ops_elec_move_ions(wo, gr, hm, st, namespace, space, ions_dyn, ions, ext_partners, mc, time, dt, save_pos)
subroutine, public propagation_ops_elec_exp_apply(te, namespace, st, mesh, hm, dt, op)
subroutine, public td_aetrs(namespace, space, hm, gr, st, tr, time, dt, ions_dyn, ions, ext_partners, mc)
Propagator with approximate enforced time-reversal symmetry.
subroutine, public td_etrs_sc(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc, sctol, scsteps)
Propagator with enforced time-reversal symmetry and self-consistency.
subroutine, public td_etrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc)
Propagator with enforced time-reversal symmetry.
subroutine, public td_caetrs(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc)
Propagator with approximate enforced time-reversal symmetry.
subroutine, public td_aetrs_sc(ks, namespace, space, hm, ext_partners, gr, st, tr, time, dt, ions_dyn, ions, mc, sctol, scsteps)
Propagator with approximate enforced time-reversal symmetry and self-consistency.
This module handles spin dimensions of the states and the k-point distribution.
This module handles groups of electronic batches and their parallel distribution.
subroutine, public states_elec_group_copy(d, group_in, group_out, copy_data, special)
make a copy of a group
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
type(type_t), public type_float
Definition: types.F90:135
subroutine, public v_ks_calc_finish(ks, hm, namespace, space, latt, st, ext_partners, force_semilocal)
Definition: v_ks.F90:1106
subroutine, public v_ks_calc_start(ks, namespace, space, hm, st, ions, latt, ext_partners, time, calc_energy, force_semilocal)
This routine starts the calculation of the Kohn-Sham potential. The routine v_ks_calc_finish must be ...
Definition: v_ks.F90:809
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:754
Definition: xc.F90:116
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:599
The calculator object.
Definition: density.F90:171
Extension of space that contains the knowledge of the spin dimension.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Stores all communicators and groups.
Definition: multicomm.F90:208
The states_elec_t class contains all electronic wave functions.
batches of electronic states
Definition: wfs_elec.F90:141
int true(void)