ADVonMises.h
Go to the documentation of this file.
1 /* ---------------------------------------------------------------------
2  * _
3  * _ __ ___ __ _ _ __ _ __ ___ ___ | |_
4  * | '_ ` _ \ / _` | '__| '_ ` _ \ / _ \| __|
5  * | | | | | | (_| | | | | | | | | (_) | |_
6  * |_| |_| |_|\__,_|_| |_| |_| |_|\___/ \__|
7  *
8  * Unit of Strength of Materials and Structural Analysis
9  * University of Innsbruck
10  * 2020 - today
11  *
12  * festigkeitslehre@uibk.ac.at
13  *
14  * Alexander Dummer alexander.dummer@uibk.ac.at
15  *
16  * This file is part of the MAteRialMOdellingToolbox (marmot).
17  *
18  * This library is free software; you can redistribute it and/or
19  * modify it under the terms of the GNU Lesser General Public
20  * License as published by the Free Software Foundation; either
21  * version 2.1 of the License, or (at your option) any later version.
22  *
23  * The full text of the license can be found in the file LICENSE.md at
24  * the top level directory of marmot.
25  * ---------------------------------------------------------------------
26  */
27 
28 #pragma once
29 #include "Marmot/MarmotConstants.h"
32 #include "Marmot/MarmotTypedefs.h"
33 #include "autodiff/forward/dual.hpp"
34 #include <Eigen/src/Core/Map.h>
35 #include <iostream>
36 #include <string>
37 #include <vector>
38 
39 using namespace Marmot;
40 
41 namespace Marmot::Materials {
48  public:
49  using MarmotMaterialHypoElasticAD::MarmotMaterialHypoElasticAD;
50 
51  const double& E;
52  const double& nu;
53  const double& yieldStress;
54  const double& HLin;
55  const double& deltaYieldStress;
56  const double& delta;
57  const double G;
58 
60 
61  protected:
62  void computeStressAD( autodiff::dual* stress,
63  const autodiff::dual* dStrain,
64  const double* timeOld,
65  const double dT,
66  double& pNewDT ) override;
67 
69 
70  public:
71  inline const static auto layout = makeLayout( {
72  { .name = "kappa", .length = 1 },
73  } );
74 
75  double& kappa;
76 
77  ADVonMisesModelStateVarManager( double* theStateVarVector )
78  : MarmotStateVarVectorManager( theStateVarVector, layout ), kappa( find( "kappa" ) ){};
79  };
80  std::unique_ptr< ADVonMisesModelStateVarManager > managedStateVars;
81 
82  int getNumberOfRequiredStateVars() override { return ADVonMisesModelStateVarManager::layout.nRequiredStateVars; }
83 
84  void assignStateVars( double* stateVars, int nStateVars ) override;
85 
86  StateView getStateView( const std::string& result ) override;
87 
88  template < typename T >
89  T fy( T kappa_ )
90  {
91  const T res = yieldStress + HLin * kappa_ + deltaYieldStress * ( 1. - exp( -delta * kappa_ ) );
92  return res;
93  }
94 
95  template < typename T >
96  T f( const T rho_, const double kappa_ )
97  {
98  return rho_ - Constants::sqrt2_3 * fy( kappa_ );
99  }
100  template < typename T >
101  T g( const T rhoTrial, const double kappa, const T deltaKappa )
102  {
103  const T kappa_ = kappa + deltaKappa;
104  return rhoTrial - Constants::sqrt6 * G * deltaKappa - Constants::sqrt2_3 * fy( kappa_ );
105  }
106  };
107 } // namespace Marmot::Materials
Marmot::Materials::ADVonMises::getNumberOfRequiredStateVars
int getNumberOfRequiredStateVars() override
Definition: ADVonMises.h:82
MarmotConstants.h
Marmot::Constants::sqrt2_3
constexpr double sqrt2_3
Definition: MarmotConstants.h:51
Marmot::Materials::ADVonMises::yieldStress
const double & yieldStress
Definition: ADVonMises.h:53
Marmot::Materials::ADVonMises::fy
T fy(T kappa_)
Definition: ADVonMises.h:89
Marmot::Math::exp
double exp(double x)
Definition: MarmotMath.cpp:13
Marmot::Materials::ADVonMises
Implementation of a isotropic J2-plasticity material for 3D stress states using automatic differentia...
Definition: ADVonMises.h:47
Marmot::Constants::sqrt6
constexpr double sqrt6
Definition: MarmotConstants.h:55
MarmotTypedefs.h
Marmot::Materials
Definition: MarmotKelvinChain.h:34
Marmot::Materials::ADVonMises::ADVonMisesModelStateVarManager
Definition: ADVonMises.h:68
Marmot::Materials::ADVonMises::ADVonMisesModelStateVarManager::layout
static const auto layout
Definition: ADVonMises.h:71
MarmotMaterialHypoElasticAD.h
MarmotStateVarVectorManager
A convenience auxiliary class for managing multiple statevars with arbitrary length in a single conse...
Definition: MarmotStateVarVectorManager.h:37
MarmotStateVarVectorManager.h
Marmot::Materials::ADVonMises::deltaYieldStress
const double & deltaYieldStress
Definition: ADVonMises.h:55
Marmot::Materials::ADVonMises::E
const double & E
Definition: ADVonMises.h:51
Marmot::Materials::ADVonMises::ADVonMisesModelStateVarManager::kappa
double & kappa
Definition: ADVonMises.h:75
StateView
Definition: MarmotUtils.h:29
Marmot
This file includes functions needed for calculations with stress and strain tensors written in voigt ...
Definition: MarmotTesting.h:37
Marmot::Materials::ADVonMises::g
T g(const T rhoTrial, const double kappa, const T deltaKappa)
Definition: ADVonMises.h:101
Marmot::Materials::ADVonMises::assignStateVars
void assignStateVars(double *stateVars, int nStateVars) override
Definition: ADVonMises.cpp:37
Marmot::Materials::ADVonMises::managedStateVars
std::unique_ptr< ADVonMisesModelStateVarManager > managedStateVars
Definition: ADVonMises.h:80
MarmotStateVarVectorManager::makeLayout
static StateVarVectorLayout makeLayout(const std::vector< StateVarEntryDefinition > &theEntries)
generate the statevar vector layout from a list of entries, defined by name and length
Definition: MarmotStateVarVectorManager.h:74
MarmotMaterialHypoElasticAD
Definition: MarmotMaterialHypoElasticAD.h:32
Marmot::Materials::ADVonMises::getStateView
StateView getStateView(const std::string &result) override
Definition: ADVonMises.cpp:46
Marmot::Materials::ADVonMises::G
const double G
Definition: ADVonMises.h:57
Marmot::Materials::ADVonMises::computeStressAD
void computeStressAD(autodiff::dual *stress, const autodiff::dual *dStrain, const double *timeOld, const double dT, double &pNewDT) override
Definition: ADVonMises.cpp:50
Marmot::Materials::ADVonMises::ADVonMises
ADVonMises(const double *materialProperties, int nMaterialProperties, int materialNumber)
Definition: ADVonMises.cpp:22
MarmotMaterial::stateVars
double * stateVars
Definition: MarmotMaterial.h:38
Marmot::Materials::ADVonMises::f
T f(const T rho_, const double kappa_)
Definition: ADVonMises.h:96
Marmot::Materials::ADVonMises::delta
const double & delta
Definition: ADVonMises.h:56
Marmot::Materials::ADVonMises::HLin
const double & HLin
Definition: ADVonMises.h:54
Marmot::Materials::ADVonMises::ADVonMisesModelStateVarManager::ADVonMisesModelStateVarManager
ADVonMisesModelStateVarManager(double *theStateVarVector)
Definition: ADVonMises.h:77
Marmot::Materials::ADVonMises::nu
const double & nu
Definition: ADVonMises.h:52
MarmotMaterial::materialNumber
const int materialNumber
Definition: MarmotMaterial.h:42
MarmotMaterial::nStateVars
int nStateVars
Definition: MarmotMaterial.h:39
MarmotMaterial::nMaterialProperties
const int nMaterialProperties
Definition: MarmotMaterial.h:36
MarmotMaterial::materialProperties
const double * materialProperties
Definition: MarmotMaterial.h:35
MarmotStateVarVectorManager::find
double & find(const std::string &name) const
get the reference to the first array element of an entry in the statevar vector
Definition: MarmotStateVarVectorManager.h:48