DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_IntrepidCell.hpp
Go to the documentation of this file.
1 //---------------------------------------------------------------------------//
2 /*
3  Copyright (c) 2012, Stuart R. Slattery
4  All rights reserved.
5 
6  Redistribution and use in source and binary forms, with or without
7  modification, are permitted provided that the following conditions are
8  met:
9 
10  *: Redistributions of source code must retain the above copyright
11  notice, this list of conditions and the following disclaimer.
12 
13  *: Redistributions in binary form must reproduce the above copyright
14  notice, this list of conditions and the following disclaimer in the
15  documentation and/or other materials provided with the distribution.
16 
17  *: Neither the name of the University of Wisconsin - Madison nor the
18  names of its contributors may be used to endorse or promote products
19  derived from this software without specific prior written permission.
20 
21  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
24  A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
25  HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
26  SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
27  LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
28  DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
29  THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30  (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
31  OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 */
33 //---------------------------------------------------------------------------//
39 //---------------------------------------------------------------------------//
40 
41 #ifndef DTK_INTREPIDCELL_HPP
42 #define DTK_INTREPIDCELL_HPP
43 
44 #include <Teuchos_RCP.hpp>
45 
46 #include <Shards_CellTopology.hpp>
47 
48 #include <Intrepid_Cubature.hpp>
49 #include <Intrepid_FieldContainer.hpp>
50 
51 namespace DataTransferKit
52 {
53 //---------------------------------------------------------------------------//
58 //---------------------------------------------------------------------------//
60 {
61  public:
63  typedef double Scalar;
64  typedef Intrepid::FieldContainer<Scalar> MDArray;
66 
67  public:
68  // Default constructor.
69  IntrepidCell() { /* ... */}
70 
71  // Constructor.
72  IntrepidCell( const shards::CellTopology &cell_topology,
73  const unsigned degree );
74 
75  // Destructor.
76  virtual ~IntrepidCell();
77 
78  // Given physical coordinates for the cell nodes (Cell,Node,Dim), assign
79  // them to the cell without allocating internal data.
80  void setCellNodeCoordinates( const MDArray &cell_node_coords );
81 
82  // Given physical coordinates for the cell nodes (Cell,Node,Dim),
83  // allocate the state of the cell object.
84  void allocateCellState( const MDArray &cell_node_coords );
85 
86  // Update the cell state of the object for the current cell node
87  // coordinates.
88  virtual void updateCellState();
89 
90  // Free function for updating the cell state for a new set of
91  // physical cells in a single call.
92  static void updateState( IntrepidCell &Intrepid_cell,
93  const MDArray &cell_node_coords );
94 
95  // Given a set of coordinates in the physical frame of the cell
96  // (Node,Dim), map them to the reference frame (Cell,Node,Dim).
97  void mapToCellReferenceFrame( const MDArray &physical_coords,
98  MDArray &reference_coords );
99 
100  // Given a set of coordinates in the reference frame of the cell
101  // (Node,Dim), map them to the physical frame (Cell,Node,Dim).
102  virtual void mapToCellPhysicalFrame( const MDArray &parametric_coords,
103  MDArray &physical_coords );
104 
105  // Determine if a point given in parametric coordinates is inside of the
106  // reference cell.
107  bool pointInReferenceCell( const MDArray &reference_point,
108  const double tolerance );
109 
110  // Determine if a point in physical coordinates is inside of the physical
111  // cell.
112  bool pointInPhysicalCell( const MDArray &point, const double tolerance );
113 
114  // Get the number of cells in the current state.
115  int getNumCells() const;
116 
117  // Get the number of cell points.
118  int getNumIntegrationPoints() const;
119 
120  // Get the spatial dimension.
121  int getSpatialDimension() const;
122 
123  // Get the cell measures (Cell). cell_measures must already be allocated.
124  void getCellMeasures( MDArray &cell_measures ) const;
125 
126  // Get the physical integration point coordinates in each cell
127  // (Cell,IP,Dim).
128  void
129  getPhysicalIntegrationCoordinates( MDArray &physical_ip_coordinates ) const;
130 
131  // Given DOFs at the quadrature points {(Cell,Node) for scalar fields,
132  // (Cell,Node,VecDim) for vector fields, and (Cell,Node,TensDim1,TensDim2)
133  // for tensor fields.} perform the numerical integration in each cell by
134  // contracting them with the weighted measures.
135  void integrate( const MDArray &dofs, MDArray &integrals ) const;
136 
137  protected:
138  // Cell topology.
139  shards::CellTopology d_topology;
140 
141  // Cell integration rule.
142  Teuchos::RCP<Intrepid::Cubature<Scalar, MDArray>> d_cubature;
143 
144  // Cubature points (IP,Dim).
145  MDArray d_cub_points;
146 
147  // Cubature weights (IP).
148  MDArray d_cub_weights;
149 
150  // Jacobian (Cell,IP,Dim,Dim).
151  MDArray d_jacobian;
152 
153  // Jacobian determinant (Cell,IP).
154  MDArray d_jacobian_det;
155 
156  // Weighted cell measures (Cell,IP).
157  MDArray d_weighted_measures;
158 
159  // Current physical cell node coordinates (Cell,Node,Dim).
160  MDArray d_cell_node_coords;
161 
162  // Current physical integration point coordinates (Cell,IP,Dim).
163  MDArray d_physical_ip_coordinates;
164 };
165 
166 //---------------------------------------------------------------------------//
167 
168 } // end namespace DataTransferKit
169 
170 //---------------------------------------------------------------------------//
171 
172 #endif // end DTK_INTREPIDCELL_HPP
173 
174 //---------------------------------------------------------------------------//
175 // end DTK_IntrepidCell.hpp
176 //---------------------------------------------------------------------------//
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 &parametric_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...
DTK_BasicEntitySet.cpp.
int getNumCells() const
Get the number of cells in the current state.