Skip to content

Commit

Permalink
Intrepid2: fixes discussed in PR trilinos#7201
Browse files Browse the repository at this point in the history
  • Loading branch information
mperego committed Apr 25, 2020
1 parent 3463bfb commit 31b0762
Show file tree
Hide file tree
Showing 2 changed files with 56 additions and 58 deletions.
53 changes: 20 additions & 33 deletions packages/intrepid2/src/Projection/Intrepid2_ProjectionStruct.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -93,7 +93,12 @@ class ProjectionStruct {
typedef typename Kokkos::Impl::is_space<SpT>::host_mirror_space::execution_space host_space_type;
typedef Kokkos::DynRankView<ValueType,host_space_type> view_type;
typedef Kokkos::View<range_type**,host_space_type> range_tag;
typedef std::array<std::array<view_type, 12>,4> view_tag;
static constexpr int numberSubCellDims = 4; //{0 for vertex, 1 for edges, 2 for faces, 3 for volumes}
//max of numVertices, numEdges, numFaces for a reference cell.
//12 is the number of edges in a Hexahderon.
//We'll need to change this if we consider generic polyhedra
static constexpr int maxSubCellsCount = 12;
typedef std::array<std::array<view_type, maxSubCellsCount>, numberSubCellDims> view_tag;
typedef Kokkos::View<unsigned**,host_space_type> key_tag;

/** \brief Returns number of basis evaluation points
Expand Down Expand Up @@ -322,23 +327,19 @@ class ProjectionStruct {
}


/** \brief Returns the range of the basis evaluation points corresponding to a subcell
\param subCellDim [in] - dimension of the subcell
\param subCellId [in] - ordinal of the subcell defined by cell topology
/** \brief Returns the range tag of the basis evaluation points subcells
\return the range of the basis evaluation points corresponding to the selected subcell
\return the range tag of the basis evaluation points on subcells
*/
const range_tag getBasisPointsRange() const {
return basisPointsRange;
}


/** \brief Returns the range of the basis evaluation points corresponding to a subcell
\param subCellDim [in] - dimension of the subcell
\param subCellId [in] - ordinal of the subcell defined by cell topology
/** \brief Returns the range tag of the basis/target evaluation points in subcells
\param evalPointType [in] - enum selecting whether the points should be computed for the basis
functions or for the target function
\return the range of the basis/target evaluation points corresponding to the selected subcell
\return the range tag of the basis/target evaluation points on subcells
*/
const range_tag getPointsRange(const EvalPointsType type) const {
if(type == BASIS)
Expand All @@ -348,23 +349,19 @@ class ProjectionStruct {
}


/** \brief Returns the range of the derivative evaluation points corresponding to a subcell
\param subCellDim [in] - dimension of the subcell
\param subCellId [in] - ordinal of the subcell defined by cell topology
/** \brief Returns the range tag of the derivative evaluation points on subcell
\return the range of the basis derivative evaluation points corresponding to the selected subcell
\return the range tag of the basis derivative evaluation points corresponding on subcell
*/
const range_tag getBasisDerivPointsRange() const {
return basisDerivPointsRange;
}

/** \brief Returns the range of the derivative evaluation points corresponding to a subcell
\param subCellDim [in] - dimension of the subcell
\param subCellId [in] - ordinal of the subcell defined by cell topology
/** \brief Returns the range tag of the basis/target derivative evaluation points on subcells
\param evalPointType [in] - enum selecting whether the points should be computed for the basis
functions or for the target function
\return the range of the basis derivative evaluation points corresponding to the selected subcell
\return the range tag of the basis/target derivative evaluation points on subcells
*/
const range_tag getDerivPointsRange(const EvalPointsType type) const {
if(type == BASIS)
Expand All @@ -374,36 +371,26 @@ class ProjectionStruct {
}


/** \brief Returns the range of the target function evaluation points corresponding to a subcell
\param subCellDim [in] - dimension of the subcell
\param subCellId [in] - ordinal of the subcell defined by cell topology
/** \brief Returns the range of the target function evaluation points on subcells
\return the range of the target function evaluation points corresponding to the selected subcell
\return the range of the target function evaluation points corresponding on subcells
*/
const range_tag getTargetPointsRange() const {
return targetPointsRange;
}


/** \brief Returns the range of the target function derivative evaluation points corresponding to a subcell
\param subCellDim [in] - dimension of the subcell
\param subCellId [in] - ordinal of the subcell defined by cell topology
/** \brief Returns the range tag of the target function derivative evaluation points on subcells
\return the range of the target function derivative evaluation points corresponding to the selected subcell
\return the range of the target function derivative evaluation points corresponding on subcells
*/
range_type getTargetDerivPointsRange(const ordinal_type subCellDim, const ordinal_type subCellId) {
return targetDerivPointsRange(subCellDim,subCellId);
}

const range_tag getTargetDerivPointsRange() const {
return targetDerivPointsRange;
}

/** \brief Returns the key of a subcell topology
\param subCellDim [in] - dimension of the subcell
\param subCellId [in] - ordinal of the subcell defined by cell topology
/** \brief Returns the key tag for subcells
\return the key of the selected subcell
\return the key tag of the selected subcells
*/
const key_tag getTopologyKey() const {
return subCellTopologyKey;
Expand Down
61 changes: 36 additions & 25 deletions packages/intrepid2/src/Projection/Intrepid2_ProjectionStructDef.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -85,15 +85,19 @@ void ProjectionStruct<SpT,ValueType>::createL2ProjectionStruct(const BasisPtrTyp

ordinal_type hasCellDofs = (cellBasis->getDofCount(dim, 0) > 0);

INTREPID2_TEST_FOR_ABORT( (numVertices > maxSubCellsCount) || (numFaces > maxSubCellsCount) || (numEdges > maxSubCellsCount),
">>> ERROR (Intrepid2::ProjectionStruct:createHDivProjectionStruct, Projections do not support a cell topology with this cub cells count");


numBasisEvalPoints += numVertices;
numTargetEvalPoints += numVertices;
view_type coord("vertex_coord", dim);

basisPointsRange = range_tag("basisPointsRange", 4,12);
targetPointsRange = range_tag("targetPointsRange", 4,12);
basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,12);
targetDerivPointsRange = range_tag("targetDerivPointsRange", 4,12);
subCellTopologyKey = key_tag("subCellTopologyKey",4,12);
basisPointsRange = range_tag("basisPointsRange", 4,maxSubCellsCount);
targetPointsRange = range_tag("targetPointsRange", 4,maxSubCellsCount);
basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,maxSubCellsCount);
targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);

maxNumBasisEvalPoints = numVertices; maxNumTargetEvalPoints = numVertices;
for(ordinal_type iv=0; iv<numVertices; ++iv) {
Expand Down Expand Up @@ -210,16 +214,20 @@ void ProjectionStruct<SpT,ValueType>::createHGradProjectionStruct(const BasisPtr
ordinal_type numFaces = (cellBasis->getDofCount(2, 0) > 0) ? cellTopo.getFaceCount(): 0;
ordinal_type numEdges = (cellBasis->getDofCount(1, 0) > 0) ? cellTopo.getEdgeCount() : 0;

INTREPID2_TEST_FOR_ABORT( (numFaces > maxSubCellsCount) || (numEdges > maxSubCellsCount),
">>> ERROR (Intrepid2::ProjectionStruct:createHDivProjectionStruct, Projections do not support a cell topology with this cub cells count");


ordinal_type hasCellDofs = (cellBasis->getDofCount(dim, 0) > 0);

maxNumBasisEvalPoints = numVertices; maxNumTargetEvalPoints = numVertices;
maxNumBasisDerivEvalPoints = 0; maxNumTargetDerivEvalPoints = 0;

basisPointsRange = range_tag("basisPointsRange", 4,12);
targetPointsRange = range_tag("targetPointsRange", 4,12);
basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,12);
targetDerivPointsRange = range_tag("targetDerivPointsRange", 4,12);
subCellTopologyKey = key_tag("subCellTopologyKey",4,12);
basisPointsRange = range_tag("basisPointsRange", numberSubCellDims,maxSubCellsCount);
targetPointsRange = range_tag("targetPointsRange", numberSubCellDims,maxSubCellsCount);
basisDerivPointsRange = range_tag("basisDerivPointsRange", numberSubCellDims,maxSubCellsCount);
targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);

numBasisEvalPoints += numVertices;
numTargetEvalPoints += numVertices;
Expand Down Expand Up @@ -325,11 +333,11 @@ void ProjectionStruct<SpT,ValueType>::createHCurlProjectionStruct(const BasisPtr
maxNumBasisEvalPoints = 0; maxNumTargetEvalPoints = 0;
maxNumBasisDerivEvalPoints = 0; maxNumTargetDerivEvalPoints = 0;

basisPointsRange = range_tag("basisPointsRange", 4,12);
targetPointsRange = range_tag("targetPointsRange", 4,12);
basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,12);
targetDerivPointsRange = range_tag("targetDerivPointsRange", 4,12);
subCellTopologyKey = key_tag("subCellTopologyKey",4,12);
basisPointsRange = range_tag("basisPointsRange", numberSubCellDims,maxSubCellsCount);
targetPointsRange = range_tag("targetPointsRange", numberSubCellDims,maxSubCellsCount);
basisDerivPointsRange = range_tag("basisDerivPointsRange", numberSubCellDims,maxSubCellsCount);
targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);

DefaultCubatureFactory cub_factory;
for(ordinal_type ie=0; ie<numEdges; ++ie) {
Expand Down Expand Up @@ -448,14 +456,17 @@ void ProjectionStruct<SpT,ValueType>::createHDivProjectionStruct(const BasisPtrT
ordinal_type numSides = cellTopo.getSideCount()*ordinal_type(cellBasis->getDofCount(sideDim, 0) > 0);
ordinal_type hasCellDofs = (cellBasis->getDofCount(dim, 0) > 0);

INTREPID2_TEST_FOR_ABORT( numSides > maxSubCellsCount,
">>> ERROR (Intrepid2::ProjectionStruct:createHDivProjectionStruct, Projections do not support a cell topology with so many sides");

maxNumBasisEvalPoints = 0; maxNumTargetEvalPoints = 0;
maxNumBasisDerivEvalPoints = 0; maxNumTargetDerivEvalPoints = 0;

basisPointsRange = range_tag("basisPointsRange", 4,12);
targetPointsRange = range_tag("targetPointsRange", 4,12);
basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,12);
targetDerivPointsRange = range_tag("targetDerivPointsRange", 4,12);
subCellTopologyKey = key_tag("subCellTopologyKey",4,12);
basisPointsRange = range_tag("basisPointsRange", numberSubCellDims,maxSubCellsCount);
targetPointsRange = range_tag("targetPointsRange", numberSubCellDims,maxSubCellsCount);
basisDerivPointsRange = range_tag("basisDerivPointsRange", numberSubCellDims,maxSubCellsCount);
targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);

Basis<host_space_type,ValueType,ValueType> *hcurlBasis = NULL;
if(cellTopo.getKey() == shards::getCellTopologyData<shards::Hexahedron<8> >()->key)
Expand Down Expand Up @@ -553,11 +564,11 @@ void ProjectionStruct<SpT,ValueType>::createHVolProjectionStruct(const BasisPtrT
numTargetEvalPoints = 0;
numTargetDerivEvalPoints = 0;

basisPointsRange = range_tag("basisPointsRange", 4,12);
targetPointsRange = range_tag("targetPointsRange", 4,12);
basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,12);
targetDerivPointsRange = range_tag("targetDerivPointsRange", 4,12);
subCellTopologyKey = key_tag("subCellTopologyKey",4,12);
basisPointsRange = range_tag("basisPointsRange", 4,maxSubCellsCount);
targetPointsRange = range_tag("targetPointsRange", 4,maxSubCellsCount);
basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,maxSubCellsCount);
targetDerivPointsRange = range_tag("targetDerivPointsRange", 4,maxSubCellsCount);
subCellTopologyKey = key_tag("subCellTopologyKey",4,maxSubCellsCount);

ordinal_type basisCubDegree = cellBasis->getDegree();
DefaultCubatureFactory cub_factory;
Expand Down

0 comments on commit 31b0762

Please sign in to comment.