45 #include <Teuchos_as.hpp> 47 #include <Intrepid_CellTools.hpp> 48 #include <Intrepid_DefaultCubatureFactory.hpp> 49 #include <Intrepid_FunctionSpaceTools.hpp> 50 #include <Intrepid_Types.hpp> 58 IntrepidCell::IntrepidCell(
const shards::CellTopology &cell_topology,
59 const unsigned degree )
60 : d_topology( cell_topology )
61 , d_cub_points( 0, 0 )
63 , d_jacobian( 0, 0, 0, 0 )
64 , d_jacobian_det( 0, 0 )
65 , d_weighted_measures( 0, 0 )
66 , d_physical_ip_coordinates( 0, 0, 0 )
68 Intrepid::DefaultCubatureFactory<Scalar, MDArray> cub_factory;
69 d_cubature = cub_factory.create( d_topology, degree );
71 unsigned num_cub_points = d_cubature->getNumPoints();
72 unsigned cub_dim = cell_topology.getDimension();
74 d_cub_points.resize( num_cub_points, cub_dim );
75 d_cub_weights.resize( num_cub_points );
76 d_cubature->getCubature( d_cub_points, d_cub_weights );
92 d_cell_node_coords = cell_node_coords;
106 int num_cells = d_cell_node_coords.dimension( 0 );
107 int num_ip = d_cub_points.dimension( 0 );
108 int space_dim = d_cub_points.dimension( 1 );
111 d_jacobian.resize( num_cells, num_ip, space_dim, space_dim );
112 d_jacobian_det.resize( num_cells, num_ip );
113 d_weighted_measures.resize( num_cells, num_ip );
114 d_physical_ip_coordinates.resize( num_cells, num_ip, space_dim );
125 Intrepid::CellTools<Scalar>::setJacobian( d_jacobian, d_cub_points,
126 d_cell_node_coords, d_topology );
129 Intrepid::CellTools<Scalar>::setJacobianDet( d_jacobian_det, d_jacobian );
132 Intrepid::FunctionSpaceTools::computeCellMeasure<Scalar>(
133 d_weighted_measures, d_jacobian_det, d_cub_weights );
136 Intrepid::CellTools<Scalar>::mapToPhysicalFrame(
137 d_physical_ip_coordinates, d_cub_points, d_cell_node_coords,
147 const MDArray &cell_node_coords )
159 MDArray &reference_coords )
161 DTK_REQUIRE( 2 == physical_coords.rank() );
162 DTK_REQUIRE( 2 == reference_coords.rank() );
163 Intrepid::CellTools<Scalar>::mapToReferenceFrame(
164 reference_coords, physical_coords, d_cell_node_coords, d_topology, 0 );
173 MDArray &physical_coords )
175 DTK_REQUIRE( 2 == parametric_coords.rank() );
176 DTK_REQUIRE( 3 == physical_coords.rank() );
177 DTK_REQUIRE( parametric_coords.dimension( 1 ) ==
178 Teuchos::as<int>( d_topology.getDimension() ) );
179 DTK_REQUIRE( physical_coords.dimension( 0 ) ==
180 d_cell_node_coords.dimension( 0 ) );
181 DTK_REQUIRE( physical_coords.dimension( 1 ) ==
182 parametric_coords.dimension( 0 ) );
183 DTK_REQUIRE( physical_coords.dimension( 2 ) ==
184 Teuchos::as<int>( d_topology.getDimension() ) );
186 Intrepid::CellTools<Scalar>::mapToPhysicalFrame(
187 physical_coords, parametric_coords, d_cell_node_coords, d_topology );
196 const double tolerance )
198 return Intrepid::CellTools<Scalar>::checkPointsetInclusion(
199 reference_point, d_topology, tolerance );
208 const double tolerance )
210 MDArray reference_point( point );
211 reference_point.initialize();
222 return d_weighted_measures.dimension( 0 );
231 return d_cub_points.dimension( 0 );
240 return d_cub_points.dimension( 1 );
250 DTK_REQUIRE( 1 == cell_measures.rank() );
251 DTK_REQUIRE( cell_measures.dimension( 0 ) ==
252 d_weighted_measures.dimension( 0 ) );
254 MDArray dofs( d_cell_node_coords.dimension( 0 ),
255 d_cub_weights.dimension( 0 ) );
256 dofs.initialize( 1.0 );
263 void IntrepidCell::getPhysicalIntegrationCoordinates(
264 MDArray &physical_ip_coordinates )
const 266 physical_ip_coordinates = d_physical_ip_coordinates;
278 Intrepid::FunctionSpaceTools::integrate<Scalar>(
279 integrals, dofs, d_weighted_measures, Intrepid::COMP_BLAS );
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.
Assertions and Design-by-Contract for error handling.
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.
Manager for Intrepid cell-level operations.
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.