ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/lj_FF.F90
Revision: 247
Committed: Mon Jan 27 18:28:11 2003 UTC (21 years, 5 months ago) by chuckv
File size: 20919 byte(s)
Log Message:
added generic atypes

File Contents

# User Rev Content
1 chuckv 239 !! Calculates Long Range forces Lennard-Jones interactions.
2     !! Corresponds to the force field defined in lj_FF.cpp
3     !! @author Charles F. Vardeman II
4     !! @author Matthew Meineke
5 chuckv 247 !! @version $Id: lj_FF.F90,v 1.11 2003-01-27 18:28:11 chuckv Exp $, $Date: 2003-01-27 18:28:11 $, $Name: not supported by cvs2svn $, $Revision: 1.11 $
6 chuckv 239
7    
8    
9 chuckv 215 module lj_ff
10     use simulation
11 chuckv 246 use definitions
12 chuckv 247 use generic_atypes
13 chuckv 230 #ifdef IS_MPI
14     use mpiSimulation
15     #endif
16 chuckv 215 implicit none
17 chuckv 230 PRIVATE
18 chuckv 215
19 chuckv 230 !! Number of lj_atypes in lj_atype_list
20 chuckv 215 integer, save :: n_lj_atypes = 0
21    
22 chuckv 247 !! Global list of lj atypes in simulation
23     type (lj_atype), pointer :: ljListHead => null()
24     type (lj_atype), pointer :: ljListTail => null()
25 chuckv 230
26 chuckv 215
27 chuckv 230 !! LJ mixing array
28 chuckv 246 type (lj_atype), dimension(:,:), pointer :: ljMixed => null()
29 chuckv 215
30 chuckv 230
31     !! Neighbor list and commom arrays
32 chuckv 239 !! This should eventally get moved to a neighbor list type
33 chuckv 230 integer, allocatable, dimension(:) :: point
34     integer, allocatable, dimension(:) :: list
35 chuckv 239 integer, parameter :: listMultiplier = 80
36     integer :: nListAllocs = 0
37     integer, parameter :: maxListAllocs = 5
38 chuckv 230
39 chuckv 247
40     !! Atype identity pointer lists
41 chuckv 230 #ifdef IS_MPI
42 chuckv 246 !! Row lj_atype pointer list
43     type (lj_atypePtr), dimension(:), pointer :: identPtrListRow => null()
44     !! Column lj_atype pointer list
45     type (lj_atypePtr), dimension(:), pointer :: identPtrListColumn => null()
46 chuckv 247 #else
47     type( lj_identPtrList ), dimension(:), pointer :: identPtrList => null()
48 chuckv 230 #endif
49    
50    
51 chuckv 246 !! Logical has lj force field module been initialized?
52     logical, save :: isljFFinit = .false.
53 chuckv 230
54    
55     !! Public methods and data
56 chuckv 215 public :: new_lj_atype
57 chuckv 230 public :: do_lj_ff
58     public :: getLjPot
59 chuckv 246 public :: init_ljFF
60 chuckv 215
61 chuckv 230
62    
63    
64 chuckv 215 contains
65    
66 chuckv 247 !! Adds a new lj_atype to the list.
67 chuckv 246 subroutine new_lj_atype(ident,mass,epsilon,sigma,status)
68 chuckv 215 real( kind = dp ), intent(in) :: mass
69 chuckv 246 real( kind = dp ), intent(in) :: epsilon
70 chuckv 230 real( kind = dp ), intent(in) :: sigma
71     integer, intent(in) :: ident
72 chuckv 215 integer, intent(out) :: status
73    
74 chuckv 247 type (lj_atype), pointer :: newLJ_atype
75 chuckv 216 integer :: alloc_error
76     integer :: atype_counter = 0
77 chuckv 230 integer :: alloc_size
78 chuckv 247 integer :: err_stat
79 chuckv 215 status = 0
80    
81 chuckv 230
82    
83     ! allocate a new atype
84 chuckv 247 allocate(newLJ_atype,stat=alloc_error)
85 chuckv 216 if (alloc_error /= 0 ) then
86     status = -1
87     return
88     end if
89    
90     ! assign our new lj_atype information
91 chuckv 247 newLJ_atype%mass = mass
92     newLJ_atype%epsilon = epsilon
93     newLJ_atype%sigma = sigma
94     newLJ_atype%sigma2 = sigma * sigma
95     newLJ_atype%sigma6 = newLJ_atype%sigma2 * newLJ_atype%sigma2 &
96     * newLJ_atype%sigma2
97 chuckv 230 ! assume that this atype will be successfully added
98 chuckv 247 newLJ_atype%atype_ident = ident
99     newLJ_atype%atype_number = n_lj_atypes + 1
100 chuckv 216
101 chuckv 247 call add_atype(newLJ_atype,ljListHead,ljListTail,err_stat)
102     if (err_stat /= 0 ) then
103     status = -1
104     return
105 chuckv 230 endif
106    
107     n_lj_atypes = n_lj_atypes + 1
108    
109    
110 chuckv 215 end subroutine new_lj_atype
111    
112    
113 chuckv 239 subroutine init_ljFF(nComponents,ident, status)
114     !! Number of components in ident array
115     integer, intent(inout) :: nComponents
116     !! Array of identities nComponents long corresponding to
117     !! ljatype ident.
118     integer, dimension(nComponents),intent(inout) :: ident
119     !! Result status, success = 0, error = -1
120     integer, intent(out) :: Status
121 chuckv 215
122 chuckv 246 integer :: alloc_stat
123    
124 chuckv 239 integer :: thisStat
125 chuckv 246 integer :: i
126 chuckv 230 #ifdef IS_MPI
127 chuckv 239 integer, allocatable, dimension(:) :: identRow
128     integer, allocatable, dimension(:) :: identCol
129     integer :: nrow
130     integer :: ncol
131     integer :: alloc_stat
132     #endif
133     status = 0
134     !! if were're not in MPI, we just update ljatypePtrList
135     #ifndef IS_MPI
136 chuckv 247 call create_IdentPtrlst(ident,ljListHead,identPtrList,thisStat)
137 chuckv 239 if ( thisStat /= 0 ) then
138     status = -1
139     return
140     endif
141 chuckv 247
142 chuckv 239 !! Allocate pointer lists
143     allocate(point(nComponents),stat=alloc_stat)
144     if (alloc_stat /=0) then
145     status = -1
146     return
147     endif
148    
149     allocate(list(nComponents*listMultiplier),stat=alloc_stat)
150     if (alloc_stat /=0) then
151     status = -1
152     return
153     endif
154 chuckv 230
155 chuckv 239 ! if were're in MPI, we also have to worry about row and col lists
156     #else
157     ! We can only set up forces if mpiSimulation has been setup.
158     if (.not. isMPISimSet()) then
159     status = -1
160     return
161     endif
162     nrow = getNrow()
163     ncol = getNcol()
164     !! Allocate temperary arrays to hold gather information
165     allocate(identRow(nrow),stat=alloc_stat)
166     if (alloc_stat /= 0 ) then
167     status = -1
168     return
169     endif
170 chuckv 215
171 chuckv 239 allocate(identCol(ncol),stat=alloc_stat)
172     if (alloc_stat /= 0 ) then
173     status = -1
174     return
175     endif
176     !! Gather idents into row and column idents
177     call gather(ident,identRow,plan_row)
178     call gather(ident,identCol,plan_col)
179    
180     !! Create row and col pointer lists
181 chuckv 247 call create_IdentPtrlst(identRow,ljListTail,identPtrListRow,thisStat)
182 chuckv 239 if (thisStat /= 0 ) then
183     status = -1
184     return
185     endif
186    
187 chuckv 247 call create_IdentPtrlst(identCol,ljListTail,identPtrListCol,thisStat)
188 chuckv 239 if (thisStat /= 0 ) then
189     status = -1
190     return
191     endif
192    
193     !! free temporary ident arrays
194     deallocate(identCol)
195     deallocate(identRow)
196    
197    
198     !! Allocate neighbor lists for mpi simulations.
199     if (.not. allocated(point)) then
200     allocate(point(nrow),stat=alloc_stat)
201     if (alloc_stat /=0) then
202     status = -1
203     return
204     endif
205    
206     allocate(list(ncol*listMultiplier),stat=alloc_stat)
207     if (alloc_stat /=0) then
208     status = -1
209     return
210     endif
211     else
212     deallocate(list)
213     deallocate(point)
214     allocate(point(nrow),stat=alloc_stat)
215     if (alloc_stat /=0) then
216     status = -1
217     return
218     endif
219    
220     allocate(list(ncol*listMultiplier),stat=alloc_stat)
221     if (alloc_stat /=0) then
222     status = -1
223     return
224     endif
225     endif
226    
227 chuckv 230 #endif
228 chuckv 239
229 chuckv 247 call createMixingList(thisStat)
230     if (thisStat /= 0) then
231     status = -1
232     return
233     endif
234 chuckv 239 isljFFinit = .true.
235    
236 chuckv 230 end subroutine init_ljFF
237    
238 chuckv 239
239    
240    
241    
242 chuckv 230
243 chuckv 247 subroutine createMixingList(status)
244     integer :: listSize
245     integer :: status
246 chuckv 230 integer :: i
247 chuckv 247 integer :: j
248 chuckv 230
249 chuckv 247 integer :: outerCounter = 0
250     integer :: innerCounter = 0
251     type (lj_atype), pointer :: tmpPtrOuter => null()
252     type (lj_atype), pointer :: tmpPtrInner => null()
253     status = 0
254 chuckv 230
255 chuckv 247 listSize = getListLen(ljListHead)
256     if (listSize == 0) then
257     status = -1
258     return
259     end if
260    
261    
262     if (.not. associated(ljMixed)) then
263     allocate(ljMixed(listSize,listSize))
264 chuckv 230 else
265 chuckv 247 status = -1
266     return
267     end if
268 chuckv 230
269 chuckv 247
270     write(*,*) "Creating mixing list"
271     tmpPtrOuter => ljListHead
272     tmpPtrInner => tmpPtrOuter%next
273     do while (associated(tmpPtrOuter))
274     outerCounter = outerCounter + 1
275     write(*,*) "Outer index is: ", outerCounter
276     write(*,*) "For atype: ", tmpPtrOuter%atype_ident
277     ! do self mixing rule
278     ljMixed(outerCounter,outerCounter)%sigma = tmpPtrOuter%sigma
279    
280     ljMixed(outerCounter,outerCounter)%sigma2 = ljMixed(outerCounter,outerCounter)%sigma &
281     * ljMixed(outerCounter,outerCounter)%sigma
282    
283     ljMixed(outerCounter,outerCounter)%sigma6 = ljMixed(outerCounter,outerCounter)%sigma2 &
284     * ljMixed(outerCounter,outerCounter)%sigma2 &
285     * ljMixed(outerCounter,outerCounter)%sigma2
286    
287     ljMixed(outerCounter,outerCounter)%epsilon = tmpPtrOuter%epsilon
288    
289     innerCounter = outerCounter + 1
290     do while (associated(tmpPtrInner))
291    
292     write(*,*) "Entered inner loop for: ", innerCounter
293     ljMixed(outerCounter,innerCounter)%sigma = &
294     calcLJMix("sigma",tmpPtrOuter%sigma, &
295     tmpPtrInner%sigma)
296    
297     ljMixed(outerCounter,innerCounter)%sigma2 = &
298     ljMixed(outerCounter,innerCounter)%sigma &
299     * ljMixed(outerCounter,innerCounter)%sigma
300    
301     ljMixed(outerCounter,innerCounter)%sigma6 = &
302     ljMixed(outerCounter,innerCounter)%sigma2 &
303     * ljMixed(outerCounter,innerCounter)%sigma2 &
304     * ljMixed(outerCounter,innerCounter)%sigma2
305    
306     ljMixed(outerCounter,innerCounter)%epsilon = &
307     calcLJMix("epsilon",tmpPtrOuter%epsilon, &
308     tmpPtrInner%epsilon)
309     ljMixed(innerCounter,outerCounter)%sigma = ljMixed(outerCounter,innerCounter)%sigma
310     ljMixed(innerCounter,outerCounter)%sigma2 = ljMixed(outerCounter,innerCounter)%sigma2
311     ljMixed(innerCounter,outerCounter)%sigma6 = ljMixed(outerCounter,innerCounter)%sigma6
312     ljMixed(innerCounter,outerCounter)%epsilon = ljMixed(outerCounter,innerCounter)%epsilon
313    
314    
315     tmpPtrInner => tmpPtrInner%next
316     innerCounter = innerCounter + 1
317     end do
318     ! advance pointers
319     tmpPtrOuter => tmpPtrOuter%next
320     if (associated(tmpPtrOuter)) then
321     tmpPtrInner => tmpPtrOuter%next
322 chuckv 230 endif
323    
324     end do
325    
326 chuckv 247 end subroutine createMixingList
327 chuckv 230
328    
329 chuckv 246
330 chuckv 230
331    
332 chuckv 246
333    
334    
335 chuckv 239 !! FORCE routine Calculates Lennard Jones forces.
336 chuckv 230 !------------------------------------------------------------->
337     subroutine do_lj_ff(q,f,potE,do_pot)
338 chuckv 246 !! Position array provided by C, dimensioned by getNlocal
339     real ( kind = dp ), dimension(3,getNlocal()) :: q
340     !! Force array provided by C, dimensioned by getNlocal
341     real ( kind = dp ), dimension(3,getNlocal()) :: f
342 chuckv 222 real ( kind = dp ) :: potE
343 mmeineke 235 logical ( kind = 2) :: do_pot
344 chuckv 239
345     type(lj_atype), pointer :: ljAtype_i
346     type(lj_atype), pointer :: ljAtype_j
347    
348 chuckv 230 #ifdef MPI
349 chuckv 246 real( kind = DP ), dimension(3,getNcol()) :: efr
350 chuckv 230 real( kind = DP ) :: pot_local
351     #else
352 chuckv 246 real( kind = DP ), dimension(3,getNlocal()) :: efr
353 chuckv 230 #endif
354    
355 chuckv 247 !! Local arrays needed for MPI
356     #ifdef IS_MPI
357     real(kind = dp), dimension(3:getNrow()) :: qRow
358     real(kind = dp), dimension(3:getNcol()) :: qCol
359    
360     real(kind = dp), dimension(3:getNrow()) :: fRow
361     real(kind = dp), dimension(3:getNcol()) :: fCol
362    
363     real(kind = dp), dimension(getNrow()) :: eRow
364     real(kind = dp), dimension(getNcol()) :: eCol
365    
366     real(kind = dp), dimension(getNlocal()) :: eTemp
367     #endif
368    
369    
370    
371 chuckv 230 real( kind = DP ) :: pe
372 chuckv 246 logical :: update_nlist
373 chuckv 216
374 chuckv 232
375 chuckv 230 integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
376     integer :: nlist
377     integer :: j_start
378     integer :: tag_i,tag_j
379     real( kind = DP ) :: r, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
380     real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq
381 chuckv 246 real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut
382 chuckv 222
383 chuckv 246 ! a rig that need to be fixed.
384     #ifdef IS_MPI
385     logical :: newtons_thrd = .true.
386     real( kind = dp ) :: pe_local
387     integer :: nlocal
388     #endif
389 chuckv 230 integer :: nrow
390     integer :: ncol
391 chuckv 246 integer :: natoms
392 chuckv 222
393 chuckv 246 ! Make sure we are properly initialized.
394 chuckv 239 if (.not. isljFFInit) then
395     write(default_error,*) "ERROR: lj_FF has not been properly initialized"
396     return
397     endif
398 chuckv 246 #ifdef IS_MPI
399     if (.not. isMPISimSet()) then
400     write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
401     return
402     endif
403     #endif
404 chuckv 222
405 chuckv 246 !! initialize local variables
406     natoms = getNlocal()
407     call getRcut(rcut,rcut2=rcutsq)
408     call getRlist(rlist,rlistsq)
409 chuckv 230 #ifndef IS_MPI
410     nrow = natoms - 1
411     ncol = natoms
412     #else
413 chuckv 239 nrow = getNrow(plan_row)
414     ncol = getNcol(plan_col)
415 chuckv 246 nlocal = natoms
416 chuckv 230 j_start = 1
417     #endif
418 chuckv 222
419 chuckv 239
420 chuckv 246 !! See if we need to update neighbor lists
421     call check(q,update_nlist)
422 chuckv 239
423     !--------------WARNING...........................
424     ! Zero variables, NOTE:::: Forces are zeroed in C
425     ! Zeroing them here could delete previously computed
426     ! Forces.
427     !------------------------------------------------
428 chuckv 230 #ifndef IS_MPI
429 chuckv 246 ! nloops = nloops + 1
430     pe = 0.0E0_DP
431    
432 chuckv 230 #else
433     f_row = 0.0E0_DP
434     f_col = 0.0E0_DP
435    
436 chuckv 246 pe_local = 0.0E0_DP
437 chuckv 230
438 chuckv 246 eRow = 0.0E0_DP
439     eCol = 0.0E0_DP
440     eTemp = 0.0E0_DP
441 chuckv 230 #endif
442     efr = 0.0E0_DP
443    
444     ! communicate MPI positions
445     #ifdef MPI
446 chuckv 239 call gather(q,qRow,plan_row3)
447     call gather(q,qCol,plan_col3)
448 chuckv 230 #endif
449    
450    
451     if (update_nlist) then
452    
453     ! save current configuration, contruct neighbor list,
454     ! and calculate forces
455 chuckv 246 call save_nlist(q)
456 chuckv 230
457     nlist = 0
458    
459    
460    
461     do i = 1, nrow
462     point(i) = nlist + 1
463     #ifdef MPI
464 chuckv 239 ljAtype_i => identPtrListRow(i)%this
465     tag_i = tagRow(i)
466     rxi = qRow(1,i)
467     ryi = qRow(2,i)
468     rzi = qRow(3,i)
469 chuckv 230 #else
470 chuckv 239 ljAtype_i => identPtrList(i)%this
471 chuckv 230 j_start = i + 1
472     rxi = q(1,i)
473     ryi = q(2,i)
474     rzi = q(3,i)
475     #endif
476    
477     inner: do j = j_start, ncol
478     #ifdef MPI
479 chuckv 239 ! Assign identity pointers and tags
480     ljAtype_j => identPtrListColumn(j)%this
481     tag_j = tagCol(j)
482 chuckv 230 if (newtons_thrd) then
483     if (tag_i <= tag_j) then
484     if (mod(tag_i + tag_j,2) == 0) cycle inner
485     else
486     if (mod(tag_i + tag_j,2) == 1) cycle inner
487     endif
488     endif
489    
490 chuckv 239 rxij = wrap(rxi - qCol(1,j), 1)
491     ryij = wrap(ryi - qCol(2,j), 2)
492     rzij = wrap(rzi - qCol(3,j), 3)
493 chuckv 230 #else
494 chuckv 239 ljAtype_j => identPtrList(j)%this
495 chuckv 230 rxij = wrap(rxi - q(1,j), 1)
496     ryij = wrap(ryi - q(2,j), 2)
497     rzij = wrap(rzi - q(3,j), 3)
498    
499     #endif
500     rijsq = rxij*rxij + ryij*ryij + rzij*rzij
501    
502     #ifdef MPI
503 chuckv 246 if (rijsq <= rlistsq .AND. &
504 chuckv 230 tag_j /= tag_i) then
505     #else
506 chuckv 246 if (rijsq < rlistsq) then
507 chuckv 230 #endif
508    
509     nlist = nlist + 1
510     if (nlist > size(list)) then
511 chuckv 246 !! "Change how nlist size is done"
512 chuckv 239 write(DEFAULT_ERROR,*) "ERROR: nlist > list size"
513 chuckv 230 endif
514     list(nlist) = j
515    
516    
517     if (rijsq < rcutsq) then
518    
519     r = dsqrt(rijsq)
520    
521 chuckv 239 call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j)
522 chuckv 230
523     #ifdef MPI
524 chuckv 246 eRow(i) = eRow(i) + pot*0.5
525     eCol(i) = eCol(i) + pot*0.5
526 chuckv 230 #else
527     pe = pe + pot
528     #endif
529    
530     efr(1,j) = -rxij
531     efr(2,j) = -ryij
532     efr(3,j) = -rzij
533    
534     do dim = 1, 3
535    
536    
537     drdx1 = efr(dim,j) / r
538     ftmp = dudr * drdx1
539    
540    
541     #ifdef MPI
542 chuckv 239 fCol(dim,j) = fCol(dim,j) - ftmp
543     fRow(dim,i) = fRow(dim,i) + ftmp
544 chuckv 230 #else
545    
546     f(dim,j) = f(dim,j) - ftmp
547     f(dim,i) = f(dim,i) + ftmp
548    
549     #endif
550     enddo
551     endif
552     endif
553     enddo inner
554     enddo
555    
556     #ifdef MPI
557     point(nrow + 1) = nlist + 1
558     #else
559     point(natoms) = nlist + 1
560     #endif
561    
562     else
563    
564     ! use the list to find the neighbors
565     do i = 1, nrow
566     JBEG = POINT(i)
567     JEND = POINT(i+1) - 1
568     ! check thiat molecule i has neighbors
569     if (jbeg .le. jend) then
570     #ifdef MPI
571 chuckv 239 ljAtype_i => identPtrListRow(i)%this
572     rxi = qRow(1,i)
573     ryi = qRow(2,i)
574     rzi = qRow(3,i)
575 chuckv 230 #else
576 chuckv 239 ljAtype_i => identPtrList(i)%this
577 chuckv 230 rxi = q(1,i)
578     ryi = q(2,i)
579     rzi = q(3,i)
580     #endif
581     do jnab = jbeg, jend
582     j = list(jnab)
583     #ifdef MPI
584 chuckv 239 ljAtype_j = identPtrListColumn(j)%this
585     rxij = wrap(rxi - q_col(1,j), 1)
586     ryij = wrap(ryi - q_col(2,j), 2)
587     rzij = wrap(rzi - q_col(3,j), 3)
588 chuckv 230 #else
589 chuckv 239 ljAtype_j = identPtrList(j)%this
590     rxij = wrap(rxi - q(1,j), 1)
591     ryij = wrap(ryi - q(2,j), 2)
592     rzij = wrap(rzi - q(3,j), 3)
593 chuckv 230 #endif
594     rijsq = rxij*rxij + ryij*ryij + rzij*rzij
595    
596     if (rijsq < rcutsq) then
597    
598     r = dsqrt(rijsq)
599    
600 chuckv 239 call getLJPot(r,pot,dudr,ljAtype_i,ljAtype_j)
601 chuckv 230 #ifdef MPI
602 chuckv 246 eRow(i) = eRow(i) + pot*0.5
603     eCol(i) = eCol(i) + pot*0.5
604 chuckv 230 #else
605 chuckv 246 pe = pe + pot
606 chuckv 230 #endif
607    
608    
609     efr(1,j) = -rxij
610     efr(2,j) = -ryij
611     efr(3,j) = -rzij
612    
613     do dim = 1, 3
614    
615     drdx1 = efr(dim,j) / r
616     ftmp = dudr * drdx1
617     #ifdef MPI
618 chuckv 239 fCol(dim,j) = fCol(dim,j) - ftmp
619     fRow(dim,i) = fRow(dim,i) + ftmp
620 chuckv 230 #else
621     f(dim,j) = f(dim,j) - ftmp
622     f(dim,i) = f(dim,i) + ftmp
623     #endif
624     enddo
625     endif
626     enddo
627     endif
628     enddo
629     endif
630    
631    
632    
633     #ifdef MPI
634     !!distribute forces
635 chuckv 239 call scatter(fRow,f,plan_row3)
636 chuckv 230
637 chuckv 239 call scatter(fCol,f_tmp,plan_col3)
638 chuckv 246
639 chuckv 230 do i = 1,nlocal
640 chuckv 246 f(1:3,i) = f(1:3,i) + f_tmp(1:3,i)
641 chuckv 230 end do
642    
643    
644    
645     if (do_pot) then
646     ! scatter/gather pot_row into the members of my column
647 chuckv 232 call scatter(e_row,e_tmp,plan_row)
648 chuckv 230
649     ! scatter/gather pot_local into all other procs
650     ! add resultant to get total pot
651     do i = 1, nlocal
652 chuckv 246 pe_local = pe_local + e_tmp(i)
653 chuckv 230 enddo
654     if (newtons_thrd) then
655     e_tmp = 0.0E0_DP
656     call scatter(e_col,e_tmp,plan_col)
657     do i = 1, nlocal
658 chuckv 246 pe_local = pe_local + e_tmp(i)
659 chuckv 230 enddo
660     endif
661     endif
662 chuckv 246 potE = pe_local
663 chuckv 230 #endif
664    
665 chuckv 246 potE = pe
666 chuckv 230
667    
668 chuckv 222 end subroutine do_lj_ff
669    
670 chuckv 239 !! Calculates the potential between two lj particles based on two lj_atype pointers, optionally returns second
671 chuckv 230 !! derivatives.
672     subroutine getLjPot(r,pot,dudr,atype1,atype2,d2,status)
673     ! arguments
674     !! Length of vector between particles
675     real( kind = dp ), intent(in) :: r
676     !! Potential Energy
677     real( kind = dp ), intent(out) :: pot
678     !! Derivatve wrt postion
679     real( kind = dp ), intent(out) :: dudr
680     !! Second Derivative, optional, used mainly for normal mode calculations.
681     real( kind = dp ), intent(out), optional :: d2
682    
683 chuckv 246 type (lj_atype), pointer :: atype1
684     type (lj_atype), pointer :: atype2
685 chuckv 222
686 chuckv 230 integer, intent(out), optional :: status
687 chuckv 222
688 chuckv 230 ! local Variables
689     real( kind = dp ) :: sigma
690     real( kind = dp ) :: sigma2
691     real( kind = dp ) :: sigma6
692 chuckv 246 real( kind = dp ) :: epsilon
693 chuckv 230
694     real( kind = dp ) :: rcut
695     real( kind = dp ) :: rcut2
696     real( kind = dp ) :: rcut6
697     real( kind = dp ) :: r2
698     real( kind = dp ) :: r6
699    
700     real( kind = dp ) :: t6
701     real( kind = dp ) :: t12
702     real( kind = dp ) :: tp6
703     real( kind = dp ) :: tp12
704     real( kind = dp ) :: delta
705    
706     logical :: doSec = .false.
707    
708     integer :: errorStat
709    
710     !! Optional Argument Checking
711     ! Check to see if we need to do second derivatives
712    
713     if (present(d2)) doSec = .true.
714     if (present(status)) status = 0
715    
716     ! Look up the correct parameters in the mixing matrix
717 chuckv 246 sigma = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma
718     sigma2 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma2
719     sigma6 = ljMixed(atype1%atype_ident,atype2%atype_ident)%sigma6
720     epsilon = ljMixed(atype1%atype_ident,atype2%atype_ident)%epsilon
721 chuckv 230
722    
723    
724     call getRcut(rcut,rcut2=rcut2,rcut6=rcut6,status=errorStat)
725    
726     r2 = r * r
727     r6 = r2 * r2 * r2
728    
729     t6 = sigma6/ r6
730     t12 = t6 * t6
731    
732    
733    
734     tp6 = sigma6 / rcut6
735     tp12 = tp6*tp6
736    
737     delta = -4.0E0_DP*epsilon * (tp12 - tp6)
738    
739     if (r.le.rcut) then
740 chuckv 246 pot = 4.0E0_DP * epsilon * (t12 - t6) + delta
741 chuckv 230 dudr = 24.0E0_DP * epsilon * (t6 - 2.0E0_DP*t12) / r
742     if(doSec) d2 = 24.0E0_DP * epsilon * (26.0E0_DP*t12 - 7.0E0_DP*t6)/r/r
743     else
744 chuckv 246 pot = 0.0E0_DP
745 chuckv 230 dudr = 0.0E0_DP
746     if(doSec) d2 = 0.0E0_DP
747     endif
748    
749     return
750    
751    
752    
753     end subroutine getLjPot
754    
755    
756 chuckv 239 !! Calculates the mixing for sigma or epslon based on character string for initialzition of mixing array.
757 chuckv 230 function calcLJMix(thisParam,param1,param2,status) result(myMixParam)
758     character(len=*) :: thisParam
759     real(kind = dp) :: param1
760     real(kind = dp) :: param2
761     real(kind = dp ) :: myMixParam
762     integer, optional :: status
763    
764    
765     myMixParam = 0.0_dp
766    
767     if (present(status)) status = 0
768    
769     select case (thisParam)
770    
771     case ("sigma")
772     myMixParam = 0.5_dp * (param1 + param2)
773 chuckv 246 case ("epsilon")
774 chuckv 230 myMixParam = sqrt(param1 * param2)
775     case default
776     status = -1
777     end select
778    
779     end function calcLJMix
780    
781    
782    
783 chuckv 215 end module lj_ff