43 #ifdef PANZER_HAVE_TEKO
45 namespace panzer_stk {
50 template <
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename Node>
51 ParameterListCallback<LocalOrdinalT,GlobalOrdinalT,Node>::ParameterListCallback(
52 const std::string & coordFieldName,
56 : coordFieldName_(coordFieldName), fieldPatterns_(fps), connManager_(connManager), ugi_(ugi), coordinatesBuilt_(false)
59 template <
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename Node>
68 for(itr=inputPL->
begin();itr!=inputPL->
end();++itr)
69 setFieldByKey(itr->first,*outputPL);
74 template <
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename Node>
75 bool ParameterListCallback<LocalOrdinalT,GlobalOrdinalT,Node>::handlesRequest(
const Teko::RequestMesg & rm)
79 if(rm.getName()==
"Parameter List")
return true;
83 template <
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename Node>
84 void ParameterListCallback<LocalOrdinalT,GlobalOrdinalT,Node>::preRequest(
const Teko::RequestMesg & rm)
93 template <
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename Node>
94 void ParameterListCallback<LocalOrdinalT,GlobalOrdinalT,Node>::setFieldByKey(
const std::string & key,
Teuchos::ParameterList & pl)
const
97 "ParameterListCallback::setFieldByKey: Coordinates have not been built!");
99 double * x =
const_cast<double *
>(&xcoords_[0]);
100 double * y =
const_cast<double *
>(&ycoords_[0]);
101 double * z =
const_cast<double *
>(&zcoords_[0]);
103 if(key==
"x-coordinates")
104 pl.
set<
double*>(key,x);
105 else if(key==
"y-coordinates")
106 pl.
set<
double*>(key,y);
107 else if(key==
"z-coordinates")
108 pl.
set<
double*>(key,z);
111 "ParameterListCallback cannot handle key=\"" << key <<
"\"");
114 template <
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename Node>
115 void ParameterListCallback<LocalOrdinalT,GlobalOrdinalT,Node>::buildArrayToVector()
117 if(arrayToVector_==Teuchos::null)
121 template <
typename LocalOrdinalT,
typename GlobalOrdinalT,
typename Node>
122 void ParameterListCallback<LocalOrdinalT,GlobalOrdinalT,Node>::buildCoordinates()
126 std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > data;
128 std::map<std::string,Teuchos::RCP<const panzer::Intrepid2FieldPattern> >::const_iterator itr;
129 for(itr=fieldPatterns_.begin();itr!=fieldPatterns_.end();++itr) {
130 std::string blockId = itr->first;
132 std::vector<std::size_t> localCellIds;
135 Kokkos::DynRankView<double,PHX::Device> & fieldData = data[blockId];
136 fieldData = Kokkos::DynRankView<double,PHX::Device>(
"fieldData",connManager_->getElementBlock(blockId).size(),fieldPattern->
numberIds());
140 connManager_->getDofCoords(blockId,*fieldPattern,localCellIds,fieldData);
144 out.setOutputToRootOnly(-1);
145 out <<
"WARNING: In ParameterListCallback::buildCoordinates(), the Intrepid2::FieldPattern in "
146 <<
"block \"" << blockId <<
"\" does not support interpolatory coordinates. "
147 <<
"This may be fine if coordinates are not actually needed. However if they are then bad things "
148 <<
"will happen. Enjoy!" << std::endl;
150 coordinatesBuilt_ =
true;
156 = arrayToVector_->template getDataVector<double>(coordFieldName_,data);
158 switch(resultVec->getNumVectors()) {
160 zcoords_.resize(resultVec->getLocalLength());
161 resultVec->getVector(2)->get1dCopy(Teuchos::arrayViewFromVector(zcoords_));
164 ycoords_.resize(resultVec->getLocalLength());
165 resultVec->getVector(1)->get1dCopy(Teuchos::arrayViewFromVector(ycoords_));
168 xcoords_.resize(resultVec->getLocalLength());
169 resultVec->getVector(0)->get1dCopy(Teuchos::arrayViewFromVector(xcoords_));
173 "ParameterListCallback::buildCoordinates: Constructed multivector has nonphysical dimensions.");
177 coordinatesBuilt_ =
true;
ConstIterator end() const
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual int numberIds() const
bool supportsInterpolatoryCoordinates() const
Does this field pattern support interpolatory coordinates?
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
params_t::ConstIterator ConstIterator
ConstIterator begin() const
#define TEUCHOS_ASSERT(assertion_test)