Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007 #include <Intrepid_MiniTensor.h>
00008 #include "Teuchos_TestForException.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010
00011 namespace LCM
00012 {
00013
00014
00015 template<typename EvalT, typename Traits>
00016 StVenantKirchhoffModel<EvalT, Traits>::
00017 StVenantKirchhoffModel(Teuchos::ParameterList* p,
00018 const Teuchos::RCP<Albany::Layouts>& dl) :
00019 LCM::ConstitutiveModel<EvalT, Traits>(p, dl)
00020 {
00021
00022 this->dep_field_map_.insert(std::make_pair("F", dl->qp_tensor));
00023 this->dep_field_map_.insert(std::make_pair("J", dl->qp_scalar));
00024 this->dep_field_map_.insert(std::make_pair("Poissons Ratio", dl->qp_scalar));
00025 this->dep_field_map_.insert(std::make_pair("Elastic Modulus", 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 template<typename EvalT, typename Traits>
00042 void StVenantKirchhoffModel<EvalT, Traits>::
00043 computeState(typename Traits::EvalData workset,
00044 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > dep_fields,
00045 std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > eval_fields)
00046 {
00047
00048 PHX::MDField<ScalarT> def_grad = *dep_fields["F"];
00049 PHX::MDField<ScalarT> J = *dep_fields["J"];
00050 PHX::MDField<ScalarT> poissons_ratio = *dep_fields["Poissons Ratio"];
00051 PHX::MDField<ScalarT> elastic_modulus = *dep_fields["Elastic Modulus"];
00052
00053 std::string cauchy = (*field_name_map_)["Cauchy_Stress"];
00054 PHX::MDField<ScalarT> stress = *eval_fields[cauchy];
00055 ScalarT lambda;
00056 ScalarT mu;
00057
00058 Intrepid::Tensor<ScalarT> F(num_dims_), C(num_dims_), sigma(num_dims_);
00059 Intrepid::Tensor<ScalarT> I(Intrepid::eye<ScalarT>(num_dims_));
00060 Intrepid::Tensor<ScalarT> S(num_dims_), E(num_dims_);
00061
00062 for (std::size_t cell(0); cell < workset.numCells; ++cell) {
00063 for (std::size_t pt(0); pt < num_pts_; ++pt) {
00064 lambda = (elastic_modulus(cell, pt) * poissons_ratio(cell, pt))
00065 / (1. + poissons_ratio(cell, pt))
00066 / (1 - 2 * poissons_ratio(cell, pt));
00067 mu = elastic_modulus(cell, pt) / (2. * (1. + poissons_ratio(cell, pt)));
00068 F.fill(&def_grad(cell, pt, 0, 0));
00069 C = F * transpose(F);
00070 E = 0.5 * ( C - I );
00071 S = lambda * Intrepid::trace(E) * I + 2.0 * mu * E;
00072 sigma = (1.0 / Intrepid::det(F) ) * F * S * Intrepid::transpose(F);
00073 for (std::size_t i = 0; i < num_dims_; ++i) {
00074 for (std::size_t j = 0; j < num_dims_; ++j) {
00075 stress(cell, pt, i, j) = sigma(i, j);
00076 }
00077 }
00078 }
00079 }
00080 }
00081
00082 }
00083