Panzer  Version of the Day
Panzer_Workset.cpp
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 #include "Panzer_Workset.hpp"
44 
45 #include "Phalanx_DataLayout.hpp"
46 #include "Phalanx_DataLayout_MDALayout.hpp"
47 
49 #include "Panzer_Workset_Builder.hpp"
50 #include "Panzer_WorksetNeeds.hpp"
51 #include "Panzer_Dimension.hpp"
52 #include "Panzer_LocalMeshInfo.hpp"
54 #include "Panzer_PointValues2.hpp"
55 #include "Panzer_HashUtils.hpp"
57 
60 
61 #include "Shards_CellTopology.hpp"
62 
63 namespace panzer {
64 
65 namespace {
66 
67 void buildLocalOrientations(const int num_cells,
68  const Kokkos::View<const panzer::LocalOrdinal*,PHX::Device> & local_cell_ids,
69  const Teuchos::RCP<const OrientationsInterface> & orientations_interface,
70  std::vector<Intrepid2::Orientation> & workset_orientations)
71 {
72  // from a list of cells in the workset, extract the subset of orientations that correspond
73 
74  const auto & local_orientations = *orientations_interface->getOrientations();
75  workset_orientations.resize(num_cells);
76 
77  // We can only apply orientations to owned and ghost cells - virtual cells are ignored (no orientations available)
78  auto local_cell_ids_host = Kokkos::create_mirror_view(local_cell_ids);
79  Kokkos::deep_copy(local_cell_ids_host, local_cell_ids);
80  for(int i=0; i<num_cells; ++i)
81  workset_orientations[i] = local_orientations[local_cell_ids_host[i]];
82 }
83 
84 void
85 applyBV2Orientations(const int num_cells,
86  BasisValues2<double> & basis_values,
87  const Kokkos::View<const panzer::LocalOrdinal*,PHX::Device> & local_cell_ids,
88  const Teuchos::RCP<const OrientationsInterface> & orientations_interface)
89 {
90  // This call exists because there is a middle man we have to go through.
91  // orientations_interface is a 'local' object, not a workset object so we need to map local cells to workset cells
92 
93  // If the object doesn't exist, feel free to skip applying orientations, they aren't needed in some cases (e.g. DG/FV)
94  if(orientations_interface.is_null())
95  return;
96 
97  // Ignore this operation if it has already been applied
98  if(basis_values.orientationsApplied())
99  return;
100 
101  // pull out the subset of orientations required for this workset
102  std::vector<Intrepid2::Orientation> workset_orientations(num_cells);
103  buildLocalOrientations(num_cells,local_cell_ids,orientations_interface, workset_orientations);
104 
105  basis_values.applyOrientations(workset_orientations,num_cells);
106 }
107 
108 }
109 
112  : num_cells(0)
113  , subcell_dim(-1)
114  , subcell_index(-1)
115  , ir_degrees(new std::vector<int>())
116  , basis_names(new std::vector<std::string>())
117  , setup_(false)
118  , num_owned_cells_(0)
119  , num_ghost_cells_(0)
120  , num_virtual_cells_(0)
121  , num_dimensions_(-1)
122 { }
123 
124 void
127  const WorksetOptions & options)
128 {
129 
130  num_cells = partition.local_cells.extent(0);
131  num_owned_cells_ = partition.num_owned_cells;
132  num_ghost_cells_ = partition.num_ghstd_cells;
134  options_ = options;
135 
137 
138  TEUCHOS_ASSERT(partition.cell_topology != Teuchos::null);
139  cell_topology_ = partition.cell_topology;
140 
141  num_dimensions_ = cell_topology_->getDimension();
142  subcell_dim = partition.subcell_dimension;
143  subcell_index = partition.subcell_index;
144  block_id = partition.element_block_name;
145  sideset_ = partition.sideset_name;
146 
147 
148  // Allocate and fill the local cell indexes for this workset
149  {
150  Kokkos::View<LocalOrdinal*, PHX::Device> cell_ids = Kokkos::View<LocalOrdinal*, PHX::Device>("cell_ids",num_cells);
151  Kokkos::deep_copy(cell_ids, partition.local_cells);
152  cell_local_ids_k = cell_ids;
153 
154  // DEPRECATED - only retain for backward compatability
155  auto local_cells_h = Kokkos::create_mirror_view(partition.local_cells);
156  Kokkos::deep_copy(local_cells_h, partition.local_cells);
157  cell_local_ids.resize(num_cells,-1);
158  for(int cell=0;cell<num_cells;++cell){
159  const int local_cell = local_cells_h(cell);
160  cell_local_ids[cell] = local_cell;
161  }
162  }
163 
164  // Allocate and fill the cell vertices
165  {
166  // Double check this
167  TEUCHOS_ASSERT(partition.cell_vertices.Rank == 3);
168 
169  // Grab the size of the cell vertices array
170  const int num_partition_cells = partition.cell_vertices.extent(0);
171  const int num_vertices_per_cell = partition.cell_vertices.extent(1);
172  const int num_dims_per_vertex = partition.cell_vertices.extent(2);
173 
174  // Make sure there isn't some strange problem going on
175  TEUCHOS_ASSERT(num_partition_cells == num_cells);
176  TEUCHOS_ASSERT(num_vertices_per_cell > 0);
177  TEUCHOS_ASSERT(num_dims_per_vertex > 0);
178 
179  // Allocate the worksets copy of the cell vertices
180  MDFieldArrayFactory af("",true);
181  cell_vertex_coordinates = af.buildStaticArray<double, Cell, NODE, Dim>("cell vertices", num_partition_cells, num_vertices_per_cell, num_dims_per_vertex);
182 
183  // Copy vertices over
184  const auto partition_vertices = partition.cell_vertices;
185  auto cvc = cell_vertex_coordinates.get_view();
186  Kokkos::parallel_for(num_cells, KOKKOS_LAMBDA (int i) {
187  for(int j=0;j<num_vertices_per_cell;++j)
188  for(int k=0;k<num_dims_per_vertex;++k)
189  cvc(i,j,k) = partition_vertices(i,j,k);
190  });
191  Kokkos::fence();
192  }
193 
194  // Add the subcell connectivity
195  if(partition.has_connectivity){
196  auto face_connectivity = Teuchos::rcp(new FaceConnectivity);
197  face_connectivity->setup(partition);
198  face_connectivity_ = face_connectivity;
199  }
200  // We have enough information to construct Basis/Point/Integration Values on the fly
201  setup_ = true;
202 
203 }
204 
205 bool
207 hasSubcellConnectivity(const unsigned int subcell_dimension) const
208 {
209  TEUCHOS_ASSERT(setup_);
210  return (subcell_dimension == (numDimensions() - 1)) and (not face_connectivity_.is_null());
211 }
212 
213 
214 const SubcellConnectivity &
216 getSubcellConnectivity(const unsigned int subcell_dimension) const
217 {
218  TEUCHOS_ASSERT(setup_);
219  TEUCHOS_TEST_FOR_EXCEPT_MSG(not hasSubcellConnectivity(subcell_dimension),
220  "Workset::getSubcellConnectivity : Requested subcell dimension "<<subcell_dimension<<" for a "<<num_dimensions_<<"D workset. This is not supported.");
221  // Right now we have only one option
222  return *face_connectivity_;
223 }
224 
227 {
228  TEUCHOS_ASSERT(face_connectivity_ != Teuchos::null);
229  return *face_connectivity_;
230 }
231 
235  const bool lazy_version) const
236 {
237  TEUCHOS_ASSERT(setup_);
238 
239  // We need unique keys for the lazy copy or else we get some weird behavior
240  size_t key = description.getKey();
241  if(lazy_version)
242  panzer::hash_combine<int>(key, 123);
243 
244  // Check if exists
245  const auto itr = integration_values_map_.find(key);
246  if(itr != integration_values_map_.end())
247  return *(itr->second);
248 
249  // Since it doesn't exist, we need to create it
250  const unsigned int subcell_dimension = numDimensions()-1;
251  int num_faces = -1;
252  if(hasSubcellConnectivity(subcell_dimension))
253  num_faces = getSubcellConnectivity(subcell_dimension).numSubcells();
254 
255  // For now, we need to make sure the descriptor lines up with the workset
257  TEUCHOS_TEST_FOR_EXCEPT_MSG(description.getSide() != getSubcellIndex(),
258  "Workset::getIntegrationValues : Attempted to build integration values for side '"<<description.getSide()<<"', but workset is constructed for side '"<<getSubcellIndex()<<"'");
259  }
260  auto ir = Teuchos::rcp(new IntegrationRule(description, cell_topology_, numCells(), num_faces));
261 
262  // Create the integration values object
263  Teuchos::RCP<IntegrationValues2<double>> iv;
264  if(lazy_version){
265  iv = Teuchos::rcp(new IntegrationValues2<double>());
266 
267  iv->setup(ir,getCellVertices(),numCells());
268 
269  // Surface integration schemes need to properly "permute" their entries to line up the surface points between cells
271  iv->setupPermutations(face_connectivity_, numVirtualCells());
272 
273  } else {
274 
275  iv = Teuchos::rcp(new IntegrationValues2<double>("",true));
276  iv->setupArrays(ir);
277  iv->evaluateValues(getCellVertices(), numCells(), face_connectivity_, numVirtualCells());
278 
279  // This is an advanced feature that requires changes to the workset construction
280  // Basically there needs to be a way to grab the side assembly for both "details" belonging to the same workset, which requires a refactor
281  // Alternatively, we can do this using a face connectivity object, but the refactor is more important atm.
282  TEUCHOS_ASSERT(not (options_.side_assembly_ and options_.align_side_points_));
283 
284  }
285 
286  integration_values_map_[key] = iv;
287  ir_degrees->push_back(iv->int_rule->cubature_degree);
288  int_rules.push_back(iv);
289 
290  return *iv;
291 
292 }
293 
297  const bool lazy_version) const
298 {
299  TEUCHOS_ASSERT(setup_);
300 
301  // Check if one exists (can be of any integration order)
302  const auto itr = basis_integration_values_map_.find(description.getKey());
303  if(itr != basis_integration_values_map_.end()){
304  for(const auto & pr : itr->second)
305  return *pr.second;
306  }
307 
308  // TODO: We currently overlap BasisIntegrationValues and BasisValues
309  // To counter this we create a fake integration rule at this point to ensure the basis values exist
310 
312 
313  // We have to have the right integrator if this is a side workset
315  TEUCHOS_ASSERT(getSubcellIndex() >= 0);
317  }
318 
319  // Now just use the other call
320  return getBasisValues(description, id, lazy_version);
321 
322 }
323 
324 
327 getBasisValues(const panzer::BasisDescriptor & basis_description,
328  const panzer::IntegrationDescriptor & integration_description,
329  const bool lazy_version) const
330 {
331  TEUCHOS_ASSERT(setup_);
332 
333  // We need unique keys for the lazy copy or else we get some weird behavior
334  size_t basis_key = basis_description.getKey();
335  if(lazy_version)
336  panzer::hash_combine<int>(basis_key, 123);
337 
338  // We need unique keys for the lazy copy or else we get some weird behavior
339  size_t integration_key = integration_description.getKey();
340  if(lazy_version)
341  panzer::hash_combine<int>(integration_key, 123);
342 
343  // Check if exists
344  const auto itr = basis_integration_values_map_.find(basis_key);
345  if(itr != basis_integration_values_map_.end()){
346  const auto & submap = itr->second;
347  const auto itr2 = submap.find(integration_key);
348  if(itr2 != submap.end())
349  return *(itr2->second);
350 
351  }
352 
353  // Get the integration values for this description
354  const auto & iv = getIntegrationValues(integration_description,lazy_version);
355  auto bir = Teuchos::rcp(new BasisIRLayout(basis_description.getType(), basis_description.getOrder(), *iv.int_rule));
356 
357  Teuchos::RCP<BasisValues2<double>> biv;
358 
359  if(lazy_version){
360 
361  // Initialized for lazy evaluation
362 
363  biv = Teuchos::rcp(new BasisValues2<double>());
364 
365  if(integration_description.getType() == IntegrationDescriptor::VOLUME)
366  biv->setupUniform(bir, iv.getUniformCubaturePointsRef(false), iv.getJacobian(false), iv.getJacobianDeterminant(false), iv.getJacobianInverse(false));
367  else
368  biv->setup(bir, iv.getCubaturePointsRef(false), iv.getJacobian(false), iv.getJacobianDeterminant(false), iv.getJacobianInverse(false));
369 
370  // pull out the subset of orientations required for this workset
371  std::vector<Intrepid2::Orientation> workset_orientations;
372  buildLocalOrientations(numOwnedCells()+numGhostCells(),getLocalCellIDs(),options_.orientations_, workset_orientations);
373 
374  biv->setOrientations(workset_orientations, numOwnedCells()+numGhostCells());
375  biv->setWeightedMeasure(iv.getWeightedMeasure(false));
376  biv->setCellVertexCoordinates(cell_vertex_coordinates);
377 
378  } else {
379 
380  // Standard, fully allocated version of BasisValues2
381 
382  biv = Teuchos::rcp(new BasisValues2<double>("", true, true));
383  biv->setupArrays(bir);
384  if((integration_description.getType() == IntegrationDescriptor::CV_BOUNDARY) or
385  (integration_description.getType() == IntegrationDescriptor::CV_SIDE) or
386  (integration_description.getType() == IntegrationDescriptor::CV_VOLUME)){
387 
388  biv->evaluateValuesCV(iv.ref_ip_coordinates,
389  iv.jac,
390  iv.jac_det,
391  iv.jac_inv,
392  getCellVertices(),
393  true,
394  numCells());
395  } else {
396 
397  if(integration_description.getType() == IntegrationDescriptor::VOLUME){
398 
399  // TODO: Eventually we will use the other call, however, that will be part of the BasisValues2 refactor
400  // The reason we don't do it now is that there are small differences (machine precision) that break EMPIRE testing
401  biv->evaluateValues(iv.cub_points,
402  iv.jac,
403  iv.jac_det,
404  iv.jac_inv,
405  iv.weighted_measure,
406  getCellVertices(),
407  true,
408  numCells());
409 
410  } else {
411 
412  biv->evaluateValues(iv.ref_ip_coordinates,
413  iv.jac,
414  iv.jac_det,
415  iv.jac_inv,
416  iv.weighted_measure,
417  getCellVertices(),
418  true,
419  numCells());
420  }
421  }
422 
423  applyBV2Orientations(numOwnedCells()+numGhostCells(),*biv,getLocalCellIDs(),options_.orientations_);
424 
425  }
426 
427  basis_integration_values_map_[basis_key][integration_key] = biv;
428  bases.push_back(biv);
429  basis_names->push_back(biv->basis_layout->name());
430 
431  return *biv;
432 
433 }
434 
435 
438 getPointValues(const panzer::PointDescriptor & description) const
439 {
440  TEUCHOS_ASSERT(setup_);
441 
442 
443  // Check if exists
444  const auto itr = point_values_map_.find(description.getKey());
445  if(itr != point_values_map_.end())
446  return *(itr->second);
447 
448  TEUCHOS_TEST_FOR_EXCEPT_MSG(not description.hasGenerator(),
449  "Point Descriptor of type '"<<description.getType()<<"' does not have associated generator.");
450 
451  auto pr = Teuchos::rcp(new PointRule(description, cell_topology_, numCells()));
452 
453  auto pv = Teuchos::rcp(new PointValues2<double>("",true));
454 
455  pv->setupArrays(pr);
456 
457  // Point values are not necessarily set at the workset level, but can be set by evaluators
458  if(description.hasGenerator())
459  if(description.getGenerator().hasPoints(*cell_topology_))
460  pv->evaluateValues(getCellVertices(), description.getGenerator().getPoints(*cell_topology_),false, numCells());
461 
462  point_values_map_[description.getKey()] = pv;
463 
464  return *pv;
465 
466 }
467 
470 getBasisValues(const panzer::BasisDescriptor & basis_description,
471  const panzer::PointDescriptor & point_description,
472  const bool lazy_version) const
473 {
474  TEUCHOS_ASSERT(setup_);
475 
476  // We need unique keys for the lazy copy or else we get some weird behavior
477  size_t basis_key = basis_description.getKey();
478  if(lazy_version)
479  panzer::hash_combine<int>(basis_key, 123);
480 
481  // Check if exists
482  const auto itr = basis_point_values_map_.find(basis_key);
483  if(itr != basis_point_values_map_.end()){
484  const auto & submap = itr->second;
485  const auto itr2 = submap.find(point_description.getKey());
486  if(itr2 != submap.end())
487  return *(itr2->second);
488 
489  }
490 
491  // Get the integration values for this description
492  const auto & pv = getPointValues(point_description);
493 
494  auto bir = Teuchos::rcp(new BasisIRLayout(basis_description.getType(), basis_description.getOrder(), *pv.point_rule));
495 
496  Teuchos::RCP<BasisValues2<double>> bpv;
497 
498  if(lazy_version){
499 
500  // Initialized for lazy evaluation
501 
502  bpv = Teuchos::rcp(new BasisValues2<double>());
503 
504  bpv->setupUniform(bir, pv.coords_ref, pv.jac, pv.jac_det, pv.jac_inv);
505 
506  // pull out the subset of orientations required for this workset
507  std::vector<Intrepid2::Orientation> workset_orientations;
508  buildLocalOrientations(numOwnedCells()+numGhostCells(),getLocalCellIDs(),options_.orientations_, workset_orientations);
509 
510  bpv->setOrientations(workset_orientations, numOwnedCells()+numGhostCells());
511  bpv->setCellVertexCoordinates(cell_vertex_coordinates);
512 
513  } else {
514 
515  // Standard fully allocated version
516 
517  bpv = Teuchos::rcp(new BasisValues2<double>("", true, false));
518  bpv->setupArrays(bir);
519  bpv->evaluateValues(pv.coords_ref,
520  pv.jac,
521  pv.jac_det,
522  pv.jac_inv,
523  numCells());
524 
525  // TODO: We call this separately due to how BasisValues2 is structured - needs to be streamlined
526  bpv->evaluateBasisCoordinates(getCellVertices(),numCells());
527 
528  applyBV2Orientations(numOwnedCells()+numGhostCells(),*bpv, getLocalCellIDs(), options_.orientations_);
529 
530  }
531 
532  basis_point_values_map_[basis_key][point_description.getKey()] = bpv;
533 
534  return *bpv;
535 
536 }
537 
541 {
542  const auto itr = _integration_rule_map.find(description.getKey());
543  if(itr == _integration_rule_map.end()){
544 
545  // Must run setup or cell topology wont be set properly
546  TEUCHOS_ASSERT(setup_);
547 
548  // Since it doesn't exist, we need to create it
549  const unsigned int subcell_dimension = numDimensions()-1;
550  int num_faces = -1;
551  if(hasSubcellConnectivity(subcell_dimension))
552  num_faces = getSubcellConnectivity(subcell_dimension).numSubcells();
553 
554  // For now, we need to make sure the descriptor lines up with the workset
556  TEUCHOS_TEST_FOR_EXCEPT_MSG(description.getSide() != getSubcellIndex(),
557  "Workset::getIntegrationValues : Attempted to build integration values for side '"<<description.getSide()<<"', but workset is constructed for side '"<<getSubcellIndex()<<"'");
558  }
559 
560  auto ir = Teuchos::rcp(new IntegrationRule(description, cell_topology_, numCells(), num_faces));
561 
562  _integration_rule_map[description.getKey()] = ir;
563 
564  return *ir;
565  }
566  return *(itr->second);
567 }
568 
569 const panzer::PureBasis &
571 getBasis(const panzer::BasisDescriptor & description) const
572 {
573  const auto itr = _pure_basis_map.find(description.getKey());
574  if(itr == _pure_basis_map.end()){
575 
576  // Must run setup or cell topology wont be set properly
577  TEUCHOS_ASSERT(setup_);
578 
579  // Create and storethe pure basis
580  Teuchos::RCP<panzer::PureBasis> basis = Teuchos::rcp(new panzer::PureBasis(description, cell_topology_, numCells()));
581  _pure_basis_map[description.getKey()] = basis;
582  return *basis;
583  }
584  return *(itr->second);
585 }
586 
587 
588 void
590 setNumberOfCells(const int o_cells,
591  const int g_cells,
592  const int v_cells)
593 {
594  num_owned_cells_ = o_cells;
595  num_ghost_cells_ = g_cells;
596  num_virtual_cells_ = v_cells;
597  num_cells = o_cells + g_cells + v_cells;
598 }
599 
600 std::ostream&
601 operator<<(std::ostream& os,
602  const panzer::Workset& w)
603 {
604  using std::endl;
605 
606  os << "Workset" << endl;
607  os << " block_id=" << w.getElementBlock() << endl;
608  os << " num_cells:" << w.num_cells << endl;
609  os << " num_owned_cells:" << w.numOwnedCells() << endl;
610  os << " num_ghost_cells:" << w.numGhostCells() << endl;
611  os << " num_virtual_cells:" << w.numVirtualCells() << endl;
612  os << " cell_local_ids (size=" << w.getLocalCellIDs().size() << ")" << endl;
613  os << " subcell_dim = " << w.getSubcellDimension() << endl;
614  os << " subcell_index = " << w.getSubcellIndex() << endl;
615 
616  os << " ir_degrees: " << endl;
617  for (std::vector<int>::const_iterator ir = w.ir_degrees->begin();
618  ir != w.ir_degrees->end(); ++ir)
619  os << " " << *ir << std::endl;
620 
621  std::vector<int>::const_iterator ir = w.ir_degrees->begin();
622  for (std::vector<Teuchos::RCP<panzer::IntegrationValues2<double> > >::const_iterator irv = w.int_rules.begin();
623  irv != w.int_rules.end(); ++irv,++ir) {
624 
625  os << " IR Values (Degree=" << *ir << "):" << endl;
626 
627  os << " cub_points:" << endl;
628  os << (*irv)->cub_points << endl;
629 
630  os << " side_cub_points:" << endl;
631  os << (*irv)->side_cub_points << endl;
632 
633  os << " cub_weights:" << endl;
634  os << (*irv)->cub_weights << endl;
635 
636  os << " node_coordinates:" << endl;
637  os << (*irv)->node_coordinates << endl;
638 
639  os << " jac:" << endl;
640  os << (*irv)->jac << endl;
641 
642  os << " jac_inv:" << endl;
643  os << (*irv)->jac_inv << endl;
644 
645  os << " jac_det:" << endl;
646  os << (*irv)->jac_det << endl;
647 
648  os << " weighted_measure:" << endl;
649  os << (*irv)->weighted_measure << endl;
650 
651  os << " covarient:" << endl;
652  os << (*irv)->covarient << endl;
653 
654  os << " contravarient:" << endl;
655  os << (*irv)->contravarient << endl;
656 
657  os << " norm_contravarient:" << endl;
658  os << (*irv)->norm_contravarient << endl;
659 
660  os << " ip_coordinates:" << endl;
661  os << (*irv)->ip_coordinates << endl;
662 
663  os << " int_rule->getName():" << (*irv)->int_rule->getName() << endl;
664  }
665 
666 
667  os << " basis_names: " << endl;
668  for (std::vector<std::string>::const_iterator b = w.basis_names->begin();
669  b != w.basis_names->end(); ++b)
670  os << " " << *b << std::endl;
671 
672  std::vector<std::string>::const_iterator b = w.basis_names->begin();
673 
674  for (std::vector<Teuchos::RCP< panzer::BasisValues2<double> > >::const_iterator bv = w.bases.begin(); bv != w.bases.end(); ++bv,++b) {
675 
676  os << " Basis Values (basis_name=" << *b << "):" << endl;
677 
678 /*
679  os << " basis_ref:" << endl;
680  os << (*bv)->basis_ref << endl;
681 
682  os << " basis:" << endl;
683  os << (*bv)->basis_scalar << endl;
684 
685  os << " grad_basis_ref:" << endl;
686  os << (*bv)->grad_basis_ref << endl;
687 
688  os << " grad_basis:" << endl;
689  os << (*bv)->grad_basis << endl;
690 
691  os << " curl_basis_ref:" << endl;
692  os << (*bv)->curl_basis_ref_vector << endl;
693 
694  os << " curl_basis:" << endl;
695  os << (*bv)->curl_basis_vector << endl;
696 
697  os << " basis_coordinates_ref:" << endl;
698  os << (*bv)->basis_coordinates_ref << endl;
699 
700  os << " basis_coordinates:" << endl;
701  os << (*bv)->basis_coordinates << endl;
702 */
703 
704  os << " basis_layout->name():" << (*bv)->basis_layout->name() << endl;
705  }
706 
707 
708 
709  return os;
710 }
711 
712 }
std::size_t getKey() const
Get unique key associated with basis of this order and type The key is used to sort through a map of ...
const panzer::IntegrationValues2< double > & getIntegrationValues(const panzer::IntegrationDescriptor &description, const bool lazy_version=false) const
Get the integration values for a given integration description.
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
int num_cells
DEPRECATED - use: numCells()
Teuchos::RCP< const shards::CellTopology > cell_topology
const panzer::IntegrationRule & getIntegrationRule(const panzer::IntegrationDescriptor &description) const
Grab the integration rule for a given integration description (throws error if integration doesn&#39;t ex...
const SubcellConnectivity & getSubcellConnectivity(const unsigned int subcell_dimension) const
Get the subcell connectivity for the workset topology.
std::map< size_t, Teuchos::RCP< const panzer::IntegrationRule > > _integration_rule_map
const int & getSide() const
Get side associated with integration - this is for backward compatibility.
PHX::View< const int * > cell_local_ids_k
bool hasGenerator() const
Check if the point descriptor has a generator for generating point values.
virtual Kokkos::DynRankView< double > getPoints(const shards::CellTopology &topo) const =0
Get the points for a particular topology.
Teuchos::RCP< panzer::SubcellConnectivity > face_connectivity_
KOKKOS_INLINE_FUNCTION int numSubcells() const
Gives number of subcells (e.g. faces) in connectivity.
int numVirtualCells() const
Number of cells not owned by any workset - these are used for boundary conditions.
std::map< size_t, std::map< size_t, Teuchos::RCP< panzer::BasisValues2< double > > > > basis_point_values_map_
const int & getType() const
Get type of integrator.
std::vector< size_t > cell_local_ids
Integral over a specific side of cells (side must be set)
Used to define options for lazy evaluation of BasisValues and IntegrationValues objects.
PHX::View< panzer::LocalOrdinal * > local_cells
Teuchos::RCP< std::vector< int > > ir_degrees
If workset corresponds to a sub cell, what is the index?
std::vector< Teuchos::RCP< panzer::BasisValues2< double > > > bases
Static basis function data, key is basis name, value is index in the static_bases vector...
std::map< size_t, Teuchos::RCP< const panzer::IntegrationValues2< double > > > integration_values_map_
Teuchos::RCP< std::vector< std::string > > basis_names
Value corresponds to basis type. Use the offest for indexing.
PHX::View< double *** > cell_vertices
int subcell_dim
DEPRECATED - use: getSubcellDimension()
unsigned int numDimensions() const
Get the cell dimension for the mesh.
std::map< size_t, Teuchos::RCP< const panzer::PureBasis > > _pure_basis_map
const std::string & getElementBlock() const
Get the element block id.
int numCells() const
Number of total cells in workset (owned, ghost, and virtual)
No integral specified - default state.
Generates a SubcellConnectivity associated with faces and cells given a partition of the local mesh...
const std::string & getType() const
Get unique string associated with the type of point descriptor. This will be used generate a hash to ...
panzer::LocalOrdinal num_owned_cells
CellCoordArray cell_vertex_coordinates
DEPRECATED - use: getCellVertices()
panzer::PointValues2< double > & getPointValues(const panzer::PointDescriptor &point_description) const
Grab the basis values for a given basis description and integration description (throws error if it d...
void setNumberOfCells(const int owned_cells, const int ghost_cells, const int virtual_cells)
Provides access to set numbers of cells (required for backwards compatibility)
panzer::LocalOrdinal num_virtual_cells
Teuchos::RCP< const shards::CellTopology > cell_topology_
Teuchos::RCP< const OrientationsInterface > orientations_
Must be set to apply orientations - if it is set, then orientations will be applied to basis values...
int getSubcellIndex() const
Get the subcell index (returns -1 if not a subcell)
bool side_assembly_
Build integration values for sides.
std::vector< Teuchos::RCP< panzer::IntegrationValues2< double > > > int_rules
std::map< size_t, Teuchos::RCP< panzer::PointValues2< double > > > point_values_map_
const PointGenerator & getGenerator() const
const panzer::BasisValues2< double > & getBasisValues(const panzer::BasisDescriptor &basis_description, const bool lazy_version=false) const
bool hasSubcellConnectivity(const unsigned int subcell_dimension) const
Check if subcell connectivity exists for a given dimension.
int getSubcellDimension() const
Get the subcell dimension.
panzer::LocalOrdinal num_ghstd_cells
Integral over all sides of cells (closed surface integral)
std::map< size_t, std::map< size_t, Teuchos::RCP< panzer::BasisValues2< double > > > > basis_integration_values_map_
Kokkos::View< const int *, PHX::Device > getLocalCellIDs() const
Get the local cell IDs for the workset.
std::string block_id
DEPRECATED - use: getElementBlock()
WorksetDetails()
Default constructor.
int subcell_index
DEPRECATED - use: getSubcellIndex()
int numOwnedCells() const
Number of cells owned by this workset.
std::size_t getKey() const
Get unique key associated with integrator of this order and type The key is used to sort through a ma...
std::ostream & operator<<(std::ostream &os, const AssemblyEngineInArgs &in)
bool align_side_points_
If workset side integration values must align with another workset, there must be a unique order assi...
std::size_t getKey() const
Get unique key associated with integrator of this order and type The key is used to sort through a ma...
const panzer::SubcellConnectivity & getFaceConnectivity() const
int numGhostCells() const
Number of cells owned by a different workset.
virtual bool hasPoints(const shards::CellTopology &topo) const =0
Check if the generator can generate points for the given topology.
Description and data layouts associated with a particular basis.
CellCoordArray getCellVertices() const
Get the vertices for the cells.
void setup(const LocalMeshPartition &partition, const WorksetOptions &options)
Constructs the workset details from a given chunk of the mesh.
int getOrder() const
Get order of basis.
const std::string & getType() const
Get type of basis.
const panzer::PureBasis & getBasis(const panzer::BasisDescriptor &description) const
Grab the pure basis (contains data layouts) for a given basis description (throws error if integratio...