OPALX (Object Oriented Parallel Accelerator Library for Exascale) MINIorX
OPALX
PreconditionedFEMPoissonSolver.h
Go to the documentation of this file.
1// Class FEMPoissonSolver
2// Solves the poisson equation using finite element methods and Conjugate
3// Gradient + a Preconditioner.
4
5#ifndef IPPL_PRECONFEMPOISSONSOLVER_H
6#define IPPL_PRECONFEMPOISSONSOLVER_H
7
8// #include "FEM/FiniteElementSpace.h"
9#include "LaplaceHelpers.h"
10#include "LinearSolvers/PCG.h"
11#include "Poisson.h"
12
13namespace ippl {
14
15 template <typename Tlhs, unsigned Dim, unsigned numElemDOFs>
16 struct EvalFunctor {
18 const Tlhs absDetDPhi;
19
23
24 KOKKOS_FUNCTION auto operator()(
25 const size_t& i, const size_t& j,
26 const Vector<Vector<Tlhs, Dim>, numElemDOFs>& grad_b_q_k) const {
27 return dot((DPhiInvT * grad_b_q_k[j]), (DPhiInvT * grad_b_q_k[i])).apply() * absDetDPhi;
28 }
29 };
30
38 template <typename FieldLHS, typename FieldRHS = FieldLHS>
39 class PreconditionedFEMPoissonSolver : public Poisson<FieldLHS, FieldRHS> {
40 constexpr static unsigned Dim = FieldLHS::dim;
41 using Tlhs = typename FieldLHS::value_type;
42
43 public:
45 using typename Base::lhs_type, typename Base::rhs_type;
46 using MeshType = typename FieldRHS::Mesh_t;
47
48 // PCG (Preconditioned Conjugate Gradient) is the solver algorithm used
51
52 // FEM Space types
54 std::conditional_t<Dim == 1, ippl::EdgeElement<Tlhs>,
55 std::conditional_t<Dim == 2, ippl::QuadrilateralElement<Tlhs>,
57
59
61
62 // default constructor (compatibility with Alpine)
64 : Base()
65 , refElement_m()
66 , quadrature_m(refElement_m, 0.0, 0.0)
67 , lagrangeSpace_m(*(new MeshType(NDIndex<Dim>(Vector<unsigned, Dim>(0)), Vector<Tlhs, Dim>(0),
69 {}
70
72 : Base(lhs, rhs)
73 , refElement_m()
74 , quadrature_m(refElement_m, 0.0, 0.0)
75 , lagrangeSpace_m(rhs.get_mesh(), refElement_m, quadrature_m, rhs.getLayout()) {
76 static_assert(std::is_floating_point<Tlhs>::value, "Not a floating point type");
78
79 // start a timer
80 static IpplTimings::TimerRef init = IpplTimings::getTimer("initFEM");
82
83 rhs.fillHalo();
84
86
87 rhs.fillHalo();
88
90 }
91
92 void setRhs(rhs_type& rhs) override {
93 Base::setRhs(rhs);
94
95 lagrangeSpace_m.initialize(rhs.get_mesh(), rhs.getLayout());
96
97 rhs.fillHalo();
98
100
101 rhs.fillHalo();
102 }
103
108 void solve() override {
109 // start a timer
112
113 const Vector<size_t, Dim> zeroNdIndex = Vector<size_t, Dim>(0);
114
115 // We can pass the zeroNdIndex here, since the transformation jacobian does not depend
116 // on translation
117 const auto firstElementVertexPoints =
119
120 // Compute Inverse Transpose Transformation Jacobian ()
121 const Vector<Tlhs, Dim> DPhiInvT =
122 refElement_m.getInverseTransposeTransformationJacobian(firstElementVertexPoints);
123
124 // Compute absolute value of the determinant of the transformation jacobian (|det D
125 // Phi_K|)
126 const Tlhs absDetDPhi = Kokkos::abs(
127 refElement_m.getDeterminantOfTransformationJacobian(firstElementVertexPoints));
128
130 DPhiInvT, absDetDPhi);
131
132 // get BC type of our RHS
133 BConds<FieldRHS, Dim>& bcField = (this->rhs_mp)->getFieldBC();
134 FieldBC bcType = bcField[0]->getBCType();
135
136 const auto algoOperator = [poissonEquationEval, &bcField, this](rhs_type field) -> lhs_type {
137 // set appropriate BCs for the field as the info gets lost in the CG iteration
138 field.setFieldBC(bcField);
139
140 field.fillHalo();
141
142 auto return_field = lagrangeSpace_m.evaluateAx(field, poissonEquationEval);
143
144 return return_field;
145 };
146
147 const auto algoOperatorL = [poissonEquationEval, &bcField, this](lhs_type field) -> lhs_type {
148 // set appropriate BCs for the field as the info gets lost in the CG iteration
149 field.setFieldBC(bcField);
150
151 field.fillHalo();
152
153 auto return_field = lagrangeSpace_m.evaluateAx_lower(field, poissonEquationEval);
154
155 return return_field;
156 };
157
158 const auto algoOperatorU = [poissonEquationEval, &bcField, this](lhs_type field) -> lhs_type {
159 // set appropriate BCs for the field as the info gets lost in the CG iteration
160 field.setFieldBC(bcField);
161
162 field.fillHalo();
163
164 auto return_field = lagrangeSpace_m.evaluateAx_upper(field, poissonEquationEval);
165
166 return return_field;
167 };
168
169 const auto algoOperatorUL = [poissonEquationEval, &bcField, this](lhs_type field) -> lhs_type {
170 // set appropriate BCs for the field as the info gets lost in the CG iteration
171 field.setFieldBC(bcField);
172
173 field.fillHalo();
174
175 auto return_field = lagrangeSpace_m.evaluateAx_upperlower(field, poissonEquationEval);
176
177 return return_field;
178 };
179
180 const auto algoOperatorInvD = [poissonEquationEval, &bcField, this](lhs_type field) -> lhs_type {
181 // set appropriate BCs for the field as the info gets lost in the CG iteration
182 field.setFieldBC(bcField);
183
184 field.fillHalo();
185
186 auto return_field = lagrangeSpace_m.evaluateAx_inversediag(field, poissonEquationEval);
187
188 return return_field;
189 };
190
191 const auto algoOperatorD = [poissonEquationEval, &bcField, this](lhs_type field) -> lhs_type {
192 // set appropriate BCs for the field as the info gets lost in the CG iteration
193 field.setFieldBC(bcField);
194
195 field.fillHalo();
196
197 auto return_field = lagrangeSpace_m.evaluateAx_diag(field, poissonEquationEval);
198
199 return return_field;
200 };
201
202 // set preconditioner for PCG
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");
212
213 pcg_algo_m.setPreconditioner(algoOperator, algoOperatorL, algoOperatorU, algoOperatorUL,
214 algoOperatorInvD, algoOperatorD, 0, 0, preconditioner_type,
215 level, degree, richardson_iterations, inner, outer, omega);
216
217 pcg_algo_m.setOperator(algoOperator);
218
219 // send boundary values to RHS (load vector) i.e. lifting (Dirichlet BCs)
220 if (bcType == CONSTANT_FACE) {
221 *(this->rhs_mp) = *(this->rhs_mp) -
222 lagrangeSpace_m.evaluateAx_lift(*(this->rhs_mp), poissonEquationEval);
223 }
224
225 // start a timer
226 static IpplTimings::TimerRef pcgTimer = IpplTimings::getTimer("pcg");
227 IpplTimings::startTimer(pcgTimer);
228
229 // run PCG -> lhs contains solution
230 pcg_algo_m(*(this->lhs_mp), *(this->rhs_mp), this->params_m);
231
232 // added for BCs to be imposed properly
233 // (they are not propagated through the preconditioner)
234 if (bcType == CONSTANT_FACE) {
235 bcField.assignGhostToPhysical(*(this->lhs_mp));
236 }
237 (this->lhs_mp)->fillHalo();
238
239 IpplTimings::stopTimer(pcgTimer);
240
241 int output = this->params_m.template get<int>("output_type");
242 if (output & Base::GRAD) {
243 *(this->grad_mp) = -grad(*(this->lhs_mp));
244 }
245
247 }
248
255
260 Tlhs getResidue() const { return pcg_algo_m.getResidue(); }
261
266 template <typename F>
267 Tlhs getL2Error(const F& analytic) {
268 Tlhs error_norm = this->lagrangeSpace_m.computeErrorL2(*(this->lhs_mp), analytic);
269 return error_norm;
270 }
271
277 Tlhs getAvg(bool Vol = false) {
278 Tlhs avg = this->lagrangeSpace_m.computeAvg(*(this->lhs_mp));
279 if (Vol) {
280 lhs_type unit((this->lhs_mp)->get_mesh(), (this->lhs_mp)->getLayout());
281 unit = 1.0;
282 Tlhs vol = this->lagrangeSpace_m.computeAvg(unit);
283 return avg/vol;
284 } else {
285 return avg;
286 }
287 }
288
289 protected:
291
292 virtual void setDefaultParameters() override {
293 this->params_m.add("max_iterations", 1000);
294 this->params_m.add("tolerance", (Tlhs)1e-13);
295 }
296
300 };
301
302} // namespace ippl
303
304#endif
ippl::UniformCartesian< double, 3 > Mesh_t
Definition: PBunchDefs.h:19
constexpr double e
The value of.
Definition: Physics.h:39
Implementations for FFT constructor/destructor and transforms.
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 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.
Definition: LagrangeSpace.h:35
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)
Definition: BConds.hpp:38
virtual T getResidue() const
Definition: PCG.h:151
virtual int getIterationCount()
Definition: PCG.h:73
virtual void setOperator(OperatorF op)
Definition: PCG.h:36
Definition: PCG.h:286
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
Definition: PCG.h:308
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.
FieldRHS rhs_type
Definition: Poisson.h:25
virtual void setRhs(rhs_type &rhs)
Definition: Poisson.h:96
FieldLHS lhs_type
Definition: Poisson.h:24
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.
Tlhs getAvg(bool Vol=false)
Query the average of the solution.
Timing::TimerRef TimerRef
Definition: IpplTimings.h:144
static TimerRef getTimer(const char *nm)
Definition: IpplTimings.h:150
static void stopTimer(TimerRef t)
Definition: IpplTimings.h:156
static void startTimer(TimerRef t)
Definition: IpplTimings.h:153
void add(const std::string &key, const T &value)
Definition: ParameterList.h:44