Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ScatterDirichletResidual_Tpetra_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_SCATTER_DIRICHLET_RESIDUAL_TPETRA_IMPL_HPP
44 #define PANZER_SCATTER_DIRICHLET_RESIDUAL_TPETRA_IMPL_HPP
45 
46 #include "Teuchos_RCP.hpp"
47 #include "Teuchos_Assert.hpp"
48 
49 #include "Phalanx_DataLayout.hpp"
50 
51 #include "Panzer_GlobalIndexer.hpp"
52 #include "Panzer_PureBasis.hpp"
57 
58 #include "Phalanx_DataLayout_MDALayout.hpp"
59 
60 #include "Teuchos_FancyOStream.hpp"
61 
62 // **********************************************************************
63 // Specialization: Residual
64 // **********************************************************************
65 
66 
67 template<typename TRAITS,typename LO,typename GO,typename NodeT>
70  const Teuchos::ParameterList& p)
71  : globalIndexer_(indexer)
72  , globalDataKey_("Residual Scatter Container")
73 {
74  std::string scatterName = p.get<std::string>("Scatter Name");
75  scatterHolder_ =
76  Teuchos::rcp(new PHX::Tag<ScalarT>(scatterName,Teuchos::rcp(new PHX::MDALayout<Dummy>(0))));
77 
78  // get names to be evaluated
79  const std::vector<std::string>& names =
80  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Dependent Names"));
81 
82  // grab map from evaluated names to field names
83  fieldMap_ = p.get< Teuchos::RCP< std::map<std::string,std::string> > >("Dependent Map");
84 
85  // determine if we are scattering an initial condition
86  scatterIC_ = p.isParameter("Scatter Initial Condition") ? p.get<bool>("Scatter Initial Condition") : false;
87 
88  Teuchos::RCP<PHX::DataLayout> dl = (!scatterIC_) ?
89  p.get< Teuchos::RCP<panzer::PureBasis> >("Basis")->functional :
90  p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis")->functional ;
91  if (!scatterIC_) {
92  side_subcell_dim_ = p.get<int>("Side Subcell Dimension");
93  local_side_id_ = p.get<int>("Local Side ID");
94  }
95 
96  // build the vector of fields that this is dependent on
97  scatterFields_.resize(names.size());
98  for (std::size_t eq = 0; eq < names.size(); ++eq) {
99  scatterFields_[eq] = PHX::MDField<const ScalarT,Cell,NODE>(names[eq],dl);
100 
101  // tell the field manager that we depend on this field
102  this->addDependentField(scatterFields_[eq]);
103  }
104 
105  checkApplyBC_ = p.isParameter("Check Apply BC") ? p.get<bool>("Check Apply BC") : false;
106  if (checkApplyBC_) {
107  applyBC_.resize(names.size());
108  for (std::size_t eq = 0; eq < names.size(); ++eq) {
109  applyBC_[eq] = PHX::MDField<const bool,Cell,NODE>(std::string("APPLY_BC_")+fieldMap_->find(names[eq])->second,dl);
110  this->addDependentField(applyBC_[eq]);
111  }
112  }
113 
114  // this is what this evaluator provides
115  this->addEvaluatedField(*scatterHolder_);
116 
117  if (p.isType<std::string>("Global Data Key"))
118  globalDataKey_ = p.get<std::string>("Global Data Key");
119 
120  this->setName(scatterName+" Scatter Residual");
121 }
122 
123 // **********************************************************************
124 template<typename TRAITS,typename LO,typename GO,typename NodeT>
126 postRegistrationSetup(typename TRAITS::SetupData /* d */,
127  PHX::FieldManager<TRAITS>& /* fm */)
128 {
129  fieldIds_.resize(scatterFields_.size());
130 
131  // load required field numbers for fast use
132  for(std::size_t fd=0;fd<scatterFields_.size();++fd) {
133  // get field ID from DOF manager
134  std::string fieldName = fieldMap_->find(scatterFields_[fd].fieldTag().name())->second;
135  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
136  }
137 
138  // get the number of nodes (Should be renamed basis)
139  num_nodes = scatterFields_[0].extent(1);
140 }
141 
142 // **********************************************************************
143 template<typename TRAITS,typename LO,typename GO,typename NodeT>
145 preEvaluate(typename TRAITS::PreEvalData d)
146 {
147  // extract linear object container
148  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject(globalDataKey_));
149 
150  if(tpetraContainer_==Teuchos::null) {
151  // extract linear object container
152  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
153  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(loc);
154 
155  dirichletCounter_ = Teuchos::null;
156  }
157  else {
158  // extract dirichlet counter from container
159  Teuchos::RCP<LOC> tpetraContainer
160  = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject("Dirichlet Counter"),true);
161 
162  dirichletCounter_ = tpetraContainer->get_f();
163  TEUCHOS_ASSERT(!Teuchos::is_null(dirichletCounter_));
164  }
165 }
166 
167 // **********************************************************************
168 template<typename TRAITS,typename LO,typename GO,typename NodeT>
170 evaluateFields(typename TRAITS::EvalData workset)
171 {
172  std::vector<GO> GIDs;
173  std::vector<LO> LIDs;
174 
175  // for convenience pull out some objects from workset
176  std::string blockId = this->wda(workset).block_id;
177  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
178 
179  Teuchos::RCP<typename LOC::VectorType> r = (!scatterIC_) ?
180  tpetraContainer_->get_f() :
181  tpetraContainer_->get_x();
182 
183  Teuchos::ArrayRCP<double> r_array = r->get1dViewNonConst();
184  Teuchos::ArrayRCP<double> dc_array = dirichletCounter_->get1dViewNonConst();
185 
186  // NOTE: A reordering of these loops will likely improve performance
187  // The "getGIDFieldOffsets may be expensive. However the
188  // "getElementGIDs" can be cheaper. However the lookup for LIDs
189  // may be more expensive!
190 
191 
192  // scatter operation for each cell in workset
193  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
194  std::size_t cellLocalId = localCellIds[worksetCellIndex];
195 
196  globalIndexer_->getElementGIDs(cellLocalId,GIDs);
197 
198  // caculate the local IDs for this element
199  LIDs.resize(GIDs.size());
200  for(std::size_t i=0;i<GIDs.size();i++)
201  LIDs[i] = r->getMap()->getLocalElement(GIDs[i]);
202 
203  // loop over each field to be scattered
204  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++) {
205  int fieldNum = fieldIds_[fieldIndex];
206 
207  if (!scatterIC_) {
208  // this call "should" get the right ordering according to the Intrepid2 basis
209  const std::pair<std::vector<int>,std::vector<int> > & indicePair
210  = globalIndexer_->getGIDFieldOffsets_closure(blockId,fieldNum, side_subcell_dim_, local_side_id_);
211  const std::vector<int> & elmtOffset = indicePair.first;
212  const std::vector<int> & basisIdMap = indicePair.second;
213 
214  // loop over basis functions
215  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
216  int offset = elmtOffset[basis];
217  LO lid = LIDs[offset];
218  if(lid<0) // not on this processor!
219  continue;
220 
221  int basisId = basisIdMap[basis];
222 
223  if (checkApplyBC_)
224  if (!applyBC_[fieldIndex](worksetCellIndex,basisId))
225  continue;
226 
227  r_array[lid] = (scatterFields_[fieldIndex])(worksetCellIndex,basisId);
228 
229  // record that you set a dirichlet condition
230  dc_array[lid] = 1.0;
231  }
232  } else {
233  // this call "should" get the right ordering according to the Intrepid2 basis
234  const std::vector<int> & elmtOffset = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
235 
236  // loop over basis functions
237  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
238  int offset = elmtOffset[basis];
239  LO lid = LIDs[offset];
240  if(lid<0) // not on this processor!
241  continue;
242 
243  r_array[lid] = (scatterFields_[fieldIndex])(worksetCellIndex,basis);
244 
245  // record that you set a dirichlet condition
246  dc_array[lid] = 1.0;
247  }
248  }
249  }
250  }
251 }
252 
253 // **********************************************************************
254 // Specialization: Tangent
255 // **********************************************************************
256 
257 
258 template<typename TRAITS,typename LO,typename GO,typename NodeT>
261  const Teuchos::ParameterList& p)
262  : globalIndexer_(indexer)
263  , globalDataKey_("Residual Scatter Container")
264 {
265  std::string scatterName = p.get<std::string>("Scatter Name");
266  scatterHolder_ =
267  Teuchos::rcp(new PHX::Tag<ScalarT>(scatterName,Teuchos::rcp(new PHX::MDALayout<Dummy>(0))));
268 
269  // get names to be evaluated
270  const std::vector<std::string>& names =
271  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Dependent Names"));
272 
273  // grab map from evaluated names to field names
274  fieldMap_ = p.get< Teuchos::RCP< std::map<std::string,std::string> > >("Dependent Map");
275 
276  // determine if we are scattering an initial condition
277  scatterIC_ = p.isParameter("Scatter Initial Condition") ? p.get<bool>("Scatter Initial Condition") : false;
278 
279  Teuchos::RCP<PHX::DataLayout> dl = (!scatterIC_) ?
280  p.get< Teuchos::RCP<panzer::PureBasis> >("Basis")->functional :
281  p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis")->functional ;
282  if (!scatterIC_) {
283  side_subcell_dim_ = p.get<int>("Side Subcell Dimension");
284  local_side_id_ = p.get<int>("Local Side ID");
285  }
286 
287  // build the vector of fields that this is dependent on
288  scatterFields_.resize(names.size());
289  for (std::size_t eq = 0; eq < names.size(); ++eq) {
290  scatterFields_[eq] = PHX::MDField<const ScalarT,Cell,NODE>(names[eq],dl);
291 
292  // tell the field manager that we depend on this field
293  this->addDependentField(scatterFields_[eq]);
294  }
295 
296  checkApplyBC_ = p.isParameter("Check Apply BC") ? p.get<bool>("Check Apply BC") : false;
297  if (checkApplyBC_) {
298  applyBC_.resize(names.size());
299  for (std::size_t eq = 0; eq < names.size(); ++eq) {
300  applyBC_[eq] = PHX::MDField<const bool,Cell,NODE>(std::string("APPLY_BC_")+fieldMap_->find(names[eq])->second,dl);
301  this->addDependentField(applyBC_[eq]);
302  }
303  }
304 
305  // this is what this evaluator provides
306  this->addEvaluatedField(*scatterHolder_);
307 
308  if (p.isType<std::string>("Global Data Key"))
309  globalDataKey_ = p.get<std::string>("Global Data Key");
310 
311  this->setName(scatterName+" Scatter Tangent");
312 }
313 
314 // **********************************************************************
315 template<typename TRAITS,typename LO,typename GO,typename NodeT>
317 postRegistrationSetup(typename TRAITS::SetupData /* d */,
318  PHX::FieldManager<TRAITS>& /* fm */)
319 {
320  fieldIds_.resize(scatterFields_.size());
321 
322  // load required field numbers for fast use
323  for(std::size_t fd=0;fd<scatterFields_.size();++fd) {
324  // get field ID from DOF manager
325  std::string fieldName = fieldMap_->find(scatterFields_[fd].fieldTag().name())->second;
326  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
327  }
328 
329  // get the number of nodes (Should be renamed basis)
330  num_nodes = scatterFields_[0].extent(1);
331 }
332 
333 // **********************************************************************
334 template<typename TRAITS,typename LO,typename GO,typename NodeT>
336 preEvaluate(typename TRAITS::PreEvalData d)
337 {
338  // extract linear object container
339  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject(globalDataKey_));
340 
341  if(tpetraContainer_==Teuchos::null) {
342  // extract linear object container
343  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
344  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(loc);
345 
346  dirichletCounter_ = Teuchos::null;
347  }
348  else {
349  // extract dirichlet counter from container
350  Teuchos::RCP<LOC> tpetraContainer
351  = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject("Dirichlet Counter"),true);
352 
353  dirichletCounter_ = tpetraContainer->get_f();
354  TEUCHOS_ASSERT(!Teuchos::is_null(dirichletCounter_));
355  }
356 
357  using Teuchos::RCP;
358  using Teuchos::rcp_dynamic_cast;
359 
360  // this is the list of parameters and their names that this scatter has to account for
361  std::vector<std::string> activeParameters =
362  rcp_dynamic_cast<ParameterList_GlobalEvaluationData>(d.gedc->getDataObject("PARAMETER_NAMES"))->getActiveParameters();
363 
364  // ETP 02/03/16: This code needs to be updated to properly handle scatterIC_
365  TEUCHOS_ASSERT(!scatterIC_);
366  dfdp_vectors_.clear();
367  for(std::size_t i=0;i<activeParameters.size();i++) {
369  rcp_dynamic_cast<LOC>(d.gedc->getDataObject(activeParameters[i]),true)->get_f();
370  Teuchos::ArrayRCP<double> vec_array = vec->get1dViewNonConst();
371  dfdp_vectors_.push_back(vec_array);
372  }
373 }
374 
375 // **********************************************************************
376 template<typename TRAITS,typename LO,typename GO,typename NodeT>
378 evaluateFields(typename TRAITS::EvalData workset)
379 {
380  std::vector<GO> GIDs;
381  std::vector<LO> LIDs;
382 
383  // for convenience pull out some objects from workset
384  std::string blockId = this->wda(workset).block_id;
385  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
386 
387  Teuchos::RCP<typename LOC::VectorType> r = (!scatterIC_) ?
388  tpetraContainer_->get_f() :
389  tpetraContainer_->get_x();
390 
391  Teuchos::ArrayRCP<double> r_array = r->get1dViewNonConst();
392  Teuchos::ArrayRCP<double> dc_array = dirichletCounter_->get1dViewNonConst();
393 
394  // NOTE: A reordering of these loops will likely improve performance
395  // The "getGIDFieldOffsets may be expensive. However the
396  // "getElementGIDs" can be cheaper. However the lookup for LIDs
397  // may be more expensive!
398 
399 
400  // scatter operation for each cell in workset
401  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
402  std::size_t cellLocalId = localCellIds[worksetCellIndex];
403 
404  globalIndexer_->getElementGIDs(cellLocalId,GIDs);
405 
406  // caculate the local IDs for this element
407  LIDs.resize(GIDs.size());
408  for(std::size_t i=0;i<GIDs.size();i++)
409  LIDs[i] = r->getMap()->getLocalElement(GIDs[i]);
410 
411  // loop over each field to be scattered
412  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++) {
413  int fieldNum = fieldIds_[fieldIndex];
414 
415  if (!scatterIC_) {
416  // this call "should" get the right ordering according to the Intrepid2 basis
417  const std::pair<std::vector<int>,std::vector<int> > & indicePair
418  = globalIndexer_->getGIDFieldOffsets_closure(blockId,fieldNum, side_subcell_dim_, local_side_id_);
419  const std::vector<int> & elmtOffset = indicePair.first;
420  const std::vector<int> & basisIdMap = indicePair.second;
421 
422  // loop over basis functions
423  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
424  int offset = elmtOffset[basis];
425  LO lid = LIDs[offset];
426  if(lid<0) // not on this processor!
427  continue;
428 
429  int basisId = basisIdMap[basis];
430 
431  if (checkApplyBC_)
432  if (!applyBC_[fieldIndex](worksetCellIndex,basisId))
433  continue;
434 
435  ScalarT value = (scatterFields_[fieldIndex])(worksetCellIndex,basisId);
436  //r_array[lid] = (scatterFields_[fieldIndex])(worksetCellIndex,basisId).val();
437 
438  // then scatter the sensitivity vectors
439  if(value.size()==0)
440  for(std::size_t d=0;d<dfdp_vectors_.size();d++)
441  dfdp_vectors_[d][lid] = 0.0;
442  else
443  for(int d=0;d<value.size();d++) {
444  dfdp_vectors_[d][lid] = value.fastAccessDx(d);
445  }
446 
447  // record that you set a dirichlet condition
448  dc_array[lid] = 1.0;
449  }
450  } else {
451  // this call "should" get the right ordering according to the Intrepid2 basis
452  const std::vector<int> & elmtOffset = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
453 
454  // loop over basis functions
455  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
456  int offset = elmtOffset[basis];
457  LO lid = LIDs[offset];
458  if(lid<0) // not on this processor!
459  continue;
460 
461  ScalarT value = (scatterFields_[fieldIndex])(worksetCellIndex,basis);
462  //r_array[lid] = (scatterFields_[fieldIndex])(worksetCellIndex,basis).val();
463 
464  // then scatter the sensitivity vectors
465  if(value.size()==0)
466  for(std::size_t d=0;d<dfdp_vectors_.size();d++)
467  dfdp_vectors_[d][lid] = 0.0;
468  else
469  for(int d=0;d<value.size();d++) {
470  dfdp_vectors_[d][lid] = value.fastAccessDx(d);
471  }
472 
473  // record that you set a dirichlet condition
474  dc_array[lid] = 1.0;
475  }
476  }
477  }
478  }
479 }
480 
481 // **********************************************************************
482 // Specialization: Jacobian
483 // **********************************************************************
484 
485 template<typename TRAITS,typename LO,typename GO,typename NodeT>
488  const Teuchos::ParameterList& p)
489  : globalIndexer_(indexer)
490  , globalDataKey_("Residual Scatter Container")
491 {
492  std::string scatterName = p.get<std::string>("Scatter Name");
493  scatterHolder_ =
494  Teuchos::rcp(new PHX::Tag<ScalarT>(scatterName,Teuchos::rcp(new PHX::MDALayout<Dummy>(0))));
495 
496  // get names to be evaluated
497  const std::vector<std::string>& names =
498  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Dependent Names"));
499 
500  // grab map from evaluated names to field names
501  fieldMap_ = p.get< Teuchos::RCP< std::map<std::string,std::string> > >("Dependent Map");
502 
504  p.get< Teuchos::RCP<panzer::PureBasis> >("Basis")->functional;
505 
506  side_subcell_dim_ = p.get<int>("Side Subcell Dimension");
507  local_side_id_ = p.get<int>("Local Side ID");
508 
509  // build the vector of fields that this is dependent on
510  scatterFields_.resize(names.size());
511  for (std::size_t eq = 0; eq < names.size(); ++eq) {
512  scatterFields_[eq] = PHX::MDField<const ScalarT,Cell,NODE>(names[eq],dl);
513 
514  // tell the field manager that we depend on this field
515  this->addDependentField(scatterFields_[eq]);
516  }
517 
518  checkApplyBC_ = p.get<bool>("Check Apply BC");
519  if (checkApplyBC_) {
520  applyBC_.resize(names.size());
521  for (std::size_t eq = 0; eq < names.size(); ++eq) {
522  applyBC_[eq] = PHX::MDField<const bool,Cell,NODE>(std::string("APPLY_BC_")+fieldMap_->find(names[eq])->second,dl);
523  this->addDependentField(applyBC_[eq]);
524  }
525  }
526 
527  // this is what this evaluator provides
528  this->addEvaluatedField(*scatterHolder_);
529 
530  if (p.isType<std::string>("Global Data Key"))
531  globalDataKey_ = p.get<std::string>("Global Data Key");
532 
533  this->setName(scatterName+" Scatter Residual (Jacobian)");
534 }
535 
536 // **********************************************************************
537 template<typename TRAITS,typename LO,typename GO,typename NodeT>
539 postRegistrationSetup(typename TRAITS::SetupData /* d */,
540  PHX::FieldManager<TRAITS>& /* fm */)
541 {
542  fieldIds_.resize(scatterFields_.size());
543  // load required field numbers for fast use
544  for(std::size_t fd=0;fd<scatterFields_.size();++fd) {
545  // get field ID from DOF manager
546  std::string fieldName = fieldMap_->find(scatterFields_[fd].fieldTag().name())->second;
547  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
548  }
549 
550  // get the number of nodes (Should be renamed basis)
551  num_nodes = scatterFields_[0].extent(1);
552  num_eq = scatterFields_.size();
553 }
554 
555 // **********************************************************************
556 template<typename TRAITS,typename LO,typename GO,typename NodeT>
558 preEvaluate(typename TRAITS::PreEvalData d)
559 {
560  // extract linear object container
561  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject(globalDataKey_));
562 
563  if(tpetraContainer_==Teuchos::null) {
564  // extract linear object container
565  Teuchos::RCP<LinearObjContainer> loc = Teuchos::rcp_dynamic_cast<LOCPair_GlobalEvaluationData>(d.gedc->getDataObject(globalDataKey_),true)->getGhostedLOC();
566  tpetraContainer_ = Teuchos::rcp_dynamic_cast<LOC>(loc);
567 
568  dirichletCounter_ = Teuchos::null;
569  }
570  else {
571  // extract dirichlet counter from container
572  Teuchos::RCP<LOC> tpetraContainer
573  = Teuchos::rcp_dynamic_cast<LOC>(d.gedc->getDataObject("Dirichlet Counter"),true);
574 
575  dirichletCounter_ = tpetraContainer->get_f();
576  TEUCHOS_ASSERT(!Teuchos::is_null(dirichletCounter_));
577  }
578 }
579 
580 // **********************************************************************
581 template<typename TRAITS,typename LO,typename GO,typename NodeT>
583 evaluateFields(typename TRAITS::EvalData workset)
584 {
585  std::vector<GO> GIDs;
586 
587  // for convenience pull out some objects from workset
588  std::string blockId = this->wda(workset).block_id;
589  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
590 
591  Teuchos::RCP<typename LOC::VectorType> r = tpetraContainer_->get_f();
592  Teuchos::RCP<typename LOC::CrsMatrixType> Jac = tpetraContainer_->get_A();
593 
594  Teuchos::ArrayRCP<double> r_array = r->get1dViewNonConst();
595  Teuchos::ArrayRCP<double> dc_array = dirichletCounter_->get1dViewNonConst();
596 
597  // NOTE: A reordering of these loops will likely improve performance
598  // The "getGIDFieldOffsets may be expensive. However the
599  // "getElementGIDs" can be cheaper. However the lookup for LIDs
600  // may be more expensive!
601 
602  // scatter operation for each cell in workset
603  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
604  std::size_t cellLocalId = localCellIds[worksetCellIndex];
605 
606  globalIndexer_->getElementGIDs(cellLocalId,GIDs);
607  auto LIDs = globalIndexer_->getElementLIDs(cellLocalId);
608 
609  // loop over each field to be scattered
610  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++) {
611  int fieldNum = fieldIds_[fieldIndex];
612 
613  // this call "should" get the right ordering according to the Intrepid2 basis
614  const std::pair<std::vector<int>,std::vector<int> > & indicePair
615  = globalIndexer_->getGIDFieldOffsets_closure(blockId,fieldNum, side_subcell_dim_, local_side_id_);
616  const std::vector<int> & elmtOffset = indicePair.first;
617  const std::vector<int> & basisIdMap = indicePair.second;
618 
619  // loop over basis functions
620  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
621  int offset = elmtOffset[basis];
622  LO lid = LIDs[offset];
623  if(lid<0) // not on this processor
624  continue;
625 
626  int basisId = basisIdMap[basis];
627 
628  if (checkApplyBC_)
629  if (!applyBC_[fieldIndex](worksetCellIndex,basisId))
630  continue;
631 
632  // zero out matrix row
633  {
634  std::size_t sz = Jac->getNumEntriesInLocalRow(lid);
635  std::size_t numEntries = 0;
636  Teuchos::Array<LO> rowIndices(sz);
637  Teuchos::Array<double> rowValues(sz);
638 
639  // Jac->getLocalRowView(lid,numEntries,rowValues,rowIndices);
640  Jac->getLocalRowCopy(lid,rowIndices,rowValues,numEntries);
641 
642  for(std::size_t i=0;i<numEntries;i++)
643  rowValues[i] = 0.0;
644 
645  Jac->replaceLocalValues(lid,rowIndices,rowValues);
646  }
647 
648  GO gid = GIDs[offset];
649  const ScalarT scatterField = (scatterFields_[fieldIndex])(worksetCellIndex,basisId);
650 
651  r_array[lid] = scatterField.val();
652  dc_array[lid] = 1.0; // mark row as dirichlet
653 
654  // loop over the sensitivity indices: all DOFs on a cell
655  std::vector<double> jacRow(scatterField.size(),0.0);
656 
657  for(int sensIndex=0;sensIndex<scatterField.size();++sensIndex)
658  jacRow[sensIndex] = scatterField.fastAccessDx(sensIndex);
659  TEUCHOS_ASSERT(jacRow.size()==GIDs.size());
660 
661  Jac->replaceGlobalValues(gid, GIDs, jacRow);
662  }
663  }
664  }
665 }
666 
667 // **********************************************************************
668 
669 #endif
T & get(const std::string &name, T def_value)
bool is_null(const std::shared_ptr< T > &p)
void resize(const size_type n, const T &val=T())
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Pushes residual values into the residual vector for a Newton-based solve.
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)