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

ConstitutiveModel_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 <Intrepid_MiniTensor.h>
00008 #include <Phalanx_DataLayout.hpp>
00009 #include <Teuchos_TestForException.hpp>
00010 
00011 namespace LCM
00012 {
00013 
00014 //------------------------------------------------------------------------------
00015 template<typename EvalT, typename Traits>
00016 ConstitutiveModel<EvalT, Traits>::
00017 ConstitutiveModel(Teuchos::ParameterList* p,
00018     const Teuchos::RCP<Albany::Layouts>& dl) :
00019     num_state_variables_(0),
00020     compute_energy_(false),
00021     compute_tangent_(false),
00022     need_integration_pt_locations_(false),
00023     have_temperature_(false),
00024     have_damage_(false)
00025 {
00026   // extract number of integration points and dimensions
00027   std::vector<PHX::DataLayout::size_type> dims;
00028   dl->qp_tensor->dimensions(dims);
00029   num_pts_  = dims[1];
00030   num_dims_ = dims[2];
00031   field_name_map_ =
00032       p->get<Teuchos::RCP<std::map<std::string, std::string> > >("Name Map");
00033 
00034   if (p->isType<bool>("Have Temperature")) {
00035     if (p->get<bool>("Have Temperature")) {
00036       have_temperature_ = true;
00037       expansion_coeff_ = p->get<RealType>("Thermal Expansion Coefficient", 0.0);
00038       ref_temperature_ = p->get<RealType>("Reference Temperature", 0.0);
00039       heat_capacity_ = p->get<RealType>("Heat Capacity", 1.0);
00040       density_ = p->get<RealType>("Density", 1.0);
00041     }
00042   }
00043 
00044   if (p->isType<bool>("Have Damage")) {
00045     if (p->get<bool>("Have Damage")) {
00046       have_damage_ = true;
00047     }
00048   }
00049 }
00050 //------------------------------------------------------------------------------
00051 template<typename EvalT, typename Traits>
00052 void ConstitutiveModel<EvalT, Traits>::
00053 computeVolumeAverage(typename Traits::EvalData workset,
00054     std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > dep_fields,
00055     std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > eval_fields)
00056 {
00057   Intrepid::Tensor<ScalarT> sig(num_dims_);
00058   Intrepid::Tensor<ScalarT> I(Intrepid::eye<ScalarT>(num_dims_));
00059 
00060   std::string cauchy = (*field_name_map_)["Cauchy_Stress"];
00061   PHX::MDField<ScalarT> stress = *eval_fields[cauchy];
00062 
00063   ScalarT volume, pbar, p;
00064 
00065   for (std::size_t cell(0); cell < workset.numCells; ++cell) {
00066     volume = pbar = 0.0;
00067     for (std::size_t pt(0); pt < num_pts_; ++pt) {
00068       sig.fill(&stress(cell,pt,0,0));
00069       pbar += weights_(cell,pt) * (1./num_dims_) * Intrepid::trace(sig);
00070       volume += weights_(cell,pt);
00071     }
00072 
00073     pbar /= volume;
00074 
00075     for (std::size_t pt(0); pt < num_pts_; ++pt) {
00076       sig.fill(&stress(cell,pt,0,0));
00077       p = (1./num_dims_) * Intrepid::trace(sig);
00078       sig += (pbar - p)*I;
00079 
00080       for (std::size_t i = 0; i < num_dims_; ++i) {
00081         stress(cell,pt,i,i) = sig(i,i);
00082       }
00083     }
00084   }
00085 }
00086 //------------------------------------------------------------------------------
00087 template<typename EvalT, typename Traits>
00088 std::string ConstitutiveModel<EvalT, Traits>::
00089 getStateVarName(int state_var)
00090 {
00091   return state_var_names_[state_var];
00092 }
00093 //------------------------------------------------------------------------------
00094 template<typename EvalT, typename Traits>
00095 Teuchos::RCP<PHX::DataLayout> ConstitutiveModel<EvalT, Traits>::
00096 getStateVarLayout(int state_var)
00097 {
00098   return state_var_layouts_[state_var];
00099 }
00100 //------------------------------------------------------------------------------
00101 template<typename EvalT, typename Traits>
00102 std::string ConstitutiveModel<EvalT, Traits>::
00103 getStateVarInitType(int state_var)
00104 {
00105   return state_var_init_types_[state_var];
00106 }
00107 //------------------------------------------------------------------------------
00108 template<typename EvalT, typename Traits>
00109 double ConstitutiveModel<EvalT, Traits>::
00110 getStateVarInitValue(int state_var)
00111 {
00112   return state_var_init_values_[state_var];
00113 }
00114 //------------------------------------------------------------------------------
00115 template<typename EvalT, typename Traits>
00116 bool ConstitutiveModel<EvalT, Traits>::
00117 getStateVarOldStateFlag(int state_var)
00118 {
00119   return state_var_old_state_flags_[state_var];
00120 }
00121 //------------------------------------------------------------------------------
00122 template<typename EvalT, typename Traits>
00123 bool ConstitutiveModel<EvalT, Traits>::
00124 getStateVarOutputFlag(int state_var)
00125 {
00126   return state_var_output_flags_[state_var];
00127 }
00128 //------------------------------------------------------------------------------
00129 }
00130 

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