Page MenuHomec4science

test_features.cpp
No OneTemporary

File Metadata

Created
Sat, May 11, 20:44

test_features.cpp

/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see <https://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "boussinesq.hh"
#include "eigenvalues.hh"
#include "isotropic_hardening.hh"
#include "kelvin.hh"
#include "mindlin.hh"
#include "model.hh"
#include "model_type.hh"
#include "volume_potential.hh"
#include "wrap.hh"
/* -------------------------------------------------------------------------- */
namespace tamaas {
namespace wrap {
using namespace py::literals;
template <template <model_type, UInt> class KOp, model_type type,
UInt tensor_order>
void wrapKOp(const std::string& basename, py::module& mod) {
constexpr UInt dim = model_type_traits<type>::dimension;
std::stringstream str;
str << basename << tensor_order;
py::class_<KOp<type, tensor_order>>(mod, str.str().c_str())
.def(py::init<Model*>())
.def("apply",
[](const KOp<type, tensor_order>& engine, Grid<Real, dim>& in,
Grid<Real, dim>& out) { engine.apply(in, out); });
}
template <model_type type>
void wrapIsotropicHardening(py::module& mod) {
constexpr UInt dim = model_type_traits<type>::dimension;
py::class_<IsotropicHardening<type>>(mod, "IsotropicHardening")
.def(py::init<Model*, Real, Real>())
.def_property("h", &IsotropicHardening<type>::getHardeningModulus,
&IsotropicHardening<type>::setHardeningModulus)
.def_property("sigma_0", &IsotropicHardening<type>::getYieldStress,
&IsotropicHardening<type>::setYieldStress)
.def("computeStress",
[](IsotropicHardening<type>& iso, Grid<Real, dim>& stress,
const Grid<Real, dim>& strain,
const Grid<Real, dim>& strain_increment) {
iso.template computeStress<false>(stress, strain,
strain_increment);
})
.def("computeStressUpdate",
[](IsotropicHardening<type>& iso, Grid<Real, dim>& stress,
const Grid<Real, dim>& strain,
const Grid<Real, dim>& strain_increment) {
iso.template computeStress<true>(stress, strain, strain_increment);
})
.def("computePlasticIncrement",
[](IsotropicHardening<type>& iso, Grid<Real, dim>& stress,
const Grid<Real, dim>& strain,
const Grid<Real, dim>& strain_increment) {
iso.template computePlasticIncrement<false>(stress, strain,
strain_increment);
})
.def("computePlasticIncrementUpdate",
[](IsotropicHardening<type>& iso, Grid<Real, dim>& stress,
const Grid<Real, dim>& strain,
const Grid<Real, dim>& strain_increment) {
iso.template computePlasticIncrement<true>(stress, strain,
strain_increment);
});
// .def("getPlasticStrain", &IsotropicHardening<type>::getPlasticStrain);
}
void wrapEigenvalues(py::module& mod) {
mod.def("eigenvalues",
[](model_type type, Grid<Real, 3>& eigs, const Grid<Real, 3>& field) {
eigenvalues(type, eigs, field);
},
"model_type"_a, "eigenvalues_out"_a, "field"_a);
}
/// Wrap temporary features for testing
void wrapTestFeatures(py::module& mod) {
auto test_module = mod.def_submodule("_test_features");
test_module.doc() =
"Module for testing new features.\n"
"DISCLAIMER: this API is subject to frequent and unannounced changes "
"and should **not** be relied upon!";
// wrapKOp<Kelvin, model_type::volume_2d, 2>("Kelvin_", test_module);
// wrapKOp<Kelvin, model_type::volume_2d, 3>("Kelvin_", test_module);
// wrapKOp<Kelvin, model_type::volume_2d, 4>("Kelvin_", test_module);
// wrapKOp<Mindlin, model_type::volume_2d, 3>("Mindlin_", test_module);
// wrapKOp<Mindlin, model_type::volume_2d, 4>("Mindlin_", test_module);
// wrapKOp<Boussinesq, model_type::volume_2d, 0>("Boussinesq_", test_module);
// wrapKOp<Boussinesq, model_type::volume_2d, 1>("Boussinesq_", test_module);
wrapIsotropicHardening<model_type::volume_2d>(test_module);
wrapEigenvalues(test_module);
}
} // namespace wrap
} // namespace tamaas

Event Timeline