59 character(len=256) :: config_str
105 class(electrons_t),
pointer :: sys
107 character(MAX_PATH_LEN) :: basename, folder, ref_name, ref_folder, folder_default
108 integer :: c_start, c_end, c_step, c_start_default, length, c_how
109 logical :: iterate_folder, subtract_file
110 integer,
parameter :: CONVERT_FORMAT = 1, fourier_transform = 2, operation = 3
133 if (basename ==
" ") basename =
""
135 length = len_trim(basename)
137 if (basename(length-3:length) ==
'.obf')
then
138 basename = trim(basename(1:length-4))
170 if (iterate_folder)
then
171 folder_default =
'output_iter/td.'
174 folder_default =
'restart'
223 if (ref_name ==
" ") ref_name =
""
225 length = len_trim(ref_name)
227 if (ref_name(length-3:length) ==
'.obf')
then
228 ref_name = trim(ref_name(1:length-4))
255 CASE(fourier_transform)
258 c_start, c_end, c_step, sys%outp, subtract_file, &
259 ref_name, ref_folder)
263 c_start, c_end, c_step, sys%outp, iterate_folder, &
264 subtract_file, ref_name, ref_folder)
267 safe_deallocate_p(sys)
275 subroutine convert_low(mesh, namespace, space, ions, psolver, mc, basename, in_folder, c_start, c_end, c_step, outp, &
276 iterate_folder, subtract_file, ref_name, ref_folder)
277 class(
mesh_t),
intent(in) :: mesh
279 class(
space_t),
intent(in) :: space
280 type(
ions_t),
intent(in) :: ions
283 character(len=*),
intent(inout) :: basename
284 character(len=*),
intent(in) :: in_folder
285 integer,
intent(in) :: c_start
286 integer,
intent(in) :: c_end
287 integer,
intent(in) :: c_step
289 logical,
intent(in) :: iterate_folder
291 logical,
intent(in) :: subtract_file
292 character(len=*),
intent(inout) :: ref_name
293 character(len=*),
intent(inout) :: ref_folder
296 integer :: ierr, ii, folder_index, output_i
297 character(MAX_PATH_LEN) :: filename, out_name, folder, frmt, restart_folder
298 real(real64),
allocatable :: read_ff(:), read_rff(:), pot(:)
302 safe_allocate(read_ff(1:mesh%np))
303 safe_allocate(read_rff(1:mesh%np))
304 safe_allocate(pot(1:mesh%np))
307 write(
message(1),
'(5a,i5,a,i5,a,i5)')
"Converting '", trim(in_folder),
"/", trim(basename), &
308 "' from ", c_start,
" to ", c_end,
" every ", c_step
311 if (subtract_file)
then
312 write(
message(1),
'(a,a,a,a)')
"Reading ref-file from ", trim(ref_folder), trim(ref_name),
".obf"
314 dir=trim(ref_folder), mesh = mesh)
317 call restart%read_mesh_function(trim(ref_name), mesh, read_rff, ierr)
320 write(
message(1),
'(2a)')
"Failed to read from ref-file ", trim(ref_name)
321 write(
message(2),
'(2a)')
"from folder ", trim(ref_folder)
329 if (iterate_folder)
then
331 folder = in_folder(1:len_trim(in_folder)-1)
332 folder_index = index(folder,
'/', .
true.)
333 restart_folder = folder(1:folder_index)
335 restart_folder = in_folder
338 dir=trim(restart_folder), mesh = mesh)
340 do ii = c_start, c_end, c_step
341 if (iterate_folder)
then
343 write(folder,
'(a,i0.7,a)') in_folder(folder_index+1:len_trim(in_folder)-1),ii,
"/"
344 write(filename,
'(a,a,a)') trim(folder), trim(basename)
345 out_name = trim(basename)
348 if (c_start /= c_end)
then
352 write(frmt,
'(a,i0,a)')
"(a,i0.",10-len_trim(basename),
")"
353 write(filename, fmt=trim(frmt)) trim(basename), ii
354 write(out_name,
'(a)') trim(filename)
357 write(filename,
'(a,a,a,a)') trim(folder),
"/", trim(basename)
359 write(out_name,
'(a)') trim(basename)
364 call restart%read_mesh_function(trim(filename), mesh, read_ff, ierr)
367 write(
message(1),
'(a,a)')
"Error reading the file ", trim(filename)
368 write(
message(2),
'(a,i4)')
"Error code: ",ierr
369 write(
message(3),
'(a)')
"Skipping...."
373 if (subtract_file)
then
374 read_ff(:) = read_ff(:) - read_rff(:)
375 write(out_name,
'(a,a)') trim(out_name),
"-ref"
378 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
379 if (outp%how(output_i) /= 0)
then
381 trim(out_name), namespace, space, mesh, read_ff,
units_out%length**(-space%dim), ierr, &
382 pos=ions%pos, atoms=ions%atom)
385 if (outp%what(option__output__potential))
then
386 write(out_name,
'(a)')
"potential"
388 call dio_function_output(outp%how(option__output__potential), trim(restart_folder)//trim(folder), &
389 trim(out_name), namespace, space, mesh, pot,
units_out%energy, ierr, pos=ions%pos, atoms=ions%atom)
397 safe_deallocate_a(read_ff)
398 safe_deallocate_a(read_rff)
399 safe_deallocate_a(pot)
406 subroutine convert_transform(mesh, namespace, space, ions, mc, kpoints, basename, in_folder, c_start, c_end, c_step, outp, &
407 subtract_file, ref_name, ref_folder)
408 class(
mesh_t),
intent(in) :: mesh
410 class(
space_t),
intent(in) :: space
411 type(
ions_t),
intent(in) :: ions
414 character(len=*),
intent(inout) :: basename
415 character(len=*),
intent(in) :: in_folder
416 integer,
intent(in) :: c_start
417 integer,
intent(in) :: c_end
418 integer,
intent(in) :: c_step
420 logical,
intent(in) :: subtract_file
421 character(len=*),
intent(inout) :: ref_name
422 character(len=*),
intent(inout) :: ref_folder
424 integer :: ierr, i_space, i_time, nn(1:3), optimize_parity(1:3), wd_info, output_i
425 integer :: i_energy, e_end, e_start, e_point, chunk_size, read_count, t_point
426 logical :: optimize(1:3)
427 integer :: folder_index
428 character(MAX_PATH_LEN) :: filename, folder, restart_folder
429 real(real64) :: fdefault, w_max
430 real(real64),
allocatable :: read_ft(:), read_rff(:), point_tmp(:,:)
432 integer,
parameter :: FAST_FOURIER = 1, standard_fourier = 2
436 type(
batch_t) :: tdrho_b, wdrho_b
437 real(real64),
allocatable :: tdrho_a(:,:,:), wdrho_a(:,:,:)
441 complex(real64),
allocatable :: out_fft(:)
443 real(real64) :: start_time
444 integer :: time_steps
447 real(real64) :: max_energy
448 real(real64) :: min_energy
456 write(
message(1),
'(a)')
'Input: TDTimeStep must be positive.'
460 call io_mkdir(
'wd.general', namespace)
483 call parse_variable(namespace,
'ConvertReadSize', mesh%np, chunk_size)
486 if (chunk_size == 0) chunk_size = mesh%np
488 if (mesh%mpi_grp%size > 1 .and. chunk_size /= mesh%np)
then
489 write(
message(1),*)
'Incompatible value for ConvertReadSize and Parallelizaion in Domains'
490 write(
message(2),*)
'Use the default value for ConvertReadSize (or set it to 0)'
495 if (c_step <= 0)
then
496 write(
message(1),
'(a)')
'Input: ConvertStep must be positive for Fourier transform.'
500 start_time = c_start * dt
502 time_steps = (c_end - c_start) / c_step
503 if (time_steps <= 0)
then
504 write(
message(1),
'(a)')
'Input: ConvertEnd must be larger than ConvertStart for Fourier transform.'
519 if (max_energy > w_max)
then
520 write(
message(1),
'(a,f12.7)')
'Impossible to set ConvertEnergyMax to ', &
522 write(
message(2),
'(a)')
'ConvertEnergyMax is too large.'
523 write(
message(3),
'(a,f12.7,a)')
'ConvertEnergyMax reset to ', &
547 safe_allocate(read_ft(0:time_steps))
549 safe_allocate(read_rff(1:mesh%np))
551 select case (ft_method)
553 nn(1) = time_steps + 1
556 safe_allocate(out_fft(0:time_steps))
561 case (standard_fourier)
571 fdefault =
m_two *
m_pi / (dt * time_steps)
577 spectrum%start_time = c_start * dt
578 spectrum%end_time = c_end * dt
579 spectrum%energy_step = dw
580 spectrum%max_energy = max_energy
581 safe_allocate(tdrho_a(0:time_steps, 1, 1))
582 safe_allocate(wdrho_a(0:time_steps, 1, 1))
588 call kick_init(kick, namespace, space, kpoints, 1)
590 e_start = nint(min_energy / dw)
591 e_end = nint(max_energy / dw)
592 write(
message(1),
'(a,1x,i0.7,a,f12.7,a,i0.7,a,f12.7,a)')
'Frequency index:',e_start,
'(',&
598 if (subtract_file)
then
599 write(
message(1),
'(a,a,a,a)')
"Reading ref-file from ", trim(ref_folder), trim(ref_name),
".obf"
602 dir=trim(ref_folder), mesh = mesh)
605 call restart%read_mesh_function(trim(ref_name), mesh, read_rff, ierr)
608 write(
message(1),
'(2a)')
"Failed to read from ref-file ", trim(ref_name)
609 write(
message(2),
'(2a)')
"from folder ", trim(ref_folder)
615 do i_energy = e_start, e_end
616 write(filename,
'(a14,i0.7,a1)')
'wd.general/wd.',i_energy,
'/'
617 write(
message(1),
'(a,a,f12.7,a,1x,i7,a)')trim(filename),
' w =', &
622 call io_mkdir(trim(filename), namespace)
626 if (mesh%parallel_in_domains)
then
628 folder = in_folder(1:len_trim(in_folder)-1)
629 folder_index = index(folder,
'/', .
true.)
630 restart_folder = folder(1:folder_index)
632 dir=trim(restart_folder), mesh = mesh)
638 safe_allocate(point_tmp(1:chunk_size, 0:time_steps))
641 do i_space = 1, mesh%np
645 do i_time = c_start, c_end, c_step
647 if (mesh%parallel_in_domains .and. i_space == 1)
then
648 write(folder,
'(a,i0.7,a)') in_folder(folder_index+1:len_trim(in_folder)-1),i_time,
"/"
649 write(filename,
'(a,a,a)') trim(folder), trim(basename)
650 call restart%read_mesh_function(trim(filename), mesh, point_tmp(:, t_point), ierr)
654 write(folder,
'(a,i0.7,a)') in_folder(1:len_trim(in_folder)-1),i_time,
"/"
655 write(filename,
'(a,a,a,a)') trim(folder), trim(basename),
".obf"
656 if (mod(i_space-1, chunk_size) == 0)
then
662 if (i_time == c_start) read_count = 0
666 if (ierr /= 0 .and. i_space == 1)
then
667 write(
message(1),
'(a,a,2i10)')
"Error reading the file ", trim(filename), i_space, i_time
668 write(
message(2),
'(a)')
"Skipping...."
669 write(
message(3),
'(a,i0)')
"Error :", ierr
674 if (i_time == c_start) read_count = read_count + 1
675 if (subtract_file)
then
676 read_ft(t_point) = point_tmp(read_count, t_point) - read_rff(i_space)
678 read_ft(t_point) = point_tmp(read_count, t_point)
681 t_point = t_point + 1
684 select case (ft_method)
691 point_tmp(read_count, 0:time_steps) = aimag(out_fft(0:time_steps)) * dt
692 case (standard_fourier)
693 tdrho_a(0:time_steps, 1, 1) = read_ft(0:time_steps)
696 call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, c_start + 1, c_start + time_steps + 1, &
697 kick%time, dt, tdrho_b)
699 c_start + 1, c_start + time_steps + 1, kick%time, dt, tdrho_b, min_energy, max_energy, &
700 spectrum%energy_step, wdrho_b)
703 do e_point = e_start, e_end
704 point_tmp(read_count, e_point) = - wdrho_a(e_point, 1, 1)
708 if (mod(i_space-1, 1000) == 0 .and.
mpi_world%is_root())
then
713 if (mesh%mpi_grp%size == 1)
then
714 if (mod(i_space, chunk_size) == 0)
then
716 write(
message(2),
'(a,i0)')
"Writing binary output: step ", i_space/chunk_size
718 do i_energy = e_start, e_end
719 write(filename,
'(a14,i0.7,a12)')
'wd.general/wd.',i_energy,
'/density.obf'
721 if (i_space == chunk_size)
then
732 call mesh%mpi_grp%barrier()
734 if (mesh%parallel_in_domains)
then
735 do i_energy = e_start, e_end
736 write(filename,
'(a14,i0.7,a1)')
'wd.general/wd.',i_energy,
'/'
737 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
738 if (outp%how(output_i) /= 0)
then
740 trim(
'density'), namespace, space, mesh, point_tmp(:, i_energy), &
741 units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
748 if (any(outp%how /= option__outputformat__binary))
then
749 do i_energy = e_start, e_end
750 write(filename,
'(a14,i0.7,a1)')
'wd.general/wd.',i_energy,
'/'
752 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
753 if ((outp%how(output_i) /= 0) .and. (outp%how(output_i) /= option__outputformat__binary))
then
755 trim(
'density'), namespace, space, mesh, read_rff, &
756 units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
763 safe_deallocate_a(point_tmp)
764 safe_deallocate_a(read_ft)
765 safe_deallocate_a(read_rff)
767 select case (ft_method)
769 safe_deallocate_a(out_fft)
770 case (standard_fourier)
772 safe_deallocate_a(tdrho_a)
773 safe_deallocate_a(wdrho_a)
782 class(
mesh_t),
intent(in) :: mesh
784 class(
space_t),
intent(in) :: space
785 type(
ions_t),
intent(in) :: ions
789 integer :: ierr, ip, i_op, length, n_operations, output_i
793 real(real64) :: f_re, f_im
794 real(real64),
allocatable :: tmp_ff(:), scalar_ff(:)
796 character(len=200) :: var, scalar_expression
797 character(len=MAX_PATH_LEN) :: folder, filename, out_folder, out_filename
813 if (
parse_block(namespace,
'ConvertScalarOperation', blk) == 0)
then
817 if (n_operations == 0)
then
818 write(
message(1),
'(a)')
'No operations found. Check the input file'
829 call parse_variable(namespace,
'ConvertOutputFolder',
"convert", out_folder)
831 call io_mkdir(out_folder, namespace)
841 call parse_variable(namespace,
'ConvertOutputFilename',
'density', out_filename)
843 safe_allocate(tmp_ff(1:mesh%np))
844 safe_allocate(scalar_ff(1:mesh%np))
847 do i_op = 1, n_operations
856 length = len_trim(filename)
858 if (filename(length-3:length) ==
'.obf')
then
859 filename = trim(filename(1:length-4))
865 dir=trim(folder), mesh = mesh, exact=.
true.)
867 call restart%read_mesh_function(trim(filename), mesh, tmp_ff, ierr)
869 write(
message(1),
'(2a)')
"Failed to read from file ", trim(filename)
870 write(
message(2),
'(2a)')
"from folder ", trim(folder)
877 call parse_expression(f_re, f_im, trim(var), real(tmp_ff(ip), real64), trim(scalar_expression))
879 scalar_ff(ip) = scalar_ff(ip) + f_re
888 call mesh%mpi_grp%barrier()
893 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
894 if (outp%how(output_i) /= 0)
then
895 call dio_function_output(outp%how(output_i), trim(out_folder), trim(out_filename), namespace, space, mesh, &
896 scalar_ff, units, ierr, pos=ions%pos, atoms=ions%atom)
900 safe_deallocate_a(tmp_ff)
901 safe_deallocate_a(scalar_ff)
program oct_convert
This utility runs in parallel and can be used for post-processing of the results of Output.
subroutine convert_low(mesh, namespace, space, ions, psolver, mc, basename, in_folder, c_start, c_end, c_step, outp, iterate_folder, subtract_file, ref_name, ref_folder)
Giving a range of input files, it writes the corresponding output files.
subroutine convert_operate(mesh, namespace, space, ions, mc, outp)
Given a set of mesh function operations it computes a a resulting mesh function from linear combinati...
subroutine convert_transform(mesh, namespace, space, ions, mc, kpoints, basename, in_folder, c_start, c_end, c_step, outp, subtract_file, ref_name, ref_folder)
Giving a range of input files, it computes the Fourier transform of the file.
initialize a batch with existing memory
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
static void convert(multi *in, multi *out, int t_in, int t_out)
This module implements batches of mesh functions.
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
subroutine, public getopt_octopus(config_str)
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
subroutine, public fft_init(this, nn, dim, type, library, optimize, optimize_parity, comm, mpi_grp, use_aligned)
subroutine, public fft_all_init(namespace)
initialize the table
subroutine, public fft_all_end()
delete all plans
integer, parameter, public fft_real
integer, parameter, public fftlib_fftw
real(real64), parameter, public m_two
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
complex(real64), parameter, public m_z0
subroutine, public global_init(communicator)
Initialise Octopus.
subroutine, public dwrite_header(fname, np_global, ierr)
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
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_end()
subroutine, public io_mkdir(fname, namespace, parents)
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
subroutine, public kick_end(kick)
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
System information (time, memory, sysname)
subroutine, public loct_progress_bar(a, maxcount)
A wrapper around the progress bar, such that it can be silenced without needing to dress the call wit...
This module defines the meshes, which are used in Octopus.
subroutine, public messages_end()
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
character(len=512), private msg
subroutine, public messages_init(output_dir)
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public print_date(str)
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_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
type(namespace_t), public global_namespace
this module contains the low-level part of the output system
this module contains the output system
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
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 dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
subroutine, public profiling_end(namespace)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
integer, parameter, public restart_undefined
integer, parameter, public restart_type_load
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)
subroutine, public add_last_slash(str)
Adds a '/' in the end of the string, only if it missing. Useful for directories.
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)
type(unit_system_t), public units_inp
the units systems for reading and writing
This module is intended to contain simple general-purpose utility functions and procedures.
character(len=256) function, public get_config_opts()
character(len=256) function, public get_optional_libraries()
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Class defining batches of mesh functions.
Class describing the electron system.
Describes mesh distribution to nodes.
Stores all communicators and groups.