CapImplicit stress response. More...
#include <CapImplicitModel.hpp>
Public Types | |
typedef EvalT::ScalarT | ScalarT |
typedef EvalT::MeshScalarT | MeshScalarT |
typedef Sacado::mpl::apply < FadType, ScalarT >::type | DFadType |
typedef Sacado::mpl::apply < FadType, DFadType >::type | D2FadType |
Public Member Functions | |
CapImplicitModel (Teuchos::ParameterList *p, const Teuchos::RCP< Albany::Layouts > &dl) | |
Constructor. | |
virtual | ~CapImplicitModel () |
Virtual Destructor. | |
virtual void | computeState (typename Traits::EvalData workset, std::map< std::string, Teuchos::RCP< PHX::MDField< ScalarT > > > dep_fields, std::map< std::string, Teuchos::RCP< PHX::MDField< ScalarT > > > eval_fields) |
Implementation of physics. | |
Private Member Functions | |
CapImplicitModel (const CapImplicitModel &) | |
Private to prohibit copying. | |
CapImplicitModel & | operator= (const CapImplicitModel &) |
Private to prohibit copying. | |
ScalarT | compute_f (Intrepid::Tensor< ScalarT > &sigma, Intrepid::Tensor< ScalarT > &alpha, ScalarT &kappa) |
std::vector< ScalarT > | initialize (Intrepid::Tensor< ScalarT > &sigmaVal, Intrepid::Tensor< ScalarT > &alphaVal, ScalarT &kappaVal, ScalarT &dgammaVal) |
void | compute_ResidJacobian (std::vector< ScalarT > const &XXVal, std::vector< ScalarT > &R, std::vector< ScalarT > &dRdX, const Intrepid::Tensor< ScalarT > &sigmaVal, const Intrepid::Tensor< ScalarT > &alphaVal, const ScalarT &kappaVal, Intrepid::Tensor4< ScalarT > const &Celastic, bool kappa_flag) |
DFadType | compute_f (Intrepid::Tensor< DFadType > &sigma, Intrepid::Tensor< DFadType > &alpha, DFadType &kappa) |
D2FadType | compute_g (Intrepid::Tensor< D2FadType > &sigma, Intrepid::Tensor< D2FadType > &alpha, D2FadType &kappa) |
Intrepid::Tensor< DFadType > | compute_dgdsigma (std::vector< DFadType > const &XX) |
DFadType | compute_Galpha (DFadType J2_alpha) |
Intrepid::Tensor< DFadType > | compute_halpha (Intrepid::Tensor< DFadType > const &dgdsigma, DFadType const J2_alpha) |
DFadType | compute_dedkappa (DFadType const kappa) |
DFadType | compute_hkappa (DFadType const I1_dgdsigma, DFadType const dedkappa) |
Private Attributes | |
RealType | A |
RealType | B |
RealType | C |
RealType | theta |
RealType | R |
RealType | kappa0 |
RealType | W |
RealType | D1 |
RealType | D2 |
RealType | calpha |
RealType | psi |
RealType | N |
RealType | L |
RealType | phi |
RealType | Q |
std::string | strainName |
std::string | stressName |
std::string | backStressName |
std::string | capParameterName |
std::string | eqpsName |
std::string | volPlasticStrainName |
CapImplicit stress response.
This evaluator computes stress based on a cap plasticity model.
Definition at line 27 of file CapImplicitModel.hpp.
typedef EvalT::ScalarT LCM::CapImplicitModel< EvalT, Traits >::ScalarT |
Reimplemented from LCM::ConstitutiveModel< EvalT, Traits >.
Definition at line 32 of file CapImplicitModel.hpp.
typedef EvalT::MeshScalarT LCM::CapImplicitModel< EvalT, Traits >::MeshScalarT |
Reimplemented from LCM::ConstitutiveModel< EvalT, Traits >.
Definition at line 33 of file CapImplicitModel.hpp.
typedef Sacado::mpl::apply<FadType,ScalarT>::type LCM::CapImplicitModel< EvalT, Traits >::DFadType |
Definition at line 34 of file CapImplicitModel.hpp.
typedef Sacado::mpl::apply<FadType,DFadType>::type LCM::CapImplicitModel< EvalT, Traits >::D2FadType |
Definition at line 35 of file CapImplicitModel.hpp.
LCM::CapImplicitModel< EvalT, Traits >::CapImplicitModel | ( | Teuchos::ParameterList * | p, | |
const Teuchos::RCP< Albany::Layouts > & | dl | |||
) |
Constructor.
Definition at line 21 of file CapImplicitModel_Def.hpp.
virtual LCM::CapImplicitModel< EvalT, Traits >::~CapImplicitModel | ( | ) | [inline, virtual] |
Virtual Destructor.
Definition at line 51 of file CapImplicitModel.hpp.
LCM::CapImplicitModel< EvalT, Traits >::CapImplicitModel | ( | const CapImplicitModel< EvalT, Traits > & | ) | [private] |
Private to prohibit copying.
void LCM::CapImplicitModel< EvalT, Traits >::computeState | ( | typename Traits::EvalData | workset, | |
std::map< std::string, Teuchos::RCP< PHX::MDField< ScalarT > > > | dep_fields, | |||
std::map< std::string, Teuchos::RCP< PHX::MDField< ScalarT > > > | eval_fields | |||
) | [virtual] |
Implementation of physics.
Definition at line 123 of file CapImplicitModel_Def.hpp.
CapImplicitModel& LCM::CapImplicitModel< EvalT, Traits >::operator= | ( | const CapImplicitModel< EvalT, Traits > & | ) | [private] |
Private to prohibit copying.
CapImplicitModel< EvalT, Traits >::ScalarT LCM::CapImplicitModel< EvalT, Traits >::compute_f | ( | Intrepid::Tensor< ScalarT > & | sigma, | |
Intrepid::Tensor< ScalarT > & | alpha, | |||
ScalarT & | kappa | |||
) | [private] |
Definition at line 356 of file CapImplicitModel_Def.hpp.
std::vector< typename CapImplicitModel< EvalT, Traits >::ScalarT > LCM::CapImplicitModel< EvalT, Traits >::initialize | ( | Intrepid::Tensor< ScalarT > & | sigmaVal, | |
Intrepid::Tensor< ScalarT > & | alphaVal, | |||
ScalarT & | kappaVal, | |||
ScalarT & | dgammaVal | |||
) | [private] |
Definition at line 470 of file CapImplicitModel_Def.hpp.
void LCM::CapImplicitModel< EvalT, Traits >::compute_ResidJacobian | ( | std::vector< ScalarT > const & | XXVal, | |
std::vector< ScalarT > & | R, | |||
std::vector< ScalarT > & | dRdX, | |||
const Intrepid::Tensor< ScalarT > & | sigmaVal, | |||
const Intrepid::Tensor< ScalarT > & | alphaVal, | |||
const ScalarT & | kappaVal, | |||
Intrepid::Tensor4< ScalarT > const & | Celastic, | |||
bool | kappa_flag | |||
) | [private] |
Definition at line 614 of file CapImplicitModel_Def.hpp.
CapImplicitModel< EvalT, Traits >::DFadType LCM::CapImplicitModel< EvalT, Traits >::compute_f | ( | Intrepid::Tensor< DFadType > & | sigma, | |
Intrepid::Tensor< DFadType > & | alpha, | |||
DFadType & | kappa | |||
) | [private] |
Definition at line 394 of file CapImplicitModel_Def.hpp.
CapImplicitModel< EvalT, Traits >::D2FadType LCM::CapImplicitModel< EvalT, Traits >::compute_g | ( | Intrepid::Tensor< D2FadType > & | sigma, | |
Intrepid::Tensor< D2FadType > & | alpha, | |||
D2FadType & | kappa | |||
) | [private] |
Definition at line 432 of file CapImplicitModel_Def.hpp.
Intrepid::Tensor< typename CapImplicitModel< EvalT, Traits >::DFadType > LCM::CapImplicitModel< EvalT, Traits >::compute_dgdsigma | ( | std::vector< DFadType > const & | XX | ) | [private] |
Definition at line 494 of file CapImplicitModel_Def.hpp.
CapImplicitModel< EvalT, Traits >::DFadType LCM::CapImplicitModel< EvalT, Traits >::compute_Galpha | ( | DFadType | J2_alpha | ) | [private] |
Definition at line 546 of file CapImplicitModel_Def.hpp.
Intrepid::Tensor< typename CapImplicitModel< EvalT, Traits >::DFadType > LCM::CapImplicitModel< EvalT, Traits >::compute_halpha | ( | Intrepid::Tensor< DFadType > const & | dgdsigma, | |
DFadType const | J2_alpha | |||
) | [private] |
Definition at line 557 of file CapImplicitModel_Def.hpp.
CapImplicitModel< EvalT, Traits >::DFadType LCM::CapImplicitModel< EvalT, Traits >::compute_dedkappa | ( | DFadType const | kappa | ) | [private] |
Definition at line 581 of file CapImplicitModel_Def.hpp.
CapImplicitModel< EvalT, Traits >::DFadType LCM::CapImplicitModel< EvalT, Traits >::compute_hkappa | ( | DFadType const | I1_dgdsigma, | |
DFadType const | dedkappa | |||
) | [private] |
Definition at line 604 of file CapImplicitModel_Def.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::A [private] |
Definition at line 117 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::B [private] |
Definition at line 118 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::C [private] |
Definition at line 119 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::theta [private] |
Definition at line 120 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::R [private] |
Definition at line 121 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::kappa0 [private] |
Definition at line 122 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::W [private] |
Definition at line 123 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::D1 [private] |
Definition at line 124 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::D2 [private] |
Definition at line 125 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::calpha [private] |
Definition at line 126 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::psi [private] |
Definition at line 127 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::N [private] |
Definition at line 128 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::L [private] |
Definition at line 129 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::phi [private] |
Definition at line 130 of file CapImplicitModel.hpp.
RealType LCM::CapImplicitModel< EvalT, Traits >::Q [private] |
Definition at line 131 of file CapImplicitModel.hpp.
std::string LCM::CapImplicitModel< EvalT, Traits >::strainName [private] |
Definition at line 133 of file CapImplicitModel.hpp.
std::string LCM::CapImplicitModel< EvalT, Traits >::stressName [private] |
Definition at line 133 of file CapImplicitModel.hpp.
std::string LCM::CapImplicitModel< EvalT, Traits >::backStressName [private] |
Definition at line 134 of file CapImplicitModel.hpp.
std::string LCM::CapImplicitModel< EvalT, Traits >::capParameterName [private] |
Definition at line 134 of file CapImplicitModel.hpp.
std::string LCM::CapImplicitModel< EvalT, Traits >::eqpsName [private] |
Definition at line 134 of file CapImplicitModel.hpp.
std::string LCM::CapImplicitModel< EvalT, Traits >::volPlasticStrainName [private] |
Definition at line 134 of file CapImplicitModel.hpp.