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 BiotCoefficient<EvalT, Traits>::
00017 BiotCoefficient(Teuchos::ParameterList& p) :
00018 biotCoefficient(p.get<std::string>("Biot Coefficient 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("Biot Coefficient 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 "Biot Coefficient", 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("Biot Coefficient 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 Biot coefficient type " << type);
00065 }
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077 isPoroElastic = true;
00078 Kskeleton_value = elmd_list->get("Skeleton Bulk Modulus Parameter Value", 10.0e5);
00079 new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00080 "Skeleton Bulk Modulus Parameter Value", this, paramLib);
00081 Kgrain_value = elmd_list->get("Grain Bulk Modulus Value", 10.0e12);
00082 new Sacado::ParameterRegistration<EvalT, SPL_Traits>(
00083 "Grain Bulk Modulus Value", this, paramLib);
00084
00085
00086
00087
00088
00089
00090
00091
00092 this->addEvaluatedField(biotCoefficient);
00093 this->setName("Biot Coefficient"+PHX::TypeString<EvalT>::value);
00094 }
00095
00096
00097 template<typename EvalT, typename Traits>
00098 void BiotCoefficient<EvalT, Traits>::
00099 postRegistrationSetup(typename Traits::SetupData d,
00100 PHX::FieldManager<Traits>& fm)
00101 {
00102 this->utils.setFieldData(biotCoefficient,fm);
00103 if (!is_constant) this->utils.setFieldData(coordVec,fm);
00104
00105 }
00106
00107
00108 template<typename EvalT, typename Traits>
00109 void BiotCoefficient<EvalT, Traits>::
00110 evaluateFields(typename Traits::EvalData workset)
00111 {
00112 std::size_t numCells = workset.numCells;
00113
00114 if (is_constant) {
00115 for (std::size_t cell=0; cell < numCells; ++cell) {
00116 for (std::size_t qp=0; qp < numQPs; ++qp) {
00117 biotCoefficient(cell,qp) = constant_value;
00118 }
00119 }
00120 }
00121 else {
00122 for (std::size_t cell=0; cell < numCells; ++cell) {
00123 for (std::size_t qp=0; qp < numQPs; ++qp) {
00124 Teuchos::Array<MeshScalarT> point(numDims);
00125 for (std::size_t i=0; i<numDims; i++)
00126 point[i] = Sacado::ScalarValue<MeshScalarT>::eval(coordVec(cell,qp,i));
00127 biotCoefficient(cell,qp) = exp_rf_kl->evaluate(point, rv);
00128 }
00129 }
00130 }
00131 if (isPoroElastic) {
00132 for (std::size_t cell=0; cell < numCells; ++cell) {
00133 for (std::size_t qp=0; qp < numQPs; ++qp) {
00134
00135 biotCoefficient(cell,qp) = 1.0 - Kskeleton_value/Kgrain_value;
00136 }
00137 }
00138 }
00139 }
00140
00141
00142 template<typename EvalT,typename Traits>
00143 typename BiotCoefficient<EvalT,Traits>::ScalarT&
00144 BiotCoefficient<EvalT,Traits>::getValue(const std::string &n)
00145 {
00146 if (n == "Biot Coefficient")
00147 return constant_value;
00148 else if (n == "Skeleton Bulk Modulus Parameter Value")
00149 return Kskeleton_value;
00150 else if (n == "Grain Bulk Modulus Value")
00151 return Kgrain_value;
00152 for (int i=0; i<rv.size(); i++) {
00153 if (n == Albany::strint("Biot Coefficient KL Random Variable",i))
00154 return rv[i];
00155 }
00156 TEUCHOS_TEST_FOR_EXCEPTION(true, Teuchos::Exceptions::InvalidParameter,
00157 std::endl <<
00158 "Error! Logic error in getting paramter " << n
00159 << " in BiotCoefficient::getValue()" << std::endl);
00160 return constant_value;
00161 }
00162
00163
00164
00165 }
00166