ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/mdtools/md_code/ssd_FF.F90
Revision: 265
Committed: Wed Feb 12 16:21:26 2003 UTC (21 years, 4 months ago) by chuckv
File size: 28684 byte(s)
Log Message:
Initial commit of ssd force module.

File Contents

# User Rev Content
1 chuckv 265 !! This Module Calculates forces due to SSD potential and VDW interactions
2     !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
3    
4     !! This module contains the Public procedures:
5    
6    
7     !! Corresponds to the force field defined in ssd_FF.cpp
8     !! @author Charles F. Vardeman II
9     !! @author Matthew Meineke
10     !! @version $Id: ssd_FF.F90,v 1.1 2003-02-12 16:21:26 chuckv Exp $, $Date: 2003-02-12 16:21:26 $, $Name: not supported by cvs2svn $, $Revision: 1.1 $
11    
12    
13    
14     module ssd_ff
15     use simulation
16     use definitions
17     use generic_atypes
18     #ifdef IS_MPI
19     use mpiSimulation
20     #endif
21     implicit none
22     PRIVATE
23    
24     !! Number of lj_atypes in lj_atype_list
25     integer, save :: n_atypes = 0
26    
27     !! Global list of lj atypes in simulation
28     type (ssd_atype), pointer :: ListHead => null()
29     type (ssd_atype), pointer :: ListTail => null()
30    
31    
32     !! LJ mixing array
33     type (ssd_atype), dimension(:,:), pointer :: Mixed => null()
34    
35    
36     !! Neighbor list and commom arrays
37     !! This should eventally get moved to a neighbor list type
38     integer, allocatable, dimension(:) :: point
39     integer, allocatable, dimension(:) :: list
40     integer, parameter :: listMultiplier = 80
41     integer :: nListAllocs = 0
42     integer, parameter :: maxListAllocs = 5
43    
44     logical, save :: firstTime = .True.
45    
46     !! Atype identity pointer lists
47     #ifdef IS_MPI
48     !! Row lj_atype pointer list
49     type (ssd_identPtrList), dimension(:), pointer :: identPtrListRow => null()
50     !! Column lj_atype pointer list
51     type (ssd_identPtrList), dimension(:), pointer :: identPtrListColumn => null()
52     #else
53     type( ssd_identPtrList ), dimension(:), pointer :: identPtrList => null()
54     #endif
55    
56    
57     !! Logical has lj force field module been initialized?
58     logical, save :: isFFinit = .false.
59    
60    
61     !! Public methods and data
62     public :: new_ssd_atype
63     public :: do_SSD_ff
64     public :: getSSDForce
65     public :: init_ssdFF
66    
67    
68    
69    
70     contains
71    
72     !! Adds a new lj_atype to the list.
73     subroutine new_ssd_atype(ident,mass,epsilon,sigma,dipoleMoment,&
74     hasDipole,v0,w0,status)
75     real( kind = dp ), intent(in) :: mass
76     real( kind = dp ), intent(in) :: epsilon
77     real( kind = dp ), intent(in) :: sigma
78     real( kind = dp ), intent(in) :: dipoleMoment
79     logical :: hasDipole
80     real( kind = dp ), intent(in) :: w0
81     real( kind = dp ), intent(in) :: v0
82     integer, intent(in) :: ident
83     integer, intent(out) :: status
84    
85     type (ssd_atype), pointer :: newSSD_atype
86     integer :: alloc_error
87     integer :: atype_counter = 0
88     integer :: alloc_size
89     integer :: err_stat
90     status = 0
91    
92    
93    
94     ! allocate a new atype
95     allocate(newSSD_atype,stat=alloc_error)
96     if (alloc_error /= 0 ) then
97     status = -1
98     return
99     end if
100    
101     ! assign our new lj_atype information
102     newSSD_atype%lj_parms%mass = mass
103     newSSD_atype%lj_parms%epsilon = epsilon
104     newSSD_atype%lj_parms%sigma = sigma
105     newSSD_atype%lj_parms%sigma2 = sigma * sigma
106     newSSD_atype%lj_parms%sigma6 = newLJ_atype%sigma2 * newLJ_atype%sigma2 &
107     * newLJ_atype%sigma2
108     newSSD_atype%dipoleMoment = dipoleMoment
109     newSSD_atype&hasDipole = hasDipole
110     newSSD_atype%w0 = w0
111     newSSD_atype%v0 = v0
112     ! assume that this atype will be successfully added
113     newSSD_atype%atype_ident = ident
114     newSSD_atype%atype_number = n_SSD_atypes + 1
115    
116     call add_atype(newSSD_atype,ListHead,ListTail,err_stat)
117     if (err_stat /= 0 ) then
118     status = -1
119     return
120     endif
121    
122     n_ssd_atypes = n_ssd_atypes + 1
123    
124    
125     end subroutine new_ssd_atype
126    
127    
128     subroutine init_ssdFF(nComponents,ident, status)
129     !! Number of components in ident array
130     integer, intent(inout) :: nComponents
131     !! Array of identities nComponents long corresponding to
132     !! ljatype ident.
133     integer, dimension(nComponents),intent(inout) :: ident
134     !! Result status, success = 0, error = -1
135     integer, intent(out) :: Status
136    
137     integer :: alloc_stat
138    
139     integer :: thisStat
140     integer :: i
141    
142     integer :: myNode
143     #ifdef IS_MPI
144     integer, allocatable, dimension(:) :: identRow
145     integer, allocatable, dimension(:) :: identCol
146     integer :: nrow
147     integer :: ncol
148     #endif
149     status = 0
150    
151    
152    
153    
154     !! if were're not in MPI, we just update ljatypePtrList
155     #ifndef IS_MPI
156     call create_IdentPtrlst(ident,ListHead,identPtrList,thisStat)
157     if ( thisStat /= 0 ) then
158     status = -1
159     return
160     endif
161    
162     !! Allocate pointer lists
163     allocate(point(nComponents),stat=alloc_stat)
164     if (alloc_stat /=0) then
165     status = -1
166     return
167     endif
168    
169     allocate(list(nComponents*listMultiplier),stat=alloc_stat)
170     if (alloc_stat /=0) then
171     status = -1
172     return
173     endif
174    
175     ! if were're in MPI, we also have to worry about row and col lists
176     #else
177    
178     ! We can only set up forces if mpiSimulation has been setup.
179     if (.not. isMPISimSet()) then
180     write(default_error,*) "MPI is not set"
181     status = -1
182     return
183     endif
184     nrow = getNrow(plan_row)
185     ncol = getNcol(plan_col)
186     mynode = getMyNode()
187     !! Allocate temperary arrays to hold gather information
188     allocate(identRow(nrow),stat=alloc_stat)
189     if (alloc_stat /= 0 ) then
190     status = -1
191     return
192     endif
193    
194     allocate(identCol(ncol),stat=alloc_stat)
195     if (alloc_stat /= 0 ) then
196     status = -1
197     return
198     endif
199    
200     !! Gather idents into row and column idents
201    
202     call gather(ident,identRow,plan_row)
203     call gather(ident,identCol,plan_col)
204    
205    
206     !! Create row and col pointer lists
207    
208     call create_IdentPtrlst(identRow,ListHead,identPtrListRow,thisStat)
209     if (thisStat /= 0 ) then
210     status = -1
211     return
212     endif
213    
214     call create_IdentPtrlst(identCol,ListHead,identPtrListColumn,thisStat)
215     if (thisStat /= 0 ) then
216     status = -1
217     return
218     endif
219    
220     !! free temporary ident arrays
221     if (allocated(identCol)) then
222     deallocate(identCol)
223     end if
224     if (allocated(identCol)) then
225     deallocate(identRow)
226     endif
227    
228     !! Allocate neighbor lists for mpi simulations.
229     if (.not. allocated(point)) then
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     else
242     deallocate(list)
243     deallocate(point)
244     allocate(point(nrow),stat=alloc_stat)
245     if (alloc_stat /=0) then
246     status = -1
247     return
248     endif
249    
250     allocate(list(ncol*listMultiplier),stat=alloc_stat)
251     if (alloc_stat /=0) then
252     status = -1
253     return
254     endif
255     endif
256    
257     #endif
258    
259     call createMixingList(thisStat)
260     if (thisStat /= 0) then
261     status = -1
262     return
263     endif
264     isFFinit = .true.
265    
266    
267     end subroutine init_ssdFF
268    
269    
270    
271    
272    
273    
274     subroutine createMixingList(status)
275     integer :: listSize
276     integer :: status
277     integer :: i
278     integer :: j
279    
280     integer :: outerCounter = 0
281     integer :: innerCounter = 0
282     type (lj_atype), pointer :: tmpPtrOuter => null()
283     type (lj_atype), pointer :: tmpPtrInner => null()
284     status = 0
285    
286     listSize = getListLen(ListHead)
287     if (listSize == 0) then
288     status = -1
289     return
290     end if
291    
292    
293     if (.not. associated(Mixed)) then
294     allocate(Mixed(listSize,listSize))
295     else
296     status = -1
297     return
298     end if
299    
300    
301    
302     tmpPtrOuter => ListHead
303     tmpPtrInner => tmpPtrOuter%next
304     do while (associated(tmpPtrOuter))
305     outerCounter = outerCounter + 1
306     ! do self mixing rule
307     Mixed(outerCounter,outerCounter)%sigma = tmpPtrOuter%sigma
308    
309     Mixed(outerCounter,outerCounter)%sigma2 = Mixed(outerCounter,outerCounter)%sigma &
310     * Mixed(outerCounter,outerCounter)%sigma
311    
312     Mixed(outerCounter,outerCounter)%sigma6 = Mixed(outerCounter,outerCounter)%sigma2 &
313     * Mixed(outerCounter,outerCounter)%sigma2 &
314     * Mixed(outerCounter,outerCounter)%sigma2
315    
316     Mixed(outerCounter,outerCounter)%epsilon = tmpPtrOuter%epsilon
317    
318     innerCounter = outerCounter + 1
319     do while (associated(tmpPtrInner))
320    
321     Mixed(outerCounter,innerCounter)%sigma = &
322     calcLJMix("sigma",tmpPtrOuter%sigma, &
323     tmpPtrInner%sigma)
324    
325     Mixed(outerCounter,innerCounter)%sigma2 = &
326     Mixed(outerCounter,innerCounter)%sigma &
327     * Mixed(outerCounter,innerCounter)%sigma
328    
329     Mixed(outerCounter,innerCounter)%sigma6 = &
330     Mixed(outerCounter,innerCounter)%sigma2 &
331     * Mixed(outerCounter,innerCounter)%sigma2 &
332     * Mixed(outerCounter,innerCounter)%sigma2
333    
334     Mixed(outerCounter,innerCounter)%epsilon = &
335     calcLJMix("epsilon",tmpPtrOuter%epsilon, &
336     tmpPtrInner%epsilon)
337     Mixed(innerCounter,outerCounter)%sigma = Mixed(outerCounter,innerCounter)%sigma
338     Mixed(innerCounter,outerCounter)%sigma2 = Mixed(outerCounter,innerCounter)%sigma2
339     Mixed(innerCounter,outerCounter)%sigma6 = Mixed(outerCounter,innerCounter)%sigma6
340     Mixed(innerCounter,outerCounter)%epsilon = Mixed(outerCounter,innerCounter)%epsilon
341    
342    
343     tmpPtrInner => tmpPtrInner%next
344     innerCounter = innerCounter + 1
345     end do
346     ! advance pointers
347     tmpPtrOuter => tmpPtrOuter%next
348     if (associated(tmpPtrOuter)) then
349     tmpPtrInner => tmpPtrOuter%next
350     endif
351    
352     end do
353    
354     end subroutine createMixingList
355    
356    
357    
358    
359    
360    
361    
362    
363     !! FORCE routine Calculates Lennard Jones forces.
364     !------------------------------------------------------------->
365     subroutine do_ssd_ff(q,f,t,u,A,potE,do_pot)
366     !! Position array provided by C, dimensioned by getNlocal
367     real ( kind = dp ), dimension(3,getNlocal()) :: q
368     !! Force array given to C, dimensioned by getNlocal
369     real ( kind = dp ), dimension(3,getNlocal()) :: f
370     !! Torsion array given to C, dimensioned by getNlocal
371     real ( kind = dp ), dimension(3,getNlocal()) :: t
372     !! Dipole unit vectors given to C, dimensioned by getNlocal
373     real ( kind = dp ), dimension(3,getNlocal()) :: u
374    
375    
376    
377     !! rotational array -- 9 element matrix
378     !! 1 - Axx
379     !! 2 - Axy
380     !! 3 - Axz
381     !! 4 - Ayx
382     !! 5 - Ayy
383     !! 6 - Ayz
384     !! 7 - Azx
385     !! 8 - Azy
386     !! 9 - Azz
387     real( kind = dp ), dimension(9,getNlocal()) :: A
388    
389     real ( kind = dp ) :: potE
390     logical ( kind = 2) :: do_pot
391    
392     type(ssd_atype), pointer :: Atype_i
393     type(ssd_atype), pointer :: Atype_j
394    
395    
396     !! Force from SSD due to (ndim,i-j(1) and j-i(2))
397     real( kind = dp ), dimension(3,2) :: fl_ij_ji = 0.0_dp
398     !! Torsion from SSD due to (ndim,i-j(1) and j-i(2))
399     real( kind = dp ), dimension(3,2) :: tl_ij_ji = 0.0_dp
400    
401     #ifdef IS_MPI
402     real( kind = DP ), dimension(3,getNcol(plan_col)) :: efr = 0.0_dp
403     real( kind = DP ) :: pot_local = 0.0_dp
404     #else
405     real( kind = DP ), dimension(3,getNlocal()) :: efr = 0.0_dp
406     #endif
407    
408     !! Local arrays needed for MPI
409     #ifdef IS_MPI
410     !! Position arrays
411     real(kind = dp), dimension(3,getNrow(plan_row)) :: qRow = 0.0_dp
412     real(kind = dp), dimension(3,getNcol(plan_col)) :: qCol = 0.0_dp
413    
414     !! Rotational Arrays -- note 9 element matrix same layout as A
415     real(kind = dp), dimension(9,getNrow(plan_row)) :: ARow = 0.0_dp
416     real(kind = dp), dimension(9,getNcol(plan_col)) :: ACol = 0.0_dp
417    
418     !! Force Arrays
419     real(kind = dp), dimension(3,getNrow(plan_row)) :: fRow = 0.0_dp
420     real(kind = dp), dimension(3,getNcol(plan_col)) :: fCol = 0.0_dp
421     real(kind = dp), dimension(3,getNlocal()) :: fMPITemp = 0.0_dp
422     !! Torsion arrays
423     real(kind = dp), dimension(3,getNrow(plan_row)) :: tRow = 0.0_dp
424     real(kind = dp), dimension(3,getNcol(plan_col)) :: tCol = 0.0_dp
425     real(kind = dp), dimension(3,getNlocal()) :: tMPITemp = 0.0_dp
426    
427    
428    
429     real(kind = dp), dimension(getNrow(plan_row)) :: eRow = 0.0_dp
430     real(kind = dp), dimension(getNcol(plan_col)) :: eCol = 0.0_dp
431    
432     real(kind = dp), dimension(getNlocal()) :: eTemp = 0.0_dp
433     #endif
434    
435     real( kind = DP ) :: pe = 0.0_dp
436     logical :: update_nlist
437    
438    
439     integer :: i, j, jbeg, jend, jnab, idim, jdim, idim2, jdim2, dim, dim2
440     integer :: nlist
441     integer :: j_start
442     integer :: tag_i,tag_j
443     real( kind = DP ) :: r, pot, ftmp, dudr, d2, drdx1, kt1, kt2, kt3, ktmp
444     real( kind = DP ) :: rxi, ryi, rzi, rxij, ryij, rzij, rijsq
445     real( kind = DP ) :: rlistsq, rcutsq,rlist,rcut
446    
447     ! a rig that need to be fixed.
448     #ifdef IS_MPI
449     logical :: newtons_thrd = .true.
450     real( kind = dp ) :: pe_local = 0.0_dp
451     #endif
452     integer :: nrow
453     integer :: ncol
454     integer :: nlocal
455    
456    
457    
458    
459     ! Make sure we are properly initialized.
460     if (.not. isljFFInit) then
461     write(default_error,*) "ERROR: lj_FF has not been properly initialized"
462     return
463     endif
464     #ifdef IS_MPI
465     if (.not. isMPISimSet()) then
466     write(default_error,*) "ERROR: mpiSimulation has not been properly initialized"
467     return
468     endif
469     #endif
470    
471     !! initialize local variables
472     nlocal = getNlocal()
473     call getRcut(rcut,rcut2=rcutsq)
474     call getRlist(rlist,rlistsq)
475    
476     #ifndef IS_MPI
477     nrow = nlocal - 1
478     ncol = nlocal
479     #else
480     nrow = getNrow(plan_row)
481     ncol = getNcol(plan_col)
482     j_start = 1
483     #endif
484    
485    
486     !! See if we need to update neighbor lists
487     call check(q,update_nlist)
488     if (firstTime) then
489     update_nlist = .true.
490     firstTime = .false.
491     endif
492    
493    
494     !--------------WARNING...........................
495     ! Zero variables, NOTE:::: Forces are zeroed in C
496     ! Zeroing them here could delete previously computed
497     ! Forces.
498     !------------------------------------------------
499     #ifndef IS_MPI
500     ! nloops = nloops + 1
501     pe = 0.0E0_DP
502    
503     #else
504     fRow = 0.0E0_DP
505     fCol = 0.0E0_DP
506    
507     pe_local = 0.0E0_DP
508    
509     eRow = 0.0E0_DP
510     eCol = 0.0E0_DP
511     eTemp = 0.0E0_DP
512     #endif
513     efr = 0.0E0_DP
514    
515    
516     #ifdef IS_MPI
517     ! communicate local MPI position arrays into row and column arrays.
518     call gather(q,qRow,plan_row3d)
519     call gather(q,qCol,plan_col3d)
520     ! communicate local MPI Rotational Array into row and column arrays.
521     call gather(A,ARow,plan_row_Rotation)
522     call gather(A,ACol,plan_col_Rotation)
523     #endif
524    
525    
526     if (update_nlist) then
527    
528     ! save current configuration, contruct neighbor list,
529     ! and calculate forces
530     call save_nlist(q)
531    
532     nlist = 0
533    
534    
535    
536     do i = 1, nrow
537     point(i) = nlist + 1
538     #ifdef IS_MPI
539     Atype_i => identPtrListRow(i)%this
540     tag_i = tagRow(i)
541     rxi = qRow(1,i)
542     ryi = qRow(2,i)
543     rzi = qRow(3,i)
544     #else
545     Atype_i => identPtrList(i)%this
546     j_start = i + 1
547     rxi = q(1,i)
548     ryi = q(2,i)
549     rzi = q(3,i)
550     #endif
551    
552     inner: do j = j_start, ncol
553     #ifdef IS_MPI
554     ! Assign identity pointers and tags
555     Atype_j => identPtrListColumn(j)%this
556     tag_j = tagColumn(j)
557     if (newtons_thrd) then
558     if (tag_i <= tag_j) then
559     if (mod(tag_i + tag_j,2) == 0) cycle inner
560     else
561     if (mod(tag_i + tag_j,2) == 1) cycle inner
562     endif
563     endif
564    
565     rxij = wrap(rxi - qCol(1,j), 1)
566     ryij = wrap(ryi - qCol(2,j), 2)
567     rzij = wrap(rzi - qCol(3,j), 3)
568     #else
569     Atype_j => identPtrList(j)%this
570     rxij = wrap(rxi - q(1,j), 1)
571     ryij = wrap(ryi - q(2,j), 2)
572     rzij = wrap(rzi - q(3,j), 3)
573    
574     #endif
575     rijsq = rxij*rxij + ryij*ryij + rzij*rzij
576    
577     #ifdef IS_MPI
578     if (rijsq <= rlistsq .AND. &
579     tag_j /= tag_i) then
580     #else
581    
582     if (rijsq < rlistsq) then
583     #endif
584    
585     nlist = nlist + 1
586     if (nlist > size(list)) then
587     !! "Change how nlist size is done"
588     write(DEFAULT_ERROR,*) "ERROR: nlist > list size"
589     endif
590     list(nlist) = j
591    
592    
593     if (rijsq < rcutsq) then
594    
595     r = dsqrt(rijsq)
596    
597     call getLJPot(r,pot,dudr,Atype_i,type_j)
598    
599     #ifdef IS_MPI
600     eRow(i) = eRow(i) + pot*0.5
601     eCol(i) = eCol(i) + pot*0.5
602     #else
603     pe = pe + pot
604     #endif
605    
606     efr(1,j) = -rxij
607     efr(2,j) = -ryij
608     efr(3,j) = -rzij
609    
610     do dim = 1, 3
611    
612    
613     drdx1 = efr(dim,j) / r
614     ftmp = dudr * drdx1
615    
616    
617     #ifdef IS_MPI
618     fCol(dim,j) = fCol(dim,j) - ftmp
619     fRow(dim,i) = fRow(dim,i) + ftmp
620     #else
621    
622     f(dim,j) = f(dim,j) - ftmp
623     f(dim,i) = f(dim,i) + ftmp
624    
625     #endif
626     enddo
627     endif
628     endif
629     enddo inner
630     enddo
631    
632     #ifdef IS_MPI
633     point(nrow + 1) = nlist + 1
634     #else
635     point(natoms) = nlist + 1
636     #endif
637    
638     else
639    
640     ! use the list to find the neighbors
641     do i = 1, nrow
642     JBEG = POINT(i)
643     JEND = POINT(i+1) - 1
644     ! check thiat molecule i has neighbors
645     if (jbeg .le. jend) then
646     #ifdef IS_MPI
647     Atype_i => identPtrListRow(i)%this
648     rxi = qRow(1,i)
649     ryi = qRow(2,i)
650     rzi = qRow(3,i)
651     #else
652     Atype_i => identPtrList(i)%this
653     rxi = q(1,i)
654     ryi = q(2,i)
655     rzi = q(3,i)
656     #endif
657     do jnab = jbeg, jend
658     j = list(jnab)
659     #ifdef IS_MPI
660     Atype_j = identPtrListColumn(j)%this
661     rxij = wrap(rxi - qCol(1,j), 1)
662     ryij = wrap(ryi - qCol(2,j), 2)
663     rzij = wrap(rzi - qCol(3,j), 3)
664     #else
665     Atype_j = identPtrList(j)%this
666     rxij = wrap(rxi - q(1,j), 1)
667     ryij = wrap(ryi - q(2,j), 2)
668     rzij = wrap(rzi - q(3,j), 3)
669     #endif
670     rijsq = rxij*rxij + ryij*ryij + rzij*rzij
671    
672     if (rijsq < rcutsq) then
673    
674     r = dsqrt(rijsq)
675    
676     call getSSDPot(r,pot,dudr,Atype_i,Atype_j)
677     #ifdef IS_MPI
678     eRow(i) = eRow(i) + pot*0.5
679     eCol(i) = eCol(i) + pot*0.5
680     #else
681     pe = pe + pot
682     #endif
683    
684    
685     efr(1,j) = -rxij
686     efr(2,j) = -ryij
687     efr(3,j) = -rzij
688    
689     do dim = 1, 3
690    
691     drdx1 = efr(dim,j) / r
692     ftmp = dudr * drdx1
693     #ifdef IS_MPI
694     fCol(dim,j) = fCol(dim,j) - ftmp
695     fRow(dim,i) = fRow(dim,i) + ftmp
696     #else
697     f(dim,j) = f(dim,j) - ftmp
698     f(dim,i) = f(dim,i) + ftmp
699     #endif
700     enddo
701     endif
702     enddo
703     endif
704     enddo
705     endif
706    
707    
708    
709     #ifdef IS_MPI
710     !!distribute forces
711    
712     call scatter(fRow,f,plan_row3d)
713    
714     call scatter(fCol,fMPITemp,plan_col3d)
715    
716     do i = 1,nlocal
717     f(1:3,i) = f(1:3,i) + fMPITemp(1:3,i)
718     end do
719    
720    
721    
722     if (do_pot) then
723     ! scatter/gather pot_row into the members of my column
724     call scatter(eRow,eTemp,plan_row)
725    
726     ! scatter/gather pot_local into all other procs
727     ! add resultant to get total pot
728     do i = 1, nlocal
729     pe_local = pe_local + eTemp(i)
730     enddo
731     if (newtons_thrd) then
732     eTemp = 0.0E0_DP
733     call scatter(eCol,eTemp,plan_col)
734     do i = 1, nlocal
735     pe_local = pe_local + eTemp(i)
736     enddo
737     endif
738     pe = pe_local
739     endif
740    
741     #endif
742    
743     potE = pe
744    
745     end subroutine do_ssd_ff
746    
747    
748    
749     !! This routine does only the sticky portion of the SSD potential
750     !! [Chandra and Ichiye, J. Chem. Phys. 111, 2701 (1999)].
751     !! The Lennard-Jones and dipolar interaction must be handled separately.
752    
753     !! We assume that the rotation matrices have already been calculated
754     !! and placed in the A(9, max_mol) array.
755    
756     !! i and j are pointers to the two SSD molecules, while atom1 and
757     !! atom2 are the pointers to the location of the atoms in the force
758     !! and position arrays. These are both necessary since the rotation
759     !! matrix is a property of the molecule, while the force is acting on
760     !! the atom. The indexing of atoms and molecules is not necessarily
761     !! the same in simulations of mixtures.
762    
763     ! subroutine getSSDSticky(natoms, i, j, atom1, atom2, dx, dy, dz, rij, pot, r2, &
764     ! flx, fly, flz, tlx, tly, tlz, a)
765     subroutine getSSDSticky(r,r_ij,r__2,a,fl_ij,tl_ij,pot)
766    
767     !! Position vector between particle i and particle j
768     real( kind = dp ) :: r
769     !! Components of position vector between particle i and particle j
770     real( kind = dp), dimension(3) :: r_ij
771     !! Position vector squared
772     real( kind = dp ) :: r__2
773     !! Rotation matrix
774     real( kind = dp ), dimension(:,:) :: a
775    
776     !! force
777    
778    
779     integer :: i, j, atom1, atom2, natoms
780     double precision dx, dy, dz, rij, pot, r2
781     double precision, dimension(natoms) ::flx, fly, flz, tlx, tly, tlz
782     double precision, dimension(9,natoms):: a
783    
784    
785     double precision xi, yi, zi, xj, yj, zj, xi2, yi2, zi2, xj2, yj2, zj2
786     double precision r3, r5, r6, s, sp, dsdr, dspdr
787     double precision wi, wj, w, wip, wjp, wp
788     double precision dwidx, dwidy, dwidz, dwjdx, dwjdy, dwjdz
789     double precision dwipdx, dwipdy, dwipdz, dwjpdx, dwjpdy, dwjpdz
790     double precision dwidux, dwiduy, dwiduz, dwjdux, dwjduy, dwjduz
791     double precision dwipdux, dwipduy, dwipduz, dwjpdux, dwjpduy, dwjpduz
792     double precision zif, zis, zjf, zjs, uglyi, uglyj
793     double precision drdx, drdy, drdz
794     double precision txi, tyi, tzi, txj, tyj, tzj
795     double precision fxii, fyii, fzii, fxjj, fyjj, fzjj
796     double precision fxij, fyij, fzij, fxji, fyji, fzji
797    
798    
799     ! Use molecular positions, since the SSD model has only one atom, and the
800     ! rotation matrix is for the molecule itself:
801    
802     r3 = r2*rij
803     r5 = r3*r2
804    
805     drdx = dx / rij
806     drdy = dy / rij
807     drdz = dz / rij
808    
809     ! rotate the inter-particle separation into the two different
810     ! body-fixed coordinate systems:
811    
812     xi = a(1,i)*dx + a(2,i)*dy + a(3,i)*dz
813     yi = a(4,i)*dx + a(5,i)*dy + a(6,i)*dz
814     zi = a(7,i)*dx + a(8,i)*dy + a(9,i)*dz
815    
816     ! negative sign because this is the vector from j to i:
817    
818     xj = -(a(1,j)*dx + a(2,j)*dy + a(3,j)*dz)
819     yj = -(a(4,j)*dx + a(5,j)*dy + a(6,j)*dz)
820     zj = -(a(7,j)*dx + a(8,j)*dy + a(9,j)*dz)
821    
822     xi2 = xi*xi
823     yi2 = yi*yi
824     zi2 = zi*zi
825    
826     xj2 = xj*xj
827     yj2 = yj*yj
828     zj2 = zj*zj
829    
830     call calc_sw_fnc(rij, s, sp, dsdr, dspdr)
831    
832     wi = 2.0d0*(xi2-yi2)*zi / r3
833     wj = 2.0d0*(xj2-yj2)*zj / r3
834     w = wi+wj
835    
836     zif = zi/rij - 0.6d0
837     zis = zi/rij + 0.8d0
838    
839     zjf = zj/rij - 0.6d0
840     zjs = zj/rij + 0.8d0
841    
842     wip = zif*zif*zis*zis - SSD_w0
843     wjp = zjf*zjf*zjs*zjs - SSD_w0
844     wp = wip + wjp
845    
846     pot = pot + 0.5d0*SSD_v0*(s*w + sp*wp)
847    
848     dwidx = 4.0d0*xi*zi/r3 - 6.0d0*xi*zi*(xi2-yi2)/r5
849     dwidy = - 4.0d0*yi*zi/r3 - 6.0d0*yi*zi*(xi2-yi2)/r5
850     dwidz = 2.0d0*(xi2-yi2)/r3 - 6.0d0*zi2*(xi2-yi2)/r5
851    
852     dwjdx = 4.0d0*xj*zj/r3 - 6.0d0*xj*zj*(xj2-yj2)/r5
853     dwjdy = - 4.0d0*yj*zj/r3 - 6.0d0*yj*zj*(xj2-yj2)/r5
854     dwjdz = 2.0d0*(xj2-yj2)/r3 - 6.0d0*zj2*(xj2-yj2)/r5
855    
856     uglyi = zif*zif*zis + zif*zis*zis
857     uglyj = zjf*zjf*zjs + zjf*zjs*zjs
858    
859     dwipdx = -2.0d0*xi*zi*uglyi/r3
860     dwipdy = -2.0d0*yi*zi*uglyi/r3
861     dwipdz = 2.0d0*(1.0d0/rij - zi2/r3)*uglyi
862    
863     dwjpdx = -2.0d0*xj*zj*uglyj/r3
864     dwjpdy = -2.0d0*yj*zj*uglyj/r3
865     dwjpdz = 2.0d0*(1.0d0/rij - zj2/r3)*uglyj
866    
867     dwidux = 4.0d0*(yi*zi2 + 0.5d0*yi*(xi2-yi2))/r3
868     dwiduy = 4.0d0*(xi*zi2 - 0.5d0*xi*(xi2-yi2))/r3
869     dwiduz = - 8.0d0*xi*yi*zi/r3
870    
871     dwjdux = 4.0d0*(yj*zj2 + 0.5d0*yj*(xj2-yj2))/r3
872     dwjduy = 4.0d0*(xj*zj2 - 0.5d0*xj*(xj2-yj2))/r3
873     dwjduz = - 8.0d0*xj*yj*zj/r3
874    
875     dwipdux = 2.0d0*yi*uglyi/rij
876     dwipduy = -2.0d0*xi*uglyi/rij
877     dwipduz = 0.0d0
878    
879     dwjpdux = 2.0d0*yj*uglyj/rij
880     dwjpduy = -2.0d0*xj*uglyj/rij
881     dwjpduz = 0.0d0
882    
883     ! do the torques first since they are easy:
884     ! remember that these are still in the body fixed axes
885    
886     txi = 0.5d0*SSD_v0*(s*dwidux + sp*dwipdux)
887     tyi = 0.5d0*SSD_v0*(s*dwiduy + sp*dwipduy)
888     tzi = 0.5d0*SSD_v0*(s*dwiduz + sp*dwipduz)
889    
890     txj = 0.5d0*SSD_v0*(s*dwjdux + sp*dwjpdux)
891     tyj = 0.5d0*SSD_v0*(s*dwjduy + sp*dwjpduy)
892     tzj = 0.5d0*SSD_v0*(s*dwjduz + sp*dwjpduz)
893    
894     ! go back to lab frame using transpose of rotation matrix:
895    
896     tlx(atom1) = tlx(atom1) + a(1,i)*txi + a(4,i)*tyi + a(7,i)*tzi
897     tly(atom1) = tly(atom1) + a(2,i)*txi + a(5,i)*tyi + a(8,i)*tzi
898     tlz(atom1) = tlz(atom1) + a(3,i)*txi + a(6,i)*tyi + a(9,i)*tzi
899    
900     tlx(atom2) = tlx(atom2) + a(1,j)*txj + a(4,j)*tyj + a(7,j)*tzj
901     tly(atom2) = tly(atom2) + a(2,j)*txj + a(5,j)*tyj + a(8,j)*tzj
902     tlz(atom2) = tlz(atom2) + a(3,j)*txj + a(6,j)*tyj + a(9,j)*tzj
903    
904     ! Now, on to the forces:
905    
906     ! first rotate the i terms back into the lab frame:
907    
908     fxii = a(1,i)*(s*dwidx+sp*dwipdx) + &
909     a(4,i)*(s*dwidy+sp*dwipdy) + &
910     a(7,i)*(s*dwidz+sp*dwipdz)
911     fyii = a(2,i)*(s*dwidx+sp*dwipdx) + &
912     a(5,i)*(s*dwidy+sp*dwipdy) + &
913     a(8,i)*(s*dwidz+sp*dwipdz)
914     fzii = a(3,i)*(s*dwidx+sp*dwipdx) + &
915     a(6,i)*(s*dwidy+sp*dwipdy) + &
916     a(9,i)*(s*dwidz+sp*dwipdz)
917    
918     fxij = -fxii
919     fyij = -fyii
920     fzij = -fzii
921    
922     fxjj = a(1,j)*(s*dwjdx+sp*dwjpdx) + &
923     a(4,j)*(s*dwjdy+sp*dwjpdy) + &
924     a(7,j)*(s*dwjdz+sp*dwjpdz)
925     fyjj = a(2,j)*(s*dwjdx+sp*dwjpdx) + &
926     a(5,j)*(s*dwjdy+sp*dwjpdy) + &
927     a(8,j)*(s*dwjdz+sp*dwjpdz)
928     fzjj = a(3,j)*(s*dwjdx+sp*dwjpdx)+ &
929     a(6,j)*(s*dwjdy+sp*dwjpdy) + &
930     a(9,j)*(s*dwjdz+sp*dwjpdz)
931    
932     fxji = -fxjj
933     fyji = -fyjj
934     fzji = -fzjj
935    
936     ! now assemble these with the radial-only terms:
937    
938     flx(atom1) = flx(atom1) + 0.5d0*SSD_v0*(dsdr*drdx*w + dspdr*drdx*wp + &
939     fxii + fxji)
940     fly(atom1) = fly(atom1) + 0.5d0*SSD_v0*(dsdr*drdy*w + dspdr*drdy*wp + &
941     fyii + fyji)
942     flz(atom1) = flz(atom1) + 0.5d0*SSD_v0*(dsdr*drdz*w + dspdr*drdz*wp + &
943     fzii + fzji)
944    
945     flx(atom2) = flx(atom2) + 0.5d0*SSD_v0*(-dsdr*drdx*w - dspdr*drdx*wp + &
946     fxjj + fxij)
947     fly(atom2) = fly(atom2) + 0.5d0*SSD_v0*(-dsdr*drdy*w - dspdr*drdy*wp + &
948     fyjj + fyij)
949     flz(atom2) = flz(atom2) + 0.5d0*SSD_v0*(-dsdr*drdz*w - dspdr*drdz*wp + &
950     fzjj + fzij)
951    
952     end subroutine getSSDSticky
953    
954     !! calculates the switching functions and their derivatives for a given
955     subroutine calc_sw_fnc(r, s, sp, dsdr, dspdr)
956    
957    
958     ! value of r
959    
960     double precision r, s, sp, dsdr, dspdr
961     double precision rl, ru, rup
962     ! distances are in angstroms
963     parameter(rl = 2.75d0, ru = 3.35d0, rup = 4.0d0)
964    
965     if (r.lt.rl) then
966     s = 1.0d0
967     sp = 1.0d0
968     dsdr = 0.0d0
969     dspdr = 0.0d0
970     elseif (r.gt.rup) then
971     s = 0.0d0
972     sp = 0.0d0
973     dsdr = 0.0d0
974     dspdr = 0.0d0
975     else
976     sp = ((rup + 2.0d0*r - 3.0d0*rl) * (rup-r)**2)/((rup - rl)**3)
977     dspdr = 6.0d0*(r-rup)*(r-rl)/((rup - rl)**3)
978    
979     if (r.gt.ru) then
980     s = 0.0d0
981     dsdr = 0.0d0
982     else
983     s = ((ru + 2.0d0*r - 3.0d0*rl) * (ru-r)**2)/((ru - rl)**3)
984     dsdr = 6.0d0*(r-ru)*(r-rl)/((ru - rl)**3)
985     endif
986     endif
987    
988     return
989     end subroutine calc_sw_fnc
990    
991    
992    
993    
994    
995    
996    
997    
998    
999     !! Calculates the mixing for sigma or epslon based on character string for initialzition of mixing array.
1000     function calcMix(thisParam,param1,param2,status) result(myMixParam)
1001     character(len=*) :: thisParam
1002     real(kind = dp) :: param1
1003     real(kind = dp) :: param2
1004     real(kind = dp ) :: myMixParam
1005     integer, optional :: status
1006    
1007    
1008     myMixParam = 0.0_dp
1009    
1010     if (present(status)) status = 0
1011    
1012     select case (thisParam)
1013    
1014     case ("sigma")
1015     myMixParam = 0.5_dp * (param1 + param2)
1016     case ("epsilon")
1017     myMixParam = sqrt(param1 * param2)
1018     case default
1019     status = -1
1020     end select
1021    
1022     end function calcMix
1023    
1024    
1025    
1026    
1027    
1028    
1029     end module lj_ff