Octopus
vibrations.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
22 use debug_oct_m
23 use global_oct_m
24 use io_oct_m
27 use mpi_oct_m
30 use space_oct_m
31 use unit_oct_m
33
34 implicit none
35
36 private
37 public :: &
51
52 type vibrations_t
53 private
54 integer, public :: num_modes
55 integer, public :: ndim
56 integer, public :: natoms
57 real(real64), allocatable, public :: dyn_matrix(:, :)
58 real(real64), allocatable, public :: infrared(:, :)
59 real(real64), allocatable, public :: normal_mode(:, :)
60 real(real64), allocatable, public :: freq(:)
61 real(real64), public :: disp
62 real(real64) :: total_mass
63 real(real64), pointer :: mass(:) => null()
64 character (len=2) :: suffix
65 character (len=80) :: filename_dynmat
66 type(unit_t) :: unit_dynmat
67 type(namespace_t), pointer, public :: namespace
68 end type vibrations_t
69
70contains
71
72 ! ---------------------------------------------------------
73 subroutine vibrations_init(this, space, natoms, mass, suffix, namespace)
74 type(vibrations_t), intent(out) :: this
75 class(space_t), intent(in) :: space
76 integer, intent(in) :: natoms
77 real(real64), target, intent(in) :: mass(:)
78 character (len=2), intent(in) :: suffix
79 type(namespace_t), target, intent(in) :: namespace
80
81 push_sub(vibrations_init)
82
83 this%ndim = space%dim
84 this%natoms = natoms
85 this%num_modes = natoms*space%dim
86 this%namespace => namespace
87 safe_allocate(this%dyn_matrix(1:this%num_modes, 1:this%num_modes))
88 safe_allocate(this%infrared(1:this%num_modes, 1:this%ndim))
89 safe_allocate(this%normal_mode(1:this%num_modes, 1:this%num_modes))
90 safe_allocate(this%freq(1:this%num_modes))
91
92 this%mass => mass
93 this%total_mass = sum(mass)
94
95 ! Since frequencies are reported as invcm, the matrix they are derived from can be expressed in invcm**2.
96 ! However, that matrix is the dynamical matrix divided by the total mass, so it has a different unit.
97 this%unit_dynmat = units_out%energy / units_out%length**2
98
99 this%suffix = suffix
100 this%filename_dynmat = vib_modes_dir//'dynamical_matrix_'//trim(this%suffix)
101 if (mpi_world%is_root()) then
102 call io_mkdir(vib_modes_dir, namespace)
103 call io_rm(this%filename_dynmat, namespace)
105 end if
106
107 pop_sub(vibrations_init)
108 end subroutine vibrations_init
109
110
111 ! ---------------------------------------------------------
112 subroutine vibrations_end(this)
113 type(vibrations_t), intent(inout) :: this
114
115 push_sub(vibrations_end)
117 safe_deallocate_a(this%dyn_matrix)
118 safe_deallocate_a(this%infrared)
119 safe_deallocate_a(this%freq)
120 safe_deallocate_a(this%normal_mode)
121
122 pop_sub(vibrations_end)
123 end subroutine vibrations_end
124
125
126 ! ---------------------------------------------------------
127 character(len=2) pure function vibrations_get_suffix(this)
128 type(vibrations_t), intent(in) :: this
129
130 vibrations_get_suffix = this%suffix
131 end function vibrations_get_suffix
132
133
134 ! ---------------------------------------------------------
136 subroutine vibrations_symmetrize_dyn_matrix(this)
137 type(vibrations_t), intent(inout) :: this
138
139 integer :: imat, jmat
140 real(real64) :: average, maxdiff
141
143
144 ! FIXME: enforce acoustic sum rule.
145
146 maxdiff = m_zero
147 do imat = 1, this%num_modes
148 do jmat = imat + 1, this%num_modes
149 average = m_half * (this%dyn_matrix(imat, jmat) + this%dyn_matrix(jmat, imat))
150 maxdiff = max(maxdiff, abs(this%dyn_matrix(imat, jmat) - this%dyn_matrix(jmat, imat)))
151 this%dyn_matrix(imat, jmat) = average
152 this%dyn_matrix(jmat, imat) = average
153 end do
154 end do
156 write(message(1),'(a)') 'Info: Symmetrizing dynamical matrix.'
157 write(message(2),'(a,es13.6,a,a)') 'Info: Maximum discrepancy from symmetry: ', &
158 units_from_atomic(this%unit_dynmat, maxdiff), &
159 " ", trim(units_abbrev(this%unit_dynmat))
160 call messages_info(2, namespace=this%namespace)
164
165 ! ---------------------------------------------------------
166 real(real64) pure function vibrations_norm_factor(this, iatom, jatom)
167 type(vibrations_t), intent(in) :: this
168 integer, intent(in) :: iatom
169 integer, intent(in) :: jatom
170
171 vibrations_norm_factor = this%total_mass / &
172 sqrt(this%mass(iatom) * this%mass(jatom))
173
174 end function vibrations_norm_factor
175
176 ! ---------------------------------------------------------
178 subroutine vibrations_out_dyn_matrix_row(this, imat)
179 type(vibrations_t), intent(in) :: this
180 integer, intent(in) :: imat
181
182 integer :: iatom, idir, jatom, jdir, jmat, iunit
183
184 if (.not. mpi_world%is_root()) return
185
187
188 iatom = vibrations_get_atom(this, imat)
189 idir = vibrations_get_dir(this, imat)
190
191 iunit = io_open(this%filename_dynmat, this%namespace, action='write', position='append')
192
193 do jmat = 1, this%num_modes
194 jatom = vibrations_get_atom(this, jmat)
195 jdir = vibrations_get_dir(this, jmat)
196
197 write(iunit, '(2(i8, i6), e25.12)') &
198 jatom, jdir, iatom, idir, units_from_atomic(this%unit_dynmat, this%dyn_matrix(jmat, imat))
199 end do
200 call io_close(iunit)
201
203 end subroutine vibrations_out_dyn_matrix_row
204
205 ! ---------------------------------------------------------
206 subroutine vibrations_out_dyn_matrix_header(this)
207 type(vibrations_t), intent(in) :: this
208
209 integer :: iunit
210
211 if (.not. mpi_world%is_root()) return
212
214
215 iunit = io_open(this%filename_dynmat, namespace=this%namespace, action='write') ! start at the beginning
216 write(iunit, '(2(a8, a6), a25)') 'atom', 'dir', 'atom', 'dir', &
217 '[' // trim(units_abbrev(this%unit_dynmat)) // ']'
218 call io_close(iunit)
219
223 ! ---------------------------------------------------------
225 subroutine vibrations_diag_dyn_matrix(this)
226 type(vibrations_t), intent(inout) :: this
227
228 integer :: imode
229
232 this%normal_mode = m_zero
233
234 this%normal_mode(:, :) = this%dyn_matrix(:, :)
235 call lalg_eigensolve(this%num_modes, this%normal_mode, this%freq)
236
237 ! FIXME: why not remove total mass here and in norm_factor?
238 this%freq(1:this%num_modes) = -this%freq(1:this%num_modes) / this%total_mass
239
240 if (any(this%freq(1:this%num_modes) < -m_epsilon)) then
241 message(1) = "There are imaginary vibrational frequencies (represented as negative)."
242 call messages_warning(1, namespace=this%namespace)
243 end if
244
245 do imode = 1, this%num_modes
246 this%freq(imode) = sign(sqrt(abs(this%freq(imode))), this%freq(imode))
247
248 ! make the largest component positive, to specify the phase
249 if (maxval(this%normal_mode(:, imode)) - abs(minval(this%normal_mode(:, imode))) < -m_epsilon) then
250 this%normal_mode(:, imode) = -this%normal_mode(:, imode)
251 end if
252 end do
253
255 end subroutine vibrations_diag_dyn_matrix
256
257
258 ! ---------------------------------------------------------
259 integer pure function vibrations_get_index(this, iatom, idim)
260 type(vibrations_t), intent(in) :: this
261 integer, intent(in) :: iatom
262 integer, intent(in) :: idim
263
264 vibrations_get_index = (iatom - 1)*this%ndim + idim
265 end function vibrations_get_index
266
267
268 ! ---------------------------------------------------------
269 integer pure function vibrations_get_atom(this, index)
270 type(vibrations_t), intent(in) :: this
271 integer, intent(in) :: index
272
273 vibrations_get_atom = 1 + (index - 1)/ this%ndim
274 end function vibrations_get_atom
275
276
277 ! ---------------------------------------------------------
278 integer pure function vibrations_get_dir(this, index)
279 type(vibrations_t), intent(in) :: this
280 integer, intent(in) :: index
281
282 vibrations_get_dir = 1 + mod(index - 1, this%ndim)
283 end function vibrations_get_dir
284
285
286 ! ---------------------------------------------------------
288 subroutine vibrations_output(this)
289 type(vibrations_t), intent(in) :: this
290
291 integer :: iunit, imat, jmat
292
293 if (.not. mpi_world%is_root()) return
294
295 push_sub(vibrations_output)
296
297 ! output frequencies and eigenvectors
298 iunit = io_open(vib_modes_dir//'normal_frequencies_'//trim(this%suffix), this%namespace, action='write')
299 write(iunit,'(a)') "# Mode Frequency [cm-1]"
300 do imat = 1, this%num_modes
301 write(iunit, '(i6,f17.8)') imat, units_from_atomic(unit_invcm, this%freq(imat))
302 end do
303 call io_close(iunit)
304
305 ! output eigenvectors
306 iunit = io_open(vib_modes_dir//'normal_modes_'//trim(this%suffix), this%namespace, action='write')
307 do imat = 1, this%num_modes
308 write(iunit, '(i6)', advance='no') imat
309 do jmat = 1, this%num_modes
310 write(iunit, '(es14.5)', advance='no') this%normal_mode(jmat, imat)
311 end do
312 write(iunit, '(1x)')
313 end do
314 call io_close(iunit)
315
316 ! output frequencies and eigenvectors for the phonon_modes_t.
317 iunit = io_open(vib_modes_dir//'vibration_modes', this%namespace, action='write')
318 write(iunit, '("Version: 0.0.0")')
319 write(iunit, '("Nmodes: ", I4)') this%num_modes
320 write(iunit, '("Np: 1")')
321 do imat=1, this%num_modes
322 write(iunit, '("# ")')
323 write(iunit, '("frequency:",E17.7)') this%freq(imat)
324 write(iunit, '(3E17.7)') this%normal_mode(:, imat)
325 end do
326 call io_close(iunit)
327
328 pop_sub(vibrations_output)
329 end subroutine vibrations_output
330
331end module vibrations_oct_m
332
333!! Local Variables:
334!! mode: f90
335!! coding: utf-8
336!! End:
double sqrt(double __x) __attribute__((__nothrow__
character(len= *), parameter, public vib_modes_dir
Definition: global.F90:274
Definition: io.F90:116
subroutine, public io_rm(fname, namespace)
Definition: io.F90:392
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:361
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
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_system_t), public units_out
character(len=2) pure function, public vibrations_get_suffix(this)
Definition: vibrations.F90:223
subroutine, public vibrations_out_dyn_matrix_header(this)
Definition: vibrations.F90:302
real(real64) pure function, public vibrations_norm_factor(this, iatom, jatom)
Definition: vibrations.F90:262
subroutine, public vibrations_diag_dyn_matrix(this)
Diagonalize the dynamical matrix.
Definition: vibrations.F90:321
subroutine, public vibrations_out_dyn_matrix_row(this, imat)
Outputs one row of the dynamical matrix.
Definition: vibrations.F90:274
subroutine, public vibrations_init(this, space, natoms, mass, suffix, namespace)
Definition: vibrations.F90:169
integer pure function, public vibrations_get_dir(this, index)
Definition: vibrations.F90:374
subroutine, public vibrations_symmetrize_dyn_matrix(this)
Symmetrize the dynamical matric, which is real symmetric matrix.
Definition: vibrations.F90:232
integer pure function, public vibrations_get_index(this, iatom, idim)
Definition: vibrations.F90:355
subroutine, public vibrations_output(this)
Outputs the eigenvectors and eigenenergies of the dynamical matrix.
Definition: vibrations.F90:384
subroutine, public vibrations_end(this)
Definition: vibrations.F90:208
integer pure function, public vibrations_get_atom(this, index)
Definition: vibrations.F90:365