Octopus
propagation_spectrum.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
23 use debug_oct_m
24 use global_oct_m
25 use io_oct_m
26 use kick_oct_m
29 use parser_oct_m
33
34 implicit none
35
36 integer :: in_file(3), out_file(3), ref_file, eq_axes, nspin, &
37 lmax, time_steps, ierr
38 logical :: calculate_tensor, reference_multipoles
39 type(spectrum_t) :: spectrum
40 type(unit_system_t) :: file_units
41 character(len=80) :: refmultipoles
42
43 ! Initialize stuff
45
46 call getopt_init(ierr)
47 refmultipoles = ""
48 if (ierr == 0) call getopt_propagation_spectrum(refmultipoles)
49 call getopt_end()
50
51 call parser_init()
52
53 call messages_init()
54
55 call io_init()
57
59
60 call spectrum_init(spectrum, global_namespace)
61
62 select case (spectrum%spectype)
64 call read_files('multipoles', refmultipoles)
65 call calculate_absorption('cross_section', global_namespace)
66 case (spectrum_p_power)
67 call calculate_dipole_power("multipoles", 'dipole_power')
69 call calculate_ftchd('ftchd', 'dynamic_structure_factor')
71 call calculate_rotatory_strength("angular", "rotatory_strength")
72 case default
73 write(message(1), '(a)') 'No PropagationSpectrumType defined,'
74 write(message(2), '(a)') 'cannot calculate the spectrum.'
75 call messages_warning(2)
76 end select
77
79 call io_end()
80 call messages_end()
81
82 call parser_end()
83
84 call global_end()
85
86contains
87
120 subroutine read_files(fname, reffname)
121 character(len=*), intent(in) :: fname
122 character(len=*), intent(in) :: reffname
123
124 type(kick_t) :: kick
125 real(real64) :: dt
126
127 push_sub(read_files)
128
129 in_file(1) = io_open(trim(fname), global_namespace, action='read', status='old', die=.false.)
130 if (in_file(1) < 0) in_file(1) = io_open('td.general/'//trim(fname), global_namespace, &
131 action='read', status='old', die=.false.)
132 if (in_file(1) >= 0) then
133 write(message(1),'(3a)') 'File "', trim(fname), '" found. This will be the only file to be processed.'
134 write(message(2),'(a)') '(If more than one file is to be used, the files should be called'
135 write(message(3),'(5a)') '"', trim(fname), '.1", "', trim(fname), '.2", etc.)'
136 write(message(4),'(a)')
137 call messages_info(4)
138
139 ! OK, so we only have one file. Now we have to see what it has inside.
140 call spectrum_mult_info(global_namespace, in_file(1), nspin, kick, time_steps, dt, file_units, lmax=lmax)
141 eq_axes = kick%pol_equiv_axes
142 if (eq_axes == 3) then
143 calculate_tensor = .true.
144 write(message(1),'(3a)') 'The file "', trim(fname), '" tells me that the system has three equivalent axes.'
145 write(message(2),'(a)') 'I will calculate the full tensor, written in file "XXXX_tensor".'
146 call messages_info(2)
147 else if (eq_axes == 2) then
148 write(message(1),'(3a)') 'The file "', trim(fname), '" tells me that the system has two equivalent axes.'
149 write(message(2),'(a)') 'However, I am only using this file; cannot calculate the full tensor.'
150 write(message(3),'(a)') 'A file "XXXX_vector" will be generated instead.'
151 call messages_warning(3)
152 calculate_tensor = .false.
153 else
154 write(message(1),'(3a)') 'The file "', trim(fname), '" tells me that the system has no usable symmetry. '
155 write(message(2),'(a)') 'However, I am only using this file; cannot calculate the full tensor.'
156 write(message(3),'(a)') 'A file "XXXX_vector" will be generated instead.'
157 call messages_warning(3)
158 calculate_tensor = .false.
159 end if
160
161 else ! We will try to load more fname.1 files...
162
163 ! In this case, we will always want the full tensor
164 calculate_tensor = .true.
165
166 in_file(1) = io_open(trim(fname)//'.1', global_namespace, action='read', &
167 status='old', die=.false.)
168 if (in_file(1) < 0) in_file(1) = io_open('td.general/'//trim(fname)//'.1', global_namespace, &
169 action='read', status='old', die=.false.)
170 if (in_file(1) < 0) then ! Could not find proper files. Die and complain.
171 write(message(1),'(5a)') 'No "', trim(fname), '" or "', trim(fname), '.1" file found. At least one of those'
172 write(message(2),'(a)') 'should be visible.'
173 call messages_fatal(2)
174 end if
175
176 call spectrum_mult_info(global_namespace, in_file(1), nspin, kick, time_steps, dt, file_units, lmax=lmax)
177 eq_axes = kick%pol_equiv_axes
178
179 if (eq_axes == 3) then
180 write(message(1),'(3a)') 'The file "', trim(fname), '.1" tells me that the system has three equivalent axes.'
181 write(message(2),'(a)') 'I will calculate the full tensor, written in file "cross_section_tensor".'
182 call messages_info(2)
183
184 else if (eq_axes == 2) then
185 in_file(2) = io_open(trim(fname)//'.2', global_namespace, action='read', status='old', die=.false.)
186 if (in_file(2) < 0) in_file(2) = io_open('td.general/'//trim(fname)//'.2', global_namespace, &
187 action='read', status='old', die=.false.)
188 if (in_file(2) < 0) then
189 write(message(1),'(3a)') 'The file "', trim(fname), '.1" tells me that the system has two equivalent axes,'
190 write(message(2),'(3a)') 'but I cannot find a "', trim(fname), '.2".'
191 call messages_fatal(2)
192 end if
193 write(message(1),'(5a)') 'Found two files, "', trim(fname), '.1" and "', trim(fname), '.2".'
194 write(message(2),'(a)') 'Two polarization axes are equivalent. I will generate the full tensor.'
195 call messages_info(2)
196
197 else ! No equivalent axes
198 in_file(2) = io_open(trim(fname)//'.2', global_namespace, action='read', &
199 status='old', die=.false.)
200 if (in_file(2) < 0) in_file(2) = io_open('td.general/'//trim(fname)//'.2', global_namespace, &
201 action='read', status='old', die=.false.)
202 if (in_file(2) < 0) then
203 write(message(1),'(3a)') 'The file "', trim(fname), '.1" tells me that the system has three inequivalent axes,'
204 write(message(2),'(3a)') 'but I cannot find a "', trim(fname), '.2".'
205 call messages_fatal(2)
206 end if
207 in_file(3) = io_open(trim(fname)//'.3', global_namespace, action='read', &
208 status='old', die=.false.)
209 if (in_file(3) < 0) in_file(3) = io_open('td.general/'//trim(fname)//'.3', global_namespace, &
210 action='read', status='old', die=.false.)
211 if (in_file(3) < 0) then
212 write(message(1),'(3a)') 'The file "', trim(fname), '.1" tells me that the system has three inequivalent axes,'
213 write(message(2),'(3a)') 'but I cannot find a "', trim(fname), '.3".'
214 call messages_fatal(2)
215 end if
216 write(message(1),'(7a)') 'Found three files, "', trim(fname), '.1", "', trim(fname), '.2" and "', trim(fname), '.3".'
217 write(message(2),'(a)') 'No symmetry information will be used.'
218 call messages_info(2)
219 end if
220
221 end if
222
223 if (reffname == "") then
224 reference_multipoles = .false.
225 else
226 reference_multipoles = .true.
227 ref_file = io_open(trim(reffname), global_namespace, action='read', status='old', die=.false.)
228 if (ref_file < 0) then
229 write(message(1),'(3a)') 'No "',trim(reffname), '" file found.'
230 call messages_fatal(1)
231 end if
232 end if
233
234 pop_sub(read_files)
235 end subroutine read_files
236
237
238 !----------------------------------------------------------------------------
239 subroutine calculate_absorption(fname, namespace)
240 character(len=*), intent(in) :: fname
241 type(namespace_t), intent(in) :: namespace
242
243 integer :: ii, jj
244 character(len=150), allocatable :: filename(:)
245
246 push_sub(calculate_absorption)
247
248 if (.not. calculate_tensor) then
249
250 out_file(1) = io_open(trim(fname)//'_vector', global_namespace, action='write')
251 if (.not. reference_multipoles) then
252 call spectrum_cross_section(spectrum, namespace, in_file(1), out_file(1))
253 else
254 call spectrum_cross_section(spectrum, namespace, in_file(1), out_file(1), ref_file)
255 end if
256 call io_close(in_file(1))
257 call io_close(out_file(1))
258
259 else
260 select case (eq_axes)
261 case (0, 1)
262 jj = 3
263 case (2)
264 jj = 2
265 case (3)
266 jj = 1
267 end select
268
269 safe_allocate(filename(1:jj))
270 do ii = 1, jj
271 write(filename(ii),'(2a,i1)') trim(fname), '_vector.',ii
272 out_file(ii) = io_open(trim(filename(ii)), global_namespace, action='write')
273 if (.not. reference_multipoles) then
274 call spectrum_cross_section(spectrum, namespace, in_file(ii), out_file(ii))
275 else
276 call spectrum_cross_section(spectrum, namespace, in_file(ii), out_file(ii), ref_file)
277 end if
278 call io_close(in_file(ii))
279 call io_close(out_file(ii))
280 in_file(ii) = io_open(trim(filename(ii)), global_namespace, action='read', status='old')
281 end do
282
283 out_file(1) = io_open(trim(fname)//'_tensor', global_namespace, action='write')
284 call spectrum_cross_section_tensor(spectrum, namespace, out_file(1), in_file(1:jj))
285 do ii = 1, jj
286 call io_close(in_file(ii))
287 end do
288 call io_close(out_file(1))
289
290 end if
291
292 pop_sub(calculate_absorption)
293 end subroutine calculate_absorption
294
295 !----------------------------------------------------------------------------
296 subroutine calculate_dipole_power(fname_in, fname_out)
297 character(len=*), intent(in) :: fname_in, fname_out
298
299 push_sub(calculate_dipole_power)
300
301 in_file(1) = io_open(trim(fname_in), global_namespace, action='read', status='old', die=.false.)
302 if (in_file(1) < 0) in_file(1) = io_open('td.general/'//trim(fname_in), global_namespace, &
303 action='read', status='old', die=.false.)
304 if (in_file(1) >= 0) then
305 write(message(1),'(3a)') 'File "', trim(fname_in), '" found.'
306 write(message(2),'(a)')
307 call messages_info(2)
308 end if
309
310 out_file(1) = io_open(trim(fname_out), global_namespace, action='write')
311 call spectrum_dipole_power(spectrum, global_namespace, in_file(1), out_file(1))
312
313 call io_close(in_file(1))
314 call io_close(out_file(1))
315
317 end subroutine calculate_dipole_power
318
319 !----------------------------------------------------------------------------
320 subroutine calculate_rotatory_strength(fname_in, fname_out)
321 character(len=*), intent(in) :: fname_in, fname_out
322
324
325 in_file(1) = io_open(trim(fname_in), global_namespace, action='read', status='old', die=.false.)
326 if (in_file(1) < 0) in_file(1) = io_open('td.general/'//trim(fname_in), global_namespace, &
327 action='read', status='old', die=.false.)
328 if (in_file(1) >= 0) then
329 write(message(1),'(3a)') 'File "', trim(fname_in), '" found.'
330 write(message(2),'(a)')
331 call messages_info(2)
332 end if
333
334 out_file(1) = io_open(trim(fname_out), global_namespace, action='write')
335 call spectrum_rotatory_strength(spectrum, global_namespace, in_file(1), out_file(1))
336
337 call io_close(in_file(1))
338 call io_close(out_file(1))
339
341 end subroutine calculate_rotatory_strength
342
343
344 !----------------------------------------------------------------------------
345 subroutine calculate_ftchd(fname_in, fname_out)
346 character(len=*), intent(in) :: fname_in, fname_out
347
348 push_sub(calculate_ftchd)
349
350 ! read files
351 in_file(1) = io_open(trim(fname_in) // '.sin', global_namespace, action='read', status='old', die=.false.)
352 in_file(2) = io_open(trim(fname_in) // '.cos', global_namespace, action='read', status='old', die=.false.)
353
354 out_file(1) = io_open(trim(fname_out), global_namespace, action='write')
355 call spectrum_dyn_structure_factor(spectrum, global_namespace, in_file(1), in_file(2), out_file(1))
356 call io_close(in_file(1))
357 call io_close(out_file(1))
358
359 pop_sub(calculate_ftchd)
360 end subroutine calculate_ftchd
361
362end program propagation_spectrum
363
364!! Local Variables:
365!! mode: f90
366!! coding: utf-8
367!! End:
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:381
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:264
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:353
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
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
subroutine, public messages_end()
Definition: messages.F90:277
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 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
type(namespace_t), public global_namespace
Definition: namespace.F90:132
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:450
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:481
subroutine, public profiling_end(namespace)
Definition: profiling.F90:413
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:255
subroutine, public spectrum_cross_section(spectrum, namespace, in_file, out_file, ref_file)
Definition: spectrum.F90:710
subroutine, public spectrum_cross_section_tensor(spectrum, namespace, out_file, in_file)
Definition: spectrum.F90:430
subroutine, public spectrum_dyn_structure_factor(spectrum, namespace, in_file_sin, in_file_cos, out_file)
Definition: spectrum.F90:1262
subroutine, public spectrum_init(spectrum, namespace, default_energy_step, default_max_energy)
Definition: spectrum.F90:213
subroutine, public spectrum_mult_info(namespace, iunit, nspin, kick, time_steps, dt, file_units, lmax)
Definition: spectrum.F90:2342
subroutine, public spectrum_dipole_power(spectrum, namespace, in_file, out_file)
Definition: spectrum.F90:1154
integer, parameter, public spectrum_energyloss
Definition: spectrum.F90:176
integer, parameter, public spectrum_rotatory
Definition: spectrum.F90:176
subroutine, public spectrum_rotatory_strength(spectrum, namespace, in_file, out_file)
Definition: spectrum.F90:1415
integer, parameter, public spectrum_absorption
Definition: spectrum.F90:176
integer, parameter, public spectrum_p_power
Definition: spectrum.F90:176
This module defines the unit system, used for input and output.
subroutine, public unit_system_init(namespace)
subroutine calculate_rotatory_strength(fname_in, fname_out)
subroutine read_files(fname, reffname)
subroutine calculate_ftchd(fname_in, fname_out)
subroutine calculate_absorption(fname, namespace)
subroutine calculate_dipole_power(fname_in, fname_out)
program propagation_spectrum
int true(void)