Octopus
casida_spectrum.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2012-2013 D. Strubbe
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22program casida_spectrum
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
31 use parser_oct_m
33 use space_oct_m
34 use unit_oct_m
36 use utils_oct_m
38
39 implicit none
40
42 real(real64) :: br, energy_step, min_energy, max_energy
43 integer :: ispin
44 end type casida_spectrum_t
45
46 integer :: ierr, idir, jdir, iatom
47 type(space_t) :: space
48 type(casida_spectrum_t) :: cs
49 real(real64), allocatable :: rotation(:,:), rot2(:,:), identity(:,:), coord(:)
50 type(block_t) :: blk
51 type(ions_t), pointer :: ions
52
53 ! Initialize stuff
55
56 call getopt_init(ierr)
57 if (ierr == 0) call getopt_casida_spectrum()
58 call getopt_end()
59
60 call parser_init()
61 call messages_init()
62 call io_init()
66
67 allocate(rotation(space%dim, space%dim))
68 allocate(rot2(space%dim, space%dim))
69 allocate(identity(space%dim, space%dim))
70 allocate(coord(space%dim))
71
72 ! Reads the spin components. This is read here, as well as in states_init.
73 call parse_variable(global_namespace, 'SpinComponents', 1, cs%ispin)
74 if (.not. varinfo_valid_option('SpinComponents', cs%ispin)) call messages_input_error(global_namespace, 'SpinComponents')
75 cs%ispin = min(2, cs%ispin)
76
77 !%Variable CasidaSpectrumBroadening
78 !%Type float
79 !%Default 0.005 Ha
80 !%Section Utilities::oct-casida_spectrum
81 !%Description
82 !% Width of the Lorentzian used to broaden the excitations.
83 !%End
84 call parse_variable(global_namespace, 'CasidaSpectrumBroadening', 0.005_real64, cs%br, units_inp%energy)
85
86 call messages_print_var_value('CasidaSpectrumBroadening', cs%br, unit = units_out%energy)
87
88 !%Variable CasidaSpectrumEnergyStep
89 !%Type float
90 !%Default 0.001 Ha
91 !%Section Utilities::oct-casida_spectrum
92 !%Description
93 !% Sampling rate for the spectrum.
94 !%End
95 call parse_variable(global_namespace, 'CasidaSpectrumEnergyStep', 0.001_real64, cs%energy_step, units_inp%energy)
96
97 call messages_print_var_value('CasidaSpectrumEnergyStep', cs%energy_step, unit = units_out%energy)
98
99 !%Variable CasidaSpectrumMinEnergy
100 !%Type float
101 !%Default 0.0
102 !%Section Utilities::oct-casida_spectrum
103 !%Description
104 !% The broadening is done for energies greater than <tt>CasidaSpectrumMinEnergy</tt>.
105 !%End
106 call parse_variable(global_namespace, 'CasidaSpectrumMinEnergy', m_zero, cs%min_energy, units_inp%energy)
107
108 call messages_print_var_value('CasidaSpectrumMinEnergy', cs%min_energy, unit = units_out%energy)
109
110 !%Variable CasidaSpectrumMaxEnergy
111 !%Type float
112 !%Default 1.0 Ha
113 !%Section Utilities::oct-casida_spectrum
114 !%Description
115 !% The broadening is done for energies smaller than <tt>CasidaSpectrumMaxEnergy</tt>.
116 !%End
117 call parse_variable(global_namespace, 'CasidaSpectrumMaxEnergy', m_one, cs%max_energy, units_inp%energy)
118
119 call messages_print_var_value('CasidaSpectrumMaxEnergy', cs%max_energy, unit = units_out%energy)
120
121 identity = m_zero
122 do idir = 1, space%dim
123 identity(idir, idir) = m_one
124 end do
125
126 !%Variable CasidaSpectrumRotationMatrix
127 !%Type block
128 !%Default identity
129 !%Section Utilities::oct-casida_spectrum
130 !%Description
131 !% Supply a rotation matrix to apply to the transition dipoles in generating the spectrum. The rotated atomic structure
132 !% will also be output. Size of matrix must be <tt>Dimensions</tt>.
133 !%End
135 if (parse_block(global_namespace, 'CasidaSpectrumRotationMatrix', blk) == 0) then
136 do idir = 1, space%dim
137 do jdir = 1, space%dim
138 call parse_block_float(blk, idir - 1, jdir - 1, rotation(idir, jdir))
139 end do
140 end do
141 call parse_block_end(blk)
142
143 message(1) = "Info: Applying rotation matrix"
144 call messages_info(1)
145 call output_tensor(rotation, space%dim, unit_one, write_average = .false., namespace=global_namespace)
146
147 ! allowing inversions is fine
148 rot2(:,:) = abs(matmul(transpose(rotation), rotation))
149 if (any(abs(rot2(:,:) - identity(:,:)) > 1e-6_real64)) then
150 write(message(1),'(a,es13.6)') "Rotation matrix is not orthogonal. max discrepancy in product = ", &
151 maxval(abs(rot2(:,:) - identity(:,:)))
152 call messages_warning(1)
153 end if
154
155 ! apply rotation to geometry
156 ions => ions_t(global_namespace)
157 do iatom = 1, ions%natoms
158 coord = ions%pos(:, iatom)
159 ions%pos(:, iatom) = matmul(rotation, coord)
160 end do
161 call ions%write_xyz(trim(casida_dir)//'rotated')
162 safe_deallocate_p(ions)
163 else
164 rotation(:,:) = identity(:,:)
165 end if
166
167 call calc_broad(cs, casida_dir, 'eps_diff', .true.)
168 call calc_broad(cs, casida_dir, 'petersilka', .false.)
169 call calc_broad(cs, casida_dir, 'tamm_dancoff', .false.)
170 call calc_broad(cs, casida_dir, 'variational', .false.)
171 call calc_broad(cs, casida_dir, 'casida', .false.)
172
173 deallocate(rotation)
174 deallocate(rot2)
175 deallocate(identity)
176 deallocate(coord)
177
179 call io_end()
180 call messages_end()
181
182 call parser_end()
183 call global_end()
184
185contains
186
187 ! ---------------------------------------------------------
188 subroutine calc_broad(cs, dir, fname, extracols)
189 type(casida_spectrum_t), intent(in) :: cs
190 character(len=*), intent(in) :: dir
191 character(len=*), intent(in) :: fname
192 logical, intent(in) :: extracols
193
194 real(real64), allocatable :: spectrum(:,:)
195 real(real64) :: omega, energy, re_tm(space%dim), im_tm(space%dim), ff(space%dim+1), tm_sq(space%dim)
196 integer :: istep, nsteps, iunit, trash(3), idir, ncols, ios
197 character(len=256) :: string
198 logical :: is_complex
199
200 push_sub(calc_broad)
201
202 nsteps = int((cs%max_energy - cs%min_energy) / cs%energy_step)
203 safe_allocate(spectrum(1:space%dim+1, 1:nsteps))
204 spectrum = m_zero
205
206 iunit = io_open(trim(dir)// fname, global_namespace, action='read', status='old', die = .false.)
207 if (iunit == -1) then
208 message(1) = 'Cannot open file "'//trim(dir)//trim(fname)//'".'
209 message(2) = 'The '//trim(fname)//' spectrum was not generated.'
210 call messages_warning(2)
211 pop_sub(calc_broad)
212 return
213 end if
214
215 ! For Casida, CV(2), Tamm-Dancoff, Petersilka: first column is the index of the excitation
216 ! For eps_diff: first two columns are occ and unocc states, then spin if spin-polarized
217 ncols = 1
218 if (extracols) then
219 ncols = ncols + cs%ispin
220 end if
221
222 read(iunit, *) ! skip header
223
224 read(iunit, '(a)') string
225 ! complex has this many columns; real lacks the imaginary columns (im_tm)
226 read(string, *, iostat = ios) trash(1:ncols), energy, (re_tm(idir), im_tm(idir), idir = 1, space%dim), ff(space%dim+1)
227 is_complex = (ios == 0)
228 ! put it back to read the data in the loop
229 backspace(iunit)
230
231 do
232 if (is_complex) then
233 read(iunit, *, iostat = ios) trash(1:ncols), energy, (re_tm(idir), im_tm(idir), idir = 1, space%dim), ff(space%dim+1)
234 im_tm = units_to_atomic(units_out%length, im_tm)
235 else
236 read(iunit, *, iostat = ios) trash(1:ncols), energy, (re_tm(idir), idir = 1, space%dim), ff(space%dim+1)
237 end if
238 re_tm = units_to_atomic(units_out%length, re_tm)
239 ! ff, the last column, is a dimensionless number
240
241 if (ios < 0) then
242 exit ! end of file
243 else if (ios > 0) then
244 message(1) = "Error parsing file " // trim(fname)
245 call messages_fatal(1)
246 end if
247
248 energy = units_to_atomic(units_out%energy, energy)
249 ! transition matrix elements by themselves are dependent on gauge in degenerate subspaces
250 ! make into oscillator strengths, as in casida_inc.F90 X(oscillator_strengths), and like the last column
251 tm_sq = (matmul(rotation, re_tm))**2
252 if (is_complex) then
253 tm_sq = tm_sq + (matmul(rotation, im_tm))**2
254 end if
255 ff(1:space%dim) = m_two * energy * tm_sq
256
257 do istep = 1, nsteps
258 omega = cs%min_energy + real(istep-1, real64) *cs%energy_step
259 spectrum(:, istep) = spectrum(:, istep) + ff*cs%br/((omega-energy)**2 + cs%br**2)/m_pi ! Lorentzian
260 end do
261 end do
262 call io_close(iunit)
263
264 ! print spectra
265 iunit = io_open(trim(dir)//"/spectrum."//fname, global_namespace, action='write')
266
267 write(iunit, '(a2,a12)', advance = 'no') '# ', 'E [' // trim(units_abbrev(units_out%energy)) // ']'
268 do idir = 1, space%dim
269 write(iunit, '(a14)', advance = 'no') '<' // index2axis(idir) // index2axis(idir) // '>'
270 end do
271 write(iunit, '(a14)') '<f>'
272
273 do istep = 1, nsteps
274 write(iunit, '(99es14.6)') units_from_atomic(units_out%energy, cs%min_energy + real(istep - 1, real64) &
275 *cs%energy_step), (units_from_atomic(unit_one/units_out%energy, spectrum(idir, istep)), idir = 1, space%dim+1)
276 end do
277
278 call io_close(iunit)
279
280 safe_deallocate_a(spectrum)
281 pop_sub(calc_broad)
282 end subroutine calc_broad
283
284end program casida_spectrum
285
286!! Local Variables:
287!! mode: f90
288!! coding: utf-8
289!! End:
subroutine calc_broad(cs, dir, fname, extracols)
program casida_spectrum
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:181
subroutine, public getopt_init(ierr)
Initializes the getopt machinery. Must be called before attempting to parse the options....
subroutine, public getopt_end
real(real64), parameter, public m_two
Definition: global.F90:190
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:382
real(real64), parameter, public m_zero
Definition: global.F90:188
character(len= *), parameter, public casida_dir
Definition: global.F90:258
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
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:265
subroutine, public init_octopus_globals(comm)
Initialise Octopus-specific global constants and files. This routine performs no initialisation calls...
Definition: global.F90:354
real(real64), parameter, public m_one
Definition: global.F90:189
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:161
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_end()
Definition: io.F90:258
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public messages_end()
Definition: messages.F90:278
subroutine, public messages_init(output_dir)
Definition: messages.F90:225
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
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
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_end(namespace)
Definition: profiling.F90:413
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:255
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
type(unit_t), public unit_one
some special units required for particular quantities
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
subroutine, public output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
Definition: utils.F90:242
character pure function, public index2axis(idir)
Definition: utils.F90:202
int true(void)