Octopus
convert.F90
Go to the documentation of this file.
1!! Copyright (C) 2013 J. Alberdi-Rodriguez, J. Jornet-Somoza
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
20
21#include "global.h"
22
23program oct_convert
24 use batch_oct_m
27 use debug_oct_m
29 use fft_oct_m
31 use global_oct_m
32 use io_oct_m
35 use ions_oct_m
36 use kick_oct_m
38 use loct_oct_m
40 use mesh_oct_m
41 use mpi_oct_m
44 use output_oct_m
46 use parser_oct_m
50 use space_oct_m
52 use string_oct_m
53 use unit_oct_m
55 use utils_oct_m
56
57 implicit none
58
59 character(len=256) :: config_str
60 integer :: ierr
61
62 call getopt_init(ierr)
63 config_str = trim(get_config_opts()) // trim(get_optional_libraries())
64 if (ierr == 0) call getopt_octopus(config_str)
65 call getopt_end()
66
67 call global_init()
68
69 call parser_init()
70
71 call messages_init()
72
73 call io_init()
75 call messages_experimental("oct-convert utility")
76
77 call print_header()
78
79 call messages_print_with_emphasis(msg="Convert mode", namespace=global_namespace)
81
85
86 call convert()
87
88 call fft_all_end()
90 call io_end()
91 call print_date("Calculation ended on ")
92 call messages_end()
93
94 call parser_end()
95
96 call global_end()
97
98contains
99
100 ! -------------
105 subroutine convert()
106 type(electrons_t), pointer :: sys
107
108 character(MAX_PATH_LEN) :: basename, folder, ref_name, ref_folder, folder_default
109 integer :: c_start, c_end, c_step, c_start_default, length, c_how
110 logical :: iterate_folder, subtract_file
111 integer, parameter :: CONVERT_FORMAT = 1, fourier_transform = 2, operation = 3
112
113 push_sub(convert)
114
115 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
116 sys => electrons_t(global_namespace, generate_epot=.false.)
117 call sys%init_parallelization(mpi_world)
118
119 message(1) = 'Info: Converting files'
120 message(2) = ''
121 call messages_info(2)
122
123 !%Variable ConvertFilename
124 !%Type string
125 !%Default "density"
126 !%Section Utilities::oct-convert
127 !%Description
128 !% Input filename. The original filename which is going to be converted in the format
129 !% specified in <tt>OutputFormat</tt>. It is going to convert various files, it should
130 !% only contain the beginning of the name. For instance, in the case of the restart
131 !% files, it should be one space ' '.
132 !%End
133 call parse_variable(global_namespace, 'ConvertFilename', 'density', basename)
134 if (basename == " ") basename = ""
135 ! Delete the extension if present
136 length = len_trim(basename)
137 if (length > 4) then
138 if (basename(length-3:length) == '.obf') then
139 basename = trim(basename(1:length-4))
140 end if
141 end if
142
143 !%Variable ConvertHow
144 !%Type integer
145 !%Default convert_format
146 !%Section Utilities::oct-convert
147 !%Description
148 !% Select how the mesh function will be converted.
149 !%Option format 1
150 !% The format of the mesh function will be convert from the binary file.obf.
151 !% The format of the output function is set by OutputHow variable.
152 !%Option fourier_transform 2
153 !% A fourier transform of the mesh function will be computed.
154 !% It requieres that ConvertStart and ConvertEnd have to be set.
155 !%Option operation 3
156 !% Convert utility will generate a new mesh function constructed by linear
157 !% combination of scalar function of different mesh functions,
158 !%End
159 call parse_variable(global_namespace, 'ConvertHow', convert_format, c_how)
160
161 !%Variable ConvertIterateFolder
162 !%Type logical
163 !%Default true
164 !%Section Utilities::oct-convert
165 !%Description
166 !% This variable decides if a folder is going to be iterated or the
167 !% filename is going to be iterated.
168 !%End
169 call parse_variable(global_namespace, 'ConvertIterateFolder', .true., iterate_folder)
170
171 if (iterate_folder) then
172 folder_default = 'td.'
173 c_start_default = 0
174 else
175 folder_default = 'restart'
176 c_start_default = 1
177 end if
178
179 !%Variable ConvertFolder
180 !%Type string
181 !%Section Utilities::oct-convert
182 !%Description
183 !% The folder name where the input files are. The default is
184 !% <tt>td.</tt> if <tt>ConvertIterateFolder = true</tt>, otherwise <tt>restart</tt>.
185 !%End
186 call parse_variable(global_namespace, 'ConvertFolder', folder_default, folder)
187 call add_last_slash(folder)
188
189 !%Variable ConvertStart
190 !%Type integer
191 !%Section Utilities::oct-convert
192 !%Description
193 !% The starting number of the filename or folder.
194 !% Default is 0 if <tt>ConvertIterateFolder = true</tt>, otherwise 1.
195 !%End
196 call parse_variable(global_namespace, 'ConvertStart', c_start_default, c_start)
197
198 !%Variable ConvertEnd
199 !%Type integer
200 !%Default 1
201 !%Section Utilities::oct-convert
202 !%Description
203 !% The last number of the filename or folder.
204 !%End
205 call parse_variable(global_namespace, 'ConvertEnd', 1, c_end)
206
207 !%Variable ConvertStep
208 !%Type integer
209 !%Default 1
210 !%Section Utilities::oct-convert
211 !%Description
212 !% The padding between the filenames or folder.
213 !%End
214 call parse_variable(global_namespace, 'ConvertStep', 1, c_step)
215
216 !%Variable ConvertSubtractFilename
217 !%Type string
218 !%Default density
219 !%Section Utilities::oct-convert
220 !%Description
221 !% Input filename. The file which is going to subtracted to rest of the files.
222 !%End
223 call parse_variable(global_namespace, 'ConvertSubtractFilename', 'density', ref_name)
224 if (ref_name == " ") ref_name = ""
225 ! Delete the extension if present
226 length = len_trim(ref_name)
227 if (length > 4) then
228 if (ref_name(length-3:length) == '.obf') then
229 ref_name = trim(ref_name(1:length-4))
230 end if
231 end if
232
233 !%Variable ConvertSubtract
234 !%Type logical
235 !%Default false
236 !%Section Utilities::oct-convert
237 !%Description
238 !% Decides if a reference file is going to be subtracted.
239 !%End
240 call parse_variable(global_namespace, 'ConvertSubtract', .false., subtract_file)
241
242 !%Variable ConvertSubtractFolder
243 !%Type string
244 !%Default .
245 !%Section Utilities::oct-convert
246 !%Description
247 !% The folder name which is going to be subtracted.
248 !%End
249 call parse_variable(global_namespace, 'ConvertSubtractFolder', '.', ref_folder)
250 call add_last_slash(folder)
251
252 select case (c_how)
253 CASE(operation)
254 call convert_operate(sys%gr, global_namespace, sys%space, sys%ions, sys%mc, sys%outp)
255
256 CASE(fourier_transform)
257 ! Compute Fourier transform
258 call convert_transform(sys%gr, global_namespace, sys%space, sys%ions, sys%mc, sys%kpoints, basename, folder, &
259 c_start, c_end, c_step, sys%outp, subtract_file, &
260 ref_name, ref_folder)
261
262 CASE(convert_format)
263 call convert_low(sys%gr, global_namespace, sys%space, sys%ions, sys%hm%psolver, sys%mc, basename, folder, &
264 c_start, c_end, c_step, sys%outp, iterate_folder, &
265 subtract_file, ref_name, ref_folder)
266 end select
267
268 safe_deallocate_p(sys)
269
270 pop_sub(convert)
271 end subroutine convert
272
273 ! ---------------------------------------------------------
276 subroutine convert_low(mesh, namespace, space, ions, psolver, mc, basename, in_folder, c_start, c_end, c_step, outp, &
277 iterate_folder, subtract_file, ref_name, ref_folder)
278 class(mesh_t), intent(in) :: mesh
279 type(namespace_t), intent(in) :: namespace
280 class(space_t), intent(in) :: space
281 type(ions_t), intent(in) :: ions
282 type(poisson_t), intent(in) :: psolver
283 type(multicomm_t), intent(in) :: mc
284 character(len=*), intent(inout) :: basename
285 character(len=*), intent(in) :: in_folder
286 integer, intent(in) :: c_start
287 integer, intent(in) :: c_end
288 integer, intent(in) :: c_step
289 type(output_t), intent(in) :: outp
290 logical, intent(in) :: iterate_folder
292 logical, intent(in) :: subtract_file
293 character(len=*), intent(inout) :: ref_name
294 character(len=*), intent(inout) :: ref_folder
295
296 type(restart_t) :: restart
297 integer :: ierr, ii, folder_index, output_i
298 character(MAX_PATH_LEN) :: filename, out_name, folder, frmt, restart_folder
299 real(real64), allocatable :: read_ff(:), read_rff(:), pot(:)
300
301 push_sub(convert_low)
302
303 safe_allocate(read_ff(1:mesh%np))
304 safe_allocate(read_rff(1:mesh%np))
305 safe_allocate(pot(1:mesh%np))
306 read_rff(:) = m_zero
307
308 write(message(1),'(5a,i5,a,i5,a,i5)') "Converting '", trim(in_folder), "/", trim(basename), &
309 "' from ", c_start, " to ", c_end, " every ", c_step
310 call messages_info(1)
311
312 if (subtract_file) then
313 write(message(1),'(a,a,a,a)') "Reading ref-file from ", trim(ref_folder), trim(ref_name),".obf"
314 call restart_init(restart, namespace, restart_undefined, restart_type_load, mc, ierr, &
315 dir=trim(ref_folder), mesh = mesh)
316 ! FIXME: why only real functions? Please generalize.
317 if (ierr == 0) then
318 call restart_read_mesh_function(restart, space, trim(ref_name), mesh, read_rff, ierr)
319 call restart_end(restart)
320 else
321 write(message(1),'(2a)') "Failed to read from ref-file ", trim(ref_name)
322 write(message(2), '(2a)') "from folder ", trim(ref_folder)
323 call messages_fatal(2)
324 end if
325 end if
326
327 ! Initialize the restart directory from <tt>ConvertFolder</tt> value.
328 ! This directory has to have the files 'grid' and 'lxyz.obf'
329 ! and the files that are going to be converged, must be inside this folder
330 if (iterate_folder) then
331 ! Delete the last / and find the previous /, if any
332 folder = in_folder(1:len_trim(in_folder)-1)
333 folder_index = index(folder, '/', .true.)
334 restart_folder = folder(1:folder_index)
335 else
336 restart_folder = in_folder
337 end if
338 call restart_init(restart, namespace, restart_undefined, restart_type_load, mc, ierr, &
339 dir=trim(restart_folder), mesh = mesh)
340 call loct_progress_bar(-1, c_end-c_start)
341 do ii = c_start, c_end, c_step
342 if (iterate_folder) then
343 ! Delete the last / and add the corresponding folder number
344 write(folder,'(a,i0.7,a)') in_folder(folder_index+1:len_trim(in_folder)-1),ii,"/"
345 write(filename, '(a,a,a)') trim(folder), trim(basename)
346 out_name = trim(basename)
347 else
348 folder = ""
349 if (c_start /= c_end) then
350 ! Here, we are only considering 10 character long filenames.
351 ! Subtract the initial part given at 'ConvertFilename' from the format and pad
352 ! with zeros.
353 write(frmt,'(a,i0,a)')"(a,i0.",10-len_trim(basename),")"
354 write(filename, fmt=trim(frmt)) trim(basename), ii
355 write(out_name, '(a)') trim(filename)
356 else
357 ! Assuming filename is given complete in the 'ConvertFilename'
358 write(filename, '(a,a,a,a)') trim(folder),"/", trim(basename)
359 filename = basename
360 write(out_name, '(a)') trim(basename)
361 end if
362 end if
363
364 ! Read the obf file
365 if (ierr == 0) then
366 call restart_read_mesh_function(restart, space, trim(filename), mesh, read_ff, ierr)
367 end if
368
369 if (ierr /= 0) then
370 write(message(1), '(a,a)') "Error reading the file ", filename
371 write(message(2), '(a,i4)') "Error code: ",ierr
372 write(message(3), '(a)') "Skipping...."
373 call messages_warning(3)
374 cycle
375 end if
376 if (subtract_file) then
377 read_ff(:) = read_ff(:) - read_rff(:)
378 write(out_name, '(a,a)') trim(out_name),"-ref"
379 end if
380 ! Write the corresponding output
381 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
382 if (outp%how(output_i) /= 0) then
383 call dio_function_output(outp%how(output_i), trim(restart_folder)//trim(folder), &
384 trim(out_name), namespace, space, mesh, read_ff, units_out%length**(-space%dim), ierr, &
385 pos=ions%pos, atoms=ions%atom)
386 end if
387 end do
388 if (outp%what(option__output__potential)) then
389 write(out_name, '(a)') "potential"
390 call dpoisson_solve(psolver, namespace, pot, read_ff)
391 call dio_function_output(outp%how(option__output__potential), trim(restart_folder)//trim(folder), &
392 trim(out_name), namespace, space, mesh, pot, units_out%energy, ierr, pos=ions%pos, atoms=ions%atom)
393 end if
394 call loct_progress_bar(ii-c_start, c_end-c_start)
395 ! It does not matter if the current write has failed for the next iteration
396 ierr = 0
397 end do
398 call restart_end(restart)
399
400 safe_deallocate_a(read_ff)
401 safe_deallocate_a(read_rff)
402 safe_deallocate_a(pot)
403 pop_sub(convert_low)
404 end subroutine convert_low
405
406 ! ---------------------------------------------------------
409 subroutine convert_transform(mesh, namespace, space, ions, mc, kpoints, basename, in_folder, c_start, c_end, c_step, outp, &
410 subtract_file, ref_name, ref_folder)
411 class(mesh_t), intent(in) :: mesh
412 type(namespace_t), intent(in) :: namespace
413 class(space_t), intent(in) :: space
414 type(ions_t), intent(in) :: ions
415 type(multicomm_t), intent(in) :: mc
416 type(kpoints_t), intent(in) :: kpoints
417 character(len=*), intent(inout) :: basename
418 character(len=*), intent(in) :: in_folder
419 integer, intent(in) :: c_start
420 integer, intent(in) :: c_end
421 integer, intent(in) :: c_step
422 type(output_t), intent(in) :: outp
423 logical, intent(in) :: subtract_file
424 character(len=*), intent(inout) :: ref_name
425 character(len=*), intent(inout) :: ref_folder
426
427 integer :: ierr, i_space, i_time, nn(1:3), optimize_parity(1:3), wd_info, output_i
428 integer :: i_energy, e_end, e_start, e_point, chunk_size, read_count, t_point
429 logical :: optimize(1:3)
430 integer :: folder_index
431 character(MAX_PATH_LEN) :: filename, folder, restart_folder
432 real(real64) :: fdefault, w_max
433 real(real64), allocatable :: read_ft(:), read_rff(:), point_tmp(:,:)
434
435 integer, parameter :: FAST_FOURIER = 1, standard_fourier = 2
436 type(kick_t) :: kick
437 integer :: ft_method
438 type(spectrum_t) :: spectrum
439 type(batch_t) :: tdrho_b, wdrho_b
440 real(real64), allocatable :: tdrho_a(:,:,:), wdrho_a(:,:,:)
441 type(fft_t) :: fft
442
443 type(restart_t) :: restart
444 complex(real64), allocatable :: out_fft(:)
445
446 real(real64) :: start_time
447 integer :: time_steps
448 real(real64) :: dt
449 real(real64) :: dw
450 real(real64) :: max_energy
451 real(real64) :: min_energy
452
453 push_sub(convert_transform)
454
455 ! set default time_step as dt from TD section
456 fdefault = m_zero
457 call parse_variable(namespace, 'TDTimeStep', fdefault, dt, unit = units_inp%time)
458 if (dt <= m_zero) then
459 write(message(1),'(a)') 'Input: TDTimeStep must be positive.'
460 write(message(2),'(a)') 'Input: TDTimeStep reset to 0. Check input file'
461 call messages_info(2)
462 end if
463
464 call io_mkdir('wd.general', namespace)
465 wd_info = io_open('wd.general/wd.info', global_namespace, action='write')
466 call messages_print_with_emphasis(msg="Fourier Transform Options", iunit=wd_info)
467
468 !%Variable ConvertEnergyMin
469 !%Type float
470 !%Default 0.0
471 !%Section Utilities::oct-convert
472 !%Description
473 !% Minimum energy to output from Fourier transform.
474 !%End
475 call parse_variable(namespace, 'ConvertEnergyMin', m_zero, min_energy, units_inp%energy)
476 call messages_print_var_value('ConvertEnergyMin', min_energy, unit = units_out%energy, iunit=wd_info)
477
478 !%Variable ConvertReadSize
479 !%Type integer
480 !%Default mesh%np
481 !%Section Utilities::oct-convert
482 !%Description
483 !% How many points are read at once. For the parallel run this has not been
484 !% yet tested, so it should be one. For the serial run, a number
485 !% of 100-1000 will speed-up the execution time by this factor.
486 !%End
487 call parse_variable(namespace, 'ConvertReadSize', mesh%np, chunk_size)
488 call messages_print_var_value('ConvertReadSize', chunk_size, iunit=wd_info)
489 !Check that default value is set when ConvertReadSize = 0
490 if (chunk_size == 0) chunk_size = mesh%np
491 ! Parallel version just work in domains and chunk_size equal to mesh%np
492 if (mesh%mpi_grp%size > 1 .and. chunk_size /= mesh%np) then
493 write(message(1),*)'Incompatible value for ConvertReadSize and Parallelizaion in Domains'
494 write(message(2),*)'Use the default value for ConvertReadSize (or set it to 0)'
495 call messages_fatal(2)
496 end if
497
498 ! Calculate the limits in frequency space.
499 start_time = c_start * dt
500 dt = dt * c_step
501 time_steps = (c_end - c_start) / c_step
502 dw = m_two * m_pi / (dt * time_steps)
503 w_max = m_two * m_pi / dt
504
505 !%Variable ConvertEnergyMax
506 !%Type float
507 !%Default w_max
508 !%Section Utilities::oct-convert
509 !%Description
510 !% Maximum energy to output from Fourier transform.
511 !%End
512 fdefault = units_from_atomic(units_inp%energy, w_max)
513 call parse_variable(namespace, 'ConvertEnergyMax',fdefault, max_energy, units_inp%energy)
514 if (max_energy > w_max) then
515 write(message(1),'(a,f12.7)')'Impossible to set ConvertEnergyMax to ', &
516 units_from_atomic(units_inp%energy, max_energy)
517 write(message(2),'(a)')'ConvertEnergyMax is too large.'
518 write(message(3),'(a,f12.7,a)')'ConvertEnergyMax reset to ', &
519 units_from_atomic(units_inp%energy, w_max),'[' // trim(units_abbrev(units_out%energy)) // ']'
520 call messages_info(3)
521 max_energy = w_max
522 end if
523 call messages_print_var_value('ConvertEnergyMax', max_energy, unit = units_out%energy, iunit=wd_info)
524
525 !%Variable ConvertFTMethod
526 !%Type integer
527 !%Default FAST_FOURIER
528 !%Section Utilities::oct-convert
529 !%Description
530 !% Describes the method used to perform the Fourier Transform
531 !%Option fast_fourier 1
532 !% Uses Fast Fourier Transform as implemented in the external library.
533 !%Option standard_fourier 2
534 !% Uses polinomial approach to the computation of discrete Fourier Transform.
535 !% It uses the same variable described in how to obtain spectrum from
536 !% a time-propagation calculation.
537 !%End
538 call parse_variable(namespace, 'ConvertFTMethod', 1, ft_method)
539 call messages_print_var_option('ConvertFTMethod', ft_method, iunit=wd_info)
540
541 !TODO: check if e_point can be used instead of e_point+1
542 safe_allocate(read_ft(0:time_steps))
543 read_ft = m_zero
544 safe_allocate(read_rff(1:mesh%np))
545
546 select case (ft_method)
547 case (fast_fourier)
548 nn(1) = time_steps + 1
549 nn(2) = 1
550 nn(3) = 1
551 safe_allocate(out_fft(0:time_steps))
552 optimize = .false.
553 optimize_parity = -1
554 call fft_init(fft, nn, 1, fft_real, fftlib_fftw, optimize, optimize_parity)
555 case (standard_fourier)
556 !%Variable ConvertEnergyStep
557 !%Type float
558 !%Default <math>2 \pi / T</math>, where <math>T</math> is the total propagation time
559 !%Section Utilities::oct-convert
560 !%Description
561 !% Energy step to output from Fourier transform.
562 !% Sampling rate for the Fourier transform. If you supply a number equal or smaller than zero, then
563 !% the sampling rate will be <math>2 \pi / T</math>, where <math>T</math> is the total propagation time.
564 !%End
565 fdefault = m_two * m_pi / (dt * time_steps)
566 call parse_variable(namespace, 'ConvertEnergyStep',fdefault, dw, units_inp%energy)
567 if (dw <= m_zero) dw = m_two * m_pi / (dt * time_steps)
568
569 call spectrum_init(spectrum, namespace, dw, w_max)
570 ! Manually setting already defined variables on spectrum.
571 spectrum%start_time = c_start * dt
572 spectrum%end_time = c_end * dt
573 spectrum%energy_step = dw
574 spectrum%max_energy = max_energy
575 safe_allocate(tdrho_a(0:time_steps, 1, 1))
576 safe_allocate(wdrho_a(0:time_steps, 1, 1))
577 end select
578 call messages_print_var_value('ConvertEnergyStep', dw, unit = units_out%energy, iunit=wd_info)
579
580 !TODO: set system variable common for all the program in
581 ! order to use call kick_init(kick, sy%st%d%nspin, sys%space%dim, sys%ions%periodic_dim)
582 call kick_init(kick, namespace, space, kpoints, 1)
583
584 e_start = nint(min_energy / dw)
585 e_end = nint(max_energy / dw)
586 write(message(1),'(a,1x,i0.7,a,f12.7,a,i0.7,a,f12.7,a)')'Frequency index:',e_start,'(',&
587 units_from_atomic(units_out%energy, e_start * dw),')-',e_end,'(',units_from_atomic(units_out%energy, e_end * dw),')'
588 write(message(2),'(a,f12.7,a)')'Frequency Step, dw: ', units_from_atomic(units_out%energy, dw), &
589 '[' // trim(units_abbrev(units_out%energy)) // ']'
590 call messages_info(2)
591
592 if (subtract_file) then
593 write(message(1),'(a,a,a,a)') "Reading ref-file from ", trim(ref_folder), trim(ref_name),".obf"
594 call messages_info(1)
595 call restart_init(restart, namespace, restart_undefined, restart_type_load, mc, ierr, &
596 dir=trim(ref_folder), mesh = mesh)
597 ! FIXME: why only real functions? Please generalize.
598 if (ierr == 0) then
599 call restart_read_mesh_function(restart, space, trim(ref_name), mesh, read_rff, ierr)
600 call restart_end(restart)
601 else
602 write(message(1),'(2a)') "Failed to read from ref-file ", trim(ref_name)
603 write(message(2), '(2a)') "from folder ", trim(ref_folder)
604 call messages_fatal(2)
605 end if
606 end if
607
608 call messages_print_with_emphasis(msg="File Information", iunit=wd_info)
609 do i_energy = e_start, e_end
610 write(filename,'(a14,i0.7,a1)')'wd.general/wd.',i_energy,'/'
611 write(message(1),'(a,a,f12.7,a,1x,i7,a)')trim(filename),' w =', &
612 units_from_atomic(units_out%energy,(i_energy) * dw), &
613 '[' // trim(units_abbrev(units_out%energy)) // ']'
614 if (mpi_world%rank == 0) write(wd_info,'(a)') message(1)
615 call messages_info(1)
616 call io_mkdir(trim(filename), namespace)
617 end do
618 call io_close(wd_info)
619
620 if (mesh%parallel_in_domains) then
621 ! Delete the last / and find the previous /, if any
622 folder = in_folder(1:len_trim(in_folder)-1)
623 folder_index = index(folder, '/', .true.)
624 restart_folder = folder(1:folder_index)
625 call restart_init(restart, namespace, restart_undefined, restart_type_load, mc, ierr, &
626 dir=trim(restart_folder), mesh = mesh)
627 end if
628
629 !For each mesh point, open density file and read corresponding point.
630 if (mpi_world%rank == 0) call loct_progress_bar(-1, mesh%np)
631
632 safe_allocate(point_tmp(1:chunk_size, 0:time_steps))
633 read_count = 0
634 ! Space
635 do i_space = 1, mesh%np
636 ! Time
637 t_point = 0
638 do i_time = c_start, c_end, c_step
639
640 if (mesh%parallel_in_domains .and. i_space == 1) then
641 write(folder,'(a,i0.7,a)') in_folder(folder_index+1:len_trim(in_folder)-1),i_time,"/"
642 write(filename, '(a,a,a)') trim(folder), trim(basename)
643 call restart_read_mesh_function(restart, space, trim(filename), mesh, point_tmp(:, t_point), ierr)
644 else
645 ! Here, we always iterate folders
646 ! Delete the last / and add the corresponding folder number
647 write(folder,'(a,i0.7,a)') in_folder(1:len_trim(in_folder)-1),i_time,"/"
648 write(filename, '(a,a,a,a)') trim(folder), trim(basename), ".obf"
649 if (mod(i_space-1, chunk_size) == 0) then
650 call profiling_in("READING")
651 !TODO: check for any error on the whole file before reading by parts.
652 call io_binary_read(trim(filename), i4_to_i8(chunk_size), point_tmp(1:chunk_size, t_point), &
653 ierr, offset=i4_to_i8(i_space-1))
654 call profiling_out("READING")
655 if (i_time == c_start) read_count = 0
656 end if
657 end if
658
659 if (ierr /= 0 .and. i_space == 1) then
660 write(message(1), '(a,a,2i10)') "Error reading the file ", trim(filename), i_space, i_time
661 write(message(2), '(a)') "Skipping...."
662 write(message(3), '(a,i0)') "Error :", ierr
663 call messages_warning(3)
664 cycle
665 end if
666
667 if (i_time == c_start) read_count = read_count + 1
668 if (subtract_file) then
669 read_ft(t_point) = point_tmp(read_count, t_point) - read_rff(i_space)
670 else
671 read_ft(t_point) = point_tmp(read_count, t_point)
672 end if
673
674 t_point = t_point + 1
675 end do ! Time
676
677 select case (ft_method)
678 case (fast_fourier)
679 call profiling_in("CONVERT_FFTW")
680 call dfft_forward(fft, read_ft, out_fft)
681 call profiling_out("CONVERT_FFTW")
682 ! Should the value be multiplied by dt ??? as in standard discrete Fourier Transform ?
683 point_tmp(read_count, 0:time_steps) = aimag(out_fft(0:time_steps)) * dt
684 case (standard_fourier)
685 tdrho_a(0:time_steps, 1, 1) = read_ft(0:time_steps)
686 call batch_init(tdrho_b, 1, 1, 1, tdrho_a)
687 call batch_init(wdrho_b, 1, 1, 1, wdrho_a)
688 call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, c_start + 1, c_start + time_steps + 1, &
689 kick%time, dt, tdrho_b)
690 call spectrum_fourier_transform(spectrum%method, spectrum%transform, spectrum%noise, &
691 c_start + 1, c_start + time_steps + 1, kick%time, dt, tdrho_b, min_energy, max_energy, &
692 spectrum%energy_step, wdrho_b)
693 call tdrho_b%end()
694 call wdrho_b%end()
695 do e_point = e_start, e_end
696 point_tmp(read_count, e_point) = - wdrho_a(e_point, 1, 1)
697 end do
698 end select
699
700 if (mod(i_space-1, 1000) == 0 .and. mpi_world%rank == 0) then
701 call loct_progress_bar(i_space-1, mesh%np)
702 end if
703
704 !print out wd densities from (ii-chunksize,ii] if running in serial
705 if (mesh%mpi_grp%size == 1) then
706 if (mod(i_space, chunk_size) == 0) then
707 write(message(1),'(a)') ""
708 write(message(2),'(a,i0)') "Writing binary output: step ", i_space/chunk_size
709 call messages_info(2)
710 do i_energy = e_start, e_end
711 write(filename,'(a14,i0.7,a12)')'wd.general/wd.',i_energy,'/density.obf'
712 ! If it is the first time entering here, write the header. But, only once
713 if (i_space == chunk_size) then
714 call dwrite_header(trim(filename), mesh%np_global, ierr)
715 end if
716 call io_binary_write(trim(filename), i4_to_i8(chunk_size), point_tmp(1:chunk_size, i_energy), ierr, &
717 nohead = .true.)
718 end do
719 end if
720 end if
721 end do ! Space
722
723 if (mpi_world%rank == 0) call loct_progress_bar(mesh%np, mesh%np)
724 call mesh%mpi_grp%barrier()
725
726 if (mesh%parallel_in_domains) then
727 do i_energy = e_start, e_end
728 write(filename,'(a14,i0.7,a1)')'wd.general/wd.',i_energy,'/'
729 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
730 if (outp%how(output_i) /= 0) then
731 call dio_function_output(0_8, trim(filename), &
732 trim('density'), namespace, space, mesh, point_tmp(:, i_energy), &
733 units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
734 end if
735 end do
736 end do
737 call restart_end(restart)
738 else
739 ! write the output files
740 if (any(outp%how /= option__outputformat__binary)) then
741 do i_energy = e_start, e_end
742 write(filename,'(a14,i0.7,a1)')'wd.general/wd.',i_energy,'/'
743 call io_binary_read(trim(filename)//'density.obf', i4_to_i8(mesh%np), read_rff, ierr)
744 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
745 if ((outp%how(output_i) /= 0) .and. (outp%how(output_i) /= option__outputformat__binary)) then
746 call dio_function_output(outp%how(output_i), trim(filename), &
747 trim('density'), namespace, space, mesh, read_rff, &
748 units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
749 end if
750 end do
751 end do
752 end if
753 end if
754
755 safe_deallocate_a(point_tmp)
756 safe_deallocate_a(read_ft)
757 safe_deallocate_a(read_rff)
758
759 select case (ft_method)
760 case (fast_fourier)
761 safe_deallocate_a(out_fft)
762 case (standard_fourier)
763 call kick_end(kick)
764 safe_deallocate_a(tdrho_a)
765 safe_deallocate_a(wdrho_a)
766 end select
767
768 pop_sub(convert_transform)
769 end subroutine convert_transform
770 ! ---------------------------------------------------------
773 subroutine convert_operate(mesh, namespace, space, ions, mc, outp)
774 class(mesh_t), intent(in) :: mesh
775 type(namespace_t), intent(in) :: namespace
776 class(space_t), intent(in) :: space
777 type(ions_t), intent(in) :: ions
778 type(multicomm_t), intent(in) :: mc
779 type(output_t) , intent(in) :: outp
780
781 integer :: ierr, ip, i_op, length, n_operations, output_i
782 type(block_t) :: blk
783 type(restart_t) :: restart
784 type(unit_t) :: units
785 real(real64) :: f_re, f_im
786 real(real64), allocatable :: tmp_ff(:), scalar_ff(:)
787
788 character(len=200) :: var, scalar_expression
789 character(len=MAX_PATH_LEN) :: folder, filename, out_folder, out_filename
790
791 push_sub(convert_operate)
792
793 !%Variable ConvertScalarOperation
794 !%Type block
795 !%Section Utilities::oct-convert
796 !%Description
797 !% This variable is used to generate a new mesh function as a linear combination
798 !% different mesh function having the same mesh. Each row defines an operation for
799 !% for a single mesh function.
800 !% The format of the block is the following: <br>
801 !% 'variable name' | 'folder' | 'file' | 'operation'
802 !%End
803 ! First, find out if there is a ConvertScalarOperation block.
804 n_operations = 0
805 if (parse_block(namespace, 'ConvertScalarOperation', blk) == 0) then
806 n_operations = parse_block_n(blk)
807 end if
808
809 if (n_operations == 0) then
810 write(message(1),'(a)')'No operations found. Check the input file'
811 call messages_fatal(1)
812 end if
813
814 !%Variable ConvertOutputFolder
815 !%Type string
816 !%Section Utilities::oct-convert
817 !%Description
818 !% The folder name where the output files will be write. The default is
819 !% <tt>convert</tt>.
820 !%End
821 call parse_variable(namespace, 'ConvertOutputFolder', "convert", out_folder)
822 call add_last_slash(out_folder)
823 call io_mkdir(out_folder, namespace)
824
825 !%Variable ConvertOutputFilename
826 !%Type string
827 !%Default "density"
828 !%Section Utilities::oct-convert
829 !%Description
830 !% Output filename. The name of the file in which the converted mesh function will be
831 !% written in the format specified in <tt>OutputFormat</tt>.
832 !%End
833 call parse_variable(namespace, 'ConvertOutputFilename', 'density', out_filename)
834
835 safe_allocate(tmp_ff(1:mesh%np))
836 safe_allocate(scalar_ff(1:mesh%np))
837 scalar_ff = m_zero
838
839 do i_op = 1, n_operations
840 !read variable name
841 call parse_block_string(blk, i_op-1, 0, var)
842 !read folder path
843 call parse_block_string(blk, i_op-1, 1, folder)
844 call add_last_slash(folder)
845 !read file
846 call parse_block_string(blk, i_op-1, 2, filename)
847 ! Delete the extension if present
848 length = len_trim(filename)
849 if (length > 4) then
850 if (filename(length-3:length) == '.obf') then
851 filename = trim(filename(1:length-4))
852 end if
853 end if
854 ! FIXME: why only real functions? Please generalize.
855 ! TODO: check if mesh function are real or complex.
856 call restart_init(restart, namespace, restart_undefined, restart_type_load, mc, ierr, &
857 dir=trim(folder), mesh = mesh, exact=.true.)
858 if (ierr == 0) then
859 call restart_read_mesh_function(restart, space, trim(filename), mesh, tmp_ff, ierr)
860 else
861 write(message(1),'(2a)') "Failed to read from file ", trim(filename)
862 write(message(2), '(2a)') "from folder ", trim(folder)
863 call messages_fatal(2)
864 end if
865 !read scalar expression
866 call parse_block_string(blk, i_op-1, 3, scalar_expression)
867
868 do ip = 1, mesh%np
869 call parse_expression(f_re, f_im, trim(var), real(tmp_ff(ip), real64), trim(scalar_expression))
870 !TODO: implement use of complex functions.
871 scalar_ff(ip) = scalar_ff(ip) + f_re
872 end do
873
874 call restart_end(restart)
875
876 end do
877
878 call parse_block_end(blk)
879
880 call mesh%mpi_grp%barrier()
881 ! Write the corresponding output
882 !TODO: add variable ConvertFunctionType to select the type(density, wfs, potential, ...)
883 ! and units of the conversion.
884 units = units_out%length**(-space%dim)
885 do output_i = lbound(outp%how, 1), ubound(outp%how, 1)
886 if (outp%how(output_i) /= 0) then
887 call dio_function_output(outp%how(output_i), trim(out_folder), trim(out_filename), namespace, space, mesh, &
888 scalar_ff, units, ierr, pos=ions%pos, atoms=ions%atom)
889 end if
890 end do
891
892 safe_deallocate_a(tmp_ff)
893 safe_deallocate_a(scalar_ff)
894
895 pop_sub(convert_operate)
896 end subroutine convert_operate
897end program
898
899!! Local Variables:
900!! mode: f90
901!! coding: utf-8
902!! End:
program oct_convert
This utility runs in parallel and can be used for post-processing of the results of Output.
Definition: convert.F90:116
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.
Definition: convert.F90:371
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...
Definition: convert.F90:867
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.
Definition: convert.F90:504
initialize a batch with existing memory
Definition: batch.F90:267
Each program/utility that needs to use the getopt features should have an interface here – the defini...
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
static void convert(multi *in, multi *out, int t_in, int t_out)
Definition: io_binary.c:5131
This module implements batches of mesh functions.
Definition: batch.F90:133
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_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 ...
Definition: fft.F90:118
subroutine, public fft_init(this, nn, dim, type, library, optimize, optimize_parity, comm, mpi_grp, use_aligned)
Definition: fft.F90:417
subroutine, public fft_all_init(namespace)
initialize the table
Definition: fft.F90:277
subroutine, public fft_all_end()
delete all plans
Definition: fft.F90:390
integer, parameter, public fft_real
Definition: fft.F90:177
integer, parameter, public fftlib_fftw
Definition: fft.F90:182
real(real64), parameter, public m_two
Definition: global.F90:189
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:381
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:324
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.
Definition: io.F90:114
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...
Definition: io.F90:166
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
subroutine, public io_end()
Definition: io.F90:265
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
subroutine, public kick_end(kick)
Definition: kick.F90:795
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
Definition: kick.F90:223
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_end()
Definition: messages.F90:277
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_init(output_dir)
Definition: messages.F90:224
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public print_date(str)
Definition: messages.F90:1017
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
type(namespace_t), public global_namespace
Definition: namespace.F90:132
this module contains the output system
Definition: output_low.F90:115
this module contains the output system
Definition: output.F90:115
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:450
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:481
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:892
subroutine, public profiling_end(namespace)
Definition: profiling.F90:413
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:255
subroutine, public restart_module_init(namespace)
Definition: restart.F90:306
integer, parameter, public restart_undefined
Definition: restart.F90:229
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:514
integer, parameter, public restart_type_load
Definition: restart.F90:225
subroutine, public restart_end(restart)
Definition: restart.F90:720
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...
Definition: spectrum.F90:2638
subroutine, public spectrum_init(spectrum, namespace, default_energy_step, default_max_energy)
Definition: spectrum.F90:213
integer default
Definition: spectrum.F90:207
subroutine, public spectrum_signal_damp(damp_type, damp_factor, time_start, time_end, t0, time_step, time_function)
Definition: spectrum.F90:2552
subroutine, public add_last_slash(str)
Adds a '/' in the end of the string, only if it missing. Useful for directories.
Definition: string.F90:161
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
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.
Definition: utils.F90:118
character(len=256) function, public get_config_opts()
Definition: utils.F90:360
character(len=256) function, public get_optional_libraries()
Definition: utils.F90:383
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Definition: utils.F90:301
Class defining batches of mesh functions.
Definition: batch.F90:159
Class describing the electron system.
Definition: electrons.F90:214
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Stores all communicators and groups.
Definition: multicomm.F90:206
output handler class
Definition: output_low.F90:163
int true(void)