| 1 | !! Fortran interface to C entry plug. | 
| 2 |  | 
| 3 | module simulation | 
| 4 | use definitions | 
| 5 | use neighborLists | 
| 6 | use force_globals | 
| 7 | use vector_class | 
| 8 | use atype_module | 
| 9 | use lj | 
| 10 | #ifdef IS_MPI | 
| 11 | use mpiSimulation | 
| 12 | #endif | 
| 13 |  | 
| 14 | implicit none | 
| 15 | PRIVATE | 
| 16 |  | 
| 17 | #define __FORTRAN90 | 
| 18 | #include "fSimulation.h" | 
| 19 |  | 
| 20 | type (simtype), public :: thisSim | 
| 21 |  | 
| 22 | logical, save :: simulation_setup_complete = .false. | 
| 23 |  | 
| 24 | integer, public, save :: natoms | 
| 25 | integer, public, save :: nExcludes_Global = 0 | 
| 26 | integer, public, save :: nExcludes_Local = 0 | 
| 27 | integer, allocatable, dimension(:,:), public :: excludesLocal | 
| 28 | integer, allocatable, dimension(:), public :: excludesGlobal | 
| 29 |  | 
| 30 | real(kind=dp), save :: rcut2 = 0.0_DP | 
| 31 | real(kind=dp), save :: rcut6 = 0.0_DP | 
| 32 | real(kind=dp), save :: rlist2 = 0.0_DP | 
| 33 | real(kind=dp), public, dimension(3), save :: box | 
| 34 |  | 
| 35 |  | 
| 36 | public :: SimulationSetup | 
| 37 | public :: getNlocal | 
| 38 | public :: setBox | 
| 39 | public :: setBox_3d | 
| 40 | public :: getBox | 
| 41 | public :: setRcut | 
| 42 | public :: getRcut | 
| 43 | public :: getRlist | 
| 44 | public :: getRrf | 
| 45 | public :: getRt | 
| 46 | public :: getDielect | 
| 47 | public :: SimUsesPBC | 
| 48 | public :: SimUsesLJ | 
| 49 | public :: SimUsesDipoles | 
| 50 | public :: SimUsesSticky | 
| 51 | public :: SimUsesRF | 
| 52 | public :: SimUsesGB | 
| 53 | public :: SimUsesEAM | 
| 54 | public :: SimRequiresPrepairCalc | 
| 55 | public :: SimRequiresPostpairCalc | 
| 56 | public :: SimUsesDirectionalAtoms | 
| 57 |  | 
| 58 | interface getBox | 
| 59 | module procedure getBox_3d | 
| 60 | module procedure getBox_1d | 
| 61 | end interface | 
| 62 |  | 
| 63 | interface setBox | 
| 64 | module procedure setBox_3d | 
| 65 | module procedure setBox_1d | 
| 66 | end interface | 
| 67 |  | 
| 68 | contains | 
| 69 |  | 
| 70 | subroutine SimulationSetup(setThisSim, nComponents, c_idents, & | 
| 71 | CnLocalExcludes, CexcludesLocal, CnGlobalExcludes, CexcludesGlobal, & | 
| 72 | status) | 
| 73 |  | 
| 74 | type (simtype) :: setThisSim | 
| 75 | integer, intent(inout) :: nComponents | 
| 76 | integer, dimension(nComponents),intent(inout) :: c_idents | 
| 77 |  | 
| 78 | integer :: CnLocalExcludes | 
| 79 | integer, dimension(2,CnLocalExcludes), intent(in) :: CexcludesLocal | 
| 80 | integer :: CnGlobalExcludes | 
| 81 | integer, dimension(CnGlobalExcludes), intent(in) :: CexcludesGlobal | 
| 82 | !!  Result status, success = 0, status = -1 | 
| 83 | integer, intent(out) :: status | 
| 84 | integer :: i, me, thisStat, alloc_stat, myNode | 
| 85 | #ifdef IS_MPI | 
| 86 | integer, allocatable, dimension(:) :: c_idents_Row | 
| 87 | integer, allocatable, dimension(:) :: c_idents_Col | 
| 88 | integer :: nrow | 
| 89 | integer :: ncol | 
| 90 | #endif | 
| 91 |  | 
| 92 | simulation_setup_complete = .false. | 
| 93 | status = 0 | 
| 94 |  | 
| 95 | ! copy C struct into fortran type | 
| 96 | thisSim = setThisSim | 
| 97 | natoms = nComponents | 
| 98 | rcut2 = thisSim%rcut * thisSim%rcut | 
| 99 | rcut6 = rcut2 * rcut2 * rcut2 | 
| 100 | rlist2 = thisSim%rlist * thisSim%rlist | 
| 101 | box = thisSim%box | 
| 102 |  | 
| 103 | nExcludes_Global = CnGlobalExcludes | 
| 104 | nExcludes_Local = CnLocalExcludes | 
| 105 |  | 
| 106 | call InitializeForceGlobals(natoms, thisStat) | 
| 107 | if (thisStat /= 0) then | 
| 108 | status = -1 | 
| 109 | return | 
| 110 | endif | 
| 111 |  | 
| 112 | call InitializeSimGlobals(thisStat) | 
| 113 | if (thisStat /= 0) then | 
| 114 | status = -1 | 
| 115 | return | 
| 116 | endif | 
| 117 |  | 
| 118 | #ifdef IS_MPI | 
| 119 | ! We can only set up forces if mpiSimulation has been setup. | 
| 120 | if (.not. isMPISimSet()) then | 
| 121 | write(default_error,*) "MPI is not set" | 
| 122 | status = -1 | 
| 123 | return | 
| 124 | endif | 
| 125 | nrow = getNrow(plan_row) | 
| 126 | ncol = getNcol(plan_col) | 
| 127 | mynode = getMyNode() | 
| 128 |  | 
| 129 | allocate(c_idents_Row(nrow),stat=alloc_stat) | 
| 130 | if (alloc_stat /= 0 ) then | 
| 131 | status = -1 | 
| 132 | return | 
| 133 | endif | 
| 134 |  | 
| 135 | allocate(c_idents_Col(ncol),stat=alloc_stat) | 
| 136 | if (alloc_stat /= 0 ) then | 
| 137 | status = -1 | 
| 138 | return | 
| 139 | endif | 
| 140 |  | 
| 141 | call gather(c_idents, c_idents_Row, plan_row) | 
| 142 | call gather(c_idents, c_idents_Col, plan_col) | 
| 143 |  | 
| 144 | do i = 1, nrow | 
| 145 | me = getFirstMatchingElement(atypes, "c_ident", c_idents_Row(i)) | 
| 146 | atid_Row(i) = me | 
| 147 | enddo | 
| 148 |  | 
| 149 | do i = 1, ncol | 
| 150 | me = getFirstMatchingElement(atypes, "c_ident", c_idents_Col(i)) | 
| 151 | atid_Col(i) = me | 
| 152 | enddo | 
| 153 |  | 
| 154 | !! free temporary ident arrays | 
| 155 | if (allocated(c_idents_Col)) then | 
| 156 | deallocate(c_idents_Col) | 
| 157 | end if | 
| 158 | if (allocated(c_idents_Row)) then | 
| 159 | deallocate(c_idents_Row) | 
| 160 | endif | 
| 161 |  | 
| 162 | #else | 
| 163 | do i = 1, nComponents | 
| 164 |  | 
| 165 | me = getFirstMatchingElement(atypes, "c_ident", c_idents(i)) | 
| 166 | atid(i) = me | 
| 167 |  | 
| 168 | enddo | 
| 169 | #endif | 
| 170 |  | 
| 171 | !! Create neighbor lists | 
| 172 | call expandNeighborList(nComponents, thisStat) | 
| 173 | if (thisStat /= 0) then | 
| 174 | status = -1 | 
| 175 | return | 
| 176 | endif | 
| 177 |  | 
| 178 | do i = 1, nExcludes_Local | 
| 179 | excludesLocal(1,i) = CexcludesLocal(1,i) | 
| 180 | excludesLocal(2,i) = CexcludesLocal(2,i) | 
| 181 | enddo | 
| 182 |  | 
| 183 | do i = 1, nExcludes_Global | 
| 184 | excludesGlobal(i) = CexcludesGlobal(i) | 
| 185 | enddo | 
| 186 |  | 
| 187 | if (status == 0) simulation_setup_complete = .true. | 
| 188 |  | 
| 189 | end subroutine SimulationSetup | 
| 190 |  | 
| 191 | subroutine setBox_3d(new_box_size) | 
| 192 | real(kind=dp), dimension(3) :: new_box_size | 
| 193 | integer :: smallest, status, i | 
| 194 |  | 
| 195 | thisSim%box = new_box_size | 
| 196 | box = thisSim%box | 
| 197 |  | 
| 198 | smallest = 1 | 
| 199 | do i = 2, 3 | 
| 200 | if (new_box_size(i) .lt. new_box_size(smallest)) smallest = i | 
| 201 | end do | 
| 202 | if (thisSim%rcut .gt. 0.5_dp * new_box_size(smallest)) & | 
| 203 | call setRcut(0.5_dp * new_box_size(smallest), status) | 
| 204 | return | 
| 205 | end subroutine setBox_3d | 
| 206 |  | 
| 207 | subroutine setBox_1d(dim, new_box_size) | 
| 208 | integer :: dim, status | 
| 209 | real(kind=dp) :: new_box_size | 
| 210 | thisSim%box(dim) = new_box_size | 
| 211 | box(dim) = thisSim%box(dim) | 
| 212 | if (thisSim%rcut .gt. 0.5_dp * new_box_size) & | 
| 213 | call setRcut(0.5_dp * new_box_size, status) | 
| 214 | end subroutine setBox_1d | 
| 215 |  | 
| 216 | subroutine setRcut(new_rcut, status) | 
| 217 | real(kind = dp) :: new_rcut | 
| 218 | integer :: myStatus, status | 
| 219 | thisSim%rcut = new_rcut | 
| 220 | rcut2 = thisSim%rcut * thisSim%rcut | 
| 221 | rcut6 = rcut2 * rcut2 * rcut2 | 
| 222 | myStatus = 0 | 
| 223 | call LJ_new_rcut(new_rcut, myStatus) | 
| 224 | if (myStatus .ne. 0) then | 
| 225 | write(default_error, *) 'LJ module refused our rcut!' | 
| 226 | status = -1 | 
| 227 | return | 
| 228 | endif | 
| 229 | status = 0 | 
| 230 | return | 
| 231 | end subroutine setRcut | 
| 232 |  | 
| 233 | function getBox_3d() result(thisBox) | 
| 234 | real( kind = dp ), dimension(3) :: thisBox | 
| 235 | thisBox = thisSim%box | 
| 236 | end function getBox_3d | 
| 237 |  | 
| 238 | function getBox_1d(dim) result(thisBox) | 
| 239 | integer, intent(in) :: dim | 
| 240 | real( kind = dp ) :: thisBox | 
| 241 |  | 
| 242 | thisBox = thisSim%box(dim) | 
| 243 | end function getBox_1d | 
| 244 |  | 
| 245 | subroutine getRcut(thisrcut,rc2,rc6,status) | 
| 246 | real( kind = dp ), intent(out) :: thisrcut | 
| 247 | real( kind = dp ), intent(out), optional :: rc2 | 
| 248 | real( kind = dp ), intent(out), optional :: rc6 | 
| 249 | integer, optional :: status | 
| 250 |  | 
| 251 | if (present(status)) status = 0 | 
| 252 |  | 
| 253 | if (.not.simulation_setup_complete ) then | 
| 254 | if (present(status)) status = -1 | 
| 255 | return | 
| 256 | end if | 
| 257 |  | 
| 258 | thisrcut = thisSim%rcut | 
| 259 | if(present(rc2)) rc2 = rcut2 | 
| 260 | if(present(rc6)) rc6 = rcut6 | 
| 261 | end subroutine getRcut | 
| 262 |  | 
| 263 | subroutine getRlist(thisrlist,rl2,status) | 
| 264 | real( kind = dp ), intent(out) :: thisrlist | 
| 265 | real( kind = dp ), intent(out), optional :: rl2 | 
| 266 |  | 
| 267 | integer, optional :: status | 
| 268 |  | 
| 269 | if (present(status)) status = 0 | 
| 270 |  | 
| 271 | if (.not.simulation_setup_complete ) then | 
| 272 | if (present(status)) status = -1 | 
| 273 | return | 
| 274 | end if | 
| 275 |  | 
| 276 | thisrlist = thisSim%rlist | 
| 277 | if(present(rl2)) rl2 = rlist2 | 
| 278 | end subroutine getRlist | 
| 279 |  | 
| 280 | function getRrf() result(rrf) | 
| 281 | real( kind = dp ) :: rrf | 
| 282 | rrf = thisSim%rrf | 
| 283 | end function getRrf | 
| 284 |  | 
| 285 | function getRt() result(rt) | 
| 286 | real( kind = dp ) :: rt | 
| 287 | rt = thisSim%rt | 
| 288 | end function getRt | 
| 289 |  | 
| 290 | function getDielect() result(dielect) | 
| 291 | real( kind = dp ) :: dielect | 
| 292 | dielect = thisSim%dielect | 
| 293 | end function getDielect | 
| 294 |  | 
| 295 | function SimUsesPBC() result(doesit) | 
| 296 | logical :: doesit | 
| 297 | doesit = thisSim%SIM_uses_PBC | 
| 298 | end function SimUsesPBC | 
| 299 |  | 
| 300 | function SimUsesLJ() result(doesit) | 
| 301 | logical :: doesit | 
| 302 | doesit = thisSim%SIM_uses_LJ | 
| 303 | end function SimUsesLJ | 
| 304 |  | 
| 305 | function SimUsesSticky() result(doesit) | 
| 306 | logical :: doesit | 
| 307 | doesit = thisSim%SIM_uses_sticky | 
| 308 | end function SimUsesSticky | 
| 309 |  | 
| 310 | function SimUsesDipoles() result(doesit) | 
| 311 | logical :: doesit | 
| 312 | doesit = thisSim%SIM_uses_dipoles | 
| 313 | end function SimUsesDipoles | 
| 314 |  | 
| 315 | function SimUsesRF() result(doesit) | 
| 316 | logical :: doesit | 
| 317 | doesit = thisSim%SIM_uses_RF | 
| 318 | end function SimUsesRF | 
| 319 |  | 
| 320 | function SimUsesGB() result(doesit) | 
| 321 | logical :: doesit | 
| 322 | doesit = thisSim%SIM_uses_GB | 
| 323 | end function SimUsesGB | 
| 324 |  | 
| 325 | function SimUsesEAM() result(doesit) | 
| 326 | logical :: doesit | 
| 327 | doesit = thisSim%SIM_uses_EAM | 
| 328 | end function SimUsesEAM | 
| 329 |  | 
| 330 | function SimUsesDirectionalAtoms() result(doesit) | 
| 331 | logical :: doesit | 
| 332 | doesit = thisSim%SIM_uses_dipoles .or. thisSim%SIM_uses_sticky .or. & | 
| 333 | thisSim%SIM_uses_GB .or. thisSim%SIM_uses_RF | 
| 334 | end function SimUsesDirectionalAtoms | 
| 335 |  | 
| 336 | function SimRequiresPrepairCalc() result(doesit) | 
| 337 | logical :: doesit | 
| 338 | doesit = thisSim%SIM_uses_EAM | 
| 339 | end function SimRequiresPrepairCalc | 
| 340 |  | 
| 341 | function SimRequiresPostpairCalc() result(doesit) | 
| 342 | logical :: doesit | 
| 343 | doesit = thisSim%SIM_uses_RF | 
| 344 | end function SimRequiresPostpairCalc | 
| 345 |  | 
| 346 | subroutine InitializeSimGlobals(thisStat) | 
| 347 | integer, intent(out) :: thisStat | 
| 348 | integer :: alloc_stat | 
| 349 |  | 
| 350 | thisStat = 0 | 
| 351 |  | 
| 352 | call FreeSimGlobals() | 
| 353 |  | 
| 354 | allocate(excludesLocal(2,nExcludes_Local), stat=alloc_stat) | 
| 355 | if (alloc_stat /= 0 ) then | 
| 356 | thisStat = -1 | 
| 357 | return | 
| 358 | endif | 
| 359 |  | 
| 360 | allocate(excludesGlobal(nExcludes_Global), stat=alloc_stat) | 
| 361 | if (alloc_stat /= 0 ) then | 
| 362 | thisStat = -1 | 
| 363 | return | 
| 364 | endif | 
| 365 |  | 
| 366 | end subroutine InitializeSimGlobals | 
| 367 |  | 
| 368 | subroutine FreeSimGlobals() | 
| 369 |  | 
| 370 | !We free in the opposite order in which we allocate in. | 
| 371 |  | 
| 372 | if (allocated(excludesGlobal)) deallocate(excludesGlobal) | 
| 373 | if (allocated(excludesLocal)) deallocate(excludesLocal) | 
| 374 |  | 
| 375 | end subroutine FreeSimGlobals | 
| 376 |  | 
| 377 | pure function getNlocal() result(nlocal) | 
| 378 | integer :: nlocal | 
| 379 | nlocal = natoms | 
| 380 | end function getNlocal | 
| 381 |  | 
| 382 |  | 
| 383 | end module simulation |