Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F86310432
npair_full_nsq_omp.cpp
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Sat, Oct 5, 18:09
Size
3 KB
Mime Type
text/x-c
Expires
Mon, Oct 7, 18:09 (2 d)
Engine
blob
Format
Raw Data
Handle
21381064
Attached To
rLAMMPS lammps
npair_full_nsq_omp.cpp
View Options
/* ----------------------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#include "npair_full_nsq_omp.h"
#include "npair_omp.h"
#include "neighbor.h"
#include "neigh_list.h"
#include "atom.h"
#include "atom_vec.h"
#include "group.h"
#include "molecule.h"
#include "domain.h"
#include "my_page.h"
#include "error.h"
using
namespace
LAMMPS_NS
;
using
namespace
NeighConst
;
/* ---------------------------------------------------------------------- */
NPairFullNsqOmp
::
NPairFullNsqOmp
(
LAMMPS
*
lmp
)
:
NPair
(
lmp
)
{}
/* ----------------------------------------------------------------------
N^2 search for all neighbors
every neighbor pair appears in list of both atoms i and j
------------------------------------------------------------------------- */
void
NPairFullNsqOmp
::
build
(
NeighList
*
list
)
{
const
int
nlocal
=
(
includegroup
)
?
atom
->
nfirst
:
atom
->
nlocal
;
const
int
bitmask
=
(
includegroup
)
?
group
->
bitmask
[
includegroup
]
:
0
;
const
int
molecular
=
atom
->
molecular
;
const
int
moltemplate
=
(
molecular
==
2
)
?
1
:
0
;
NPAIR_OMP_INIT
;
#if defined(_OPENMP)
#pragma omp parallel default(none) shared(list)
#endif
NPAIR_OMP_SETUP
(
nlocal
);
int
i
,
j
,
n
,
itype
,
jtype
,
which
,
imol
,
iatom
;
tagint
tagprev
;
double
xtmp
,
ytmp
,
ztmp
,
delx
,
dely
,
delz
,
rsq
;
int
*
neighptr
;
double
**
x
=
atom
->
x
;
int
*
type
=
atom
->
type
;
int
*
mask
=
atom
->
mask
;
tagint
*
tag
=
atom
->
tag
;
tagint
*
molecule
=
atom
->
molecule
;
tagint
**
special
=
atom
->
special
;
int
**
nspecial
=
atom
->
nspecial
;
int
nall
=
atom
->
nlocal
+
atom
->
nghost
;
int
*
molindex
=
atom
->
molindex
;
int
*
molatom
=
atom
->
molatom
;
Molecule
**
onemols
=
atom
->
avec
->
onemols
;
int
*
ilist
=
list
->
ilist
;
int
*
numneigh
=
list
->
numneigh
;
int
**
firstneigh
=
list
->
firstneigh
;
// each thread has its own page allocator
MyPage
<
int
>
&
ipage
=
list
->
ipage
[
tid
];
ipage
.
reset
();
// loop over owned atoms, storing neighbors
for
(
i
=
ifrom
;
i
<
ito
;
i
++
)
{
n
=
0
;
neighptr
=
ipage
.
vget
();
itype
=
type
[
i
];
xtmp
=
x
[
i
][
0
];
ytmp
=
x
[
i
][
1
];
ztmp
=
x
[
i
][
2
];
if
(
moltemplate
)
{
imol
=
molindex
[
i
];
iatom
=
molatom
[
i
];
tagprev
=
tag
[
i
]
-
iatom
-
1
;
}
// loop over all atoms, owned and ghost
// skip i = j
for
(
j
=
0
;
j
<
nall
;
j
++
)
{
if
(
includegroup
&&
!
(
mask
[
j
]
&
bitmask
))
continue
;
if
(
i
==
j
)
continue
;
jtype
=
type
[
j
];
if
(
exclude
&&
exclusion
(
i
,
j
,
itype
,
jtype
,
mask
,
molecule
))
continue
;
delx
=
xtmp
-
x
[
j
][
0
];
dely
=
ytmp
-
x
[
j
][
1
];
delz
=
ztmp
-
x
[
j
][
2
];
rsq
=
delx
*
delx
+
dely
*
dely
+
delz
*
delz
;
if
(
rsq
<=
cutneighsq
[
itype
][
jtype
])
{
if
(
molecular
)
{
if
(
!
moltemplate
)
which
=
find_special
(
special
[
i
],
nspecial
[
i
],
tag
[
j
]);
else
if
(
imol
>=
0
)
which
=
find_special
(
onemols
[
imol
]
->
special
[
iatom
],
onemols
[
imol
]
->
nspecial
[
iatom
],
tag
[
j
]
-
tagprev
);
else
which
=
0
;
if
(
which
==
0
)
neighptr
[
n
++
]
=
j
;
else
if
(
domain
->
minimum_image_check
(
delx
,
dely
,
delz
))
neighptr
[
n
++
]
=
j
;
else
if
(
which
>
0
)
neighptr
[
n
++
]
=
j
^
(
which
<<
SBBITS
);
}
else
neighptr
[
n
++
]
=
j
;
}
}
ilist
[
i
]
=
i
;
firstneigh
[
i
]
=
neighptr
;
numneigh
[
i
]
=
n
;
ipage
.
vgot
(
n
);
if
(
ipage
.
status
())
error
->
one
(
FLERR
,
"Neighbor list overflow, boost neigh_modify one"
);
}
NPAIR_OMP_CLOSE
;
list
->
inum
=
nlocal
;
list
->
gnum
=
0
;
}
Event Timeline
Log In to Comment