dcdfort
Modern Fortran library for analyzing DCD trajectories
Functions/Subroutines
dcdfort_utils Module Reference

Module that contains some useful utilities. More...

Functions/Subroutines

real(8) function, dimension(3) pbc (a, box)
 Corrects for the periodic boundary condition. More...
 
real(8) function, dimension(3) cross (a, b)
 Performs cross product between two vectors. More...
 
real(8) function distance2 (a, b, box)
 Calculates the distance squared between two points. More...
 
real(8) function distance (a, b, box)
 Calculates the distance between two points. More...
 
real(8) function, dimension(3) bond_vector (a, b, box)
 Calculates the bond vector between two points. More...
 
real(8) function magnitude (a)
 Calculates the magnitude of a vector. More...
 
real(8) function bond_angle (a, b, c, box)
 Calculates the bond angle between two vectors. More...
 
real(8) function dihedral_angle (i, j, k, l, box)
 Calculates the dihedral angle between two planes formed by four atoms. More...
 

Detailed Description

Module that contains some useful utilities.

Function/Subroutine Documentation

◆ bond_angle()

real(8) function dcdfort_utils::bond_angle ( real(8), dimension(3), intent(in)  a,
real(8), dimension(3), intent(in)  b,
real(8), dimension(3), intent(in)  c,
real(8), dimension(6), intent(in)  box 
)

Calculates the bond angle between two vectors.

Calculates the angle between the vector formed by a-b and b-c.

Parameters
[in]afirst point
[in]bmiddle point
[in]cthird point
[in]boxbox, if pbc to be accounted for
Returns
bond angle between a-b-c
162 
163  implicit none
164  real(8) :: bond_angle
165  real(8), intent(in), dimension(3) :: a, b, c
166  real(8), intent(in) :: box(6)
167  real(8), dimension(3) :: bond1, bond2
168 
169  bond1 = bond_vector(b, a, box)
170  bond2 = bond_vector(b, c, box)
171 
172  bond_angle = dacos(dot_product(bond1, bond2)/(magnitude(bond1)*magnitude(bond2)))
173 

◆ bond_vector()

real(8) function, dimension(3) dcdfort_utils::bond_vector ( real(8), dimension(3), intent(in)  a,
real(8), dimension(3), intent(in)  b,
real(8), dimension(6), intent(in)  box 
)

Calculates the bond vector between two points.

Parameters
[in]afirst point
[in]bsecond point
[in]boxbox, if pbc to be accounted for
Returns
bond vector pointing from a to b
131 
132  implicit none
133  real(8) :: bond_vector(3)
134  real(8), intent(in), dimension(3) :: a, b
135  real(8), intent(in) :: box(6)
136 
137  bond_vector = pbc(a-b, box)
138 

◆ cross()

real(8) function, dimension(3) dcdfort_utils::cross ( real(8), dimension(3), intent(in)  a,
real(8), dimension(3), intent(in)  b 
)

Performs cross product between two vectors.

Parameters
[in]afirst vector
[in]bsecond vector
Returns
resulting cross product
72 
73  implicit none
74  real(8) :: cross(3)
75  real(8), intent(in), dimension(3) :: a, b
76 
77  cross(1) = a(2) * b(3) - a(3) * b(2)
78  cross(2) = a(3) * b(1) - a(1) * b(3)
79  cross(3) = a(1) * b(2) - a(2) * b(1)
80 

◆ dihedral_angle()

real(8) function dcdfort_utils::dihedral_angle ( real(8), dimension(3), intent(in)  i,
real(8), dimension(3), intent(in)  j,
real(8), dimension(3), intent(in)  k,
real(8), dimension(3), intent(in)  l,
real(8), dimension(6), intent(in)  box 
)

Calculates the dihedral angle between two planes formed by four atoms.

Calculates the dihedral angle between the vectors formed by i-j, j-k, k-l

