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

HyperelasticDamageModel_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 "Teuchos_TestForException.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010 
00011 namespace LCM
00012 {
00013 
00014 //------------------------------------------------------------------------------
00015 template<typename EvalT, typename Traits>
00016 HyperelasticDamageModel<EvalT, Traits>::
00017 HyperelasticDamageModel(Teuchos::ParameterList* p,
00018     const Teuchos::RCP<Albany::Layouts>& dl) :
00019     LCM::ConstitutiveModel<EvalT, Traits>(p, dl),
00020       max_damage_(p->get<RealType>("Maximum Damage", 1.0)),
00021       damage_saturation_(p->get<RealType>("Damage Saturation", 1.0))
00022 {
00023   // define the dependent fields
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   this->dep_field_map_.insert(std::make_pair("Poissons Ratio", dl->qp_scalar));
00027   this->dep_field_map_.insert(std::make_pair("Elastic Modulus", dl->qp_scalar));
00028   this->dep_field_map_.insert(std::make_pair("Delta Time", dl->workset_scalar));
00029 
00030   // define the evaluated fields
00031   std::string cauchy = (*field_name_map_)["Cauchy_Stress"];
00032   this->eval_field_map_.insert(std::make_pair(cauchy, dl->qp_tensor));
00033   this->eval_field_map_.insert(std::make_pair("Damage_Source", dl->qp_scalar));
00034   this->eval_field_map_.insert(std::make_pair("alpha", dl->qp_scalar));
00035 
00036   // deal with damage
00037   if (have_damage_) {
00038     this->dep_field_map_.insert(std::make_pair("damage", dl->qp_scalar));
00039   } else {
00040     this->eval_field_map_.insert(std::make_pair("local damage", dl->qp_scalar));
00041   }
00042 
00043   // define the state variables
00044   //
00045   this->num_state_variables_++;
00046   this->state_var_names_.push_back(cauchy);
00047   this->state_var_layouts_.push_back(dl->qp_tensor);
00048   this->state_var_init_types_.push_back("scalar");
00049   this->state_var_init_values_.push_back(0.0);
00050   this->state_var_old_state_flags_.push_back(false);
00051   this->state_var_output_flags_.push_back(true);
00052   //
00053   //
00054   this->num_state_variables_++;
00055   this->state_var_names_.push_back("alpha");
00056   this->state_var_layouts_.push_back(dl->qp_scalar);
00057   this->state_var_init_types_.push_back("scalar");
00058   this->state_var_init_values_.push_back(0.0);
00059   this->state_var_old_state_flags_.push_back(true);
00060   this->state_var_output_flags_.push_back(true);
00061 }
00062 //------------------------------------------------------------------------------
00063 template<typename EvalT, typename Traits>
00064 void HyperelasticDamageModel<EvalT, Traits>::
00065 computeState(typename Traits::EvalData workset,
00066     std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > dep_fields,
00067     std::map<std::string, Teuchos::RCP<PHX::MDField<ScalarT> > > eval_fields)
00068 {
00069   // extract dependent MDFields
00070   PHX::MDField<ScalarT> def_grad = *dep_fields["F"];
00071   PHX::MDField<ScalarT> J = *dep_fields["J"];
00072   PHX::MDField<ScalarT> poissons_ratio = *dep_fields["Poissons Ratio"];
00073   PHX::MDField<ScalarT> elastic_modulus = *dep_fields["Elastic Modulus"];
00074   PHX::MDField<ScalarT> delta_time = *dep_fields["Delta Time"];
00075   // extract evaluated MDFields
00076   std::string cauchy = (*field_name_map_)["Cauchy_Stress"];
00077   PHX::MDField<ScalarT> stress = *eval_fields[cauchy];
00078   PHX::MDField<ScalarT> alpha = *eval_fields["alpha"];
00079   PHX::MDField<ScalarT> damage;
00080   if (!have_damage_) {
00081     damage = *eval_fields["local damage"];
00082   }
00083   PHX::MDField<ScalarT> source = *eval_fields["Damage_Source"];
00084   //get state variables
00085   Albany::MDArray alpha_old = (*workset.stateArrayPtr)["alpha_old"];
00086   ScalarT kappa;
00087   ScalarT mu, mubar;
00088   ScalarT Jm53, Jm23;
00089   ScalarT smag;
00090   ScalarT energy;
00091   ScalarT dt = delta_time(0);
00092 
00093   Intrepid::Tensor<ScalarT> F(num_dims_), b(num_dims_), sigma(num_dims_);
00094   Intrepid::Tensor<ScalarT> I(Intrepid::eye<ScalarT>(num_dims_));
00095   Intrepid::Tensor<ScalarT> s(num_dims_), n(num_dims_);
00096 
00097   for (std::size_t cell(0); cell < workset.numCells; ++cell) {
00098     for (std::size_t pt(0); pt < num_pts_; ++pt) {
00099       kappa = elastic_modulus(cell, pt)
00100           / (3. * (1. - 2. * poissons_ratio(cell, pt)));
00101       mu = elastic_modulus(cell, pt) / (2. * (1. + poissons_ratio(cell, pt)));
00102       Jm53 = std::pow(J(cell, pt), -5. / 3.);
00103       Jm23 = Jm53 * J(cell, pt);
00104 
00105       F.fill(&def_grad(cell, pt, 0, 0));
00106       b = F * transpose(F);
00107       mubar = (1.0 / 3.0) * mu * Jm23 * Intrepid::trace(b);
00108       sigma = 0.5 * kappa * (J(cell, pt) - 1. / J(cell, pt)) * I
00109           + mu * Jm53 * Intrepid::dev(b);
00110 
00111       energy = 0.5 * kappa
00112           * (0.5 * (J(cell, pt) * J(cell, pt) - 1.0) - std::log(J(cell, pt)))
00113           + 0.5 * mu * (Jm23 * Intrepid::trace(b) - 3.0);
00114 
00115       if (have_temperature_) {
00116         ScalarT delta_temp = temperature_(cell, pt) - ref_temperature_;
00117         energy += heat_capacity_
00118             * ((delta_temp) - temperature_(cell, pt)
00119                 * std::log(temperature_(cell, pt) / ref_temperature_))
00120             - 3.0 * expansion_coeff_ * (J(cell, pt) - 1.0 / J(cell, pt))
00121                 * delta_temp;
00122         sigma -= expansion_coeff_ * (1.0 + 1.0 / (J(cell, pt) * J(cell, pt)))
00123             * delta_temp * I;
00124       }
00125 
00126       alpha(cell, pt) = std::max((ScalarT) alpha_old(cell, pt), energy);
00127 
00128       source(cell, pt) = (max_damage_ / damage_saturation_)
00129           * std::exp(-alpha(cell, pt) / damage_saturation_)
00130           * (alpha(cell, pt) - alpha_old(cell, pt)) / dt;
00131 
00132       if (!have_damage_) {
00133         damage(cell, pt) = max_damage_
00134             * (1.0 - std::exp(-alpha(cell, pt) / damage_saturation_));
00135       }
00136 
00137       for (std::size_t i = 0; i < num_dims_; ++i) {
00138         for (std::size_t j = 0; j < num_dims_; ++j) {
00139           stress(cell, pt, i, j) = (1.0 - damage(cell, pt)) * sigma(i, j);
00140         }
00141       }
00142     }
00143   }
00144 }
00145 //------------------------------------------------------------------------------
00146 }
00147 

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