41 #ifndef DTK_PROJECTIONPRIMITIVES_HPP 42 #define DTK_PROJECTIONPRIMITIVES_HPP 48 #include <Teuchos_ParameterList.hpp> 50 #include <Shards_CellTopology.hpp> 52 #include <Intrepid_Basis.hpp> 53 #include <Intrepid_FieldContainer.hpp> 68 Intrepid::FieldContainer<double> ¢er );
72 const Teuchos::ParameterList ¶meters,
73 const Intrepid::FieldContainer<double> &point,
74 const Intrepid::FieldContainer<double> &face_nodes,
75 const Intrepid::FieldContainer<double> &face_node_normals,
76 const shards::CellTopology &face_topology );
82 const Teuchos::ParameterList ¶meters,
83 const Intrepid::FieldContainer<double> &point,
84 const Intrepid::FieldContainer<double> &face_nodes,
85 const Intrepid::FieldContainer<double> &face_node_normals,
86 const shards::CellTopology &face_topology,
87 Intrepid::FieldContainer<double> ¶metric_point,
88 Intrepid::FieldContainer<double> &physical_point,
int &face_edge_id,
92 static bool projectPointFeatureToEdgeFeature(
93 const Teuchos::ParameterList ¶meters,
94 const Intrepid::FieldContainer<double> &point,
95 const Intrepid::FieldContainer<double> &point_normal,
96 const Intrepid::FieldContainer<double> &edge_nodes,
97 const Intrepid::FieldContainer<double> &edge_node_normals,
98 Intrepid::FieldContainer<double> &projected_point,
int &edge_node_id );
103 const Teuchos::ParameterList ¶meters,
104 const Intrepid::FieldContainer<double> &edge_1,
105 const Intrepid::FieldContainer<double> &edge_2,
106 const Intrepid::FieldContainer<double> &edge_2_node_normals,
107 Intrepid::FieldContainer<double> &edge_1_intersection,
108 Intrepid::FieldContainer<double> &edge_2_intersection,
109 int &edge_1_node_id,
int &edge_2_node_id );
114 static double distanceToFaceBilinearSurface(
115 const Teuchos::ParameterList ¶meters,
116 const Intrepid::FieldContainer<double> &point,
117 const Intrepid::FieldContainer<double> &face_edge_nodes,
118 const Intrepid::FieldContainer<double> &face_edge_node_normals );
127 #endif // end DTK_PROJECTIONPRIMITIVES_HPP static void referenceCellCenter(const shards::CellTopology &cell_topo, Intrepid::FieldContainer< double > ¢er)
Get the center of the reference cell of the given topology.
Assertions and Design-by-Contract for error handling.
static void projectPointToFace(const Teuchos::ParameterList ¶meters, const Intrepid::FieldContainer< double > &point, const Intrepid::FieldContainer< double > &face_nodes, const Intrepid::FieldContainer< double > &face_node_normals, const shards::CellTopology &face_topology, Intrepid::FieldContainer< double > ¶metric_point, Intrepid::FieldContainer< double > &physical_point, int &face_edge_id, int &face_node_id)
Project a point onto a face and return the physical coordinates of the projected point on that face...
A stateless class of projection primitive operations.
static bool edgeEdgeIntersection(const Teuchos::ParameterList ¶meters, const Intrepid::FieldContainer< double > &edge_1, const Intrepid::FieldContainer< double > &edge_2, const Intrepid::FieldContainer< double > &edge_2_node_normals, Intrepid::FieldContainer< double > &edge_1_intersection, Intrepid::FieldContainer< double > &edge_2_intersection, int &edge_1_node_id, int &edge_2_node_id)
Intersect two edges in 3 dimensions and return their intersection point realized on both edges...
static bool pointInFaceVolumeOfInfluence(const Teuchos::ParameterList ¶meters, const Intrepid::FieldContainer< double > &point, const Intrepid::FieldContainer< double > &face_nodes, const Intrepid::FieldContainer< double > &face_node_normals, const shards::CellTopology &face_topology)
Determine if a point is within the volume of influence of a face.