• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

ConstitutiveModelInterface_Def.hpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 
00007 #include "Teuchos_TestForException.hpp"
00008 #include "Teuchos_RCP.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010 
00011 #include "AnisotropicDamageModel.hpp"
00012 #include "AnisotropicHyperelasticDamageModel.hpp"
00013 #include "ElasticDamageModel.hpp"
00014 #include "GursonHMRModel.hpp"
00015 #include "GursonModel.hpp"
00016 #include "J2FiberModel.hpp"
00017 #include "J2Model.hpp"
00018 #include "CreepModel.hpp"
00019 #include "MooneyRivlinModel.hpp"
00020 #include "NeohookeanModel.hpp"
00021 #include "RIHMRModel.hpp"
00022 #include "StVenantKirchhoffModel.hpp"
00023 #include "AAAModel.hpp"
00024 #include "LinearElasticModel.hpp"
00025 #include "HyperelasticDamageModel.hpp"
00026 #include "CapExplicitModel.hpp"
00027 #include "CapImplicitModel.hpp"
00028 #include "DruckerPragerModel.hpp"
00029 #include "CrystalPlasticityModel.hpp"
00030 #include "TvergaardHutchinsonModel.hpp"
00031 
00032 namespace LCM
00033 {
00034 
00035 //------------------------------------------------------------------------------
00036 template<typename EvalT, typename Traits>
00037 ConstitutiveModelInterface<EvalT, Traits>::
00038 ConstitutiveModelInterface(Teuchos::ParameterList& p,
00039                            const Teuchos::RCP<Albany::Layouts>& dl):
00040   have_temperature_(false),
00041   have_damage_(false),
00042   volume_average_pressure_(p.get<bool>("Volume Average Pressure", false))
00043 {
00044   Teuchos::ParameterList* plist = p.get<Teuchos::ParameterList*>("Material Parameters");
00045   plist->set<bool>("Volume Average Pressure", volume_average_pressure_);
00046   this->initializeModel(plist,dl);
00047 
00048   // construct the dependent fields
00049   std::map<std::string, Teuchos::RCP<PHX::DataLayout> >
00050   dependent_map = model_->getDependentFieldMap();
00051   typename std::map<std::string, Teuchos::RCP<PHX::DataLayout> >::iterator miter;
00052   for (miter = dependent_map.begin();
00053       miter != dependent_map.end();
00054       ++miter) {
00055     Teuchos::RCP<PHX::MDField<ScalarT> > temp_field =
00056         Teuchos::rcp(new PHX::MDField<ScalarT>(miter->first, miter->second));
00057     dep_fields_map_.insert(std::make_pair(miter->first, temp_field));
00058   }
00059 
00060   // register dependent fields
00061   typename std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > >::iterator it;
00062   for (it = dep_fields_map_.begin();
00063       it != dep_fields_map_.end();
00064       ++it) {
00065     this->addDependentField(*(it->second));
00066   }
00067 
00068   // optionally deal with integration point locations
00069   if (model_->getIntegrationPointLocationFlag()) {
00070     PHX::MDField<MeshScalarT, Cell, QuadPoint, Dim> cv("Coord Vec",
00071         dl->qp_vector);
00072     coord_vec_ = cv;
00073     this->addDependentField(coord_vec_);
00074   }
00075 
00076   // optionally deal with temperature
00077   if (p.isType<std::string>("Temperature Name")) {
00078     have_temperature_ = true;
00079     PHX::MDField<ScalarT, Cell, QuadPoint> t(p.get<std::string>("Temperature Name"),
00080         dl->qp_scalar);
00081     temperature_ = t;
00082     this->addDependentField(temperature_);
00083   }
00084 
00085   // optionally deal with damage
00086   if (p.isType<std::string>("Damage Name")) {
00087     have_damage_ = true;
00088     PHX::MDField<ScalarT, Cell, QuadPoint> d(p.get<std::string>("Damage Name"),
00089         dl->qp_scalar);
00090     damage_ = d;
00091     this->addDependentField(damage_);
00092   }
00093 
00094   // optional volume averaging needs integration weights
00095   if (volume_average_pressure_) {
00096     PHX::MDField<MeshScalarT, Cell, QuadPoint> w(p.get<std::string>("Weights Name"),
00097         dl->qp_scalar);
00098     weights_ = w;
00099     this->addDependentField(weights_);
00100   }
00101 
00102   // construct the evaluated fields
00103   std::map<std::string, Teuchos::RCP<PHX::DataLayout> >
00104   eval_map = model_->getEvaluatedFieldMap();
00105   for (miter = eval_map.begin();
00106       miter != eval_map.end();
00107       ++miter) {
00108     Teuchos::RCP<PHX::MDField<ScalarT> > temp_field =
00109         Teuchos::rcp(new PHX::MDField<ScalarT>(miter->first, miter->second));
00110     eval_fields_map_.insert(std::make_pair(miter->first, temp_field));
00111   }
00112 
00113   // register evaluated fields
00114   for (it = eval_fields_map_.begin();
00115       it != eval_fields_map_.end();
00116       ++it) {
00117     this->addEvaluatedField(*(it->second));
00118   }
00119 
00120   this->setName("ConstitutiveModelInterface" + PHX::TypeString<EvalT>::value);
00121 }
00122 
00123 //------------------------------------------------------------------------------
00124 template<typename EvalT, typename Traits>
00125 void ConstitutiveModelInterface<EvalT, Traits>::
00126 postRegistrationSetup(typename Traits::SetupData d,
00127     PHX::FieldManager<Traits>& fm)
00128 {
00129   TEUCHOS_TEST_FOR_EXCEPTION(dep_fields_map_.size() == 0, std::logic_error,
00130       "something is wrong in the LCM::CMI");
00131   TEUCHOS_TEST_FOR_EXCEPTION(eval_fields_map_.size() == 0, std::logic_error,
00132       "something is wrong in the LCM::CMI");
00133   // dependent fields
00134   typename std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > >::iterator it;
00135   for (it = dep_fields_map_.begin();
00136       it != dep_fields_map_.end();
00137       ++it) {
00138     this->utils.setFieldData(*(it->second), fm);
00139   }
00140 
00141   // optionally deal with integration point locations
00142   if (model_->getIntegrationPointLocationFlag()) {
00143     this->utils.setFieldData(coord_vec_, fm);
00144     model_->setCoordVecField(coord_vec_);
00145   }
00146 
00147   // optionally deal with temperature
00148   if (have_temperature_) {
00149     this->utils.setFieldData(temperature_, fm);
00150     model_->setTemperatureField(temperature_);
00151   }
00152 
00153   // optionally deal with damage
00154   if (have_damage_) {
00155     this->utils.setFieldData(damage_, fm);
00156     model_->setDamageField(damage_);
00157   }
00158 
00159   // optionally deal with volume averaging
00160   if (volume_average_pressure_) {
00161     this->utils.setFieldData(weights_, fm);
00162     model_->setWeightsField(weights_);
00163   }
00164 
00165   // evaluated fields
00166   for (it = eval_fields_map_.begin();
00167       it != eval_fields_map_.end();
00168       ++it) {
00169     this->utils.setFieldData(*(it->second), fm);
00170   }
00171 }
00172 
00173 //------------------------------------------------------------------------------
00174 template<typename EvalT, typename Traits>
00175 void ConstitutiveModelInterface<EvalT, Traits>::
00176 evaluateFields(typename Traits::EvalData workset)
00177 {
00178   model_->computeState(workset, dep_fields_map_, eval_fields_map_);
00179   if (volume_average_pressure_) {
00180     model_->computeVolumeAverage(workset,dep_fields_map_, eval_fields_map_);
00181   }
00182 }
00183 
00184 //------------------------------------------------------------------------------
00185 template<typename EvalT, typename Traits>
00186 void ConstitutiveModelInterface<EvalT, Traits>::
00187 fillStateVariableStruct(int state_var)
00188 {
00189   sv_struct_.name = model_->getStateVarName(state_var);
00190   sv_struct_.data_layout = model_->getStateVarLayout(state_var);
00191   sv_struct_.init_type = model_->getStateVarInitType(state_var);
00192   sv_struct_.init_value = model_->getStateVarInitValue(state_var);
00193   sv_struct_.register_old_state = model_->getStateVarOldStateFlag(state_var);
00194   sv_struct_.output_to_exodus = model_->getStateVarOutputFlag(state_var);
00195 }
00196 
00197 //------------------------------------------------------------------------------
00198 template<typename EvalT, typename Traits>
00199 void ConstitutiveModelInterface<EvalT, Traits>::
00200 initializeModel(Teuchos::ParameterList* p,
00201     const Teuchos::RCP<Albany::Layouts>& dl)
00202 {
00203   std::string model_name =
00204       p->sublist("Material Model").get<std::string>("Model Name");
00205 
00206   if (model_name == "Neohookean") {
00207     this->model_ = Teuchos::rcp(new LCM::NeohookeanModel<EvalT, Traits>(p, dl));
00208   } else if (model_name == "Creep") {
00209     this->model_ = Teuchos::rcp(new LCM::CreepModel<EvalT, Traits>(p, dl));
00210   } else if (model_name == "J2") {
00211     this->model_ = Teuchos::rcp(new LCM::J2Model<EvalT, Traits>(p, dl));
00212   } else if (model_name == "CrystalPlasticity") {
00213     this->model_ = Teuchos::rcp(new LCM::CrystalPlasticityModel<EvalT, Traits>(p, dl));
00214   } else if (model_name == "AHD") {
00215     this->model_ = Teuchos::rcp(
00216         new LCM::AnisotropicHyperelasticDamageModel<EvalT, Traits>(p, dl));
00217   } else if (model_name == "Gurson") {
00218     this->model_ = Teuchos::rcp(new LCM::GursonModel<EvalT, Traits>(p, dl));
00219   } else if (model_name == "GursonHMR") {
00220     this->model_ = Teuchos::rcp(new LCM::GursonHMRModel<EvalT, Traits>(p, dl));
00221   } else if (model_name == "Mooney Rivlin") {
00222     this->model_ = Teuchos::rcp(
00223         new LCM::MooneyRivlinModel<EvalT, Traits>(p, dl));
00224   } else if (model_name == "RIHMR") {
00225     this->model_ = Teuchos::rcp(new LCM::RIHMRModel<EvalT, Traits>(p, dl));
00226   } else if (model_name == "J2Fiber") {
00227     this->model_ = Teuchos::rcp(new LCM::J2FiberModel<EvalT, Traits>(p, dl));
00228   } else if (model_name == "Anisotropic Damage") {
00229     this->model_ = Teuchos::rcp(
00230         new LCM::AnisotropicDamageModel<EvalT, Traits>(p, dl));
00231   } else if (model_name == "Elastic Damage") {
00232     this->model_ = Teuchos::rcp(
00233         new LCM::ElasticDamageModel<EvalT, Traits>(p, dl));
00234   } else if (model_name == "Saint Venant Kirchhoff") {
00235     this->model_ = Teuchos::rcp(
00236         new LCM::StVenantKirchhoffModel<EvalT, Traits>(p, dl));
00237   } else if (model_name == "AAA") {
00238     this->model_ = Teuchos::rcp(new LCM::AAAModel<EvalT, Traits>(p, dl));
00239   } else if (model_name == "Linear Elastic") {
00240     this->model_ = Teuchos::rcp(
00241         new LCM::LinearElasticModel<EvalT, Traits>(p, dl));
00242   } else if (model_name == "Hyperelastic Damage") {
00243     this->model_ = Teuchos::rcp(
00244         new LCM::HyperelasticDamageModel<EvalT, Traits>(p, dl));
00245   } else if (model_name == "Cap Explicit") {
00246     this->model_ = Teuchos::rcp(
00247         new LCM::CapExplicitModel<EvalT, Traits>(p, dl));
00248   } else if (model_name == "Cap Implicit") {
00249       this->model_ = Teuchos::rcp(
00250         new LCM::CapImplicitModel<EvalT, Traits>(p, dl));
00251   } else if (model_name == "Drucker Prager") {
00252       this->model_ = Teuchos::rcp(
00253         new LCM::DruckerPragerModel<EvalT, Traits>(p, dl));
00254   } else if (model_name == "Tvergaard Hutchinson") {
00255       this->model_ = Teuchos::rcp(
00256         new LCM::TvergaardHutchinsonModel<EvalT, Traits>(p, dl));
00257   } else {
00258     TEUCHOS_TEST_FOR_EXCEPTION(true,
00259         std::logic_error,
00260         "Undefined material model name");
00261   }
00262 }
00263 
00264 //------------------------------------------------------------------------------
00265 }
00266 

Generated on Wed Mar 26 2014 18:36:37 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1