00001
00002
00003
00004
00005
00006
00007 #include "Teuchos_TestForException.hpp"
00008 #include "Phalanx_DataLayout.hpp"
00009
00010 #include "Intrepid_FunctionSpaceTools.hpp"
00011
00012 namespace PHAL {
00013
00014
00015 template<typename EvalT, typename Traits>
00016 DOFGradInterpolation<EvalT, Traits>::
00017 DOFGradInterpolation(const Teuchos::ParameterList& p,
00018 const Teuchos::RCP<Albany::Layouts>& dl) :
00019 val_node (p.get<std::string> ("Variable Name"), dl->node_scalar),
00020 GradBF (p.get<std::string> ("Gradient BF Name"), dl->node_qp_gradient),
00021 grad_val_qp (p.get<std::string> ("Gradient Variable Name"), dl->qp_gradient)
00022 {
00023 this->addDependentField(val_node);
00024 this->addDependentField(GradBF);
00025 this->addEvaluatedField(grad_val_qp);
00026
00027 this->setName("DOFGradInterpolation"+PHX::TypeString<EvalT>::value);
00028
00029 std::vector<PHX::DataLayout::size_type> dims;
00030 GradBF.fieldTag().dataLayout().dimensions(dims);
00031 numNodes = dims[1];
00032 numQPs = dims[2];
00033 numDims = dims[3];
00034 }
00035
00036
00037 template<typename EvalT, typename Traits>
00038 void DOFGradInterpolation<EvalT, Traits>::
00039 postRegistrationSetup(typename Traits::SetupData d,
00040 PHX::FieldManager<Traits>& fm)
00041 {
00042 this->utils.setFieldData(val_node,fm);
00043 this->utils.setFieldData(GradBF,fm);
00044 this->utils.setFieldData(grad_val_qp,fm);
00045 }
00046
00047
00048 template<typename EvalT, typename Traits>
00049 void DOFGradInterpolation<EvalT, Traits>::
00050 evaluateFields(typename Traits::EvalData workset)
00051 {
00052
00053
00054
00055
00056 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00057 for (std::size_t qp=0; qp < numQPs; ++qp) {
00058 for (std::size_t dim=0; dim<numDims; dim++) {
00059 ScalarT& gvqp = grad_val_qp(cell,qp,dim);
00060 gvqp = val_node(cell, 0) * GradBF(cell, 0, qp, dim);
00061 for (std::size_t node= 1 ; node < numNodes; ++node) {
00062 gvqp += val_node(cell, node) * GradBF(cell, node, qp, dim);
00063 }
00064 }
00065 }
00066 }
00067 }
00068
00069
00070 template<typename Traits>
00071 DOFGradInterpolation<PHAL::AlbanyTraits::Jacobian, Traits>::
00072 DOFGradInterpolation(const Teuchos::ParameterList& p,
00073 const Teuchos::RCP<Albany::Layouts>& dl) :
00074 val_node (p.get<std::string> ("Variable Name"), dl->node_scalar),
00075 GradBF (p.get<std::string> ("Gradient BF Name"), dl->node_qp_gradient),
00076 grad_val_qp (p.get<std::string> ("Gradient Variable Name"), dl->qp_gradient)
00077 {
00078 this->addDependentField(val_node);
00079 this->addDependentField(GradBF);
00080 this->addEvaluatedField(grad_val_qp);
00081
00082 this->setName("DOFGradInterpolation"+PHX::TypeString<PHAL::AlbanyTraits::Jacobian>::value);
00083
00084 std::vector<PHX::DataLayout::size_type> dims;
00085 GradBF.fieldTag().dataLayout().dimensions(dims);
00086 numNodes = dims[1];
00087 numQPs = dims[2];
00088 numDims = dims[3];
00089
00090 offset = p.get<int>("Offset of First DOF");
00091 }
00092
00093
00094 template<typename Traits>
00095 void DOFGradInterpolation<PHAL::AlbanyTraits::Jacobian, Traits>::
00096 postRegistrationSetup(typename Traits::SetupData d,
00097 PHX::FieldManager<Traits>& fm)
00098 {
00099 this->utils.setFieldData(val_node,fm);
00100 this->utils.setFieldData(GradBF,fm);
00101 this->utils.setFieldData(grad_val_qp,fm);
00102 }
00103
00104
00105 template<typename Traits>
00106 void DOFGradInterpolation<PHAL::AlbanyTraits::Jacobian, Traits>::
00107 evaluateFields(typename Traits::EvalData workset)
00108 {
00109
00110
00111
00112
00113 int num_dof = val_node(0,0,0).size();
00114 int neq = num_dof / numNodes;
00115
00116 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00117 for (std::size_t qp=0; qp < numQPs; ++qp) {
00118 for (std::size_t dim=0; dim<numDims; dim++) {
00119 ScalarT& gvqp = grad_val_qp(cell,qp,dim);
00120 gvqp = FadType(num_dof, val_node(cell, 0).val() * GradBF(cell, 0, qp, dim));
00121 gvqp.fastAccessDx(offset) = val_node(cell, 0).fastAccessDx(offset) * GradBF(cell, 0, qp, dim);
00122 for (std::size_t node= 1 ; node < numNodes; ++node) {
00123 gvqp.val() += val_node(cell, node).val() * GradBF(cell, node, qp, dim);
00124 gvqp.fastAccessDx(neq*node+offset) += val_node(cell, node).fastAccessDx(neq*node+offset) * GradBF(cell, node, qp, dim);
00125 }
00126 }
00127 }
00128 }
00129 }
00130
00131
00132
00133
00134 template<typename EvalT, typename Traits>
00135 DOFGradInterpolation_noDeriv<EvalT, Traits>::
00136 DOFGradInterpolation_noDeriv(const Teuchos::ParameterList& p,
00137 const Teuchos::RCP<Albany::Layouts>& dl) :
00138 val_node (p.get<std::string> ("Variable Name"), dl->node_scalar),
00139 GradBF (p.get<std::string> ("Gradient BF Name"), dl->node_qp_gradient),
00140 grad_val_qp (p.get<std::string> ("Gradient Variable Name"), dl->qp_gradient)
00141 {
00142 this->addDependentField(val_node);
00143 this->addDependentField(GradBF);
00144 this->addEvaluatedField(grad_val_qp);
00145
00146 this->setName("DOFGradInterpolation_noDeriv"+PHX::TypeString<EvalT>::value);
00147
00148 std::vector<PHX::DataLayout::size_type> dims;
00149 GradBF.fieldTag().dataLayout().dimensions(dims);
00150 numNodes = dims[1];
00151 numQPs = dims[2];
00152 numDims = dims[3];
00153 }
00154
00155
00156 template<typename EvalT, typename Traits>
00157 void DOFGradInterpolation_noDeriv<EvalT, Traits>::
00158 postRegistrationSetup(typename Traits::SetupData d,
00159 PHX::FieldManager<Traits>& fm)
00160 {
00161 this->utils.setFieldData(val_node,fm);
00162 this->utils.setFieldData(GradBF,fm);
00163 this->utils.setFieldData(grad_val_qp,fm);
00164 }
00165
00166
00167 template<typename EvalT, typename Traits>
00168 void DOFGradInterpolation_noDeriv<EvalT, Traits>::
00169 evaluateFields(typename Traits::EvalData workset)
00170 {
00171
00172
00173
00174
00175 for (std::size_t cell=0; cell < workset.numCells; ++cell) {
00176 for (std::size_t qp=0; qp < numQPs; ++qp) {
00177 for (std::size_t dim=0; dim<numDims; dim++) {
00178 MeshScalarT& gvqp = grad_val_qp(cell,qp,dim);
00179 gvqp = val_node(cell, 0) * GradBF(cell, 0, qp, dim);
00180 for (std::size_t node= 1 ; node < numNodes; ++node) {
00181 gvqp += val_node(cell, node) * GradBF(cell, node, qp, dim);
00182 }
00183 }
00184 }
00185 }
00186 }
00187
00188
00189
00190 }
00191