ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/lj_FF.F90
Revision: 270
Committed: Fri Feb 14 17:08:46 2003 UTC (21 years, 6 months ago) by mmeineke
File size: 21237 byte(s)
Log Message:
added libmdCode and a couple help scripts

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