ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/UseTheForce/DarkSide/simulation.F90
Revision: 1530
Committed: Tue Dec 28 21:47:55 2010 UTC (14 years, 4 months ago) by gezelter
File size: 23593 byte(s)
Log Message:
Moved switching functions and force options over to the C++ side, and
removed them from Fortran.

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 #ifdef IS_MPI
53 use mpiSimulation
54 #endif
55
56 implicit none
57 PRIVATE
58
59 #define __FORTRAN90
60 #include "brains/fSimulation.h"
61
62 type (simtype), public, save :: thisSim
63
64 logical, save :: simulation_setup_complete = .false.
65
66 integer, public, save :: nLocal, nGlobal
67 integer, public, save :: nGroups, nGroupGlobal
68 integer, public, save :: nExcludes = 0
69 integer, public, save :: nOneTwo = 0
70 integer, public, save :: nOneThree = 0
71 integer, public, save :: nOneFour = 0
72
73 integer, allocatable, dimension(:,:), public :: excludes
74 integer, allocatable, dimension(:), public :: molMembershipList
75 integer, allocatable, dimension(:), public :: groupListRow
76 integer, allocatable, dimension(:), public :: groupStartRow
77 integer, allocatable, dimension(:), public :: groupListCol
78 integer, allocatable, dimension(:), public :: groupStartCol
79 integer, allocatable, dimension(:), public :: groupListLocal
80 integer, allocatable, dimension(:), public :: groupStartLocal
81 integer, allocatable, dimension(:), public :: nSkipsForLocalAtom
82 integer, allocatable, dimension(:,:), public :: skipsForLocalAtom
83 integer, allocatable, dimension(:), public :: nSkipsForRowAtom
84 integer, allocatable, dimension(:,:), public :: skipsForRowAtom
85 integer, allocatable, dimension(:), public :: nTopoPairsForAtom
86 integer, allocatable, dimension(:,:), public :: toposForAtom
87 integer, allocatable, dimension(:,:), public :: topoDistance
88 real(kind=dp), allocatable, dimension(:), public :: mfactRow
89 real(kind=dp), allocatable, dimension(:), public :: mfactCol
90 real(kind=dp), allocatable, dimension(:), public :: mfactLocal
91
92 real(kind=dp), public, dimension(3,3), save :: Hmat, HmatInv
93 real(kind=dp), save :: DangerRcut
94 logical, public, save :: boxIsOrthorhombic
95
96 public :: SimulationSetup
97 public :: getNlocal
98 public :: setBox
99 public :: checkBox
100 public :: SimUsesPBC
101 public :: SimUsesAtomicVirial
102 public :: SimUsesDirectionalAtoms
103 public :: SimUsesMetallicAtoms
104 public :: SimRequiresSkipCorrection
105 public :: SimRequiresSelfCorrection
106 public :: setHmatDangerousRcutValue
107
108 contains
109
110 subroutine SimulationSetup(setThisSim, CnGlobal, CnLocal, c_idents, &
111 CnExcludes, Cexcludes, CnOneTwo, ConeTwo, CnOneThree, ConeThree, &
112 CnOneFour, ConeFour, CmolMembership, Cmfact, CnGroups, &
113 CglobalGroupMembership, status)
114
115 type (simtype) :: setThisSim
116 integer, intent(inout) :: CnGlobal, CnLocal
117 integer, dimension(CnLocal), intent(inout) :: c_idents
118
119 integer :: CnExcludes
120 integer, dimension(2,CnExcludes), intent(in) :: Cexcludes
121 integer :: CnOneTwo
122 integer, dimension(2,CnOneTwo), intent(in) :: ConeTwo
123 integer :: CnOneThree
124 integer, dimension(2,CnOneThree), intent(in) :: ConeThree
125 integer :: CnOneFour
126 integer, dimension(2,CnOneFour), intent(in) :: ConeFour
127
128 integer, dimension(CnGlobal),intent(in) :: CmolMembership
129 !! Result status, success = 0, status = -1
130 integer, intent(out) :: status
131 integer :: i, j, me, thisStat, alloc_stat, myNode, id1, id2
132 integer :: ia, jend
133
134 !! mass factors used for molecular cutoffs
135 real ( kind = dp ), dimension(CnLocal) :: Cmfact
136 integer, intent(in):: CnGroups
137 integer, dimension(CnGlobal), intent(in):: CglobalGroupMembership
138 integer :: maxSkipsForLocalAtom, maxToposForAtom, glPointer
139 integer :: maxSkipsForRowAtom
140
141 #ifdef IS_MPI
142 integer, allocatable, dimension(:) :: c_idents_Row
143 integer, allocatable, dimension(:) :: c_idents_Col
144 integer :: nAtomsInRow, nGroupsInRow, aid
145 integer :: nAtomsInCol, nGroupsInCol, gid
146 #endif
147
148 simulation_setup_complete = .false.
149 status = 0
150
151 ! copy C struct into fortran type
152
153 nLocal = CnLocal
154 nGlobal = CnGlobal
155 nGroups = CnGroups
156
157 thisSim = setThisSim
158
159 nExcludes = CnExcludes
160
161 call InitializeForceGlobals(nLocal, thisStat)
162 if (thisStat /= 0) then
163 write(default_error,*) "SimSetup: InitializeForceGlobals error"
164 status = -1
165 return
166 endif
167
168 call InitializeSimGlobals(thisStat)
169 if (thisStat /= 0) then
170 write(default_error,*) "SimSetup: InitializeSimGlobals error"
171 status = -1
172 return
173 endif
174
175 #ifdef IS_MPI
176 ! We can only set up forces if mpiSimulation has been setup.
177 if (.not. isMPISimSet()) then
178 write(default_error,*) "MPI is not set"
179 status = -1
180 return
181 endif
182 nAtomsInRow = getNatomsInRow(plan_atom_row)
183 nAtomsInCol = getNatomsInCol(plan_atom_col)
184 nGroupsInRow = getNgroupsInRow(plan_group_row)
185 nGroupsInCol = getNgroupsInCol(plan_group_col)
186 mynode = getMyNode()
187
188 call gather(c_idents, c_idents_Row, plan_atom_row)
189 call gather(c_idents, c_idents_Col, plan_atom_col)
190
191 do i = 1, nAtomsInRow
192 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i))
193 atid_Row(i) = me
194 enddo
195
196 do i = 1, nAtomsInCol
197 me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i))
198 atid_Col(i) = me
199 enddo
200
201 #endif
202
203 #ifdef IS_MPI
204 allocate(groupStartRow(nGroupsInRow+1),stat=alloc_stat)
205 if (alloc_stat /= 0 ) then
206 status = -1
207 return
208 endif
209 allocate(groupStartCol(nGroupsInCol+1),stat=alloc_stat)
210 if (alloc_stat /= 0 ) then
211 status = -1
212 return
213 endif
214 allocate(groupListRow(nAtomsInRow),stat=alloc_stat)
215 if (alloc_stat /= 0 ) then
216 status = -1
217 return
218 endif
219 allocate(groupListCol(nAtomsInCol),stat=alloc_stat)
220 if (alloc_stat /= 0 ) then
221 status = -1
222 return
223 endif
224 allocate(mfactRow(nAtomsInRow),stat=alloc_stat)
225 if (alloc_stat /= 0 ) then
226 status = -1
227 return
228 endif
229 allocate(mfactCol(nAtomsInCol),stat=alloc_stat)
230 if (alloc_stat /= 0 ) then
231 status = -1
232 return
233 endif
234 allocate(mfactLocal(nLocal),stat=alloc_stat)
235 if (alloc_stat /= 0 ) then
236 status = -1
237 return
238 endif
239
240 glPointer = 1
241
242 do i = 1, nGroupsInRow
243
244 gid = GroupRowToGlobal(i)
245 groupStartRow(i) = glPointer
246
247 do j = 1, nAtomsInRow
248 aid = AtomRowToGlobal(j)
249 if (CglobalGroupMembership(aid) .eq. gid) then
250 groupListRow(glPointer) = j
251 glPointer = glPointer + 1
252 endif
253 enddo
254 enddo
255 groupStartRow(nGroupsInRow+1) = nAtomsInRow + 1
256
257 glPointer = 1
258
259 do i = 1, nGroupsInCol
260
261 gid = GroupColToGlobal(i)
262 groupStartCol(i) = glPointer
263
264 do j = 1, nAtomsInCol
265 aid = AtomColToGlobal(j)
266 if (CglobalGroupMembership(aid) .eq. gid) then
267 groupListCol(glPointer) = j
268 glPointer = glPointer + 1
269 endif
270 enddo
271 enddo
272 groupStartCol(nGroupsInCol+1) = nAtomsInCol + 1
273
274 mfactLocal = Cmfact
275
276 call gather(mfactLocal, mfactRow, plan_atom_row)
277 call gather(mfactLocal, mfactCol, plan_atom_col)
278
279 if (allocated(mfactLocal)) then
280 deallocate(mfactLocal)
281 end if
282 #else
283 allocate(groupStartRow(nGroups+1),stat=alloc_stat)
284 if (alloc_stat /= 0 ) then
285 status = -1
286 return
287 endif
288 allocate(groupStartCol(nGroups+1),stat=alloc_stat)
289 if (alloc_stat /= 0 ) then
290 status = -1
291 return
292 endif
293 allocate(groupListRow(nLocal),stat=alloc_stat)
294 if (alloc_stat /= 0 ) then
295 status = -1
296 return
297 endif
298 allocate(groupListCol(nLocal),stat=alloc_stat)
299 if (alloc_stat /= 0 ) then
300 status = -1
301 return
302 endif
303 allocate(mfactRow(nLocal),stat=alloc_stat)
304 if (alloc_stat /= 0 ) then
305 status = -1
306 return
307 endif
308 allocate(mfactCol(nLocal),stat=alloc_stat)
309 if (alloc_stat /= 0 ) then
310 status = -1
311 return
312 endif
313 allocate(mfactLocal(nLocal),stat=alloc_stat)
314 if (alloc_stat /= 0 ) then
315 status = -1
316 return
317 endif
318
319 glPointer = 1
320 do i = 1, nGroups
321 groupStartRow(i) = glPointer
322 groupStartCol(i) = glPointer
323 do j = 1, nLocal
324 if (CglobalGroupMembership(j) .eq. i) then
325 groupListRow(glPointer) = j
326 groupListCol(glPointer) = j
327 glPointer = glPointer + 1
328 endif
329 enddo
330 enddo
331 groupStartRow(nGroups+1) = nLocal + 1
332 groupStartCol(nGroups+1) = nLocal + 1
333
334 do i = 1, nLocal
335 mfactRow(i) = Cmfact(i)
336 mfactCol(i) = Cmfact(i)
337 end do
338
339 #endif
340
341 ! We build the local atid's for both mpi and nonmpi
342 do i = 1, nLocal
343 me = getFirstMatchingElement(atypes, "c_ident", c_idents(i))
344 atid(i) = me
345 c_idents_local(i) = c_idents(i)
346 enddo
347
348 do i = 1, nExcludes
349 excludes(1,i) = Cexcludes(1,i)
350 excludes(2,i) = Cexcludes(2,i)
351 enddo
352
353 #ifdef IS_MPI
354 allocate(nSkipsForRowAtom(nAtomsInRow), stat=alloc_stat)
355 #endif
356
357 allocate(nSkipsForLocalAtom(nLocal), stat=alloc_stat)
358
359 if (alloc_stat /= 0 ) then
360 thisStat = -1
361 write(*,*) 'Could not allocate nSkipsForAtom array'
362 return
363 endif
364
365 #ifdef IS_MPI
366 maxSkipsForRowAtom = 0
367 do j = 1, nAtomsInRow
368 nSkipsForRowAtom(j) = 0
369 id1 = AtomRowToGlobal(j)
370 do i = 1, nExcludes
371 if (excludes(1,i) .eq. id1 ) then
372 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
373 if (nSkipsForRowAtom(j) .gt. maxSkipsForRowAtom) then
374 maxSkipsForRowAtom = nSkipsForRowAtom(j)
375 endif
376 endif
377 if (excludes(2,i) .eq. id1 ) then
378 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
379 if (nSkipsForRowAtom(j) .gt. maxSkipsForRowAtom) then
380 maxSkipsForRowAtom = nSkipsForRowAtom(j)
381 endif
382 endif
383 end do
384 enddo
385 #endif
386 maxSkipsForLocalAtom = 0
387 do j = 1, nLocal
388 nSkipsForLocalAtom(j) = 0
389 #ifdef IS_MPI
390 id1 = AtomLocalToGlobal(j)
391 #else
392 id1 = j
393 #endif
394 do i = 1, nExcludes
395 if (excludes(1,i) .eq. id1 ) then
396 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
397 if (nSkipsForLocalAtom(j) .gt. maxSkipsForLocalAtom) then
398 maxSkipsForLocalAtom = nSkipsForLocalAtom(j)
399 endif
400 endif
401 if (excludes(2,i) .eq. id1 ) then
402 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
403 if (nSkipsForLocalAtom(j) .gt. maxSkipsForLocalAtom) then
404 maxSkipsForLocalAtom = nSkipsForLocalAtom(j)
405 endif
406 endif
407 end do
408 enddo
409
410 #ifdef IS_MPI
411 allocate(skipsForRowAtom(nAtomsInRow, maxSkipsForRowAtom), stat=alloc_stat)
412 #endif
413 allocate(skipsForLocalAtom(nLocal, maxSkipsForLocalAtom), stat=alloc_stat)
414
415 if (alloc_stat /= 0 ) then
416 write(*,*) 'Could not allocate skipsForAtom arrays'
417 return
418 endif
419
420 #ifdef IS_MPI
421 do j = 1, nAtomsInRow
422 nSkipsForRowAtom(j) = 0
423 id1 = AtomRowToGlobal(j)
424 do i = 1, nExcludes
425 if (excludes(1,i) .eq. id1 ) then
426 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
427 ! exclude lists have global ID's
428 id2 = excludes(2,i)
429 skipsForRowAtom(j, nSkipsForRowAtom(j)) = id2
430 endif
431 if (excludes(2, i) .eq. id1 ) then
432 nSkipsForRowAtom(j) = nSkipsForRowAtom(j) + 1
433 ! exclude lists have global ID's
434 id2 = excludes(1,i)
435 skipsForRowAtom(j, nSkipsForRowAtom(j)) = id2
436 endif
437 end do
438 enddo
439 #endif
440 do j = 1, nLocal
441 nSkipsForLocalAtom(j) = 0
442 #ifdef IS_MPI
443 id1 = AtomLocalToGlobal(j)
444 #else
445 id1 = j
446 #endif
447 do i = 1, nExcludes
448 if (excludes(1,i) .eq. id1 ) then
449 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
450 ! exclude lists have global ID's
451 #ifdef IS_MPI
452 id2 = AtomGlobalToLocal(excludes(2,i))
453 #else
454 id2 = excludes(2,i)
455 #endif
456 skipsForLocalAtom(j, nSkipsForLocalAtom(j)) = id2
457 endif
458 if (excludes(2, i) .eq. id1 ) then
459 nSkipsForLocalAtom(j) = nSkipsForLocalAtom(j) + 1
460 ! exclude lists have global ID's
461 #ifdef IS_MPI
462 id2 = AtomGlobalToLocal(excludes(1,i))
463 #else
464 id2 = excludes(1,i)
465 #endif
466 skipsForLocalAtom(j, nSkipsForLocalAtom(j)) = id2
467 endif
468 end do
469 enddo
470
471 do i = 1, nGlobal
472 molMemberShipList(i) = CmolMembership(i)
473 enddo
474
475 #ifdef IS_MPI
476 allocate(nTopoPairsForAtom(nAtomsInRow), stat=alloc_stat)
477 #else
478 allocate(nTopoPairsForAtom(nLocal), stat=alloc_stat)
479 #endif
480 if (alloc_stat /= 0 ) then
481 thisStat = -1
482 write(*,*) 'Could not allocate nTopoPairsForAtom array'
483 return
484 endif
485
486 #ifdef IS_MPI
487 jend = nAtomsInRow
488 #else
489 jend = nLocal
490 #endif
491
492 do j = 1, jend
493 nTopoPairsForAtom(j) = 0
494 #ifdef IS_MPI
495 id1 = AtomRowToGlobal(j)
496 #else
497 id1 = j
498 #endif
499 do i = 1, CnOneTwo
500 if (ConeTwo(1,i) .eq. id1 ) then
501 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
502 endif
503 if (ConeTwo(2,i) .eq. id1 ) then
504 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
505 endif
506 end do
507
508 do i = 1, CnOneThree
509 if (ConeThree(1,i) .eq. id1 ) then
510 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
511 endif
512 if (ConeThree(2,i) .eq. id1 ) then
513 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
514 endif
515 end do
516
517 do i = 1, CnOneFour
518 if (ConeFour(1,i) .eq. id1 ) then
519 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
520 endif
521 if (ConeFour(2,i) .eq. id1 ) then
522 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
523 endif
524 end do
525 enddo
526
527 maxToposForAtom = maxval(nTopoPairsForAtom)
528 #ifdef IS_MPI
529 allocate(toposForAtom(nAtomsInRow, maxToposForAtom), stat=alloc_stat)
530 allocate(topoDistance(nAtomsInRow, maxToposForAtom), stat=alloc_stat)
531 #else
532 allocate(toposForAtom(nLocal, maxToposForAtom), stat=alloc_stat)
533 allocate(topoDistance(nLocal, maxToposForAtom), stat=alloc_stat)
534 #endif
535 if (alloc_stat /= 0 ) then
536 write(*,*) 'Could not allocate topoDistance array'
537 return
538 endif
539
540 #ifdef IS_MPI
541 jend = nAtomsInRow
542 #else
543 jend = nLocal
544 #endif
545 do j = 1, jend
546 nTopoPairsForAtom(j) = 0
547 #ifdef IS_MPI
548 id1 = AtomRowToGlobal(j)
549 #else
550 id1 = j
551 #endif
552 do i = 1, CnOneTwo
553 if (ConeTwo(1,i) .eq. id1 ) then
554 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
555 id2 = ConeTwo(2,i)
556 toposForAtom(j, nTopoPairsForAtom(j)) = id2
557 topoDistance(j, nTopoPairsForAtom(j)) = 1
558 endif
559 if (ConeTwo(2, i) .eq. id1 ) then
560 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
561 id2 = ConeTwo(1,i)
562 toposForAtom(j, nTopoPairsForAtom(j)) = id2
563 topoDistance(j, nTopoPairsForAtom(j)) = 1
564 endif
565 end do
566
567 do i = 1, CnOneThree
568 if (ConeThree(1,i) .eq. id1 ) then
569 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
570 id2 = ConeThree(2,i)
571 toposForAtom(j, nTopoPairsForAtom(j)) = id2
572 topoDistance(j, nTopoPairsForAtom(j)) = 2
573 endif
574 if (ConeThree(2, i) .eq. id1 ) then
575 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
576 id2 = ConeThree(1,i)
577 toposForAtom(j, nTopoPairsForAtom(j)) = id2
578 topoDistance(j, nTopoPairsForAtom(j)) = 2
579 endif
580 end do
581
582 do i = 1, CnOneFour
583 if (ConeFour(1,i) .eq. id1 ) then
584 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
585 id2 = ConeFour(2,i)
586 toposForAtom(j, nTopoPairsForAtom(j)) = id2
587 topoDistance(j, nTopoPairsForAtom(j)) = 3
588 endif
589 if (ConeFour(2, i) .eq. id1 ) then
590 nTopoPairsForAtom(j) = nTopoPairsForAtom(j) + 1
591 id2 = ConeFour(1,i)
592 toposForAtom(j, nTopoPairsForAtom(j)) = id2
593 topoDistance(j, nTopoPairsForAtom(j)) = 3
594 endif
595 end do
596 enddo
597
598 if (status == 0) simulation_setup_complete = .true.
599
600 end subroutine SimulationSetup
601
602 subroutine setBox(cHmat, cHmatInv, cBoxIsOrthorhombic)
603 real(kind=dp), dimension(3,3) :: cHmat, cHmatInv
604 integer :: cBoxIsOrthorhombic
605 integer :: smallest, status
606
607 Hmat = cHmat
608 HmatInv = cHmatInv
609 if (cBoxIsOrthorhombic .eq. 0 ) then
610 boxIsOrthorhombic = .false.
611 else
612 boxIsOrthorhombic = .true.
613 endif
614
615 call checkBox()
616 return
617 end subroutine setBox
618
619 subroutine checkBox()
620
621 integer :: i
622 real(kind=dp), dimension(3) :: hx, hy, hz, ax, ay, az, piped
623 character(len = statusMsgSize) :: errMsg
624
625 hx = Hmat(1,:)
626 hy = Hmat(2,:)
627 hz = Hmat(3,:)
628
629 ax = cross_product(hy, hz)
630 ay = cross_product(hx, hz)
631 az = cross_product(hx, hy)
632
633 ax = ax / length(ax)
634 ay = ay / length(ay)
635 az = az / length(az)
636
637 piped(1) = abs(dot_product(ax, hx))
638 piped(2) = abs(dot_product(ay, hy))
639 piped(3) = abs(dot_product(az, hz))
640
641 do i = 1, 3
642 if ((0.5_dp * piped(i)).lt.DangerRcut) then
643 write(errMsg, '(a94,f9.4,a1)') 'One of the dimensions of the Periodic ' &
644 // 'Box is smaller than ' // newline // tab // &
645 'the largest cutoff radius' // &
646 ' (rCut = ', DangerRcut, ')'
647 call handleError("checkBox", errMsg)
648
649 end if
650 enddo
651 return
652 end subroutine checkBox
653
654 function SimUsesPBC() result(doesit)
655 logical :: doesit
656 doesit = thisSim%SIM_uses_PBC
657 end function SimUsesPBC
658
659 function SimUsesAtomicVirial() result(doesit)
660 logical :: doesit
661 doesit = thisSim%SIM_uses_AtomicVirial
662 end function SimUsesAtomicVirial
663
664 function SimUsesDirectionalAtoms() result(doesit)
665 logical :: doesit
666 doesit = thisSim%SIM_uses_DirectionalAtoms
667 end function SimUsesDirectionalAtoms
668
669 function SimUsesMetallicAtoms() result(doesit)
670 logical :: doesit
671 doesit = thisSim%SIM_uses_MetallicAtoms
672 end function SimUsesMetallicAtoms
673
674 function SimRequiresSkipCorrection() result(doesit)
675 logical :: doesit
676 doesit = thisSim%SIM_requires_SkipCorrection
677 end function SimRequiresSkipCorrection
678
679 function SimRequiresSelfCorrection() result(doesit)
680 logical :: doesit
681 doesit = thisSim%SIM_requires_SelfCorrection
682 end function SimRequiresSelfCorrection
683
684 subroutine InitializeSimGlobals(thisStat)
685 integer, intent(out) :: thisStat
686 integer :: alloc_stat
687
688 thisStat = 0
689
690 call FreeSimGlobals()
691
692 allocate(excludes(2,nExcludes), stat=alloc_stat)
693 if (alloc_stat /= 0 ) then
694 thisStat = -1
695 return
696 endif
697
698 allocate(molMembershipList(nGlobal), stat=alloc_stat)
699 if (alloc_stat /= 0 ) then
700 thisStat = -1
701 return
702 endif
703
704 end subroutine InitializeSimGlobals
705
706 subroutine FreeSimGlobals()
707
708 !We free in the opposite order in which we allocate in.
709 if (allocated(topoDistance)) deallocate(topoDistance)
710 if (allocated(toposForAtom)) deallocate(toposForAtom)
711 if (allocated(nTopoPairsForAtom)) deallocate(nTopoPairsForAtom)
712 if (allocated(skipsForLocalAtom)) deallocate(skipsForLocalAtom)
713 if (allocated(nSkipsForLocalAtom)) deallocate(nSkipsForLocalAtom)
714 if (allocated(skipsForRowAtom)) deallocate(skipsForRowAtom)
715 if (allocated(nSkipsForRowAtom)) deallocate(nSkipsForRowAtom)
716 if (allocated(mfactLocal)) deallocate(mfactLocal)
717 if (allocated(mfactCol)) deallocate(mfactCol)
718 if (allocated(mfactRow)) deallocate(mfactRow)
719 if (allocated(groupListCol)) deallocate(groupListCol)
720 if (allocated(groupListRow)) deallocate(groupListRow)
721 if (allocated(groupStartCol)) deallocate(groupStartCol)
722 if (allocated(groupStartRow)) deallocate(groupStartRow)
723 if (allocated(molMembershipList)) deallocate(molMembershipList)
724 if (allocated(excludes)) deallocate(excludes)
725
726 end subroutine FreeSimGlobals
727
728 pure function getNlocal() result(n)
729 integer :: n
730 n = nLocal
731 end function getNlocal
732
733 subroutine setHmatDangerousRcutValue(dangerWillRobinson)
734 real(kind=dp), intent(in) :: dangerWillRobinson
735 DangerRcut = dangerWillRobinson
736
737 call checkBox()
738
739 return
740 end subroutine setHmatDangerousRcutValue
741
742 end module simulation
743

Properties

Name Value
svn:keywords Author Id Revision Date