Go to the documentation of this file.00001
00002
00003
00004
00005
00006 #include <Intrepid_MiniTensor.h>
00007 #include <Teuchos_TestForException.hpp>
00008 #include <Phalanx_DataLayout.hpp>
00009
00010 namespace LCM
00011 {
00012
00013
00014 template<typename EvalT, typename Traits>
00015 AAAModel<EvalT, Traits>::
00016 AAAModel(Teuchos::ParameterList* p,
00017 const Teuchos::RCP<Albany::Layouts>& dl) :
00018 LCM::ConstitutiveModel<EvalT, Traits>(p, dl),
00019 alpha_(p->get<RealType>("alpha", 0.0)),
00020 beta_(p->get<RealType>("beta", 0.0)),
00021 mult_(p->get<RealType>("mult", 0.0))
00022 {
00023
00024 this->dep_field_map_.insert(std::make_pair("F", dl->qp_tensor));
00025 this->dep_field_map_.insert(std::make_pair("J", dl->qp_scalar));
00026
00027
00028 std::string cauchy = (*field_name_map_)["Cauchy_Stress"];
00029 this->eval_field_map_.insert(std::make_pair(cauchy, dl->qp_tensor));
00030
00031
00032 this->num_state_variables_++;
00033 this->state_var_names_.push_back(cauchy);
00034 this->state_var_layouts_.push_back(dl->qp_tensor);
00035 this->state_var_init_types_.push_back("scalar");
00036 this->state_var_init_values_.push_back(0.0);
00037 this->state_var_old_state_flags_.push_back(false);
00038 this->state_var_output_flags_.push_back(true);
00039 }
00040
00041
00042 template<typename EvalT, typename Traits>
00043 void AAAModel<EvalT, Traits>::
00044 computeState(typename Traits::EvalData workset,
00045 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > dep_fields,
00046 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > eval_fields)
00047 {
00048
00049 PHX::MDField<ScalarT> defGrad = *dep_fields["F"];
00050 PHX::MDField<ScalarT> J = *dep_fields["J"];
00051
00052 std::string cauchy = (*field_name_map_)["Cauchy_Stress"];
00053 PHX::MDField<ScalarT> stress = *eval_fields[cauchy];
00054
00055 Intrepid::Tensor<ScalarT> F(num_dims_);
00056 Intrepid::Tensor<ScalarT> S(num_dims_);
00057 Intrepid::Tensor<ScalarT> B(num_dims_);
00058 Intrepid::Tensor<ScalarT> Id = Intrepid::identity<ScalarT>(num_dims_);
00059
00060
00061 ScalarT mu = 2.0 * (alpha_);
00062
00063
00064 ScalarT kappa = mult_ * mu;
00065
00066 for (std::size_t cell(0); cell < workset.numCells; ++cell) {
00067 for (std::size_t pt(0); pt < num_pts_; ++pt) {
00068 F.fill(&defGrad(cell, pt, 0, 0));
00069 B = F * Intrepid::transpose(F);
00070
00071 ScalarT pressure = kappa * (J(cell, pt) - 1.0);
00072
00073
00074 S = -pressure * Id
00075 + 2.0 * (alpha_ + 2.0 * beta_ * (Intrepid::I1(B) - 3.0)) * B;
00076
00077 for (std::size_t i(0); i < num_dims_; ++i) {
00078 for (std::size_t j(0); j < num_dims_; ++j) {
00079 stress(cell, pt, i, j) = S(i, j);
00080 }
00081 }
00082 }
00083 }
00084 }
00085
00086 }