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

ConcurrentMultiscaleProblem.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 #if !defined(LCM_ConcurrentMultiscaleProblem_hpp)
00008 #define LCM_ConcurrentMultiscaleProblem_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 
00025 class ConcurrentMultiscaleProblem : public Albany::AbstractProblem {
00026 public:
00027 
00028   typedef Intrepid::FieldContainer<RealType> FC;
00029 
00033   ConcurrentMultiscaleProblem(
00034       Teuchos::RCP<Teuchos::ParameterList> const & params,
00035       Teuchos::RCP<ParamLib> const & param_lib,
00036       int const num_dims,
00037       Teuchos::RCP<const Epetra_Comm> const & comm);
00038 
00042   virtual
00043   ~ConcurrentMultiscaleProblem();
00044 
00046   Teuchos::RCP<std::map<std::string, std::string> >
00047   constructFieldNameMap(bool surface_flag);
00048 
00052   virtual
00053   int
00054   spatialDimension() const {return num_dims_;}
00055 
00059   virtual
00060   void
00061   buildProblem(
00062       Teuchos::ArrayRCP<Teuchos::RCP<Albany::MeshSpecsStruct> > mesh_specs,
00063       StateManager & state_mgr);
00064 
00068   virtual
00069   Teuchos::Array< Teuchos::RCP<const PHX::FieldTag> >
00070   buildEvaluators(
00071       PHX::FieldManager<PHAL::AlbanyTraits> & fm0,
00072       Albany::MeshSpecsStruct const & mesh_specs,
00073       Albany::StateManager & state_mgr,
00074       Albany::FieldManagerChoice fm_choice,
00075       Teuchos::RCP<Teuchos::ParameterList> const & response_list);
00076 
00080   Teuchos::RCP<Teuchos::ParameterList const>
00081   getValidProblemParameters() const;
00082 
00086   void
00087   getAllocatedStates(
00088       Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<FC> > > old_state,
00089       Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<FC> > > new_state
00090   ) const;
00091 
00092 private:
00093 
00097   ConcurrentMultiscaleProblem(ConcurrentMultiscaleProblem const &);
00098 
00102   ConcurrentMultiscaleProblem &
00103   operator=(ConcurrentMultiscaleProblem const &);
00104 
00105   QCAD::MaterialDatabase &
00106   matDB()
00107   {return *material_db_;}
00108 
00109 public:
00110 
00115   template <typename EvalT>
00116   Teuchos::RCP<const PHX::FieldTag>
00117   constructEvaluators(
00118       PHX::FieldManager<PHAL::AlbanyTraits> & fm0,
00119       Albany::MeshSpecsStruct const & mesh_specs,
00120       Albany::StateManager & state_mgr,
00121       Albany::FieldManagerChoice fm_choice,
00122       Teuchos::RCP<Teuchos::ParameterList> & response_list);
00123 
00127   void
00128   constructDirichletEvaluators(Albany::MeshSpecsStruct const & mesh_specs);
00129 
00130 protected:
00131 
00135   bool have_source_;
00136 
00140   int num_dims_;
00141 
00145   int num_pts_;
00146 
00150   int num_nodes_;
00151 
00155   int num_vertices_;
00156 
00160   std::map< std::string, bool > coarse_overlap_map_;
00161 
00165   std::map< std::string, bool > fine_overlap_map_;
00166 
00170   std::map< std::string, bool > lm_overlap_map_;
00171 
00175   Teuchos::RCP<QCAD::MaterialDatabase> material_db_;
00176 
00180   Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<FC> > > old_state_;
00181 
00185   Teuchos::ArrayRCP<Teuchos::ArrayRCP<Teuchos::RCP<FC> > > new_state_;
00186 
00187 };
00188 
00189 } // namespace Albany
00190 
00191 
00192 #include "Albany_Utils.hpp"
00193 #include "Albany_ProblemUtils.hpp"
00194 #include "Albany_ResponseUtilities.hpp"
00195 #include "Albany_EvaluatorUtils.hpp"
00196 #include "PHAL_NSMaterialProperty.hpp"
00197 #include "PHAL_Source.hpp"
00198 #include "PHAL_SaveStateField.hpp"
00199 
00200 #include "FieldNameMap.hpp"
00201 
00202 #include "MechanicsResidual.hpp"
00203 #include "Time.hpp"
00204 
00205 // Constitutive Model Interface and parameters
00206 #include "Kinematics.hpp"
00207 #include "ConstitutiveModelInterface.hpp"
00208 #include "ConstitutiveModelParameters.hpp"
00209 
00210 //
00211 //
00212 //
00213 template <typename EvalT>
00214 Teuchos::RCP<const PHX::FieldTag>
00215 Albany::ConcurrentMultiscaleProblem::
00216 constructEvaluators(
00217     PHX::FieldManager<PHAL::AlbanyTraits> & fm0,
00218     Albany::MeshSpecsStruct const & mesh_specs,
00219     Albany::StateManager & state_mgr,
00220     Albany::FieldManagerChoice fm_choice,
00221     Teuchos::RCP<Teuchos::ParameterList> & response_list)
00222 {
00223   using Teuchos::RCP;
00224   using Teuchos::rcp;
00225   using Teuchos::ParameterList;
00226   using PHX::DataLayout;
00227   using PHX::MDALayout;
00228   using std::vector;
00229   using PHAL::AlbanyTraits;
00230   using shards::CellTopology;
00231   using shards::getCellTopologyData;
00232 
00233   // get the name of the current element block
00234   std::string
00235   eb_name = mesh_specs.ebName;
00236 
00237   // get the name of the material model to be used (and make sure there is one)
00238   Teuchos::ParameterList &
00239   material_model_sublist =
00240       matDB().getElementBlockSublist(eb_name,"Material Model");
00241 
00242   std::string
00243   material_model_name = material_model_sublist.get<std::string>("Model Name");
00244 
00245   TEUCHOS_TEST_FOR_EXCEPTION(
00246       material_model_name.length() == 0, std::logic_error,
00247       "A material model must be defined for block: " + eb_name);
00248 
00249 #ifdef ALBANY_VERBOSE
00250   *out << "In ConcurrentMultiscaleProblem::constructEvaluators" << '\n';
00251   *out << "element block name: " << eb_name << '\n';
00252   *out << "material model name: " << material_model_name << '\n';
00253 #endif
00254 
00255   // define cell topologies
00256   RCP<CellTopology>
00257   composite_cell_type =
00258     rcp(new CellTopology(getCellTopologyData<shards::Tetrahedron<11> >()));
00259 
00260   RCP<shards::CellTopology>
00261   cell_type = rcp(new CellTopology (&mesh_specs.ctd));
00262 
00263   // Check if we are setting the composite tet flag
00264   bool
00265   is_composite = false;
00266 
00267   bool const
00268   is_composite_block_present =
00269       matDB().isElementBlockParam(eb_name, "Use Composite Tet 10");
00270 
00271   if (is_composite_block_present == true) {
00272     is_composite =
00273         matDB().getElementBlockParam<bool>(eb_name, "Use Composite Tet 10");
00274   }
00275 
00276   // get the intrepid basis for the given cell topology
00277   RCP<Intrepid::Basis<RealType, Intrepid::FieldContainer<RealType> > >
00278   intrepid_basis = Albany::getIntrepidBasis(mesh_specs.ctd, is_composite);
00279 
00280   bool const
00281   is_composite_cell_type =
00282       is_composite &&
00283       mesh_specs.ctd.dimension == 3 &&
00284       mesh_specs.ctd.node_count == 10;
00285 
00286   if (is_composite_cell_type == true) {
00287     cell_type = composite_cell_type;
00288   }
00289 
00290   Intrepid::DefaultCubatureFactory<RealType>
00291   cubature_factory;
00292 
00293   RCP <Intrepid::Cubature<RealType> >
00294   cubature = cubature_factory.create(*cell_type, mesh_specs.cubatureDegree);
00295 
00296   // Note that these are the volume element quantities
00297   num_nodes_ = intrepid_basis->getCardinality();
00298   int const
00299   workset_size = mesh_specs.worksetSize;
00300 
00301   num_dims_ = cubature->getDimension();
00302   num_pts_ = cubature->getNumPoints();
00303   num_vertices_ = num_nodes_;
00304 
00305 #ifdef ALBANY_VERBOSE
00306   *out << "Field Dimensions: Workset=" << workset_size
00307        << ", Vertices= " << num_vertices_
00308        << ", Nodes= " << num_nodes_
00309        << ", QuadPts= " << num_pts_
00310        << ", Dim= " << num_dims_ << '\n';
00311 #endif
00312 
00313   // Construct standard FEM evaluators with standard field names
00314   RCP<Albany::Layouts>
00315   dl =
00316       rcp(
00317           new Albany::Layouts(
00318               workset_size, num_vertices_, num_nodes_, num_pts_, num_dims_
00319           )
00320       );
00321 
00322   std::string const
00323   msg = "Data Layout Usage in Mechanics problems assume vecDim = num_dims_";
00324 
00325   TEUCHOS_TEST_FOR_EXCEPTION(
00326       dl->vectorAndGradientLayoutsAreEquivalent == false,
00327       std::logic_error,
00328       msg);
00329 
00330   Albany::EvaluatorUtils<EvalT, PHAL::AlbanyTraits>
00331   eu(dl);
00332 
00333   int
00334   offset = 0;
00335 
00336   // Temporary variable used numerous times below
00337   RCP<PHX::Evaluator<AlbanyTraits> > ev;
00338 
00339   // Define Field Names
00340   {
00341     Teuchos::ArrayRCP<std::string>
00342     dof_names(1);
00343 
00344     Teuchos::ArrayRCP<std::string>
00345     resid_names(1);
00346 
00347     dof_names[0] = "Displacement";
00348     resid_names[0] = dof_names[0] + " Residual";
00349 
00350     fm0.template registerEvaluator<EvalT>(
00351         eu.constructGatherSolutionEvaluator_noTransient(true, dof_names)
00352     );
00353 
00354     fm0.template registerEvaluator<EvalT>(
00355         eu.constructGatherCoordinateVectorEvaluator()
00356     );
00357 
00358     fm0.template registerEvaluator<EvalT>(
00359         eu.constructDOFVecInterpolationEvaluator(dof_names[0])
00360     );
00361 
00362     fm0.template registerEvaluator<EvalT>(
00363         eu.constructDOFVecGradInterpolationEvaluator(dof_names[0])
00364     );
00365 
00366     fm0.template registerEvaluator<EvalT>(
00367         eu.constructMapToPhysicalFrameEvaluator(cell_type, cubature)
00368     );
00369 
00370     fm0.template registerEvaluator<EvalT>(
00371         eu.constructComputeBasisFunctionsEvaluator(
00372             cell_type, intrepid_basis, cubature
00373         )
00374     );
00375 
00376     fm0.template registerEvaluator<EvalT>(
00377         eu.constructScatterResidualEvaluator(true, resid_names)
00378     );
00379 
00380     offset += num_dims_;
00381   }
00382 
00383   bool const
00384   have_lagrange_multiplier = lm_overlap_map_[eb_name];
00385 
00386   if (have_lagrange_multiplier == true) { // add Lagrange multiplier field
00387     Teuchos::ArrayRCP<std::string>
00388     dof_names(1);
00389 
00390     Teuchos::ArrayRCP<std::string>
00391     resid_names(1);
00392 
00393     dof_names[0] = "lagrange_multiplier";
00394     resid_names[0] = dof_names[0] + "_residual";
00395 
00396     fm0.template registerEvaluator<EvalT>(
00397         eu.constructGatherSolutionEvaluator_noTransient(
00398             true, dof_names, offset
00399         )
00400     );
00401 
00402     fm0.template registerEvaluator<EvalT>(
00403         eu.constructGatherCoordinateVectorEvaluator()
00404     );
00405 
00406     fm0.template registerEvaluator<EvalT>(
00407         eu.constructDOFVecInterpolationEvaluator(dof_names[0], offset)
00408     );
00409 
00410     fm0.template registerEvaluator<EvalT>(
00411         eu.constructDOFVecGradInterpolationEvaluator(dof_names[0], offset)
00412     );
00413 
00414     fm0.template registerEvaluator<EvalT>(
00415         eu.constructMapToPhysicalFrameEvaluator(cell_type, cubature)
00416     );
00417 
00418     fm0.template registerEvaluator<EvalT>(
00419         eu.constructComputeBasisFunctionsEvaluator(
00420             cell_type, intrepid_basis, cubature
00421         )
00422     );
00423 
00424     offset += num_dims_;
00425   }
00426 
00427   // generate the field name map to deal with outputting surface element info
00428   LCM::FieldNameMap
00429   field_name_map(false);
00430 
00431   RCP<std::map<std::string, std::string> >
00432   fnm = field_name_map.getMap();
00433 
00434   std::string cauchy       = (*fnm)["Cauchy_Stress"];
00435   std::string Fp           = (*fnm)["Fp"];
00436   std::string eqps         = (*fnm)["eqps"];
00437 
00438   { // Time
00439     RCP<ParameterList>
00440     p = rcp(new ParameterList("Time"));
00441     p->set<std::string>("Time Name", "Time");
00442     p->set<std::string>("Delta Time Name", "Delta Time");
00443     p->set<RCP<DataLayout> >("Workset Scalar Data Layout", dl->workset_scalar);
00444     p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00445     p->set<bool>("Disable Transient", true);
00446     ev = rcp(new LCM::Time<EvalT,AlbanyTraits>(*p));
00447     fm0.template registerEvaluator<EvalT>(ev);
00448     p = state_mgr.registerStateVariable(
00449         "Time",
00450         dl->workset_scalar,
00451         dl->dummy,
00452         eb_name,
00453         "scalar",
00454         0.0,
00455         true
00456     );
00457     ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00458     fm0.template registerEvaluator<EvalT>(ev);
00459   }
00460 
00461   if (have_lagrange_multiplier == true) {
00462     RCP<ParameterList>
00463     p = rcp(new ParameterList("Save Lagrange Multiplier"));
00464     p = state_mgr.registerStateVariable(
00465         "Lagrange Multiplier",
00466         dl->qp_scalar,
00467         dl->dummy,
00468         eb_name,
00469         "scalar",
00470         0.0,
00471         true,
00472         false
00473     );
00474     ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00475     fm0.template registerEvaluator<EvalT>(ev);
00476   }
00477 
00478   if (have_source_ == true) { // Source
00479     RCP<ParameterList>
00480     p = rcp(new ParameterList);
00481     p->set<std::string>("Source Name", "Source");
00482     p->set<std::string>("Variable Name", "Displacement");
00483     p->set< RCP<DataLayout> >("QP Scalar Data Layout", dl->qp_scalar);
00484 
00485     p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00486     Teuchos::ParameterList &
00487     param_list = params->sublist("Source Functions");
00488     p->set<Teuchos::ParameterList*>("Parameter List", &param_list);
00489 
00490     ev = rcp(new PHAL::Source<EvalT,AlbanyTraits>(*p));
00491     fm0.template registerEvaluator<EvalT>(ev);
00492   }
00493 
00494   { // Constitutive Model Parameters
00495     RCP<ParameterList>
00496     p = rcp(new ParameterList("Constitutive Model Parameters"));
00497 
00498     std::string const
00499     name = matDB().getElementBlockParam<std::string>(eb_name,"material");
00500 
00501     Teuchos::ParameterList &
00502     param_list = matDB().getElementBlockSublist(eb_name, name);
00503 
00504     p->set<Teuchos::ParameterList*>("Material Parameters", &param_list);
00505 
00506     RCP<LCM::ConstitutiveModelParameters<EvalT,AlbanyTraits> >
00507     cmp_ev = rcp(
00508         new LCM::ConstitutiveModelParameters<EvalT, AlbanyTraits>(*p, dl)
00509     );
00510     fm0.template registerEvaluator<EvalT>(cmp_ev);
00511   }
00512 
00513   {
00514     RCP<ParameterList>
00515     p = rcp(new ParameterList("Constitutive Model Interface"));
00516 
00517     std::string const
00518     name = matDB().getElementBlockParam<std::string>(eb_name, "material");
00519 
00520     Teuchos::ParameterList &
00521     param_list = matDB().getElementBlockSublist(eb_name, name);
00522 
00523     param_list.set<RCP<std::map<std::string, std::string> > >("Name Map", fnm);
00524     p->set<Teuchos::ParameterList*>("Material Parameters", &param_list);
00525 
00526     RCP<LCM::ConstitutiveModelInterface<EvalT,AlbanyTraits> >
00527     cmi_ev = rcp(
00528         new LCM::ConstitutiveModelInterface<EvalT, AlbanyTraits>(*p, dl)
00529     );
00530     fm0.template registerEvaluator<EvalT>(cmi_ev);
00531 
00532     // register state variables
00533     for (int sv(0); sv < cmi_ev->getNumStateVars(); ++sv) {
00534       cmi_ev->fillStateVariableStruct(sv);
00535       p = state_mgr.registerStateVariable(
00536           cmi_ev->getName(),
00537           cmi_ev->getLayout(),
00538           dl->dummy,
00539           eb_name,
00540           cmi_ev->getInitType(),
00541           cmi_ev->getInitValue(),
00542           cmi_ev->getStateFlag(),
00543           cmi_ev->getOutputFlag()
00544       );
00545       ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00546       fm0.template registerEvaluator<EvalT>(ev);
00547     }
00548   }
00549 
00550 
00551   { // Kinematics quantities
00552     RCP<ParameterList>
00553     p = rcp(new ParameterList("Kinematics"));
00554 
00555     std::string const
00556     wva_str = "Weighted Volume Average J";
00557 
00558     bool const
00559     is_wva = matDB().isElementBlockParam(eb_name, wva_str);
00560 
00561     // set flags to optionally volume average J with a weighted average
00562     if (is_wva == true) {
00563       bool const
00564       ebp_wva = matDB().getElementBlockParam<bool>(eb_name, wva_str);
00565       p->set<bool>(wva_str, ebp_wva);
00566     }
00567 
00568     std::string const
00569     asp_str = "Average J Stabilization Parameter";
00570 
00571     bool const
00572     is_asp = matDB().isElementBlockParam(eb_name, asp_str);
00573 
00574     if (is_asp == true) {
00575       bool const
00576       ebp_asp = matDB().getElementBlockParam<RealType>(eb_name, asp_str);
00577       p->set<RealType>(asp_str, ebp_asp);
00578     }
00579 
00580     // set flag for return strain and velocity gradient
00581     bool have_strain(false), have_velocity_gradient(false);
00582 
00583     std::string const
00584     str_str = "Strain Flag";
00585 
00586     bool const
00587     is_str = matDB().isElementBlockParam(eb_name, str_str);
00588 
00589     if (is_str == true) {
00590       bool const
00591       ebp_str = matDB().getElementBlockParam<bool>(eb_name, str_str);
00592       p->set<bool>(str_str, ebp_str);
00593 
00594       if (ebp_str == true) {
00595         p->set<std::string>("Strain Name", "Strain");
00596       }
00597     }
00598 
00599     std::string const
00600     vgf_str = "Velocity Gradient Flag";
00601 
00602     bool const
00603     is_vgf = matDB().isElementBlockParam(eb_name, vgf_str);
00604 
00605     if (is_vgf == true) {
00606       bool const
00607       ebp_vgf = matDB().getElementBlockParam<bool>(eb_name, vgf_str);
00608       p->set<bool>(vgf_str, ebp_vgf);
00609 
00610       if (ebp_vgf == true) {
00611         p->set<std::string>("Velocity Gradient Name", "Velocity Gradient");
00612       }
00613     }
00614 
00615     // send in integration weights and the displacement gradient
00616     p->set<std::string>("Weights Name","Weights");
00617     p->set<std::string>("Gradient QP Variable Name", "Displacement Gradient");
00618 
00619     //Outputs: F, J
00620     p->set<std::string>("DefGrad Name", "F"); //dl->qp_tensor also
00621     p->set<std::string>("DetDefGrad Name", "J");
00622 
00623     ev = rcp(new LCM::Kinematics<EvalT,AlbanyTraits>(*p, dl));
00624     fm0.template registerEvaluator<EvalT>(ev);
00625 
00626 
00627     // optional output
00628     bool output_flag(false);
00629     if (matDB().isElementBlockParam(eb_name, "Output Deformation Gradient")) {
00630       output_flag = matDB().getElementBlockParam<bool>(
00631           eb_name,
00632           "Output Deformation Gradient"
00633       );
00634     }
00635 
00636     p = state_mgr.registerStateVariable(
00637         "F",
00638         dl->qp_tensor,
00639         dl->dummy,
00640         eb_name,
00641         "identity",
00642         1.0,
00643         output_flag
00644     );
00645     ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00646     fm0.template registerEvaluator<EvalT>(ev);
00647 
00648     // need J and J_old to perform time integration for poromechanics problem
00649     output_flag = false;
00650 
00651     if (matDB().isElementBlockParam(eb_name, "Output J")) {
00652       output_flag = matDB().getElementBlockParam<bool>(
00653           eb_name,
00654           "Output J"
00655       );
00656     }
00657 
00658     if (output_flag == true) {
00659       p = state_mgr.registerStateVariable(
00660           "J",
00661           dl->qp_scalar,
00662           dl->dummy,
00663           eb_name,
00664           "scalar",
00665           1.0,
00666           true
00667       );
00668       ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00669       fm0.template registerEvaluator<EvalT>(ev);
00670     }
00671 
00672     // Optional output: strain
00673     if (have_strain == true) {
00674       output_flag = false;
00675       if (matDB(). isElementBlockParam(eb_name, "Output Strain")) {
00676         output_flag =
00677           matDB().getElementBlockParam<bool>(eb_name, "Output Strain");
00678       }
00679 
00680       p = state_mgr.registerStateVariable(
00681           "Strain",
00682           dl->qp_tensor,
00683           dl->dummy,
00684           eb_name,
00685           "scalar",
00686           0.0,
00687           output_flag
00688       );
00689       ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00690       fm0.template registerEvaluator<EvalT>(ev);
00691     }
00692 
00693     // Optional output: velocity gradient
00694     if (have_velocity_gradient == true) {
00695       output_flag = false;
00696       if(matDB().isElementBlockParam(eb_name, "Output Velocity Gradient")) {
00697         output_flag = matDB().getElementBlockParam<bool>(
00698             eb_name,
00699             "Output Velocity Gradient"
00700         );
00701       }
00702 
00703       p = state_mgr.registerStateVariable(
00704           "Velocity Gradient",
00705           dl->qp_tensor,
00706           dl->dummy,
00707           eb_name,
00708           "scalar",
00709           0.0,
00710           output_flag
00711       );
00712       ev = rcp(new PHAL::SaveStateField<EvalT,AlbanyTraits>(*p));
00713       fm0.template registerEvaluator<EvalT>(ev);
00714     }
00715 
00716   }
00717 
00718   { // Residual
00719     RCP<ParameterList>
00720     p = rcp(new ParameterList("Displacement Residual"));
00721     //Input
00722     p->set<std::string>("Stress Name", cauchy);
00723     p->set<std::string>("DefGrad Name", "F");
00724     p->set<std::string>("DetDefGrad Name", "J");
00725     p->set<std::string>("Weighted Gradient BF Name", "wGrad BF");
00726     p->set<std::string>("Weighted BF Name", "wBF");
00727 
00728     // Strain flag for small deformation problem
00729     if (matDB().isElementBlockParam(eb_name,"Strain Flag")) {
00730       p->set<bool>("Strain Flag","Strain Flag");
00731     }
00732 
00733     p->set<RCP<ParamLib> >("Parameter Library", paramLib);
00734     //Output
00735     p->set<std::string>("Residual Name", "Displacement Residual");
00736     ev = rcp(new LCM::MechanicsResidual<EvalT,AlbanyTraits>(*p,dl));
00737     fm0.template registerEvaluator<EvalT>(ev);
00738   }
00739 
00740   if (fm_choice == Albany::BUILD_RESID_FM)  {
00741     Teuchos::RCP<const PHX::FieldTag>
00742     ret_tag;
00743     {
00744       PHX::Tag<typename EvalT::ScalarT>
00745       res_tag("Scatter", dl->dummy);
00746       fm0.requireField<EvalT>(res_tag);
00747       ret_tag = res_tag.clone();
00748     }
00749     if (have_lagrange_multiplier == true) {
00750       PHX::Tag<typename EvalT::ScalarT>
00751       lagrange_multiplier_tag("Scatter Lagrange Multiplier", dl->dummy);
00752       fm0.requireField<EvalT>(lagrange_multiplier_tag);
00753       ret_tag = lagrange_multiplier_tag.clone();
00754     }
00755     return ret_tag;
00756   }
00757   else if (fm_choice == Albany::BUILD_RESPONSE_FM) {
00758     Albany::ResponseUtilities<EvalT, PHAL::AlbanyTraits> respUtils(dl);
00759     return respUtils.constructResponses(fm0, *response_list, state_mgr);
00760   }
00761 
00762   return Teuchos::null;
00763 }
00764 
00765 #endif // LCM_ConcurrentMultiscaleProblem_hpp

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