IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
FEMMaxwellDiffusionSolver.h
Go to the documentation of this file.
1// Class FEMMaxwellDifussionSolver
2// Solves the electric diffusion probelm given by curl(curl(E)) + E = f in
3// the domain and n x E = 0 on the boundary.
4
5#ifndef IPPL_FEM_MAXWELL_DIFFUSION_SOLVER_H
6#define IPPL_FEM_MAXWELL_DIFFUSION_SOLVER_H
7
8#include "LinearSolvers/PCG.h"
9#include "Maxwell.h"
10#include <iomanip>
11#include <iostream>
12#include <fstream>
13
14namespace ippl {
15
26 template <typename T, unsigned Dim, unsigned numElementDOFs>
27 struct EvalFunctor {
35
43
50
69 KOKKOS_FUNCTION auto operator()(size_t i, size_t j,
70 const ippl::Vector<ippl::Vector<T, Dim>, numElementDOFs>& curl_b_q_k,
71 const ippl::Vector<ippl::Vector<T, Dim>, numElementDOFs>& val_b_q_k) const {
72
73 T curlTerm = dot(DPhiInvT*curl_b_q_k[j], DPhiInvT*curl_b_q_k[i]).apply();
74 T massTerm = dot(val_b_q_k[j], val_b_q_k[i]).apply();
75 return (curlTerm + massTerm)*absDetDPhi;
76 }
77 };
78
87 template <typename FieldType>
88 class FEMMaxwellDiffusionSolver : public Maxwell<FieldType, FieldType> {
89 constexpr static unsigned Dim = FieldType::dim;
90
91 // we call value_type twice as in theory we expect the field to store
92 // vector data represented by an ippl::Vector
93 using T = typename FieldType::value_type::value_type;
94
96
97 public:
99 using MeshType = typename FieldType::Mesh_t;
100
101 // PCG (Preconditioned Conjugate Gradient) is the solver algorithm used
104
105 // FEM Space types
106 using ElementType = std::conditional_t<Dim == 2, ippl::QuadrilateralElement<T>,
108
110
112
113 // default constructor (compatibility with Alpine)
115 : Base()
116 , rhsVector_m(nullptr)
117 , refElement_m()
118 , quadrature_m(refElement_m, 0.0, 0.0)
119 , nedelecSpace_m(*(new MeshType(NDIndex<Dim>(Vector<unsigned, Dim>(0)), Vector<T, Dim>(0),
121 {}
122
123 FEMMaxwellDiffusionSolver(FieldType& lhs, FieldType& rhs, const FEMVector<point_t>& rhsVector)
124 : Base(lhs, lhs, rhs)
125 , rhsVector_m(nullptr)
126 , refElement_m()
127 , quadrature_m(refElement_m, 0.0, 0.0)
128 , nedelecSpace_m(rhs.get_mesh(), refElement_m, quadrature_m, rhs.getLayout()) {
129
130 static_assert(std::is_floating_point<T>::value, "Not a floating point type");
132
133 // Calcualte the rhs, using the Nedelec space
135 std::make_unique<FEMVector<T>>(nedelecSpace_m.evaluateLoadVector(rhsVector));
136
137 rhsVector_m->accumulateHalo();
138 rhsVector_m->fillHalo();
139
140 }
141
142 void setRhs(FieldType& rhs, const FEMVector<point_t>& rhsVector) {
143
144 Base::setRhs(rhs);
145
146 // Calcualte the rhs, using the Nedelec space
148 std::make_unique<FEMVector<T>>(nedelecSpace_m.evaluateLoadVector(rhsVector));
149
150 rhsVector_m->accumulateHalo();
151 rhsVector_m->fillHalo();
152 }
153
157 void solve() override {
158
159 const Vector<size_t, Dim> zeroNdIndex = Vector<size_t, Dim>(0);
160
161 // We can pass the zeroNdIndex here, since the transformation
162 // jacobian does not depend on translation
163 const auto firstElementVertexPoints =
164 nedelecSpace_m.getElementMeshVertexPoints(zeroNdIndex);
165
166 // Compute Inverse Transpose Transformation Jacobian ()
167 const Vector<T, Dim> DPhiInvT =
168 refElement_m.getInverseTransposeTransformationJacobian(firstElementVertexPoints);
169
170 // Compute absolute value of the determinant of the transformation
171 // jacobian (|det D Phi_K|)
172 const T absDetDPhi = Kokkos::abs(
173 refElement_m.getDeterminantOfTransformationJacobian(firstElementVertexPoints));
174
175 // Create the functor object which stores the function we have to
176 // solve for the lhs
178 DPhiInvT, absDetDPhi);
179
180 // The Ax operator
181 const auto algoOperator = [maxwellDiffusionEval, this](FEMVector<T> vector)
182 -> FEMVector<T> {
183
184 vector.fillHalo();
185
186 FEMVector<T> return_vector = nedelecSpace_m.evaluateAx(vector,maxwellDiffusionEval);
187
188 return_vector.accumulateHalo();
189
190 return return_vector;
191 };
192
193 // setup the CG solver
194 pcg_algo_m.setOperator(algoOperator);
195
196 // Create the coefficient vector for the solution
197 FEMVector<T> lhsVector = rhsVector_m->deepCopy();
198
199 // Solve the system using CG
200 try {
201 pcg_algo_m(lhsVector, *rhsVector_m, this->params_m);
202 } catch (IpplException& e) {
203 std::string msg = e.where() + ": " + e.what() + "\n";
204 Kokkos::abort(msg.c_str());
205 }
206
207 // store solution.
208 lhsVector_m = std::make_unique<FEMVector<T>>(lhsVector);
209
210 // set the boundary values to the correct values.
211 lhsVector.fillHalo();
212
213 }
214
220 int getIterationCount() { return pcg_algo_m.getIterationCount(); }
221
226 T getResidue() const { return pcg_algo_m.getResidue(); }
227
251 Kokkos::View<point_t*> reconstructToPoints(const Kokkos::View<point_t*>& positions) {
252 return this->nedelecSpace_m.reconstructToPoints(positions, *lhsVector_m);
253 }
254
255
256
264 template <typename F>
265 T getL2Error(const F& analytic) {
266 T error_norm = this->nedelecSpace_m.computeError(*lhsVector_m, analytic);
267 return error_norm;
268 }
269
270
271 protected:
272
277
282 virtual void setDefaultParameters() {
283 this->params_m.add("max_iterations", 10);
284 this->params_m.add("tolerance", (T)1e-13);
285 }
286
292 std::unique_ptr<FEMVector<T>> rhsVector_m;
293
299 std::unique_ptr<FEMVector<T>> lhsVector_m;
300
305
310
318 };
319
320} // namespace ippl
321
322
323
324#endif // IPPL_FEM_MAXWELL_DIFFUSION_SOLVER_H
Definition Archive.h:20
KOKKOS_INLINE_FUNCTION detail::meta_dot< E1, E2 > dot(const detail::Expression< E1, N1 > &u, const detail::Expression< E2, N2 > &v)
1D vector used in the context of FEM.
Definition FEMVector.h:32
void fillHalo()
Copy values from neighboring ranks into local halo.
Definition FEMVector.hpp:38
void accumulateHalo()
Accumulate halo values in neighbor.
Definition FEMVector.hpp:95
A class representing a Nedelec 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.
KOKKOS_FUNCTION auto operator()(size_t i, size_t j, const ippl::Vector< ippl::Vector< T, Dim >, numElementDOFs > &curl_b_q_k, const ippl::Vector< ippl::Vector< T, Dim >, numElementDOFs > &val_b_q_k) const
Returns the evaluation of (curl(b_i)*curl(b_j) + b_i*b_j)*absDetDPhi.
const T absDetDPhi
The determinant of the Jacobian.
const Vector< T, Dim > DPhiInvT
The inverse transpose Jacobian.
EvalFunctor(Vector< T, Dim > DPhiInvT, T absDetDPhi)
Constructor.
NedelecSpace< T, Dim, 1, ElementType, QuadratureType, FieldType > NedelecType
Maxwell< FieldType, FieldType > Base
std::unique_ptr< FEMVector< T > > rhsVector_m
FEM represenation of the rhs We use this to store the rhs b of the System Ax = b used in the Galerkin...
void setRhs(FieldType &rhs, const FEMVector< point_t > &rhsVector)
virtual void setDefaultParameters()
Sets the default values for the CG solver. Defaults are: max Iterations = 10, tolerance = 1e-13.
Kokkos::View< point_t * > reconstructToPoints(const Kokkos::View< point_t * > &positions)
Reconstructs function values at arbitrary points in the mesh.
void solve() override
Solve the equation using finite element methods.
std::conditional_t< Dim==2, ippl::QuadrilateralElement< T >, ippl::HexahedralElement< T > > ElementType
GaussJacobiQuadrature< T, 5, ElementType > QuadratureType
T getL2Error(const F &analytic)
Given an analytical solution computes the L2 norm error.
QuadratureType quadrature_m
The quadrature rule we use.
FEMMaxwellDiffusionSolver(FieldType &lhs, FieldType &rhs, const FEMVector< point_t > &rhsVector)
CG< FEMVector< T >, FEMVector< T >, FEMVector< T >, FEMVector< T >, FEMVector< T >, FEMVector< T >, FEMVector< T > > PCGSolverAlgorithm_t
typename FieldType::value_type::value_type T
PCGSolverAlgorithm_t pcg_algo_m
The CG Solver we use.
ElementType refElement_m
the reference element we have.
NedelecType nedelecSpace_m
The Nedelec Space object.
std::unique_ptr< FEMVector< T > > lhsVector_m
FEM represenation of the solution vector We use this to store the solution x of the System Ax = b use...
virtual const char * what() const
virtual const std::string & where() const