Skip to content

Commit

Permalink
[NL] Add a create function for common FE type.
Browse files Browse the repository at this point in the history
  • Loading branch information
endJunction committed May 27, 2019
1 parent e27b61e commit b63c384
Show file tree
Hide file tree
Showing 11 changed files with 41 additions and 51 deletions.
12 changes: 12 additions & 0 deletions NumLib/Fem/FiniteElement/TemplateIsoparametric.h
Original file line number Diff line number Diff line change
Expand Up @@ -157,4 +157,16 @@ class TemplateIsoparametric
const MeshElementType* _ele;
};

/// Creates a TemplateIsoparametric element for the given shape functions and
/// the underlying mesh element.
template <typename ShapeFunction, typename ShapeMatricesType>
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>
createIsoparametricFiniteElement(MeshLib::Element const& e)
{
using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;

return FemType{
*static_cast<const typename ShapeFunction::MeshElement*>(&e)};
}
} // namespace NumLib
5 changes: 3 additions & 2 deletions NumLib/Function/Interpolation.h
Original file line number Diff line number Diff line change
Expand Up @@ -108,9 +108,10 @@ void interpolateToHigherOrderNodes(
using SF = LowerOrderShapeFunction;
using ShapeMatricesType = ShapeMatrixPolicyType<SF, GlobalDim>;
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices;
using FemType = TemplateIsoparametric<SF, ShapeMatricesType>;

FemType fe(*static_cast<const typename SF::MeshElement*>(&element));
auto const fe =
NumLib::createIsoparametricFiniteElement<SF, ShapeMatricesType>(
element);
int const number_base_nodes = element.getNumberOfBaseNodes();
int const number_all_nodes = element.getNumberOfNodes();

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -89,12 +89,8 @@ class ConstraintDirichletBoundaryConditionLocalAssembler final
{
(void)local_matrix_size; // unused, but needed for the interface


using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(
&_surface_element));
auto const fe = NumLib::createIsoparametricFiniteElement<
ShapeFunction, ShapeMatricesType>(_surface_element);

auto const n_integration_points =
_integration_method.getNumberOfPoints();
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -52,10 +52,8 @@ class PythonBoundaryConditionLocalAssembler final
{
using ShapeMatricesType =
ShapeMatrixPolicyType<ShapeFunction, GlobalDim>;
using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;
FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(
&this->_element));
auto const fe = NumLib::createIsoparametricFiniteElement<
ShapeFunction, ShapeMatricesType>(Base::_element);

unsigned const num_integration_points =
Base::_integration_method.getNumberOfPoints();
Expand Down
7 changes: 2 additions & 5 deletions ProcessLib/ComponentTransport/ComponentTransportFEM.h
Original file line number Diff line number Diff line change
Expand Up @@ -809,11 +809,8 @@ class LocalAssemblerData : public ComponentTransportLocalAssemblerInterface
{
// eval shape matrices at given point
auto const shape_matrices = [&]() {
using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(
&_element));
auto const fe = NumLib::createIsoparametricFiniteElement<
ShapeFunction, ShapeMatricesType>(_element);

typename ShapeMatricesType::ShapeMatrices shape_matrices(
ShapeFunction::DIM, GlobalDim, ShapeFunction::NPOINTS);
Expand Down
7 changes: 2 additions & 5 deletions ProcessLib/GroundwaterFlow/GroundwaterFlowFEM.h
Original file line number Diff line number Diff line change
Expand Up @@ -138,11 +138,8 @@ class LocalAssemblerData : public GroundwaterFlowLocalAssemblerInterface
std::vector<double> const& local_x) const override
{
// eval dNdx and invJ at p
using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(
&_element));
auto const fe = NumLib::createIsoparametricFiniteElement<
ShapeFunction, ShapeMatricesType>(_element);

typename ShapeMatricesType::ShapeMatrices shape_matrices(
ShapeFunction::DIM, GlobalDim, ShapeFunction::NPOINTS);
Expand Down
7 changes: 2 additions & 5 deletions ProcessLib/HT/HTFEM.h
Original file line number Diff line number Diff line change
Expand Up @@ -132,11 +132,8 @@ class HTFEM : public HTLocalAssemblerInterface
std::vector<double> const& local_x) const override
{
// eval dNdx and invJ at given point
using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(
&this->_element));
auto const fe = NumLib::createIsoparametricFiniteElement<
ShapeFunction, ShapeMatricesType>(_element);

