Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F92232647
py_aka_solid_mechanics_model.cc
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, Nov 18, 14:24
Size
5 KB
Mime Type
text/x-c
Expires
Wed, Nov 20, 14:24 (2 d)
Engine
blob
Format
Raw Data
Handle
22378362
Attached To
rAKA akantu
py_aka_solid_mechanics_model.cc
View Options
/* -------------------------------------------------------------------------- */
#include "py_aka_array.hh"
/* -------------------------------------------------------------------------- */
#include <non_linear_solver.hh>
#include <solid_mechanics_model.hh>
/* -------------------------------------------------------------------------- */
#include <pybind11/operators.h>
#include <pybind11/pybind11.h>
#include <pybind11/stl.h>
/* -------------------------------------------------------------------------- */
namespace py = pybind11;
/* -------------------------------------------------------------------------- */
namespace akantu {
/* -------------------------------------------------------------------------- */
#define def_deprecated(func_name, mesg) \
def(func_name, [](py::args, py::kwargs) { AKANTU_ERROR(mesg); })
#define def_function_nocopy(func_name) \
def(#func_name, \
[](SolidMechanicsModel & self) -> decltype(auto) { \
return self.func_name(); \
}, \
py::return_value_policy::reference)
#define def_function(func_name) \
def(#func_name, [](SolidMechanicsModel & self) -> decltype(auto) { \
return self.func_name(); \
})
/* -------------------------------------------------------------------------- */
[[gnu::visibility("default")]] void
register_solid_mechanics_model(py::module & mod) {
py::class_<SolidMechanicsModelOptions>(mod, "SolidMechanicsModelOptions")
.def(py::init<AnalysisMethod>(),
py::arg("analysis_method") = _explicit_lumped_mass);
py::class_<SolidMechanicsModel, Model, ModelSolver>(mod,
"SolidMechanicsModel")
.def(py::init<Mesh &, UInt, const ID &, const MemoryID &,
const ModelType>(),
py::arg("mesh"), py::arg("spatial_dimension") = _all_dimensions,
py::arg("id") = "solid_mechanics_model", py::arg("memory_id") = 0,
py::arg("model_type") = ModelType::_solid_mechanics_model)
.def("initFull",
[](SolidMechanicsModel & self,
const SolidMechanicsModelOptions & options) {
self.initFull(options);
},
py::arg("_analysis_method") = SolidMechanicsModelOptions())
.def("initFull",
[](SolidMechanicsModel & self,
const AnalysisMethod & _analysis_method) {
self.initFull(SolidMechanicsModelOptions(_analysis_method));
},
py::arg("_analysis_method"))
.def_deprecated("applyDirichletBC", "Deprecated: use applyBC")
.def("applyBC",
[](SolidMechanicsModel & self,
BC::Dirichlet::DirichletFunctor & func,
const std::string & element_group) {
self.applyBC(func, element_group);
})
.def("applyBC",
[](SolidMechanicsModel & self, BC::Neumann::NeumannFunctor & func,
const std::string & element_group) {
self.applyBC(func, element_group);
})
.def("setTimeStep", &SolidMechanicsModel::setTimeStep,
py::arg("time_step"), py::arg("solver_id") = "")
.def("getEnergy",
py::overload_cast<const std::string &>(
&SolidMechanicsModel::getEnergy),
py::arg("energy_id"))
.def_function(assembleStiffnessMatrix)
.def_function(assembleInternalForces)
.def_function(assembleMass)
.def_function(assembleMassLumped)
.def_function(getStableTimeStep)
.def_function_nocopy(getExternalForce)
.def_function_nocopy(getDisplacement)
.def_function_nocopy(getPreviousDisplacement)
.def_function_nocopy(getIncrement)
.def_function_nocopy(getInternalForce)
.def_function_nocopy(getMass)
.def_function_nocopy(getVelocity)
.def_function_nocopy(getAcceleration)
.def_function_nocopy(getInternalForce)
.def_function_nocopy(getBlockedDOFs)
.def_function_nocopy(getMesh)
.def("dump", py::overload_cast<>(&SolidMechanicsModel::dump))
.def("dump",
py::overload_cast<const std::string &>(&SolidMechanicsModel::dump))
.def("dump", py::overload_cast<const std::string &, UInt>(
&SolidMechanicsModel::dump))
.def("dump", py::overload_cast<const std::string &, Real, UInt>(
&SolidMechanicsModel::dump))
.def("getMaterial",
py::overload_cast<UInt>(&SolidMechanicsModel::getMaterial),
py::return_value_policy::reference)
.def("getMaterial",
py::overload_cast<const std::string &>(
&SolidMechanicsModel::getMaterial),
py::return_value_policy::reference)
.def("getMaterialIndex", &SolidMechanicsModel::getMaterialIndex);
}
} // namespace akantu
Event Timeline
Log In to Comment