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
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
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))
93 this%total_mass = sum(mass)
100 this%filename_dynmat =
vib_modes_dir//
'dynamical_matrix_'//trim(this%suffix)
103 call io_rm(this%filename_dynmat, namespace)
113 type(vibrations_t),
intent(inout) :: this
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)
127 character(len=2) pure function vibrations_get_suffix(this)
128 type(vibrations_t),
intent(in) :: this
130 vibrations_get_suffix = this%suffix
137 type(vibrations_t),
intent(inout) :: this
139 integer :: imat, jmat
140 real(real64) :: average, maxdiff
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
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)
168 integer,
intent(in) :: iatom
169 integer,
intent(in) :: jatom
172 sqrt(this%mass(iatom) * this%mass(jatom))
180 integer,
intent(in) :: imat
182 integer :: iatom, idir, jatom, jdir, jmat, iunit
184 if (.not. mpi_world%is_root())
return
191 iunit = io_open(this%filename_dynmat, this%namespace, action=
'write', position=
'append')
193 do jmat = 1, this%num_modes
197 write(iunit,
'(2(i8, i6), e25.12)') &
198 jatom, jdir, iatom, idir, units_from_atomic(this%unit_dynmat, this%dyn_matrix(jmat, imat))
211 if (.not. mpi_world%is_root())
return
215 iunit = io_open(this%filename_dynmat, namespace=this%namespace, action=
'write')
216 write(iunit,
'(2(a8, a6), a25)')
'atom',
'dir',
'atom',
'dir', &
217 '[' // trim(units_abbrev(this%unit_dynmat)) //
']'
232 this%normal_mode = m_zero
234 this%normal_mode(:, :) = this%dyn_matrix(:, :)
235 call lalg_eigensolve(this%num_modes, this%normal_mode, this%freq)
238 this%freq(1:this%num_modes) = -this%freq(1:this%num_modes) / this%total_mass
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)
245 do imode = 1, this%num_modes
246 this%freq(imode) = sign(
sqrt(abs(this%freq(imode))), this%freq(imode))
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)
261 integer,
intent(in) :: iatom
262 integer,
intent(in) :: idim
271 integer,
intent(in) :: index
280 integer,
intent(in) :: index
291 integer :: iunit, imat, jmat
293 if (.not. mpi_world%is_root())
return
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))
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)
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)
double sqrt(double __x) __attribute__((__nothrow__
character(len= *), parameter, public vib_modes_dir
subroutine, public io_rm(fname, namespace)
subroutine, public io_mkdir(fname, namespace, parents)
type(mpi_grp_t), public mpi_world
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
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)
subroutine, public vibrations_out_dyn_matrix_header(this)
real(real64) pure function, public vibrations_norm_factor(this, iatom, jatom)
subroutine, public vibrations_diag_dyn_matrix(this)
Diagonalize the dynamical matrix.
subroutine, public vibrations_out_dyn_matrix_row(this, imat)
Outputs one row of the dynamical matrix.
subroutine, public vibrations_init(this, space, natoms, mass, suffix, namespace)
integer pure function, public vibrations_get_dir(this, index)
subroutine, public vibrations_symmetrize_dyn_matrix(this)
Symmetrize the dynamical matric, which is real symmetric matrix.
integer pure function, public vibrations_get_index(this, iatom, idim)
subroutine, public vibrations_output(this)
Outputs the eigenvectors and eigenenergies of the dynamical matrix.
subroutine, public vibrations_end(this)
integer pure function, public vibrations_get_atom(this, index)