Octopus
vibrational.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
19#include "global.h"
20
21program vibrational
22 use batch_oct_m
24 use debug_oct_m
25 use global_oct_m
26 use io_oct_m
27 use ions_oct_m
28 use, intrinsic :: iso_fortran_env
30 use mpi_oct_m
32 use parser_oct_m
35 use unit_oct_m
37
38 implicit none
39
40 integer :: iunit, ierr, ii, jj, iter, read_iter, ntime, nvaf, nvel, ivel
41 real(real64), allocatable :: vaf(:), time(:), velocities(:, :), ftvaf(:)
42 type(ions_t), pointer :: ions
43 type(spectrum_t) :: spectrum
44 type(batch_t) :: vafb, ftvafb
45 real(real64) :: ww, curtime, vaftime, deltat
46 integer :: ifreq, max_freq
47 integer :: skip
48
49 ! Initialize stuff
51
52 call getopt_init(ierr)
53 call getopt_end()
54
55 call parser_init()
56
57 call messages_init()
58
59 call io_init()
61
63
64 call spectrum_init(spectrum, global_namespace, &
65 default_energy_step = units_to_atomic(unit_invcm, 0.2_real64), &
66 default_max_energy = units_to_atomic(unit_invcm, 5000.0_real64))
67
68 !%Variable VibrationalSpectrumTimeStepFactor
69 !%Type integer
70 !%Default 10
71 !%Section Utilities::oct-vibrational_spectrum
72 !%Description
73 !% In the calculation of the vibrational spectrum, it is not necessary
74 !% to read the velocity at every time step. This variable controls
75 !% the integer factor between the simulation time step and the
76 !% time step used to calculate the vibrational spectrum.
77 !%End
78
79 call messages_obsolete_variable(global_namespace, 'PropagationSpectrumTimeStepFactor', 'VibrationalSpectrumTimeStepFactor')
80 call parse_variable(global_namespace, 'VibrationalSpectrumTimeStepFactor', 10, skip)
81 if (skip <= 0) call messages_input_error(global_namespace, 'VibrationalSpectrumTimeStepFactor')
82
83 max_freq = spectrum_nenergy_steps(spectrum)
84
85 if (spectrum%end_time < m_zero) spectrum%end_time = huge(spectrum%end_time)
86
88
89 ! Opens the coordinates files.
90 iunit = io_open('td.general/coordinates', global_namespace, action='read')
91
92 call io_skip_header(iunit)
93
94 ntime = 1
95 iter = 1
96
97 ! check the number of time steps we will read
98 do
99 read(unit = iunit, iostat = ierr, fmt = *) read_iter, curtime, &
100 ((ions%pos(jj, ii), jj = 1, 3), ii = 1, ions%natoms), &
101 ((ions%vel(jj, ii), jj = 1, 3), ii = 1, ions%natoms)
102
103 curtime = units_to_atomic(units_out%time, curtime)
104
105 if (ierr /= 0 .or. curtime >= spectrum%end_time) then
106 iter = iter - 1 ! last iteration is not valid
107 ntime = ntime - 1
108 exit
109 end if
110
111 if (iter /= read_iter + 1) then
112 call messages_write("Error while reading file 'td.general/coordinates',", new_line = .true.)
113 call messages_write('expected iteration ')
114 call messages_write(iter - 1)
115 call messages_write(', got iteration ')
116 call messages_write(read_iter)
117 call messages_write('.')
118 call messages_fatal()
119 end if
120
121 ! ntime counts how many steps are gonna be used
122 if (curtime >= spectrum%start_time .and. mod(iter, skip) == 0) ntime = ntime + 1
123
124 iter = iter + 1 !counts number of timesteps (with time larger than zero up to SpecEndTime)
125 end do
126
127 call io_close(iunit)
128
129 nvel = ions%natoms*ions%space%dim
130
131 safe_allocate(time(1:ntime))
132 safe_allocate(velocities(1:nvel, 1:ntime))
133
134 ! Opens the coordinates files.
135 iunit = io_open('td.general/coordinates', global_namespace, action='read')
136
137 call io_skip_header(iunit)
138
139 ntime = 1
140 iter = 1
141
142 do
143 read(unit = iunit, iostat = ierr, fmt = *) read_iter, curtime, &
144 ((ions%pos(jj, ii), jj = 1, 3), ii = 1, ions%natoms), &
145 ((ions%vel(jj, ii), jj = 1, 3), ii = 1, ions%natoms)
146
147 curtime = units_to_atomic(units_out%time, curtime)
148
149 if (ierr /= 0 .or. curtime >= spectrum%end_time) then
150 iter = iter - 1 ! last iteration is not valid
151 ntime = ntime - 1
152 exit
153 end if
154
155 assert(iter == read_iter + 1)
156
157 if (curtime >= spectrum%start_time .and. mod(iter, skip) == 0) then
158
159 time(ntime) = curtime
160 ivel = 1
161 do ii = 1, ions%natoms
162 do jj = 1, ions%space%dim
163 velocities(ivel, ntime) = units_to_atomic(units_out%velocity, ions%vel(jj, ii))
164 ivel = ivel + 1
165 end do
166 end do
167
168 ntime = ntime + 1
169 end if
170
171 iter = iter + 1 !counts number of timesteps (with time larger than zero up to SpecEndTime)
172 end do
173
174 call io_close(iunit)
175
176 deltat = time(2) - time(1)
177
178 !%Variable VibrationalSpectrumTime
179 !%Type integer
180 !%Section Utilities::oct-vibrational_spectrum
181 !%Description
182 !% This variable controls the maximum time for the calculation of
183 !% the velocity autocorrelation function. The default is the total
184 !% propagation time.
185 !%End
186 call parse_variable(global_namespace, 'VibrationalSpectrumTime', ntime*deltat, vaftime)
187
188 nvaf = int(vaftime/deltat)
189
190 safe_allocate(vaf(1:ntime))
191
192 call messages_write('Time step = ')
193 call messages_write(deltat, units = units_out%time)
194 call messages_info()
195
196 call messages_new_line()
197 call messages_write('Calculating the velocity autocorrelation function')
198 call messages_info()
199
200 call calculate_vaf(vaf)
201
202 !print the vaf
203 iunit = io_open('td.general/velocity_autocorrelation', global_namespace, action='write')
204
205800 FORMAT(80('#'))
206 write(unit = iunit, iostat = ierr, fmt = 800)
207 write(unit = iunit, iostat = ierr, fmt = '(8a)') '# HEADER'
208 write(unit = iunit, iostat = ierr, fmt = '(a,4x,a6,a7,a1,10x,a10)') &
209 '# Iter', 'time [',units_out%time%abbrev,']', 'VAF [a.u.]'
210 write(unit = iunit, iostat = ierr, fmt = 800)
211
212 do jj = 1, nvaf
213 write(unit = iunit, iostat = ierr, fmt = *) jj, units_from_atomic(units_out%time, (jj - 1)*deltat), vaf(jj)
214 end do
215
216 call io_close(iunit)
217
218
219 safe_allocate(ftvaf(1:max_freq))
220
221 ftvaf = m_one
222
223 call batch_init(vafb, vaf)
224
225 call spectrum_signal_damp(spectrum%damp, spectrum%damp_factor, 1, nvaf, m_zero, deltat, vafb)
226
227 call batch_init(ftvafb, ftvaf)
228
229 call spectrum_fourier_transform(spectrum%method, spectrum_transform_cos, spectrum%noise, &
230 1, nvaf, m_zero, deltat, vafb, spectrum%min_energy, spectrum%max_energy, spectrum%energy_step, ftvafb)
231
232 call vafb%end()
233 call ftvafb%end()
234
235
236 !and print the spectrum
237 iunit = io_open('td.general/vibrational_spectrum', global_namespace, action='write')
238
239 write(unit = iunit, iostat = ierr, fmt = 800)
240 write(unit = iunit, iostat = ierr, fmt = '(8a)') '# HEADER'
241 write(unit = iunit, iostat = ierr, fmt = '(a17,6x,a15)') '# Energy [1/cm]', 'Spectrum [a.u.]'
242 write(unit = iunit, iostat = ierr, fmt = 800)
243
244 do ifreq = 1, max_freq
245 ww = spectrum%energy_step*(ifreq - 1) + spectrum%min_energy
246 write(unit = iunit, iostat = ierr, fmt = '(2e20.10)') units_from_atomic(unit_invcm, ww), ftvaf(ifreq)
247 end do
248
249 call io_close(iunit)
250
251 safe_deallocate_a(vaf)
252
253 safe_deallocate_a(ftvaf)
254
255 safe_deallocate_p(ions)
256
257 safe_deallocate_a(time)
258
259 safe_deallocate_a(velocities)
260
262 call io_end()
263 call messages_end()
264
265 call parser_end()
266 call global_end()
267
268contains
269
270 subroutine calculate_vaf(vaf)
271 real(real64), intent(out) :: vaf(:)
272 integer :: itm, itn
273
274 push_sub(calculate_vaf)
275
276 write (message(1), '(a)') "Read velocities from '"// &
277 trim(io_workpath('td.general/coordinates', global_namespace))//"'"
278 call messages_info(1)
279
280 !calculating the vaf, formula from
281 !
282 ! http://www.timteatro.net/2010/09/29/velocity-autocorrelation-and-vibrational-spectrum-calculation/
283
284 vaf = m_zero
285
286 do itm = 1, ntime
287 vaf(itm) = m_zero
288
289 do itn = 1, ntime - itm + 1
290
291 do ivel = 1, nvel
292 vaf(itm) = vaf(itm) + velocities(ivel, itm + itn - 1)*velocities(ivel, itn)
293 end do
294 end do
295
296 vaf(itm) = vaf(itm)/real(ntime - itm + 1, real64)
297
298 end do
299
300 do itm = 2, ntime
301 vaf(itm) = vaf(itm)/vaf(1)
302 end do
303 vaf(1) = m_one
304
305 pop_sub(calculate_vaf)
306 end subroutine calculate_vaf
307
308end program vibrational
309
310!! Local Variables:
311!! mode: f90
312!! coding: utf-8
313!! End:
initialize a batch with existing memory
Definition: batch.F90:277
This module implements batches of mesh functions.
Definition: batch.F90:135
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:494
real(real64), parameter, public m_zero
Definition: global.F90:200
type(mpi_comm), parameter, public serial_dummy_comm
Alias MPI_COMM_UNDEFINED for the specific use case of initialising Octopus utilities with no MPI supp...
Definition: global.F90:294
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:432
real(real64), parameter, public m_one
Definition: global.F90:201
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_skip_header(iunit)
Definition: io.F90:646
subroutine, public io_end()
Definition: io.F90:271
character(len=max_path_len) function, public io_workpath(path, namespace)
construct path name from given name and namespace
Definition: io.F90:318
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
subroutine, public messages_end()
Definition: messages.F90:273
subroutine, public messages_init(output_dir)
Definition: messages.F90:220
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1000
subroutine, public messages_new_line()
Definition: messages.F90:1089
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_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
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
type(namespace_t), public global_namespace
Definition: namespace.F90:135
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:410
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:442
subroutine, public profiling_end(namespace)
Definition: profiling.F90:415
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:257
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
subroutine, public spectrum_signal_damp(damp_type, damp_factor, time_start, time_end, t0, time_step, time_function)
Definition: spectrum.F90:2559
integer, parameter, public spectrum_transform_cos
Definition: spectrum.F90:173
pure integer function, public spectrum_nenergy_steps(spectrum)
Definition: spectrum.F90:2956
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_t), public unit_invcm
For vibrational frequencies.
type(unit_system_t), public units_out
subroutine, public unit_system_init(namespace)
int true(void)
program vibrational
subroutine calculate_vaf(vaf)