Skip to content

Commit

Permalink
Merge Pull Request trilinos#6482 from CamelliaDPG/Trilinos/intrepid2-…
Browse files Browse the repository at this point in the history
…add-hierarchical-bases-simplices

Automatically Merged using Trilinos Pull Request AutoTester
PR Title: Intrepid2: add support for hierarchical H(grad) bases on triangle and tetrahedron.
PR Author: CamelliaDPG
  • Loading branch information
trilinos-autotester authored Apr 13, 2020
2 parents 2827cd5 + ef5553b commit a9bd7f9
Show file tree
Hide file tree
Showing 21 changed files with 2,579 additions and 119 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -359,7 +359,7 @@ namespace Intrepid2 {
const PointViewType /* inputPoints */,
const EOperator /* operatorType */ = OPERATOR_VALUE ) const {
INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
">>> ERROR (Basis::getValues): this method (FEM) is not supported or should be over-riden accordingly by derived classes.");
">>> ERROR (Basis::getValues): this method (FEM) is not supported or should be overridden accordingly by derived classes.");
}

/** \brief Evaluation of an FVD basis evaluation on a <strong>physical cell</strong>.
Expand Down Expand Up @@ -388,7 +388,7 @@ namespace Intrepid2 {
const PointViewType /* cellVertices */,
const EOperator /* operatorType */ = OPERATOR_VALUE ) const {
INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
">>> ERROR (Basis::getValues): this method (FVM) is not supported or should be over-riden accordingly by derived classes.");
">>> ERROR (Basis::getValues): this method (FVM) is not supported or should be overridden accordingly by derived classes.");
}


Expand All @@ -399,7 +399,7 @@ namespace Intrepid2 {
void
getDofCoords( ScalarViewType /* dofCoords */ ) const {
INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
">>> ERROR (Basis::getDofCoords): this method is not supported or should be over-riden accordingly by derived classes.");
">>> ERROR (Basis::getDofCoords): this method is not supported or should be overridden accordingly by derived classes.");
}

/** \brief Coefficients for computing degrees of freedom for Lagrangian basis
Expand All @@ -414,7 +414,7 @@ namespace Intrepid2 {
void
getDofCoeffs( ScalarViewType /* dofCoeffs */ ) const {
INTREPID2_TEST_FOR_EXCEPTION( true, std::logic_error,
">>> ERROR (Basis::getDofCoeffs): this method is not supported or should be over-riden accordingly by derived classes.");
">>> ERROR (Basis::getDofCoeffs): this method is not supported or should be overridden accordingly by derived classes.");
}

/** \brief For hierarchical bases, returns the field ordinals that have at most the specified degree in each dimension.
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -63,12 +63,22 @@

namespace Intrepid2
{
//! \brief EmptyBasisFamily allows us to set a default void family for a given topology
class EmptyBasisFamily
{
public:
using HGRAD = void;
using HCURL = void;
using HDIV = void;
using HVOL = void;
};

/** \class Intrepid2::DerivedBasisFamily
\brief A family of basis functions, constructed from H(vol) and H(grad) bases on the line.
At present, only hypercube topologies (line, quadrilateral, hexahedron) are supported, but other topologies will be supported in the future.
*/
template<class LineBasisHGRAD, class LineBasisHVOL>
template<class LineBasisHGRAD, class LineBasisHVOL, class TriangleBasisFamily = EmptyBasisFamily, class TetrahedronBasisFamily = EmptyBasisFamily>
class DerivedBasisFamily
{
public:
Expand All @@ -89,11 +99,17 @@ namespace Intrepid2
using HDIV_QUAD = Basis_Derived_HDIV_QUAD <HGRAD_LINE, HVOL_LINE>;
using HVOL_QUAD = Basis_Derived_HVOL_QUAD <HVOL_LINE>;

// hexahedral bases
// hexahedron bases
using HGRAD_HEX = Basis_Derived_HGRAD_HEX<HGRAD_LINE>;
using HCURL_HEX = Basis_Derived_HCURL_HEX<HGRAD_LINE, HVOL_LINE>;
using HDIV_HEX = Basis_Derived_HDIV_HEX <HGRAD_LINE, HVOL_LINE>;
using HVOL_HEX = Basis_Derived_HVOL_HEX <HVOL_LINE>;

// triangle bases
using HGRAD_TRI = typename TriangleBasisFamily::HGRAD;

// tetrahedron bases
using HGRAD_TET = typename TetrahedronBasisFamily::HGRAD;
};

