11 #include "PanzerAdaptersSTK_config.hpp"
12 #ifdef PANZER_HAVE_TEKO
16 namespace panzer_stk {
21 ParameterListCallbackBlocked::ParameterListCallbackBlocked(
25 : connManager_(connManager), blocked_ugi_(blocked_ugi), aux_blocked_ugi_(aux_blocked_ugi)
29 ParameterListCallbackBlocked::request(
const Teko::RequestMesg & rm)
37 for(itr=inputPL->
begin();itr!=inputPL->
end();++itr) {
38 std::string * str_ptr = 0;
40 setFieldByKey(itr->first,field,*outputPL);
46 bool ParameterListCallbackBlocked::handlesRequest(
const Teko::RequestMesg & rm)
50 if(rm.getName()==
"Parameter List") {
51 bool isHandled =
true;
54 if(pl->
isType<std::string>(
"x-coordinates")) {
55 field = pl->
get<std::string>(
"x-coordinates");
60 if(pl->
isType<std::string>(
"y-coordinates")) {
62 if(field != pl->
get<std::string>(
"y-coordinates")) {
66 if(pl->
isType<std::string>(
"z-coordinates")) {
68 if(field != pl->
get<std::string>(
"z-coordinates")) {
72 if(pl->
isType<std::string>(
"Coordinates")){
73 field = pl->
get<std::string>(
"Coordinates");
75 #ifdef PANZER_HAVE_EPETRA_STACK
76 if(pl->
isType<std::string>(
"Coordinates-Epetra")){
77 field = pl->
get<std::string>(
"Coordinates-Epetra");
86 void ParameterListCallbackBlocked::preRequest(
const Teko::RequestMesg & rm)
90 const std::string&
field(getHandledField(*rm.getParameterList()));
95 std::vector<Teuchos::RCP<panzer::GlobalIndexer>>
96 fieldDOFMngrs = blocked_ugi_->getFieldDOFManagers();
97 for (
int b(0); b < static_cast<int>(fieldDOFMngrs.size()); ++b)
99 for (
int f(0); f < fieldDOFMngrs[b]->getNumFields(); ++f)
101 if (fieldDOFMngrs[b]->getFieldString(f) ==
field)
109 aux_blocked_ugi_->getFieldBlock(aux_blocked_ugi_->getFieldNum(field));
111 block = blocked_ugi_->getFieldBlock(blocked_ugi_->getFieldNum(field));
114 #ifdef PANZER_HAVE_EPETRA_STACK
115 if (rm.getParameterList()->isType<std::string>(
"Coordinates-Epetra")) {
116 buildArrayToVectorEpetra(block, field, useAux);
117 buildCoordinatesEpetra(field, useAux);
121 buildArrayToVectorTpetra(block, field, useAux);
122 buildCoordinatesTpetra(field, useAux);
126 void ParameterListCallbackBlocked::setFieldByKey(
const std::string & key,
const std::string & field,
Teuchos::ParameterList & pl)
const
129 if(key==
"x-coordinates") {
130 double * x =
const_cast<double *
>(&getCoordinateByField(0,field)[0]);
131 pl.
set<
double*>(key,x);
133 else if(key==
"y-coordinates") {
134 double * y =
const_cast<double *
>(&getCoordinateByField(1,field)[0]);
135 pl.
set<
double*>(key,y);
137 else if(key==
"z-coordinates") {
138 double * z =
const_cast<double *
>(&getCoordinateByField(2,field)[0]);
139 pl.
set<
double*>(key,z);
140 }
else if(key ==
"Coordinates") {
142 #ifdef PANZER_HAVE_EPETRA_STACK
143 }
else if(key ==
"Coordinates-Epetra") {
150 "ParameterListCallback cannot handle key=\"" << key <<
"\"");
153 void ParameterListCallbackBlocked::buildArrayToVectorTpetra(
int block,
const std::string & field,
const bool useAux)
155 if(arrayToVectorTpetra_[field]==Teuchos::null) {
158 ugi = aux_blocked_ugi_->getFieldDOFManagers()[block];
160 ugi = blocked_ugi_->getFieldDOFManagers()[block];
165 #ifdef PANZER_HAVE_EPETRA_STACK
166 void ParameterListCallbackBlocked::buildArrayToVectorEpetra(
int block,
const std::string & field,
const bool useAux)
168 if(arrayToVectorEpetra_[field]==Teuchos::null) {
171 ugi = aux_blocked_ugi_->getFieldDOFManagers()[block];
173 ugi = blocked_ugi_->getFieldDOFManagers()[block];
179 void ParameterListCallbackBlocked::buildCoordinatesTpetra(
const std::string & field,
const bool useAux)
181 std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > data;
185 std::vector<std::string> elementBlocks;
187 aux_blocked_ugi_->getElementBlockIds(elementBlocks);
189 blocked_ugi_->getElementBlockIds(elementBlocks);
190 for(std::size_t i=0;i<elementBlocks.size();++i) {
191 std::string blockId = elementBlocks[i];
192 std::vector<std::size_t> localCellIds;
195 Kokkos::DynRankView<double,PHX::Device> & fieldData = data[blockId];
196 fieldData = Kokkos::DynRankView<double,PHX::Device>(
"fieldData",connManager_->getElementBlock(blockId).size(),fieldPattern->
numberIds());
200 connManager_->getDofCoords(blockId,*fieldPattern,localCellIds,fieldData);
204 out.setOutputToRootOnly(-1);
205 out <<
"WARNING: In ParameterListCallback::buildCoordinates(), the Intrepid2::FieldPattern in "
206 <<
"block \"" << blockId <<
"\" does not support interpolatory coordinates. "
207 <<
"This may be fine if coordinates are not actually needed. However if they are then bad things "
208 <<
"will happen. Enjoy!" << std::endl;
214 coordsVecTp_ = arrayToVectorTpetra_[
field]->template getDataVector<double>(
field,data);
216 switch(coordsVecTp_->getNumVectors()) {
218 zcoords_[
field].resize(coordsVecTp_->getLocalLength());
219 coordsVecTp_->getVector(2)->get1dCopy(Teuchos::arrayViewFromVector(zcoords_[field]));
222 ycoords_[
field].resize(coordsVecTp_->getLocalLength());
223 coordsVecTp_->getVector(1)->get1dCopy(Teuchos::arrayViewFromVector(ycoords_[field]));
226 xcoords_[
field].resize(coordsVecTp_->getLocalLength());
227 coordsVecTp_->getVector(0)->get1dCopy(Teuchos::arrayViewFromVector(xcoords_[field]));
231 "ParameterListCallback::buildCoordinates: Constructed multivector has nonphysical dimensions.");
236 #ifdef PANZER_HAVE_EPETRA_STACK
237 void ParameterListCallbackBlocked::buildCoordinatesEpetra(
const std::string & field,
const bool useAux)
239 std::map<std::string,Kokkos::DynRankView<double,PHX::Device> > data;
243 std::vector<std::string> elementBlocks;
245 aux_blocked_ugi_->getElementBlockIds(elementBlocks);
247 blocked_ugi_->getElementBlockIds(elementBlocks);
248 for(std::size_t i=0;i<elementBlocks.size();++i) {
249 std::string blockId = elementBlocks[i];
250 std::vector<std::size_t> localCellIds;
253 Kokkos::DynRankView<double,PHX::Device> & fieldData = data[blockId];
254 fieldData = Kokkos::DynRankView<double,PHX::Device>(
"fieldData",connManager_->getElementBlock(blockId).size(),fieldPattern->
numberIds());
258 connManager_->getDofCoords(blockId,*fieldPattern,localCellIds,fieldData);
262 out.setOutputToRootOnly(-1);
263 out <<
"WARNING: In ParameterListCallback::buildCoordinates(), the Intrepid2::FieldPattern in "
264 <<
"block \"" << blockId <<
"\" does not support interpolatory coordinates. "
265 <<
"This may be fine if coordinates are not actually needed. However if they are then bad things "
266 <<
"will happen. Enjoy!" << std::endl;
272 coordsVecEp_ = arrayToVectorEpetra_[
field]->template getDataVector<double>(
field,data);
276 std::string ParameterListCallbackBlocked::
280 if(pl.
isType<std::string>(
"x-coordinates"))
281 return pl.
get<std::string>(
"x-coordinates");
282 else if(pl.
isType<std::string>(
"Coordinates"))
283 return pl.
get<std::string>(
"Coordinates");
284 #ifdef PANZER_HAVE_EPETRA_STACK
285 else if(pl.
isType<std::string>(
"Coordinates-Epetra"))
286 return pl.
get<std::string>(
"Coordinates-Epetra");
289 #ifdef PANZER_HAVE_EPETRA_STACK
290 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::logic_error,
"Neither x-coordinates nor Coordinates or Coordinates-Epetra field provided.");
296 const std::vector<double> & ParameterListCallbackBlocked::
297 getCoordinateByField(
int dim,
const std::string & field)
const
303 const std::map<std::string,std::vector<double> > * coord =
nullptr;
304 if(dim==0) coord = &xcoords_;
305 else if(dim==1) coord = &ycoords_;
306 else if(dim==2) coord = &zcoords_;
309 "ParameterListCallbackBlocked::getCoordinateByField: dim is not in range of [0,2] for field \""
313 std::map<std::string,std::vector<double> >::const_iterator itr;
314 itr = coord->find(field);
317 "ParameterListCallbackBlocked::getCoordinateByField: Coordinates for field \"" + field +
318 "\" dimension " << dim <<
" have not been built!");
324 ::getFieldPattern(
const std::string & fieldName,
const bool useAux)
const
326 std::vector<std::string> elementBlocks;
328 aux_blocked_ugi_->getElementBlockIds(elementBlocks);
330 blocked_ugi_->getElementBlockIds(elementBlocks);
332 for(std::size_t e=0;e<elementBlocks.size();e++) {
333 std::string blockId = elementBlocks[e];
335 if(blocked_ugi_->fieldInBlock(fieldName,blockId))
336 return Teuchos::rcp_dynamic_cast<const panzer::Intrepid2FieldPattern>(blocked_ugi_->getFieldPattern(blockId,fieldName),
true);
338 if(aux_blocked_ugi_->fieldInBlock(fieldName,blockId))
339 return Teuchos::rcp_dynamic_cast<const panzer::Intrepid2FieldPattern>(aux_blocked_ugi_->getFieldPattern(blockId,fieldName),
true);
342 return Teuchos::null;
ConstIterator end() const
T & get(const std::string &name, T def_value)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
T & getValue(T *ptr) const
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
const ParameterEntry & entry(ConstIterator i) const
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral...
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)