typename ShapeMatricesType::ShapeMatrices shape_matrices(
ShapeFunction::DIM, GlobalDim, ShapeFunction::NPOINTS);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -92,10 +92,8 @@ class PythonSourceTermLocalAssembler final
{
using ShapeMatricesType =
ShapeMatrixPolicyType<ShapeFunction, GlobalDim>;
using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;
FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(
&_element));
auto const fe = NumLib::createIsoparametricFiniteElement<
ShapeFunction, ShapeMatricesType>(_element);

unsigned const num_integration_points =
_integration_method.getNumberOfPoints();
Expand Down
7 changes: 2 additions & 5 deletions ProcessLib/SurfaceFlux/SurfaceFluxLocalAssembler.h
Original file line number Diff line number Diff line change
Expand Up @@ -71,11 +71,8 @@ class SurfaceFluxLocalAssembler final
_bulk_element_id(bulk_element_ids[surface_element.getID()]),
_bulk_face_id(bulk_face_ids[surface_element.getID()])
{
using FemType =
NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(
&surface_element));
auto const fe = NumLib::createIsoparametricFiniteElement<
ShapeFunction, ShapeMatricesType>(_surface_element);

std::size_t const n_integration_points =
_integration_method.getNumberOfPoints();
Expand Down
21 changes: 9 additions & 12 deletions ProcessLib/Utils/InitShapeMatrices.h
Original file line number Diff line number Diff line change
Expand Up @@ -28,10 +28,9 @@ initShapeMatrices(MeshLib::Element const& e, bool is_axially_symmetric,
Eigen::aligned_allocator<typename ShapeMatricesType::ShapeMatrices>>
shape_matrices;

using FemType = NumLib::TemplateIsoparametric<
ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(&e));
auto const fe =
NumLib::createIsoparametricFiniteElement<ShapeFunction,
ShapeMatricesType>(e);

unsigned const n_integration_points = integration_method.getNumberOfPoints();

Expand All @@ -52,10 +51,9 @@ double interpolateXCoordinate(
MeshLib::Element const& e,
typename ShapeMatricesType::ShapeMatrices::ShapeType const& N)
{
using FemType = NumLib::TemplateIsoparametric<
ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(&e));
auto const fe =
NumLib::createIsoparametricFiniteElement<ShapeFunction,
ShapeMatricesType>(e);

return fe.interpolateZerothCoordinate(N);
}
Expand All @@ -65,10 +63,9 @@ std::array<double, 3> interpolateCoordinates(
MeshLib::Element const& e,
typename ShapeMatricesType::ShapeMatrices::ShapeType const& N)
{
using FemType = NumLib::TemplateIsoparametric<
ShapeFunction, ShapeMatricesType>;

FemType fe(*static_cast<const typename ShapeFunction::MeshElement*>(&e));
auto const fe =
NumLib::createIsoparametricFiniteElement<ShapeFunction,
ShapeMatricesType>(e);

return fe.interpolateCoordinates(N);
}
Expand Down
6 changes: 3 additions & 3 deletions Tests/NumLib/TestFunctionInterpolation.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -47,8 +47,6 @@ TEST(NumLibFunctionInterpolationTest, Linear1DElement)
using ShapeFunction = NumLib::ShapeLine2;
using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, 1u>;

using FemType = NumLib::TemplateIsoparametric<ShapeFunction, ShapeMatricesType>;

using IntegrationMethod = NumLib::GaussLegendreIntegrationPolicy<
ShapeFunction::MeshElement>::IntegrationMethod;

Expand All @@ -60,7 +58,9 @@ TEST(NumLibFunctionInterpolationTest, Linear1DElement)
MeshLib::Line(std::array<MeshLib::Node*, 2>{{&pt_a, &pt_b}});

// set up shape function
FemType finite_element(*static_cast<const ShapeFunction::MeshElement*>(&element));
auto const finite_element =
NumLib::createIsoparametricFiniteElement<ShapeFunction,
ShapeMatricesType>(element);

const unsigned integration_order = 2;
IntegrationMethod integration_method(integration_order);
Expand Down

0 comments on commit b63c384

Please sign in to comment.