IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
FEMPoissonSolver.h
Go to the documentation of this file.
1// Class FEMPoissonSolver
2// Solves the poisson equation using finite element methods and Conjugate
3// Gradient
4
5#ifndef IPPL_FEMPOISSONSOLVER_H
6#define IPPL_FEMPOISSONSOLVER_H
7
8#include "LinearSolvers/PCG.h"
9#include "Poisson.h"
10
11namespace ippl {
12
13 template <typename Tlhs, unsigned Dim, unsigned numElemDOFs>
14 struct EvalFunctor {
16 const Tlhs absDetDPhi;
17
21
22 KOKKOS_FUNCTION auto operator()(
23 const size_t& i, const size_t& j,
24 const Vector<Vector<Tlhs, Dim>, numElemDOFs>& grad_b_q_k) const {
25 return dot((DPhiInvT * grad_b_q_k[j]), (DPhiInvT * grad_b_q_k[i])).apply() * absDetDPhi;
26 }
27 };
28
36 template <typename FieldLHS, typename FieldRHS = FieldLHS, unsigned Order = 1, unsigned QuadNumNodes = 5>
37 class FEMPoissonSolver : public Poisson<FieldLHS, FieldRHS> {
38 constexpr static unsigned Dim = FieldLHS::dim;
39 using Tlhs = typename FieldLHS::value_type;
40
41 public:
43 using typename Base::lhs_type, typename Base::rhs_type;
44 using MeshType = typename FieldRHS::Mesh_t;
45
46 // PCG (Preconditioned Conjugate Gradient) is the solver algorithm used
49
50 // FEM Space types
52 std::conditional_t<Dim == 1, ippl::EdgeElement<Tlhs>,
53 std::conditional_t<Dim == 2, ippl::QuadrilateralElement<Tlhs>,
55
57
59
60 // default constructor (compatibility with Alpine)
62 : Base()
63 , refElement_m()
64 , quadrature_m(refElement_m, 0.0, 0.0)
65 , lagrangeSpace_m(*(new MeshType(NDIndex<Dim>(Vector<unsigned, Dim>(0)), Vector<Tlhs, Dim>(0),
67 {}
68
70 : Base(lhs, rhs)
71 , refElement_m()
72 , quadrature_m(refElement_m, 0.0, 0.0)
73 , lagrangeSpace_m(rhs.get_mesh(), refElement_m, quadrature_m, rhs.getLayout()) {
74 static_assert(std::is_floating_point<Tlhs>::value, "Not a floating point type");
76
77 // start a timer
78 static IpplTimings::TimerRef init = IpplTimings::getTimer("initFEM");
80
81 rhs.fillHalo();
82
83 lagrangeSpace_m.evaluateLoadVector(rhs);
84
85 rhs.fillHalo();
86
88 }
89
90 void setRhs(rhs_type& rhs) override {
91 Base::setRhs(rhs);
92
93 lagrangeSpace_m.initialize(rhs.get_mesh(), rhs.getLayout());
94
95 rhs.fillHalo();
96
97 lagrangeSpace_m.evaluateLoadVector(rhs);
98
99 rhs.fillHalo();
100 }
101
106 void solve() override {
107 // start a timer
110
111 const Vector<size_t, Dim> zeroNdIndex = Vector<size_t, Dim>(0);
112
113 // We can pass the zeroNdIndex here, since the transformation jacobian does not depend
114 // on translation
115 const auto firstElementVertexPoints =
116 lagrangeSpace_m.getElementMeshVertexPoints(zeroNdIndex);
117
118 // Compute Inverse Transpose Transformation Jacobian ()
119 const Vector<Tlhs, Dim> DPhiInvT =
120 refElement_m.getInverseTransposeTransformationJacobian(firstElementVertexPoints);
121
122 // Compute absolute value of the determinant of the transformation jacobian (|det D
123 // Phi_K|)
124 const Tlhs absDetDPhi = Kokkos::abs(
125 refElement_m.getDeterminantOfTransformationJacobian(firstElementVertexPoints));
126
128 DPhiInvT, absDetDPhi);
129
130 // get BC type of our RHS
131 BConds<FieldRHS, Dim>& bcField = (this->rhs_mp)->getFieldBC();
132 FieldBC bcType = bcField[0]->getBCType();
133
134 const auto algoOperator = [poissonEquationEval, &bcField, this](rhs_type field) -> lhs_type {
135 // set appropriate BCs for the field as the info gets lost in the CG iteration
136 field.setFieldBC(bcField);
137
138 field.fillHalo();
139
140 auto return_field = lagrangeSpace_m.evaluateAx(field, poissonEquationEval);
141
142 return return_field;
143 };
144
145 pcg_algo_m.setOperator(algoOperator);
146
147 // send boundary values to RHS (load vector) i.e. lifting (Dirichlet BCs)
148 if (bcType == CONSTANT_FACE) {
149 *(this->rhs_mp) = *(this->rhs_mp) -
150 lagrangeSpace_m.evaluateAx_lift(*(this->rhs_mp), poissonEquationEval);
151 }
152
153 // start a timer
154 static IpplTimings::TimerRef pcgTimer = IpplTimings::getTimer("pcg");
155 IpplTimings::startTimer(pcgTimer);
156
157 pcg_algo_m(*(this->lhs_mp), *(this->rhs_mp), this->params_m);
158
159 (this->lhs_mp)->fillHalo();
160
161 IpplTimings::stopTimer(pcgTimer);
162
163 int output = this->params_m.template get<int>("output_type");
164 if (output & Base::GRAD) {
165 *(this->grad_mp) = -grad(*(this->lhs_mp));
166 }
167
169 }
170
176 int getIterationCount() { return pcg_algo_m.getIterationCount(); }
177
182 Tlhs getResidue() const { return pcg_algo_m.getResidue(); }
183
188 template <typename F>
189 Tlhs getL2Error(const F& analytic) {
190 Tlhs error_norm = this->lagrangeSpace_m.computeErrorL2(*(this->lhs_mp), analytic);
191 return error_norm;
192 }
193
199 Tlhs getAvg(bool Vol = false) {
200 Tlhs avg = this->lagrangeSpace_m.computeAvg(*(this->lhs_mp));
201 if (Vol) {
202 lhs_type unit((this->lhs_mp)->get_mesh(), (this->lhs_mp)->getLayout());
203 unit = 1.0;
204 Tlhs vol = this->lagrangeSpace_m.computeAvg(unit);
205 return avg/vol;
206 } else {
207 return avg;
208 }
209 }
210
211 protected:
213
214 virtual void setDefaultParameters() override {
215 this->params_m.add("max_iterations", 1000);
216 this->params_m.add("tolerance", (Tlhs)1e-13);
217 }
218
222 };
223
224} // namespace ippl
225
226#endif
Definition Archive.h:20
detail::meta_grad< Field > grad(Field &u)
FieldBC
Definition BcTypes.h:33
@ CONSTANT_FACE
Definition BcTypes.h:35
KOKKOS_INLINE_FUNCTION auto & get(Tuple< Ts... > &t)
Accessor function to get an element mutable reference at a specific index from a Tuple.
Definition Tuple.h:314
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.
Definition PCG.h:17
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
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
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
FieldRHS rhs_type
Definition Poisson.h:25
virtual void setRhs(rhs_type &rhs)
Definition Poisson.h:96
FieldLHS lhs_type
Definition Poisson.h:24
Timing::TimerRef TimerRef
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)