Parameters
[in]ifirst point
[in]jmiddle point
[in]kthird point
[in]lfourth point
[in]boxbox, if pbc to be accounted for
Returns
dihedral angle forms by i-j-k-l
185 
186  implicit none
187  real(8) :: dihedral_angle
188  real(8), intent(in), dimension(3) :: i, j, k, l
189  real(8), intent(in), dimension(6) :: box
190  real(8) :: a_mag, b_mag, g_mag
191  real(8), dimension(3) :: h, g, f, a, b, cross_ba
192 
193  h = bond_vector(k, l, box)
194  g = bond_vector(k, j, box)
195  f = bond_vector(j, i, box)
196  a = cross(f, g)
197  b = cross(h, g)
198  cross_ba = cross(b, a)
199  a_mag = magnitude(a)
200  b_mag = magnitude(b)
201  g_mag = magnitude(g)
202 
203  dihedral_angle = atan2(dot_product(cross_ba, g) / (a_mag*b_mag*g_mag), dot_product(a, b) / (a_mag*b_mag))
204 

◆ distance()

real(8) function dcdfort_utils::distance ( real(8), dimension(3), intent(in)  a,
real(8), dimension(3), intent(in)  b,
real(8), dimension(6), intent(in), optional  box 
)

Calculates the distance between two points.

Parameters
[in]afirst point
[in]bsecond point
[in]boxbox, if pbc to be accounted for
Returns
distance
111 
112  implicit none
113  real(8) :: distance
114  real(8), intent(in), dimension(3) :: a, b
115  real(8), intent(in), optional :: box(6)
116 
117  if (present(box)) then
118  distance = dsqrt(distance2(a, b, box))
119  else
120  distance = dsqrt(distance2(a, b))
121  end if
122 

◆ distance2()

real(8) function dcdfort_utils::distance2 ( real(8), dimension(3), intent(in)  a,
real(8), dimension(3), intent(in)  b,
real(8), dimension(6), intent(in), optional  box 
)

Calculates the distance squared between two points.

Parameters
[in]afirst point
[in]bsecond point
[in]boxbox, if pbc to be accounted for
Returns
distance squared
89 
90  implicit none
91  real(8) :: distance2
92  real(8), intent(in), dimension(3) :: a, b
93  real(8) :: c(3)
94  real(8), intent(in), optional :: box(6)
95 
96  if (present(box)) then
97  c = pbc(a - b, box)
98  else
99  c = a - b
100  end if
101  distance2 = dot_product(c, c)
102 

◆ magnitude()

real(8) function dcdfort_utils::magnitude ( real(8), dimension(3), intent(in)  a)

Calculates the magnitude of a vector.

Parameters
[in]avector
Returns
magnitude of a
145 
146  implicit none
147  real(8) :: magnitude
148  real(8), intent(in) :: a(3)
149 
150  magnitude = dsqrt(dot_product(a, a))
151 

◆ pbc()

real(8) function, dimension(3) dcdfort_utils::pbc ( real(8), dimension(3), intent(in)  a,
real(8), dimension(6), intent(in)  box 
)

Corrects for the periodic boundary condition.

Moves particle (or vector) the distance of half the box if it is more than half the distance of the box

Parameters
[in]aoriginal coordinates
[in]boxsimulation box
Returns
the shifted coordinates
41 
42  implicit none
43  real(8), intent(in) :: a(3), box(6)
44  real(8) :: pbc(3), tbox(3,3) = 0.0
45  integer :: i
46 
47  ! A = box(1), B = box(2), C = box(3)
48  ! alpha = box(4), beta = box(5), gamma = box(6)
49 
50  ! convert angles to box vectors
51  tbox(1,1) = box(1)
52 
53  tbox(1,2) = box(2)*cos(box(6)*degreestoradians)
54  tbox(2,2) = sqrt(box(2)**2-tbox(1,2)**2)
55 
56  tbox(1,3) = box(3)*cos(box(5)*degreestoradians)
57  tbox(2,3) = (box(2)*box(3)*cos(box(4)*degreestoradians) - tbox(1,2)*tbox(1,3))/tbox(2,2)
58  tbox(3,3) = sqrt(box(2)**2-tbox(1,3)**2-tbox(2,3)**2)
59 
60  pbc = a
61  do i = 3, 1, -1
62  pbc(1:i) = pbc(1:i) - tbox(1:i,i) * nint(pbc(i) / tbox(i,i))
63  end do
64