5#ifndef IPPL_FEM_MAXWELL_DIFFUSION_SOLVER_H
6#define IPPL_FEM_MAXWELL_DIFFUSION_SOLVER_H
26 template <
typename T,
unsigned Dim,
unsigned numElementDOFs>
74 T massTerm =
dot(val_b_q_k[j], val_b_q_k[i]).apply();
87 template <
typename FieldType>
89 constexpr static unsigned Dim = FieldType::dim;
93 using T =
typename FieldType::value_type::value_type;
106 using ElementType = std::conditional_t<Dim == 2, ippl::QuadrilateralElement<T>,
124 :
Base(lhs, lhs, rhs)
130 static_assert(std::is_floating_point<T>::value,
"Not a floating point type");
135 std::make_unique<FEMVector<T>>(
nedelecSpace_m.evaluateLoadVector(rhsVector));
148 std::make_unique<FEMVector<T>>(
nedelecSpace_m.evaluateLoadVector(rhsVector));
163 const auto firstElementVertexPoints =
168 refElement_m.getInverseTransposeTransformationJacobian(firstElementVertexPoints);
172 const T absDetDPhi = Kokkos::abs(
173 refElement_m.getDeterminantOfTransformationJacobian(firstElementVertexPoints));
178 DPhiInvT, absDetDPhi);
181 const auto algoOperator = [maxwellDiffusionEval,
this](
FEMVector<T> vector)
190 return return_vector;
203 std::string msg = e.
where() +
": " + e.
what() +
"\n";
204 Kokkos::abort(msg.c_str());
208 lhsVector_m = std::make_unique<FEMVector<T>>(lhsVector);
264 template <
typename F>
283 this->
params_m.add(
"max_iterations", 10);
284 this->
params_m.add(
"tolerance", (
T)1e-13);
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.
void fillHalo()
Copy values from neighboring ranks into local halo.
void accumulateHalo()
Accumulate halo values in neighbor.
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.
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...
typename FieldType::Mesh_t MeshType
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
FEMMaxwellDiffusionSolver()
T getL2Error(const F &analytic)
Given an analytical solution computes the L2 norm error.
static constexpr unsigned Dim
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