Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_Workset_Builder.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #include "PanzerDiscFE_config.hpp"
12 
14 
17 
18 template
20 panzer::buildWorksets(const WorksetNeeds & needs,
21  const std::string & elementBlock,
22  const std::vector<std::size_t>& local_cell_ids,
23  const Kokkos::DynRankView<double,PHX::Device>& node_coordinates);
24 
25 template
27 panzer::buildBCWorkset(const WorksetNeeds& needs,
28  const std::string& elementBlock,
29  const std::vector<std::size_t>& local_cell_ids,
30  const std::vector<std::size_t>& local_side_ids,
31  const Kokkos::DynRankView<double,PHX::Device>& node_coordinates,
32  const bool populate_value_arrays);
33 
34 template
36 panzer::buildBCWorkset(const WorksetNeeds & needs_a,
37  const std::string & blockid_a,
38  const std::vector<std::size_t>& local_cell_ids_a,
39  const std::vector<std::size_t>& local_side_ids_a,
40  const Kokkos::DynRankView<double,PHX::Device> & node_coordinates_a,
41  const panzer::WorksetNeeds & needs_b,
42  const std::string & blockid_b,
43  const std::vector<std::size_t>& local_cell_ids_b,
44  const std::vector<std::size_t>& local_side_ids_b,
45  const Kokkos::DynRankView<double,PHX::Device> & node_coordinates_b);
46 
47 template
50  const std::string & eblock_a,
51  const std::vector<std::size_t>& local_cell_ids_a,
52  const std::vector<std::size_t>& local_side_ids_a,
53  const Kokkos::DynRankView<double,PHX::Device> & node_coordinates_a,
54  const panzer::WorksetNeeds & needs_b,
55  const std::string & eblock_b,
56  const std::vector<std::size_t>& local_cell_ids_b,
57  const std::vector<std::size_t>& local_side_ids_b,
58  const Kokkos::DynRankView<double,PHX::Device> & node_coordinates_b);
59 
60 namespace panzer {
61 
62 void populateValueArrays(std::size_t num_cells,bool isSide,const WorksetNeeds & needs,
63  WorksetDetails & details,const Teuchos::RCP<WorksetDetails> other_details)
64 {
65  using Teuchos::RCP;
66  using Teuchos::rcp;
67 
68  // check to see if the number of owned, ghosted and virtual cells are needed
69  // this is required for backwards compatibility
71  if(details.numOwnedCells()==-1)
72  details.setNumberOfCells(num_cells,0,0);
73  if(other_details!=Teuchos::null and other_details->numOwnedCells()==-1)
74  other_details->setNumberOfCells(num_cells,0,0);
75 
76  panzer::MDFieldArrayFactory mdArrayFactory("",true);
77 
78  // setup the integration rules and bases
79 
80  std::vector<RCP<const panzer::PureBasis> > bases;
81  std::vector<RCP<const panzer::IntegrationRule> > int_rules;
82  if(isSide) {
83  const panzer::CellData side_cell_data(num_cells,
84  details.subcell_index,
85  needs.cellData.getCellTopology());
86 
87  for(std::size_t i=0;i<needs.int_rules.size();i++)
88  int_rules.push_back(rcp(new IntegrationRule(needs.int_rules[i]->cubature_degree,side_cell_data)));
89 
90  for(std::size_t i=0;i<needs.bases.size();i++)
91  bases.push_back(rcp(new PureBasis(needs.bases[i]->type(),needs.bases[i]->order(),side_cell_data)));
92  }
93  else {
94  int_rules = needs.int_rules;
95  bases = needs.bases;
96  }
97 
98  details.ir_degrees = rcp(new std::vector<int>(0));
99  details.basis_names = rcp(new std::vector<std::string>(0));
100 
101  for(std::size_t i=0;i<int_rules.size();i++) {
102 
103  details.ir_degrees->push_back(int_rules[i]->cubature_degree);
104 
107  iv2->setupArrays(int_rules[i]);
108  if (Teuchos::nonnull(other_details))
109  iv2->evaluateValues(details.cell_node_coordinates, other_details->int_rules[i]->ip_coordinates,num_cells);
110  else
111  iv2->evaluateValues(details.cell_node_coordinates,num_cells);
112 
113  details.int_rules.push_back(iv2);
114 
115  // Need to create all combinations of basis/ir pairings
116  for(std::size_t b=0;b<bases.size();b++) {
117  RCP<panzer::BasisIRLayout> b_layout = rcp(new panzer::BasisIRLayout(bases[b],*int_rules[i]));
118  details.basis_names->push_back(b_layout->name());
119 
120  std::size_t int_degree_index = std::distance(details.ir_degrees->begin(),
121  std::find(details.ir_degrees->begin(),
122  details.ir_degrees->end(),
123  int_rules[i]->order()));
125  rcp(new panzer::BasisValues2<double>("",true,true));
126  bv2->setupArrays(b_layout);
127  bv2->evaluateValues(details.int_rules[int_degree_index]->cub_points,
128  details.int_rules[int_degree_index]->jac,
129  details.int_rules[int_degree_index]->jac_det,
130  details.int_rules[int_degree_index]->jac_inv,
131  details.int_rules[int_degree_index]->weighted_measure,
132  details.cell_node_coordinates,
133  true,
134  num_cells);
135 
136  details.bases.push_back(bv2);
137  }
138  }
139 
140 }
141 
142 }
std::vector< Teuchos::RCP< const PureBasis > > bases
bool nonnull(const std::shared_ptr< T > &p)
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...
Teuchos::RCP< std::vector< std::string > > basis_names
Value corresponds to basis type. Use the offest for indexing.
Teuchos::RCP< std::vector< Workset > > buildWorksets(const WorksetNeeds &needs, const std::string &elementBlock, const std::vector< std::size_t > &local_cell_ids, const ArrayT &node_coordinates)
int numOwnedCells() const
Number of cells owned by this workset.
Teuchos::RCP< const shards::CellTopology > getCellTopology() const
Get CellTopology for the base cell.
std::vector< Teuchos::RCP< const IntegrationRule > > int_rules
Teuchos::RCP< std::vector< Workset > > buildEdgeWorksets(const WorksetNeeds &needs_a, const std::string &eblock_a, const std::vector< std::size_t > &local_cell_ids_a, const std::vector< std::size_t > &local_side_ids_a, const ArrayT &node_coordinates_a, const WorksetNeeds &needs_b, const std::string &eblock_b, const std::vector< std::size_t > &local_cell_ids_b, const std::vector< std::size_t > &local_side_ids_b, const ArrayT &node_coordinates_b)
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)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Data for determining cell topology and dimensionality.
std::vector< Teuchos::RCP< panzer::IntegrationValues2< double > > > int_rules
Teuchos::RCP< std::map< unsigned, Workset > > buildBCWorkset(const WorksetNeeds &needs, const std::string &elementBlock, const std::vector< std::size_t > &local_cell_ids, const std::vector< std::size_t > &local_side_ids, const ArrayT &node_coordinates, const bool populate_value_arrays=true)
std::string name() const
Unique key for workset indexing composed of basis name and point rule name.
int subcell_index
DEPRECATED - use: getSubcellIndex()
CellCoordArray cell_node_coordinates
DEPRECATED - use: getCellNodes()
Description and data layouts associated with a particular basis.
void populateValueArrays(std::size_t num_cells, bool isSide, const WorksetNeeds &needs, WorksetDetails &details, const Teuchos::RCP< WorksetDetails > other_details)