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