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);
159 (this->
lhs_mp)->fillHalo();
163 int output = this->
params_m.template get<int>(
"output_type");
188 template <
typename F>
ippl::UniformCartesian< double, 3 > Mesh_t
constexpr double e
The value of.
Implementations for FFT constructor/destructor and transforms.
detail::meta_grad< Field > grad(Field &u)
KOKKOS_INLINE_FUNCTION detail::meta_dot< E1, E2 > dot(const detail::Expression< E1, N1 > &u, const detail::Expression< E2, N2 > &v)
KOKKOS_FUNCTION vertex_points_t getElementMeshVertexPoints(const indices_t &elementNDIndex) const
Get all the global vertex points of an element (given by its NDIndex).
A class representing a Lagrange space for finite element methods on a structured, rectilinear grid.
T computeAvg(const FieldLHS &u_h) const
Given a field, compute the average.
FieldLHS evaluateAx_lift(FieldLHS &field, F &evalFunction) const
Assemble the left stiffness matrix A of the system but only for the boundary values,...
void evaluateLoadVector(FieldRHS &field) const
Assemble the load vector b of the system Ax = b.
void initialize(UniformCartesian< T, Dim > &mesh, const Layout_t &layout)
Initialize a LagrangeSpace object created with the default constructor.
FieldLHS evaluateAx(FieldLHS &field, F &evalFunction) const
Assembly operations ///////////////////////////////////////////////.
T computeErrorL2(const FieldLHS &u_h, const F &u_sol) const
Error norm computations ///////////////////////////////////////////.
This is class represents the Gauss-Jacobi quadrature rule on a reference element.
virtual T getResidue() const
virtual int getIterationCount()
virtual void setOperator(OperatorF op)
Representation of the lhs of the problem we are trying to solve.
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.
const Vector< Tlhs, Dim > DPhiInvT
A solver for the poisson equation using finite element methods and Conjugate Gradient (CG)
typename FieldRHS::Mesh_t MeshType
FEMPoissonSolver(lhs_type &lhs, rhs_type &rhs)
void setRhs(rhs_type &rhs) override
virtual void setDefaultParameters() override
Tlhs getResidue() const
Query the residue.
static constexpr unsigned Dim
QuadratureType quadrature_m
int getIterationCount()
Query how many iterations were required to obtain the solution the last time this solver was used.
Tlhs getL2Error(const F &analytic)
Query the L2-norm error compared to a given (analytical) sol.
LagrangeType lagrangeSpace_m
typename FieldLHS::value_type Tlhs
std::conditional_t< Dim==1, ippl::EdgeElement< Tlhs >, std::conditional_t< Dim==2, ippl::QuadrilateralElement< Tlhs >, ippl::HexahedralElement< Tlhs > > > ElementType
Tlhs getAvg(bool Vol=false)
Query the average of the solution.
void solve() override
Solve the poisson equation using finite element methods.
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)
void add(const std::string &key, const T &value)