Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F90984938
local_material_damage.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
Wed, Nov 6, 16:08
Size
3 KB
Mime Type
text/x-c
Expires
Fri, Nov 8, 16:08 (2 d)
Engine
blob
Format
Raw Data
Handle
22088057
Attached To
rAKA akantu
local_material_damage.cc
View Options
/**
* @file local_material_damage.cc
*
* @author Guillaume Anciaux <guillaume.anciaux@epfl.ch>
* @author Marion Estelle Chambart <marion.chambart@epfl.ch>
* @author Nicolas Richart <nicolas.richart@epfl.ch>
*
* @date creation: Mon Jan 18 2016
*
* @brief Specialization of the material class for the damage material
*
*
* Copyright (©) 2015 EPFL (Ecole Polytechnique Fédérale de Lausanne) Laboratory
* (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* Akantu is free software: you can redistribute it and/or modify it under the
* terms of the GNU Lesser General Public License as published by the Free
* Software Foundation, either version 3 of the License, or (at your option) any
* later version.
*
* Akantu 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 Lesser General Public License for more
* details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with Akantu. If not, see <http://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "local_material_damage.hh"
#include "solid_mechanics_model.hh"
namespace akantu {
/* -------------------------------------------------------------------------- */
LocalMaterialDamage::LocalMaterialDamage(SolidMechanicsModel & model,
const ID & id)
: Material(model, id), damage("damage", *this) {
AKANTU_DEBUG_IN();
this->registerParam("E", E, 0., _pat_parsable, "Young's modulus");
this->registerParam("nu", nu, 0.5, _pat_parsable, "Poisson's ratio");
this->registerParam("lambda", lambda, _pat_readable,
"First Lamé coefficient");
this->registerParam("mu", mu, _pat_readable, "Second Lamé coefficient");
this->registerParam("kapa", kpa, _pat_readable, "Bulk coefficient");
this->registerParam("Yd", Yd, 50., _pat_parsmod);
this->registerParam("Sd", Sd, 5000., _pat_parsmod);
damage.initialize(1);
AKANTU_DEBUG_OUT();
}
/* -------------------------------------------------------------------------- */
void LocalMaterialDamage::initMaterial() {
AKANTU_DEBUG_IN();
Material::initMaterial();
lambda = nu * E / ((1 + nu) * (1 - 2 * nu));
mu = E / (2 * (1 + nu));
kpa = lambda + 2. / 3. * mu;
AKANTU_DEBUG_OUT();
}
/* -------------------------------------------------------------------------- */
void LocalMaterialDamage::computeStress(ElementType el_type,
GhostType ghost_type) {
AKANTU_DEBUG_IN();
Real * dam = damage(el_type, ghost_type).storage();
MATERIAL_STRESS_QUADRATURE_POINT_LOOP_BEGIN(el_type, ghost_type);
computeStressOnQuad(grad_u, sigma, *dam);
++dam;
MATERIAL_STRESS_QUADRATURE_POINT_LOOP_END;
AKANTU_DEBUG_OUT();
}
/* -------------------------------------------------------------------------- */
void LocalMaterialDamage::computePotentialEnergy(ElementType el_type) {
AKANTU_DEBUG_IN();
Real * epot = potential_energy(el_type).storage();
MATERIAL_STRESS_QUADRATURE_POINT_LOOP_BEGIN(el_type, _not_ghost);
computePotentialEnergyOnQuad(grad_u, sigma, *epot);
epot++;
MATERIAL_STRESS_QUADRATURE_POINT_LOOP_END;
AKANTU_DEBUG_OUT();
}
static bool material_is_alocated_local_damage [[gnu::unused]] =
MaterialFactory::getInstance().registerAllocator(
"local_damage",
[](UInt, const ID &, SolidMechanicsModel & model,
const ID & id) -> std::unique_ptr<Material> {
return std::make_unique<LocalMaterialDamage>(model, id);
});
} // namespace akantu
Event Timeline
Log In to Comment