Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F66373514
lm_python_containers.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
Mon, Jun 10, 03:24
Size
5 KB
Mime Type
text/x-c++
Expires
Wed, Jun 12, 03:24 (2 d)
Engine
blob
Format
Raw Data
Handle
18212768
Attached To
rLIBMULTISCALE LibMultiScale
lm_python_containers.hh
View Options
#ifndef LM_PYTHON_CONTAINERS_HH__
#define LM_PYTHON_CONTAINERS_HH__
/* -------------------------------------------------------------------------- */
#include <pybind11/eigen.h>
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
/* -------------------------------------------------------------------------- */
#include "container.hh"
#include "container_array.hh"
#include "ref_node_continuum.hh"
/* -------------------------------------------------------------------------- */
namespace
py
=
pybind11
;
__BEGIN_LIBMULTISCALE__
/* -------------------------------------------------------------------------- */
inline
void
declare_containers
(
py
::
module
&
m
)
{
py
::
class_
<
ContainerInterface
,
LMObject
,
std
::
shared_ptr
<
ContainerInterface
>>
(
m
,
"ContainerInterface"
,
py
::
multiple_inheritance
())
.
def
(
"size"
,
&
ContainerInterface
::
size
);
py
::
class_
<
ContainerArray
<
Real
>
,
ContainerInterface
,
std
::
shared_ptr
<
ContainerArray
<
Real
>>>
(
m
,
"ContainerArrayReal"
,
py
::
multiple_inheritance
())
.
def
(
"getDim"
,
&
ContainerArray
<
Real
>::
getDim
)
.
def
(
"size"
,
&
ContainerArray
<
Real
>::
size
)
.
def
(
"array"
,
&
ContainerArray
<
Real
>::
array
)
.
def
(
py
::
init
<
UInt
,
UInt
>
(),
py
::
arg
(
"nrows"
)
=
0
,
py
::
arg
(
"ncols"
)
=
1
)
.
def
(
py
::
init
<
const
LMID
&
,
UInt
,
UInt
>
(),
py
::
arg
(
"id"
),
py
::
arg
(
"nrows"
)
=
0
,
py
::
arg
(
"ncols"
)
=
1
)
.
def
(
"__getitem__"
,
[](
ContainerArray
<
Real
>
&
a
,
py
::
object
slices
)
->
decltype
(
auto
)
{
py
::
object
obj
=
py
::
cast
(
a
.
array
());
return
obj
[
slices
];
})
.
def
(
"__setitem__"
,
[](
ContainerArray
<
Real
>
&
a
,
py
::
object
slices
,
py
::
object
val
)
{
py
::
object
obj
=
py
::
cast
(
a
.
array
());
obj
[
slices
]
=
val
;
});
}
/* -------------------------------------------------------------------------- */
// references
template
<
typename
Ref
>
struct
declare_reference
{
static
inline
void
exec
(
py
::
module
&
m
,
const
char
*
class_name
)
{
py
::
class_
<
Ref
>
(
m
,
class_name
);
}
};
// node-like references
template
<
typename
Ref
,
UInt
Dim
>
inline
void
declare_elem_reference
(
py
::
module
&
m
,
const
char
*
class_name
)
{
py
::
class_
<
Ref
>
(
m
,
class_name
).
def
(
"contains"
,
[](
Ref
&
self
,
py
::
array
x
)
{
auto
X
=
py
::
cast
<
Array1D
<
Dim
>>
(
x
);
return
self
.
contains
(
X
);
});
}
// references specialization
template
<
UInt
Dim
>
struct
declare_reference
<
RefElemAkantu
<
Dim
>>
{
using
Ref
=
RefElemAkantu
<
Dim
>
;
static
inline
void
exec
(
py
::
module
&
m
,
const
char
*
class_name
)
{
if
constexpr
(
Dim
>
1
)
{
declare_elem_reference
<
Ref
,
Dim
>
(
m
,
class_name
);
}
else
py
::
class_
<
Ref
>
(
m
,
class_name
);
}
};
/* -------------------------------------------------------------------------- */
// python iterator on top of std iterator
template
<
typename
Cont
>
struct
PythonIter
{
PythonIter
(
Cont
&
c
)
:
it
(
c
.
begin
()),
end
(
c
.
end
())
{}
using
Ref
=
std
::
decay_t
<
decltype
(
*
std
::
declval
<
typename
Cont
::
iterator
>
())
>
;
auto
next
()
{
if
(
not
first
)
{
first
=
false
;
++
it
;
}
if
(
it
==
end
)
throw
py
::
stop_iteration
();
auto
&
val
=
*
it
;
return
val
;
}
bool
first
=
true
;
typename
Cont
::
iterator
it
;
typename
Cont
::
iterator
end
;
};
/* -------------------------------------------------------------------------- */
// Point-wise container
template
<
typename
container_class
>
struct
declare_container
{
static
inline
void
exec
(
py
::
module
&
m
,
const
char
*
class_name
)
{
std
::
string
iter_name
=
std
::
string
(
class_name
)
+
"Iterator"
;
py
::
class_
<
PythonIter
<
container_class
>>
(
m
,
iter_name
.
c_str
())
.
def
(
"__next__"
,
&
PythonIter
<
container_class
>::
next
,
py
::
return_value_policy
::
reference
);
std
::
string
container_name
=
std
::
string
(
class_name
)
+
"Container"
;
py
::
class_
<
container_class
,
ContainerInterface
,
std
::
shared_ptr
<
container_class
>>
(
m
,
container_name
.
c_str
())
.
def
(
"__iter__"
,
[](
container_class
&
self
)
{
return
PythonIter
<
container_class
>
(
self
);
});
std
::
string
ref_name
=
std
::
string
(
class_name
)
+
"Ref"
;
declare_reference
<
typename
PythonIter
<
container_class
>::
Ref
>::
exec
(
m
,
ref_name
.
c_str
());
}
};
/* -------------------------------------------------------------------------- */
// Mesh-wise container
template
<
typename
ContainerNodes
,
typename
ContainerElems
>
struct
declare_container
<
ContainerMesh
<
ContainerNodes
,
ContainerElems
>>
{
using
container_class
=
ContainerMesh
<
ContainerNodes
,
ContainerElems
>
;
static
inline
void
exec
(
py
::
module
&
m
,
const
char
*
class_name
)
{
std
::
string
container_name
=
std
::
string
(
class_name
)
+
"MeshContainer"
;
py
::
class_
<
container_class
,
ContainerInterface
,
std
::
shared_ptr
<
container_class
>>
(
m
,
container_name
.
c_str
())
.
def
(
"getContainerNodes"
,
&
ContainerMesh
<
ContainerNodes
,
ContainerElems
>::
getContainerNodes
)
.
def
(
"getContainerElems"
,
&
ContainerMesh
<
ContainerNodes
,
ContainerElems
>::
getContainerElems
);
using
container_nodes
=
std
::
decay_t
<
decltype
(
std
::
declval
<
container_class
>
().
getContainerNodes
())
>
;
std
::
string
class_nodes
=
std
::
string
(
class_name
)
+
"Nodes"
;
declare_container
<
container_nodes
>::
exec
(
m
,
class_nodes
.
c_str
());
using
container_elems
=
std
::
decay_t
<
decltype
(
std
::
declval
<
container_class
>
().
getContainerElems
())
>
;
std
::
string
class_elems
=
std
::
string
(
class_name
)
+
"Elems"
;
declare_container
<
container_elems
>::
exec
(
m
,
class_elems
.
c_str
());
}
};
__END_LIBMULTISCALE__
#endif
// LM_PYTHON_CONTAINERS_HH__
Event Timeline
Log In to Comment