IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
ippl::FEMMaxwellDiffusionSolver< FieldType > Class Template Reference

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>

Inheritance diagram for ippl::FEMMaxwellDiffusionSolver< FieldType >:
Collaboration diagram for ippl::FEMMaxwellDiffusionSolver< FieldType >:

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 &params)

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, Dimpoint_t

Static Private Attributes

static constexpr unsigned Dim = FieldType::dim

Detailed Description

template<typename FieldType>
class ippl::FEMMaxwellDiffusionSolver< FieldType >

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.

Template Parameters
FieldTypeThe type used to represent a field on a mesh.

Definition at line 88 of file FEMMaxwellDiffusionSolver.h.

Member Typedef Documentation

◆ Base

template<typename FieldType>
using ippl::FEMMaxwellDiffusionSolver< FieldType >::Base = Maxwell<FieldType, FieldType>

Definition at line 98 of file FEMMaxwellDiffusionSolver.h.

◆ ElementType

template<typename FieldType>
using ippl::FEMMaxwellDiffusionSolver< FieldType >::ElementType
Initial value:

Definition at line 106 of file FEMMaxwellDiffusionSolver.h.

◆ MeshType

template<typename FieldType>
using ippl::FEMMaxwellDiffusionSolver< FieldType >::MeshType = typename FieldType::Mesh_t

Definition at line 99 of file FEMMaxwellDiffusionSolver.h.

◆ NedelecType

template<typename FieldType>
using ippl::FEMMaxwellDiffusionSolver< FieldType >::NedelecType = NedelecSpace<T, Dim, 1, ElementType, QuadratureType, FieldType>

Definition at line 111 of file FEMMaxwellDiffusionSolver.h.

◆ PCGSolverAlgorithm_t

template<typename FieldType>
using ippl::FEMMaxwellDiffusionSolver< FieldType >::PCGSolverAlgorithm_t
Initial value:

Definition at line 102 of file FEMMaxwellDiffusionSolver.h.

◆ point_t

template<typename FieldType>
typedef Vector<T,Dim> ippl::FEMMaxwellDiffusionSolver< FieldType >::point_t
private

Definition at line 95 of file FEMMaxwellDiffusionSolver.h.

◆ QuadratureType

template<typename FieldType>
using ippl::FEMMaxwellDiffusionSolver< FieldType >::QuadratureType = GaussJacobiQuadrature<T, 5, ElementType>

Definition at line 109 of file FEMMaxwellDiffusionSolver.h.

◆ T

template<typename FieldType>
using ippl::FEMMaxwellDiffusionSolver< FieldType >::T = typename FieldType::value_type::value_type
private

Definition at line 93 of file FEMMaxwellDiffusionSolver.h.

Constructor & Destructor Documentation

◆ FEMMaxwellDiffusionSolver() [1/2]

template<typename FieldType>
ippl::FEMMaxwellDiffusionSolver< FieldType >::FEMMaxwellDiffusionSolver ( )
inline

Definition at line 114 of file FEMMaxwellDiffusionSolver.h.

References Dim, nedelecSpace_m, quadrature_m, refElement_m, and rhsVector_m.

◆ FEMMaxwellDiffusionSolver() [2/2]

template<typename FieldType>
ippl::FEMMaxwellDiffusionSolver< FieldType >::FEMMaxwellDiffusionSolver ( FieldType & lhs,
FieldType & rhs,
const FEMVector< point_t > & rhsVector )
inline

Definition at line 123 of file FEMMaxwellDiffusionSolver.h.

References nedelecSpace_m, quadrature_m, refElement_m, rhsVector_m, and setDefaultParameters().

Here is the call graph for this function:

Member Function Documentation

◆ getIterationCount()

template<typename FieldType>
int ippl::FEMMaxwellDiffusionSolver< FieldType >::getIterationCount ( )
inline

Query how many iterations were required to obtain the solution the last time this solver was used

Returns
Iteration count of last solve

Definition at line 220 of file FEMMaxwellDiffusionSolver.h.

References pcg_algo_m.

◆ getL2Error()

template<typename FieldType>
template<typename F>
T ippl::FEMMaxwellDiffusionSolver< FieldType >::getL2Error ( const F & analytic)
inline

Given an analytical solution computes the L2 norm error.

 

Parameters
analyticalThe analytical solution (functor)
Returns
error - The error ||u - analytical||_L2

Definition at line 265 of file FEMMaxwellDiffusionSolver.h.

References lhsVector_m, and nedelecSpace_m.

◆ getResidue()

template<typename FieldType>
T ippl::FEMMaxwellDiffusionSolver< FieldType >::getResidue ( ) const
inline

Query the residue

Returns
Residue norm from last solve

Definition at line 226 of file FEMMaxwellDiffusionSolver.h.

References pcg_algo_m.

