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