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

CMResidualCoarse_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_Mechanics.h>
00008 #include <Teuchos_TestForException.hpp>
00009 #include <Phalanx_DataLayout.hpp>
00010 #include <Sacado_ParameterRegistration.hpp>
00011 
00012 using std::string;
00013 using std::size_t;
00014 
00015 namespace LCM
00016 {
00017 
00018 //
00019 //
00020 //
00021 template<typename EvalT, typename Traits>
00022 CMResidualCoarse<EvalT, Traits>::
00023 CMResidualCoarse(Teuchos::ParameterList & p,
00024     Teuchos::RCP<Albany::Layouts> const & dl) :
00025       stress_(p.get<string>("Stress Name"), dl->qp_tensor),
00026       def_grad_(p.get<string>("DefGrad Name"), dl->qp_tensor),
00027       w_grad_bf_(p.get<string>("Weighted Gradient BF Name"),dl->node_qp_vector),
00028       w_bf_(p.get<string>("Weighted BF Name"), dl->node_qp_scalar),
00029       residual_(p.get<string>("Residual Name"), dl->node_vector),
00030       have_body_force_(p.get<bool>("Have Body Force", false))
00031 {
00032   this->addDependentField(stress_);
00033   this->addDependentField(def_grad_);
00034   this->addDependentField(w_grad_bf_);
00035   this->addDependentField(w_bf_);
00036 
00037   this->addEvaluatedField(residual_);
00038 
00039   this->setName("CMResidualCoarse" + PHX::TypeString<EvalT>::value);
00040 
00041   std::vector<PHX::DataLayout::size_type>
00042   dims;
00043 
00044   w_grad_bf_.fieldTag().dataLayout().dimensions(dims);
00045   num_nodes_ = dims[1];
00046   num_pts_ = dims[2];
00047   num_dims_ = dims[3];
00048 
00049   Teuchos::RCP<ParamLib>
00050   paramLib = p.get<Teuchos::RCP<ParamLib> >("Parameter Library");
00051 }
00052 
00053 //
00054 //
00055 //
00056 template<typename EvalT, typename Traits>
00057 void CMResidualCoarse<EvalT, Traits>::
00058 postRegistrationSetup(typename Traits::SetupData d,
00059     PHX::FieldManager<Traits>& fm)
00060 {
00061   this->utils.setFieldData(stress_, fm);
00062   this->utils.setFieldData(def_grad_, fm);
00063   this->utils.setFieldData(w_grad_bf_, fm);
00064   this->utils.setFieldData(w_bf_, fm);
00065   this->utils.setFieldData(residual_, fm);
00066   if (have_body_force_ == true) {
00067     this->utils.setFieldData(body_force_, fm);
00068   }
00069 }
00070 
00071 //
00072 //
00073 //
00074 template<typename EvalT, typename Traits>
00075 void CMResidualCoarse<EvalT, Traits>::
00076 evaluateFields(typename Traits::EvalData workset)
00077 {
00078   Intrepid::Tensor<ScalarT>
00079   F(num_dims_), P(num_dims_), sig(num_dims_);
00080 
00081   Intrepid::Tensor<ScalarT>
00082   I(Intrepid::eye<ScalarT>(num_dims_));
00083 
00084   // initialize residual
00085   for (size_t cell = 0; cell < workset.numCells; ++cell) {
00086      for (size_t node = 0; node < num_nodes_; ++node) {
00087        for (size_t dim = 0; dim < num_dims_; ++dim) {
00088          residual_(cell, node, dim) = 0.0;
00089        }
00090      }
00091      for (size_t pt = 0; pt < num_pts_; ++pt) {
00092        F.fill(&def_grad_(cell, pt, 0, 0));
00093        sig.fill(&stress_(cell, pt, 0, 0));
00094 
00095        // map Cauchy stress to 1st PK
00096        P = Intrepid::piola(F, sig);
00097 
00098        for (size_t node = 0; node < num_nodes_; ++node) {
00099          for (size_t i = 0; i < num_dims_; ++i) {
00100            for (size_t j = 0; j < num_dims_; ++j) {
00101              residual_(cell, node, i) +=
00102                  P(i, j) * w_grad_bf_(cell, node, pt, j);
00103            }
00104          }
00105        }
00106      }
00107    }
00108 
00109   // optional body force
00110   if (have_body_force_) {
00111     for (size_t cell = 0; cell < workset.numCells; ++cell) {
00112       for (size_t node = 0; node < num_nodes_; ++node) {
00113         for (size_t pt = 0; pt < num_pts_; ++pt) {
00114           for (size_t dim = 0; dim < num_dims_; ++dim) {
00115             residual_(cell, node, dim) +=
00116                 w_bf_(cell, node, pt) * body_force_(cell, pt, dim);
00117           }
00118         }
00119       }
00120     }
00121   }
00122 }
00123 
00124 }
00125 

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