|
IPPL (Independent Parallel Particle Layer)
IPPL
|
A solver for the electric diffusion equation given by \( \nabla \times \nabla \times E + E = f \text{ in } \Omega\) and \( n \times E = 0 \text{ on } \partial \Omega\) using the Nédélec basis functions. More...
#include <FEMMaxwellDiffusionSolver.h>
Public Types | |
| using | Base = Maxwell<FieldType, FieldType> |
| using | MeshType = typename FieldType::Mesh_t |
| using | PCGSolverAlgorithm_t |
| using | ElementType |
| using | QuadratureType = GaussJacobiQuadrature<T, 5, ElementType> |
| using | NedelecType = NedelecSpace<T, Dim, 1, ElementType, QuadratureType, FieldType> |
Public Member Functions | |
| FEMMaxwellDiffusionSolver () | |
| FEMMaxwellDiffusionSolver (FieldType &lhs, FieldType &rhs, const FEMVector< point_t > &rhsVector) | |
| void | setRhs (FieldType &rhs, const FEMVector< point_t > &rhsVector) |
| void | solve () override |
| Solve the equation using finite element methods. | |
| int | getIterationCount () |
| T | getResidue () const |
| Kokkos::View< point_t * > | reconstructToPoints (const Kokkos::View< point_t * > &positions) |
| Reconstructs function values at arbitrary points in the mesh. | |
| template<typename F> | |
| T | getL2Error (const F &analytic) |
| Given an analytical solution computes the L2 norm error. | |
| virtual void | setSources (FieldType &four_current) |
| void | setEMFields (FieldType &E, FieldType &B) |
| void | mergeParameters (const ParameterList ¶ms) |
Protected Member Functions | |
| virtual void | setDefaultParameters () |
| Sets the default values for the CG solver. Defaults are: max Iterations = 10, tolerance = 1e-13. | |
Protected Attributes | |
| PCGSolverAlgorithm_t | pcg_algo_m |
| The CG Solver we use. | |
| 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 FEM scheme. | |
| 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 used in the Galerkin FEM scheme. | |
| ElementType | refElement_m |
| the reference element we have. | |
| QuadratureType | quadrature_m |
| The quadrature rule we use. | |
| NedelecType | nedelecSpace_m |
| The Nedelec Space object. | |
| ParameterList | params_m |
| FieldType * | JN_mp |
| FieldType * | En_mp |
| FieldType * | Bn_mp |
Private Types | |
| using | T = typename FieldType::value_type::value_type |
| typedef Vector< T, Dim > | point_t |
Static Private Attributes | |
| static constexpr unsigned | Dim = FieldType::dim |
A solver for the electric diffusion equation given by \( \nabla \times \nabla \times E + E = f \text{ in } \Omega\) and \( n \times E = 0 \text{ on } \partial \Omega\) using the Nédélec basis functions.
| FieldType | The type used to represent a field on a mesh. |
Definition at line 88 of file FEMMaxwellDiffusionSolver.h.
| using ippl::FEMMaxwellDiffusionSolver< FieldType >::Base = Maxwell<FieldType, FieldType> |
Definition at line 98 of file FEMMaxwellDiffusionSolver.h.
| using ippl::FEMMaxwellDiffusionSolver< FieldType >::ElementType |
Definition at line 106 of file FEMMaxwellDiffusionSolver.h.
| using ippl::FEMMaxwellDiffusionSolver< FieldType >::MeshType = typename FieldType::Mesh_t |
Definition at line 99 of file FEMMaxwellDiffusionSolver.h.
| using ippl::FEMMaxwellDiffusionSolver< FieldType >::NedelecType = NedelecSpace<T, Dim, 1, ElementType, QuadratureType, FieldType> |
Definition at line 111 of file FEMMaxwellDiffusionSolver.h.
| using ippl::FEMMaxwellDiffusionSolver< FieldType >::PCGSolverAlgorithm_t |
Definition at line 102 of file FEMMaxwellDiffusionSolver.h.
|
private |
Definition at line 95 of file FEMMaxwellDiffusionSolver.h.
| using ippl::FEMMaxwellDiffusionSolver< FieldType >::QuadratureType = GaussJacobiQuadrature<T, 5, ElementType> |
Definition at line 109 of file FEMMaxwellDiffusionSolver.h.
|
private |
Definition at line 93 of file FEMMaxwellDiffusionSolver.h.
|
inline |
Definition at line 114 of file FEMMaxwellDiffusionSolver.h.
References Dim, nedelecSpace_m, quadrature_m, refElement_m, and rhsVector_m.
|
inline |
Definition at line 123 of file FEMMaxwellDiffusionSolver.h.
References nedelecSpace_m, quadrature_m, refElement_m, rhsVector_m, and setDefaultParameters().
|
inline |
Query how many iterations were required to obtain the solution the last time this solver was used
Definition at line 220 of file FEMMaxwellDiffusionSolver.h.
References pcg_algo_m.
|
inline |
Given an analytical solution computes the L2 norm error.
| analytical | The analytical solution (functor) |
Definition at line 265 of file FEMMaxwellDiffusionSolver.h.
References lhsVector_m, and nedelecSpace_m.
|
inline |
Query the residue
Definition at line 226 of file FEMMaxwellDiffusionSolver.h.
References pcg_algo_m.
|
inlineinherited |
|
inline |
Reconstructs function values at arbitrary points in the mesh.
This function can be used to retrieve the values of a solution function at arbitrary points inside of the mesh.
positions only contains positions which are inside of local domain of this MPI rank (i.e. each rank gets its own unique position ) and where positions contains the positions of all ranks (i.e. positions is the same for all ranks). If in the future it can be guaranteed, that each rank will get its own positions then certain parts of the function implementation can be removed. Instructions for this are given in the implementation itself.| positions | The points at which the function should be evaluated. A Kokkos::View which stores in each element a 2D/3D point. |
Kokkos::View where each element corresponts to the function value at the point described by the same element inside of positions. Definition at line 251 of file FEMMaxwellDiffusionSolver.h.
References lhsVector_m, and nedelecSpace_m.
|
inlineprotectedvirtual |
Sets the default values for the CG solver. Defaults are: max Iterations = 10, tolerance = 1e-13.
Definition at line 282 of file FEMMaxwellDiffusionSolver.h.
References ippl::Maxwell< FieldType, FieldType >::params_m.
Referenced by FEMMaxwellDiffusionSolver().
|
inlineinherited |
|
inline |
Definition at line 142 of file FEMMaxwellDiffusionSolver.h.
References nedelecSpace_m, and rhsVector_m.
|
inlinevirtualinherited |
|
inlineoverridevirtual |
Solve the equation using finite element methods.
Implements ippl::Maxwell< FieldType, FieldType >.
Definition at line 157 of file FEMMaxwellDiffusionSolver.h.
References ippl::FEMVector< T >::accumulateHalo(), ippl::FEMVector< T >::fillHalo(), lhsVector_m, nedelecSpace_m, ippl::Maxwell< FieldType, FieldType >::params_m, pcg_algo_m, refElement_m, rhsVector_m, IpplException::what(), and IpplException::where().
|
protectedinherited |
Definition at line 81 of file Maxwell.h.
Referenced by setEMFields().
|
staticconstexprprivate |
Definition at line 89 of file FEMMaxwellDiffusionSolver.h.
Referenced by FEMMaxwellDiffusionSolver().
|
protectedinherited |
Definition at line 80 of file Maxwell.h.
Referenced by setEMFields().
|
protectedinherited |
Definition at line 77 of file Maxwell.h.
Referenced by setSources().
|
protected |
FEM represenation of the solution vector We use this to store the solution x of the System Ax = b used in the Galerkin FEM scheme.
Definition at line 299 of file FEMMaxwellDiffusionSolver.h.
Referenced by getL2Error(), reconstructToPoints(), and solve().
|
protected |
The Nedelec Space object.
This is the representation of the Nedelec space that we have and which we use to interact with all the Nedelec stuff.
Definition at line 317 of file FEMMaxwellDiffusionSolver.h.
Referenced by FEMMaxwellDiffusionSolver(), FEMMaxwellDiffusionSolver(), getL2Error(), reconstructToPoints(), setRhs(), and solve().
|
protectedinherited |
Definition at line 74 of file Maxwell.h.
Referenced by mergeParameters(), ippl::FEMMaxwellDiffusionSolver< FieldType >::setDefaultParameters(), and ippl::FEMMaxwellDiffusionSolver< FieldType >::solve().
|
protected |
The CG Solver we use.
Definition at line 276 of file FEMMaxwellDiffusionSolver.h.
Referenced by getIterationCount(), getResidue(), and solve().
|
protected |
The quadrature rule we use.
Definition at line 309 of file FEMMaxwellDiffusionSolver.h.
Referenced by FEMMaxwellDiffusionSolver(), and FEMMaxwellDiffusionSolver().
|
protected |
the reference element we have.
Definition at line 304 of file FEMMaxwellDiffusionSolver.h.
Referenced by FEMMaxwellDiffusionSolver(), FEMMaxwellDiffusionSolver(), and solve().
|
protected |
FEM represenation of the rhs We use this to store the rhs b of the System Ax = b used in the Galerkin FEM scheme.
Definition at line 292 of file FEMMaxwellDiffusionSolver.h.
Referenced by FEMMaxwellDiffusionSolver(), FEMMaxwellDiffusionSolver(), setRhs(), and solve().