5#ifndef IPPL_PRECONFEMPOISSONSOLVER_H
6#define IPPL_PRECONFEMPOISSONSOLVER_H
15 template <
typename Tlhs,
unsigned Dim,
unsigned numElemDOFs>
25 const size_t& i,
const size_t& j,
38 template <
typename FieldLHS,
typename FieldRHS = FieldLHS>
40 constexpr static unsigned Dim = FieldLHS::dim;
41 using Tlhs =
typename FieldLHS::value_type;
54 std::conditional_t<Dim == 1, ippl::EdgeElement<Tlhs>,
55 std::conditional_t<Dim == 2, ippl::QuadrilateralElement<Tlhs>,
76 static_assert(std::is_floating_point<Tlhs>::value,
"Not a floating point type");
117 const auto firstElementVertexPoints =
122 refElement_m.getInverseTransposeTransformationJacobian(firstElementVertexPoints);
126 const Tlhs absDetDPhi = Kokkos::abs(
127 refElement_m.getDeterminantOfTransformationJacobian(firstElementVertexPoints));
130 DPhiInvT, absDetDPhi);
134 FieldBC bcType = bcField[0]->getBCType();
136 const auto algoOperator = [poissonEquationEval, &bcField,
this](
rhs_type field) ->
lhs_type {
138 field.setFieldBC(bcField);
147 const auto algoOperatorL = [poissonEquationEval, &bcField,
this](
lhs_type field) ->
lhs_type {
149 field.setFieldBC(bcField);
158 const auto algoOperatorU = [poissonEquationEval, &bcField,
this](
lhs_type field) ->
lhs_type {
160 field.setFieldBC(bcField);
169 const auto algoOperatorUL = [poissonEquationEval, &bcField,
this](
lhs_type field) ->
lhs_type {
171 field.setFieldBC(bcField);
180 const auto algoOperatorInvD = [poissonEquationEval, &bcField,
this](
lhs_type field) ->
lhs_type {
182 field.setFieldBC(bcField);
191 const auto algoOperatorD = [poissonEquationEval, &bcField,
this](
lhs_type field) ->
lhs_type {
193 field.setFieldBC(bcField);
203 std::string preconditioner_type =
204 this->
params_m.template get<std::string>(
"preconditioner_type");
205 int level = this->
params_m.template get<int>(
"newton_level");
206 int degree = this->
params_m.template get<int>(
"chebyshev_degree");
207 int inner = this->
params_m.template get<int>(
"gauss_seidel_inner_iterations");
208 int outer = this->
params_m.template get<int>(
"gauss_seidel_outer_iterations");
209 double omega = this->
params_m.template get<double>(
"ssor_omega");
210 int richardson_iterations =
211 this->
params_m.template get<int>(
"richardson_iterations");
214 algoOperatorInvD, algoOperatorD, 0, 0, preconditioner_type,
215 level, degree, richardson_iterations, inner, outer, omega);
237 (this->
lhs_mp)->fillHalo();
241 int output = this->
params_m.template get<int>(
"output_type");
266 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_upperlower(FieldLHS &field, F &evalFunction) const
FieldLHS evaluateAx_inversediag(FieldLHS &field, F &evalFunction) const
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 ///////////////////////////////////////////.
FieldLHS evaluateAx_diag(FieldLHS &field, F &evalFunction) const
FieldLHS evaluateAx_upper(FieldLHS &field, F &evalFunction) const
FieldLHS evaluateAx_lower(FieldLHS &field, F &evalFunction) const
void assignGhostToPhysical(Field &field)
virtual T getResidue() const
virtual int getIterationCount()
virtual void setOperator(OperatorF op)
void setPreconditioner(OperatorF &&op, LowerF &&lower, UpperF &&upper, UpperLowerF &&upper_and_lower, InverseDiagF &&inverse_diagonal, DiagF &&diagonal, double alpha, double beta, std::string preconditioner_type="", int level=5, int degree=31, int richardson_iterations=4, int inner=2, int outer=2, double omega=1.57079632679) override
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.
virtual void setRhs(rhs_type &rhs)
A solver for the poisson equation using finite element methods and Conjugate Gradient (CG)
void solve() override
Solve the poisson equation using finite element methods.
int getIterationCount()
Query how many iterations were required to obtain the solution the last time this solver was used.
PreconditionedFEMPoissonSolver(lhs_type &lhs, rhs_type &rhs)
std::conditional_t< Dim==1, ippl::EdgeElement< Tlhs >, std::conditional_t< Dim==2, ippl::QuadrilateralElement< Tlhs >, ippl::HexahedralElement< Tlhs > > > ElementType
Tlhs getL2Error(const F &analytic)
Query the L2-norm error compared to a given (analytical) sol.
LagrangeType lagrangeSpace_m
static constexpr unsigned Dim
virtual void setDefaultParameters() override
typename FieldLHS::value_type Tlhs
PreconditionedFEMPoissonSolver()
void setRhs(rhs_type &rhs) override
typename FieldRHS::Mesh_t MeshType
PCGSolverAlgorithm_t pcg_algo_m
QuadratureType quadrature_m
Tlhs getResidue() const
Query the residue.
Tlhs getAvg(bool Vol=false)
Query the average of the solution.
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)