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
117 message(1) =
'Info: Converting files'
132 if (basename ==
" ") basename =
""
134 length = len_trim(basename)
136 if (basename(length-3:length) ==
'.obf')
then
137 basename = trim(basename(1:length-4))
169 if (iterate_folder)
then
170 folder_default =
'output_iter/td.'
173 folder_default =
'restart'
222 if (ref_name ==
" ") ref_name =
""
224 length = len_trim(ref_name)
226 if (ref_name(length-3:length) ==
'.obf')
then
227 ref_name = trim(ref_name(1:length-4))
254 CASE(fourier_transform)
257 c_start, c_end, c_step, sys%outp, subtract_file, &
258 ref_name, ref_folder)
262 c_start, c_end, c_step, sys%outp, iterate_folder, &
263 subtract_file, ref_name, ref_folder)
266 safe_deallocate_p(sys)
274 subroutine convert_low(mesh, namespace, space, ions, psolver, mc, basename, in_folder, c_start, c_end, c_step, outp, &
275 iterate_folder, subtract_file, ref_name, ref_folder)
276 class(
mesh_t),
intent(in) :: mesh
278 class(
space_t),
intent(in) :: space
279 type(
ions_t),
intent(in) :: ions
282 character(len=*),
intent(inout) :: basename
283 character(len=*),
intent(in) :: in_folder
284 integer,
intent(in) :: c_start
285 integer,
intent(in) :: c_end
286 integer,
intent(in) :: c_step
288 logical,
intent(in) :: iterate_folder
290 logical,
intent(in) :: subtract_file
291 character(len=*),
intent(inout) :: ref_name
292 character(len=*),
intent(inout) :: ref_folder
295 integer :: ierr, ii, folder_index, output_i
296 character(MAX_PATH_LEN) :: filename, out_name, folder, frmt, restart_folder
297 real(real64),
allocatable :: read_ff(:), read_rff(:), pot(:)
301 safe_allocate(read_ff(1:mesh%np))
302 safe_allocate(read_rff(1:mesh%np))
303 safe_allocate(pot(1:mesh%np))
306 write(
message(1),
'(5a,i5,a,i5,a,i5)')
"Converting '", trim(in_folder),
"/", trim(basename), &
307 "' from ", c_start,
" to ", c_end,
" every ", c_step
310 if (subtract_file)
then
311 write(
message(1),
'(a,a,a,a)')
"Reading ref-file from ", trim(ref_folder), trim(ref_name),
".obf"
313 dir=trim(ref_folder), mesh = mesh)
316 call restart%read_mesh_function(trim(ref_name), mesh, read_rff, ierr)
319 write(
message(1),
'(2a)')
"Failed to read from ref-file ", trim(ref_name)
320 write(
message(2),
'(2a)')
"from folder ", trim(ref_folder)
328 if (iterate_folder)
then
330 folder = in_folder(1:len_trim(in_folder)-1)
331 folder_index = index(folder,
'/', .
true.)
332 if(folder_index > 0) restart_folder = folder(1:folder_index)
334 restart_folder = in_folder
337 dir=trim(restart_folder), mesh = mesh)
339 do ii = c_start, c_end, c_step
340 if (iterate_folder)
then
342 write(folder,
'(a,i0.7,a)') in_folder(folder_index+1:len_trim(in_folder)-1),ii,
"/"
343 write(filename,
'(a,a,a)') trim(folder), trim(basename)
344 out_name = trim(basename)
347 if (c_start /= c_end)
then
351 write(frmt,
'(a,i0,a)')
"(a,i0.",10-len_trim(basename),
")"
352 write(filename, fmt=trim(frmt)) trim(basename), ii
353 write(out_name,
'(a)') trim(filename)
356 write(filename,
'(a,a,a,a)') trim(folder),
"/", trim(basename)
358 write(out_name,
'(a)') trim(basename)
363 call restart%read_mesh_function(trim(filename), mesh, read_ff, ierr)
366 write(
message(1),
'(a,a)')
"Error reading the file ", trim(filename)
367 write(
message(2),
'(a,i4)')
"Error code: ",ierr
368 write(
message(3),
'(a)')
"Skipping...."
372 if (subtract_file)
then
373 read_ff(:) = read_ff(:) - read_rff(:)
374 write(out_name,
'(a,a)') trim(out_name),
"-ref"
377 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
378 if (outp%how(output_i) /= 0)
then
380 trim(out_name), namespace, space, mesh, read_ff,
units_out%length**(-space%dim), ierr, &
381 pos=ions%pos, atoms=ions%atom)
384 if (outp%what(option__output__potential))
then
385 write(out_name,
'(a)')
"potential"
387 call dio_function_output(outp%how(option__output__potential), trim(restart_folder)//trim(folder), &
388 trim(out_name), namespace, space, mesh, pot,
units_out%energy, ierr, pos=ions%pos, atoms=ions%atom)
396 safe_deallocate_a(read_ff)
397 safe_deallocate_a(read_rff)
398 safe_deallocate_a(pot)
405 subroutine convert_transform(mesh, namespace, space, ions, mc, kpoints, basename, in_folder, c_start, c_end, c_step, outp, &
406 subtract_file, ref_name, ref_folder)
407 class(
mesh_t),
intent(in) :: mesh
409 class(
space_t),
intent(in) :: space
410 type(
ions_t),
intent(in) :: ions
413 character(len=*),
intent(inout) :: basename
414 character(len=*),
intent(in) :: in_folder
415 integer,
intent(in) :: c_start
416 integer,
intent(in) :: c_end
417 integer,
intent(in) :: c_step
419 logical,
intent(in) :: subtract_file
420 character(len=*),
intent(inout) :: ref_name
421 character(len=*),
intent(inout) :: ref_folder
423 integer :: ierr, i_space, i_time, nn(1:3), optimize_parity(1:3), wd_info, output_i
424 integer :: i_energy, e_end, e_start, e_point, chunk_size, read_count, t_point
425 logical :: optimize(1:3)
426 integer :: folder_index
427 character(MAX_PATH_LEN) :: filename, folder, restart_folder
428 real(real64) :: fdefault, w_max
429 real(real64),
allocatable :: read_ft(:), read_rff(:), point_tmp(:,:)
431 integer,
parameter :: FAST_FOURIER = 1, standard_fourier = 2
435 type(
batch_t) :: tdrho_b, wdrho_b
436 real(real64),
allocatable :: tdrho_a(:,:,:), wdrho_a(:,:,:)
440 complex(real64),
allocatable :: out_fft(:)
442 real(real64) :: start_time
443 integer :: time_steps
446 real(real64) :: max_energy
447 real(real64) :: min_energy
455 write(
message(1),
'(a)')
'Input: TDTimeStep must be positive.'
459 call io_mkdir(
'wd.general', namespace)
482 call parse_variable(namespace,
'ConvertReadSize', mesh%np, chunk_size)
485 if (chunk_size == 0) chunk_size = mesh%np
487 if (mesh%mpi_grp%size > 1 .and. chunk_size /= mesh%np)
then
488 write(
message(1),*)
'Incompatible value for ConvertReadSize and Parallelizaion in Domains'
489 write(
message(2),*)
'Use the default value for ConvertReadSize (or set it to 0)'
494 if (c_step <= 0)
then
495 write(
message(1),
'(a)')
'Input: ConvertStep must be positive for Fourier transform.'
499 start_time = c_start * dt
501 time_steps = (c_end - c_start) / c_step
502 if (time_steps <= 0)
then
503 write(
message(1),
'(a)')
'Input: ConvertEnd must be larger than ConvertStart for Fourier transform.'
518 if (max_energy > w_max)
then
519 write(
message(1),
'(a,f12.7)')
'Impossible to set ConvertEnergyMax to ', &
521 write(
message(2),
'(a)')
'ConvertEnergyMax is too large.'
522 write(
message(3),
'(a,f12.7,a)')
'ConvertEnergyMax reset to ', &
546 safe_allocate(read_ft(0:time_steps))
548 safe_allocate(read_rff(1:mesh%np))
550 select case (ft_method)
552 nn(1) = time_steps + 1
555 safe_allocate(out_fft(0:time_steps))
560 case (standard_fourier)
570 fdefault =
m_two *
m_pi / (dt * time_steps)
576 spectrum%start_time = c_start * dt
577 spectrum%end_time = c_end * dt
578 spectrum%energy_step = dw
579 spectrum%max_energy = max_energy
580 safe_allocate(tdrho_a(0:time_steps, 1, 1))
581 safe_allocate(wdrho_a(0:time_steps, 1, 1))
587 call kick_init(kick, namespace, space, kpoints, 1)
589 e_start = nint(min_energy / dw)
590 e_end = nint(max_energy / dw)
591 write(
message(1),
'(a,1x,i0.7,a,f12.7,a,i0.7,a,f12.7,a)')
'Frequency index:',e_start,
'(',&
597 if (subtract_file)
then
598 write(
message(1),
'(a,a,a,a)')
"Reading ref-file from ", trim(ref_folder), trim(ref_name),
".obf"
601 dir=trim(ref_folder), mesh = mesh)
604 call restart%read_mesh_function(trim(ref_name), mesh, read_rff, ierr)
607 write(
message(1),
'(2a)')
"Failed to read from ref-file ", trim(ref_name)
608 write(
message(2),
'(2a)')
"from folder ", trim(ref_folder)
614 do i_energy = e_start, e_end
615 write(filename,
'(a14,i0.7,a1)')
'wd.general/wd.',i_energy,
'/'
616 write(
message(1),
'(a,a,f12.7,a,1x,i7,a)')trim(filename),
' w =', &
621 call io_mkdir(trim(filename), namespace)
625 if (mesh%parallel_in_domains)
then
627 folder = in_folder(1:len_trim(in_folder)-1)
628 folder_index = index(folder,
'/', .
true.)
629 restart_folder = folder(1:folder_index)
631 dir=trim(restart_folder), mesh = mesh)
637 safe_allocate(point_tmp(1:chunk_size, 0:time_steps))
640 do i_space = 1, mesh%np
644 do i_time = c_start, c_end, c_step
646 if (mesh%parallel_in_domains .and. i_space == 1)
then
647 write(folder,
'(a,i0.7,a)') in_folder(folder_index+1:len_trim(in_folder)-1),i_time,
"/"
648 write(filename,
'(a,a,a)') trim(folder), trim(basename)
649 call restart%read_mesh_function(trim(filename), mesh, point_tmp(:, t_point), ierr)
653 write(folder,
'(a,i0.7,a)') in_folder(1:len_trim(in_folder)-1),i_time,
"/"
654 write(filename,
'(a,a,a,a)') trim(folder), trim(basename),
".obf"
655 if (mod(i_space-1, chunk_size) == 0)
then
661 if (i_time == c_start) read_count = 0
665 if (ierr /= 0 .and. i_space == 1)
then
666 write(
message(1),
'(a,a,2i10)')
"Error reading the file ", trim(filename), i_space, i_time
667 write(
message(2),
'(a)')
"Skipping...."
668 write(
message(3),
'(a,i0)')
"Error :", ierr
673 if (i_time == c_start) read_count = read_count + 1
674 if (subtract_file)
then
675 read_ft(t_point) = point_tmp(read_count, t_point) - read_rff(i_space)
677 read_ft(t_point) = point_tmp(read_count, t_point)
680 t_point = t_point + 1
683 select case (ft_method)
690 point_tmp(read_count, 0:time_steps) = aimag(out_fft(0:time_steps)) * dt
691 case (standard_fourier)
692 tdrho_a(0:time_steps, 1, 1) = read_ft(0:time_steps)
695 call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, c_start + 1, c_start + time_steps + 1, &
696 kick%time, dt, tdrho_b)
698 c_start + 1, c_start + time_steps + 1, kick%time, dt, tdrho_b, min_energy, max_energy, &
699 spectrum%energy_step, wdrho_b)
702 do e_point = e_start, e_end
703 point_tmp(read_count, e_point) = - wdrho_a(e_point, 1, 1)
707 if (mod(i_space-1, 1000) == 0 .and.
mpi_world%is_root())
then
712 if (mesh%mpi_grp%size == 1)
then
713 if (mod(i_space, chunk_size) == 0)
then
715 write(
message(2),
'(a,i0)')
"Writing binary output: step ", i_space/chunk_size
717 do i_energy = e_start, e_end
718 write(filename,
'(a14,i0.7,a12)')
'wd.general/wd.',i_energy,
'/density.obf'
720 if (i_space == chunk_size)
then
731 if (mesh%mpi_grp%size == 1)
then
732 if (mod(mesh%np, chunk_size) /= 0)
then
733 do i_energy = e_start, e_end
734 write(filename,
'(a14,i0.7,a12)')
'wd.general/wd.',i_energy,
'/density.obf'
735 if (mesh%np < chunk_size)
call dwrite_header(trim(filename), mesh%np_global, ierr)
737 point_tmp(1:mod(mesh%np, chunk_size), i_energy), ierr, nohead=.
true.)
743 call mesh%mpi_grp%barrier()
745 if (mesh%parallel_in_domains)
then
746 do i_energy = e_start, e_end
747 write(filename,
'(a14,i0.7,a1)')
'wd.general/wd.',i_energy,
'/'
748 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
749 if (outp%how(output_i) /= 0)
then
751 trim(
'density'), namespace, space, mesh, point_tmp(:, i_energy), &
752 units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
759 if (any(outp%how /= option__outputformat__binary))
then
760 do i_energy = e_start, e_end
761 write(filename,
'(a14,i0.7,a1)')
'wd.general/wd.',i_energy,
'/'
763 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
764 if ((outp%how(output_i) /= 0) .and. (outp%how(output_i) /= option__outputformat__binary))
then
766 trim(
'density'), namespace, space, mesh, read_rff, &
767 units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
774 safe_deallocate_a(point_tmp)
775 safe_deallocate_a(read_ft)
776 safe_deallocate_a(read_rff)
778 select case (ft_method)
780 safe_deallocate_a(out_fft)
781 case (standard_fourier)
782 safe_deallocate_a(tdrho_a)
783 safe_deallocate_a(wdrho_a)
794 class(
mesh_t),
intent(in) :: mesh
796 class(
space_t),
intent(in) :: space
797 type(
ions_t),
intent(in) :: ions
801 integer :: ierr, ip, i_op, length, n_operations, output_i
805 real(real64) :: f_re, f_im
806 real(real64),
allocatable :: tmp_ff(:), scalar_ff(:)
808 character(len=200) :: var, scalar_expression
809 character(len=MAX_PATH_LEN) :: folder, filename, out_folder, out_filename
825 if (
parse_block(namespace,
'ConvertScalarOperation', blk) == 0)
then
829 if (n_operations == 0)
then
830 write(
message(1),
'(a)')
'No operations found. Check the input file'
841 call parse_variable(namespace,
'ConvertOutputFolder',
"convert", out_folder)
843 call io_mkdir(out_folder, namespace)
853 call parse_variable(namespace,
'ConvertOutputFilename',
'density', out_filename)
855 safe_allocate(tmp_ff(1:mesh%np))
856 safe_allocate(scalar_ff(1:mesh%np))
859 do i_op = 1, n_operations
868 length = len_trim(filename)
870 if (filename(length-3:length) ==
'.obf')
then
871 filename = trim(filename(1:length-4))
877 dir=trim(folder), mesh = mesh, exact=.
true.)
879 call restart%read_mesh_function(trim(filename), mesh, tmp_ff, ierr)
881 write(
message(1),
'(2a)')
"Failed to read mesh function ", trim(filename)
882 write(
message(2),
'(a,i0)')
"Error code: ", ierr
886 write(
message(1),
'(2a)')
"Failed to read from file ", trim(filename)
887 write(
message(2),
'(2a)')
"from folder ", trim(folder)
894 call parse_expression(f_re, f_im, trim(var), real(tmp_ff(ip), real64), trim(scalar_expression))
896 scalar_ff(ip) = scalar_ff(ip) + f_re
905 call mesh%mpi_grp%barrier()
910 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
911 if (outp%how(output_i) /= 0)
then
912 call dio_function_output(outp%how(output_i), trim(out_folder), trim(out_filename), namespace, space, mesh, &
913 scalar_ff, units, ierr, pos=ions%pos, atoms=ions%atom)
917 safe_deallocate_a(tmp_ff)
918 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.