ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/simulation.F90
Revision: 1473
Committed: Tue Jul 20 15:43:00 2010 UTC (14 years, 9 months ago) by gezelter
File size: 28646 byte(s)
Log Message:
fixing c/fortran bugs

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. Redistributions of source code must retain the above copyright
10 !! notice, this list of conditions and the following disclaimer.
11 !!
12 !! 2. Redistributions in binary form must reproduce the above copyright
13 !! notice, this list of conditions and the following disclaimer in the
14 !! documentation and/or other materials provided with the
15 !! distribution.
16 !!
17 !! This software is provided "AS IS," without a warranty of any
18 !! kind. All express or implied conditions, representations and
19 !! warranties, including any implied warranty of merchantability,
20 !! fitness for a particular purpose or non-infringement, are hereby
21 !! excluded. The University of Notre Dame and its licensors shall not
22 !! be liable for any damages suffered by licensee as a result of
23 !! using, modifying or distributing the software or its
24 !! derivatives. In no event will the University of Notre Dame or its
25 !! licensors be liable for any lost revenue, profit or data, or for
26 !! direct, indirect, special, consequential, incidental or punitive
27 !! damages, however caused and regardless of the theory of liability,
28 !! arising out of the use of or inability to use software, even if the
29 !! University of Notre Dame has been advised of the possibility of
30 !! such damages.
31 !!
32 !! SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 !! research, please cite the appropriate papers when you publish your
34 !! work. Good starting points are:
35 !!
36 !! [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 !! [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 !! [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 !! [4] Vardeman & Gezelter, in progress (2009).
40 !!
41
42 !! Fortran interface to C entry plug.
43
44 module simulation
45 use definitions
46 use status
47 use linearAlgebra
48 use neighborLists
49 use force_globals
50 use vector_class
51 use atype_module
52 use switcheroo
53 #ifdef IS_MPI
54 use mpiSimulation
55 #endif
56
57 implicit none
58 PRIVATE
59
60 #define __FORTRAN90
61 #include "brains/fSimulation.h"
62 #include "UseTheForce/fSwitchingFunction.h"
63 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
64
65 type (simtype), public, save :: thisSim
66
67 logical, save :: simulation_setup_complete = .false.
68
69 integer, public, save :: nLocal, nGlobal
70 integer, public, save :: nGroups, nGroupGlobal
71 integer, public, save :: nExcludes = 0
72 integer, public, save :: nOneTwo = 0
73 integer, public, save :: nOneThree = 0
74 integer, public, save :: nOneFour = 0
75
76 integer, allocatable, dimension(:,:), public :: excludes
77 integer, allocatable, dimension(:), public :: molMembershipList
78 integer, allocatable, dimension(:), public :: groupListRow
79 integer, allocatable, dimension(:), public :: groupStartRow
80 integer, allocatable, dimension(:), public :: groupListCol
81 integer, allocatable, dimension(:), public :: groupStartCol
82 integer, allocatable, dimension(:), public :: groupListLocal
83 integer, allocatable, dimension(:), public :: groupStartLocal
84 integer, allocatable, dimension(:), public :: nSkipsForLocalAtom
85 integer, allocatable, dimension(:,:), public :: skipsForLocalAtom
86 integer, allocatable, dimension(:), public :: nSkipsForRowAtom
87 integer, allocatable, dimension(:,:), public :: skipsForRowAtom
88 integer, allocatable, dimension(:), public :: nTopoPairsForAtom
89 integer, allocatable, dimension(:,:), public :: toposForAtom
90 integer, allocatable, dimension(:,:), public :: topoDistance
91 real(kind=dp), allocatable, dimension(:), public :: mfactRow
92 real(kind=dp), allocatable, dimension(:), public :: mfactCol
93 real(kind=dp), allocatable, dimension(:), public :: mfactLocal
94
95 logical, allocatable, dimension(:) :: simHasAtypeMap
96 #ifdef IS_MPI
97 logical, allocatable, dimension(:) :: simHasAtypeMapTemp
98 #endif
99
100 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
101 real(kind=dp), save :: DangerRcut
102 logical, public, save :: boxIsOrthorhombic
103
104 public :: SimulationSetup
105 public :: getNlocal
106 public :: setBox
107 public :: checkBox
108 public :: SimUsesPBC
109 public :: SimUsesAtomicVirial
110
111 public :: SimUsesDirectionalAtoms
112 public :: SimUsesLennardJones
113 public :: SimUsesElectrostatics
114 public :: SimUsesCharges
115 public :: SimUsesDipoles
116 public :: SimUsesSticky
117 public :: SimUsesStickyPower
118 public :: SimUsesGayBerne
119 public :: SimUsesEAM
120 public :: SimUsesShapes
121 public :: SimUsesFLARB
122 public :: SimUsesRF
123 public :: SimUsesSF
124 public :: SimUsesSP
125 public :: SimUsesBoxDipole
126 public :: SimRequiresPrepairCalc
127 public :: SimRequiresPostpairCalc
128 public :: SimHasAtype
129 public :: SimUsesSC
130 public :: SimUsesMNM
131 public :: setHmatDangerousRcutValue
132
133 contains
134
135 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
136 CnExcludes, Cexcludes, CnOneTwo, ConeTwo, CnOneThree, ConeThree, &
137 CnOneFour, ConeFour, CmolMembership, Cmfact, CnGroups, &
138 CglobalGroupMembership, status)
139
140 type (simtype) :: setThisSim
141 integer, intent(inout) :: CnGlobal, CnLocal
142 integer, dimension(CnLocal), intent(inout) :: c_idents
143
144 integer :: CnExcludes
145 integer, dimension(2,CnExcludes), intent(in) :: Cexcludes
146 integer :: CnOneTwo
147 integer, dimension(2,CnOneTwo), intent(in) :: ConeTwo
148 integer :: CnOneThree
149 integer, dimension(2,CnOneThree), intent(in) :: ConeThree
150 integer :: CnOneFour
151 integer, dimension(2,CnOneFour), intent(in) :: ConeFour
152
153 integer, dimension(CnGlobal),intent(in) :: CmolMembership
154 !! Result status, success = 0, status = -1
155 integer, intent(out) :: status
156 integer :: i, j, me, thisStat, alloc_stat, myNode, id1, id2
157 integer :: ia, jend
158
159 !! mass factors used for molecular cutoffs
160 real ( kind = dp ), dimension(CnLocal) :: Cmfact
161 integer, intent(in):: CnGroups
162 integer, dimension(CnGlobal), intent(in):: CglobalGroupMembership
163 integer :: maxSkipsForLocalAtom, maxToposForAtom, glPointer
164 integer :: maxSkipsForRowAtom
165
166 #ifdef IS_MPI
167 integer, allocatable, dimension(:) :: c_idents_Row
168 integer, allocatable, dimension(:) :: c_idents_Col
169 integer :: nAtomsInRow, nGroupsInRow, aid
170 integer :: nAtomsInCol, nGroupsInCol, gid
171 #endif
172
173 simulation_setup_complete = .false.
174 status = 0
175
176 ! copy C struct into fortran type
177
178 nLocal = CnLocal
179 nGlobal = CnGlobal
180 nGroups = CnGroups
181
182 thisSim = setThisSim
183
184 nExcludes = CnExcludes
185
186 call InitializeForceGlobals(nLocal, thisStat)
187 if (thisStat /= 0) then
188 write(default_error,*) "SimSetup: InitializeForceGlobals error"
189 status = -1
190 return
191 endif
192
193 call InitializeSimGlobals(thisStat)
194 if (thisStat /= 0) then
195 write(default_error,*) "SimSetup: InitializeSimGlobals error"
196 status = -1
197 return
198 endif
199
200 #ifdef IS_MPI
201 ! We can only set up forces if mpiSimulation has been setup.
202 if (.not. isMPISimSet()) then
203 write(default_error,*) "MPI is not set"
204 status = -1
205 return
206 endif
207 nAtomsInRow = getNatomsInRow(plan_atom_row)
208 nAtomsInCol = getNatomsInCol(plan_atom_col)
209 nGroupsInRow = getNgroupsInRow(plan_group_row)
210 nGroupsInCol = getNgroupsInCol(plan_group_col)
211 mynode = getMyNode()
212
213 call gather(c_idents, c_idents_Row, plan_atom_row)
214 call gather(c_idents, c_idents_Col, plan_atom_col)
215
216 do i = 1, nAtomsInRow
217 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
218 atid_Row(i) = me
219 enddo
220
221 do i = 1, nAtomsInCol
222 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
223 atid_Col(i) = me
224 enddo
225
226 #endif
227
228 #ifdef IS_MPI
229 allocate(groupStartRow(nGroupsInRow+1),stat=alloc_stat)
230 if (alloc_stat /= 0 ) then
231 status = -1
232 return
233 endif
234 allocate(groupStartCol(nGroupsInCol+1),stat=alloc_stat)
235 if (alloc_stat /= 0 ) then
236 status = -1
237 return
238 endif
239 allocate(groupListRow(nAtomsInRow),stat=alloc_stat)
240 if (alloc_stat /= 0 ) then
241 status = -1
242 return
243 endif
244 allocate(groupListCol(nAtomsInCol),stat=alloc_stat)
245 if (alloc_stat /= 0 ) then
246 status = -1
247 return
248 endif
249 allocate(mfactRow(nAtomsInRow),stat=alloc_stat)
250 if (alloc_stat /= 0 ) then
251 status = -1
252 return
253 endif
254 allocate(mfactCol(nAtomsInCol),stat=alloc_stat)
255 if (alloc_stat /= 0 ) then
256 status = -1
257 return
258 endif
259 allocate(mfactLocal(nLocal),stat=alloc_stat)
260 if (alloc_stat /= 0 ) then
261 status = -1
262 return
263 endif
264
265 glPointer = 1
266
267 do i = 1, nGroupsInRow
268
269 gid = GroupRowToGlobal(i)
270 groupStartRow(i) = glPointer
271
272 do j = 1, nAtomsInRow
273 aid = AtomRowToGlobal(j)
274 if (CglobalGroupMembership(aid) .eq. gid) then
275 groupListRow(glPointer) = j
276 glPointer = glPointer + 1
277 endif
278 enddo
279 enddo
280 groupStartRow(nGroupsInRow+1) = nAtomsInRow + 1
281
282 glPointer = 1
283
284 do i = 1, nGroupsInCol
285
286 gid = GroupColToGlobal(i)
287 groupStartCol(i) = glPointer
288
289 do j = 1, nAtomsInCol
290 aid = AtomColToGlobal(j)
291 if (CglobalGroupMembership(aid) .eq. gid) then
292 groupListCol(glPointer) = j
293 glPointer = glPointer + 1
294 endif
295 enddo
296 enddo
297 groupStartCol(nGroupsInCol+1) = nAtomsInCol + 1
298
299 mfactLocal = Cmfact
300
301 call gather(mfactLocal, mfactRow, plan_atom_row)
302 call gather(mfactLocal, mfactCol, plan_atom_col)
303
304 if (allocated(mfactLocal)) then
305 deallocate(mfactLocal)
306 end if
307 #else
308 allocate(groupStartRow(nGroups+1),stat=alloc_stat)
309 if (alloc_stat /= 0 ) then
310 status = -1
311 return
312 endif
313 allocate(groupStartCol(nGroups+1),stat=alloc_stat)
314 if (alloc_stat /= 0 ) then
315 status = -1
316 return
317 endif
318 allocate(groupListRow(nLocal),stat=alloc_stat)
319 if (alloc_stat /= 0 ) then
320 status = -1
321 return
322 endif
323 allocate(groupListCol(nLocal),stat=alloc_stat)
324 if (alloc_stat /= 0 ) then
325 status = -1
326 return
327 endif
328 allocate(mfactRow(nLocal),stat=alloc_stat)
329 if (alloc_stat /= 0 ) then
330 status = -1
331 return
332 endif
333 allocate(mfactCol(nLocal),stat=alloc_stat)
334 if (alloc_stat /= 0 ) then
335 status = -1
336 return
337 endif
338 allocate(mfactLocal(nLocal),stat=alloc_stat)
339 if (alloc_stat /= 0 ) then
340 status = -1
341 return
342 endif
343
344 glPointer = 1
345 do i = 1, nGroups
346 groupStartRow(i) = glPointer
347 groupStartCol(i) = glPointer
348 do j = 1, nLocal
349 if (CglobalGroupMembership(j) .eq. i) then
350 groupListRow(glPointer) = j
351 groupListCol(glPointer) = j
352 glPointer = glPointer + 1
353 endif
354 enddo
355 enddo
356 groupStartRow(nGroups+1) = nLocal + 1
357 groupStartCol(nGroups+1) = nLocal + 1
358
359 do i = 1, nLocal
360 mfactRow(i) = Cmfact(i)
361 mfactCol(i) = Cmfact(i)
362 end do
363
364 #endif
365
366 ! We build the local atid's for both mpi and nonmpi
367 do i = 1, nLocal
368 me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
369 atid(i) = me
370 c_idents_local(i) = c_idents(i)
371 enddo
372
373 do i = 1, nExcludes
374 excludes(1,i) = Cexcludes(1,i)
375 excludes(2,i) = Cexcludes(2,i)
376 enddo
377
378 #ifdef IS_MPI
379 allocate(nSkipsForRowAtom(nAtomsInRow), stat=alloc_stat)
380 #endif
381
382 allocate(nSkipsForLocalAtom(nLocal), stat=alloc_stat)
383
384 if (alloc_stat /= 0 ) then
385 thisStat = -1
386 write(*,*) 'Could not allocate nSkipsForAtom array'
387 return
388 endif
389
390 #ifdef IS_MPI
391 maxSkipsForRowAtom = 0
392 do j = 1, nAtomsInRow
393 nSkipsForRowAtom(j) = 0
394 id1 = AtomRowToGlobal(j)
395 do i = 1, nExcludes
396 if (excludes(1,i) .eq. id1 ) then
397 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
398 if (nSkipsForRowAtom(j) .gt. maxSkipsForRowAtom) then
399 maxSkipsForRowAtom = nSkipsForRowAtom(j)
400 endif
401 endif
402 if (excludes(2,i) .eq. id1 ) then
403 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
404 if (nSkipsForRowAtom(j) .gt. maxSkipsForRowAtom) then
405 maxSkipsForRowAtom = nSkipsForRowAtom(j)
406 endif
407 endif
408 end do
409 enddo
410 #endif
411 maxSkipsForLocalAtom = 0
412 do j = 1, nLocal
413 nSkipsForLocalAtom(j) = 0
414 #ifdef IS_MPI
415 id1 = AtomLocalToGlobal(j)
416 #else
417 id1 = j
418 #endif
419 do i = 1, nExcludes
420 if (excludes(1,i) .eq. id1 ) then
421 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
422 if (nSkipsForLocalAtom(j) .gt. maxSkipsForLocalAtom) then
423 maxSkipsForLocalAtom = nSkipsForLocalAtom(j)
424 endif
425 endif
426 if (excludes(2,i) .eq. id1 ) then
427 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
428 if (nSkipsForLocalAtom(j) .gt. maxSkipsForLocalAtom) then
429 maxSkipsForLocalAtom = nSkipsForLocalAtom(j)
430 endif
431 endif
432 end do
433 enddo
434
435 #ifdef IS_MPI
436 allocate(skipsForRowAtom(nAtomsInRow, maxSkipsForRowAtom), stat=alloc_stat)
437 #endif
438 allocate(skipsForLocalAtom(nLocal, maxSkipsForLocalAtom), stat=alloc_stat)
439
440 if (alloc_stat /= 0 ) then
441 write(*,*) 'Could not allocate skipsForAtom arrays'
442 return
443 endif
444
445 #ifdef IS_MPI
446 do j = 1, nAtomsInRow
447 nSkipsForRowAtom(j) = 0
448 id1 = AtomRowToGlobal(j)
449 do i = 1, nExcludes
450 if (excludes(1,i) .eq. id1 ) then
451 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
452 ! exclude lists have global ID's
453 id2 = excludes(2,i)
454 skipsForRowAtom(j, nSkipsForRowAtom(j)) = id2
455 endif
456 if (excludes(2, i) .eq. id1 ) then
457 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
458 ! exclude lists have global ID's
459 id2 = excludes(1,i)
460 skipsForRowAtom(j, nSkipsForRowAtom(j)) = id2
461 endif
462 end do
463 enddo
464 #endif
465 do j = 1, nLocal
466 nSkipsForLocalAtom(j) = 0
467 #ifdef IS_MPI
468 id1 = AtomLocalToGlobal(j)
469 #else
470 id1 = j
471 #endif
472 do i = 1, nExcludes
473 if (excludes(1,i) .eq. id1 ) then
474 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
475 ! exclude lists have global ID's
476 #ifdef IS_MPI
477 id2 = AtomGlobalToLocal(excludes(2,i))
478 #else
479 id2 = excludes(2,i)
480 #endif
481 skipsForLocalAtom(j, nSkipsForLocalAtom(j)) = id2
482 endif
483 if (excludes(2, i) .eq. id1 ) then
484 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
485 ! exclude lists have global ID's
486 #ifdef IS_MPI
487 id2 = AtomGlobalToLocal(excludes(1,i))
488 #else
489 id2 = excludes(1,i)
490 #endif
491 skipsForLocalAtom(j, nSkipsForLocalAtom(j)) = id2
492 endif
493 end do
494 enddo
495
496 do i = 1, nGlobal
497 molMemberShipList(i) = CmolMembership(i)
498 enddo
499
500 #ifdef IS_MPI
501 allocate(nTopoPairsForAtom(nAtomsInRow), stat=alloc_stat)
502 #else
503 allocate(nTopoPairsForAtom(nLocal), stat=alloc_stat)
504 #endif
505 if (alloc_stat /= 0 ) then
506 thisStat = -1
507 write(*,*) 'Could not allocate nTopoPairsForAtom array'
508 return
509 endif
510
511 #ifdef IS_MPI
512 jend = nAtomsInRow
513 #else
514 jend = nLocal
515 #endif
516
517 do j = 1, jend
518 nTopoPairsForAtom(j) = 0
519 #ifdef IS_MPI
520 id1 = AtomRowToGlobal(j)
521 #else
522 id1 = j
523 #endif
524 do i = 1, CnOneTwo
525 if (ConeTwo(1,i) .eq. id1 ) then
526 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
527 endif
528 if (ConeTwo(2,i) .eq. id1 ) then
529 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
530 endif
531 end do
532
533 do i = 1, CnOneThree
534 if (ConeThree(1,i) .eq. id1 ) then
535 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
536 endif
537 if (ConeThree(2,i) .eq. id1 ) then
538 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
539 endif
540 end do
541
542 do i = 1, CnOneFour
543 if (ConeFour(1,i) .eq. id1 ) then
544 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
545 endif
546 if (ConeFour(2,i) .eq. id1 ) then
547 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
548 endif
549 end do
550 enddo
551
552 maxToposForAtom = maxval(nTopoPairsForAtom)
553 #ifdef IS_MPI
554 allocate(toposForAtom(nAtomsInRow, maxToposForAtom), stat=alloc_stat)
555 allocate(topoDistance(nAtomsInRow, maxToposForAtom), stat=alloc_stat)
556 #else
557 allocate(toposForAtom(nLocal, maxToposForAtom), stat=alloc_stat)
558 allocate(topoDistance(nLocal, maxToposForAtom), stat=alloc_stat)
559 #endif
560 if (alloc_stat /= 0 ) then
561 write(*,*) 'Could not allocate topoDistance array'
562 return
563 endif
564
565 #ifdef IS_MPI
566 jend = nAtomsInRow
567 #else
568 jend = nLocal
569 #endif
570 do j = 1, jend
571 nTopoPairsForAtom(j) = 0
572 #ifdef IS_MPI
573 id1 = AtomRowToGlobal(j)
574 #else
575 id1 = j
576 #endif
577 do i = 1, CnOneTwo
578 if (ConeTwo(1,i) .eq. id1 ) then
579 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
580 id2 = ConeTwo(2,i)
581 toposForAtom(j, nTopoPairsForAtom(j)) = id2
582 topoDistance(j, nTopoPairsForAtom(j)) = 1
583 endif
584 if (ConeTwo(2, i) .eq. id1 ) then
585 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
586 id2 = ConeTwo(1,i)
587 toposForAtom(j, nTopoPairsForAtom(j)) = id2
588 topoDistance(j, nTopoPairsForAtom(j)) = 1
589 endif
590 end do
591
592 do i = 1, CnOneThree
593 if (ConeThree(1,i) .eq. id1 ) then
594 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
595 id2 = ConeThree(2,i)
596 toposForAtom(j, nTopoPairsForAtom(j)) = id2
597 topoDistance(j, nTopoPairsForAtom(j)) = 2
598 endif
599 if (ConeThree(2, i) .eq. id1 ) then
600 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
601 id2 = ConeThree(1,i)
602 toposForAtom(j, nTopoPairsForAtom(j)) = id2
603 topoDistance(j, nTopoPairsForAtom(j)) = 2
604 endif
605 end do
606
607 do i = 1, CnOneFour
608 if (ConeFour(1,i) .eq. id1 ) then
609 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
610 id2 = ConeFour(2,i)
611 toposForAtom(j, nTopoPairsForAtom(j)) = id2
612 topoDistance(j, nTopoPairsForAtom(j)) = 3
613 endif
614 if (ConeFour(2, i) .eq. id1 ) then
615 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
616 id2 = ConeFour(1,i)
617 toposForAtom(j, nTopoPairsForAtom(j)) = id2
618 topoDistance(j, nTopoPairsForAtom(j)) = 3
619 endif
620 end do
621 enddo
622
623 call createSimHasAtype(alloc_stat)
624 if (alloc_stat /= 0) then
625 status = -1
626 end if
627
628 if (status == 0) simulation_setup_complete = .true.
629
630 end subroutine SimulationSetup
631
632 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
633 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
634 integer :: cBoxIsOrthorhombic
635 integer :: smallest, status
636
637 Hmat = cHmat
638 HmatInv = cHmatInv
639 if (cBoxIsOrthorhombic .eq. 0 ) then
640 boxIsOrthorhombic = .false.
641 else
642 boxIsOrthorhombic = .true.
643 endif
644
645 call checkBox()
646 return
647 end subroutine setBox
648
649 subroutine checkBox()
650
651 integer :: i
652 real(kind=dp), dimension(3) :: hx, hy, hz, ax, ay, az, piped
653 character(len = statusMsgSize) :: errMsg
654
655 hx = Hmat(1,:)
656 hy = Hmat(2,:)
657 hz = Hmat(3,:)
658
659 ax = cross_product(hy, hz)
660 ay = cross_product(hx, hz)
661 az = cross_product(hx, hy)
662
663 ax = ax / length(ax)
664 ay = ay / length(ay)
665 az = az / length(az)
666
667 piped(1) = abs(dot_product(ax, hx))
668 piped(2) = abs(dot_product(ay, hy))
669 piped(3) = abs(dot_product(az, hz))
670
671 do i = 1, 3
672 if ((0.5_dp * piped(i)).lt.DangerRcut) then
673 write(errMsg, '(a94,f9.4,a1)') 'One of the dimensions of the Periodic ' &
674 // 'Box is smaller than ' // newline // tab // &
675 'the largest cutoff radius' // &
676 ' (rCut = ', DangerRcut, ')'
677 call handleError("checkBox", errMsg)
678
679 end if
680 enddo
681 return
682 end subroutine checkBox
683
684 function SimUsesPBC() result(doesit)
685 logical :: doesit
686 doesit = thisSim%SIM_uses_PBC
687 end function SimUsesPBC
688
689 function SimUsesAtomicVirial() result(doesit)
690 logical :: doesit
691 doesit = thisSim%SIM_uses_AtomicVirial
692 end function SimUsesAtomicVirial
693
694 function SimUsesDirectionalAtoms() result(doesit)
695 logical :: doesit
696 doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_Sticky .or. &
697 thisSim%SIM_uses_StickyPower .or. &
698 thisSim%SIM_uses_GayBerne .or. thisSim%SIM_uses_Shapes
699 end function SimUsesDirectionalAtoms
700
701 function SimUsesLennardJones() result(doesit)
702 logical :: doesit
703 doesit = thisSim%SIM_uses_LennardJones
704 end function SimUsesLennardJones
705
706 function SimUsesElectrostatics() result(doesit)
707 logical :: doesit
708 doesit = thisSim%SIM_uses_Electrostatics
709 end function SimUsesElectrostatics
710
711 function SimUsesCharges() result(doesit)
712 logical :: doesit
713 doesit = thisSim%SIM_uses_Charges
714 end function SimUsesCharges
715
716 function SimUsesDipoles() result(doesit)
717 logical :: doesit
718 doesit = thisSim%SIM_uses_Dipoles
719 end function SimUsesDipoles
720
721 function SimUsesSticky() result(doesit)
722 logical :: doesit
723 doesit = thisSim%SIM_uses_Sticky
724 end function SimUsesSticky
725
726 function SimUsesStickyPower() result(doesit)
727 logical :: doesit
728 doesit = thisSim%SIM_uses_StickyPower
729 end function SimUsesStickyPower
730
731 function SimUsesGayBerne() result(doesit)
732 logical :: doesit
733 doesit = thisSim%SIM_uses_GayBerne
734 end function SimUsesGayBerne
735
736 function SimUsesEAM() result(doesit)
737 logical :: doesit
738 doesit = thisSim%SIM_uses_EAM
739 end function SimUsesEAM
740
741
742 function SimUsesSC() result(doesit)
743 logical :: doesit
744 doesit = thisSim%SIM_uses_SC
745 end function SimUsesSC
746
747 function SimUsesMNM() result(doesit)
748 logical :: doesit
749 doesit = thisSim%SIM_uses_MNM
750 end function SimUsesMNM
751
752
753 function SimUsesShapes() result(doesit)
754 logical :: doesit
755 doesit = thisSim%SIM_uses_Shapes
756 end function SimUsesShapes
757
758 function SimUsesFLARB() result(doesit)
759 logical :: doesit
760 doesit = thisSim%SIM_uses_FLARB
761 end function SimUsesFLARB
762
763 function SimUsesRF() result(doesit)
764 logical :: doesit
765 doesit = thisSim%SIM_uses_RF
766 end function SimUsesRF
767
768 function SimUsesSF() result(doesit)
769 logical :: doesit
770 doesit = thisSim%SIM_uses_SF
771 end function SimUsesSF
772
773 function SimUsesSP() result(doesit)
774 logical :: doesit
775 doesit = thisSim%SIM_uses_SP
776 end function SimUsesSP
777
778 function SimUsesBoxDipole() result(doesit)
779 logical :: doesit
780 doesit = thisSim%SIM_uses_BoxDipole
781 end function SimUsesBoxDipole
782
783 function SimRequiresPrepairCalc() result(doesit)
784 logical :: doesit
785 doesit = thisSim%SIM_uses_EAM .or. thisSim%SIM_uses_SC
786 end function SimRequiresPrepairCalc
787
788 function SimRequiresPostpairCalc() result(doesit)
789 logical :: doesit
790 doesit = thisSim%SIM_uses_RF .or. thisSim%SIM_uses_SF &
791 .or. thisSim%SIM_uses_SP .or. thisSim%SIM_uses_BoxDipole
792 end function SimRequiresPostpairCalc
793
794 ! Function returns true if the simulation has this atype
795 function SimHasAtype(thisAtype) result(doesit)
796 logical :: doesit
797 integer :: thisAtype
798 doesit = .false.
799 if(.not.allocated(SimHasAtypeMap)) return
800
801 doesit = SimHasAtypeMap(thisAtype)
802
803 end function SimHasAtype
804
805 subroutine createSimHasAtype(status)
806 integer, intent(out) :: status
807 integer :: alloc_stat
808 integer :: me_i
809 integer :: mpiErrors
810 integer :: nAtypes
811 status = 0
812
813 nAtypes = getSize(atypes)
814 ! Setup logical map for atypes in simulation
815 if (.not.allocated(SimHasAtypeMap)) then
816 allocate(SimHasAtypeMap(nAtypes),stat=alloc_stat)
817 if (alloc_stat /= 0 ) then
818 status = -1
819 return
820 end if
821 SimHasAtypeMap = .false.
822 end if
823
824 ! Loop through the local atoms and grab the atypes present
825 do me_i = 1,nLocal
826 SimHasAtypeMap(atid(me_i)) = .true.
827 end do
828 ! For MPI, we need to know all possible atypes present in
829 ! simulation on all processors. Use LOR operation to set map.
830 #ifdef IS_MPI
831 if (.not.allocated(SimHasAtypeMapTemp)) then
832 allocate(SimHasAtypeMapTemp(nAtypes),stat=alloc_stat)
833 if (alloc_stat /= 0 ) then
834 status = -1
835 return
836 end if
837 end if
838 call mpi_allreduce(SimHasAtypeMap, SimHasAtypeMaptemp, nAtypes, &
839 mpi_logical, MPI_LOR, mpi_comm_world, mpiErrors)
840 simHasAtypeMap = simHasAtypeMapTemp
841 deallocate(simHasAtypeMapTemp)
842 #endif
843 end subroutine createSimHasAtype
844
845 subroutine InitializeSimGlobals(thisStat)
846 integer, intent(out) :: thisStat
847 integer :: alloc_stat
848
849 thisStat = 0
850
851 call FreeSimGlobals()
852
853 allocate(excludes(2,nExcludes), stat=alloc_stat)
854 if (alloc_stat /= 0 ) then
855 thisStat = -1
856 return
857 endif
858
859 allocate(molMembershipList(nGlobal), stat=alloc_stat)
860 if (alloc_stat /= 0 ) then
861 thisStat = -1
862 return
863 endif
864
865 end subroutine InitializeSimGlobals
866
867 subroutine FreeSimGlobals()
868
869 !We free in the opposite order in which we allocate in.
870 if (allocated(topoDistance)) deallocate(topoDistance)
871 if (allocated(toposForAtom)) deallocate(toposForAtom)
872 if (allocated(nTopoPairsForAtom)) deallocate(nTopoPairsForAtom)
873 if (allocated(skipsForLocalAtom)) deallocate(skipsForLocalAtom)
874 if (allocated(nSkipsForLocalAtom)) deallocate(nSkipsForLocalAtom)
875 if (allocated(skipsForRowAtom)) deallocate(skipsForRowAtom)
876 if (allocated(nSkipsForRowAtom)) deallocate(nSkipsForRowAtom)
877 if (allocated(mfactLocal)) deallocate(mfactLocal)
878 if (allocated(mfactCol)) deallocate(mfactCol)
879 if (allocated(mfactRow)) deallocate(mfactRow)
880 if (allocated(groupListCol)) deallocate(groupListCol)
881 if (allocated(groupListRow)) deallocate(groupListRow)
882 if (allocated(groupStartCol)) deallocate(groupStartCol)
883 if (allocated(groupStartRow)) deallocate(groupStartRow)
884 if (allocated(molMembershipList)) deallocate(molMembershipList)
885 if (allocated(excludes)) deallocate(excludes)
886
887 end subroutine FreeSimGlobals
888
889 pure function getNlocal() result(n)
890 integer :: n
891 n = nLocal
892 end function getNlocal
893
894 subroutine setHmatDangerousRcutValue(dangerWillRobinson)
895 real(kind=dp), intent(in) :: dangerWillRobinson
896 DangerRcut = dangerWillRobinson
897
898 call checkBox()
899
900 return
901 end subroutine setHmatDangerousRcutValue
902
903
904
905 end module simulation

Properties

Name Value
svn:keywords Author Id Revision Date