forked from lammps/lammps
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathTubePotTrue.f90
443 lines (415 loc) · 24.2 KB
/
TubePotTrue.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
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
! ------------ ----------------------------------------------------------
! LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
! http://lammps.sandia.gov, Sandia National Laboratories
! Steve Plimpton, [email protected]
!
! Copyright (2003) Sandia Corporation. Under the terms of Contract
! DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
! certain rights in this software. This software is distributed under
! the GNU General Public License.
!
! See the README file in the top-level LAMMPS directory.
!
! Contributing author: Alexey N. Volkov, UA, [email protected]
!-------------------------------------------------------------------------
module TubePotTrue !********************************************************************************
!
! TMD Library: True tubular potential and transfer function
!
!---------------------------------------------------------------------------------------------------
!
! Intel Fortran
!
! Alexey N. Volkov, University of Alabama, [email protected], Version 09.01, 2017
!
!---------------------------------------------------------------------------------------------------
!
! This module implements calculation of the true potential and transfer functions for interaction
! between two cylinder segments of nanotubes by direct integration over the surfaces of both
! segments.
!
!***************************************************************************************************
use TPMGeom
use TubePotBase
use iso_c_binding, only : c_int, c_double, c_char
implicit none
!---------------------------------------------------------------------------------------------------
! Constants
!---------------------------------------------------------------------------------------------------
integer(c_int), parameter :: TPTNXMAX = 257
integer(c_int), parameter :: TPTNEMAX = 128
!---------------------------------------------------------------------------------------------------
! Types
!---------------------------------------------------------------------------------------------------
type TPTSEG !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
real(c_double) :: X, Y, Z
real(c_double) :: Psi, Theta, Phi ! Euler's angles
real(c_double) :: R ! Segment radius
real(c_double) :: L ! Segment length
integer(c_int) :: NX, NE ! Number of nodes for numerical integration
real(c_double) :: DX, DE ! Spacings
real(c_double), dimension(0:2,0:2) :: M ! Transformation matrix
real(c_double), dimension(0:TPTNXMAX-1,0:TPTNXMAX-1,0:2) :: Rtab! Node coordinates
end type TPTSEG !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!---------------------------------------------------------------------------------------------------
! Global variables
!---------------------------------------------------------------------------------------------------
type(TPTSEG) :: TPTSeg1, TPTSeg2 ! Two segments
contains !******************************************************************************************
subroutine TPTSegAxisVector ( S, Laxis ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(in) :: S
real(c_double), dimension(0:2), intent(out) :: Laxis
!-------------------------------------------------------------------------------------------
Laxis(0:2) = S%M(2,0:2)
end subroutine TPTSegAxisVector !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine TPTSegRadVector ( S, Lrad, Eps ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(in) :: S
real(c_double), dimension(0:2), intent(out) :: Lrad
real(c_double), intent(in) :: Eps
!-------------------------------------------------------------------------------------------
real(c_double) :: Ce, Se
!-------------------------------------------------------------------------------------------
Ce = cos ( Eps )
Se = sin ( Eps )
Lrad(0) = Ce * S%M(0,0) + Se * S%M(1,0)
Lrad(1) = Ce * S%M(0,1) + Se * S%M(1,1)
Lrad(2) = Ce * S%M(0,2) + Se * S%M(1,2)
end subroutine TPTSegRadVector !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine TPTRadiusVector ( S, R, X, Eps ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(in) :: S
real(c_double), dimension(0:2), intent(out) :: R
real(c_double), intent(in) :: X, Eps
!-------------------------------------------------------------------------------------------
real(c_double), dimension(0:2) :: Laxis, Lrad
!-------------------------------------------------------------------------------------------
call TPTSegAxisVector ( S, Laxis )
call TPTSegRadVector ( S, Lrad, Eps )
R(0) = S%X + X * Laxis(0) + S%R * Lrad(0)
R(1) = S%Y + X * Laxis(1) + S%R * Lrad(1)
R(2) = S%Z + X * Laxis(2) + S%R * Lrad(2)
end subroutine TPTRadiusVector !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine TPTCalcSegNodeTable ( S ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(inout) :: S
!-------------------------------------------------------------------------------------------
real(c_double) :: X, Eps
integer(c_int) :: i, j
!-------------------------------------------------------------------------------------------
X = - S%L / 2.0
call RotationMatrix3 ( S%M, S%Psi, S%Theta, S%Phi )
do i = 0, S%NX - 1
Eps = 0.0d+00
do j = 0, S%NE - 1
call TPTRadiusVector ( S, S%Rtab(i,j,0:2), X, Eps )
Eps = Eps + S%DE
end do
X = X + S%DX
end do
end subroutine TPTCalcSegNodeTable !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine TPTSetSegPosition1 ( S, Rcenter, Laxis, L ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(inout) :: S
real(c_double), dimension(0:2), intent(in) :: Rcenter, Laxis
real(c_double), intent(in) :: L
!-------------------------------------------------------------------------------------------
S%L = L
S%DX = L / ( S%NX - 1 )
call EulerAngles ( S%Psi, S%Theta, Laxis )
S%Phi= 0.0d+00
S%X = Rcenter(0)
S%Y = Rcenter(1)
S%Z = Rcenter(2)
call TPTCalcSegNodeTable ( S )
end subroutine TPTSetSegPosition1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine TPTSetSegPosition2 ( S, R1, R2 ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(inout) :: S
real(c_double), dimension(0:2), intent(in) :: R1, R2
!-------------------------------------------------------------------------------------------
real(c_double), dimension(0:2) :: R, Laxis
real(c_double) :: L
!-------------------------------------------------------------------------------------------
R = 0.5 * ( R1 + R2 )
Laxis = R2 - R1
L = S_V3norm3 ( Laxis )
Laxis = Laxis / L
call TPTSetSegPosition1 ( S, R, Laxis, L )
end subroutine TPTSetSegPosition2 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
integer(c_int) function TPTCheckIntersection ( S1, S2 ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(in) :: S1, S2
!-------------------------------------------------------------------------------------------
integer(c_int) :: i, j
real(c_double) :: L1, L2, Displacement, D
real(c_double), dimension(0:2) :: Laxis, Q, R
!-------------------------------------------------------------------------------------------
L2 = S1%L / 2.0
L1 = - L2
call TPTSegAxisVector ( S1, Laxis )
R(0) = S1%X
R(1) = S1%Y
R(2) = S1%Z
do i = 0, S2%NX - 1
do j = 0, S2%NE - 1
call LinePoint ( Displacement, Q, R, Laxis, S2%Rtab(i,j,0:2) )
D = sqrt ( sqr ( Q(0) - S2%Rtab(i,j,0) ) + sqr ( Q(1) - S2%Rtab(i,j,1) ) &
+ sqr ( Q(2) - S2%Rtab(i,j,2) ) )
if ( Displacement > L1 .and. Displacement < L2 .and. D < S1%R ) then
TPTCheckIntersection = 1
return
end if
end do
end do
TPTCheckIntersection = 0
end function TPTCheckIntersection !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
integer(c_int) function TPTCalcPointRange ( S, Xmin, Xmax, Re ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!
type(TPTSEG), intent(in) :: S
real(c_double), intent(out) :: Xmin, Xmax
real(c_double), dimension(0:2), intent(in) :: Re
!-------------------------------------------------------------------------------------------
real(c_double) :: Displacement, Distance
real(c_double), dimension(0:2) :: Laxis, Q, R
!-------------------------------------------------------------------------------------------
call TPTSegAxisVector ( S, Laxis )
R(0) = S%X
R(1) = S%Y
R(2) = S%Z
call LinePoint ( Displacement, Q, R, Laxis, Re )
Distance = sqrt ( sqr ( Q(0) - Re(0) ) + sqr ( Q(1) - Re(1) ) + sqr ( Q(2) - Re(2) ) ) - S%R
if ( TPBRcutoff < Distance ) then
Xmin = 0.0d+00
Xmax = 0.0d+00
TPTCalcPointRange = 0
return
end if
Distance = sqrt ( TPBRcutoff * TPBRcutoff - Distance * Distance )
Xmin = Displacement - Distance
Xmax = Displacement + Distance
TPTCalcPointRange = 1
end function TPTCalcPointRange !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine TPTGetEnds ( R1_1, R1_2, R2_1, R2_2, X1_1, X1_2, X2_1, X2_2, H, A ) !!!!!!!!!!!!!
real(c_double), dimension(0:2), intent(out) :: R1_1, R1_2, R2_1, R2_2
real(c_double), intent(in) :: X1_1, X1_2, X2_1, X2_2, H, A
!-------------------------------------------------------------------------------------------
R1_1(0) = 0.0d+00
R1_1(1) = 0.0d+00
R1_1(2) = X1_1
R1_2(0) = 0.0d+00
R1_2(1) = 0.0d+00
R1_2(2) = X1_2
R2_1(0) = H
R2_1(1) = - X2_1 * sin ( A )
R2_1(2) = X2_1 * cos ( A )
R2_2(0) = H
R2_2(1) = - X2_2 * sin ( A )
R2_2(2) = X2_2 * cos ( A )
end subroutine TPTGetEnds !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!---------------------------------------------------------------------------------------------------
! Tubular potential
!---------------------------------------------------------------------------------------------------
integer(c_int) function TPTPointPotential ( Q, U, F, R, S ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! This function returns the potential U and force F applied to an atom in position R and
! produced by the segment S.
!-------------------------------------------------------------------------------------------
real(c_double), intent(out) :: Q, U
real(c_double), dimension(0:2), intent(out) :: F
real(c_double), dimension(0:2), intent(in) :: R
type(TPTSEG), intent(in) :: S
!-------------------------------------------------------------------------------------------
integer(c_int) :: i, j
real(c_double), dimension(0:2) :: RR, FF
real(c_double) :: QQ, UU, UUU, FFF, Rabs
real(c_double) :: Coeff, Xmin, Xmax, X
!-------------------------------------------------------------------------------------------
TPTPointPotential = 0
Q = 0.0d+00
U = 0.0d+00
F = 0.0d+00
if ( TPTCalcPointRange ( S, Xmin, Xmax, R ) == 0 ) return
X = - S%L / 2.0
do i = 0, S%NX - 1
if ( X > Xmin .and. X < Xmax ) then
QQ = 0.0d+00
UU = 0.0d+00
FF = 0.0d+00
do j = 0, S%NE - 1
RR(0:2) = S%Rtab(i,j,0:2) - R(0:2)
Rabs = S_V3norm3 ( RR )
if ( Rabs < TPBRcutoff ) then
QQ = QQ + TPBQCalc0 ( Rabs )
call TPBUCalc1 ( UUU, FFF, Rabs )
UU = UU + UUU
FFF = FFF / Rabs
FF = FF + FFF * RR
TPTPointPotential = 1
end if
end do
if ( i == 0 .or. i == S%NX - 1 ) then
Q = Q + 0.5d+00 * QQ
U = U + 0.5d+00 * UU
F = F + 0.5d+00 * FF
else
Q = Q + QQ
U = U + UU
F = F + FF
end if
end if
X = X + S%DX
end do
Coeff = TPBD * S%DX * S%R * S%DE
Q = Q * S%DX * S%R * S%DE
U = U * Coeff
F = F * Coeff
end function TPTPointPotential !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
integer(c_int) function TPTSectionPotential ( Q, U, F, M, S, i, Ssource ) !!!!!!!!!!!!!!!!!!
! This function returns the potential U, force F and torque M produced by the segment Ssource
! and applied to the i-th circular cross-section of the segment S.
!-------------------------------------------------------------------------------------------
real(c_double), intent(out) :: Q, U
real(c_double), dimension(0:2), intent(out) :: F, M
type(TPTSEG), intent(in) :: S, Ssource
integer(c_int), intent(in) :: i
!-------------------------------------------------------------------------------------------
integer(c_int) :: j
real(c_double), dimension(0:2) :: R, Fp, Mp, Lrad
real(c_double) :: Qp, Up, Eps
real(c_double) :: Coeff
!-------------------------------------------------------------------------------------------
TPTSectionPotential = 0
Q = 0.0d+00
U = 0.0d+00
F = 0.0d+00
M = 0.0d+00
Eps = 0.0d+00
do j = 0, S%NE - 1
call TPTSegRadVector ( S, Lrad, Eps )
if ( TPTPointPotential ( Qp, Up, Fp, S%Rtab(i,j,0:2), Ssource ) == 1 ) then
Q = Q + Qp
U = U + Up
F = F + Fp
R(0) = S%Rtab(i,j,0) - S%X
R(1) = S%Rtab(i,j,1) - S%Y
R(2) = S%Rtab(i,j,2) - S%Z
call V3_V3xxV3 ( Mp, R, Fp )
M = M + Mp
TPTSectionPotential = 1
end if
Eps = Eps + S%DE
end do
Coeff = TPBD * S%R * S%DE
Q = Q * S%R * S%DE
U = U * Coeff
F = F * Coeff
M = M * Coeff
end function TPTSectionPotential !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
integer(c_int) function TPTSegmentPotential ( Q, U, F, M, S, Ssource ) !!!!!!!!!!!!!!!!!!!!!!
! This function returns the potential U, force F and torque M produced by the segment
! Ssource and applied to the segment S.
!-------------------------------------------------------------------------------------------
real(c_double), intent(out) :: Q, U
real(c_double), dimension(0:2), intent(out) :: F, M
type(TPTSEG), intent(in) :: S, Ssource
integer(c_int) :: i
real(c_double), dimension(0:2) :: Fc, Mc
real(c_double) :: Qc, Uc
!-------------------------------------------------------------------------------------------
TPTSegmentPotential = 0
Q = 0.0d+00
U = 0.0d+00
F = 0.0d+00
M = 0.0d+00
if ( TPTCheckIntersection ( S, Ssource ) == 1 ) then
TPTSegmentPotential = 2
return
end if
do i = 0, S%NX - 1
if ( TPTSectionPotential ( Qc, Uc, Fc, Mc, S, i, Ssource ) == 1 ) then
if ( i == 0 .or. i == S%NX - 1 ) then
Q = Q + 0.5d+00 * Qc
U = U + 0.5d+00 * Uc
F = F + 0.5d+00 * Fc
M = M + 0.5d+00 * Mc
else
Q = Q + Qc
U = U + Uc
F = F + Fc
M = M + Mc
end if
TPTSegmentPotential = 1
end if
end do
Q = Q * S%DX
U = U * S%DX
F = F * S%DX
M = M * S%DX
end function TPTSegmentPotential !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!---------------------------------------------------------------------------------------------------
! Forces
!---------------------------------------------------------------------------------------------------
subroutine TPTSegmentForces ( F1, F2, F, M, Laxis, L ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
real(c_double), dimension(0:2), intent(out) :: F1, F2
real(c_double), dimension(0:2), intent(in) :: F, M, Laxis
real(c_double), intent(in) :: L
!-------------------------------------------------------------------------------------------
real(c_double), dimension(0:2) :: MM, FF, FFF
!-------------------------------------------------------------------------------------------
FF = 0.5d+00 * F
MM = M / L
call V3_V3xxV3 ( FFF, MM, Laxis )
F1 = FF - FFF
F2 = FF + FFF
end subroutine TPTSegmentForces !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
integer(c_int) function TPTInteractionF ( Q, U, F1_1, F1_2, F2_1, F2_2, R1_1, R1_2, R2_1, R2_2 )
! This function returns the potential and forces applied to the ends of segments.
!-------------------------------------------------------------------------------------------
real(c_double), intent(out) :: Q, U
real(c_double), dimension(0:2), intent(out) :: F1_1, F1_2, F2_1, F2_2
real(c_double), dimension(0:2), intent(in) :: R1_1, R1_2, R2_1, R2_2
!-------------------------------------------------------------------------------------------
real(c_double), dimension(0:2) :: R1, R2, Laxis1, Laxis2, DR, F1, M1, F2, M2
real(c_double) :: L1, L2
!-------------------------------------------------------------------------------------------
R1 = 0.5 * ( R1_1 + R1_2 )
R2 = 0.5 * ( R2_1 + R2_2 )
Laxis1 = R1_2 - R1_1
Laxis2 = R2_2 - R2_1
L1 = S_V3norm3 ( Laxis1 )
L2 = S_V3norm3 ( Laxis2 )
Laxis1 = Laxis1 / L1
Laxis2 = Laxis2 / L2
DR = R2 - R1
call TPTSetSegPosition1 ( TPTSeg1, R1, Laxis1, L1 )
call TPTSetSegPosition1 ( TPTSeg2, R2, Laxis2, L2 )
TPTInteractionF = TPTSegmentPotential ( Q, U, F1, M1, TPTSeg1, TPTSeg2 )
if ( TPTInteractionF .ne. 1 ) return
call V3_V3xxV3 ( M2, DR, F1 )
F2 = - F1
M2 = - M1 - M2
call TPTSegmentForces ( F1_1, F1_2, F1, M1, Laxis1, L1 )
call TPTSegmentForces ( F2_1, F2_2, F2, M2, Laxis2, L2 )
end function TPTInteractionF !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!---------------------------------------------------------------------------------------------------
! Initialization
!---------------------------------------------------------------------------------------------------
subroutine TPTInit ( R1, R2, NX, NE ) !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
real(c_double), intent(in) :: R1, R2
integer(c_int), intent(in) :: NX, NE
!-------------------------------------------------------------------------------------------
TPTSeg1%X = 0.0d+00
TPTSeg1%Y = 0.0d+00
TPTSeg1%Z = 0.0d+00
TPTSeg1%Psi = 0.0d+00
TPTSeg1%Theta = 0.0d+00
TPTSeg1%Phi = 0.0d+00
TPTSeg1%R = R1
TPTSeg1%NX = NX
TPTSeg1%NE = NE
TPTSeg1%DE = M_2PI / NE
TPTSeg2%X = 0.0d+00
TPTSeg2%Y = 0.0d+00
TPTSeg2%Z = 0.0d+00
TPTSeg2%Psi = 0.0d+00
TPTSeg2%Theta = 0.0d+00
TPTSeg2%Phi = 0.0d+00
TPTSeg2%R = R2
TPTSeg2%NX = NX
TPTSeg2%NE = NE
TPTSeg2%DE = M_2PI / NE
end subroutine TPTInit !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
end module TubePotTrue !****************************************************************************