◆ mergeParameters()

void ippl::Maxwell< FieldType, FieldType >::mergeParameters ( const ParameterList & params)
inlineinherited

Merges another parameter set into the solver's parameters, overwriting existing parameters in case of conflict

Parameters
paramsParameter list with desired values

Definition at line 63 of file Maxwell.h.

References params_m.

◆ reconstructToPoints()

template<typename FieldType>
Kokkos::View< point_t * > ippl::FEMMaxwellDiffusionSolver< FieldType >::reconstructToPoints ( const Kokkos::View< point_t * > & positions)
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.

Note
Currently the function is able to handle both cases, where we have that 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.
Parameters
positionsThe points at which the function should be evaluated. A Kokkos::View which stores in each element a 2D/3D point.
Returns
The function evaluated at the given points, stored inside of 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.

◆ setDefaultParameters()

template<typename FieldType>
virtual void ippl::FEMMaxwellDiffusionSolver< FieldType >::setDefaultParameters ( )
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().

Here is the caller graph for this function:

◆ setEMFields()

void ippl::Maxwell< FieldType, FieldType >::setEMFields ( FieldType & E,
FieldType & B )
inlineinherited

Set the problem LHS (electromagnetic fields)

Parameters
EThe electric field
BThe magnetic field

Definition at line 52 of file Maxwell.h.

References Bn_mp, and En_mp.

◆ setRhs()

template<typename FieldType>
void ippl::FEMMaxwellDiffusionSolver< FieldType >::setRhs ( FieldType & rhs,
const FEMVector< point_t > & rhsVector )
inline

Definition at line 142 of file FEMMaxwellDiffusionSolver.h.

References nedelecSpace_m, and rhsVector_m.

◆ setSources()

virtual void ippl::Maxwell< FieldType, FieldType >::setSources ( FieldType & four_current)
inlinevirtualinherited

Set the problem RHS (charge & current densities)

Parameters
four_currentThe four current field (rho, J)

Definition at line 45 of file Maxwell.h.

References JN_mp.

◆ solve()

template<typename FieldType>
void ippl::FEMMaxwellDiffusionSolver< FieldType >::solve ( )
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().

Here is the call graph for this function:

Member Data Documentation

◆ Bn_mp

FieldType* ippl::Maxwell< FieldType, FieldType >::Bn_mp
protectedinherited

Definition at line 81 of file Maxwell.h.

Referenced by setEMFields().

◆ Dim

template<typename FieldType>
unsigned ippl::FEMMaxwellDiffusionSolver< FieldType >::Dim = FieldType::dim
staticconstexprprivate

Definition at line 89 of file FEMMaxwellDiffusionSolver.h.

Referenced by FEMMaxwellDiffusionSolver().

◆ En_mp

FieldType* ippl::Maxwell< FieldType, FieldType >::En_mp
protectedinherited

Definition at line 80 of file Maxwell.h.

Referenced by setEMFields().

◆ JN_mp

FieldType* ippl::Maxwell< FieldType, FieldType >::JN_mp
protectedinherited

Definition at line 77 of file Maxwell.h.

Referenced by setSources().

◆ lhsVector_m

template<typename FieldType>
std::unique_ptr<FEMVector<T> > ippl::FEMMaxwellDiffusionSolver< FieldType >::lhsVector_m
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().

◆ nedelecSpace_m

template<typename FieldType>
NedelecType ippl::FEMMaxwellDiffusionSolver< FieldType >::nedelecSpace_m
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().

◆ params_m

◆ pcg_algo_m

template<typename FieldType>
PCGSolverAlgorithm_t ippl::FEMMaxwellDiffusionSolver< FieldType >::pcg_algo_m
protected

The CG Solver we use.

Definition at line 276 of file FEMMaxwellDiffusionSolver.h.

Referenced by getIterationCount(), getResidue(), and solve().

◆ quadrature_m

template<typename FieldType>
QuadratureType ippl::FEMMaxwellDiffusionSolver< FieldType >::quadrature_m
protected

The quadrature rule we use.

Definition at line 309 of file FEMMaxwellDiffusionSolver.h.

Referenced by FEMMaxwellDiffusionSolver(), and FEMMaxwellDiffusionSolver().

◆ refElement_m

template<typename FieldType>
ElementType ippl::FEMMaxwellDiffusionSolver< FieldType >::refElement_m
protected

the reference element we have.

Definition at line 304 of file FEMMaxwellDiffusionSolver.h.

Referenced by FEMMaxwellDiffusionSolver(), FEMMaxwellDiffusionSolver(), and solve().

◆ rhsVector_m

template<typename FieldType>
std::unique_ptr<FEMVector<T> > ippl::FEMMaxwellDiffusionSolver< FieldType >::rhsVector_m
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().


The documentation for this class was generated from the following file: