ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/DarkSide/suttonchen.F90
Revision: 2680
Committed: Thu Mar 30 21:08:48 2006 UTC (18 years, 3 months ago) by chuckv
File size: 19342 byte(s)
Log Message:
Fixed memory leak bug where SCList capacity was not reset to zero on destruction.

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 #ifdef IS_MPI
54 use mpiSimulation
55 #endif
56 implicit none
57 PRIVATE
58 #define __FORTRAN90
59 #include "UseTheForce/DarkSide/fInteractionMap.h"
60
61 INTEGER, PARAMETER :: DP = selected_real_kind(15)
62
63 logical, save :: SC_FF_initialized = .false.
64 integer, save :: SC_Mixing_Policy
65 real(kind = dp), save :: SC_rcut
66 logical, save :: haveRcut = .false.
67 logical, save :: haveMixingMap = .false.
68 logical, save :: useGeometricDistanceMixing = .false.
69 logical, save :: cleanArrays = .true.
70 logical, save :: arraysAllocated = .false.
71
72
73 character(len = statusMsgSize) :: errMesg
74 integer :: sc_err
75
76 character(len = 200) :: errMsg
77 character(len=*), parameter :: RoutineName = "Sutton-Chen MODULE"
78 !! Logical that determines if eam arrays should be zeroed
79 logical :: nmflag = .false.
80
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
99
100 !! Arrays for MPI storage
101 #ifdef IS_MPI
102 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_col
103 real( kind = dp),save, dimension(:), allocatable :: dfrhodrho_row
104 real( kind = dp),save, dimension(:), allocatable :: frho_row
105 real( kind = dp),save, dimension(:), allocatable :: frho_col
106 real( kind = dp),save, dimension(:), allocatable :: rho_row
107 real( kind = dp),save, dimension(:), allocatable :: rho_col
108 real( kind = dp),save, dimension(:), allocatable :: rho_tmp
109 #endif
110
111 type, private :: SCTypeList
112 integer :: nSCTypes = 0
113 integer :: currentSCtype = 0
114
115 type (SCtype), pointer :: SCtypes(:) => null()
116 integer, pointer :: atidToSCtype(:) => null()
117 end type SCTypeList
118
119 type (SCTypeList), save :: SCList
120
121
122
123
124 type :: MixParameters
125 real(kind=DP) :: alpha
126 real(kind=DP) :: epsilon
127 real(kind=DP) :: m
128 real(Kind=DP) :: n
129 real(kind=DP) :: vpair_pot
130 real(kind=dp) :: rCut
131 logical :: rCutWasSet = .false.
132
133 end type MixParameters
134
135 type(MixParameters), dimension(:,:), allocatable :: MixingMap
136
137
138
139 public :: setCutoffSC
140 public :: do_SC_pair
141 public :: newSCtype
142 public :: calc_SC_prepair_rho
143 public :: calc_SC_preforce_Frho
144 public :: clean_SC
145 public :: destroySCtypes
146 public :: getSCCut
147 ! public :: setSCDefaultCutoff
148 ! public :: setSCUniformCutoff
149
150
151 contains
152
153
154 subroutine newSCtype(c_ident,c,m,n,alpha,epsilon,status)
155 real (kind = dp ) :: c ! Density Scaling
156 real (kind = dp ) :: m ! Density Exponent
157 real (kind = dp ) :: n ! Pair Potential Exponent
158 real (kind = dp ) :: alpha ! Length Scaling
159 real (kind = dp ) :: epsilon ! Energy Scaling
160
161
162 integer :: c_ident
163 integer :: status
164
165 integer :: nAtypes,nSCTypes,myATID
166 integer :: maxVals
167 integer :: alloc_stat
168 integer :: current
169 integer,pointer :: Matchlist(:) => null()
170
171 status = 0
172
173
174 !! Assume that atypes has already been set and get the total number of types in atypes
175
176
177 ! check to see if this is the first time into
178 if (.not.associated(SCList%SCTypes)) then
179 call getMatchingElementList(atypes, "is_SC", .true., nSCtypes, MatchList)
180 SCList%nSCtypes = nSCtypes
181 allocate(SCList%SCTypes(nSCTypes))
182 nAtypes = getSize(atypes)
183 allocate(SCList%atidToSCType(nAtypes))
184 end if
185
186 SCList%currentSCType = SCList%currentSCType + 1
187 current = SCList%currentSCType
188
189 myATID = getFirstMatchingElement(atypes, "c_ident", c_ident)
190 SCList%atidToSCType(myATID) = current
191
192
193
194 SCList%SCTypes(current)%atid = c_ident
195 SCList%SCTypes(current)%alpha = alpha
196 SCList%SCTypes(current)%c = c
197 SCList%SCTypes(current)%m = m
198 SCList%SCTypes(current)%n = n
199 SCList%SCTypes(current)%epsilon = epsilon
200 end subroutine newSCtype
201
202
203 subroutine destroySCTypes()
204 if (associated(SCList%SCtypes)) then
205 deallocate(SCList%SCTypes)
206 SCList%SCTypes=>null()
207 end if
208 if (associated(SCList%atidToSCtype)) then
209 deallocate(SCList%atidToSCtype)
210 SCList%atidToSCtype=>null()
211 end if
212 ! Reset Capacity
213 SCList% nSCTypes = 0
214 SCList%currentSCtype=0
215
216 end subroutine destroySCTypes
217
218
219
220 function getSCCut(atomID) result(cutValue)
221 integer, intent(in) :: atomID
222 integer :: scID
223 real(kind=dp) :: cutValue
224
225 scID = SCList%atidToSCType(atomID)
226 cutValue = 2.0_dp * SCList%SCTypes(scID)%alpha
227 end function getSCCut
228
229
230
231
232 subroutine createMixingMap()
233 integer :: nSCtypes, i, j
234 real ( kind = dp ) :: e1, e2,m1,m2,alpha1,alpha2,n1,n2
235 real ( kind = dp ) :: rcut6, tp6, tp12
236 logical :: isSoftCore1, isSoftCore2, doShift
237
238 if (SCList%currentSCtype == 0) then
239 call handleError("SuttonChen", "No members in SCMap")
240 return
241 end if
242
243 nSCtypes = SCList%nSCtypes
244
245 if (.not. allocated(MixingMap)) then
246 allocate(MixingMap(nSCtypes, nSCtypes))
247 endif
248 useGeometricDistanceMixing = usesGeometricDistanceMixing()
249 do i = 1, nSCtypes
250
251 e1 = SCList%SCtypes(i)%epsilon
252 m1 = SCList%SCtypes(i)%m
253 n1 = SCList%SCtypes(i)%n
254 alpha1 = SCList%SCtypes(i)%alpha
255
256 do j = i, nSCtypes
257
258
259 e2 = SCList%SCtypes(j)%epsilon
260 m2 = SCList%SCtypes(j)%m
261 n2 = SCList%SCtypes(j)%n
262 alpha2 = SCList%SCtypes(j)%alpha
263
264 if (useGeometricDistanceMixing) then
265 MixingMap(i,j)%alpha = sqrt(alpha1 * alpha2) !SC formulation
266 else
267 MixingMap(i,j)%alpha = 0.5_dp * (alpha1 + alpha2) ! Goddard formulation
268 endif
269
270 MixingMap(i,j)%epsilon = sqrt(e1 * e2)
271 MixingMap(i,j)%m = 0.5_dp*(m1+m2)
272 MixingMap(i,j)%n = 0.5_dp*(n1+n2)
273 MixingMap(i,j)%alpha = 0.5_dp*(alpha1+alpha2)
274 MixingMap(i,j)%rcut = 2.0_dp *MixingMap(i,j)%alpha
275 MixingMap(i,j)%vpair_pot = MixingMap(i,j)%epsilon* &
276 (MixingMap(i,j)%alpha/MixingMap(i,j)%rcut)**MixingMap(i,j)%n
277
278 if (i.ne.j) then
279 MixingMap(j,i)%epsilon = MixingMap(i,j)%epsilon
280 MixingMap(j,i)%m = MixingMap(i,j)%m
281 MixingMap(j,i)%n = MixingMap(i,j)%n
282 MixingMap(j,i)%alpha = MixingMap(i,j)%alpha
283 MixingMap(j,i)%rcut = MixingMap(i,j)%rcut
284 MixingMap(j,i)%vpair_pot = MixingMap(i,j)%vpair_pot
285 endif
286 enddo
287 enddo
288
289 haveMixingMap = .true.
290
291 end subroutine createMixingMap
292
293
294 !! routine checks to see if array is allocated, deallocates array if allocated
295 !! and then creates the array to the required size
296 subroutine allocateSC()
297 integer :: status
298
299 #ifdef IS_MPI
300 integer :: nAtomsInRow
301 integer :: nAtomsInCol
302 #endif
303 integer :: alloc_stat
304
305
306 status = 0
307 #ifdef IS_MPI
308 nAtomsInRow = getNatomsInRow(plan_atom_row)
309 nAtomsInCol = getNatomsInCol(plan_atom_col)
310 #endif
311
312
313
314 if (allocated(frho)) deallocate(frho)
315 allocate(frho(nlocal),stat=alloc_stat)
316 if (alloc_stat /= 0) then
317 status = -1
318 end if
319
320 if (allocated(rho)) deallocate(rho)
321 allocate(rho(nlocal),stat=alloc_stat)
322 if (alloc_stat /= 0) then
323 status = -1
324 end if
325
326 if (allocated(dfrhodrho)) deallocate(dfrhodrho)
327 allocate(dfrhodrho(nlocal),stat=alloc_stat)
328 if (alloc_stat /= 0) then
329 status = -1
330 end if
331
332 #ifdef IS_MPI
333
334 if (allocated(rho_tmp)) deallocate(rho_tmp)
335 allocate(rho_tmp(nlocal),stat=alloc_stat)
336 if (alloc_stat /= 0) then
337 status = -1
338 end if
339
340
341 if (allocated(frho_row)) deallocate(frho_row)
342 allocate(frho_row(nAtomsInRow),stat=alloc_stat)
343 if (alloc_stat /= 0) then
344 status = -1
345 end if
346 if (allocated(rho_row)) deallocate(rho_row)
347 allocate(rho_row(nAtomsInRow),stat=alloc_stat)
348 if (alloc_stat /= 0) then
349 status = -1
350 end if
351 if (allocated(dfrhodrho_row)) deallocate(dfrhodrho_row)
352 allocate(dfrhodrho_row(nAtomsInRow),stat=alloc_stat)
353 if (alloc_stat /= 0) then
354 status = -1
355 end if
356
357
358 ! Now do column arrays
359
360 if (allocated(frho_col)) deallocate(frho_col)
361 allocate(frho_col(nAtomsInCol),stat=alloc_stat)
362 if (alloc_stat /= 0) then
363 status = -1
364 end if
365 if (allocated(rho_col)) deallocate(rho_col)
366 allocate(rho_col(nAtomsInCol),stat=alloc_stat)
367 if (alloc_stat /= 0) then
368 status = -1
369 end if
370 if (allocated(dfrhodrho_col)) deallocate(dfrhodrho_col)
371 allocate(dfrhodrho_col(nAtomsInCol),stat=alloc_stat)
372 if (alloc_stat /= 0) then
373 status = -1
374 end if
375
376 #endif
377 if (status == -1) then
378 call handleError("SuttonChen:allocateSC","Error in allocating SC arrays")
379 end if
380 arraysAllocated = .true.
381 end subroutine allocateSC
382
383 !! C sets rcut to be the largest cutoff of any atype
384 !! present in this simulation. Doesn't include all atypes
385 !! sim knows about, just those in the simulation.
386 subroutine setCutoffSC(rcut, status)
387 real(kind=dp) :: rcut
388 integer :: status
389 status = 0
390
391 SC_rcut = rcut
392
393 end subroutine setCutoffSC
394
395 !! This array allocates module arrays if needed and builds mixing map.
396 subroutine clean_SC()
397 if (.not.arraysAllocated) call allocateSC()
398 if (.not.haveMixingMap) call createMixingMap()
399 ! clean non-IS_MPI first
400 frho = 0.0_dp
401 rho = 0.0_dp
402 dfrhodrho = 0.0_dp
403 ! clean MPI if needed
404 #ifdef IS_MPI
405 frho_row = 0.0_dp
406 frho_col = 0.0_dp
407 rho_row = 0.0_dp
408 rho_col = 0.0_dp
409 rho_tmp = 0.0_dp
410 dfrhodrho_row = 0.0_dp
411 dfrhodrho_col = 0.0_dp
412 #endif
413 end subroutine clean_SC
414
415
416
417 !! Calculates rho_r
418 subroutine calc_sc_prepair_rho(atom1,atom2,d,r,rijsq, rcut)
419 integer :: atom1,atom2
420 real(kind = dp), dimension(3) :: d
421 real(kind = dp), intent(inout) :: r
422 real(kind = dp), intent(inout) :: rijsq, rcut
423 ! value of electron density rho do to atom i at atom j
424 real(kind = dp) :: rho_i_at_j
425 ! value of electron density rho do to atom j at atom i
426 real(kind = dp) :: rho_j_at_i
427
428 ! we don't use the derivatives, dummy variables
429 real( kind = dp) :: drho
430 integer :: sc_err
431
432 integer :: atid1,atid2 ! Global atid
433 integer :: myid_atom1 ! EAM atid
434 integer :: myid_atom2
435
436
437 ! check to see if we need to be cleaned at the start of a force loop
438
439 if (cleanArrays) call clean_SC()
440 cleanArrays = .false.
441
442 #ifdef IS_MPI
443 Atid1 = Atid_row(Atom1)
444 Atid2 = Atid_col(Atom2)
445 #else
446 Atid1 = Atid(Atom1)
447 Atid2 = Atid(Atom2)
448 #endif
449
450 Myid_atom1 = SCList%atidtoSCtype(Atid1)
451 Myid_atom2 = SCList%atidtoSCtype(Atid2)
452
453
454
455 rho_i_at_j = (MixingMap(Myid_atom1,Myid_atom2)%alpha/r)&
456 **MixingMap(Myid_atom1,Myid_atom2)%m
457 rho_j_at_i = rho_i_at_j
458
459 #ifdef IS_MPI
460 rho_col(atom2) = rho_col(atom2) + rho_i_at_j
461 rho_row(atom1) = rho_row(atom1) + rho_j_at_i
462 #else
463 rho(atom2) = rho(atom2) + rho_i_at_j
464 rho(atom1) = rho(atom1) + rho_j_at_i
465 #endif
466
467 end subroutine calc_sc_prepair_rho
468
469
470 !! Calculate the rho_a for all local atoms
471 subroutine calc_sc_preforce_Frho(nlocal,pot)
472 integer :: nlocal
473 real(kind=dp) :: pot
474 integer :: i,j
475 integer :: atom
476 integer :: atype1
477 integer :: atid1
478 integer :: myid
479
480
481 !! Scatter the electron density from pre-pair calculation back to local atoms
482 #ifdef IS_MPI
483 call scatter(rho_row,rho,plan_atom_row,sc_err)
484 if (sc_err /= 0 ) then
485 write(errMsg,*) " Error scattering rho_row into rho"
486 call handleError(RoutineName,errMesg)
487 endif
488 call scatter(rho_col,rho_tmp,plan_atom_col,sc_err)
489 if (sc_err /= 0 ) then
490 write(errMsg,*) " Error scattering rho_col into rho"
491 call handleError(RoutineName,errMesg)
492
493 endif
494
495 rho(1:nlocal) = rho(1:nlocal) + rho_tmp(1:nlocal)
496 #endif
497
498
499
500 !! Calculate F(rho) and derivative
501 do atom = 1, nlocal
502 Myid = SCList%atidtoSctype(Atid(atom))
503 frho(atom) = -SCList%SCTypes(Myid)%c * &
504 SCList%SCTypes(Myid)%epsilon * sqrt(rho(atom))
505
506 dfrhodrho(atom) = 0.5_dp*frho(atom)/rho(atom)
507
508 pot = pot + frho(atom)
509 enddo
510
511 #ifdef IS_MPI
512 !! communicate f(rho) and derivatives back into row and column arrays
513 call gather(frho,frho_row,plan_atom_row, sc_err)
514 if (sc_err /= 0) then
515 call handleError("cal_eam_forces()","MPI gather frho_row failure")
516 endif
517 call gather(dfrhodrho,dfrhodrho_row,plan_atom_row, sc_err)
518 if (sc_err /= 0) then
519 call handleError("cal_eam_forces()","MPI gather dfrhodrho_row failure")
520 endif
521 call gather(frho,frho_col,plan_atom_col, sc_err)
522 if (sc_err /= 0) then
523 call handleError("cal_eam_forces()","MPI gather frho_col failure")
524 endif
525 call gather(dfrhodrho,dfrhodrho_col,plan_atom_col, sc_err)
526 if (sc_err /= 0) then
527 call handleError("cal_eam_forces()","MPI gather dfrhodrho_col failure")
528 endif
529 #endif
530
531
532 end subroutine calc_sc_preforce_Frho
533
534
535 !! Does Sutton-Chen pairwise Force calculation.
536 subroutine do_sc_pair(atom1, atom2, d, rij, r2, rcut, sw, vpair, fpair, &
537 pot, f, do_pot)
538 !Arguments
539 integer, intent(in) :: atom1, atom2
540 real( kind = dp ), intent(in) :: rij, r2, rcut
541 real( kind = dp ) :: pot, sw, vpair
542 real( kind = dp ), dimension(3,nLocal) :: f
543 real( kind = dp ), intent(in), dimension(3) :: d
544 real( kind = dp ), intent(inout), dimension(3) :: fpair
545
546 logical, intent(in) :: do_pot
547
548 real( kind = dp ) :: drdx,drdy,drdz
549 real( kind = dp ) :: dvpdr
550 real( kind = dp ) :: drhodr
551 real( kind = dp ) :: dudr
552 real( kind = dp ) :: drhoidr,drhojdr
553 real( kind = dp ) :: Fx,Fy,Fz
554 real( kind = dp ) :: d2pha,phb,d2phb
555 real( kind = dp ) :: pot_temp,vptmp
556 real( kind = dp ) :: epsilonij,aij,nij,mij,vcij
557 integer :: id1,id2
558 integer :: mytype_atom1
559 integer :: mytype_atom2
560 integer :: atid1,atid2
561 !Local Variables
562
563 ! write(*,*) "Frho: ", Frho(atom1)
564
565 cleanArrays = .true.
566
567 dvpdr = 0.0E0_DP
568
569
570 #ifdef IS_MPI
571 atid1 = atid_row(atom1)
572 atid2 = atid_col(atom2)
573 #else
574 atid1 = atid(atom1)
575 atid2 = atid(atom2)
576 #endif
577
578 mytype_atom1 = SCList%atidToSCType(atid1)
579 mytype_atom2 = SCList%atidTOSCType(atid2)
580
581 drdx = d(1)/rij
582 drdy = d(2)/rij
583 drdz = d(3)/rij
584
585
586 epsilonij = MixingMap(mytype_atom1,mytype_atom2)%epsilon
587 aij = MixingMap(mytype_atom1,mytype_atom2)%alpha
588 nij = MixingMap(mytype_atom1,mytype_atom2)%n
589 mij = MixingMap(mytype_atom1,mytype_atom2)%m
590 vcij = MixingMap(mytype_atom1,mytype_atom2)%vpair_pot
591
592 vptmp = epsilonij*((aij/rij)**nij)
593
594
595 dvpdr = -nij*vptmp/rij
596 drhodr = -mij*((aij/rij)**mij)/rij
597
598
599 dudr = drhodr*(dfrhodrho(atom1)+dfrhodrho(atom2)) &
600 + dvpdr
601
602 pot_temp = vptmp + vcij
603
604
605 #ifdef IS_MPI
606 dudr = drhodr*(dfrhodrho_row(atom1)+dfrhodrho_col(atom2)) &
607 + dvpdr
608
609 #else
610 dudr = drhodr*(dfrhodrho(atom1)+dfrhodrho(atom2)) &
611 + dvpdr
612 #endif
613
614
615 fx = dudr * drdx
616 fy = dudr * drdy
617 fz = dudr * drdz
618
619
620 #ifdef IS_MPI
621 if (do_pot) then
622 pot_Row(METALLIC_POT,atom1) = pot_Row(METALLIC_POT,atom1) + (pot_temp)*0.5
623 pot_Col(METALLIC_POT,atom2) = pot_Col(METALLIC_POT,atom2) + (pot_temp)*0.5
624 end if
625
626 f_Row(1,atom1) = f_Row(1,atom1) + fx
627 f_Row(2,atom1) = f_Row(2,atom1) + fy
628 f_Row(3,atom1) = f_Row(3,atom1) + fz
629
630 f_Col(1,atom2) = f_Col(1,atom2) - fx
631 f_Col(2,atom2) = f_Col(2,atom2) - fy
632 f_Col(3,atom2) = f_Col(3,atom2) - fz
633 #else
634
635 if(do_pot) then
636 pot = pot + pot_temp
637 end if
638
639 f(1,atom1) = f(1,atom1) + fx
640 f(2,atom1) = f(2,atom1) + fy
641 f(3,atom1) = f(3,atom1) + fz
642
643 f(1,atom2) = f(1,atom2) - fx
644 f(2,atom2) = f(2,atom2) - fy
645 f(3,atom2) = f(3,atom2) - fz
646 #endif
647
648
649 #ifdef IS_MPI
650 id1 = AtomRowToGlobal(atom1)
651 id2 = AtomColToGlobal(atom2)
652 #else
653 id1 = atom1
654 id2 = atom2
655 #endif
656
657 if (molMembershipList(id1) .ne. molMembershipList(id2)) then
658
659 fpair(1) = fpair(1) + fx
660 fpair(2) = fpair(2) + fy
661 fpair(3) = fpair(3) + fz
662
663 endif
664
665
666 end subroutine do_sc_pair
667
668
669
670 end module suttonchen