Intrepid2
Intrepid2_CubatureControlVolumeDef.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
16#ifndef __INTREPID2_CUBATURE_CONTROLVOLUME_DEF_HPP__
17#define __INTREPID2_CUBATURE_CONTROLVOLUME_DEF_HPP__
18
19namespace Intrepid2 {
20
21 template <typename DT, typename PT, typename WT>
23 CubatureControlVolume(const shards::CellTopology cellTopology) {
24
25 // define primary cell topology with given one
26 primaryCellTopo_ = cellTopology;
27
28 // subcell is defined either hex or quad according to dimension
29 switch (primaryCellTopo_.getDimension()) {
30 case 2:
31 subcvCellTopo_ = shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >());
32 break;
33 case 3:
34 subcvCellTopo_ = shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<8> >());
35 break;
36 }
37
38 // computation order is always one;
39 degree_ = 1;
40
41 // create subcell cubature points and weights and cache them
42 const ordinal_type subcvDegree = 2;
43 auto subcvCubature = DefaultCubatureFactory::create<DT,PT,WT>(subcvCellTopo_, subcvDegree);
44
45 const ordinal_type numSubcvPoints = subcvCubature->getNumPoints();
46 const ordinal_type subcvDim = subcvCubature->getDimension();
47
48 subcvCubaturePoints_ = Kokkos::DynRankView<PT,DT>("CubatureControlVolume::subcvCubaturePoints_",
49 numSubcvPoints, subcvDim);
50 subcvCubatureWeights_ = Kokkos::DynRankView<WT,DT>("CubatureControlVolume::subcvCubatureWeights_",
51 numSubcvPoints);
52
53 subcvCubature->getCubature(subcvCubaturePoints_,
54 subcvCubatureWeights_);
55 }
56
57 template <typename DT, typename PT, typename WT>
58 void
60 getCubature( PointViewType cubPoints,
61 weightViewType cubWeights,
62 PointViewType cellCoords ) const {
63#ifdef HAVE_INTREPID2_DEBUG
64 INTREPID2_TEST_FOR_EXCEPTION( cubPoints.rank() != 3, std::invalid_argument,
65 ">>> ERROR (CubatureControlVolume): cubPoints must have rank 3 (C,P,D).");
66 INTREPID2_TEST_FOR_EXCEPTION( cubWeights.rank() != 2, std::invalid_argument,
67 ">>> ERROR (CubatureControlVolume): cubWeights must have rank 2 (C,P).");
68 INTREPID2_TEST_FOR_EXCEPTION( cellCoords.rank() != 3, std::invalid_argument,
69 ">>> ERROR (CubatureControlVolume): cellCoords must have rank 3 of (C,P,D).");
70
71 INTREPID2_TEST_FOR_EXCEPTION( cubPoints.extent(0) != cellCoords.extent(0) ||
72 cubPoints.extent(0) != cubWeights.extent(0), std::invalid_argument,
73 ">>> ERROR (CubatureControlVolume): cubPoints, cubWeights and cellCoords dimension(0) are not consistent, numCells");
74
75 INTREPID2_TEST_FOR_EXCEPTION( cubPoints.extent(1) != cellCoords.extent(1) ||
76 cubPoints.extent(1) != cubWeights.extent(1), std::invalid_argument,
77 ">>> ERROR (CubatureControlVolume): cubPoints, cubWeights and cellCoords dimension(1) are not consistent, numNodesPerCell");
78
79 INTREPID2_TEST_FOR_EXCEPTION( cubPoints.extent(2) != cellCoords.extent(2) ||
80 static_cast<ordinal_type>(cubPoints.extent(2)) != getDimension(), std::invalid_argument,
81 ">>> ERROR (CubatureControlVolume): cubPoints, cellCoords, this->getDimension() are not consistent, spaceDim.");
82#endif
83 typedef Kokkos::DynRankView<PT,DT> tempPointViewType;
84
85 // get array dimensions
86 const ordinal_type numCells = cellCoords.extent(0);
87 const ordinal_type numNodesPerCell = cellCoords.extent(1);
88 const ordinal_type spaceDim = cellCoords.extent(2);
89
90 const ordinal_type numNodesPerSubcv = subcvCellTopo_.getNodeCount();
91 tempPointViewType subcvCoords("CubatureControlVolume::subcvCoords_",
92 numCells, numNodesPerCell, numNodesPerSubcv, spaceDim);
94 cellCoords,
95 primaryCellTopo_);
96
97 const ordinal_type numSubcvPoints = subcvCubaturePoints_.extent(0);
98 tempPointViewType subcvJacobian("CubatureControlVolume::subcvJacobian_",
99 numCells, numNodesPerCell, numSubcvPoints, spaceDim, spaceDim);
100
101 tempPointViewType subcvJacobianDet("CubatureControlVolume::subcvJacobDet_",
102 numCells, numNodesPerCell, numSubcvPoints);
103
104 // numNodesPerCell is maximum 8; this repeated run is necessary because of cell tools input consideration
105 for (ordinal_type node=0;node<numNodesPerCell;++node) {
106 auto subcvJacobianNode = Kokkos::subdynrankview(subcvJacobian, Kokkos::ALL(), node, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL());
107 auto subcvCoordsNode = Kokkos::subdynrankview(subcvCoords, Kokkos::ALL(), node, Kokkos::ALL(), Kokkos::ALL());
108 auto subcvJacobianDetNode = Kokkos::subdynrankview(subcvJacobianDet, Kokkos::ALL(), node, Kokkos::ALL());
109
110 CellTools<DT>::setJacobian(subcvJacobianNode, // C, P, D, D
111 subcvCubaturePoints_, // P, D
112 subcvCoordsNode, // C, N, D
113 subcvCellTopo_);
114
115 CellTools<DT>::setJacobianDet(subcvJacobianDetNode, // C, P
116 subcvJacobianNode); // C, P, D, D
117 }
118
119 //typedef typename ExecSpace<typename PointViewType::execution_space,DT>::ExecSpaceType ExecSpaceType;
120
121 const auto loopSize = numCells;
122 Kokkos::RangePolicy<typename DT::execution_space,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
123
125 Kokkos::parallel_for( policy, FunctorType(cubPoints,
126 cubWeights,
127 subcvCoords,
128 subcvCubatureWeights_,
129 subcvJacobianDet) );
130 }
131
132}
133
134#endif
135
static void setJacobianDet(JacobianDetViewType jacobianDet, const JacobianViewType jacobian)
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F.
static void getSubcvCoords(Kokkos::DynRankView< subcvCoordValueType, subcvCoordProperties... > subcvCoords, const Kokkos::DynRankView< cellCoordValueType, cellCoordProperties... > cellCoords, const shards::CellTopology primaryCell)
Computes coordinates of sub-control volumes in each primary cell.
static void setJacobian(JacobianViewType jacobian, const PointViewType points, const WorksetType worksetCell, const Teuchos::RCP< HGradBasisType > basis, const int startCell=0, const int endCell=-1)
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
virtual void getCubature(PointViewType cubPoints, weightViewType cubWeights, PointViewType cellCoords) const override
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
CubatureControlVolume(const shards::CellTopology cellTopology)