/** \brief Factory method for line bases in the given family.
Expand Down Expand Up @@ -152,12 +168,12 @@ namespace Intrepid2
}
}

/** \brief Factory method for isotropic hexahedral bases in the given family.
/** \brief Factory method for isotropic bases on the hexahedron in the given family.
\param [in] fs - the function space for the basis.
\param [in] polyOrder - the polynomial order of the basis.
*/
template<class BasisFamily>
static typename BasisFamily::BasisPtr getHexahedralBasis(Intrepid2::EFunctionSpace fs, int polyOrder)
static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs, int polyOrder)
{
using Teuchos::rcp;
switch (fs)
Expand All @@ -171,14 +187,14 @@ namespace Intrepid2
}
}

/** \brief Factory method for potentially anisotropic hexahedral bases in the given family.
/** \brief Factory method for potentially anisotropic hexahedron bases in the given family.
\param [in] fs - the function space for the basis.
\param [in] polyOrder_x - the polynomial order of the basis in the x dimension.
\param [in] polyOrder_y - the polynomial order of the basis in the y dimension.
\param [in] polyOrder_z - the polynomial order of the basis in the z dimension.
*/
template<class BasisFamily>
static typename BasisFamily::BasisPtr getHexahedralBasis(Intrepid2::EFunctionSpace fs, int polyOrder_x, int polyOrder_y, int polyOrder_z)
static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs, int polyOrder_x, int polyOrder_y, int polyOrder_z)
{
using Teuchos::rcp;
switch (fs)
Expand All @@ -192,6 +208,44 @@ namespace Intrepid2
}
}

/** \brief Factory method for isotropic tetrahedron bases in the given family.
\param [in] fs - the function space for the basis.
\param [in] polyOrder - the polynomial order of the basis.
*/
template<class BasisFamily>
static typename BasisFamily::BasisPtr getTetrahedronBasis(Intrepid2::EFunctionSpace fs, int polyOrder)
{
using Teuchos::rcp;
switch (fs)
{
// case FUNCTION_SPACE_HVOL: return rcp(new typename BasisFamily::HVOL_TET (polyOrder));
// case FUNCTION_SPACE_HCURL: return rcp(new typename BasisFamily::HCURL_TET(polyOrder));
// case FUNCTION_SPACE_HDIV: return rcp(new typename BasisFamily::HDIV_TET (polyOrder));
case FUNCTION_SPACE_HGRAD: return rcp(new typename BasisFamily::HGRAD_TET(polyOrder));
default:
INTREPID2_TEST_FOR_EXCEPTION(true, std::invalid_argument, "Unsupported function space");
}
}

/** \brief Factory method for isotropic triangle bases in the given family.
\param [in] fs - the function space for the basis.
\param [in] polyOrder - the polynomial order of the basis.
*/
template<class BasisFamily>
static typename BasisFamily::BasisPtr getTriangleBasis(Intrepid2::EFunctionSpace fs, int polyOrder)
{
using Teuchos::rcp;
switch (fs)
{
// case FUNCTION_SPACE_HVOL: return rcp(new typename BasisFamily::HVOL_TRI (polyOrder));
// case FUNCTION_SPACE_HCURL: return rcp(new typename BasisFamily::HCURL_TRI(polyOrder));
// case FUNCTION_SPACE_HDIV: return rcp(new typename BasisFamily::HDIV_TRI (polyOrder));
case FUNCTION_SPACE_HGRAD: return rcp(new typename BasisFamily::HGRAD_TRI(polyOrder));
default:
INTREPID2_TEST_FOR_EXCEPTION(true, std::invalid_argument, "Unsupported function space");
}
}

