Panzer  Version of the Day
Panzer_ScatterDirichletResidual_BlockedTpetra_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_BLOCEDEPETRA_IMPL_HPP
44 #define PANZER_SCATTER_DIRICHLET_RESIDUAL_BLOCEDEPETRA_IMPL_HPP
45 
46 #include "Teuchos_RCP.hpp"
47 #include "Teuchos_Assert.hpp"
48 
49 #include "Phalanx_DataLayout.hpp"
50 
51 #include "Epetra_Map.h"
52 #include "Epetra_Vector.h"
53 #include "Epetra_CrsMatrix.h"
54 
57 #include "Panzer_PureBasis.hpp"
59 
60 #include "Phalanx_DataLayout_MDALayout.hpp"
61 
62 #include "Thyra_SpmdVectorBase.hpp"
63 #include "Thyra_ProductVectorBase.hpp"
64 #include "Thyra_BlockedLinearOpBase.hpp"
65 #include "Thyra_get_Epetra_Operator.hpp"
66 
67 #include "Teuchos_FancyOStream.hpp"
68 
69 #include <unordered_map>
70 
71 template <typename EvalT,typename TRAITS,typename LO,typename GO,typename NodeT>
74  const Teuchos::ParameterList& p)
75 {
76  std::string scatterName = p.get<std::string>("Scatter Name");
77  Teuchos::RCP<PHX::FieldTag> scatterHolder =
78  Teuchos::rcp(new PHX::Tag<ScalarT>(scatterName,Teuchos::rcp(new PHX::MDALayout<Dummy>(0))));
79 
80  // get names to be evaluated
81  const std::vector<std::string>& names =
82  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Dependent Names"));
83 
84  Teuchos::RCP<PHX::DataLayout> dl =
85  p.get< Teuchos::RCP<panzer::PureBasis> >("Basis")->functional;
86 
87  // build the vector of fields that this is dependent on
88  for (std::size_t eq = 0; eq < names.size(); ++eq) {
89  PHX::MDField<const ScalarT,Cell,NODE> scatterField = PHX::MDField<const ScalarT,Cell,NODE>(names[eq],dl);
90 
91  // tell the field manager that we depend on this field
92  this->addDependentField(scatterField.fieldTag());
93  }
94 
95  // this is what this evaluator provides
96  this->addEvaluatedField(*scatterHolder);
97 
98  this->setName(scatterName+" Scatter Residual");
99 }
100 
101 // **********************************************************************
102 // Specialization: Residual
103 // **********************************************************************
104 
105 
106 template <typename TRAITS,typename LO,typename GO,typename NodeT>
109  const Teuchos::ParameterList& p)
110  : globalIndexer_(indexer)
111  , globalDataKey_("Residual Scatter Container")
112 {
113  std::string scatterName = p.get<std::string>("Scatter Name");
114  scatterHolder_ =
115  Teuchos::rcp(new PHX::Tag<ScalarT>(scatterName,Teuchos::rcp(new PHX::MDALayout<Dummy>(0))));
116 
117  // get names to be evaluated
118  const std::vector<std::string>& names =
119  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Dependent Names"));
120 
121  // grab map from evaluated names to field names
122  fieldMap_ = p.get< Teuchos::RCP< std::map<std::string,std::string> > >("Dependent Map");
123 
124  // determine if we are scattering an initial condition
125  scatterIC_ = p.isParameter("Scatter Initial Condition") ? p.get<bool>("Scatter Initial Condition") : false;
126 
127  Teuchos::RCP<PHX::DataLayout> dl = (!scatterIC_) ?
128  p.get< Teuchos::RCP<panzer::PureBasis> >("Basis")->functional :
129  p.get< Teuchos::RCP<const panzer::PureBasis> >("Basis")->functional;
130  if (!scatterIC_) {
131  side_subcell_dim_ = p.get<int>("Side Subcell Dimension");
132  local_side_id_ = p.get<int>("Local Side ID");
133  }
134 
135  // build the vector of fields that this is dependent on
136  scatterFields_.resize(names.size());
137  for (std::size_t eq = 0; eq < names.size(); ++eq) {
138  scatterFields_[eq] = PHX::MDField<const ScalarT,Cell,NODE>(names[eq],dl);
139 
140  // tell the field manager that we depend on this field
141  this->addDependentField(scatterFields_[eq]);
142  }
143 
144  checkApplyBC_ = p.isParameter("Check Apply BC") ? p.get<bool>("Check Apply BC") : false;
145  if (checkApplyBC_) {
146  applyBC_.resize(names.size());
147  for (std::size_t eq = 0; eq < names.size(); ++eq) {
148  applyBC_[eq] = PHX::MDField<const bool,Cell,NODE>(std::string("APPLY_BC_")+fieldMap_->find(names[eq])->second,dl);
149  this->addDependentField(applyBC_[eq]);
150  }
151  }
152 
153  // this is what this evaluator provides
154  this->addEvaluatedField(*scatterHolder_);
155 
156  if (p.isType<std::string>("Global Data Key"))
157  globalDataKey_ = p.get<std::string>("Global Data Key");
158 
159  this->setName(scatterName+" Scatter Residual");
160 }
161 
162 // **********************************************************************
163 template <typename TRAITS,typename LO,typename GO,typename NodeT>
165 postRegistrationSetup(typename TRAITS::SetupData d,
167 {
168  fieldIds_.resize(scatterFields_.size());
169  // load required field numbers for fast use
170  for(std::size_t fd=0;fd<scatterFields_.size();++fd) {
171  // get field ID from DOF manager
172  std::string fieldName = fieldMap_->find(scatterFields_[fd].fieldTag().name())->second;
173  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
174 
175  // fill field data object
176  this->utils.setFieldData(scatterFields_[fd],fm);
177 
178  if (checkApplyBC_)
179  this->utils.setFieldData(applyBC_[fd],fm);
180  }
181 
182  // get the number of nodes (Should be renamed basis)
183  num_nodes = scatterFields_[0].dimension(1);
184 }
185 
186 // **********************************************************************
187 template <typename TRAITS,typename LO,typename GO,typename NodeT>
189 preEvaluate(typename TRAITS::PreEvalData d)
190 {
191  // extract dirichlet counter from container
192  Teuchos::RCP<const ContainerType> blockContainer
193  = Teuchos::rcp_dynamic_cast<ContainerType>(d.gedc.getDataObject("Dirichlet Counter"),true);
194 
195  dirichletCounter_ = Teuchos::rcp_dynamic_cast<Thyra::ProductVectorBase<double> >(blockContainer->get_f(),true);
196  TEUCHOS_ASSERT(!Teuchos::is_null(dirichletCounter_));
197 
198  // extract linear object container
199  blockedContainer_ = Teuchos::rcp_dynamic_cast<const ContainerType>(d.gedc.getDataObject(globalDataKey_),true);
200  TEUCHOS_ASSERT(!Teuchos::is_null(blockedContainer_));
201 }
202 
203 // **********************************************************************
204 template <typename TRAITS,typename LO,typename GO,typename NodeT>
206 evaluateFields(typename TRAITS::EvalData workset)
207 {
208  using Teuchos::RCP;
209  using Teuchos::ArrayRCP;
210  using Teuchos::ptrFromRef;
211  using Teuchos::rcp_dynamic_cast;
212 
213  using Thyra::VectorBase;
214  using Thyra::SpmdVectorBase;
216 
217  Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
218  out.setShowProcRank(true);
219  out.setOutputToRootOnly(-1);
220 
221  std::vector<std::pair<int,GO> > GIDs;
222  std::vector<LO> LIDs;
223 
224  // for convenience pull out some objects from workset
225  std::string blockId = this->wda(workset).block_id;
226  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
227 
228  RCP<ProductVectorBase<double> > r = (!scatterIC_) ?
229  rcp_dynamic_cast<ProductVectorBase<double> >(blockedContainer_->get_f(),true) :
230  rcp_dynamic_cast<ProductVectorBase<double> >(blockedContainer_->get_x(),true);
231 
232  // NOTE: A reordering of these loops will likely improve performance
233  // The "getGIDFieldOffsets may be expensive. However the
234  // "getElementGIDs" can be cheaper. However the lookup for LIDs
235  // may be more expensive!
236 
237  // scatter operation for each cell in workset
238  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
239  std::size_t cellLocalId = localCellIds[worksetCellIndex];
240 
241  globalIndexer_->getElementGIDs(cellLocalId,GIDs);
242 
243  // caculate the local IDs for this element
244  LIDs.resize(GIDs.size());
245  for(std::size_t i=0;i<GIDs.size();i++) {
246  // used for doing local ID lookups
247  RCP<const MapType> r_map = blockedContainer_->getMapForBlock(GIDs[i].first);
248 
249  LIDs[i] = r_map->getLocalElement(GIDs[i].second);
250  }
251 
252  // loop over each field to be scattered
253  Teuchos::ArrayRCP<double> local_r;
254  Teuchos::ArrayRCP<double> local_dc;
255  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++) {
256  int fieldNum = fieldIds_[fieldIndex];
257  int indexerId = globalIndexer_->getFieldBlock(fieldNum);
258 
259  RCP<SpmdVectorBase<double> > dc = rcp_dynamic_cast<SpmdVectorBase<double> >(dirichletCounter_->getNonconstVectorBlock(indexerId));
260  dc->getNonconstLocalData(ptrFromRef(local_dc));
261 
262  // grab local data for inputing
263  RCP<SpmdVectorBase<double> > block_r = rcp_dynamic_cast<SpmdVectorBase<double> >(r->getNonconstVectorBlock(indexerId));
264  block_r->getNonconstLocalData(ptrFromRef(local_r));
265 
266  if (!scatterIC_) {
267  // this call "should" get the right ordering according to the Intrepid2 basis
268  const std::pair<std::vector<int>,std::vector<int> > & indicePair
269  = globalIndexer_->getGIDFieldOffsets_closure(blockId,fieldNum, side_subcell_dim_, local_side_id_);
270  const std::vector<int> & elmtOffset = indicePair.first;
271  const std::vector<int> & basisIdMap = indicePair.second;
272 
273  // loop over basis functions
274  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
275  int offset = elmtOffset[basis];
276  int lid = LIDs[offset];
277  if(lid<0) // not on this processor!
278  continue;
279 
280  int basisId = basisIdMap[basis];
281 
282  if (checkApplyBC_)
283  if (!applyBC_[fieldIndex](worksetCellIndex,basisId))
284  continue;
285 
286  local_r[lid] = (scatterFields_[fieldIndex])(worksetCellIndex,basisId);
287 
288  // record that you set a dirichlet condition
289  local_dc[lid] = 1.0;
290  }
291  } else {
292  // this call "should" get the right ordering according to the Intrepid2 basis
293  const std::vector<int> & elmtOffset = globalIndexer_->getGIDFieldOffsets(blockId,fieldNum);
294 
295  // loop over basis functions
296  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
297  int offset = elmtOffset[basis];
298  int lid = LIDs[offset];
299  if(lid<0) // not on this processor!
300  continue;
301 
302  local_r[lid] = (scatterFields_[fieldIndex])(worksetCellIndex,basis);
303 
304  // record that you set a dirichlet condition
305  local_dc[lid] = 1.0;
306  }
307  }
308  }
309  }
310 }
311 
312 // **********************************************************************
313 // Specialization: Jacobian
314 // **********************************************************************
315 
316 template <typename TRAITS,typename LO,typename GO,typename NodeT>
319  const Teuchos::ParameterList& p)
320  : globalIndexer_(indexer)
321  , globalDataKey_("Residual Scatter Container")
322 {
323  std::string scatterName = p.get<std::string>("Scatter Name");
324  scatterHolder_ =
325  Teuchos::rcp(new PHX::Tag<ScalarT>(scatterName,Teuchos::rcp(new PHX::MDALayout<Dummy>(0))));
326 
327  // get names to be evaluated
328  const std::vector<std::string>& names =
329  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Dependent Names"));
330 
331  // grab map from evaluated names to field names
332  fieldMap_ = p.get< Teuchos::RCP< std::map<std::string,std::string> > >("Dependent Map");
333 
334  Teuchos::RCP<PHX::DataLayout> dl =
335  p.get< Teuchos::RCP<panzer::PureBasis> >("Basis")->functional;
336 
337  side_subcell_dim_ = p.get<int>("Side Subcell Dimension");
338  local_side_id_ = p.get<int>("Local Side ID");
339 
340  // build the vector of fields that this is dependent on
341  scatterFields_.resize(names.size());
342  for (std::size_t eq = 0; eq < names.size(); ++eq) {
343  scatterFields_[eq] = PHX::MDField<const ScalarT,Cell,NODE>(names[eq],dl);
344 
345  // tell the field manager that we depend on this field
346  this->addDependentField(scatterFields_[eq]);
347  }
348 
349  checkApplyBC_ = p.get<bool>("Check Apply BC");
350  if (checkApplyBC_) {
351  applyBC_.resize(names.size());
352  for (std::size_t eq = 0; eq < names.size(); ++eq) {
353  applyBC_[eq] = PHX::MDField<const bool,Cell,NODE>(std::string("APPLY_BC_")+fieldMap_->find(names[eq])->second,dl);
354  this->addDependentField(applyBC_[eq]);
355  }
356  }
357 
358  // this is what this evaluator provides
359  this->addEvaluatedField(*scatterHolder_);
360 
361  if (p.isType<std::string>("Global Data Key"))
362  globalDataKey_ = p.get<std::string>("Global Data Key");
363 
364  this->setName(scatterName+" Scatter Residual (Jacobian)");
365 }
366 
367 // **********************************************************************
368 template <typename TRAITS,typename LO,typename GO,typename NodeT>
370 postRegistrationSetup(typename TRAITS::SetupData d,
372 {
373  fieldIds_.resize(scatterFields_.size());
374  // load required field numbers for fast use
375  for(std::size_t fd=0;fd<scatterFields_.size();++fd) {
376  // get field ID from DOF manager
377  std::string fieldName = fieldMap_->find(scatterFields_[fd].fieldTag().name())->second;
378  fieldIds_[fd] = globalIndexer_->getFieldNum(fieldName);
379 
380  // fill field data object
381  this->utils.setFieldData(scatterFields_[fd],fm);
382 
383  if (checkApplyBC_)
384  this->utils.setFieldData(applyBC_[fd],fm);
385  }
386 
387  // get the number of nodes (Should be renamed basis)
388  num_nodes = scatterFields_[0].dimension(1);
389  num_eq = scatterFields_.size();
390 }
391 
392 // **********************************************************************
393 template <typename TRAITS,typename LO,typename GO,typename NodeT>
395 preEvaluate(typename TRAITS::PreEvalData d)
396 {
397  // extract dirichlet counter from container
398  Teuchos::RCP<const ContainerType> blockContainer
399  = Teuchos::rcp_dynamic_cast<const ContainerType>(d.gedc.getDataObject("Dirichlet Counter"),true);
400 
401  dirichletCounter_ = Teuchos::rcp_dynamic_cast<Thyra::ProductVectorBase<double> >(blockContainer->get_f(),true);
402  TEUCHOS_ASSERT(!Teuchos::is_null(dirichletCounter_));
403 
404  // extract linear object container
405  blockedContainer_ = Teuchos::rcp_dynamic_cast<const ContainerType>(d.gedc.getDataObject(globalDataKey_),true);
406  TEUCHOS_ASSERT(!Teuchos::is_null(blockedContainer_));
407 }
408 
409 // **********************************************************************
410 template <typename TRAITS,typename LO,typename GO,typename NodeT>
412 evaluateFields(typename TRAITS::EvalData workset)
413 {
414  using Teuchos::RCP;
415  using Teuchos::ArrayRCP;
416  using Teuchos::ptrFromRef;
417  using Teuchos::rcp_dynamic_cast;
418 
419  using Thyra::VectorBase;
420  using Thyra::SpmdVectorBase;
423 
424  std::vector<std::pair<int,GO> > GIDs;
425  std::vector<LO> LIDs;
426 
427  // for convenience pull out some objects from workset
428  std::string blockId = this->wda(workset).block_id;
429  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
430 
431  RCP<ProductVectorBase<double> > r = rcp_dynamic_cast<ProductVectorBase<double> >(blockedContainer_->get_f());
432  Teuchos::RCP<BlockedLinearOpBase<double> > Jac = rcp_dynamic_cast<BlockedLinearOpBase<double> >(blockedContainer_->get_A());
433 
434  int numFieldBlocks = globalIndexer_->getNumFieldBlocks();
435  std::vector<int> blockOffsets(numFieldBlocks+1); // number of fields, plus a sentinnel
436  for(int blk=0;blk<numFieldBlocks;blk++) {
437  int blockOffset = globalIndexer_->getBlockGIDOffset(blockId,blk);
438  blockOffsets[blk] = blockOffset;
439  }
440 
441  std::unordered_map<std::pair<int,int>,Teuchos::RCP<CrsMatrixType>,panzer::pair_hash> jacTpetraBlocks;
442 
443  // NOTE: A reordering of these loops will likely improve performance
444  // The "getGIDFieldOffsets may be expensive. However the
445  // "getElementGIDs" can be cheaper. However the lookup for LIDs
446  // may be more expensive!
447 
448  // scatter operation for each cell in workset
449  for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
450  std::size_t cellLocalId = localCellIds[worksetCellIndex];
451 
452  globalIndexer_->getElementGIDs(cellLocalId,GIDs);
453  blockOffsets[numFieldBlocks] = GIDs.size();
454 
455  // caculate the local IDs for this element
456  LIDs.resize(GIDs.size());
457  for(std::size_t i=0;i<GIDs.size();i++) {
458  // used for doing local ID lookups
459  RCP<const MapType> r_map = blockedContainer_->getMapForBlock(GIDs[i].first);
460 
461  LIDs[i] = r_map->getLocalElement(GIDs[i].second);
462  }
463 
464  // loop over each field to be scattered
465  Teuchos::ArrayRCP<double> local_r, local_dc;
466  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++) {
467  int fieldNum = fieldIds_[fieldIndex];
468  int blockRowIndex = globalIndexer_->getFieldBlock(fieldNum);
469 
470  RCP<SpmdVectorBase<double> > dc = rcp_dynamic_cast<SpmdVectorBase<double> >(dirichletCounter_->getNonconstVectorBlock(blockRowIndex));
471  dc->getNonconstLocalData(ptrFromRef(local_dc));
472 
473  // grab local data for inputing
474  RCP<SpmdVectorBase<double> > block_r = rcp_dynamic_cast<SpmdVectorBase<double> >(r->getNonconstVectorBlock(blockRowIndex));
475  block_r->getNonconstLocalData(ptrFromRef(local_r));
476 
477  // this call "should" get the right ordering according to the Intrepid2 basis
478  const std::pair<std::vector<int>,std::vector<int> > & indicePair
479  = globalIndexer_->getGIDFieldOffsets_closure(blockId,fieldNum, side_subcell_dim_, local_side_id_);
480  const std::vector<int> & elmtOffset = indicePair.first;
481  const std::vector<int> & basisIdMap = indicePair.second;
482 
483  // loop over basis functions
484  for(std::size_t basis=0;basis<elmtOffset.size();basis++) {
485  int offset = elmtOffset[basis];
486  int lid = LIDs[offset];
487  if(lid<0) // not on this processor
488  continue;
489 
490  int basisId = basisIdMap[basis];
491 
492  if (checkApplyBC_)
493  if (!applyBC_[fieldIndex](worksetCellIndex,basisId))
494  continue;
495 
496  // zero out matrix row
497  for(int blockColIndex=0;blockColIndex<numFieldBlocks;blockColIndex++) {
498  int start = blockOffsets[blockColIndex];
499  int end = blockOffsets[blockColIndex+1];
500 
501  if(end-start<=0)
502  continue;
503 
504  // check hash table for jacobian sub block
505  std::pair<int,int> blockIndex = std::make_pair(blockRowIndex,blockColIndex);
506  Teuchos::RCP<CrsMatrixType> subJac = jacTpetraBlocks[blockIndex];
507 
508  // if you didn't find one before, add it to the hash table
509  if(subJac==Teuchos::null) {
510  Teuchos::RCP<Thyra::LinearOpBase<double> > tOp = Jac->getNonconstBlock(blockIndex.first,blockIndex.second);
511 
512  // block operator is null, don't do anything (it is excluded)
513  if(Teuchos::is_null(tOp))
514  continue;
515 
516  Teuchos::RCP<OperatorType> tpetra_Op = rcp_dynamic_cast<ThyraLinearOp>(tOp)->getTpetraOperator();
517  subJac = rcp_dynamic_cast<CrsMatrixType>(tpetra_Op,true);
518  jacTpetraBlocks[blockIndex] = subJac;
519  }
520 
521  std::size_t sz = subJac->getNumEntriesInLocalRow(lid);
522  std::size_t numEntries = 0;
523  Teuchos::Array<LO> rowIndices(sz);
524  Teuchos::Array<double> rowValues(sz);
525 
526  subJac->getLocalRowCopy(lid,rowIndices,rowValues,numEntries);
527 
528  for(std::size_t i=0;i<numEntries;i++)
529  rowValues[i] = 0.0;
530 
531  subJac->replaceLocalValues(lid,rowIndices,rowValues);
532  }
533 
534  const ScalarT scatterField = (scatterFields_[fieldIndex])(worksetCellIndex,basisId);
535 
536  local_r[lid] = scatterField.val();
537  local_dc[lid] = 1.0; // mark row as dirichlet
538 
539  // loop over the sensitivity indices: all DOFs on a cell
540  std::vector<double> jacRow(scatterField.size(),0.0);
541 
542  for(int sensIndex=0;sensIndex<scatterField.size();++sensIndex)
543  jacRow[sensIndex] = scatterField.fastAccessDx(sensIndex);
544  TEUCHOS_ASSERT(jacRow.size()==GIDs.size());
545 
546  for(int blockColIndex=0;blockColIndex<numFieldBlocks;blockColIndex++) {
547  int start = blockOffsets[blockColIndex];
548  int end = blockOffsets[blockColIndex+1];
549 
550  if(end-start<=0)
551  continue;
552 
553  // check hash table for jacobian sub block
554  std::pair<int,int> blockIndex = std::make_pair(blockRowIndex,blockColIndex);
555  Teuchos::RCP<CrsMatrixType> subJac = jacTpetraBlocks[blockIndex];
556 
557  // if you didn't find one before, add it to the hash table
558  if(subJac==Teuchos::null) {
559  Teuchos::RCP<Thyra::LinearOpBase<double> > tOp = Jac->getNonconstBlock(blockIndex.first,blockIndex.second);
560 
561  // block operator is null, don't do anything (it is excluded)
562  if(Teuchos::is_null(tOp))
563  continue;
564 
565  Teuchos::RCP<OperatorType> tpetra_Op = rcp_dynamic_cast<ThyraLinearOp>(tOp)->getTpetraOperator();
566  subJac = rcp_dynamic_cast<CrsMatrixType>(tpetra_Op,true);
567  jacTpetraBlocks[blockIndex] = subJac;
568  }
569 
570  // Sum Jacobian
571  subJac->replaceLocalValues(lid, Teuchos::arrayViewFromVector(LIDs).view(start,end-start),
572  Teuchos::arrayViewFromVector(jacRow).view(start,end-start));
573  }
574  }
575  }
576  }
577 }
578 
579 // **********************************************************************
580 
581 #endif
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::Point > > scatterFields_
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
Pushes residual values into the residual vector for a Newton-based solve.
ScatterDirichletResidual_BlockedTpetra(const Teuchos::RCP< const BlockedDOFManager< LO, GO > > &indexer)
void postRegistrationSetup(typename TRAITS::SetupData d, PHX::FieldManager< TRAITS > &vm)