ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 2717
Committed: Mon Apr 17 21:49:12 2006 UTC (18 years, 3 months ago) by gezelter
File size: 18887 byte(s)
Log Message:
Many performance improvements

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 simulation
48 use force_globals
49 use status
50 use atype_module
51 use vector_class
52 use fForceOptions
53 use interpolation
54 #ifdef IS_MPI
55 use mpiSimulation
56 #endif
57 implicit none
58 PRIVATE
59 #define __FORTRAN90
60 #include "UseTheForce/DarkSide/fInteractionMap.h"
61
62 INTEGER, PARAMETER :: DP = selected_real_kind(15)
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.0d0
258 vvals(1) = 0.0d0
259 phivals(1) = 0.0d0
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, &
271 0.0d0, 0.0d0, .true.)
272 call newSpline(MixingMap(i,j)%phi, rvals, phivals, &
273 0.0d0, 0.0d0, .true.)
274
275 MixingMap(i,j)%epsilon = epsilon
276 MixingMap(i,j)%m = m
277 MixingMap(i,j)%n = n
278 MixingMap(i,j)%alpha = alpha
279 MixingMap(i,j)%rCut = rcut
280 MixingMap(i,j)%vCut = vCut
281 enddo
282 enddo
283
284 haveMixingMap = .true.
285
286 end subroutine createMixingMap
287
288
289 !! routine checks to see if array is allocated, deallocates array if allocated
290 !! and then creates the array to the required size
291 subroutine allocateSC()
292 integer :: status
293
294 #ifdef IS_MPI
295 integer :: nAtomsInRow
296 integer :: nAtomsInCol
297 #endif
298 integer :: alloc_stat
299
300
301 status = 0
302 #ifdef IS_MPI
303 nAtomsInRow = getNatomsInRow(plan_atom_row)
304 nAtomsInCol = getNatomsInCol(plan_atom_col)
305 #endif
306
307 if (allocated(frho)) deallocate(frho)
308 allocate(frho(nlocal),stat=alloc_stat)
309 if (alloc_stat /= 0) then
310 status = -1
311 end if
312
313 if (allocated(rho)) deallocate(rho)
314 allocate(rho(nlocal),stat=alloc_stat)
315 if (alloc_stat /= 0) then
316 status = -1
317 end if
318
319 if (allocated(dfrhodrho)) deallocate(dfrhodrho)
320 allocate(dfrhodrho(nlocal),stat=alloc_stat)
321 if (alloc_stat /= 0) then
322 status = -1
323 end if
324
325 #ifdef IS_MPI
326
327 if (allocated(rho_tmp)) deallocate(rho_tmp)
328 allocate(rho_tmp(nlocal),stat=alloc_stat)
329 if (alloc_stat /= 0) then
330 status = -1
331 end if
332
333
334 if (allocated(frho_row)) deallocate(frho_row)
335 allocate(frho_row(nAtomsInRow),stat=alloc_stat)
336 if (alloc_stat /= 0) then
337 status = -1
338 end if
339 if (allocated(rho_row)) deallocate(rho_row)
340 allocate(rho_row(nAtomsInRow),stat=alloc_stat)
341 if (alloc_stat /= 0) then
342 status = -1
343 end if
344 if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
345 allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
346 if (alloc_stat /= 0) then
347 status = -1
348 end if
349
350
351 ! Now do column arrays
352
353 if (allocated(frho_col)) deallocate(frho_col)
354 allocate(frho_col(nAtomsInCol),stat=alloc_stat)
355 if (alloc_stat /= 0) then
356 status = -1
357 end if
358 if (allocated(rho_col)) deallocate(rho_col)
359 allocate(rho_col(nAtomsInCol),stat=alloc_stat)
360 if (alloc_stat /= 0) then
361 status = -1
362 end if
363 if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
364 allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
365 if (alloc_stat /= 0) then
366 status = -1
367 end if
368
369 #endif
370 if (status == -1) then
371 call handleError("SuttonChen:allocateSC","Error in allocating SC arrays")
372 end if
373 arraysAllocated = .true.
374 end subroutine allocateSC
375
376 subroutine setCutoffSC(rcut)
377 real(kind=dp) :: rcut
378 SC_rcut = rcut
379 end subroutine setCutoffSC
380
381 !! This array allocates module arrays if needed and builds mixing map.
382 subroutine clean_SC()
383 if (.not.arraysAllocated) call allocateSC()
384 if (.not.haveMixingMap) call createMixingMap()
385 ! clean non-IS_MPI first
386 frho = 0.0_dp
387 rho = 0.0_dp
388 dfrhodrho = 0.0_dp
389 ! clean MPI if needed
390 #ifdef IS_MPI
391 frho_row = 0.0_dp
392 frho_col = 0.0_dp
393 rho_row = 0.0_dp
394 rho_col = 0.0_dp
395 rho_tmp = 0.0_dp
396 dfrhodrho_row = 0.0_dp
397 dfrhodrho_col = 0.0_dp
398 #endif
399 end subroutine clean_SC
400
401 !! Calculates rho_r
402 subroutine calc_sc_prepair_rho(atom1, atom2, d, r, rijsq, rcut)
403 integer :: atom1,atom2
404 real(kind = dp), dimension(3) :: d
405 real(kind = dp), intent(inout) :: r
406 real(kind = dp), intent(inout) :: rijsq, rcut
407 ! value of electron density rho do to atom i at atom j
408 real(kind = dp) :: rho_i_at_j
409 ! value of electron density rho do to atom j at atom i
410 real(kind = dp) :: rho_j_at_i
411
412 ! we don't use the derivatives, dummy variables
413 real( kind = dp) :: drho
414 integer :: sc_err
415
416 integer :: atid1,atid2 ! Global atid
417 integer :: myid_atom1 ! EAM atid
418 integer :: myid_atom2
419
420
421 ! check to see if we need to be cleaned at the start of a force loop
422
423 if (cleanArrays) call clean_SC()
424 cleanArrays = .false.
425
426 #ifdef IS_MPI
427 Atid1 = Atid_row(Atom1)
428 Atid2 = Atid_col(Atom2)
429 #else
430 Atid1 = Atid(Atom1)
431 Atid2 = Atid(Atom2)
432 #endif
433
434 Myid_atom1 = SCList%atidtoSCtype(Atid1)
435 Myid_atom2 = SCList%atidtoSCtype(Atid2)
436
437 call lookupUniformSpline(MixingMap(myid_atom1,myid_atom2)%phi, r, &
438 rho_i_at_j)
439 rho_j_at_i = rho_i_at_j
440
441 #ifdef IS_MPI
442 rho_col(atom2) = rho_col(atom2) + rho_i_at_j
443 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
444 #else
445 rho(atom2) = rho(atom2) + rho_i_at_j
446 rho(atom1) = rho(atom1) + rho_j_at_i
447 #endif
448
449 end subroutine calc_sc_prepair_rho
450
451
452 !! Calculate the rho_a for all local atoms
453 subroutine calc_sc_preforce_Frho(nlocal,pot)
454 integer :: nlocal
455 real(kind=dp) :: pot
456 integer :: i,j
457 integer :: atom
458 integer :: atype1
459 integer :: atid1
460 integer :: myid
461
462 !! Scatter the electron density from pre-pair calculation back to
463 !! local atoms
464 #ifdef IS_MPI
465 call scatter(rho_row,rho,plan_atom_row,sc_err)
466 if (sc_err /= 0 ) then
467 write(errMsg,*) " Error scattering rho_row into rho"
468 call handleError(RoutineName,errMesg)
469 endif
470 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
471 if (sc_err /= 0 ) then
472 write(errMsg,*) " Error scattering rho_col into rho"
473 call handleError(RoutineName,errMesg)
474
475 endif
476
477 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
478 #endif
479
480 !! Calculate F(rho) and derivative
481 do atom = 1, nlocal
482 Myid = SCList%atidtoSctype(Atid(atom))
483 frho(atom) = - SCList%SCTypes(Myid)%c * &
484 SCList%SCTypes(Myid)%epsilon * sqrt(rho(atom))
485
486 dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom)
487
488 pot = pot + frho(atom)
489 enddo
490
491 #ifdef IS_MPI
492 !! communicate f(rho) and derivatives back into row and column arrays
493 call gather(frho,frho_row,plan_atom_row, sc_err)
494 if (sc_err /= 0) then
495 call handleError("cal_eam_forces()","MPI gather frho_row failure")
496 endif
497 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
498 if (sc_err /= 0) then
499 call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
500 endif
501 call gather(frho,frho_col,plan_atom_col, sc_err)
502 if (sc_err /= 0) then
503 call handleError("cal_eam_forces()","MPI gather frho_col failure")
504 endif
505 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
506 if (sc_err /= 0) then
507 call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
508 endif
509 #endif
510
511
512 end subroutine calc_sc_preforce_Frho
513
514 !! Does Sutton-Chen pairwise Force calculation.
515 subroutine do_sc_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
516 pot, f, do_pot)
517 !Arguments
518 integer, intent(in) :: atom1, atom2
519 real( kind = dp ), intent(in) :: rij, r2, rcut
520 real( kind = dp ) :: pot, sw, vpair
521 real( kind = dp ), dimension(3,nLocal) :: f
522 real( kind = dp ), intent(in), dimension(3) :: d
523 real( kind = dp ), intent(inout), dimension(3) :: fpair
524
525 logical, intent(in) :: do_pot
526
527 real( kind = dp ) :: drdx, drdy, drdz
528 real( kind = dp ) :: dvpdr
529 real( kind = dp ) :: rhtmp, drhodr
530 real( kind = dp ) :: dudr
531 real( kind = dp ) :: Fx,Fy,Fz
532 real( kind = dp ) :: pot_temp, vptmp
533 real( kind = dp ) :: rcij, vcij
534 integer :: id1, id2
535 integer :: mytype_atom1
536 integer :: mytype_atom2
537 integer :: atid1, atid2
538 !Local Variables
539
540 cleanArrays = .true.
541
542 #ifdef IS_MPI
543 atid1 = atid_row(atom1)
544 atid2 = atid_col(atom2)
545 #else
546 atid1 = atid(atom1)
547 atid2 = atid(atom2)
548 #endif
549
550 mytype_atom1 = SCList%atidToSCType(atid1)
551 mytype_atom2 = SCList%atidTOSCType(atid2)
552
553 drdx = d(1)/rij
554 drdy = d(2)/rij
555 drdz = d(3)/rij
556
557 rcij = MixingMap(mytype_atom1,mytype_atom2)%rCut
558 vcij = MixingMap(mytype_atom1,mytype_atom2)%vCut
559
560 call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%phi, &
561 rij, rhtmp, drhodr)
562
563 call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%V, &
564 rij, vptmp, dvpdr)
565
566 #ifdef IS_MPI
567 dudr = drhodr*(dfrhodrho_row(atom1) + dfrhodrho_col(atom2)) + dvpdr
568 #else
569 dudr = drhodr*(dfrhodrho(atom1)+ dfrhodrho(atom2)) + dvpdr
570 #endif
571
572 pot_temp = vptmp - vcij
573
574 fx = dudr * drdx
575 fy = dudr * drdy
576 fz = dudr * drdz
577
578 #ifdef IS_MPI
579 if (do_pot) then
580 pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + (pot_temp)*0.5
581 pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + (pot_temp)*0.5
582 end if
583
584 f_Row(1,atom1) = f_Row(1,atom1) + fx
585 f_Row(2,atom1) = f_Row(2,atom1) + fy
586 f_Row(3,atom1) = f_Row(3,atom1) + fz
587
588 f_Col(1,atom2) = f_Col(1,atom2) - fx
589 f_Col(2,atom2) = f_Col(2,atom2) - fy
590 f_Col(3,atom2) = f_Col(3,atom2) - fz
591 #else
592
593 if(do_pot) then
594 pot = pot + pot_temp
595 end if
596
597 f(1,atom1) = f(1,atom1) + fx
598 f(2,atom1) = f(2,atom1) + fy
599 f(3,atom1) = f(3,atom1) + fz
600
601 f(1,atom2) = f(1,atom2) - fx
602 f(2,atom2) = f(2,atom2) - fy
603 f(3,atom2) = f(3,atom2) - fz
604 #endif
605
606
607 #ifdef IS_MPI
608 id1 = AtomRowToGlobal(atom1)
609 id2 = AtomColToGlobal(atom2)
610 #else
611 id1 = atom1
612 id2 = atom2
613 #endif
614
615 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
616
617 fpair(1) = fpair(1) + fx
618 fpair(2) = fpair(2) + fy
619 fpair(3) = fpair(3) + fz
620
621 endif
622 end subroutine do_sc_pair
623 end module suttonchen