DataTransferKit - Multiphysics Solution Transfer Services  2.0
DTK_IntrepidSideCell.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_INTREPIDSIDECELL_HPP
42 #define DTK_INTREPIDSIDECELL_HPP
43 
44 #include "DTK_IntrepidCell.hpp"
45 
46 #include <Teuchos_RCP.hpp>
47 
48 #include <Shards_CellTopology.hpp>
49 
50 namespace DataTransferKit
51 {
52 //---------------------------------------------------------------------------//
57 //---------------------------------------------------------------------------//
59 {
60  public:
62  typedef IntrepidCell Base;
64  typedef Base::Scalar Scalar;
65  typedef Base::MDArray MDArray;
67 
68  public:
69  // Constructor.
70  IntrepidSideCell( const shards::CellTopology &side_topology,
71  const unsigned side_id,
72  const shards::CellTopology &parent_topology,
73  const unsigned degree );
74 
75  // Update the cell state of the object for the current cell node
76  // coordinates.
77  void updateCellState();
78 
79  // Given a set of coordinates in the reference frame of the side of the
80  // parent cell, map them to the physical frame.
81  void mapToCellPhysicalFrame( const MDArray &parametric_coords,
82  MDArray &physical_coords );
83 
84  // Compute the physical normals of the side at the integration points
85  // (IP,DIM).
86  void getPhysicalSideNormalsAtIntegrationPoints( MDArray &side_normals );
87 
88  // Compute the physical normals of the side at a given reference point.
89  void
90  getPhysicalSideNormalsAtReferencePoint( const MDArray &parametric_coords,
91  MDArray &side_normals );
92 
93  private:
94  // Reference cell side id on the parent cell.
95  unsigned d_side_id;
96 
97  // Parent cell topology.
98  shards::CellTopology d_parent_topology;
99 };
100 
101 //---------------------------------------------------------------------------//
102 
103 } // end namespace DataTransferKit
104 
105 //---------------------------------------------------------------------------//
106 
107 #endif // end DTK_INTREPIDSIDECELL_HPP
108 
109 //---------------------------------------------------------------------------//
110 // end DTK_IntrepidSideCell.hpp
111 //---------------------------------------------------------------------------//
void getPhysicalSideNormalsAtReferencePoint(const MDArray &parametric_coords, MDArray &side_normals)
Compute the physical normals of the side at a given reference point.
void getPhysicalSideNormalsAtIntegrationPoints(MDArray &side_normals)
Compute the physical normals of the side.
Manager for Intrepid cell-level operations.
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...
Manager for Intrepid cell-level operations.
void updateCellState()
Update the cell state of the object for the current cell node coordinates.
Manager for Intrepid cell-level operations on cell sides.
DTK_BasicEntitySet.cpp.
IntrepidSideCell(const shards::CellTopology &side_topology, const unsigned side_id, const shards::CellTopology &parent_topology, const unsigned degree)
Constructor.