28 use,
intrinsic :: iso_fortran_env
44 integer :: iunit, ierr, ii, jj, iter, read_iter, ntime, nvel, ivel
45 integer :: istart, iend, energy_steps, out_file
46 real(real64),
allocatable :: time(:), velocities(:, :)
47 real(real64),
allocatable :: total_current(:, :), ftcurr(:, :, :), curr(:, :)
48 real(real64),
allocatable :: heat_current(:,:), ftheatcurr(:,:,:), heatcurr(:,:)
49 complex(real64),
allocatable :: invdielectric(:, :)
50 type(ions_t),
pointer :: ions
51 type(space_t) :: space
52 type(spectrum_t) :: spectrum
54 type(batch_t) :: currb, ftcurrb, heatcurrb, ftheatcurrb
55 real(real64) :: ww, curtime, deltat, v0
56 character(len=MAX_PATH_LEN) :: ref_filename
57 integer :: ref_file, time_steps_ref, kk
58 real(real64),
allocatable :: velcm(:), vel0(:), current(:), current_ref(:, :)
59 real(real64) :: dt_ref, tt, start_time
60 integer :: ifreq, max_freq
61 integer :: skip, smearing
62 logical :: from_forces
63 character(len=120) :: header
64 real(real64) :: excess_charge, qtot, javerage
114 if (spectrum%end_time <
m_zero) spectrum%end_time = huge(spectrum%end_time)
120 qtot = -(ions%val_charge() + excess_charge)
122 if (from_forces)
then
124 call messages_write(
'Info: Reading coordinates from td.general/coordinates')
134 nvel = ions%natoms*ions%space%dim
136 safe_allocate(time(1:ntime))
137 safe_allocate(velocities(1:nvel, 1:ntime))
148 read(unit = iunit, iostat = ierr, fmt = *) read_iter, curtime, &
149 ((ions%pos(jj, ii), jj = 1, 3), ii = 1, ions%natoms), &
150 ((ions%vel(jj, ii), jj = 1, 3), ii = 1, ions%natoms)
154 if (ierr /= 0 .or. curtime >= spectrum%end_time)
then
160 assert(iter == read_iter + 1)
162 if (curtime >= spectrum%start_time .and. mod(iter, skip) == 0)
then
164 time(ntime) = curtime
166 do ii = 1, ions%natoms
167 do jj = 1, ions%space%dim
186 call messages_write(
'Info: Reading total current from td.general/total_current')
191 if (iunit /= -1)
then
200 safe_allocate(total_current(1:space%dim, 1:ntime))
201 safe_allocate(heat_current(1:space%dim, 1:ntime))
202 safe_allocate(time(1:ntime))
205 read(iunit, *) read_iter, time(iter), (total_current(ii, iter), ii = 1, space%dim)
216 ref_file =
io_open(trim(ref_filename)//
'/total_current', action=
'read', status=
'old')
220 if (time_steps_ref < ntime)
then
221 message(1) =
"The reference calculation does not contain enought time steps"
225 if (.not.
is_close(dt_ref, time(2)-time(1)))
then
226 message(1) =
"The time step of the reference calculation is different from the current calculation"
231 time_steps_ref = time_steps_ref + 1
232 safe_allocate(current_ref(1:space%dim, 1:time_steps_ref))
234 do ii = 1, time_steps_ref
235 read(ref_file, *) jj, tt, (current_ref(kk, ii), kk = 1, space%dim)
240 total_current(kk, iter) = total_current(kk, iter) - current_ref(kk, iter)
243 safe_deallocate_a(current_ref)
245 start_time = spectrum%start_time
256 javerage = javerage + total_current(kk, iter)
258 javerage = javerage/ntime
260 total_current(kk, iter) = total_current(kk, iter) - javerage
267 call messages_write(
"Cannot find the 'td.general/total_current' file.")
273 safe_allocate(total_current(1:space%dim, 1:ntime))
274 safe_allocate(heat_current(1:space%dim, 1:ntime))
275 safe_allocate(time(1:ntime))
277 total_current(1:space%dim, 1:ntime) =
m_zero
283 if (iunit /= -1)
then
294 read(iunit, *) read_iter, time(iter), (heat_current(ii, iter), ii = 1, space%dim)
301 call messages_write(
"Cannot find the 'td.general/heat_current' file.")
305 call messages_write(
"Check ConductivityFromForces to possibly compute it from forces.")
308 heat_current(1:space%dim, 1:ntime) =
m_zero
314 safe_allocate(curr(ntime, 1:space%dim))
315 safe_allocate(heatcurr(ntime, 1:space%dim))
316 safe_allocate(vel0(1:space%dim))
318 safe_allocate(velcm(1:space%dim))
319 safe_allocate(current(1:space%dim))
325 if (from_forces)
then
329 do ii = 1, ions%natoms
331 vel0(jj) = vel0(jj) + velocities(ivel, iter)/real(ions%natoms, real64)
341 do ii = 1, ions%natoms
343 velcm(jj) = velcm(jj) + velocities(ivel, iter)/real(ions%natoms, real64)
344 current(jj) = current(jj) + &
345 ions%mass(ii)/ions%latt%rcell_volume*(velocities(ivel, iter) - vel0(jj))
350 curr(iter, 1:space%dim) = vel0(1:space%dim)*qtot/ions%latt%rcell_volume + current(1:space%dim)
353 curr(iter, 1:space%dim) = total_current(1:space%dim, iter)/ions%latt%rcell_volume
354 heatcurr(iter,1:space%dim) = heat_current(1:space%dim, iter)/ions%latt%rcell_volume
357 if (from_forces)
write(iunit,*) iter, iter*deltat, curr(iter, 1:space%dim)
361 safe_deallocate_a(velcm)
362 safe_deallocate_a(current)
365 if (from_forces)
call io_close(iunit)
369 istart = max(1, istart)
372 safe_allocate(ftcurr(1:energy_steps, 1:space%dim, 1:2))
375 call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, istart, iend, spectrum%start_time, deltat, currb)
377 call batch_init(ftcurrb, 1, 1, space%dim, ftcurr(:, :, 1))
379 istart, iend, spectrum%start_time, deltat, currb, spectrum%min_energy, spectrum%max_energy, &
380 spectrum%energy_step, ftcurrb)
383 call batch_init(ftcurrb, 1, 1, space%dim, ftcurr(:, :, 2))
385 istart, iend, spectrum%start_time, deltat, currb, spectrum%min_energy, spectrum%max_energy, &
386 spectrum%energy_step, ftcurrb)
396 write(unit = iunit, iostat = ierr, fmt =
'(a)') &
397 '###########################################################################################################################'
398 write(unit = iunit, iostat = ierr, fmt =
'(8a)')
'# HEADER'
399 write(unit = iunit, iostat = ierr, fmt =
'(a,a,a)') &
401 write(unit = iunit, iostat = ierr, fmt =
'(a)') &
402 '###########################################################################################################################'
404 v0 = norm2(vel0(1:space%dim))
405 if (.not. from_forces .or. v0 < epsilon(v0)) v0 =
m_one
414 vel0(1:space%dim) = vel0(1:space%dim) /
p_c
415 v0 = norm2(vel0(1:space%dim))
419 do ifreq = 1, energy_steps
420 ww = spectrum%energy_step*(ifreq - 1) + spectrum%min_energy
422 transpose(ftcurr(ifreq, 1:space%dim, 1:2)/v0)
431 safe_allocate(invdielectric(1:space%dim, 1:energy_steps))
432 do ifreq = 1, energy_steps
433 ww = max((ifreq-1)*spectrum%energy_step + spectrum%min_energy,
m_epsilon)
435 invdielectric(1:space%dim, ifreq) = (vel0(1:space%dim) +
m_four *
m_pi * &
436 cmplx(ftcurr(ifreq, 1:space%dim, 2),-ftcurr(ifreq, 1:space%dim, 1), real64) / ww) / v0
441 select case (space%dim)
443 write(header,
'(3a15)')
'# energy',
'Re x',
'Im x'
445 write(header,
'(5a15)')
'# energy',
'Re x',
'Im x',
'Re y',
'Im y'
447 write(header,
'(7a15)')
'# energy',
'Re x',
'Im x',
'Re y',
'Im y',
'Re z',
'Im z'
449 write(out_file,
'(a)') trim(header)
450 do ifreq = 1, energy_steps
451 ww = (ifreq-1)*spectrum%energy_step + spectrum%min_energy
452 select case (space%dim)
454 write(out_file,
'(3e15.6)') ww, &
455 real(invdielectric(1, ifreq), real64), aimag(invdielectric(1, ifreq))
457 write(out_file,
'(5e15.6)') ww, &
458 real(invdielectric(1, ifreq), real64), aimag(invdielectric(1, ifreq)), &
459 real(invdielectric(2, ifreq), real64), aimag(invdielectric(2, ifreq))
461 write(out_file,
'(7e15.6)') ww, &
462 real(invdielectric(1, ifreq), real64), aimag(invdielectric(1, ifreq)), &
463 real(invdielectric(2, ifreq), real64), aimag(invdielectric(2, ifreq)), &
464 real(invdielectric(3, ifreq), real64), aimag(invdielectric(3, ifreq))
469 safe_deallocate_a(ftcurr)
470 safe_deallocate_a(invdielectric)
474 if (from_forces)
then
475 safe_allocate(ftheatcurr(1:energy_steps, 1:space%dim, 1:2))
479 call batch_init(heatcurrb, 1, 1, space%dim, heatcurr)
483 call batch_init(ftheatcurrb, 1, 1, space%dim, ftheatcurr(:, :, 1))
485 1, ntime,
m_zero, deltat, heatcurrb, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftheatcurrb)
486 call ftheatcurrb%end()
488 call batch_init(ftheatcurrb, 1, 1, space%dim, ftheatcurr(:, :, 2))
490 1, ntime,
m_zero, deltat, heatcurrb, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftheatcurrb)
491 call ftheatcurrb%end()
500 write(unit = iunit, iostat = ierr, fmt =
'(a)') &
501 '###########################################################################################################################'
502 write(unit = iunit, iostat = ierr, fmt =
'(8a)')
'# HEADER'
503 write(unit = iunit, iostat = ierr, fmt =
'(a,a,a)') &
505 write(unit = iunit, iostat = ierr, fmt =
'(a)') &
506 '###########################################################################################################################'
508 v0 = norm2(vel0(1:space%dim))
509 if (.not. from_forces .or. v0 < epsilon(v0)) v0 =
m_one
511 do ifreq = 1, energy_steps
512 ww = spectrum%energy_step*(ifreq - 1) + spectrum%min_energy
514 transpose(ftheatcurr(ifreq, 1:space%dim, 1:2)/v0)
520 safe_deallocate_a(ftheatcurr)
525 safe_deallocate_a(vel0)
526 safe_deallocate_p(ions)
528 safe_deallocate_a(total_current)
529 safe_deallocate_a(heat_current)
530 safe_deallocate_a(curr)
531 safe_deallocate_a(heatcurr)
532 safe_deallocate_a(time)
initialize a batch with existing memory
This module implements batches of mesh functions.
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
type(mpi_comm), parameter, public serial_dummy_comm
Alias MPI_COMM_UNDEFINED for the specific use case of initialising Octopus utilities with no MPI supp...
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
real(real64), parameter, public m_epsilon
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
subroutine, public io_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
subroutine, public io_close(iunit, grp)
subroutine, public io_skip_header(iunit)
subroutine, public io_end()
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
This module is intended to contain "only mathematical" functions and procedures.
subroutine, public messages_end()
subroutine, public messages_init(output_dir)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_grp_t), public mpi_world
type(namespace_t), public global_namespace
logical function, public parse_is_defined(namespace, name)
subroutine, public parser_init()
Initialise the Octopus parser.
subroutine, public parser_end()
End the Octopus parser.
integer function, public parse_block(namespace, name, blk, check_varinfo_)
subroutine, public profiling_end(namespace)
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
integer, parameter, public smear_semiconductor
subroutine, public spectrum_fix_time_limits(spectrum, time_steps, dt, istart, iend, ntiter)
subroutine, public spectrum_fourier_transform(method, transform, noise, time_start, time_end, t0, time_step, time_function, energy_start, energy_end, energy_step, energy_function)
Computes the sine, cosine, (or "exponential") Fourier transform of the real function given in the tim...
subroutine, public spectrum_init(spectrum, namespace, default_energy_step, default_max_energy)
subroutine, public spectrum_signal_damp(damp_type, damp_factor, time_start, time_end, t0, time_step, time_function)
integer, parameter, public spectrum_transform_cos
integer, parameter, public spectrum_transform_sin
subroutine, public spectrum_count_time_steps(namespace, iunit, time_steps, dt)
pure integer function, public spectrum_nenergy_steps(spectrum)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
subroutine, public unit_system_init(namespace)