44 #include <Teuchos_as.hpp> 46 #include <Intrepid_CellTools.hpp> 47 #include <Intrepid_DefaultCubatureFactory.hpp> 48 #include <Intrepid_FunctionSpaceTools.hpp> 57 const unsigned side_id,
58 const shards::CellTopology &parent_topology,
59 const unsigned degree )
60 :
Base( side_topology, degree )
61 , d_side_id( side_id )
62 , d_parent_topology( parent_topology )
65 MDArray mapped_cub_points( this->d_cubature->getNumPoints(),
66 parent_topology.getDimension() );
67 Intrepid::CellTools<Scalar>::mapToReferenceSubcell(
68 mapped_cub_points, this->d_cub_points, side_topology.getDimension(),
69 d_side_id, d_parent_topology );
70 this->d_cub_points = mapped_cub_points;
80 unsigned space_dim = d_parent_topology.getDimension();
83 Intrepid::CellTools<Scalar>::setJacobian(
84 this->d_jacobian, this->d_cub_points, this->d_cell_node_coords,
93 Intrepid::FunctionSpaceTools::computeFaceMeasure<Scalar>(
94 this->d_weighted_measures, this->d_jacobian, this->d_cub_weights,
95 d_side_id, d_parent_topology );
102 Intrepid::FunctionSpaceTools::computeEdgeMeasure<Scalar>(
103 this->d_weighted_measures, this->d_jacobian, this->d_cub_weights,
104 d_side_id, d_parent_topology );
109 DTK_INSIST( 3 == space_dim || 2 == space_dim );
114 Intrepid::CellTools<Scalar>::mapToPhysicalFrame(
115 this->d_physical_ip_coordinates, this->d_cub_points,
116 this->d_cell_node_coords, d_parent_topology );
125 MDArray &physical_coords )
127 DTK_REQUIRE( 2 == parametric_coords.rank() );
128 DTK_REQUIRE( 3 == physical_coords.rank() );
129 DTK_REQUIRE( parametric_coords.dimension( 1 ) ==
130 Teuchos::as<int>( this->d_topology.getDimension() ) );
131 DTK_REQUIRE( physical_coords.dimension( 0 ) ==
132 this->d_cell_node_coords.dimension( 0 ) );
133 DTK_REQUIRE( physical_coords.dimension( 1 ) ==
134 parametric_coords.dimension( 0 ) );
135 DTK_REQUIRE( physical_coords.dimension( 2 ) ==
136 Teuchos::as<int>( d_parent_topology.getDimension() ) );
138 MDArray mapped_coords( parametric_coords.dimension( 0 ),
139 d_parent_topology.getDimension() );
141 Intrepid::CellTools<Scalar>::mapToReferenceSubcell(
142 mapped_coords, parametric_coords, this->d_topology.getDimension(),
143 d_side_id, d_parent_topology );
145 Intrepid::CellTools<Scalar>::mapToPhysicalFrame(
146 physical_coords, mapped_coords, this->d_cell_node_coords,
155 MDArray &side_normals )
158 Intrepid::CellTools<Scalar>::getPhysicalSideNormals(
159 side_normals, this->d_jacobian, d_side_id, d_parent_topology );
168 const MDArray ¶metric_coords, MDArray &side_normals )
170 int num_cells = d_cell_node_coords.dimension( 0 );
171 int num_point = parametric_coords.dimension( 0 );
172 int space_dim = parametric_coords.dimension( 1 );
173 MDArray jacobian( num_cells, num_point, space_dim, space_dim );
175 Intrepid::CellTools<Scalar>::setJacobian( d_jacobian, parametric_coords,
176 this->d_cell_node_coords,
179 Intrepid::CellTools<Scalar>::getPhysicalSideNormals(
180 side_normals, d_jacobian, d_side_id, d_parent_topology );
void getPhysicalSideNormalsAtReferencePoint(const MDArray ¶metric_coords, MDArray &side_normals)
Compute the physical normals of the side at a given reference point.
Manager for Intrepid cell-level operations on cell sides.
void getPhysicalSideNormalsAtIntegrationPoints(MDArray &side_normals)
Compute the physical normals of the side.
Manager for Intrepid cell-level operations.
Assertions and Design-by-Contract for error handling.
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 updateCellState()
Update the cell state of the object for the current cell node coordinates.
IntrepidSideCell(const shards::CellTopology &side_topology, const unsigned side_id, const shards::CellTopology &parent_topology, const unsigned degree)
Constructor.