Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_CellData.hpp
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#ifndef PANZER_CELL_DATA_HPP
12#define PANZER_CELL_DATA_HPP
13
14#include "PanzerDiscFE_config.hpp"
15
16#include "Teuchos_Assert.hpp"
17#include "Teuchos_RCP.hpp"
18
19#include "Shards_CellTopology.hpp"
20#include "Shards_BasicTopologies.hpp"
21
22namespace panzer {
23
33 class CellData {
34
35 public:
36
38 explicit CellData()
39 { }
40
44 explicit CellData(std::size_t num_cells,
45 const Teuchos::RCP<const shards::CellTopology> & ct) :
46 m_num_cells(num_cells),
47 m_is_side(false),
48 m_side(-1),
49 m_cell_topo(ct)
50 { }
51
55 explicit CellData(std::size_t num_cells,
56 int local_side_id,const Teuchos::RCP<const shards::CellTopology> & ct) :
57 m_num_cells(num_cells),
58 m_is_side(local_side_id >= 0),
59 m_side(local_side_id),
60 m_cell_topo(ct)
61 { }
62
63 bool isSide() const
64 {return m_is_side;}
65
66 int side() const
67 {
68 TEUCHOS_TEST_FOR_EXCEPTION(!m_is_side, std::logic_error,
69 "Cannot return side index, CellData is not a side!");
70 return m_side;
71 }
72
73 std::size_t numCells() const
74 {return m_num_cells;}
75
77 int baseCellDimension() const
78 {return m_cell_topo->getDimension();}
79
81 // TODO BWR this terminology is a bit confusing... this is not the base cell topology
82 Teuchos::RCP<const shards::CellTopology> getCellTopology() const
83 { return m_cell_topo; }
84
85 private:
86 std::size_t m_num_cells;
88 int m_side;
89
90
91 Teuchos::RCP<const shards::CellTopology> m_cell_topo;
92 };
93
94}
95
96#endif
Data for determining cell topology and dimensionality.
Teuchos::RCP< const shards::CellTopology > m_cell_topo
int baseCellDimension() const
Dimension of the base cell. NOT the dimension of the local side, even if the side() method returns tr...
CellData(std::size_t num_cells, const Teuchos::RCP< const shards::CellTopology > &ct)
Teuchos::RCP< const shards::CellTopology > getCellTopology() const
Get CellTopology for the base cell.
std::size_t m_num_cells
std::size_t numCells() const
CellData(std::size_t num_cells, int local_side_id, const Teuchos::RCP< const shards::CellTopology > &ct)