ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 1245
Committed: Tue May 27 16:39:06 2008 UTC (17 years, 2 months ago) by chuckv
File size: 19179 byte(s)
Log Message:
Checking in changes for Hefland moment calculations

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 SCList%atidToSCType = -1
175 end if
176
177 SCList%currentSCType = SCList%currentSCType + 1
178 current = SCList%currentSCType
179
180 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
181 SCList%atidToSCType(myATID) = current
182
183
184 SCList%SCTypes(current)%atid = c_ident
185 SCList%SCTypes(current)%alpha = alpha
186 SCList%SCTypes(current)%c = c
187 SCList%SCTypes(current)%m = m
188 SCList%SCTypes(current)%n = n
189 SCList%SCTypes(current)%epsilon = epsilon
190 end subroutine newSCtype
191
192
193 subroutine destroySCTypes()
194 if (associated(SCList%SCtypes)) then
195 deallocate(SCList%SCTypes)
196 SCList%SCTypes=>null()
197 end if
198 if (associated(SCList%atidToSCtype)) then
199 deallocate(SCList%atidToSCtype)
200 SCList%atidToSCtype=>null()
201 end if
202 ! Reset Capacity
203 SCList%nSCTypes = 0
204 SCList%currentSCtype=0
205
206 end subroutine destroySCTypes
207
208 function getSCCut(atomID) result(cutValue)
209 integer, intent(in) :: atomID
210 integer :: scID
211 real(kind=dp) :: cutValue
212
213 scID = SCList%atidToSCType(atomID)
214 cutValue = 2.0_dp * SCList%SCTypes(scID)%alpha
215 end function getSCCut
216
217 subroutine createMixingMap()
218 integer :: nSCtypes, i, j, k
219 real ( kind = dp ) :: e1, e2, m1, m2, alpha1, alpha2, n1, n2
220 real ( kind = dp ) :: epsilon, m, n, alpha, rCut, vCut, dr, r
221 real ( kind = dp ), dimension(np) :: rvals, vvals, phivals
222
223 if (SCList%currentSCtype == 0) then
224 call handleError("SuttonChen", "No members in SCMap")
225 return
226 end if
227
228 nSCtypes = SCList%nSCtypes
229
230 if (.not. allocated(MixingMap)) then
231 allocate(MixingMap(nSCtypes, nSCtypes))
232 endif
233 useGeometricDistanceMixing = usesGeometricDistanceMixing()
234 do i = 1, nSCtypes
235
236 e1 = SCList%SCtypes(i)%epsilon
237 m1 = SCList%SCtypes(i)%m
238 n1 = SCList%SCtypes(i)%n
239 alpha1 = SCList%SCtypes(i)%alpha
240
241 do j = 1, nSCtypes
242
243 e2 = SCList%SCtypes(j)%epsilon
244 m2 = SCList%SCtypes(j)%m
245 n2 = SCList%SCtypes(j)%n
246 alpha2 = SCList%SCtypes(j)%alpha
247
248 if (useGeometricDistanceMixing) then
249 alpha = sqrt(alpha1 * alpha2) !SC formulation
250 else
251 alpha = 0.5_dp * (alpha1 + alpha2) ! Goddard formulation
252 endif
253 rcut = 2.0_dp * alpha
254 epsilon = sqrt(e1 * e2)
255 m = 0.5_dp*(m1+m2)
256 n = 0.5_dp*(n1+n2)
257
258 dr = (rCut) / dble(np-1)
259 rvals(1) = 0.0_dp
260 vvals(1) = 0.0_dp
261 phivals(1) = 0.0_dp
262
263 do k = 2, np
264 r = dble(k-1)*dr
265 rvals(k) = r
266 vvals(k) = epsilon*((alpha/r)**n)
267 phivals(k) = (alpha/r)**m
268 enddo
269
270 vCut = epsilon*((alpha/rCut)**n)
271
272 call newSpline(MixingMap(i,j)%V, rvals, vvals, .true.)
273 call newSpline(MixingMap(i,j)%phi, rvals, phivals, .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
427
428 #ifdef IS_MPI
429 Atid1 = Atid_row(Atom1)
430 Atid2 = Atid_col(Atom2)
431 #else
432 Atid1 = Atid(Atom1)
433 Atid2 = Atid(Atom2)
434 #endif
435
436 Myid_atom1 = SCList%atidtoSCtype(Atid1)
437 Myid_atom2 = SCList%atidtoSCtype(Atid2)
438
439 call lookupUniformSpline(MixingMap(myid_atom1,myid_atom2)%phi, r, &
440 rho_i_at_j)
441 rho_j_at_i = rho_i_at_j
442
443 #ifdef IS_MPI
444 rho_col(atom2) = rho_col(atom2) + rho_i_at_j
445 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
446 #else
447 rho(atom2) = rho(atom2) + rho_i_at_j
448 rho(atom1) = rho(atom1) + rho_j_at_i
449 #endif
450
451 end subroutine calc_sc_prepair_rho
452
453
454 !! Calculate the rho_a for all local atoms
455 subroutine calc_sc_preforce_Frho(nlocal,pot)
456 integer :: nlocal
457 real(kind=dp) :: pot
458 integer :: i,j
459 integer :: atom
460 integer :: atype1
461 integer :: atid1
462 integer :: myid
463
464 !! Scatter the electron density from pre-pair calculation back to
465 !! local atoms
466
467
468 if (cleanArrays) call clean_SC()
469 cleanArrays = .false.
470
471
472 #ifdef IS_MPI
473 call scatter(rho_row,rho,plan_atom_row,sc_err)
474 if (sc_err /= 0 ) then
475 write(errMsg,*) " Error scattering rho_row into rho"
476 call handleError(RoutineName,errMesg)
477 endif
478 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
479 if (sc_err /= 0 ) then
480 write(errMsg,*) " Error scattering rho_col into rho"
481 call handleError(RoutineName,errMesg)
482
483 endif
484
485 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
486 #endif
487
488 !! Calculate F(rho) and derivative
489 do atom = 1, nlocal
490 Myid = SCList%atidtoSctype(Atid(atom))
491 ! Myid is set to -1 for non SC atoms.
492 ! Punt if we are a non-SC atom type.
493 if (Myid == -1) then
494 frho(atom) = 0.0_dp
495 dfrhodrho(atom) = 0.0_dp
496 else
497 frho(atom) = - SCList%SCTypes(Myid)%c * &
498 SCList%SCTypes(Myid)%epsilon * sqrt(rho(atom))
499
500 dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom)
501 end if
502
503 pot = pot + frho(atom)
504 enddo
505
506 #ifdef IS_MPI
507 !! communicate f(rho) and derivatives back into row and column arrays
508 call gather(frho,frho_row,plan_atom_row, sc_err)
509 if (sc_err /= 0) then
510 call handleError("cal_eam_forces()","MPI gather frho_row failure")
511 endif
512 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
513 if (sc_err /= 0) then
514 call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
515 endif
516 call gather(frho,frho_col,plan_atom_col, sc_err)
517 if (sc_err /= 0) then
518 call handleError("cal_eam_forces()","MPI gather frho_col failure")
519 endif
520 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
521 if (sc_err /= 0) then
522 call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
523 endif
524 #endif
525
526
527 end subroutine calc_sc_preforce_Frho
528
529 !! Does Sutton-Chen pairwise Force calculation.
530 subroutine do_sc_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
531 pot, f, do_pot)
532 !Arguments
533 integer, intent(in) :: atom1, atom2
534 real( kind = dp ), intent(in) :: rij, r2, rcut
535 real( kind = dp ) :: pot, sw, vpair
536 real( kind = dp ), dimension(3,nLocal) :: f
537 real( kind = dp ), intent(in), dimension(3) :: d
538 real( kind = dp ), intent(inout), dimension(3) :: fpair
539
540 logical, intent(in) :: do_pot
541
542 real( kind = dp ) :: drdx, drdy, drdz
543 real( kind = dp ) :: dvpdr
544 real( kind = dp ) :: rhtmp, drhodr
545 real( kind = dp ) :: dudr
546 real( kind = dp ) :: Fx,Fy,Fz
547 real( kind = dp ) :: pot_temp, vptmp
548 real( kind = dp ) :: rcij, vcij
549 integer :: id1, id2
550 integer :: mytype_atom1
551 integer :: mytype_atom2
552 integer :: atid1, atid2
553 !Local Variables
554
555 cleanArrays = .true.
556
557 #ifdef IS_MPI
558 atid1 = atid_row(atom1)
559 atid2 = atid_col(atom2)
560 #else
561 atid1 = atid(atom1)
562 atid2 = atid(atom2)
563 #endif
564
565 mytype_atom1 = SCList%atidToSCType(atid1)
566 mytype_atom2 = SCList%atidTOSCType(atid2)
567
568 drdx = d(1)/rij
569 drdy = d(2)/rij
570 drdz = d(3)/rij
571
572 rcij = MixingMap(mytype_atom1,mytype_atom2)%rCut
573 vcij = MixingMap(mytype_atom1,mytype_atom2)%vCut
574
575 call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%phi, &
576 rij, rhtmp, drhodr)
577
578 call lookupUniformSpline1d(MixingMap(mytype_atom1, mytype_atom2)%V, &
579 rij, vptmp, dvpdr)
580
581 #ifdef IS_MPI
582 dudr = drhodr*(dfrhodrho_row(atom1) + dfrhodrho_col(atom2)) + dvpdr
583 #else
584 dudr = drhodr*(dfrhodrho(atom1)+ dfrhodrho(atom2)) + dvpdr
585 #endif
586
587 pot_temp = vptmp - vcij
588
589 vpair = vpair + pot_temp
590
591 fx = dudr * drdx
592 fy = dudr * drdy
593 fz = dudr * drdz
594
595 #ifdef IS_MPI
596 if (do_pot) then
597 pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + (pot_temp)*0.5
598 pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + (pot_temp)*0.5
599 end if
600
601 f_Row(1,atom1) = f_Row(1,atom1) + fx
602 f_Row(2,atom1) = f_Row(2,atom1) + fy
603 f_Row(3,atom1) = f_Row(3,atom1) + fz
604
605 f_Col(1,atom2) = f_Col(1,atom2) - fx
606 f_Col(2,atom2) = f_Col(2,atom2) - fy
607 f_Col(3,atom2) = f_Col(3,atom2) - fz
608 #else
609
610 if(do_pot) then
611 pot = pot + pot_temp
612 end if
613
614 f(1,atom1) = f(1,atom1) + fx
615 f(2,atom1) = f(2,atom1) + fy
616 f(3,atom1) = f(3,atom1) + fz
617
618 f(1,atom2) = f(1,atom2) - fx
619 f(2,atom2) = f(2,atom2) - fy
620 f(3,atom2) = f(3,atom2) - fz
621 #endif
622
623
624 #ifdef IS_MPI
625 id1 = AtomRowToGlobal(atom1)
626 id2 = AtomColToGlobal(atom2)
627 #else
628 id1 = atom1
629 id2 = atom2
630 #endif
631
632 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
633
634 fpair(1) = fpair(1) + fx
635 fpair(2) = fpair(2) + fy
636 fpair(3) = fpair(3) + fz
637
638 endif
639 end subroutine do_sc_pair
640 end module suttonchen