/** \brief Factory method for isotropic bases in the given family on the specified cell topology.
\param [in] cellTopo - the cell topology on which the basis is defined.
\param [in] fs - the function space for the basis.
Expand All @@ -208,7 +262,9 @@ namespace Intrepid2
{
case shards::Line<>::key: return getLineBasis<BasisFamily>(fs,polyOrder);
case shards::Quadrilateral<>::key: return getQuadrilateralBasis<BasisFamily>(fs,polyOrder);
case shards::Hexahedron<>::key: return getHexahedralBasis<BasisFamily>(fs,polyOrder);
case shards::Triangle<>::key: return getTriangleBasis<BasisFamily>(fs,polyOrder);
case shards::Hexahedron<>::key: return getHexahedronBasis<BasisFamily>(fs,polyOrder);
case shards::Tetrahedron<>::key: return getTetrahedronBasis<BasisFamily>(fs,polyOrder);
default:
INTREPID2_TEST_FOR_EXCEPTION(true, std::invalid_argument, "Unsupported cell topology");
}
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -106,13 +106,13 @@ namespace Intrepid2
*/
virtual
const char*
getName() const {
getName() const override {
return "Intrepid2_DerivedBasis_HGRAD_HEX";
}

/** \brief True if orientation is required
*/
virtual bool requireOrientation() const {
virtual bool requireOrientation() const override {
return (this->getDegree() > 2);
}

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -94,13 +94,13 @@ namespace Intrepid2
*/
virtual
const char*
getName() const {
getName() const override {
return "Intrepid2_DerivedBasis_HGRAD_QUAD";
}

/** \brief True if orientation is required
*/
virtual bool requireOrientation() const {
virtual bool requireOrientation() const override {
return (this->getDegree() > 2);
}

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -169,9 +169,9 @@ namespace Intrepid2 {
}

/** \class Intrepid2::Basis_HCURL_HEX_In_FEM
\brief Implementation of the default H(curl)-compatible FEM basis on Hexahedral cell
\brief Implementation of the default H(curl)-compatible FEM basis on Hexahedron cell
Implements Nedelec basis of degree n on the reference Hexahedral cell. The basis has
Implements Nedelec basis of degree n on the reference Hexahedron cell. The basis has
cardinality 3n (n+1)^2 and spans a INCOMPLETE polynomial space.
*/
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -158,9 +158,9 @@ namespace Intrepid2 {
}

/** \class Intrepid2::Basis_HDIV_HEX_In_FEM
\brief Implementation of the default H(div)-compatible FEM basis on Hexahedral cell
\brief Implementation of the default H(div)-compatible FEM basis on Hexahedron cell
Implements Raviart-Thomas basis of degree n on the reference Hexahedral cell. The basis has
Implements Raviart-Thomas basis of degree n on the reference Hexahedron cell. The basis has
cardinality 3(n+1)n^2 and spans a INCOMPLETE polynomial space.
*/
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -241,7 +241,7 @@ namespace Intrepid2 {
dofCoords("dofCoordsHost", this->basisCardinality_,this->basisCellTopology_.getDimension());

dofCoords(0,0) = 1.0/3.0; dofCoords(0,1) = 0.0; dofCoords(0,2) = 1.0/3.0;
dofCoords(1,0) = 1.0/3.0; dofCoords(1,1) = 1.0/3.0; dofCoords(1,2) = 1.0/3.0;;
dofCoords(1,0) = 1.0/3.0; dofCoords(1,1) = 1.0/3.0; dofCoords(1,2) = 1.0/3.0;
dofCoords(2,0) = 0.0; dofCoords(2,1) = 1.0/3.0; dofCoords(2,2) = 1.0/3.0;
dofCoords(3,0) = 1.0/3.0; dofCoords(3,1) = 1.0/3.0; dofCoords(3,2) = 0.0;

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -204,7 +204,7 @@ namespace Intrepid2 {

output.access(4, 0) = 2.*z*(1. + z);
output.access(4, 1) = 0.;
output.access(4, 2) = ((-1. + 4.*x)*(1. + 2.*z))/2.;;
output.access(4, 2) = ((-1. + 4.*x)*(1. + 2.*z))/2.;
output.access(4, 3) = 0.;
output.access(4, 4) = 0.;
output.access(4, 5) = x*(-1. + 2.*x);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -124,7 +124,7 @@ namespace Intrepid2 {
}

/** \class Intrepid2::Basis_HVOL_C0_FEM
\brief Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedral and tetrahedral cells.
\brief Implementation of the default HVOL-compatible FEM contstant basis on triangle, quadrilateral, hexahedron and tetrahedron cells.
*/
template<typename ExecSpaceType = void,
typename outputValueType = double,
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -52,12 +52,42 @@
#include "Intrepid2_DerivedBasisFamily.hpp"

#include "Intrepid2_IntegratedLegendreBasis_HGRAD_LINE.hpp"
#include "Intrepid2_IntegratedLegendreBasis_HGRAD_TRI.hpp"
#include "Intrepid2_IntegratedLegendreBasis_HGRAD_TET.hpp"
#include "Intrepid2_LegendreBasis_HVOL_LINE.hpp"

namespace Intrepid2 {
// the following defines a family of hierarchical basis functions that matches the unpermuted ESEAS basis functions
// each basis member is associated with appropriate subcell topologies, making this suitable for continuous Galerkin finite elements.

template<typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
typename OutputScalar = double,
typename PointScalar = double,
bool defineVertexFunctions = true>
class HierarchicalTriangleBasisFamily
{
public:
// we will fill these in as we implement them
using HGRAD = IntegratedLegendreBasis_HGRAD_TRI<ExecutionSpace,OutputScalar,PointScalar,defineVertexFunctions>;
using HCURL = void;
using HDIV = void;
using HVOL = void;
};

template<typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
typename OutputScalar = double,
typename PointScalar = double,
bool defineVertexFunctions = true>
class HierarchicalTetrahedronBasisFamily
{
public:
// we will fill these in as we implement them
using HGRAD = IntegratedLegendreBasis_HGRAD_TET<ExecutionSpace,OutputScalar,PointScalar,defineVertexFunctions>;
using HCURL = void;
using HDIV = void;
using HVOL = void;
};

/** \class Intrepid2::HierarchicalBasisFamily
\brief A family of hierarchical basis functions, constructed in a way that follows work by Fuentes et al.
Expand All @@ -82,7 +112,10 @@ namespace Intrepid2 {
typename OutputScalar = double,
typename PointScalar = double>
using HierarchicalBasisFamily = DerivedBasisFamily< IntegratedLegendreBasis_HGRAD_LINE<ExecutionSpace,OutputScalar,PointScalar,true>,
LegendreBasis_HVOL_LINE<ExecutionSpace,OutputScalar,PointScalar> >;
LegendreBasis_HVOL_LINE<ExecutionSpace,OutputScalar,PointScalar>,
HierarchicalTriangleBasisFamily<ExecutionSpace,OutputScalar,PointScalar>,
HierarchicalTetrahedronBasisFamily<ExecutionSpace,OutputScalar,PointScalar>
>;

/** \class Intrepid2::HierarchicalBasisFamily
\brief A family of hierarchical basis functions, constructed in a way that follows work by Fuentes et al., suitable for use in DG contexts.
Expand All @@ -96,7 +129,10 @@ namespace Intrepid2 {
typename OutputScalar = double,
typename PointScalar = double>
using DGHierarchicalBasisFamily = DerivedBasisFamily< IntegratedLegendreBasis_HGRAD_LINE<ExecutionSpace,OutputScalar,PointScalar,false>,
LegendreBasis_HVOL_LINE<ExecutionSpace,OutputScalar,PointScalar> >;
LegendreBasis_HVOL_LINE<ExecutionSpace,OutputScalar,PointScalar>,
HierarchicalTriangleBasisFamily<ExecutionSpace,OutputScalar,PointScalar,false>,
HierarchicalTetrahedronBasisFamily<ExecutionSpace,OutputScalar,PointScalar,false>
>;

}

Expand Down
Loading

0 comments on commit a9bd7f9

Please sign in to comment.