00001
00002
00003
00004
00005
00006
00007 #include <fstream>
00008 #include "Teuchos_TestForException.hpp"
00009 #include "Phalanx_DataLayout.hpp"
00010 #include "Sacado_ParameterRegistration.hpp"
00011 #include "Albany_Utils.hpp"
00012
00013 namespace LCM {
00014
00015 template<typename EvalT, typename Traits>
00016 HardeningModulus<EvalT, Traits>::
00017 HardeningModulus(Teuchos::ParameterList& p) :
00018 hardeningModulus(p.get<std::string>("QP Variable Name"),
00019 p.get<Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout"))
00020 {
00021 Teuchos::ParameterList* elmd_list =
00022 p.get<Teuchos::ParameterList*>("Parameter List");
00023
00024 Teuchos::RCP<PHX::DataLayout> vector_dl =
00025 p.get< Teuchos::RCP<PHX::DataLayout> >("QP Vector Data Layout");
00026 std::vector<PHX::DataLayout::size_type> dims;
00027 vector_dl->dimensions(dims);
00028 numQPs = dims[1];
00029 numDims = dims[2];
00030
00031 Teuchos::RCP<ParamLib> paramLib =
00032 p.get< Teuchos::RCP<ParamLib> >("Parameter Library", Teuchos::null);
00033
00034 std::string type = elmd_list->get("Hardening Modulus Type", "Constant");
00035 if (type == "Constant") {
00036 is_constant = true;
00037 constant_value = elmd_list->get("Value", 1.0);
00038
00039
00040 new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00041 "Hardening Modulus", this, paramLib);
00042 }
00043 else if (type == "Truncated KL Expansion") {
00044 is_constant = false;
00045 PHX::MDField<MeshScalarT,Cell,QuadPoint,Dim>
00046 fx(p.get<std::string>("QP Coordinate Vector Name"), vector_dl);
00047 coordVec = fx;
00048 this->addDependentField(coordVec);
00049
00050 exp_rf_kl =
00051 Teuchos::rcp(new Stokhos::KL::ExponentialRandomField<MeshScalarT>(*elmd_list));
00052 int num_KL = exp_rf_kl->stochasticDimension();
00053
00054
00055 rv.resize(num_KL);
00056 for (int i=0; i<num_KL; i++) {
00057 std::string ss = Albany::strint("Hardening Modulus KL Random Variable",i);
00058 new Sacado::ParameterRegistration<EvalT, SPL_Traits>(ss, this, paramLib);
00059 rv[i] = elmd_list->get(ss, 0.0);
00060 }
00061 }
00062 else {
00063 TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00064 "Invalid hardening modulus type " << type);
00065 }
00066
00067
00068
00069 if ( p.isType<std::string>("QP Temperature Name") ) {
00070 Teuchos::RCP<PHX::DataLayout> scalar_dl =
00071 p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
00072 PHX::MDField<ScalarT,Cell,QuadPoint>
00073 tmp(p.get<std::string>("QP Temperature Name"), scalar_dl);
00074 Temperature = tmp;
00075 this->addDependentField(Temperature);
00076 isThermoElastic = true;
00077 dHdT_value = elmd_list->get("dHdT Value", 0.0);
00078 refTemp = p.get<RealType>("Reference Temperature", 0.0);
00079 new Sacado::ParameterRegistration<EvalT, SPL_Traits>("dHdT Value", this, paramLib);
00080 }
00081 else {
00082 isThermoElastic=false;
00083 dHdT_value=0.0;
00084 }
00085
00086
00087 this->addEvaluatedField(hardeningModulus);
00088 this->setName("Hardening Modulus"+PHX::TypeString<EvalT>::value);
00089 }
00090
00091
00092 template<typename EvalT, typename Traits>
00093 void HardeningModulus<EvalT, Traits>::
00094 postRegistrationSetup(typename Traits::SetupData d,
00095 PHX::FieldManager<Traits>& fm)
00096 {
00097 this->utils.setFieldData(hardeningModulus,fm);
00098 if (!is_constant) this->utils.setFieldData(coordVec,fm);
00099 if (isThermoElastic) this->utils.setFieldData(Temperature,fm);
00100 }
00101
00102
00103 template<typename EvalT, typename Traits>
00104 void HardeningModulus<EvalT, Traits>::
00105 evaluateFields(typename Traits::EvalData workset)
00106 {
00107 std::size_t numCells = workset.numCells;
00108
00109 if (is_constant) {
00110 for (std::size_t cell=0; cell < numCells; ++cell) {
00111 for (std::size_t qp=0; qp < numQPs; ++qp) {
00112 hardeningModulus(cell,qp) = constant_value;
00113 }
00114 }
00115 }
00116 else {
00117 for (std::size_t cell=0; cell < numCells; ++cell) {
00118 for (std::size_t qp=0; qp < numQPs; ++qp) {
00119 Teuchos::Array<MeshScalarT> point(numDims);
00120 for (std::size_t i=0; i<numDims; i++)
00121 point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00122 hardeningModulus(cell,qp) = exp_rf_kl->evaluate(point, rv);
00123 }
00124 }
00125 }
00126 if (isThermoElastic) {
00127 for (std::size_t cell=0; cell < numCells; ++cell) {
00128 for (std::size_t qp=0; qp < numQPs; ++qp) {
00129 hardeningModulus(cell,qp) -= dHdT_value * (Temperature(cell,qp) - refTemp);
00130 }
00131 }
00132 }
00133 }
00134
00135
00136 template<typename EvalT,typename Traits>
00137 typename HardeningModulus<EvalT,Traits>::ScalarT&
00138 HardeningModulus<EvalT,Traits>::getValue(const std::string &n)
00139 {
00140 if (n == "Hardening Modulus")
00141 return constant_value;
00142 else if (n == "dHdT Value")
00143 return dHdT_value;
00144 for (int i=0; i<rv.size(); i++) {
00145 if (n == Albany::strint("Hardening Modulus KL Random Variable",i))
00146 return rv[i];
00147 }
00148 TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00149 std::endl <<
00150 "Error! Logic error in getting paramter " << n
00151 << " in HardeningModulus::getValue()" << std::endl);
00152 return constant_value;
00153 }
00154
00155
00156
00157 }
00158