DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_IntrepidCell.cpp
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 #include "DTK_DBC.hpp"
42 
43 #include "DTK_IntrepidCell.hpp"
44 
45 #include <Teuchos_as.hpp>
46 
47 #include <Intrepid_CellTools.hpp>
48 #include <Intrepid_DefaultCubatureFactory.hpp>
49 #include <Intrepid_FunctionSpaceTools.hpp>
50 #include <Intrepid_Types.hpp>
51 
52 namespace DataTransferKit
53 {
54 //---------------------------------------------------------------------------//
58 IntrepidCell::IntrepidCell( const shards::CellTopology &cell_topology,
59  const unsigned degree )
60  : d_topology( cell_topology )
61  , d_cub_points( 0, 0 )
62  , d_cub_weights( 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 )
67 {
68  Intrepid::DefaultCubatureFactory<Scalar, MDArray> cub_factory;
69  d_cubature = cub_factory.create( d_topology, degree );
70 
71  unsigned num_cub_points = d_cubature->getNumPoints();
72  unsigned cub_dim = cell_topology.getDimension();
73 
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 );
77 }
78 
79 //---------------------------------------------------------------------------//
84 
85 //---------------------------------------------------------------------------//
90 void IntrepidCell::setCellNodeCoordinates( const MDArray &cell_node_coords )
91 {
92  d_cell_node_coords = cell_node_coords;
93 }
94 
95 //---------------------------------------------------------------------------//
100 void IntrepidCell::allocateCellState( const MDArray &cell_node_coords )
101 {
102  // Store the cell node coords as the current state.
103  setCellNodeCoordinates( cell_node_coords );
104 
105  // Get required dimensions.
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 );
109 
110  // Resize arrays.
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 );
115 }
116 
117 //---------------------------------------------------------------------------//
123 {
124  // Compute the Jacobian.
125  Intrepid::CellTools<Scalar>::setJacobian( d_jacobian, d_cub_points,
126  d_cell_node_coords, d_topology );
127 
128  // Compute the determinant of the Jacobian.
129  Intrepid::CellTools<Scalar>::setJacobianDet( d_jacobian_det, d_jacobian );
130 
131  // Compute the cell measures.
132  Intrepid::FunctionSpaceTools::computeCellMeasure<Scalar>(
133  d_weighted_measures, d_jacobian_det, d_cub_weights );
134 
135  // Compute physical frame integration point coordinates.
136  Intrepid::CellTools<Scalar>::mapToPhysicalFrame(
137  d_physical_ip_coordinates, d_cub_points, d_cell_node_coords,
138  d_topology );
139 }
140 
141 //---------------------------------------------------------------------------//
147  const MDArray &cell_node_coords )
148 {
149  intrepid_cell.allocateCellState( cell_node_coords );
150  intrepid_cell.updateCellState();
151 }
152 
153 //---------------------------------------------------------------------------//
158 void IntrepidCell::mapToCellReferenceFrame( const MDArray &physical_coords,
159  MDArray &reference_coords )
160 {
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 );
165 }
166 
167 //---------------------------------------------------------------------------//
172 void IntrepidCell::mapToCellPhysicalFrame( const MDArray &parametric_coords,
173  MDArray &physical_coords )
174 {
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() ) );
185 
186  Intrepid::CellTools<Scalar>::mapToPhysicalFrame(
187  physical_coords, parametric_coords, d_cell_node_coords, d_topology );
188 }
189 
190 //---------------------------------------------------------------------------//
195 bool IntrepidCell::pointInReferenceCell( const MDArray &reference_point,
196  const double tolerance )
197 {
198  return Intrepid::CellTools<Scalar>::checkPointsetInclusion(
199  reference_point, d_topology, tolerance );
200 }
201 
202 //---------------------------------------------------------------------------//
207 bool IntrepidCell::pointInPhysicalCell( const MDArray &point,
208  const double tolerance )
209 {
210  MDArray reference_point( point );
211  reference_point.initialize();
212  mapToCellReferenceFrame( point, reference_point );
213  return pointInReferenceCell( reference_point, tolerance );
214 }
215 
216 //---------------------------------------------------------------------------//
221 {
222  return d_weighted_measures.dimension( 0 );
223 }
224 
225 //---------------------------------------------------------------------------//
230 {
231  return d_cub_points.dimension( 0 );
232 }
233 
234 //---------------------------------------------------------------------------//
239 {
240  return d_cub_points.dimension( 1 );
241 }
242 
243 //---------------------------------------------------------------------------//
248 void IntrepidCell::getCellMeasures( MDArray &cell_measures ) const
249 {
250  DTK_REQUIRE( 1 == cell_measures.rank() );
251  DTK_REQUIRE( cell_measures.dimension( 0 ) ==
252  d_weighted_measures.dimension( 0 ) );
253 
254  MDArray dofs( d_cell_node_coords.dimension( 0 ),
255  d_cub_weights.dimension( 0 ) );
256  dofs.initialize( 1.0 );
257  integrate( dofs, cell_measures );
258 }
259 
260 //---------------------------------------------------------------------------//
261 // Get the physical cell point coordinates in each cell
262 // (Cell,IP,Dim).
263 void IntrepidCell::getPhysicalIntegrationCoordinates(
264  MDArray &physical_ip_coordinates ) const
265 {
266  physical_ip_coordinates = d_physical_ip_coordinates;
267 }
268 
269 //---------------------------------------------------------------------------//
276 void IntrepidCell::integrate( const MDArray &dofs, MDArray &integrals ) const
277 {
278  Intrepid::FunctionSpaceTools::integrate<Scalar>(
279  integrals, dofs, d_weighted_measures, Intrepid::COMP_BLAS );
280 }
281 
282 //---------------------------------------------------------------------------//
283 
284 } // end namespace DataTransferKit
285 
286 //---------------------------------------------------------------------------//
287 // end DTK_IntrepidCell.cpp
288 //---------------------------------------------------------------------------//
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 &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.
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...
DTK_BasicEntitySet.cpp.
int getNumCells() const
Get the number of cells in the current state.