ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/lj_FF.F90
Revision: 288
Committed: Thu Feb 27 18:42:52 2003 UTC (21 years, 6 months ago) by chuckv
File size: 23050 byte(s)
Log Message:
Changed lj_FF to use new neighbor list module.

File Contents

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