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

ThermoPoroPlasticityProblem.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 THERMOPOROPLASTICITYPROBLEM_HPP
00008 #define THERMOPOROPLASTICITYPROBLEM_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 "Albany_ProblemUtils.hpp"
00019 #include "Albany_ResponseUtilities.hpp"
00020 #include "Albany_EvaluatorUtils.hpp"
00021 #include "PHAL_AlbanyTraits.hpp"
00022 
00023 namespace Albany {
00024 
00028   class ThermoPoroPlasticityProblem : public Albany::AbstractProblem {
00029   public:
00030   
00032     ThermoPoroPlasticityProblem(const Teuchos::RCP<Teuchos::ParameterList>& params,
00033         const Teuchos::RCP<ParamLib>& paramLib,
00034         const int numEq);
00035 
00037     virtual ~ThermoPoroPlasticityProblem();
00038 
00040     virtual int spatialDimension() const { return numDim; }
00041 
00043     virtual void buildProblem(
00044       Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> >  meshSpecs,
00045       StateManager& stateMgr);
00046 
00047     // Build evaluators
00048     virtual Teuchos::Array< Teuchos::RCP<const PHX::FieldTag> >
00049     buildEvaluators(
00050       PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00051       const Albany::MeshSpecsStruct& meshSpecs,
00052       Albany::StateManager& stateMgr,
00053       Albany::FieldManagerChoice fmchoice,
00054       const Teuchos::RCP<Teuchos::ParameterList>& responseList);
00055 
00057     Teuchos::RCP<const Teuchos::ParameterList> getValidProblemParameters() const;
00058 
00059     void getAllocatedStates(
00060          Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > oldState_,
00061    Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > newState_
00062    ) const;
00063 
00064   private:
00065 
00067     ThermoPoroPlasticityProblem(const ThermoPoroPlasticityProblem&);
00068     
00070     ThermoPoroPlasticityProblem& operator=(const ThermoPoroPlasticityProblem&);
00071 
00072   public:
00073 
00075     template <typename EvalT> 
00076     Teuchos::RCP<const PHX::FieldTag>
00077     constructEvaluators(
00078       PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00079       const Albany::MeshSpecsStruct& meshSpecs,
00080       Albany::StateManager& stateMgr,
00081       Albany::FieldManagerChoice fmchoice,
00082       const Teuchos::RCP<Teuchos::ParameterList>& responseList);
00083 
00084     void constructDirichletEvaluators(const Albany::MeshSpecsStruct& meshSpecs);
00085 
00086   protected:
00087 
00089     bool haveSource;
00090     int T_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 TEMP_offset; // Position of TEMP unknown in nodal DOFs
00093     int numDim;    //Number of spatial dimensions and displacement variable 
00094 
00095     std::string matModel;
00096 
00097     Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > oldState;
00098     Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<Intrepid::FieldContainer<RealType> > > > newState;
00099   };
00100 }
00101 
00102 #include "Teuchos_RCP.hpp"
00103 #include "Teuchos_ParameterList.hpp"
00104 
00105 #include "Albany_AbstractProblem.hpp"
00106 
00107 #include "Phalanx.hpp"
00108 #include "PHAL_Workset.hpp"
00109 #include "PHAL_Dimension.hpp"
00110 #include "PHAL_AlbanyTraits.hpp"
00111 #include "Albany_Utils.hpp"
00112 #include "Albany_ProblemUtils.hpp"
00113 #include "Albany_ResponseUtilities.hpp"
00114 #include "Albany_EvaluatorUtils.hpp"
00115 
00116 #include "Time.hpp"
00117 // #include "Strain.hpp"
00118 #include "DefGrad.hpp"
00119 #include "PHAL_SaveStateField.hpp"
00120 #include "Porosity.hpp"
00121 #include "BiotCoefficient.hpp"
00122 #include "BiotModulus.hpp"
00123 #include "PHAL_ThermalConductivity.hpp"
00124 #include "KCPermeability.hpp"
00125 #include "ElasticModulus.hpp"
00126 #include "PoissonsRatio.hpp"
00127 #include "ShearModulus.hpp"
00128 #include "BulkModulus.hpp"
00129 // #include "TotalStress.hpp"
00130 #include "PHAL_Source.hpp"
00131 #include "ThermoPoroPlasticityResidMass.hpp"
00132 #include "ThermoPoroPlasticityResidMomentum.hpp"
00133 #include "ThermoPoroPlasticityResidEnergy.hpp"
00134 #include "PHAL_NSMaterialProperty.hpp"
00135 
00136 #include "MixtureThermalExpansion.hpp"
00137 #include "MixtureSpecificHeat.hpp"
00138 
00139 #include "J2Stress.hpp"
00140 #include "Neohookean.hpp"
00141 #include "PisdWdF.hpp"
00142 #include "HardeningModulus.hpp"
00143 #include "YieldStrength.hpp"
00144 #include "SaturationModulus.hpp"
00145 #include "SaturationExponent.hpp"
00146 #include "DislocationDensity.hpp"
00147 #include "TLPoroStress.hpp"
00148 
00149 
00150 template <typename EvalT>
00151 Teuchos::RCP<const PHX::FieldTag>
00152 Albany::ThermoPoroPlasticityProblem::constructEvaluators(
00153   PHX::FieldManager<PHAL::AlbanyTraits>& fm0,
00154   const Albany::MeshSpecsStruct& meshSpecs,
00155   Albany::StateManager& stateMgr,
00156   Albany::FieldManagerChoice fieldManagerChoice,
00157   const Teuchos::RCP<Teuchos::ParameterList>& responseList)
00158 {
00159    using Teuchos::RCP;
00160    using Teuchos::rcp;
00161    using Teuchos::ParameterList;
00162    using PHX::DataLayout;
00163    using PHX::MDALayout;
00164    using std::vector;
00165    using PHAL::AlbanyTraits;
00166 
00167    // get the name of the current element block
00168    std::string elementBlockName = meshSpecs.ebName;
00169 
00170 
00171    RCP<shards::CellTopology> cellType = rcp(new shards::CellTopology (&meshSpecs.ctd));
00172    RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > >
00173      intrepidBasis = Albany::getIntrepidBasis(meshSpecs.ctd);
00174 
00175    const int numNodes = intrepidBasis->getCardinality();
00176    const int worksetSize = meshSpecs.worksetSize;
00177 
00178    Intrepid::DefaultCubatureFactory<RealType> cubFactory;
00179    RCP <Intrepid::Cubature<RealType> > cubature = cubFactory.create(*cellType, meshSpecs.cubatureDegree);
00180 
00181    const int numQPts = cubature->getNumPoints();
00182    const int numVertices = cellType->getNodeCount();
00183 
00184    *out << "Field Dimensions: Workset=" << worksetSize 
00185         << ", Vertices= " << numVertices
00186         << ", Nodes= " << numNodes
00187         << ", QuadPts= " << numQPts
00188         << ", Dim= " << numDim << std::endl;
00189 
00190 
00191    // Construct standard FEM evaluators with standard field names                              
00192    RCP<Albany::Layouts> dl = rcp(new Albany::Layouts(worksetSize,numVertices,numNodes,numQPts,numDim));
00193    TEUCHOS_TEST_FOR_EXCEPTION(dl->vectorAndGradientLayoutsAreEquivalent==false, std::logic_error,
00194                               "Data Layout Usage in Mechanics problems assume vecDim = numDim");
00195    Albany::EvaluatorUtils<EvalT, PHAL::AlbanyTraits> evalUtils(dl);
00196    std::string scatterName="Scatter PoreFluid";
00197 
00198 
00199    // ----------------------setup the solution field ---------------//
00200 
00201    // Displacement Variable
00202    Teuchos::ArrayRCP<std::string> dof_names(1);
00203    dof_names[0] = "Displacement";
00204    Teuchos::ArrayRCP<std::string> resid_names(1);
00205    resid_names[0] = dof_names[0]+" Residual";
00206 
00207    fm0.template registerEvaluator<EvalT>
00208      (evalUtils.constructDOFVecInterpolationEvaluator(dof_names[0], X_offset));
00209 
00210    fm0.template registerEvaluator<EvalT>
00211      (evalUtils.constructDOFVecGradInterpolationEvaluator(dof_names[0], X_offset));
00212 
00213    fm0.template registerEvaluator<EvalT>
00214      (evalUtils.constructGatherSolutionEvaluator_noTransient(true, dof_names, X_offset));
00215 
00216    fm0.template registerEvaluator<EvalT>
00217      (evalUtils.constructScatterResidualEvaluator(true, resid_names, X_offset));
00218 
00219    // Pore Pressure Variable
00220    Teuchos::ArrayRCP<std::string> tdof_names(1);
00221    tdof_names[0] = "Pore Pressure";
00222    Teuchos::ArrayRCP<std::string> tdof_names_dot(1);
00223    tdof_names_dot[0] = tdof_names[0]+"_dot";
00224    Teuchos::ArrayRCP<std::string> tresid_names(1);
00225    tresid_names[0] = tdof_names[0]+" Residual";
00226 
00227    fm0.template registerEvaluator<EvalT>
00228      (evalUtils.constructDOFInterpolationEvaluator(tdof_names[0], T_offset));
00229 
00230    (evalUtils.constructDOFInterpolationEvaluator(tdof_names_dot[0], T_offset));
00231 
00232    fm0.template registerEvaluator<EvalT>
00233      (evalUtils.constructDOFGradInterpolationEvaluator(tdof_names[0], T_offset));
00234 
00235    fm0.template registerEvaluator<EvalT>
00236      (evalUtils.constructGatherSolutionEvaluator(false, tdof_names, tdof_names_dot, T_offset));
00237 
00238    fm0.template registerEvaluator<EvalT>
00239      (evalUtils.constructScatterResidualEvaluator(false, tresid_names, T_offset, scatterName));
00240 
00241    // Temperature Variable
00242    Teuchos::ArrayRCP<std::string> tempdof_names(1);
00243    tempdof_names[0] = "Temperature";
00244    Teuchos::ArrayRCP<std::string> tempdof_names_dot(1);
00245    tempdof_names_dot[0] = tempdof_names[0]+"_dot";
00246    Teuchos::ArrayRCP<std::string> tempresid_names(1);
00247    tempresid_names[0] = tempdof_names[0]+" Residual";
00248 
00249    fm0.template registerEvaluator<EvalT>
00250      (evalUtils.constructDOFInterpolationEvaluator(tempdof_names[0], TEMP_offset));
00251 
00252    (evalUtils.constructDOFInterpolationEvaluator(tempdof_names_dot[0], TEMP_offset));
00253 
00254    fm0.template registerEvaluator<EvalT>
00255      (evalUtils.constructDOFGradInterpolationEvaluator(tempdof_names[0], TEMP_offset));
00256 
00257    fm0.template registerEvaluator<EvalT>
00258      (evalUtils.constructGatherSolutionEvaluator(false, tempdof_names, tempdof_names_dot, TEMP_offset));
00259 
00260    fm0.template registerEvaluator<EvalT>
00261      (evalUtils.constructScatterResidualEvaluator(false, tempresid_names, TEMP_offset, "Scatter Temperature"));
00262 
00263    // ----------------------setup the solution field ---------------//
00264 
00265    // General FEM stuff
00266    fm0.template registerEvaluator<EvalT>
00267      (evalUtils.constructGatherCoordinateVectorEvaluator());
00268 
00269    fm0.template registerEvaluator<EvalT>
00270      (evalUtils.constructMapToPhysicalFrameEvaluator(cellType, cubature));
00271 
00272    fm0.template registerEvaluator<EvalT>
00273      (evalUtils.constructComputeBasisFunctionsEvaluator(cellType, intrepidBasis, cubature));
00274 
00275    // Temporary variable used numerous times below
00276    Teuchos::RCP<PHX::Evaluator<AlbanyTraits> > ev;
00277 
00278    { // Time
00279      RCP<ParameterList> p = rcp(new ParameterList);
00280 
00281      p->set<std::string>("Time Name", "Time");
00282      p->set<std::string>("Delta Time Name", " Delta Time");
00283      p->set< RCP<DataLayout> >("Workset Scalar Data Layout", dl->workset_scalar);
00284      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00285      p->set<bool>("Disable Transient", true);
00286 
00287      ev = rcp(new LCM::Time<EvalT,AlbanyTraits>(*p));
00288      fm0.template registerEvaluator<EvalT>(ev);
00289      p = stateMgr.registerStateVariable("Time",dl->workset_scalar, dl->dummy, elementBlockName, "scalar", 0.0, true);
00290      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00291      fm0.template registerEvaluator<EvalT>(ev);
00292    }
00293 
00294    { // Constant Stabilization Parameter
00295      RCP<ParameterList> p = rcp(new ParameterList);
00296 
00297      p->set<std::string>("Material Property Name", "Stabilization Parameter");
00298      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00299      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00300      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00301      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00302      Teuchos::ParameterList& paramList = params->sublist("Stabilization Parameter");
00303      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00304 
00305      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00306      fm0.template registerEvaluator<EvalT>(ev);
00307    }
00308 
00309 
00310    { // Constant Reference Temperature
00311      RCP<ParameterList> p = rcp(new ParameterList);
00312 
00313      p->set<std::string>("Material Property Name", "Reference Temperature");
00314      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00315      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00316      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00317      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00318      Teuchos::ParameterList& paramList = params->sublist("Reference Temperature");
00319      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00320      double refT = paramList.get("Value", 0.0);
00321 
00322      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00323      fm0.template registerEvaluator<EvalT>(ev);
00324 
00325    }
00326 
00327    { // Skeleton Thermal Expansion
00328      RCP<ParameterList> p = rcp(new ParameterList);
00329 
00330      p->set<std::string>("Material Property Name", "Skeleton Thermal Expansion");
00331      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00332      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00333      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00334      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00335      Teuchos::ParameterList& paramList = params->sublist("Skeleton Thermal Expansion");
00336      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00337      double skAlpha = paramList.get("Value", 0.0);
00338 
00339      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00340      fm0.template registerEvaluator<EvalT>(ev);
00341 
00342      p = stateMgr.registerStateVariable("Skeleton Thermal Expansion",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", skAlpha, true);
00343      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00344      fm0.template registerEvaluator<EvalT>(ev);
00345    }
00346 
00347    { // Pore-Fluid Thermal Expansion
00348      RCP<ParameterList> p = rcp(new ParameterList);
00349 
00350      p->set<std::string>("Material Property Name", "Pore-Fluid Thermal Expansion");
00351      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00352      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00353      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00354      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00355      Teuchos::ParameterList& paramList = params->sublist("Pore-Fluid Thermal Expansion");
00356      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00357      double fAlpha = paramList.get("Value", 0.0);
00358 
00359      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00360      fm0.template registerEvaluator<EvalT>(ev);
00361 
00362      p = stateMgr.registerStateVariable("Pore-Fluid Thermal Expansion",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", fAlpha, true);
00363      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00364      fm0.template registerEvaluator<EvalT>(ev);
00365    }
00366 
00367    { // Skeleton Density
00368      RCP<ParameterList> p = rcp(new ParameterList);
00369 
00370      p->set<std::string>("Material Property Name", "Skeleton Density");
00371      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00372      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00373      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00374      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00375      Teuchos::ParameterList& paramList = params->sublist("Skeleton Density");
00376      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00377 
00378      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00379      fm0.template registerEvaluator<EvalT>(ev);
00380    }
00381 
00382    { // Pore-Fluid Density
00383      RCP<ParameterList> p = rcp(new ParameterList);
00384 
00385      p->set<std::string>("Material Property Name", "Pore-Fluid Density");
00386      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00387      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00388      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00389      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00390      Teuchos::ParameterList& paramList = params->sublist("Pore-Fluid Density");
00391      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00392 
00393      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00394      fm0.template registerEvaluator<EvalT>(ev);
00395    }
00396 
00397    { // Skeleton Specific Heat
00398      RCP<ParameterList> p = rcp(new ParameterList);
00399 
00400      p->set<std::string>("Material Property Name", "Skeleton Specific Heat");
00401      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00402      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00403      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00404      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00405      Teuchos::ParameterList& paramList = params->sublist("Skeleton Specific Heat");
00406      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00407 
00408      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00409      fm0.template registerEvaluator<EvalT>(ev);
00410    }
00411 
00412    { // Pore-Fluid Specific Heat
00413      RCP<ParameterList> p = rcp(new ParameterList);
00414 
00415      p->set<std::string>("Material Property Name", "Pore-Fluid Specific Heat");
00416      p->set< RCP<DataLayout> >("Data Layout", dl->qp_scalar);
00417      p->set<std::string>("Coordinate Vector Name", "Coord Vec");
00418      p->set< RCP<DataLayout> >("Coordinate Vector Data Layout", dl->qp_vector);
00419      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00420      Teuchos::ParameterList& paramList = params->sublist("Pore-Fluid Specific Heat");
00421      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00422 
00423      ev = rcp(new PHAL::NSMaterialProperty<EvalT,AlbanyTraits>(*p));
00424      fm0.template registerEvaluator<EvalT>(ev);
00425    }
00426 
00427    { // Mixture Specific Heat
00428      RCP<ParameterList> p = rcp(new ParameterList("Mixture Specific Heat"));
00429 
00430      //Input
00431      p->set<std::string>("Pore-Fluid Specific Heat Name", "Pore-Fluid Specific Heat");
00432      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00433      p->set<std::string>("Skeleton Specific Heat Name", "Skeleton Specific Heat");
00434 
00435      p->set<std::string>("Pore-Fluid Density Name", "Pore-Fluid Density");
00436      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00437      p->set<std::string>("Skeleton Density Name", "Skeleton Density");
00438      p->set<std::string>("DetDefGrad Name", "Jacobian");
00439 
00440      p->set<std::string>("Porosity Name", "Porosity");
00441 
00442      //Output
00443      p->set<std::string>("Mixture Specific Heat Name", "Mixture Specific Heat");
00444 
00445      ev = rcp(new LCM::MixtureSpecificHeat<EvalT,AlbanyTraits>(*p));
00446      fm0.template registerEvaluator<EvalT>(ev);
00447      p = stateMgr.registerStateVariable("Mixture Specific Heat",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0);
00448      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00449      fm0.template registerEvaluator<EvalT>(ev);
00450    }
00451 
00452    { // Mixture Thermal Expansion
00453      RCP<ParameterList> p = rcp(new ParameterList("Mixture Thermal Expansion"));
00454 
00455      //Input
00456      p->set<std::string>("Pore-Fluid Thermal Expansion Name", "Pore-Fluid Thermal Expansion");
00457      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00458      p->set<std::string>("Skeleton Thermal Expansion Name", "Skeleton Thermal Expansion");
00459 
00460      p->set<std::string>("Porosity Name", "Porosity");
00461      p->set<std::string>("Biot Coefficient Name", "Biot Coefficient");
00462      p->set<std::string>("DetDefGrad Name", "Jacobian");
00463 
00464      //Output
00465      p->set<std::string>("Mixture Thermal Expansion Name", "Mixture Thermal Expansion");
00466 
00467      ev = rcp(new LCM::MixtureThermalExpansion<EvalT,AlbanyTraits>(*p));
00468      fm0.template registerEvaluator<EvalT>(ev);
00469      p = stateMgr.registerStateVariable("Mixture Thermal Expansion",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0);
00470      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00471      fm0.template registerEvaluator<EvalT>(ev);
00472    }
00473 
00474 
00475 
00476 /*
00477    { // Strain
00478      RCP<ParameterList> p = rcp(new ParameterList("Strain"));
00479 
00480      //Input
00481      p->set<std::string>("Gradient QP Variable Name", "Displacement Gradient");
00482 
00483      //Output
00484      p->set<std::string>("Strain Name", "Strain"); //dl->qp_tensor also
00485 
00486      ev = rcp(new LCM::Strain<EvalT,AlbanyTraits>(*p,dl));
00487      fm0.template registerEvaluator<EvalT>(ev);
00488      p = stateMgr.registerStateVariable("Strain",dl->qp_tensor, dl->dummy, elementBlockName, "scalar", 0.0, true);
00489      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00490      fm0.template registerEvaluator<EvalT>(ev);
00491    }
00492    */
00493 
00494    {  // Porosity
00495      RCP<ParameterList> p = rcp(new ParameterList);
00496 
00497 
00498      p->set<std::string>("Porosity Name", "Porosity");
00499      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00500      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00501      Teuchos::ParameterList& paramList = params->sublist("Porosity");
00502      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00503      double initPorosity = paramList.get("Value", 0.0);
00504 
00505      // Setting this turns on dependence of strain and pore pressure)
00506      p->set<std::string>("DetDefGrad Name", "Jacobian");
00507 
00508      // porosity update based on Coussy's poromechanics (see p.79)
00509      p->set<std::string>("QP Pore Pressure Name", "Pore Pressure");
00510      p->set<std::string>("Biot Coefficient Name", "Biot Coefficient");
00511 
00512      p->set<std::string>("QP Temperature Name", "Temperature");
00513      p->set<std::string>("Skeleton Thermal Expansion Name", "Skeleton Thermal Expansion");
00514      p->set<std::string>("Reference Temperature Name", "Reference Temperature");
00515 
00516      ev = rcp(new LCM::Porosity<EvalT,AlbanyTraits>(*p,dl));
00517      fm0.template registerEvaluator<EvalT>(ev);
00518      p = stateMgr.registerStateVariable("Porosity",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", initPorosity, true);
00519      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00520      fm0.template registerEvaluator<EvalT>(ev);
00521    }
00522 
00523 
00524 
00525    { // Biot Coefficient
00526      RCP<ParameterList> p = rcp(new ParameterList);
00527 
00528      p->set<std::string>("Biot Coefficient Name", "Biot Coefficient");
00529      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00530      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00531      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00532      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00533 
00534      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00535      Teuchos::ParameterList& paramList = params->sublist("Biot Coefficient");
00536      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00537      double initBiotCoeff = paramList.get("Value", 1.0);
00538 
00539      // Setting this turns on linear dependence on porosity
00540      p->set<std::string>("Porosity Name", "Porosity");
00541 
00542      ev = rcp(new LCM::BiotCoefficient<EvalT,AlbanyTraits>(*p));
00543      fm0.template registerEvaluator<EvalT>(ev);
00544      p = stateMgr.registerStateVariable("Biot Coefficient",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", initBiotCoeff, true);
00545      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00546      fm0.template registerEvaluator<EvalT>(ev);
00547    }
00548 
00549    { // Biot Modulus
00550      RCP<ParameterList> p = rcp(new ParameterList);
00551 
00552      p->set<std::string>("Biot Modulus Name", "Biot Modulus");
00553      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00554      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00555      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00556      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00557 
00558      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00559      Teuchos::ParameterList& paramList = params->sublist("Biot Modulus");
00560      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00561 
00562      // Setting this turns on linear dependence on porosity and Biot's coeffcient
00563      p->set<std::string>("Porosity Name", "Porosity");
00564      p->set<std::string>("Biot Coefficient Name", "Biot Coefficient");
00565 
00566      ev = rcp(new LCM::BiotModulus<EvalT,AlbanyTraits>(*p));
00567      fm0.template registerEvaluator<EvalT>(ev);
00568      p = stateMgr.registerStateVariable("Biot Modulus",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 1.0e12);
00569      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00570      fm0.template registerEvaluator<EvalT>(ev);
00571    }
00572 
00573    { // Thermal conductivity
00574      RCP<ParameterList> p = rcp(new ParameterList);
00575 
00576      p->set<std::string>("QP Variable Name", "Thermal Conductivity");
00577      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00578      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00579      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00580      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00581 
00582      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00583      Teuchos::ParameterList& paramList = params->sublist("Thermal Conductivity");
00584      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00585 
00586      ev = rcp(new PHAL::ThermalConductivity<EvalT,AlbanyTraits>(*p));
00587      fm0.template registerEvaluator<EvalT>(ev);
00588    }
00589 
00590    { // Kozeny-Carman Permeaiblity
00591      RCP<ParameterList> p = rcp(new ParameterList);
00592 
00593      p->set<std::string>("Kozeny-Carman Permeability Name", "Kozeny-Carman Permeability");
00594      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00595      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00596      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00597      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00598 
00599      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00600      Teuchos::ParameterList& paramList = params->sublist("Kozeny-Carman Permeability");
00601      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00602 
00603      // Setting this turns on Kozeny-Carman relation
00604      p->set<std::string>("Porosity Name", "Porosity");
00605      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00606 
00607      ev = rcp(new LCM::KCPermeability<EvalT,AlbanyTraits>(*p));
00608      fm0.template registerEvaluator<EvalT>(ev);
00609      p = stateMgr.registerStateVariable("Kozeny-Carman Permeability",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0);
00610      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00611      fm0.template registerEvaluator<EvalT>(ev);
00612    }
00613 
00614    // Skeleton parameter
00615 
00616    { // Bulk Modulus
00617      RCP<ParameterList> p = rcp(new ParameterList);
00618 
00619      p->set<std::string>("QP Variable Name", "Bulk Modulus");
00620      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00621      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00622      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00623      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00624 
00625      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00626      Teuchos::ParameterList& paramList = params->sublist("Bulk Modulus");
00627      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00628 
00629 
00630      ev = rcp(new LCM::BulkModulus<EvalT,AlbanyTraits>(*p));
00631      fm0.template registerEvaluator<EvalT>(ev);
00632    }
00633 
00634    { // Elastic Modulus
00635      RCP<ParameterList> p = rcp(new ParameterList);
00636 
00637      p->set<std::string>("QP Variable Name", "Elastic Modulus");
00638      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00639      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00640      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00641      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00642 
00643      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00644      Teuchos::ParameterList& paramList = params->sublist("Elastic Modulus");
00645      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00646 
00647      p->set<std::string>("Porosity Name", "Porosity"); // porosity is defined at Cubature points
00648      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00649 
00650      ev = rcp(new LCM::ElasticModulus<EvalT,AlbanyTraits>(*p));
00651      fm0.template registerEvaluator<EvalT>(ev);
00652    }
00653 
00654    { // Shear Modulus
00655      RCP<ParameterList> p = rcp(new ParameterList);
00656 
00657      p->set<std::string>("QP Variable Name", "Shear Modulus");
00658      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00659      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00660      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00661      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00662 
00663      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00664      Teuchos::ParameterList& paramList = params->sublist("Shear Modulus");
00665      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00666 
00667      ev = rcp(new LCM::ShearModulus<EvalT,AlbanyTraits>(*p));
00668      fm0.template registerEvaluator<EvalT>(ev);
00669    }
00670 
00671    { // Poissons Ratio 
00672      RCP<ParameterList> p = rcp(new ParameterList);
00673 
00674      p->set<std::string>("QP Variable Name", "Poissons Ratio");
00675      p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00676      p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00677      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00678      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00679 
00680      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00681      Teuchos::ParameterList& paramList = params->sublist("Poissons Ratio");
00682      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00683 
00684      // Setting this turns on linear dependence of nu on T, nu = nu_ + dnudT*T)
00685      //p->set<std::string>("QP Pore Pressure Name", "Pore Pressure");
00686 
00687      ev = rcp(new LCM::PoissonsRatio<EvalT,AlbanyTraits>(*p));
00688      fm0.template registerEvaluator<EvalT>(ev);
00689    }
00690 
00691    if (matModel == "Neohookean")
00692    {
00693      { // Stress
00694        RCP<ParameterList> p = rcp(new ParameterList("Stress"));
00695 
00696        //Input
00697        p->set<std::string>("DefGrad Name", "Deformation Gradient");
00698        p->set<std::string>("Elastic Modulus Name", "Elastic Modulus");
00699        p->set<std::string>("Poissons Ratio Name", "Poissons Ratio");
00700        p->set<std::string>("DetDefGrad Name", "Jacobian");
00701 
00702        //Output
00703        p->set<std::string>("Stress Name", matModel);
00704 
00705        ev = rcp(new LCM::Neohookean<EvalT,AlbanyTraits>(*p,dl));
00706        fm0.template registerEvaluator<EvalT>(ev);
00707        p = stateMgr.registerStateVariable(matModel,dl->qp_tensor, dl->dummy, elementBlockName, "scalar", 0.0);
00708        ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00709        fm0.template registerEvaluator<EvalT>(ev);
00710      }
00711    }
00712    else if (matModel == "Neohookean AD")
00713    {
00714      RCP<ParameterList> p = rcp(new ParameterList("Stress"));
00715 
00716      //Input
00717      p->set<std::string>("Elastic Modulus Name", "Elastic Modulus");
00718      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00719      p->set<std::string>("Poissons Ratio Name", "Poissons Ratio");  // dl->qp_scalar also
00720 
00721      p->set<std::string>("DefGrad Name", "Deformation Gradient");
00722      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00723 
00724      //Output
00725      p->set<std::string>("Stress Name", matModel); //dl->qp_tensor also
00726 
00727      ev = rcp(new LCM::PisdWdF<EvalT,AlbanyTraits>(*p));
00728      fm0.template registerEvaluator<EvalT>(ev);
00729      p = stateMgr.registerStateVariable(matModel,dl->qp_tensor, dl->dummy, elementBlockName, "scalar", 0.0);
00730      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00731      fm0.template registerEvaluator<EvalT>(ev);
00732    }
00733    else if (matModel == "J2")
00734    {
00735      { // Hardening Modulus
00736        RCP<ParameterList> p = rcp(new ParameterList);
00737 
00738        p->set<std::string>("QP Variable Name", "Hardening Modulus");
00739        p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00740        p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00741        p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00742        p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00743 
00744        p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00745        Teuchos::ParameterList& paramList = params->sublist("Hardening Modulus");
00746        p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00747 
00748        ev = rcp(new LCM::HardeningModulus<EvalT,AlbanyTraits>(*p));
00749        fm0.template registerEvaluator<EvalT>(ev);
00750      }
00751 
00752      { // Yield Strength
00753        RCP<ParameterList> p = rcp(new ParameterList);
00754 
00755        p->set<std::string>("QP Variable Name", "Yield Strength");
00756        p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00757        p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00758        p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00759        p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00760 
00761        p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00762        Teuchos::ParameterList& paramList = params->sublist("Yield Strength");
00763        p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00764 
00765        ev = rcp(new LCM::YieldStrength<EvalT,AlbanyTraits>(*p));
00766        fm0.template registerEvaluator<EvalT>(ev);
00767      }
00768 
00769      { // Saturation Modulus
00770        RCP<ParameterList> p = rcp(new ParameterList);
00771 
00772        p->set<std::string>("Saturation Modulus Name", "Saturation Modulus");
00773        p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00774        p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00775        p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00776        p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00777 
00778        p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00779        Teuchos::ParameterList& paramList = params->sublist("Saturation Modulus");
00780        p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00781 
00782        ev = rcp(new LCM::SaturationModulus<EvalT,AlbanyTraits>(*p));
00783        fm0.template registerEvaluator<EvalT>(ev);
00784      }
00785 
00786      { // Saturation Exponent
00787        RCP<ParameterList> p = rcp(new ParameterList);
00788 
00789        p->set<std::string>("Saturation Exponent Name", "Saturation Exponent");
00790        p->set<std::string>("QP Coordinate Vector Name", "Coord Vec");
00791        p->set< RCP<DataLayout> >("Node Data Layout", dl->node_scalar);
00792        p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00793        p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
00794 
00795        p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00796        Teuchos::ParameterList& paramList = params->sublist("Saturation Exponent");
00797        p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00798 
00799        ev = rcp(new LCM::SaturationExponent<EvalT,AlbanyTraits>(*p));
00800        fm0.template registerEvaluator<EvalT>(ev);
00801      }
00802 
00803      if ( numDim == 3 && params->get("Compute Dislocation Density Tensor", false) )
00804      { // Dislocation Density Tensor
00805        RCP<ParameterList> p = rcp(new ParameterList("Dislocation Density"));
00806 
00807        //Input
00808        p->set<std::string>("Fp Name", "Fp");
00809        p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00810        p->set<std::string>("BF Name", "BF");
00811        p->set< RCP<DataLayout> >("Node QP Scalar Data Layout", dl->node_qp_scalar);
00812        p->set<std::string>("Gradient BF Name", "Grad BF");
00813        p->set< RCP<DataLayout> >("Node QP Vector Data Layout", dl->node_qp_vector);
00814 
00815        //Output
00816        p->set<std::string>("Dislocation Density Name", "G"); //dl->qp_tensor also
00817 
00818        //Declare what state data will need to be saved (name, layout, init_type)
00819        ev = rcp(new LCM::DislocationDensity<EvalT,AlbanyTraits>(*p));
00820        fm0.template registerEvaluator<EvalT>(ev);
00821        p = stateMgr.registerStateVariable("G",dl->qp_tensor, dl->dummy, elementBlockName, "scalar", 0.0);
00822        ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00823        fm0.template registerEvaluator<EvalT>(ev);
00824      }
00825 
00826      {// Stress
00827            RCP<ParameterList> p = rcp(new ParameterList("Stress"));
00828 
00829            //Input
00830            p->set<std::string>("DefGrad Name", "Deformation Gradient");
00831            p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00832 
00833            p->set<std::string>("Elastic Modulus Name", "Elastic Modulus");
00834            p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00835 
00836            p->set<std::string>("Poissons Ratio Name", "Poissons Ratio");  // dl->qp_scalar also
00837            p->set<std::string>("Hardening Modulus Name", "Hardening Modulus"); // dl->qp_scalar also
00838            p->set<std::string>("Saturation Modulus Name", "Saturation Modulus"); // dl->qp_scalar also
00839            p->set<std::string>("Saturation Exponent Name", "Saturation Exponent"); // dl->qp_scalar also
00840            p->set<std::string>("Yield Strength Name", "Yield Strength"); // dl->qp_scalar also
00841            p->set<std::string>("DetDefGrad Name", "Jacobian");  // dl->qp_scalar also
00842 
00843            //Output
00844            p->set<std::string>("Stress Name", matModel); //dl->qp_tensor also
00845            p->set<std::string>("Fp Name", "Fp");  // dl->qp_tensor also
00846            p->set<std::string>("Eqps Name", "eqps");  // dl->qp_scalar also
00847 
00848            //Declare what state data will need to be saved (name, layout, init_type)
00849 
00850            ev = rcp(new LCM::J2Stress<EvalT,AlbanyTraits>(*p));
00851            fm0.template registerEvaluator<EvalT>(ev);
00852            p = stateMgr.registerStateVariable(matModel,dl->qp_tensor, dl->dummy, elementBlockName, "scalar", 0.0);
00853            ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00854            fm0.template registerEvaluator<EvalT>(ev);
00855            p = stateMgr.registerStateVariable("Fp",dl->qp_tensor, dl->dummy, elementBlockName, "identity", 1.0, true);
00856            ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00857            fm0.template registerEvaluator<EvalT>(ev);
00858            p = stateMgr.registerStateVariable("eqps",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0, true);
00859            ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00860            fm0.template registerEvaluator<EvalT>(ev);
00861          }
00862    }
00863 
00864 
00865    { // Total Stress
00866      RCP<ParameterList> p = rcp(new ParameterList("Total Stress"));
00867 
00868      //Input
00869      p->set<std::string>("Stress Name", matModel);
00870      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00871 
00872      p->set<std::string>("DefGrad Name", "Deformation Gradient");
00873      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00874 
00875 
00876      p->set<std::string>("Biot Coefficient Name", "Biot Coefficient");  // dl->qp_scalar also
00877      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00878 
00879      p->set<std::string>("QP Variable Name", "Pore Pressure");
00880      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00881 
00882 
00883      p->set<std::string>("DetDefGrad Name", "Jacobian");  // dl->qp_scalar also
00884 
00885      //Output
00886      p->set<std::string>("Total Stress Name", "Total Stress"); //dl->qp_tensor also
00887 
00888 
00889      ev = rcp(new LCM::TLPoroStress<EvalT,AlbanyTraits>(*p));
00890      fm0.template registerEvaluator<EvalT>(ev);
00891      p = stateMgr.registerStateVariable("Total Stress",dl->qp_tensor, dl->dummy, elementBlockName, "scalar", 0.0);
00892      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00893      fm0.template registerEvaluator<EvalT>(ev);
00894      p = stateMgr.registerStateVariable("Pore Pressure",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0, true);
00895      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00896      fm0.template registerEvaluator<EvalT>(ev);
00897 
00898 
00899    }
00900 
00901    if (haveSource) { // Source
00902      RCP<ParameterList> p = rcp(new ParameterList);
00903 
00904      p->set<std::string>("Source Name", "Source");
00905      p->set<std::string>("QP Variable Name", "Pore Pressure");
00906      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00907 
00908      p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00909      Teuchos::ParameterList& paramList = params->sublist("Source Functions");
00910      p->set<Teuchos::ParameterList*>("Parameter List", &paramList);
00911 
00912      ev = rcp(new PHAL::Source<EvalT,AlbanyTraits>(*p));
00913      fm0.template registerEvaluator<EvalT>(ev);
00914    }
00915 
00916    { // Deformation Gradient
00917      RCP<ParameterList> p = rcp(new ParameterList("Deformation Gradient"));
00918 
00919      //Inputs: flags, weights, GradU
00920      const bool weighted_Volume_Averaged_J = params->get("weighted_Volume_Averaged_J", false);
00921      p->set<bool>("weighted_Volume_Averaged_J Name", weighted_Volume_Averaged_J);
00922      p->set<std::string>("Weights Name","Weights");
00923      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00924      p->set<std::string>("Gradient QP Variable Name", "Displacement Gradient");
00925      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00926 
00927      //Outputs: F, J
00928      p->set<std::string>("DefGrad Name", "Deformation Gradient"); //dl->qp_tensor also
00929      p->set<std::string>("DetDefGrad Name", "Jacobian");
00930      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00931 
00932      ev = rcp(new LCM::DefGrad<EvalT,AlbanyTraits>(*p));
00933      fm0.template registerEvaluator<EvalT>(ev);
00934 
00935      p = stateMgr.registerStateVariable("Displacement Gradient",dl->qp_tensor,
00936                                             dl->dummy, elementBlockName, "identity",1.0,true);
00937          ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00938          fm0.template registerEvaluator<EvalT>(ev);
00939 
00940      p = stateMgr.registerStateVariable("Jacobian",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 1.0, true);
00941            ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00942            fm0.template registerEvaluator<EvalT>(ev);
00943 
00944 
00945 
00946    }
00947 
00948    { // Displacement Resid
00949      RCP<ParameterList> p = rcp(new ParameterList("Displacement Residual"));
00950 
00951      //Input
00952      p->set<std::string>("Total Stress Name", "Total Stress");
00953      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00954 
00955      p->set<std::string>("DefGrad Name", "Deformation Gradient");
00956      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
00957 
00958      p->set<std::string>("DetDefGrad Name", "Jacobian");
00959      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00960 
00961      p->set<std::string>("Skeleton Thermal Expansion Name", "Skeleton Thermal Expansion");
00962      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00963 
00964      p->set<std::string>("Bulk Modulus Name", "Bulk Modulus");
00965      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00966 
00967      p->set<std::string>("Temperature Name", "Temperature");
00968      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00969 
00970      p->set<std::string>("Reference Temperature Name", "Reference Temperature");
00971      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00972 
00973      p->set<std::string>("Weighted BF Name", "wBF");
00974      p->set< RCP<DataLayout> >("Node QP Scalar Data Layout", dl->node_qp_scalar);
00975 
00976 
00977      p->set<std::string>("Weighted Gradient BF Name", "wGrad BF");
00978      p->set< RCP<DataLayout> >("Node QP Vector Data Layout", dl->node_qp_vector);
00979 
00980      p->set<bool>("Disable Transient", true);
00981 
00982 
00983      //Output
00984      p->set<std::string>("Residual Name", "Displacement Residual");
00985      p->set< RCP<DataLayout> >("Node Vector Data Layout", dl->node_vector);
00986 
00987      ev = rcp(new LCM::ThermoPoroPlasticityResidMomentum<EvalT,AlbanyTraits>(*p));
00988      fm0.template registerEvaluator<EvalT>(ev);
00989 
00990 
00991 
00992    }
00993 
00994 
00995 
00996    { // Pore Pressure Resid
00997      RCP<ParameterList> p = rcp(new ParameterList("Pore Pressure Residual"));
00998 
00999      //Input
01000 
01001      // Input from nodal points
01002      p->set<std::string>("Weighted BF Name", "wBF");
01003      p->set< RCP<DataLayout> >("Node QP Scalar Data Layout", dl->node_qp_scalar);
01004 
01005      p->set<bool>("Have Source", false);
01006      p->set<std::string>("Source Name", "Source");
01007 
01008      p->set<bool>("Have Absorption", false);
01009 
01010      // Input from cubature points
01011      p->set<std::string>("QP Pore Pressure Name", "Pore Pressure");
01012      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01013 
01014      p->set<std::string>("Pore-Fluid Density Name", "Pore-Fluid Density");
01015      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01016 
01017      p->set<std::string>("Elastic Modulus Name", "Elastic Modulus");
01018      p->set<std::string>("Poissons Ratio Name", "Poissons Ratio");
01019 
01020      p->set<std::string>("QP Temperature Name", "Temperature");
01021 
01022      p->set<std::string>("Reference Temperature Name", "Reference Temperature");
01023      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01024 
01025      p->set<std::string>("Material Property Name", "Stabilization Parameter");
01026      p->set<std::string>("Thermal Conductivity Name", "Thermal Conductivity");
01027      p->set<std::string>("Porosity Name", "Porosity");
01028      p->set<std::string>("Kozeny-Carman Permeability Name", "Kozeny-Carman Permeability");
01029      p->set<std::string>("Biot Coefficient Name", "Biot Coefficient");
01030      p->set<std::string>("Biot Modulus Name", "Biot Modulus");
01031 
01032      p->set<std::string>("Mixture Thermal Expansion Name", "Mixture Thermal Expansion");
01033      p->set<std::string>("Skeleton Thermal Expansion Name", "Skeleton Thermal Expansion");
01034      p->set<std::string>("Pore-Fluid Thermal Expansion Name", "Pore-Fluid Thermal Expansion");
01035 
01036      p->set<std::string>("Gradient QP Variable Name", "Pore Pressure Gradient");
01037      p->set<std::string>("Temperature Gradient Name", "Temperature Gradient");
01038      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
01039 
01040      p->set<std::string>("Weighted Gradient BF Name", "wGrad BF");
01041      p->set< RCP<DataLayout> >("Node QP Vector Data Layout", dl->node_qp_vector);
01042 
01043      p->set<std::string>("Weights Name","Weights");
01044 
01045      p->set<std::string>("Delta Time Name", " Delta Time");
01046      p->set< RCP<DataLayout> >("Workset Scalar Data Layout", dl->workset_scalar);
01047 
01048      p->set<std::string>("DefGrad Name", "Deformation Gradient");
01049      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
01050 
01051      p->set<std::string>("DetDefGrad Name", "Jacobian");
01052      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01053 
01054      //Output
01055      p->set<std::string>("Residual Name", "Pore Pressure Residual");
01056      p->set< RCP<DataLayout> >("Node Scalar Data Layout", dl->node_scalar);
01057 
01058      ev = rcp(new LCM::ThermoPoroPlasticityResidMass<EvalT,AlbanyTraits>(*p));
01059      fm0.template registerEvaluator<EvalT>(ev);
01060 
01061    }
01062 
01063    { // Temperature Resid
01064      RCP<ParameterList> p = rcp(new ParameterList("Temperature Residual"));
01065 
01066      //Input
01067 
01068      // Input from nodal points
01069      p->set<std::string>("Weighted BF Name", "wBF");
01070      p->set< RCP<DataLayout> >("Node QP Scalar Data Layout", dl->node_qp_scalar);
01071 
01072      p->set<bool>("Have Source", false);
01073      p->set<std::string>("Source Name", "Source");
01074 
01075      p->set<bool>("Have Absorption", false);
01076 
01077      // Input from cubature points
01078      p->set<std::string>("QP Pore Pressure Name", "Pore Pressure");
01079      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01080 
01081      p->set<std::string>("Bulk Modulus Name", "Bulk Modulus");
01082      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01083 
01084      p->set<std::string>("QP Temperature Name", "Temperature");
01085      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01086 
01087      p->set<std::string>("Mixture Specific Heat Name", "Mixture Specific Heat");
01088      p->set<std::string>("Pore-Fluid Specific Heat Name", "Pore-Fluid Specific Heat");
01089      p->set<std::string>("Skeleton Thermal Expansion Name", "Skeleton Thermal Expansion");
01090      p->set<std::string>("Reference Temperature Name", "Reference Temperature");
01091 
01092      p->set<std::string>("Material Property Name", "Stabilization Parameter");
01093      p->set<std::string>("Thermal Conductivity Name", "Thermal Conductivity");
01094      p->set<std::string>("Porosity Name", "Porosity");
01095      p->set<std::string>("Kozeny-Carman Permeability Name", "Kozeny-Carman Permeability");
01096      p->set<std::string>("Biot Coefficient Name", "Biot Coefficient");
01097      p->set<std::string>("Biot Modulus Name", "Biot Modulus");
01098 
01099      p->set<std::string>("Mixture Thermal Expansion Name", "Mixture Thermal Expansion");
01100 
01101      p->set<std::string>("Gradient QP Variable Name", "Temperature Gradient");
01102      p->set<std::string>("Pore Pressure Gradient Name", "Pore Pressure Gradient");
01103      p->set< RCP<DataLayout> >("QP Vector Data Layout", dl->qp_vector);
01104 
01105      p->set<std::string>("Weighted Gradient BF Name", "wGrad BF");
01106      p->set< RCP<DataLayout> >("Node QP Vector Data Layout", dl->node_qp_vector);
01107 
01108      // Inputs: X, Y at nodes, Cubature, and Basis
01109      p->set<std::string>("Coordinate Vector Name","Coord Vec");
01110      p->set< RCP<DataLayout> >("Coordinate Data Layout", dl->vertices_vector);
01111      p->set< RCP<Intrepid::Cubature<RealType> > >("Cubature", cubature);
01112      p->set<RCP<shards::CellTopology> >("Cell Type", cellType);
01113 
01114      p->set<std::string>("Weights Name","Weights");
01115 
01116      p->set<std::string>("Delta Time Name", " Delta Time");
01117      p->set< RCP<DataLayout> >("Workset Scalar Data Layout", dl->workset_scalar);
01118 
01119      p->set<std::string>("DefGrad Name", "Deformation Gradient");
01120      p->set< RCP<DataLayout> >("QP Tensor Data Layout", dl->qp_tensor);
01121 
01122      p->set<std::string>("DetDefGrad Name", "Jacobian");
01123      p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
01124 
01125      p->set<std::string>("Elastic Modulus Name", "Elastic Modulus");
01126      p->set<std::string>("Poissons Ratio Name", "Poissons Ratio");
01127 
01128      //Output
01129      p->set<std::string>("Residual Name", "Temperature Residual");
01130      p->set< RCP<DataLayout> >("Node Scalar Data Layout", dl->node_scalar);
01131 
01132      ev = rcp(new LCM::ThermoPoroPlasticityResidEnergy<EvalT,AlbanyTraits>(*p));
01133      fm0.template registerEvaluator<EvalT>(ev);
01134 
01135      p = stateMgr.registerStateVariable("Temperature",dl->qp_scalar, dl->dummy, elementBlockName, "scalar", 0.0, true);
01136      ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
01137      fm0.template registerEvaluator<EvalT>(ev);
01138 
01139 
01140    }
01141 
01142 
01143    if (fieldManagerChoice == Albany::BUILD_RESID_FM)  {
01144      PHX::Tag<typename EvalT::ScalarT> res_tag("Scatter", dl->dummy);
01145      fm0.requireField<EvalT>(res_tag);
01146 
01147      PHX::Tag<typename EvalT::ScalarT> res_tag2(scatterName, dl->dummy);
01148      fm0.requireField<EvalT>(res_tag2);
01149 
01150      PHX::Tag<typename EvalT::ScalarT> res_tag3("Scatter Temperature", dl->dummy);
01151      fm0.requireField<EvalT>(res_tag3);
01152 
01153      return res_tag.clone();
01154    }
01155    else if (fieldManagerChoice == Albany::BUILD_RESPONSE_FM) {
01156      Albany::ResponseUtilities<EvalT, PHAL::AlbanyTraits> respUtils(dl);
01157      return respUtils.constructResponses(fm0, *responseList, stateMgr);
01158    }
01159 
01160    return Teuchos::null;
01161 }
01162 
01163 #endif // THERMOPOROPLASTICITYPROBLEM_HPP

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