Page MenuHomec4science

isotropic_hardening.hh
No OneTemporary

File Metadata

Created
Tue, May 7, 19:24

isotropic_hardening.hh

/*
* SPDX-License-Indentifier: AGPL-3.0-or-later
*
* Copyright (©) 2016-2022 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/>.
*
*/
/* -------------------------------------------------------------------------- */
#ifndef ISOTROPIC_HARDENING_HH
#define ISOTROPIC_HARDENING_HH
/* -------------------------------------------------------------------------- */
#include "grid.hh"
#include "influence.hh"
#include "internal.hh"
#include "material.hh"
#include "model.hh"
#include "model_type.hh"
/* -------------------------------------------------------------------------- */
namespace tamaas {
/* -------------------------------------------------------------------------- */
template <model_type type>
class IsotropicHardening : Material<type> {
using parent = Material<type>;
using trait = typename parent::trait;
static constexpr UInt dim = parent::trait::dimension;
using Field = typename parent::Field;
public:
/// Constructor
IsotropicHardening(Model* model, Real sigma_0, Real h);
/// Compute plastic strain increment with radial return algorithm
void
computeInelasticDeformationIncrement(Field& increment, const Field& strain,
const Field& strain_increment) override;
/// Compute stress
void computeStress(Field& stress, const Field& strain,
const Field& strain_increment) override;
/// Update internal variables
void update() override;
void applyTangentIncrement(Field& output, const Field& input,
const Field& strain,
const Field& strain_increment) const;
/// Linear hardening function
static __device__ __host__ Real hardening(Real p, Real h, Real sigma_0) {
return sigma_0 + h * p;
}
Real getHardeningModulus() const { return h; }
Real getYieldStress() const { return sigma_0; }
const GridBase<Real>& getPlasticStrain() const { return *plastic_strain; }
GridBase<Real>& getPlasticStrain() { return *plastic_strain; }
void setHardeningModulus(Real h_) {
if (h_ < 0)
TAMAAS_EXCEPTION("Hardening modulus should be positive");
h = h_;
}
void setYieldStress(Real sigma_0_) {
if (sigma_0_ < 0)
TAMAAS_EXCEPTION("Yield stress should be positive");
sigma_0 = sigma_0_;
}
protected:
Real sigma_0; /// < initial yield stress
Real h; /// < hardening modulus
Internal<Real, 3> plastic_strain, cumulated_plastic_strain;
};
} // namespace tamaas
/* -------------------------------------------------------------------------- */
#endif // ISOTROPIC_HARDENING_HH

Event Timeline