1 #ifndef MAXWELL_HPP_INCLUDED
2 #define MAXWELL_HPP_INCLUDED
4 #include <Eigen/Sparse>
11 #include <unsupported/Eigen/src/IterativeSolvers/Scaling.h>
36 typedef Eigen::SparseLU<SystemMatrixType,Eigen::COLAMDOrdering<int> >
SolverType;
38 static constexpr std::string_view
SolverName =
"SparseLU";
43 bool use_threads =
true ,
44 std::array<int,2>
const* dqr =
nullptr ,
45 std::ostream &output = std::cout );
72 void assembleRHS(Eigen::Ref<const Eigen::VectorXd>
const & rho ,
73 Eigen::Ref<const Eigen::VectorXd>
const & J ,
74 Eigen::Ref<const Eigen::VectorXd>
const & EOld ,
75 Eigen::Ref<const Eigen::VectorXd>
const & BOld );
82 Eigen::VectorXd
solve();
84 Eigen::VectorXd
solve(
const Eigen::VectorXd &rhs);
90 template<
typename Derived>
91 double norm(Eigen::MatrixBase<Derived>
const &E ,
97 template<
typename Derived>
98 double normd(Eigen::MatrixBase<Derived>
const &E ,
114 void assembleLocalContribution(Eigen::Ref<const Eigen::MatrixXd>
const & A,
size_t iT,
size_t kL,
size_t kR, std::forward_list<Eigen::Triplet<double>> * triplets)
const;
115 void assembleLocalContributionH(Eigen::Ref<const Eigen::MatrixXd>
const & A,
size_t iT,
size_t k, std::forward_list<Eigen::Triplet<double>> * triplets)
const;
116 void assembleLocalContribution(Eigen::Ref<const Eigen::VectorXd>
const & R,
size_t iT,
size_t k);
118 void assembleSystem();
123 std::ostream & _output;
125 Eigen::VectorXd _rhs;
127 Eigen::VectorXd _interpOne;
128 std::vector<Eigen::MatrixXd> _ALoc00;
129 std::vector<Eigen::MatrixXd> _ALoc11;
130 std::vector<Eigen::MatrixXd> _ALoc22;
131 Eigen::VectorXd _scalingL, _scalingR;
137 :
u(maxwell.dimensionSystem()),
sizeG(maxwell.dimensionG()),
sizeE(maxwell.dimensionE()),
sizeB(maxwell.dimensionB()),
150 Eigen::VectorBlock<Eigen::VectorXd>
h()
155 Eigen::VectorBlock<Eigen::VectorXd>
G()
160 Eigen::VectorBlock<Eigen::VectorXd>
E()
165 Eigen::VectorBlock<Eigen::VectorXd>
B()
176 : _ddrcore(mesh,degree,use_threads,dqr, output), _dt(dt), _use_threads(use_threads), _output(output), _rhs(dimensionSystem())
178 auto oneFunc = [](size_t ,
const Eigen::Vector<double,2> &)->Eigen::Vector<double,1> {
179 return Eigen::Vector<double,1>{1.};
181 _interpOne = _ddrcore.template interpolate<0>(oneFunc,dqr);
186 void MaxwellProblem::assembleLocalContribution(Eigen::Ref<const Eigen::MatrixXd>
const & A,
size_t iT,
size_t kL,
size_t kR,
187 std::forward_list<Eigen::Triplet<double>> * triplets)
const
207 for (
size_t iDL = 0; iDL <= 2; ++iDL) {
208 const size_t nbLDofsL = _ddrcore.
dofspace(kL).numLocalDofs(iDL);
209 if (nbLDofsL == 0)
continue;
210 std::vector<size_t>
const & boundariesL = _ddrcore.
mesh()->get_boundary(iDL,2,iT);
211 for (
size_t iDR = 0; iDR <= 2; ++iDR) {
212 const size_t nbLDofsR = _ddrcore.
dofspace(kR).numLocalDofs(iDR);
213 if (nbLDofsR == 0)
continue;
214 std::vector<size_t>
const & boundariesR = _ddrcore.
mesh()->get_boundary(iDR,2,iT);
215 for (
size_t iFL = 0; iFL < boundariesL.size(); ++iFL) {
216 const size_t offsetL = gOffsetL + _ddrcore.
dofspace(kL).globalOffset(iDL,boundariesL[iFL]);
217 const size_t lOffsetL = _ddrcore.
dofspace(kL).localOffset(iDL,2,iFL,iT);
218 for (
size_t iFR = 0; iFR < boundariesR.size(); ++iFR) {
219 const size_t offsetR = gOffsetR + _ddrcore.
dofspace(kR).globalOffset(iDR,boundariesR[iFR]);
220 const size_t lOffsetR = _ddrcore.
dofspace(kR).localOffset(iDR,2,iFR,iT);
221 for (
size_t iRow = 0; iRow < nbLDofsL; ++iRow) {
222 for (
size_t iCol = 0; iCol < nbLDofsR; ++iCol) {
223 triplets->emplace_front(offsetL + iRow,offsetR + iCol,A(lOffsetL + iRow,lOffsetR + iCol));
232 void MaxwellProblem::assembleLocalContributionH(Eigen::Ref<const Eigen::MatrixXd>
const & A,
size_t iT,
size_t k, std::forward_list<Eigen::Triplet<double>> * triplets)
const
234 assert(A.cols() == 1 &&
"Expected only 1 harmonic form");
245 for (
size_t iD = 0; iD <= 2; ++iD) {
246 const size_t nbLDofs = _ddrcore.
dofspace(k).numLocalDofs(iD);
247 if (nbLDofs == 0)
continue;
248 std::vector<size_t>
const & boundaries = _ddrcore.
mesh()->get_boundary(iD,2,iT);
249 for (
size_t iF = 0; iF < boundaries.size(); ++iF) {
250 const size_t offset = gOffset + _ddrcore.
dofspace(k).globalOffset(iD,boundaries[iF]);
251 const size_t lOffset = _ddrcore.
dofspace(k).localOffset(iD,2,iF,iT);
252 for (
size_t iRow = 0; iRow < nbLDofs; ++iRow) {
253 triplets->emplace_front(offset+iRow,0,A(lOffset+iRow,0));
254 triplets->emplace_front(0,offset+iRow,A(lOffset+iRow,0));
260 void MaxwellProblem::assembleLocalContribution(Eigen::Ref<const Eigen::VectorXd>
const & R,
size_t iT,
size_t k)
272 for (
size_t iD = 0; iD <= 2; ++iD) {
273 const size_t nbLDofs = _ddrcore.
dofspace(k).numLocalDofs(iD);
274 if (nbLDofs == 0)
continue;
275 std::vector<size_t>
const & boundaries = _ddrcore.
mesh()->get_boundary(iD,2,iT);
276 for (
size_t iF = 0; iF < boundaries.size(); ++iF) {
277 const size_t offset = gOffset + _ddrcore.
dofspace(k).globalOffset(iD,boundaries[iF]);
278 const size_t lOffset = _ddrcore.
dofspace(k).localOffset(iD,2,iF,iT);
279 for (
size_t iRow = 0; iRow < nbLDofs; ++iRow) {
280 std::atomic_ref<double> atomicRHS(_rhs[offset+iRow]);
281 atomicRHS.fetch_add(R(lOffset+iRow),std::memory_order_relaxed);
287 void MaxwellProblem::assembleSystem()
289 const size_t nb_cell = _ddrcore.
mesh()->n_cells(2);
291 std::vector<Eigen::MatrixXd> ALoc0h(nb_cell);
292 std::vector<Eigen::MatrixXd> ALoc01(nb_cell);
294 std::vector<Eigen::MatrixXd> ALoc12(nb_cell);
295 _ALoc00.resize(nb_cell);
296 _ALoc11.resize(nb_cell);
297 _ALoc22.resize(nb_cell);
298 std::function<void(
size_t start,
size_t end)> assemble_local = [&](
size_t start,
size_t end)->
void {
299 for (
size_t iT = start; iT < end; iT++) {
304 Eigen::MatrixXd loc0h = loc00*_ddrcore.
dofspace(0).restrict(2,iT,_interpOne);
305 Eigen::MatrixXd loc01 = _ddrcore.
compose_diff(0,2,iT).transpose()*loc11;
307 Eigen::MatrixXd loc12 = _ddrcore.
compose_diff(1,2,iT).transpose()*loc22;
318 _output<<
"[MaxwellProblem] Assembling local contributions..."<<std::flush;
320 _output<<
"\r[MaxwellProblem] Assembled local contributions "<<std::endl;
322 auto batch_local_assembly = [&](
size_t start,
size_t end, std::forward_list<Eigen::Triplet<double>> * triplets)->
void {
323 for (
size_t iT = start; iT < end; iT++) {
325 assembleLocalContributionH(ALoc0h[iT],iT,0,triplets);
326 assembleLocalContribution(ALoc01[iT],iT,0,1,triplets);
327 assembleLocalContribution(ALoc01[iT].transpose(),iT,1,0,triplets);
329 assembleLocalContribution(0.5*_dt*ALoc12[iT],iT,1,2,triplets);
330 assembleLocalContribution(0.5*_dt*ALoc12[iT].transpose(),iT,2,1,triplets);
331 assembleLocalContribution(-1.*_ALoc11[iT],iT,1,1,triplets);
332 assembleLocalContribution(_ALoc22[iT],iT,2,2,triplets);
335 _output<<
"[MaxwellProblem] Assembling global system from local contributions..."<<std::flush;
337 batch_local_assembly,_use_threads);
338 _output<<
"\r[MaxwellProblem] Assembled global system from local contributions "<<std::endl;
341 void MaxwellProblem::assembleRHS(Eigen::Ref<const Eigen::VectorXd>
const & rho, Eigen::Ref<const Eigen::VectorXd>
const & J,Eigen::Ref<const Eigen::VectorXd>
const & EOld, Eigen::Ref<const Eigen::VectorXd>
const & BOld)
343 auto batch_local_assembly = [&](
size_t start,
size_t end)->
void {
344 for (
size_t iT = start; iT < end; ++iT) {
348 assert(
static_cast<size_t>(rho.size()) ==
dimensionG() &&
"Incorrect size of rho");
349 loc = -_ALoc00[iT]*_ddrcore.
dofspace(0).restrict(2,iT,rho);
350 assembleLocalContribution(loc,iT,0);
353 assert(
static_cast<size_t>(J.size()) ==
dimensionE() &&
"Incorrect size of J");
354 assert(
static_cast<size_t>(EOld.size()) ==
dimensionE() &&
"Incorrect size of EOld");
355 loc = _ALoc11[iT]*_ddrcore.
dofspace(1).restrict(2,iT,_dt*J - EOld);
356 loc -= 0.5*_dt*_ddrcore.
compose_diff(1,2,iT).transpose()*_ALoc22[iT]*_ddrcore.
dofspace(2).restrict(2,iT,BOld);
357 assembleLocalContribution(loc,iT,1);
359 assert(
static_cast<size_t>(BOld.size()) ==
dimensionB() &&
"Incorrect size of BOld");
360 loc = _ALoc22[iT]*_ddrcore.
dofspace(2).restrict(2,iT,BOld);
361 loc -= 0.5*_dt*_ALoc22[iT]*_ddrcore.
compose_diff(1,2,iT)*_ddrcore.
dofspace(1).restrict(2,iT,EOld);
362 assembleLocalContribution(loc,iT,2);
365 _output<<
"[MaxwellProblem] Assembling RHS..."<<std::flush;
367 parallel_for(_ddrcore.
mesh()->n_cells(2),batch_local_assembly,_use_threads);
368 _output<<
"\r[MaxwellProblem] Assembled RHS "<<std::endl;
373 double const err = (_system*u - _rhs).
norm()/_rhs.norm();
375 std::cerr<<
"[MaxwellProblem] The relative error on the solution is "<<err<<
". This could indicate that the solver silently failed"<<std::endl;
381 _output <<
"[MaxwellProblem] Setting solver "<<
SolverName<<
" with "<<
dimensionSystem()<<
" degrees of freedom"<<std::endl;
383 Eigen::IterScaling<SystemMatrixType> iterScaling;
384 iterScaling.computeRef(systemScaled);
385 _scalingL = iterScaling.LeftScaling();
386 _scalingR = iterScaling.RightScaling();
387 _solver.compute(systemScaled);
388 if (_solver.info() != Eigen::Success) {
389 std::cerr <<
"[MaxwellProblem] Failed to factorize the system" << std::endl;
390 throw std::runtime_error(
"Factorization failed");
394 Eigen::VectorXd scaledRhs = _rhs.cwiseProduct(_scalingR);
395 Eigen::VectorXd u = _solver.solve(scaledRhs);
396 u = u.cwiseProduct(_scalingL);
397 if (_solver.info() != Eigen::Success) {
398 std::cerr <<
"[MaxwellProblem] Failed to solve the system" << std::endl;
399 throw std::runtime_error(
"Solver failed");
406 Eigen::VectorXd scaledRhs = rhs.cwiseProduct(_scalingR);
407 Eigen::VectorXd u = _solver.solve(scaledRhs);
408 u = u.cwiseProduct(_scalingL);
409 if (_solver.info() != Eigen::Success) {
410 std::cerr <<
"[MaxwellProblem] Failed to solve the system" << std::endl;
411 throw std::runtime_error(
"Solver failed");
417 template<
typename Derived>
420 assert(
static_cast<size_t>(E.size()) == _ddrcore.
dofspace(k).dimensionMesh() &&
"Mismatched dimension, wrong form degree?");
421 const size_t nb_cell = _ddrcore.
mesh()->n_cells(2);
422 std::vector<double> accErr(nb_cell);
423 std::function<void(
size_t start,
size_t end)> compute_local = [&](
size_t start,
size_t end)->
void {
424 for (
size_t iT = start; iT < end; ++iT){
425 Eigen::VectorXd locE = _ddrcore.
dofspace(k).restrict(2,iT,E);
426 Eigen::MatrixXd matL2;
440 accErr[iT] = locE.dot(matL2*locE);
445 for (
double err : accErr) {
448 return std::sqrt(acc);
451 template<
typename Derived>
454 if (k == 2)
return 0.;
455 assert(
static_cast<size_t>(E.size()) == _ddrcore.
dofspace(k).dimensionMesh() &&
"Mismatched dimension, wrong form degree?");
456 const size_t nb_cell = _ddrcore.
mesh()->n_cells(2);
457 std::vector<double> accErr(nb_cell);
458 std::function<void(
size_t start,
size_t end)> compute_local = [&](
size_t start,
size_t end)->
void {
459 for (
size_t iT = start; iT < end; ++iT){
461 Eigen::MatrixXd matL2;
472 accErr[iT] = locE.dot(matL2*locE);
477 for (
double err : accErr) {
480 return std::sqrt(acc);
486 std::forward_list<Eigen::Triplet<double>> triplets;
487 for (
size_t iD = 0; iD < 3; ++iD) {
488 const size_t nDofsL = _ddrcore.
dofspace(k+1).numLocalDofs(iD);
489 if (nDofsL == 0)
continue;
490 for (
size_t iT = 0; iT < _ddrcore.
mesh()->n_cells(iD); ++iT) {
492 const size_t offsetL = _ddrcore.
dofspace(k+1).globalOffset(iD,iT);
493 const size_t lOffsetL = _ddrcore.
dofspace(k+1).localOffset(iD,iT);
494 for (
size_t iDR = 0; iDR<=iD;++iDR) {
495 const size_t nDofsR = _ddrcore.
dofspace(k).numLocalDofs(iDR);
496 if (nDofsR == 0)
continue;
497 std::vector<size_t>
const & boundariesR = _ddrcore.
mesh()->get_boundary(iDR,iD,iT);
498 for (
size_t iFR = 0; iFR< boundariesR.size(); ++iFR) {
499 const size_t offsetR = _ddrcore.
dofspace(k).globalOffset(iDR,boundariesR[iFR]);
500 const size_t lOffsetR = _ddrcore.
dofspace(k).localOffset(iDR,iD,iFR,iT);
501 for (
size_t iRow = 0; iRow < nDofsL; ++iRow) {
502 for (
size_t iCol = 0; iCol < nDofsR; ++iCol) {
503 triplets.emplace_front(offsetL+iRow,offsetR+iCol,A(lOffsetL+iRow,lOffsetR+iCol));
510 rv.setFromTriplets(triplets.begin(),triplets.end());
517 std::forward_list<Eigen::Triplet<double>> triplets;
518 std::vector<Eigen::MatrixXd>
const & A = (k == 0)? _ALoc00 : ((k == 1)? _ALoc11 : _ALoc22);
519 for (
size_t iT = 0; iT<_ddrcore.
mesh()->n_cells(2); ++iT) {
520 for (
size_t iDL = 0; iDL <= 2; ++iDL) {
521 const size_t nbLDofsL = _ddrcore.
dofspace(k).numLocalDofs(iDL);
522 if (nbLDofsL == 0)
continue;
523 std::vector<size_t>
const & boundariesL = _ddrcore.
mesh()->get_boundary(iDL,2,iT);
524 for (
size_t iDR = 0; iDR <= 2; ++iDR) {
525 const size_t nbLDofsR = _ddrcore.
dofspace(k).numLocalDofs(iDR);
526 if (nbLDofsR == 0)
continue;
527 std::vector<size_t>
const & boundariesR = _ddrcore.
mesh()->get_boundary(iDR,2,iT);
528 for (
size_t iFL = 0; iFL < boundariesL.size(); ++iFL) {
529 const size_t offsetL = _ddrcore.
dofspace(k).globalOffset(iDL,boundariesL[iFL]);
530 const size_t lOffsetL = _ddrcore.
dofspace(k).localOffset(iDL,2,iFL,iT);
531 for (
size_t iFR = 0; iFR < boundariesR.size(); ++iFR) {
532 const size_t offsetR = _ddrcore.
dofspace(k).globalOffset(iDR,boundariesR[iFR]);
533 const size_t lOffsetR = _ddrcore.
dofspace(k).localOffset(iDR,2,iFR,iT);
534 for (
size_t iRow = 0; iRow < nbLDofsL; ++iRow) {
535 for (
size_t iCol = 0; iCol < nbLDofsR; ++iCol) {
536 triplets.emplace_front(offsetL + iRow,offsetR + iCol,A[iT](lOffsetL + iRow,lOffsetR + iCol));
544 rv.setFromTriplets(triplets.begin(),triplets.end());
Eigen::MatrixXd computeL2Product(size_t k, size_t i_cell) const
Compute the local product (including the cell boundary)
Definition: ddr_spaces.cpp:272
const Mesh< dimension > * mesh() const
Return the Mesh associated with the class.
Definition: ddr_spaces.hpp:72
DOFSpace< dimension > const & dofspace(size_t k) const
Return DOFSpace associated to k forms.
Definition: ddr_spaces.hpp:68
Eigen::MatrixXd compose_diff(size_t k, size_t d, size_t i_cell) const
Compute the projected differential operator.
Definition: ddr_spaces.cpp:248
2 dimensional Maxwell equation on a manifold without boundary
Definition: maxwell.hpp:31
void compute()
Setup the solver.
Definition: maxwell.hpp:380
double normd(Eigen::MatrixBase< Derived > const &E, size_t k) const
Compute the discrete norm of the differential of a -form.
Definition: maxwell.hpp:452
static constexpr std::string_view SolverName
Name of the solver.
Definition: maxwell.hpp:38
SystemMatrixType globalDiff(size_t k) const
Compute the global discrete differential.
Definition: maxwell.hpp:483
size_t dimensionG() const
Return the number of unknowns attached to the space of 0 forms.
Definition: maxwell.hpp:57
Eigen::VectorXd solve()
Solve the system and return the solution.
Definition: maxwell.hpp:393
void assembleRHS(Eigen::Ref< const Eigen::VectorXd > const &rho, Eigen::Ref< const Eigen::VectorXd > const &J, Eigen::Ref< const Eigen::VectorXd > const &EOld, Eigen::Ref< const Eigen::VectorXd > const &BOld)
Assemble the Right-Hand side from the given interpolates.
Definition: maxwell.hpp:341
size_t dimensionH() const
Return the number of unknowns attached to the space of harmonic 0 forms.
Definition: maxwell.hpp:51
Eigen::SparseLU< SystemMatrixType, Eigen::COLAMDOrdering< int > > SolverType
Type used to store the solver.
Definition: maxwell.hpp:36
Eigen::SparseMatrix< double > SystemMatrixType
Type used to store the system.
Definition: maxwell.hpp:34
size_t dimensionE() const
Return the number of unknowns attached to the space of 1 forms.
Definition: maxwell.hpp:60
size_t dimensionB() const
Return the number of unknowns attached to the space of 2 forms.
Definition: maxwell.hpp:62
SystemMatrixType globalMass(size_t k) const
Compute the global mass matrix.
Definition: maxwell.hpp:514
size_t dimensionSystem() const
Return the total number of unknowns.
Definition: maxwell.hpp:64
double timeStep() const
Return the time-step.
Definition: maxwell.hpp:69
double norm(Eigen::MatrixBase< Derived > const &E, size_t k) const
Compute the discrete norm of a -form.
Definition: maxwell.hpp:418
DDR_Spaces< 2 > const & ddrcore() const
Return the associated DDR_Spaces.
Definition: maxwell.hpp:67
bool validateSolution(Eigen::Ref< const Eigen::VectorXd > const &u) const
Check if the vector u if a solution up to a given relative accuracy.
Definition: maxwell.hpp:371
MaxwellProblem(const Mesh< 2 > &mesh, int degree, double timestep, bool use_threads=true, std::array< int, 2 > const *dqr=nullptr, std::ostream &output=std::cout)
Definition: maxwell.hpp:175
Main data structure for the mesh.
Definition: mesh.hpp:54
Discrete operators for DDR-PEC.
Definition: maxwell.hpp:23
static void parallel_for(unsigned nb_elements, std::function< void(size_t start, size_t end)> functor, bool use_threads=true)
Generic function to execute threaded processes.
Definition: parallel_for.hpp:42
Eigen::SparseMatrix< double > parallel_assembly(size_t nb_elements, std::pair< size_t, size_t > systemSize, FType localAssembly, bool use_threads=true)
Function to assemble a global sparse matrix from a procedure that compute local contributions.
Definition: parallel_for.hpp:74
Helper to parallelize operations using pthreads.
Construct a convenient wrapper around a vector of global unknowns to manipulate individual components...
Definition: maxwell.hpp:135
int offsetG
Definition: maxwell.hpp:171
Eigen::VectorBlock< Eigen::VectorXd > E()
Return the electric part.
Definition: maxwell.hpp:160
MaxwellVector & operator=(const MaxwellVector &other)
Copy another MaxwellVector.
Definition: maxwell.hpp:140
int offsetE
Definition: maxwell.hpp:171
int sizeE
Definition: maxwell.hpp:170
Eigen::VectorXd u
Definition: maxwell.hpp:169
int sizeB
Definition: maxwell.hpp:170
int sizeG
Definition: maxwell.hpp:170
Eigen::VectorBlock< Eigen::VectorXd > B()
Return the magnetic part.
Definition: maxwell.hpp:165
Eigen::VectorBlock< Eigen::VectorXd > h()
Return the harmonic part.
Definition: maxwell.hpp:150
int offsetB
Definition: maxwell.hpp:171
MaxwellVector(MaxwellProblem const &maxwell)
Definition: maxwell.hpp:136
Eigen::VectorBlock< Eigen::VectorXd > G()
Return the ghost part.
Definition: maxwell.hpp:155