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

GradientDamageProblem.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 #ifndef GRADIENTDAMAGEPROBLEM_HPP
00008 #define GRADIENTDAMAGEPROBLEM_HPP
00009 
00010 #include "Teuchos_RCP.hpp"
00011 #include "Teuchos_ParameterList.hpp"
00012 
00013 #include "Albany_AbstractProblem.hpp"
00014 
00015 #include "Phalanx.hpp"
00016 #include "PHAL_Workset.hpp"
00017 #include "PHAL_Dimension.hpp"
00018 #include "PHAL_AlbanyTraits.hpp"
00019 
00020 namespace Albany {
00021 
00022   class GradientDamageProblem : public Albany::AbstractProblem {
00023   public:
00024   
00026     GradientDamageProblem( const Teuchos::RCP<Teuchos::ParameterList>& params,
00027                            const Teuchos::RCP<ParamLib>& paramLib,
00028                            const int numEq);
00029 
00031     virtual ~GradientDamageProblem();
00032 
00034     virtual int spatialDimension() const { return numDim; }
00035 
00037     virtual void buildProblem(
00038       Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >  meshSpecs,
00039       StateManager& stateMgr);
00040 
00041     // Build evaluators
00042     virtual Teuchos::Array< Teuchos::RCP<const PHX::FieldTag> >
00043     buildEvaluators(
00044       PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00045       const Albany::MeshSpecsStruct& meshSpecs,
00046       Albany::StateManager& stateMgr,
00047       Albany::FieldManagerChoice fmchoice,
00048       const Teuchos::RCP<Teuchos::ParameterList>& responseList);
00049 
00051     Teuchos::RCP<const Teuchos::ParameterList> getValidProblemParameters() const;
00052 
00053     void getAllocatedStates(
00054          Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > oldState_,
00055          Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > newState_
00056          ) const;
00057 
00058   private:
00059 
00061     GradientDamageProblem(const GradientDamageProblem&);
00062     
00064     GradientDamageProblem& operator=(const GradientDamageProblem&);
00065 
00066   public:
00067 
00069     template <typename EvalT> 
00070     Teuchos::RCP<const PHX::FieldTag>
00071     constructEvaluators(
00072       PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00073       const Albany::MeshSpecsStruct& meshSpecs,
00074       Albany::StateManager& stateMgr,
00075       Albany::FieldManagerChoice fmchoice,
00076       const Teuchos::RCP<Teuchos::ParameterList>& responseList);
00077 
00078     void constructDirichletEvaluators(const Albany::MeshSpecsStruct& meshSpecs);
00079 
00080   protected:
00081 
00083     bool haveSource;
00084 
00085     // counting helpers
00086     int numQPts;
00087     int numNodes;
00088     int numVertices;
00089 
00090     int D_offset;  //Position of T unknown in nodal DOFs
00091     int X_offset;  //Position of X unknown in nodal DOFs, followed by Y,Z
00092     int numDim;    //Number of spatial dimensions and displacement variable 
00093 
00094     // string to store material model name
00095     std::string matModel;
00096 
00097     // state containers
00098     Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > oldState;
00099     Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > newState;
00100   };
00101 
00102 }
00103 
00104 #include "Albany_Utils.hpp"
00105 #include "Albany_ProblemUtils.hpp"
00106 #include "Albany_ResponseUtilities.hpp"
00107 #include "Albany_EvaluatorUtils.hpp"
00108 
00109 #include "BulkModulus.hpp"
00110 #include "ShearModulus.hpp"
00111 #include "PHAL_Source.hpp"
00112 #include "DefGrad.hpp"
00113 #include "HardeningModulus.hpp"
00114 #include "YieldStrength.hpp"
00115 #include "SaturationModulus.hpp"
00116 #include "SaturationExponent.hpp"
00117 #include "PHAL_SaveStateField.hpp"
00118 #include "TLElasResid.hpp"
00119 #include "J2Damage.hpp"
00120 #include "DamageLS.hpp"
00121 #include "DamageSource.hpp"
00122 #include "DamageResid.hpp"
00123 
00124 template <typename EvalT>
00125 Teuchos::RCP<const PHX::FieldTag>
00126 Albany::GradientDamageProblem::constructEvaluators(
00127   PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00128   const Albany::MeshSpecsStruct& meshSpecs,
00129   Albany::StateManager& stateMgr,
00130   Albany::FieldManagerChoice fieldManagerChoice,
00131   const Teuchos::RCP<Teuchos::ParameterList>& responseList)
00132 {
00133    using Teuchos::RCP;
00134    using Teuchos::rcp;
00135    using Teuchos::ParameterList;
00136    using PHX::DataLayout;
00137    using PHX::MDALayout;
00138    using std::vector;
00139    using std::map;
00140    using PHAL::AlbanyTraits;
00141 
00142    // get the name of the current element block
00143    std::string elementBlockName = meshSpecs.ebName;
00144 
00145    RCP<shards::CellTopology> cellType = rcp(new shards::CellTopology (&meshSpecs.ctd));
00146    RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > >
00147      intrepidBasis = Albany::getIntrepidBasis(meshSpecs.ctd);
00148 
00149    numNodes = intrepidBasis->getCardinality();
00150    const int worksetSize = meshSpecs.worksetSize;
00151 
00152    Intrepid::DefaultCubatureFactory<RealType> cubFactory;
00153    RCP <Intrepid::Cubature<RealType> > cubature = cubFactory.create(*cellType, meshSpecs.cubatureDegree);
00154 
00155    const int numQPts = cubature->getNumPoints();
00156    const int numVertices = cellType->getNodeCount();
00157 
00158    *out << "Field Dimensions: Workset=" << worksetSize 
00159         << ", Vertices= " << numVertices
00160         << ", Nodes= " << numNodes
00161         << ", QuadPts= " << numQPts
00162         << ", Dim= " << numDim << std::endl;
00163 
00164 
00165    // Construct standard FEM evaluators with standard field names                              
00166    RCP<Albany::Layouts> dl = rcp(new Albany::Layouts(worksetSize,numVertices,numNodes,numQPts,numDim));
00167    TEUCHOS_TEST_FOR_EXCEPTION(dl->vectorAndGradientLayoutsAreEquivalent==false, std::logic_error,
00168                               "Data Layout Usage in Mechanics problems assume vecDim = numDim");
00169    Albany::EvaluatorUtils<EvalT, PHAL::AlbanyTraits> evalUtils(dl);
00170    std::string scatterName="Scatter Damage";
00171 
00172    // Displacement Variable
00173    Teuchos::ArrayRCP<std::string> dof_names(1);
00174    dof_names[0] = "Displacement";
00175    Teuchos::ArrayRCP<std::string> resid_names(1);
00176    resid_names[0] = "Mechanical Residual";
00177 
00178    fm0.template registerEvaluator<EvalT>
00179      (evalUtils.constructDOFVecInterpolationEvaluator(dof_names[0], X_offset));
00180 
00181    fm0.template registerEvaluator<EvalT>
00182      (evalUtils.constructDOFVecGradInterpolationEvaluator(dof_names[0], X_offset));
00183 
00184    fm0.template registerEvaluator<EvalT>
00185      (evalUtils.constructGatherSolutionEvaluator_noTransient(true, dof_names, X_offset)); 
00186    fm0.template registerEvaluator<EvalT>
00187      (evalUtils.constructScatterResidualEvaluator(true, resid_names, X_offset));
00188 
00189    // Damage Variable
00190    Teuchos::ArrayRCP<std::string> ddof_names(1);
00191    ddof_names[0] = "Damage";
00192    Teuchos::ArrayRCP<std::string> ddof_names_dot(1);
00193    ddof_names_dot[0] = ddof_names[0]+"_dot";
00194    Teuchos::ArrayRCP<std::string> dresid_names(1);
00195    dresid_names[0] = ddof_names[0]+" Residual";
00196 
00197    fm0.template registerEvaluator<EvalT>
00198      (evalUtils.constructDOFInterpolationEvaluator(ddof_names[0], D_offset));
00199 
00200    fm0.template registerEvaluator<EvalT>
00201      (evalUtils.constructDOFInterpolationEvaluator(ddof_names_dot[0], D_offset));
00202 
00203    fm0.template registerEvaluator<EvalT>
00204      (evalUtils.constructDOFGradInterpolationEvaluator(ddof_names[0], D_offset));
00205 
00206    fm0.template registerEvaluator<EvalT>
00207      (evalUtils.constructGatherSolutionEvaluator(false, ddof_names, ddof_names_dot, D_offset));
00208 
00209    fm0.template registerEvaluator<EvalT>
00210      (evalUtils.constructScatterResidualEvaluator(false, dresid_names, D_offset, scatterName));
00211 
00212    // General FEM Stuff
00213    fm0.template registerEvaluator<EvalT>
00214      (evalUtils.constructGatherCoordinateVectorEvaluator());
00215 
00216    fm0.template registerEvaluator<EvalT>
00217      (evalUtils.constructMapToPhysicalFrameEvaluator(cellType, cubature));
00218 
00219    fm0.template registerEvaluator<EvalT>
00220      (evalUtils.constructComputeBasisFunctionsEvaluator(cellType, intrepidBasis, cubature));
00221 
00222    // Temporary variable used numerous times below
00223    Teuchos::RCP<PHX::Evaluator<AlbanyTraits> > ev;
00224 
00225    { // Bulk Modulus
00226      RCP<ParameterList> p = rcp(new ParameterList);
00227 
00228      p->set<std::string>("QP Variable Name", "Bulk Modulus");
00229      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00230      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00231      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00232      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00233 
00234      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00235      Teuchos::ParameterList& paramList = params->sublist("Bulk Modulus");
00236      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00237 
00238      ev = rcp(new LCM::BulkModulus<EvalT,AlbanyTraits>(*p));
00239      fm0.template registerEvaluator<EvalT>(ev);
00240    }
00241 
00242    { // Shear Modulus 
00243      RCP<ParameterList> p = rcp(new ParameterList);
00244 
00245      p->set<std::string>("QP Variable Name", "Shear Modulus");
00246      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00247      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00248      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00249      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00250 
00251      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00252      Teuchos::ParameterList& paramList = params->sublist("Shear Modulus");
00253      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00254 
00255      ev = rcp(new LCM::ShearModulus<EvalT,AlbanyTraits>(*p));
00256      fm0.template registerEvaluator<EvalT>(ev);
00257    }
00258 
00259    if (haveSource) { // Source
00260      TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00261                                 "Error!  Sources not implemented in Elasticity yet!");
00262 
00263      RCP<ParameterList> p = rcp(new ParameterList);
00264 
00265      p->set<std::string>("Source Name", "Source");
00266      p->set<std::string>("Variable Name", "Displacement");
00267      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00268 
00269      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00270      Teuchos::ParameterList& paramList = params->sublist("Source Functions");
00271      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00272 
00273      ev = rcp(new PHAL::Source<EvalT,AlbanyTraits>(*p));
00274      fm0.template registerEvaluator<EvalT>(ev);
00275    }
00276 
00277    { // Deformation Gradient
00278      RCP<ParameterList> p = rcp(new ParameterList("DefGrad"));
00279 
00280      //Inputs: flags, weights, GradU
00281      const bool avgJ = params->get("avgJ", false);
00282      p->set<bool>("avgJ Name", avgJ);
00283      const bool volavgJ = params->get("volavgJ", false);
00284      p->set<bool>("volavgJ Name", volavgJ);
00285      const bool weighted_Volume_Averaged_J = params->get("weighted_Volume_Averaged_J", false);
00286      p->set<bool>("weighted_Volume_Averaged_J Name", weighted_Volume_Averaged_J);
00287      p->set<std::string>("Weights Name","Weights");
00288      p->set<std::string>("Gradient QP Variable Name", "Displacement Gradient");
00289      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00290 
00291      //Outputs: F, J
00292      p->set<std::string>("DefGrad Name", "Deformation Gradient"); //dl->qp_tensor also
00293      p->set<std::string>("DetDefGrad Name", "Determinant of Deformation Gradient"); 
00294      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00295 
00296      ev = rcp(new LCM::DefGrad<EvalT,AlbanyTraits>(*p));
00297      fm0.template registerEvaluator<EvalT>(ev);
00298    }
00299    { // Hardening Modulus
00300      RCP<ParameterList> p = rcp(new ParameterList);
00301 
00302      p->set<std::string>("QP Variable Name", "Hardening Modulus");
00303      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00304      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00305      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00306      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00307 
00308      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00309      Teuchos::ParameterList& paramList = params->sublist("Hardening Modulus");
00310      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00311 
00312      ev = rcp(new LCM::HardeningModulus<EvalT,AlbanyTraits>(*p));
00313      fm0.template registerEvaluator<EvalT>(ev);
00314    }
00315 
00316    { // Yield Strength
00317      RCP<ParameterList> p = rcp(new ParameterList);
00318 
00319      p->set<std::string>("QP Variable Name", "Yield Strength");
00320      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00321      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00322      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00323      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00324 
00325      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00326      Teuchos::ParameterList& paramList = params->sublist("Yield Strength");
00327      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00328 
00329      ev = rcp(new LCM::YieldStrength<EvalT,AlbanyTraits>(*p));
00330      fm0.template registerEvaluator<EvalT>(ev);
00331    }
00332 
00333    { // Saturation Modulus
00334      RCP<ParameterList> p = rcp(new ParameterList);
00335 
00336      p->set<std::string>("Saturation Modulus Name", "Saturation Modulus");
00337      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00338      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00339      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00340      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00341 
00342      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00343      Teuchos::ParameterList& paramList = params->sublist("Saturation Modulus");
00344      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00345 
00346      ev = rcp(new LCM::SaturationModulus<EvalT,AlbanyTraits>(*p));
00347      fm0.template registerEvaluator<EvalT>(ev);
00348    }
00349 
00350    { // Saturation Exponent
00351      RCP<ParameterList> p = rcp(new ParameterList);
00352 
00353      p->set<std::string>("Saturation Exponent Name", "Saturation Exponent");
00354      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00355      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00356      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00357      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00358 
00359      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00360      Teuchos::ParameterList& paramList = params->sublist("Saturation Exponent");
00361      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00362 
00363      ev = rcp(new LCM::SaturationExponent<EvalT,AlbanyTraits>(*p));
00364      fm0.template registerEvaluator<EvalT>(ev);
00365    }
00366 
00367    {// Stress
00368      RCP<ParameterList> p = rcp(new ParameterList("Stress"));
00369 
00370      //Input
00371      p->set<std::string>("DefGrad Name", "Deformation Gradient");
00372      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00373 
00374      p->set<std::string>("Bulk Modulus Name", "Bulk Modulus");
00375      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00376 
00377      p->set<std::string>("Shear Modulus Name", "Shear Modulus");  // dl->qp_scalar also
00378      p->set<std::string>("Hardening Modulus Name", "Hardening Modulus"); // dl->qp_scalar also
00379      p->set<std::string>("Yield Strength Name", "Yield Strength"); // dl->qp_scalar also
00380      p->set<std::string>("Saturation Modulus Name", "Saturation Modulus"); // dl->qp_scalar also
00381      p->set<std::string>("Saturation Exponent Name", "Saturation Exponent"); // dl->qp_scalar also
00382      p->set<std::string>("DetDefGrad Name", "Determinant of Deformation Gradient");  // dl->qp_scalar also
00383      p->set<std::string>("Damage Name", "Damage");
00384 
00385      //Output
00386      p->set<std::string>("Stress Name", "Stress"); //dl->qp_tensor also
00387      p->set<std::string>("DP Name", "DP"); // dl->qp_scalar also
00388      p->set<std::string>("Effective Stress Name", "Effective Stress"); // dl->qp_scalar also
00389      p->set<std::string>("Energy Name", "Energy"); // dl->qp_scalar also
00390 
00391      p->set<std::string>("Fp Name", "Fp");  // dl->qp_tensor also
00392      p->set<std::string>("Eqps Name", "eqps");  // dl->qp_scalar also
00393 
00394  
00395      //Declare what state data will need to be saved (name, layout, init_type)
00396      // A :true: as 5th argument declares that the previous state needs to be saved
00397 
00398      ev = rcp(new LCM::J2Damage<EvalT,AlbanyTraits>(*p));
00399      fm0.template registerEvaluator<EvalT>(ev);
00400      p = stateMgr.registerStateVariable("Stress",dl->qp_tensor, dl->dummy, elementBlockName, "scalar", 0.0);
00401      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00402      fm0.template registerEvaluator<EvalT>(ev);
00403      p = stateMgr.registerStateVariable("Fp",dl->qp_tensor, dl->dummy, elementBlockName, "identity",1.0,true);
00404      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00405      fm0.template registerEvaluator<EvalT>(ev);
00406      p = stateMgr.registerStateVariable("eqps",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0,true);
00407      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00408      fm0.template registerEvaluator<EvalT>(ev);
00409    }
00410 
00411    { // Displacement Resid
00412      RCP<ParameterList> p = rcp(new ParameterList("Mechanical Residual"));
00413 
00414      //Input
00415      p->set<std::string>("Stress Name", "Stress");
00416      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00417 
00418      p->set<std::string>("DefGrad Name", "Deformation Gradient"); //dl->qp_tensor also
00419 
00420      p->set<std::string>("DetDefGrad Name", "Determinant of Deformation Gradient");
00421      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00422 
00423      p->set<std::string>("Weighted Gradient BF Name", "wGrad BF");
00424      p->set< RCP<DataLayout> >("Node QP Vector Data Layout", dl->node_qp_vector);
00425 
00426      p->set<std::string>("Weighted BF Name", "wBF");
00427      p->set< RCP<DataLayout> >("Node QP Scalar Data Layout", dl->node_qp_scalar);
00428      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00429 
00430      //Output
00431      p->set<std::string>("Residual Name", "Mechanical Residual");
00432      p->set< RCP<DataLayout> >("Node Vector Data Layout", dl->node_vector);
00433 
00434      ev = rcp(new LCM::TLElasResid<EvalT,AlbanyTraits>(*p));
00435      fm0.template registerEvaluator<EvalT>(ev);
00436    }
00437 
00438    { // Damage length scale
00439      RCP<ParameterList> p = rcp(new ParameterList);
00440 
00441      p->set<std::string>("QP Variable Name", "Damage Length Scale");
00442      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00443      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00444      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00445      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00446 
00447      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00448      Teuchos::ParameterList& paramList = params->sublist("Damage Length Scale");
00449      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00450 
00451      ev = rcp(new LCM::DamageLS<EvalT,AlbanyTraits>(*p));
00452      fm0.template registerEvaluator<EvalT>(ev);
00453    }
00454 
00455    { // Damage Source
00456      RCP<ParameterList> p = rcp(new ParameterList("Damage Source"));
00457 
00458      //Input
00459      RealType gc = params->get("gc", 1.0);
00460      p->set<RealType>("gc Name", gc);
00461      p->set<std::string>("Bulk Modulus Name", "Bulk Modulus");
00462      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00463      p->set<std::string>("Damage Name", "Damage");
00464      p->set<std::string>("DP Name", "DP");
00465      p->set<std::string>("Effective Stress Name", "Effective Stress");
00466      p->set<std::string>("Energy Name", "Energy");
00467      p->set<std::string>("DetDefGrad Name", "Determinant of Deformation Gradient");
00468      p->set<std::string>("Damage Length Scale Name", "Damage Length Scale");
00469 
00470      //Output
00471      p->set<std::string>("Damage Source Name", "Damage Source");
00472 
00473      ev = rcp(new LCM::DamageSource<EvalT,AlbanyTraits>(*p));
00474      fm0.template registerEvaluator<EvalT>(ev);
00475      p = stateMgr.registerStateVariable("Damage Source",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0,true);
00476      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00477      fm0.template registerEvaluator<EvalT>(ev);
00478      p = stateMgr.registerStateVariable("Damage",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0,true);
00479      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00480      fm0.template registerEvaluator<EvalT>(ev);
00481    }
00482 
00483    { // Damage Resid
00484      RCP<ParameterList> p = rcp(new ParameterList("Damage Resid"));
00485 
00486      //Input
00487      RealType gc = params->get("gc", 0.0);
00488      p->set<RealType>("gc Name", gc);
00489      p->set<std::string>("Weighted BF Name", "wBF");
00490      p->set< RCP<DataLayout> >("Node QP Scalar Data Layout", dl->node_qp_scalar);
00491      p->set<std::string>("QP Variable Name", "Damage");
00492 
00493      p->set<std::string>("QP Time Derivative Variable Name", "Damage_dot");
00494 
00495      p->set<std::string>("Damage Source Name", "Damage Source");  //dl->qp_scalar
00496 
00497      p->set<std::string>("Damage Length Scale Name", "Damage Length Scale");
00498      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00499 
00500      p->set<std::string>("Gradient QP Variable Name", "Damage Gradient");
00501      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00502 
00503      p->set<std::string>("Weighted Gradient BF Name", "wGrad BF");
00504      p->set< RCP<DataLayout> >("Node QP Vector Data Layout", dl->node_qp_vector);
00505 
00506      //Output
00507      p->set<std::string>("Residual Name", "Damage Residual");
00508      p->set< RCP<DataLayout> >("Node Scalar Data Layout", dl->node_scalar);
00509 
00510      ev = rcp(new LCM::DamageResid<EvalT,AlbanyTraits>(*p));
00511      fm0.template registerEvaluator<EvalT>(ev);
00512    }
00513 
00514    if (fieldManagerChoice == Albany::BUILD_RESID_FM)  {
00515      PHX::Tag<typename EvalT::ScalarT> res_tag("Scatter", dl->dummy);
00516      fm0.requireField<EvalT>(res_tag);
00517 
00518      PHX::Tag<typename EvalT::ScalarT> res_tag2(scatterName, dl->dummy);
00519      fm0.requireField<EvalT>(res_tag2);
00520 
00521      return res_tag.clone();
00522    }
00523    else if (fieldManagerChoice == Albany::BUILD_RESPONSE_FM) {
00524      Albany::ResponseUtilities<EvalT, PHAL::AlbanyTraits> respUtils(dl);
00525      return respUtils.constructResponses(fm0, *responseList, stateMgr);
00526    }
00527 
00528    return Teuchos::null;
00529 }
00530 
00531 #endif // ALBANY_GRADIENTDAMAGEPROBLEM_HPP

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