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

GradientDamageProblem.cpp

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 #include "GradientDamageProblem.hpp"
00007 #include "Albany_Utils.hpp"
00008 #include "Albany_ProblemUtils.hpp"
00009 
00010 Albany::GradientDamageProblem::
00011 GradientDamageProblem(
00012           const Teuchos::RCP<Teuchos::ParameterList>& params_,
00013           const Teuchos::RCP<ParamLib>& paramLib_,
00014           const int numDim_) :
00015   Albany::AbstractProblem(params_, paramLib_, numDim_ + 1),
00016   haveSource(false),
00017   numDim(numDim_)
00018 {
00019  
00020   std::string& method = params->get("Name", "GradientDamage ");
00021   *out << "Problem Name = " << method << std::endl;
00022   
00023   haveSource =  params->isSublist("Source Functions");
00024 
00025 // Changing this ifdef changes ordering from  (X,Y,D) to (D,X,Y)
00026 //#define NUMBER_D_FIRST
00027 #ifdef NUMBER_D_FIRST
00028   D_offset=0;
00029   X_offset=1;
00030 #else
00031   X_offset=0;
00032   D_offset=numDim;
00033 #endif
00034 
00035 // the following function returns the problem information required for setting the rigid body modes (RBMs) for elasticity problems
00036 //written by IK, Feb. 2012
00037 
00038   int numScalar = 1;
00039   int nullSpaceDim = 0;
00040   if (numDim == 1) {nullSpaceDim = 0; }
00041   else {
00042     if (numDim == 2) {nullSpaceDim = 3; }
00043     if (numDim == 3) {nullSpaceDim = 6; }
00044   }
00045 
00046   rigidBodyModes->setParameters(numDim + 1, numDim, numScalar, nullSpaceDim);
00047 
00048 }
00049 
00050 Albany::GradientDamageProblem::
00051 ~GradientDamageProblem()
00052 {
00053 }
00054 
00055 void
00056 Albany::GradientDamageProblem::
00057 buildProblem(
00058   Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >  meshSpecs,
00059   Albany::StateManager& stateMgr)
00060 {
00061   /* Construct All Phalanx Evaluators */
00062   TEUCHOS_TEST_FOR_EXCEPTION(meshSpecs.size()!=1,std::logic_error,"Problem supports one Material Block");
00063   fm.resize(1);
00064   fm[0]  = Teuchos::rcp(new PHX::FieldManager<PHAL::AlbanyTraits>);
00065   buildEvaluators(*fm[0], *meshSpecs[0], stateMgr, BUILD_RESID_FM, 
00066       Teuchos::null);
00067   constructDirichletEvaluators(*meshSpecs[0]);
00068 }
00069 
00070 Teuchos::Array<Teuchos::RCP<const PHX::FieldTag> >
00071 Albany::GradientDamageProblem::
00072 buildEvaluators(
00073   PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00074   const Albany::MeshSpecsStruct& meshSpecs,
00075   Albany::StateManager& stateMgr,
00076   Albany::FieldManagerChoice fmchoice,
00077   const Teuchos::RCP<Teuchos::ParameterList>& responseList)
00078 {
00079   // Call constructeEvaluators<EvalT>(*rfm[0], *meshSpecs[0], stateMgr);
00080   // for each EvalT in PHAL::AlbanyTraits::BEvalTypes
00081   ConstructEvaluatorsOp<GradientDamageProblem> op(
00082     *this, fm0, meshSpecs, stateMgr, fmchoice, responseList);
00083   boost::mpl::for_each<PHAL::AlbanyTraits::BEvalTypes>(op);
00084   return *op.tags;
00085 }
00086 
00087 void
00088 Albany::GradientDamageProblem::constructDirichletEvaluators(
00089         const Albany::MeshSpecsStruct& meshSpecs)
00090 {
00091   // Construct Dirichlet evaluators for all nodesets and names
00092   std::vector<std::string> dirichletNames(neq);
00093   dirichletNames[X_offset] = "X";
00094   if (numDim>1) dirichletNames[X_offset+1] = "Y";
00095   if (numDim>2) dirichletNames[X_offset+2] = "Z";
00096   dirichletNames[D_offset] = "D";
00097   Albany::BCUtils<Albany::DirichletTraits> dirUtils;
00098   dfm = dirUtils.constructBCEvaluators(meshSpecs.nsNames, dirichletNames,
00099                                        this->params, this->paramLib);
00100 }
00101 
00102 Teuchos::RCP<const Teuchos::ParameterList>
00103 Albany::GradientDamageProblem::getValidProblemParameters() const
00104 {
00105   Teuchos::RCP<Teuchos::ParameterList> validPL =
00106     this->getGenericProblemParams("ValidGradientDamageProblemParams");
00107 
00108   validPL->sublist("Bulk Modulus", false, "");
00109   validPL->sublist("Shear Modulus", false, "");
00110   validPL->sublist("Poissons Ratio", false, "");
00111   validPL->sublist("Damage Length Scale", false, "");
00112   validPL->sublist("Hardening Modulus", false, "");
00113   validPL->sublist("Yield Strength", false, "");
00114   validPL->sublist("Saturation Modulus", false, "");
00115   validPL->sublist("Saturation Exponent", false, "");
00116   validPL->set<double>("gc", false, "");
00117   validPL->set<bool>("avgJ", false, "Flag to indicate the J should be averaged");
00118   validPL->set<bool>("volavgJ", false, "Flag to indicate the J should be volume averaged");
00119   validPL->set<bool>("Use Composite Tet 10", false, "Flag to use the compostie tet 10 basis in Intrepid");
00120 
00121   return validPL;
00122 }
00123 
00124 void
00125 Albany::GradientDamageProblem::getAllocatedStates(
00126    Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > oldState_,
00127    Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > newState_
00128    ) const
00129 {
00130   oldState_ = oldState;
00131   newState_ = newState;
00132 }

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