ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 2756
Committed: Wed May 17 15:37:15 2006 UTC (18 years, 2 months ago) by gezelter
File size: 18793 byte(s)
Log Message:
Getting fortran side prepped for single precision...

File Contents

# Content
1 !!
2 !! Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 !!
4 !! The University of Notre Dame grants you ("Licensee") a
5 !! non-exclusive, royalty free, license to use, modify and
6 !! redistribute this software in source and binary code form, provided
7 !! that the following conditions are met:
8 !!
9 !! 1. Acknowledgement of the program authors must be made in any
10 !! publication of scientific results based in part on use of the
11 !! program. An acceptable form of acknowledgement is citation of
12 !! the article in which the program was described (Matthew
13 !! A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 !! J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 !! Parallel Simulation Engine for Molecular Dynamics,"
16 !! J. Comput. Chem. 26, pp. 252-271 (2005))
17 !!
18 !! 2. Redistributions of source code must retain the above copyright
19 !! notice, this list of conditions and the following disclaimer.
20 !!
21 !! 3. Redistributions in binary form must reproduce the above copyright
22 !! notice, this list of conditions and the following disclaimer in the
23 !! documentation and/or other materials provided with the
24 !! distribution.
25 !!
26 !! This software is provided "AS IS," without a warranty of any
27 !! kind. All express or implied conditions, representations and
28 !! warranties, including any implied warranty of merchantability,
29 !! fitness for a particular purpose or non-infringement, are hereby
30 !! excluded. The University of Notre Dame and its licensors shall not
31 !! be liable for any damages suffered by licensee as a result of
32 !! using, modifying or distributing the software or its
33 !! derivatives. In no event will the University of Notre Dame or its
34 !! licensors be liable for any lost revenue, profit or data, or for
35 !! direct, indirect, special, consequential, incidental or punitive
36 !! damages, however caused and regardless of the theory of liability,
37 !! arising out of the use of or inability to use software, even if the
38 !! University of Notre Dame has been advised of the possibility of
39 !! such damages.
40 !!
41
42 !! Impliments Sutton-Chen Metallic Potential
43 !! See A.P.SUTTON and J.CHEN,PHIL MAG LETT 61,139-146,1990
44
45
46 module suttonchen
47 use definitions
48 use simulation
49 use force_globals
50 use status
51 use atype_module
52 use vector_class
53 use fForceOptions
54 use interpolation
55 #ifdef IS_MPI
56 use mpiSimulation
57 #endif
58 implicit none
59 PRIVATE
60 #define __FORTRAN90
61 #include "UseTheForce/DarkSide/fInteractionMap.h"
62
63 !! number of points for the spline approximations
64 INTEGER, PARAMETER :: np = 3000
65
66 logical, save :: SC_FF_initialized = .false.
67 integer, save :: SC_Mixing_Policy
68 real(kind = dp), save :: SC_rcut
69 logical, save :: haveRcut = .false.
70 logical, save :: haveMixingMap = .false.
71 logical, save :: useGeometricDistanceMixing = .false.
72 logical, save :: cleanArrays = .true.
73 logical, save :: arraysAllocated = .false.
74
75
76 character(len = statusMsgSize) :: errMesg
77 integer :: sc_err
78
79 character(len = 200) :: errMsg
80 character(len=*), parameter :: RoutineName = "Sutton-Chen MODULE"
81
82 type, private :: SCtype
83 integer :: atid
84 real(kind=dp) :: c
85 real(kind=dp) :: m
86 real(kind=dp) :: n
87 real(kind=dp) :: alpha
88 real(kind=dp) :: epsilon
89 real(kind=dp) :: sc_rcut
90 end type SCtype
91
92
93 !! Arrays for derivatives used in force calculation
94 real( kind = dp), dimension(:), allocatable :: rho
95 real( kind = dp), dimension(:), allocatable :: frho
96 real( kind = dp), dimension(:), allocatable :: dfrhodrho
97
98 !! Arrays for MPI storage
99 #ifdef IS_MPI
100 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col
101 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row
102 real( kind = dp),save, dimension(:), allocatable :: frho_row
103 real( kind = dp),save, dimension(:), allocatable :: frho_col
104 real( kind = dp),save, dimension(:), allocatable :: rho_row
105 real( kind = dp),save, dimension(:), allocatable :: rho_col
106 real( kind = dp),save, dimension(:), allocatable :: rho_tmp
107 #endif
108
109 type, private :: SCTypeList
110 integer :: nSCTypes = 0
111 integer :: currentSCtype = 0
112 type (SCtype), pointer :: SCtypes(:) => null()
113 integer, pointer :: atidToSCtype(:) => null()
114 end type SCTypeList
115
116 type (SCTypeList), save :: SCList
117
118 type:: MixParameters
119 real(kind=DP) :: alpha
120 real(kind=DP) :: epsilon
121 real(kind=DP) :: m
122 real(Kind=DP) :: n
123 real(kind=dp) :: rCut
124 real(kind=dp) :: vCut
125 logical :: rCutWasSet = .false.
126 type(cubicSpline) :: V
127 type(cubicSpline) :: phi
128 end type MixParameters
129
130 type(MixParameters), dimension(:,:), allocatable :: MixingMap
131
132 public :: setCutoffSC
133 public :: do_SC_pair
134 public :: newSCtype
135 public :: calc_SC_prepair_rho
136 public :: calc_SC_preforce_Frho
137 public :: clean_SC
138 public :: destroySCtypes
139 public :: getSCCut
140 ! public :: setSCDefaultCutoff
141 ! public :: setSCUniformCutoff
142
143
144 contains
145
146
147 subroutine newSCtype(c_ident,c,m,n,alpha,epsilon,status)
148 real (kind = dp ) :: c ! Density Scaling
149 real (kind = dp ) :: m ! Density Exponent
150 real (kind = dp ) :: n ! Pair Potential Exponent
151 real (kind = dp ) :: alpha ! Length Scaling
152 real (kind = dp ) :: epsilon ! Energy Scaling
153 integer :: c_ident
154 integer :: status
155 integer :: nAtypes,nSCTypes,myATID
156 integer :: maxVals
157 integer :: alloc_stat
158 integer :: current
159 integer,pointer :: Matchlist(:) => null()
160
161 status = 0
162
163
164 !! Assume that atypes has already been set and get the total number of types in atypes
165
166
167 ! check to see if this is the first time into
168 if (.not.associated(SCList%SCTypes)) then
169 call getMatchingElementList(atypes, "is_SC", .true., nSCtypes, MatchList)
170 SCList%nSCtypes = nSCtypes
171 allocate(SCList%SCTypes(nSCTypes))
172 nAtypes = getSize(atypes)
173 allocate(SCList%atidToSCType(nAtypes))
174 end if
175
176 SCList%currentSCType = SCList%currentSCType + 1
177 current = SCList%currentSCType
178
179 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
180 SCList%atidToSCType(myATID) = current
181
182 SCList%SCTypes(current)%atid = c_ident
183 SCList%SCTypes(current)%alpha = alpha
184 SCList%SCTypes(current)%c = c
185 SCList%SCTypes(current)%m = m
186 SCList%SCTypes(current)%n = n
187 SCList%SCTypes(current)%epsilon = epsilon
188 end subroutine newSCtype
189
190
191 subroutine destroySCTypes()
192 if (associated(SCList%SCtypes)) then
193 deallocate(SCList%SCTypes)
194 SCList%SCTypes=>null()
195 end if
196 if (associated(SCList%atidToSCtype)) then
197 deallocate(SCList%atidToSCtype)
198 SCList%atidToSCtype=>null()
199 end if
200 ! Reset Capacity
201 SCList%nSCTypes = 0
202 SCList%currentSCtype=0
203
204 end subroutine destroySCTypes
205
206 function getSCCut(atomID) result(cutValue)
207 integer, intent(in) :: atomID
208 integer :: scID
209 real(kind=dp) :: cutValue
210
211 scID = SCList%atidToSCType(atomID)
212 cutValue = 2.0_dp * SCList%SCTypes(scID)%alpha
213 end function getSCCut
214
215 subroutine createMixingMap()
216 integer :: nSCtypes, i, j, k
217 real ( kind = dp ) :: e1, e2, m1, m2, alpha1, alpha2, n1, n2
218 real ( kind = dp ) :: epsilon, m, n, alpha, rCut, vCut, dr, r
219 real ( kind = dp ), dimension(np) :: rvals, vvals, phivals
220
221 if (SCList%currentSCtype == 0) then
222 call handleError("SuttonChen", "No members in SCMap")
223 return
224 end if
225
226 nSCtypes = SCList%nSCtypes
227
228 if (.not. allocated(MixingMap)) then
229 allocate(MixingMap(nSCtypes, nSCtypes))
230 endif
231 useGeometricDistanceMixing = usesGeometricDistanceMixing()
232 do i = 1, nSCtypes
233
234 e1 = SCList%SCtypes(i)%epsilon
235 m1 = SCList%SCtypes(i)%m
236 n1 = SCList%SCtypes(i)%n
237 alpha1 = SCList%SCtypes(i)%alpha
238
239 do j = 1, nSCtypes
240
241 e2 = SCList%SCtypes(j)%epsilon
242 m2 = SCList%SCtypes(j)%m
243 n2 = SCList%SCtypes(j)%n
244 alpha2 = SCList%SCtypes(j)%alpha
245
246 if (useGeometricDistanceMixing) then
247 alpha = sqrt(alpha1 * alpha2) !SC formulation
248 else
249 alpha = 0.5_dp * (alpha1 + alpha2) ! Goddard formulation
250 endif
251 rcut = 2.0_dp * alpha
252 epsilon = sqrt(e1 * e2)
253 m = 0.5_dp*(m1+m2)
254 n = 0.5_dp*(n1+n2)
255
256 dr = (rCut) / dble(np-1)
257 rvals(1) = 0.0_dp
258 vvals(1) = 0.0_dp
259 phivals(1) = 0.0_dp
260
261 do k = 2, np
262 r = dble(k-1)*dr
263 rvals(k) = r
264 vvals(k) = epsilon*((alpha/r)**n)
265 phivals(k) = (alpha/r)**m
266 enddo
267
268 vCut = epsilon*((alpha/rCut)**n)
269
270 call newSpline(MixingMap(i,j)%V, rvals, vvals, .true.)
271 call newSpline(MixingMap(i,j)%phi, rvals, phivals, .true.)
272
273 MixingMap(i,j)%epsilon = epsilon
274 MixingMap(i,j)%m = m
275 MixingMap(i,j)%n = n
276 MixingMap(i,j)%alpha = alpha
277 MixingMap(i,j)%rCut = rcut
278 MixingMap(i,j)%vCut = vCut
279 enddo
280 enddo
281
282 haveMixingMap = .true.
283
284 end subroutine createMixingMap
285
286
287 !! routine checks to see if array is allocated, deallocates array if allocated
288 !! and then creates the array to the required size
289 subroutine allocateSC()
290 integer :: status
291
292 #ifdef IS_MPI
293 integer :: nAtomsInRow
294 integer :: nAtomsInCol
295 #endif
296 integer :: alloc_stat
297
298
299 status = 0
300 #ifdef IS_MPI
301 nAtomsInRow = getNatomsInRow(plan_atom_row)
302 nAtomsInCol = getNatomsInCol(plan_atom_col)
303 #endif
304
305 if (allocated(frho)) deallocate(frho)
306 allocate(frho(nlocal),stat=alloc_stat)
307 if (alloc_stat /= 0) then
308 status = -1
309 end if
310
311 if (allocated(rho)) deallocate(rho)
312 allocate(rho(nlocal),stat=alloc_stat)
313 if (alloc_stat /= 0) then
314 status = -1
315 end if
316
317 if (allocated(dfrhodrho)) deallocate(dfrhodrho)
318 allocate(dfrhodrho(nlocal),stat=alloc_stat)
319 if (alloc_stat /= 0) then
320 status = -1
321 end if
322
323 #ifdef IS_MPI
324
325 if (allocated(rho_tmp)) deallocate(rho_tmp)
326 allocate(rho_tmp(nlocal),stat=alloc_stat)
327 if (alloc_stat /= 0) then
328 status = -1
329 end if
330
331
332 if (allocated(frho_row)) deallocate(frho_row)
333 allocate(frho_row(nAtomsInRow),stat=alloc_stat)
334 if (alloc_stat /= 0) then
335 status = -1
336 end if
337 if (allocated(rho_row)) deallocate(rho_row)
338 allocate(rho_row(nAtomsInRow),stat=alloc_stat)
339 if (alloc_stat /= 0) then
340 status = -1
341 end if
342 if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
343 allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
344 if (alloc_stat /= 0) then
345 status = -1
346 end if
347
348
349 ! Now do column arrays
350
351 if (allocated(frho_col)) deallocate(frho_col)
352 allocate(frho_col(nAtomsInCol),stat=alloc_stat)
353 if (alloc_stat /= 0) then
354 status = -1
355 end if
356 if (allocated(rho_col)) deallocate(rho_col)
357 allocate(rho_col(nAtomsInCol),stat=alloc_stat)
358 if (alloc_stat /= 0) then
359 status = -1
360 end if
361 if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
362 allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
363 if (alloc_stat /= 0) then
364 status = -1
365 end if
366
367 #endif
368 if (status == -1) then
369 call handleError("SuttonChen:allocateSC","Error in allocating SC arrays")
370 end if
371 arraysAllocated = .true.
372 end subroutine allocateSC
373
374 subroutine setCutoffSC(rcut)
375 real(kind=dp) :: rcut
376 SC_rcut = rcut
377 end subroutine setCutoffSC
378
379 !! This array allocates module arrays if needed and builds mixing map.
380 subroutine clean_SC()
381 if (.not.arraysAllocated) call allocateSC()
382 if (.not.haveMixingMap) call createMixingMap()
383 ! clean non-IS_MPI first
384 frho = 0.0_dp
385 rho = 0.0_dp
386 dfrhodrho = 0.0_dp
387 ! clean MPI if needed
388 #ifdef IS_MPI
389 frho_row = 0.0_dp
390 frho_col = 0.0_dp
391 rho_row = 0.0_dp
392 rho_col = 0.0_dp
393 rho_tmp = 0.0_dp
394 dfrhodrho_row = 0.0_dp
395 dfrhodrho_col = 0.0_dp
396 #endif
397 end subroutine clean_SC
398
399 !! Calculates rho_r
400 subroutine calc_sc_prepair_rho(atom1, atom2, d, r, rijsq, rcut)
401 integer :: atom1,atom2
402 real(kind = dp), dimension(3) :: d
403 real(kind = dp), intent(inout) :: r
404 real(kind = dp), intent(inout) :: rijsq, rcut
405 ! value of electron density rho do to atom i at atom j
406 real(kind = dp) :: rho_i_at_j
407 ! value of electron density rho do to atom j at atom i
408 real(kind = dp) :: rho_j_at_i
409
410 ! we don't use the derivatives, dummy variables
411 real( kind = dp) :: drho
412 integer :: sc_err
413
414 integer :: atid1,atid2 ! Global atid
415 integer :: myid_atom1 ! EAM atid
416 integer :: myid_atom2
417
418
419 ! check to see if we need to be cleaned at the start of a force loop
420
421 if (cleanArrays) call clean_SC()
422 cleanArrays = .false.
423
424 #ifdef IS_MPI
425 Atid1 = Atid_row(Atom1)
426 Atid2 = Atid_col(Atom2)
427 #else
428 Atid1 = Atid(Atom1)
429 Atid2 = Atid(Atom2)
430 #endif
431
432 Myid_atom1 = SCList%atidtoSCtype(Atid1)
433 Myid_atom2 = SCList%atidtoSCtype(Atid2)
434
435 call lookupUniformSpline(MixingMap(myid_atom1,myid_atom2)%phi, r, &
436 rho_i_at_j)
437 rho_j_at_i = rho_i_at_j
438
439 #ifdef IS_MPI
440 rho_col(atom2) = rho_col(atom2) + rho_i_at_j
441 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
442 #else
443 rho(atom2) = rho(atom2) + rho_i_at_j
444 rho(atom1) = rho(atom1) + rho_j_at_i
445 #endif
446
447 end subroutine calc_sc_prepair_rho
448
449
450 !! Calculate the rho_a for all local atoms
451 subroutine calc_sc_preforce_Frho(nlocal,pot)
452 integer :: nlocal
453 real(kind=dp) :: pot
454 integer :: i,j
455 integer :: atom
456 integer :: atype1
457 integer :: atid1
458 integer :: myid
459
460 !! Scatter the electron density from pre-pair calculation back to
461 !! local atoms
462 #ifdef IS_MPI
463 call scatter(rho_row,rho,plan_atom_row,sc_err)
464 if (sc_err /= 0 ) then
465 write(errMsg,*) " Error scattering rho_row into rho"
466 call handleError(RoutineName,errMesg)
467 endif
468 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
469 if (sc_err /= 0 ) then
470 write(errMsg,*) " Error scattering rho_col into rho"
471 call handleError(RoutineName,errMesg)
472
473 endif
474
475 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
476 #endif
477
478 !! Calculate F(rho) and derivative
479 do atom = 1, nlocal
480 Myid = SCList%atidtoSctype(Atid(atom))
481 frho(atom) = - SCList%SCTypes(Myid)%c * &
482 SCList%SCTypes(Myid)%epsilon * sqrt(rho(atom))
483
484 dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom)
485
486 pot = pot + frho(atom)
487 enddo
488
489 #ifdef IS_MPI
490 !! communicate f(rho) and derivatives back into row and column arrays
491 call gather(frho,frho_row,plan_atom_row, sc_err)
492 if (sc_err /= 0) then
493 call handleError("cal_eam_forces()","MPI gather frho_row failure")
494 endif
495 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
496 if (sc_err /= 0) then
497 call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
498 endif
499 call gather(frho,frho_col,plan_atom_col, sc_err)
500 if (sc_err /= 0) then
501 call handleError("cal_eam_forces()","MPI gather frho_col failure")
502 endif
503 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
504 if (sc_err /= 0) then
505 call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
506 endif
507 #endif
508
509
510 end subroutine calc_sc_preforce_Frho
511
512 !! Does Sutton-Chen pairwise Force calculation.
513 subroutine do_sc_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
514 pot, f, do_pot)
515 !Arguments
516 integer, intent(in) :: atom1, atom2
517 real( kind = dp ), intent(in) :: rij, r2, rcut
518 real( kind = dp ) :: pot, sw, vpair
519 real( kind = dp ), dimension(3,nLocal) :: f
520 real( kind = dp ), intent(in), dimension(3) :: d
521 real( kind = dp ), intent(inout), dimension(3) :: fpair
522
523 logical, intent(in) :: do_pot
524
525 real( kind = dp ) :: drdx, drdy, drdz
526 real( kind = dp ) :: dvpdr
527 real( kind = dp ) :: rhtmp, drhodr
528 real( kind = dp ) :: dudr
529 real( kind = dp ) :: Fx,Fy,Fz
530 real( kind = dp ) :: pot_temp, vptmp
531 real( kind = dp ) :: rcij, vcij
532 integer :: id1, id2
533 integer :: mytype_atom1
534 integer :: mytype_atom2
535 integer :: atid1, atid2
536 !Local Variables
537
538 cleanArrays = .true.
539
540 #ifdef IS_MPI
541 atid1 = atid_row(atom1)
542 atid2 = atid_col(atom2)
543 #else
544 atid1 = atid(atom1)
545 atid2 = atid(atom2)
546 #endif
547
548 mytype_atom1 = SCList%atidToSCType(atid1)
549 mytype_atom2 = SCList%atidTOSCType(atid2)
550
551 drdx = d(1)/rij
552 drdy = d(2)/rij
553 drdz = d(3)/rij
554
555 rcij = MixingMap(mytype_atom1,mytype_atom2)%rCut
556 vcij = MixingMap(mytype_atom1,mytype_atom2)%vCut
557
558 call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%phi, &
559 rij, rhtmp, drhodr)
560
561 call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%V, &
562 rij, vptmp, dvpdr)
563
564 #ifdef IS_MPI
565 dudr = drhodr*(dfrhodrho_row(atom1) + dfrhodrho_col(atom2)) + dvpdr
566 #else
567 dudr = drhodr*(dfrhodrho(atom1)+ dfrhodrho(atom2)) + dvpdr
568 #endif
569
570 pot_temp = vptmp - vcij
571
572 fx = dudr * drdx
573 fy = dudr * drdy
574 fz = dudr * drdz
575
576 #ifdef IS_MPI
577 if (do_pot) then
578 pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + (pot_temp)*0.5
579 pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + (pot_temp)*0.5
580 end if
581
582 f_Row(1,atom1) = f_Row(1,atom1) + fx
583 f_Row(2,atom1) = f_Row(2,atom1) + fy
584 f_Row(3,atom1) = f_Row(3,atom1) + fz
585
586 f_Col(1,atom2) = f_Col(1,atom2) - fx
587 f_Col(2,atom2) = f_Col(2,atom2) - fy
588 f_Col(3,atom2) = f_Col(3,atom2) - fz
589 #else
590
591 if(do_pot) then
592 pot = pot + pot_temp
593 end if
594
595 f(1,atom1) = f(1,atom1) + fx
596 f(2,atom1) = f(2,atom1) + fy
597 f(3,atom1) = f(3,atom1) + fz
598
599 f(1,atom2) = f(1,atom2) - fx
600 f(2,atom2) = f(2,atom2) - fy
601 f(3,atom2) = f(3,atom2) - fz
602 #endif
603
604
605 #ifdef IS_MPI
606 id1 = AtomRowToGlobal(atom1)
607 id2 = AtomColToGlobal(atom2)
608 #else
609 id1 = atom1
610 id2 = atom2
611 #endif
612
613 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
614
615 fpair(1) = fpair(1) + fx
616 fpair(2) = fpair(2) + fy
617 fpair(3) = fpair(3) + fz
618
619 endif
620 end subroutine do_sc_pair
621 end module suttonchen