82    integer :: in_file, out_file, ref_file, ii, jj, kk
 
   83    integer :: time_steps, time_steps_ref, energy_steps, istart, iend, ntiter
 
   84    real(real64)   :: dt, dt_ref, tt, ww, norm_Aext
 
   85    real(real64), 
allocatable :: vecpot(:, :), Aext(:), Eind_real(:, :), Eind_imag(:, :)
 
   86    complex(real64), 
allocatable :: dielectric(:), chi(:), invdielectric(:)
 
   87    real(real64), 
allocatable :: vecpot_ref(:, :)
 
   88    type(spectrum_t)  :: spectrum
 
   90    type(space_t)     :: space
 
   91    type(lattice_vectors_t) :: latt
 
   92    type(batch_t)     :: vecpotb, Eind_realb, Eind_imagb
 
   93    character(len=120) :: header
 
   94    real(real64) :: start_time
 
   95    character(len=MAX_PATH_LEN) :: ref_filename
 
   96    complex(real64),
allocatable :: Eind(:)
 
  103    safe_allocate(aext(1:space%dim))
 
  104    safe_allocate(eind(1:space%dim))
 
  116      message(1) = 
"Cannot find the GaugeVectorField in the input file" 
  121    if(space%dim > 1) 
then 
  122      message(1) = 
"This program assumes that the gauge field is along a  " 
  123      message(2) = 
"direction for which the dielectric tensor has no off-diagonal terms." 
  124      message(3) = 
"If this is not the case the dielectric function and the" 
  125      message(4) = 
"susceptibility will be wrong." 
  129    start_time = spectrum%start_time
 
  141    spectrum%start_time = ceiling(spectrum%start_time/dt)*dt
 
  161      if (time_steps_ref < time_steps) 
then 
  162        message(1) = 
"The reference calculation does not contain enought time steps" 
  166      if (.not. 
is_close(dt_ref, dt)) 
then 
  167        message(1) = 
"The time step of the reference calculation is different from the current calculation" 
  174    time_steps = time_steps + 1
 
  176    safe_allocate(vecpot(1:time_steps, 1:space%dim*3))
 
  180    do ii = 1, time_steps
 
  181      read(in_file, *) jj, tt, (vecpot(ii, kk), kk = 1, space%dim*3)
 
  188      time_steps_ref = time_steps_ref + 1
 
  189      safe_allocate(vecpot_ref(1:time_steps_ref, 1:space%dim*3))
 
  191      do ii = 1, time_steps_ref
 
  192        read(ref_file, *) jj, tt, (vecpot_ref(ii, kk), kk = 1, space%dim*3)
 
  195      do ii = 1, time_steps
 
  196        do kk = 1, space%dim*3
 
  197          vecpot(ii, kk) = vecpot(ii, kk) - vecpot_ref(ii, kk)
 
  202    write(
message(1), 
'(a, i7, a)') 
"Info: Read ", time_steps, 
" steps from file '"// &
 
  216    norm_aext = norm2(aext(1:space%dim))
 
  218    safe_allocate(eind_real(1:energy_steps, 1:space%dim))
 
  219    safe_allocate(eind_imag(1:energy_steps, 1:space%dim))
 
  222    call batch_init(vecpotb, 1, 1, space%dim, vecpot(:, space%dim+1:space%dim*2))
 
  223    call batch_init(eind_realb, 1, 1, space%dim, eind_real)
 
  224    call batch_init(eind_imagb, 1, 1, space%dim, eind_imag)
 
  226    call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, istart, iend, spectrum%start_time, dt, vecpotb)
 
  229      istart, iend, spectrum%start_time, dt, vecpotb, spectrum%min_energy, &
 
  230      spectrum%max_energy, spectrum%energy_step, eind_realb)
 
  233      istart, iend, spectrum%start_time, dt, vecpotb, spectrum%min_energy, &
 
  234      spectrum%max_energy, spectrum%energy_step, eind_imagb)
 
  238    call eind_realb%end()
 
  239    call eind_imagb%end()
 
  241    safe_allocate(invdielectric(1:energy_steps))
 
  242    safe_allocate(dielectric(1:energy_steps))
 
  243    safe_allocate(chi(1:energy_steps))
 
  245    do kk = 1, energy_steps
 
  246      ww = (kk-1)*spectrum%energy_step + spectrum%min_energy
 
  253      eind = cmplx(eind_real(kk, 1:space%dim), eind_imag(kk, 1:space%dim), real64)
 
  254      invdielectric(kk) = dot_product(aext, aext + eind) /norm_aext**2
 
  256      dielectric(kk) = 
m_one / invdielectric(kk)
 
  262    write(header, 
'(7a15)') 
'#        energy', 
'Re', 
'Im' 
  264    write(out_file,
'(a)') trim(header)
 
  265    do kk = 1, energy_steps
 
  266      ww = (kk-1)*spectrum%energy_step + spectrum%min_energy
 
  267      write(out_file, 
'(e15.6)', advance=
'no') ww
 
  268      write(out_file, 
'(2e15.6)', advance=
'no') real(invdielectric(kk), real64), aimag(invdielectric(kk))
 
  269      write(out_file, 
'()')
 
  274    write(out_file,
'(a)') trim(header)
 
  275    do kk = 1, energy_steps
 
  276      ww = (kk-1)*spectrum%energy_step + spectrum%min_energy
 
  277      write(out_file, 
'(e15.6)', advance=
'no') ww
 
  278      write(out_file, 
'(2e15.6)', advance=
'no') real(dielectric(kk), real64), aimag(dielectric(kk))
 
  279      write(out_file, 
'()')
 
  284    write(out_file,
'(a)') trim(header)
 
  285    do kk = 1, energy_steps
 
  286      ww = (kk-1)*spectrum%energy_step + spectrum%min_energy
 
  287      write(out_file, 
'(e15.6)', advance=
'no') ww
 
  288      write(out_file, 
'(2e15.6)', advance=
'no') real(chi(kk), real64), aimag(chi(kk))
 
  289      write(out_file, 
'()')
 
  293    safe_deallocate_a(dielectric)
 
  294    safe_deallocate_a(invdielectric)
 
  295    safe_deallocate_a(chi)
 
  296    safe_deallocate_a(vecpot)
 
  297    safe_deallocate_a(vecpot_ref)
 
  298    safe_deallocate_a(aext)
 
  299    safe_deallocate_a(eind)
 
  300    safe_deallocate_a(eind_real)
 
  301    safe_deallocate_a(eind_imag)
 
subroutine dielectric_function_compute()
program dielectric_function
A utility used to obtain the dielectric function from a kick calculation using the gauge field approa...
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_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_one
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()
character(len=max_path_len) function, public io_workpath(path, namespace)
construct path name from given name and namespace
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)
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
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.
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.
This module defines the unit system, used for input and output.
subroutine, public unit_system_init(namespace)