Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F85228911
container_neighbor_atoms.hh
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
Fri, Sep 27, 16:05
Size
9 KB
Mime Type
text/x-c++
Expires
Sun, Sep 29, 16:05 (2 d)
Engine
blob
Format
Raw Data
Handle
21145569
Attached To
rLIBMULTISCALE LibMultiScale
container_neighbor_atoms.hh
View Options
/**
* @file container_neighbor_atoms.hh
*
* @author Guillaume Anciaux <guillaume.anciaux@epfl.ch>
* @author Srinivasa Babu Ramisetti <srinivasa.ramisetti@epfl.ch>
*
* @date Wed Jan 23 15:53:40 2013
*
* @brief This class virtualize a container of neighbors for a collection of atoms
*
* @section LICENSE
*
* Copyright (©) 2010-2011 EPFL (Ecole Polytechnique Fédérale de Lausanne)
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* LibMultiScale is free software: you can redistribute it and/or modify it under the
* terms of the GNU Lesser General Public License as published by the Free
* Software Foundation, either version 3 of the License, or (at your option) any
* later version.
*
* LibMultiScale is distributed in the hope that it will be useful, but WITHOUT ANY
* WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
* A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
* details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with LibMultiScale. If not, see <http://www.gnu.org/licenses/>.
*
*/
#ifndef __LIBMULTISCALE_CONTAINER_NEIGHBOR_ATOMS_HH__
#define __LIBMULTISCALE_CONTAINER_NEIGHBOR_ATOMS_HH__
/* -------------------------------------------------------------------------- */
#include "spatial_grid_libmultiscale.hh"
#include "container_neighbor_atoms_interface.hh"
#include "container.hh"
/* -------------------------------------------------------------------------- */
__BEGIN_LIBMULTISCALE__
template
<
typename
T
,
UInt
Dim
>
class
ContainerNeighborAtoms
:
public
SpatialGridLibMultiScale
<
T
,
Dim
>
,
public
ContainerNeighborAtomsInterface
,
public
Container_base
<
T
>
{
/* ------------------------------------------------------------------------ */
/* Typedefs */
/* ------------------------------------------------------------------------ */
public
:
class
iterator
;
/* ------------------------------------------------------------------------ */
/* Constructors/Destructors */
/* ------------------------------------------------------------------------ */
public
:
ContainerNeighborAtoms
(
Cube
&
c
,
Real
p
);
ContainerNeighborAtoms
(
Cube
&
c
,
Real
x
,
Real
y
,
Real
z
);
ContainerNeighborAtoms
(
Cube
&
c
,
UInt
*
p
);
virtual
~
ContainerNeighborAtoms
();
/* ------------------------------------------------------------------------ */
/* Methods */
/* ------------------------------------------------------------------------ */
public
:
template
<
typename
Cont
>
void
fillGrid
(
Cont
&
c
);
template
<
typename
Cont
>
void
fillGridWithInitialPositions
(
Cont
&
c
);
iterator
getIterator
();
void
setPBCFlag
(
bool
*
pbc
);
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
private
:
UInt
pbc_flag
[
3
];
};
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
ContainerNeighborAtoms
<
T
,
Dim
>::
ContainerNeighborAtoms
(
Cube
&
c
,
Real
p
)
:
SpatialGridLibMultiScale
<
T
,
Dim
>
(
c
,
p
){
for
(
UInt
i
=
0
;
i
<
Dim
;
++
i
)
{
pbc_flag
[
i
]
=
0
;
}
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
ContainerNeighborAtoms
<
T
,
Dim
>::
ContainerNeighborAtoms
(
Cube
&
c
,
Real
x
,
Real
y
,
Real
z
)
:
SpatialGridLibMultiScale
<
T
,
Dim
>
(
c
,
x
,
y
,
z
){
for
(
UInt
i
=
0
;
i
<
Dim
;
++
i
)
{
pbc_flag
[
i
]
=
0
;
}
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
ContainerNeighborAtoms
<
T
,
Dim
>::
ContainerNeighborAtoms
(
Cube
&
c
,
UInt
*
p
)
:
SpatialGridLibMultiScale
<
T
,
Dim
>
(
c
,
p
){
for
(
UInt
i
=
0
;
i
<
Dim
;
++
i
)
{
pbc_flag
[
i
]
=
0
;
}
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
ContainerNeighborAtoms
<
T
,
Dim
>::~
ContainerNeighborAtoms
(){
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
typename
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
ContainerNeighborAtoms
<
T
,
Dim
>::
getIterator
(){
typename
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
it
(
*
this
);
return
it
;
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
template
<
typename
Cont
>
void
ContainerNeighborAtoms
<
T
,
Dim
>::
fillGrid
(
Cont
&
c
){
typename
Cont
::
iterator
it
=
c
.
getIterator
();
it
.
includeGhosts
(
true
);
typename
Cont
::
Ref
at
;
for
(
at
=
it
.
getFirst
()
;
!
it
.
end
()
;
at
=
it
.
getNext
())
{
Real
x
[
Dim
];
at
.
getPositions
(
x
);
this
->
addElement
(
at
,
x
);
}
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
template
<
typename
Cont
>
void
ContainerNeighborAtoms
<
T
,
Dim
>::
fillGridWithInitialPositions
(
Cont
&
c
){
typename
Cont
::
iterator
it
=
c
.
getIterator
();
typename
Cont
::
Ref
at
;
for
(
at
=
it
.
getFirst
()
;
!
it
.
end
()
;
at
=
it
.
getNext
())
{
Real
x
[
Dim
];
at
.
getPositions0
(
x
);
this
->
addElement
(
at
,
x
);
}
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
void
ContainerNeighborAtoms
<
T
,
Dim
>::
setPBCFlag
(
bool
*
pbc
){
for
(
UInt
i
=
0
;
i
<
Dim
;
++
i
)
{
pbc_flag
[
i
]
=
pbc
[
i
];
}
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
class
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
:
public
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator_base
{
/* ------------------------------------------------------------------------ */
/* Constructors/Destructors */
/* ------------------------------------------------------------------------ */
public
:
iterator
(
ContainerNeighborAtoms
<
T
,
Dim
>
&
c
);
~
iterator
();
/* ------------------------------------------------------------------------ */
/* Methods */
/* ------------------------------------------------------------------------ */
template
<
typename
R
>
T
&
getFirst
(
R
&
dof
);
template
<
typename
R
>
T
&
getFirst0
(
R
&
dof
);
T
&
getFirst
(
Real
*
x
);
T
&
getNext
();
void
findNextNeighborBlock
();
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
protected
:
UInt
central_block
;
std
::
vector
<
T
>
*
current_block
;
std
::
vector
<
UInt
>
neighborhood_indexes
;
UInt
current_el
;
T
NULL_ELEM
;
};
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
::
iterator
(
ContainerNeighborAtoms
<
T
,
Dim
>
&
c
)
:
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator_base
(
c
){
if
(
Dim
==
1
)
neighborhood_indexes
.
reserve
(
3
);
if
(
Dim
==
2
)
neighborhood_indexes
.
reserve
(
9
);
if
(
Dim
==
3
)
neighborhood_indexes
.
reserve
(
27
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
::~
iterator
()
{
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
T
&
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
::
getFirst
(
Real
*
x
)
{
ContainerNeighborAtoms
<
T
,
Dim
>
&
c
=
static_cast
<
ContainerNeighborAtoms
<
T
,
Dim
>
&>
(
this
->
my_cont
);
c
.
getNeighborhood
(
x
,
neighborhood_indexes
,
c
.
pbc_flag
);
findNextNeighborBlock
();
current_el
=
0
;
if
(
!
current_block
)
{
this
->
end_flag
=
true
;
return
NULL_ELEM
;
}
this
->
end_flag
=
false
;
return
current_block
->
at
(
0
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
template
<
typename
R
>
T
&
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
::
getFirst
(
R
&
dof
)
{
Real
x
[
3
];
dof
.
getPositions
(
x
);
return
getFirst
(
x
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
template
<
typename
R
>
T
&
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
::
getFirst0
(
R
&
dof
)
{
Real
x
[
3
];
dof
.
getPositions0
(
x
);
return
getFirst
(
x
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
T
&
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
::
getNext
()
{
if
(
current_el
+
1
>=
current_block
->
size
())
{
findNextNeighborBlock
();
if
(
!
current_block
)
{
this
->
end_flag
=
true
;
return
NULL_ELEM
;
}
else
current_el
=
0
;
}
else
++
current_el
;
return
current_block
->
at
(
current_el
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
,
UInt
Dim
>
void
ContainerNeighborAtoms
<
T
,
Dim
>::
iterator
::
findNextNeighborBlock
(){
ContainerNeighborAtoms
<
T
,
Dim
>
&
c
=
static_cast
<
ContainerNeighborAtoms
<
T
,
Dim
>
&>
(
this
->
my_cont
);
if
(
neighborhood_indexes
.
size
()
==
0
){
current_block
=
NULL
;
return
;
}
do
{
UInt
neigh
=
neighborhood_indexes
.
back
();
LM_ASSERT
(
neigh
!=
UINT_MAX
,
"invalid neighbor block number"
);
neighborhood_indexes
.
pop_back
();
current_block
=
(
c
.
getBlocks
())[
neigh
];
}
while
(
!
current_block
&&
neighborhood_indexes
.
size
()
>
0
);
}
/* -------------------------------------------------------------------------- */
__END_LIBMULTISCALE__
#endif
/* __LIBMULTISCALE_CONTAINER_NEIGHBOR_ATOMS_HH__ */
Event Timeline
Log In to Comment