Page MenuHomec4science

material_anisotropic.hh
No OneTemporary

File Metadata

Created
Sun, May 5, 05:01

material_anisotropic.hh

/**
* @file material_anisotropic.hh
*
* @author Ali Falsafi<ali.falsafi@epfl.ch>
*
* @date 9 Jul 2018
*
* @brief Base class for materials (constitutive models)
*
* Copyright © 2017 Till Junge
*
* µSpectre is free software; you can redistribute it and/or
* modify it under the terms of the GNU General Public License as
* published by the Free Software Foundation, either version 3, or (at
* your option) any later version.
*
* µSpectre 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
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with GNU Emacs; see the file COPYING. If not, write to the
* Free Software Foundation, Inc., 59 Temple Place - Suite 330,
* Boston, MA 02111-1307, USA.
*/
#ifndef SRC_MATERIALS_MATERIAL_ANISOTROPIC_HH_
#define SRC_MATERIALS_MATERIAL_ANISOTROPIC_HH_
#include "materials/material_base.hh"
#include "materials/material_muSpectre_base.hh"
#include "materials/stress_transformations_PK2.hh"
#include "common/common.hh"
#include "common/T4_map_proxy.hh"
namespace muSpectre {
template <Dim_t DimS, Dim_t DimM>
class MaterialAnisotropic;
// traits for anisotropic material
template <Dim_t DimS, Dim_t DimM>
struct MaterialMuSpectre_traits<MaterialAnisotropic<DimS, DimM>> {
using Parent = MaterialMuSpectre_traits<void>; //!< base for elasticity
//! global field collection
using GFieldCollection_t =
typename MaterialBase<DimS, DimM>::GFieldCollection_t;
//! expected map type for strain fields
using StrainMap_t =
MatrixFieldMap<GFieldCollection_t, Real, DimM, DimM, true>;
//! expected map type for stress fields
using StressMap_t = MatrixFieldMap<GFieldCollection_t, Real, DimM, DimM>;
//! expected map type for tangent stiffness fields
using TangentMap_t = T4MatrixFieldMap<GFieldCollection_t, Real, DimM>;
//! declare what type of strain measure your law takes as input
constexpr static auto strain_measure{StrainMeasure::GreenLagrange};
//! declare what type of stress measure your law yields as output
constexpr static auto stress_measure{StressMeasure::PK2};
//! anisotropicity without internal variables
using InternalVariables = std::tuple<>;
};
/**
* Material implementation for anisotropic constitutive law
*/
template <Dim_t DimS, Dim_t DimM = DimS>
class MaterialAnisotropic
: public MaterialMuSpectre<MaterialAnisotropic<DimS, DimM>, DimS, DimM> {
//! base class
using Parent = MaterialMuSpectre<MaterialAnisotropic, DimS, DimM>;
using Stiffness_t = T4Mat<Real, DimM>;
//! traits of this material
using traits = MaterialMuSpectre_traits<MaterialAnisotropic>;
//! this law does not have any internal variables
using InternalVariables = typename traits::InternalVariables;
public:
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
using parent =
MaterialMuSpectre<MaterialAnisotropic<DimS, DimM>, DimS, DimM>;
//! global field collection
using GFieldCollection_t =
typename MaterialBase<DimS, DimM>::GFieldCollection_t;
//! expected map type for tangent stiffness fields
using Tangent_t = T4MatrixFieldMap<GFieldCollection_t, Real, DimM>;
//! Default constructor
MaterialAnisotropic() = delete;
// constructor
// a std::vector is utilized as the input of the constructor to
// enable us to check its length so to prevent user mistake
MaterialAnisotropic(std::string name, std::vector<Real> input_c);
//! Copy constructor
MaterialAnisotropic(const MaterialAnisotropic & other) = delete;
//! Move constructor
MaterialAnisotropic(MaterialAnisotropic && other) = delete;
//! Destructor
virtual ~MaterialAnisotropic() = default;
template <class s_t>
inline decltype(auto) evaluate_stress(s_t && E);
/**
* evaluates both second Piola-Kirchhoff stress and stiffness given
* the Green-Lagrange strain (or Cauchy stress and stiffness if
* called with a small strain tensor) and the local stiffness tensor.
*/
template <class s_t>
inline decltype(auto) evaluate_stress_tangent(s_t && E);
/**
* return the empty internals tuple
*/
InternalVariables & get_internals() { return this->internal_variables; }
// takes the elements of the C and makes it:
static auto c_maker(std::vector<Real> input) -> Stiffness_t;
auto get_C() { return this->C; }
protected:
Stiffness_t C; //!< stiffness tensor
//! empty tuple
InternalVariables internal_variables{};
};
/* ---------------------------------------------------------------------- */
template <Dim_t DimS, Dim_t DimM>
template <class s_t>
decltype(auto) MaterialAnisotropic<DimS, DimM>::evaluate_stress(s_t && E) {
return Matrices::tensmult(this->C, E);
}
/* ---------------------------------------------------------------------- */
template <Dim_t DimS, Dim_t DimM>
template <class s_t>
decltype(auto)
MaterialAnisotropic<DimS, DimM>::evaluate_stress_tangent(s_t && E) {
return std::make_tuple(evaluate_stress(E), this->C);
}
} // namespace muSpectre
#endif // SRC_MATERIALS_MATERIAL_ANISOTROPIC_HH_

Event Timeline