-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathm_interaction.f90
163 lines (140 loc) · 3.35 KB
/
m_interaction.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
!----------------------------------------------------------------------
! MELQUIADES : Metropolis Monte Carlo Program
!----------------------------------------------------------------------
!bop
!
! !Module: m_interaction
!
! !Description: This module contains a routine for
! interaction energy calculation.
!\\
!\\
! !Interface:
!
module m_interaction
!
! !Uses:
use m_kind
use m_simtype
use m_boxtype
use m_precells
use m_constants, only : c_kcalmol, pi
use m_zeros
use m_pbc
use m_ljones
!
! !Public member functions:
!
private
!
public :: r_interactions
!
! !Revision history:
! 06Aug 2015 Asdrubal Lozada
!
!eop
!----------------------------------------------------------------------
contains
!
!bop
!
! !Iroutine: r_interactions
!
! !Description: This routine uses a neighbors list for intermolecular
!energy calculation. Through a scheme of case selection different types
!of potential can be used. In the current version this potential functions
!are: The Lennard-Jones, Buckingham, Yukawa and Coulombic potentials.
!
!\\
!\\
! !Interface:
subroutine r_interactions( com, edge, i, jb, ids, energy, virial, t, y, x )
!
implicit none
!--------------------
!
! !Input arguments:
type(simulation), intent(inout) :: y
type(box), pointer :: x
type(temporary), intent(inout) :: t
real(rkind), dimension(:), intent(in) :: com
real(rkind), dimension(:), intent(in) :: edge
real(rkind), intent(inout) :: energy, virial
integer, intent(in) :: i, jb, ids
!-------------------------------------
!
! !Revision history:
! 06Aug 2015 Asdrubal Lozada
!
!eop
!----------------------------------------------------------------
! Local variables
real(rkind) :: rix, riy, riz
real(rkind) :: rjx, rjy, rjz
real(rkind), dimension(3) :: rij
real(rkind), dimension(6) :: v
!---------------------------------
integer :: icell
integer :: k, j, iter
real(rkind) :: rijsq, rm, rp
energy = 0.0_rkind
virial = 0.0_rkind
v = 0.0_rkind
rix = com(1)
riy = com(2)
riz = com(3)
icell = f_idcell( com, y, x )
call r_idneighs( icell, y, x )
do k = 1, y%m_viz
if( .not. y%m_solute ) then
iter = x%m_ncell(k)
j = x%m_head(iter)
else
iter = x%m_ncelsa(k)
j = x%m_hedsb(iter)
end if
do while(j /= 0)
if( j /= i .and. j >= jb ) then
rjx = x%m_cmass(1,j); rjy = x%m_cmass(2,j); rjz = x%m_cmass(3,j)
rij(1) = rix - rjx
rij(2) = riy - rjy
rij(3) = riz - rjz
if( y%m_pbcs ) then
call r_pbc(rij,edge)
end if
rijsq = rij(1)*rij(1) + rij(2)*rij(2) + rij(3)*rij(3)
rm = rijsq
if( .not.y%m_solute ) then
rp = y%m_rpair
else
rp = y%m_rsol
end if
if( rm<=rp ) then
!
!boc
select case(y%m_potens)
case('m_ljc')
call p_ljones(rij, v, i, j, ids, t, x)
!! case('m_bc')
!! call p_ljones(rij, v, i, j, ids, t, x)
!! call p_buckhm(rij, v, i, j, ids, t, x)
!! case('m_ljy')
!! call p_ljones(rij, v, i, j, ids, t, x)
!! call p_yukawa(rij, v, i, j, ids, t, x)
case default
write(*,*) 'Test without potential'
end select
!eoc
end if
end if
if( .not.y%m_solute ) then
j = x%m_list(j)
else
j = x%m_listb(j)
end if
end do ! While
end do ! k
energy = v(1) + v(2) + v(3)
virial = v(4) + v(5) + v(6)
end subroutine r_interactions
end module m_interaction