5#ifndef IPPL_FEMPOISSONSOLVER_H
6#define IPPL_FEMPOISSONSOLVER_H
13 template <
typename Tlhs,
unsigned Dim,
unsigned numElemDOFs>
23 const size_t& i,
const size_t& j,
36 template <
typename FieldLHS,
typename FieldRHS = FieldLHS,
unsigned Order = 1,
unsigned QuadNumNodes = 5>
38 constexpr static unsigned Dim = FieldLHS::dim;
39 using Tlhs =
typename FieldLHS::value_type;
52 std::conditional_t<Dim == 1, ippl::EdgeElement<Tlhs>,
53 std::conditional_t<Dim == 2, ippl::QuadrilateralElement<Tlhs>,
74 static_assert(std::is_floating_point<Tlhs>::value,
"Not a floating point type");
115 const auto firstElementVertexPoints =
120 refElement_m.getInverseTransposeTransformationJacobian(firstElementVertexPoints);
124 const Tlhs absDetDPhi = Kokkos::abs(
125 refElement_m.getDeterminantOfTransformationJacobian(firstElementVertexPoints));
128 DPhiInvT, absDetDPhi);
132 FieldBC bcType = bcField[0]->getBCType();
134 const auto algoOperator = [poissonEquationEval, &bcField,
this](
rhs_type field) ->
lhs_type {
136 field.setFieldBC(bcField);
140 auto return_field =
lagrangeSpace_m.evaluateAx(field, poissonEquationEval);
159 (this->
lhs_mp)->fillHalo();
188 template <
typename F>
215 this->
params_m.add(
"max_iterations", 1000);
detail::meta_grad< Field > grad(Field &u)
KOKKOS_INLINE_FUNCTION auto & get(Tuple< Ts... > &t)
Accessor function to get an element mutable reference at a specific index from a Tuple.
KOKKOS_INLINE_FUNCTION detail::meta_dot< E1, E2 > dot(const detail::Expression< E1, N1 > &u, const detail::Expression< E2, N2 > &v)
A class representing a Lagrange space for finite element methods on a structured, rectilinear grid.
This is class represents the Gauss-Jacobi quadrature rule on a reference element.
Representation of the lhs of the problem we are trying to solve.
const T absDetDPhi
The determinant of the Jacobian.
KOKKOS_FUNCTION auto operator()(const size_t &i, const size_t &j, const Vector< Vector< Tlhs, Dim >, numElemDOFs > &grad_b_q_k) const
EvalFunctor(Vector< Tlhs, Dim > DPhiInvT, Tlhs absDetDPhi)
const Vector< T, Dim > DPhiInvT
The inverse transpose Jacobian.
typename FieldRHS::Mesh_t MeshType
FEMPoissonSolver(lhs_type &lhs, rhs_type &rhs)
void setRhs(rhs_type &rhs) override
virtual void setDefaultParameters() override
Poisson< FieldLHS, FieldRHS > Base
GaussJacobiQuadrature< Tlhs, QuadNumNodes, ElementType > QuadratureType
static constexpr unsigned Dim
CG< lhs_type, lhs_type, lhs_type, lhs_type, lhs_type, FieldLHS, FieldRHS > PCGSolverAlgorithm_t
QuadratureType quadrature_m
LagrangeSpace< Tlhs, Dim, Order, ElementType, QuadratureType, FieldLHS, FieldRHS > LagrangeType
Tlhs getL2Error(const F &analytic)
std::conditional_t< Dim==1, ippl::EdgeElement< Tlhs >, std::conditional_t< Dim==2, ippl::QuadrilateralElement< Tlhs >, ippl::HexahedralElement< Tlhs > > > ElementType
LagrangeType lagrangeSpace_m
typename FieldLHS::value_type Tlhs
Tlhs getAvg(bool Vol=false)
void solve() override
Solve the poisson equation using finite element methods. The problem is described by -laplace(lhs) = ...
PCGSolverAlgorithm_t pcg_algo_m
virtual void setRhs(rhs_type &rhs)
Timing::TimerRef TimerRef
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)