41 #ifndef DTK_INTREPIDCELL_HPP 42 #define DTK_INTREPIDCELL_HPP 44 #include <Teuchos_RCP.hpp> 46 #include <Shards_CellTopology.hpp> 48 #include <Intrepid_Cubature.hpp> 49 #include <Intrepid_FieldContainer.hpp> 64 typedef Intrepid::FieldContainer<Scalar> MDArray;
73 const unsigned degree );
93 const MDArray &cell_node_coords );
98 MDArray &reference_coords );
103 MDArray &physical_coords );
108 const double tolerance );
129 getPhysicalIntegrationCoordinates( MDArray &physical_ip_coordinates )
const;
135 void integrate(
const MDArray &dofs, MDArray &integrals )
const;
139 shards::CellTopology d_topology;
142 Teuchos::RCP<Intrepid::Cubature<Scalar, MDArray>> d_cubature;
145 MDArray d_cub_points;
148 MDArray d_cub_weights;
154 MDArray d_jacobian_det;
157 MDArray d_weighted_measures;
160 MDArray d_cell_node_coords;
163 MDArray d_physical_ip_coordinates;
172 #endif // end DTK_INTREPIDCELL_HPP int getNumIntegrationPoints() const
Get the number of cell points.
bool pointInPhysicalCell(const MDArray &point, const double tolerance)
Determine if a point given in physical coordinates is inside of the phyiscal cell.
static void updateState(IntrepidCell &Intrepid_cell, const MDArray &cell_node_coords)
Free function for updating the cell state for a new set of physical cells in a single call...
virtual ~IntrepidCell()
Destructor.
int getSpatialDimension() const
Get the spatial dimension.
Manager for Intrepid cell-level operations.
void getCellMeasures(MDArray &cell_measures) const
Get the cell measures (Cell). cell_measures must all ready be allocated.
bool pointInReferenceCell(const MDArray &reference_point, const double tolerance)
Determine if a point given in parametric coordinates is inside of the reference cell.
virtual void mapToCellPhysicalFrame(const MDArray ¶metric_coords, MDArray &physical_coords)
Given a set of coordinates in the reference frame of the cell, map them to the physical frame...
void allocateCellState(const MDArray &cell_node_coords)
Given physical coordinates for the cell nodes (Cell,Node,Dim), allocate the state of the cell object...
void integrate(const MDArray &dofs, MDArray &integrals) const
Given DOFs at the quadrature points {(Cell,Node) for scalar fields, (Cell,Node,VecDim) for vector fie...
virtual void updateCellState()
Update the state of the cell object for the current cell node coordinates.
void mapToCellReferenceFrame(const MDArray &physical_coords, MDArray &reference_coords)
Given a set of coordinates in the physical frame of the cell, map them to the reference frame of the ...
void setCellNodeCoordinates(const MDArray &cell_node_coords)
Given physical coordinates for the cell nodes (Cell,Node,Dim), assign them to the cell without alloca...
int getNumCells() const
Get the number of cells in the current state.