HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
lasddr-yangmills.hpp
Go to the documentation of this file.
1 #ifndef SDDR_YANGMILLS_HPP
2 #define SDDR_YANGMILLS_HPP
3 
4 #include <iostream>
5 
6 #include <boost/math/constants/constants.hpp>
7 
8 #include <Eigen/Sparse>
9 #include <unsupported/Eigen/SparseExtra>
11 
12 #include <mesh.hpp>
13 #include <mesh_builder.hpp>
14 #include <GMpoly_cell.hpp>
15 #include <parallel_for.hpp>
16 
17 #include <laddrcore.hpp>
18 #include <lasxgrad.hpp>
19 #include <lasxcurl.hpp>
20 #include <lasxdiv.hpp>
21 
27 namespace HArDCore3D
28 {
29 
36  struct YangMillsNorms
37  {
39  YangMillsNorms(double norm_E, double norm_A, double norm_lambda):
40  E(norm_E),
41  A(norm_A),
42  lambda(norm_lambda)
43  {
44  // do nothing
45  };
46 
47  double E;
48  double A;
49  double lambda;
50  };
51 
53  struct YangMills
54  {
55  typedef Eigen::SparseMatrix<double> SystemMatrixType;
57  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType;
58  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> ElectricFieldType;
59  typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> MagneticFieldType;
61  typedef std::function<Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TForcingTermType;
62  typedef std::function<Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TElectricFieldType;
63  typedef std::function<Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TMagneticFieldType;
65  typedef std::vector<ForcingTermType> LAForcingTermType;
66  typedef std::vector<ElectricFieldType> LAElectricFieldType;
67  typedef std::vector<MagneticFieldType> LAMagneticFieldType;
69  typedef std::vector<TForcingTermType> TLAForcingTermType;
70  typedef std::vector<TElectricFieldType> TLAElectricFieldType;
71  typedef std::vector<TMagneticFieldType> TLAMagneticFieldType;
72 
74  YangMills(
75  const DDRCore & ddrcore,
76  const LADDRCore & laddrcore,
77  const LieAlgebra & liealgebra,
78  size_t nonlinear_discretisation,
79  bool use_threads,
80  std::ostream & output = std::cout
81  );
82 
85  double dt
86  );
89  const Eigen::VectorXd & interp_f,
90  const Eigen::VectorXd & interp_dE,
91  const Eigen::VectorXd & interp_A,
92  const Eigen::VectorXd & E_i,
93  const Eigen::VectorXd & A_i,
94  double dt,
95  double theta,
96  double nonlinear_coeff
97  );
100  const Eigen::VectorXd & E_i,
101  const Eigen::VectorXd & A_i,
102  const Eigen::VectorXd & Elambda_k,
103  double dt,
104  double theta,
105  double nonlinear_coeff
106  );
108  double stoppingCrit(
109  const Eigen::VectorXd & v,
110  const Eigen::VectorXd & u
111  );
113  void setNonlinearRes(
114  const Eigen::VectorXd & Elambda_k,
115  const Eigen::VectorXd & E_i,
116  const Eigen::VectorXd & A_i,
117  double dt,
118  double theta,
119  double nonlinear_coeff
120  );
121 
124  const LAMagneticFieldType & curl_A,
125  double dt
126  );
128  Eigen::VectorXd computeInitialConditions(const Eigen::MatrixXd & Eh,
129  const Eigen::MatrixXd & Ah,
130  const double nonlinear_coeff,
131  const size_t solver
132  );
133 
137  std::vector<Eigen::MatrixXd> epsBkt_v(size_t iT,
138  boost::multi_array<double, 3> & crossij_Pot_T,
139  const std::vector<Eigen::VectorXd> & vec_list
140  ) const;
141 
143  std::vector<Eigen::MatrixXd> L2v_Bkt(size_t iT,
144  boost::multi_array<double, 3> & intPciPcjPgk,
145  const std::vector<Eigen::VectorXd> & vec_list,
146  const size_t & entry
147  ) const;
148 
152  std::vector<Eigen::MatrixXd> L2v_epsBkt(size_t iT,
153  boost::multi_array<double, 3> & crossij_Pot_T,
154  const std::vector<Eigen::VectorXd> & v_L2prod
155  ) const;
156 
158  inline size_t dimensionSpace() const
159  {
160  return m_lasxcurl.dimension() + m_lasxgrad.dimension();
161  }
162 
164  inline size_t nbSCDOFs_E() const
165  {
166  return m_nloc_sc_E.sum();
167  }
168 
170  inline size_t nbSCDOFs_lambda() const
171  {
172  return m_nloc_sc_lambda.sum();
173  }
174 
176  inline size_t nbSCDOFs() const
177  {
178  return nbSCDOFs_E() + nbSCDOFs_lambda();
179  }
180 
182  inline size_t sizeSystem() const
183  {
184  return dimensionSpace() - nbSCDOFs();
185  }
186 
188  inline const LieAlgebra & lieAlg() const
189  {
190  return m_liealg;
191  }
192 
194  inline const LASXGrad & laSXGrad() const
195  {
196  return m_lasxgrad;
197  }
198 
200  inline const LASXCurl & laSXCurl() const
201  {
202  return m_lasxcurl;
203  }
204 
206  inline const LASXDiv & laSXDiv() const
207  {
208  return m_lasxdiv;
209  }
210 
212  inline const SXGrad & xSGrad() const
213  {
214  return m_sxgrad;
215  }
216 
218  inline const SXCurl & xSCurl() const
219  {
220  return m_sxcurl;
221  }
222 
224  inline const SXDiv & xSDiv() const
225  {
226  return m_sxdiv;
227  }
228 
230  inline const SystemMatrixType & systemMatrix() const {
231  return m_A;
232  }
233 
236  return m_A;
237  }
238 
240  inline const Eigen::VectorXd & systemVectorNewton() const {
241  return m_b;
242  }
243 
245  inline Eigen::VectorXd & systemVectorNewton() {
246  return m_b;
247  }
248 
250  inline const Eigen::VectorXd & systemVector() const {
251  return m_b_i;
252  }
253 
255  inline Eigen::VectorXd & systemVector() {
256  return m_b_i;
257  }
258 
260  inline const Eigen::VectorXd & nonlinearRes() const {
261  return m_b_k;
262  }
263 
265  inline Eigen::VectorXd & nonlinearRes() {
266  return m_b_k;
267  }
268 
270  inline const SystemMatrixType & scMatrix() const {
271  return m_sc_A;
272  }
273 
275  inline Eigen::VectorXd & scVector() {
276  return m_sc_b;
277  }
278 
280  inline const double & stabilizationParameter() const {
281  return m_stab_par;
282  }
283 
285  inline double & stabilizationParameter() {
286  return m_stab_par;
287  }
288 
290  std::vector<YangMillsNorms> computeYangMillsNorms(
291  const std::vector<Eigen::VectorXd> & list_dofs
292  ) const;
293 
295  template<typename outValue, typename TFct>
296  std::function<outValue(const Eigen::Vector3d &)> contractPara(const TFct &F, double t) const
297  {
298  std::function<outValue(const Eigen::Vector3d &)> f = [&F, t](const Eigen::Vector3d &x) -> outValue {return F(t, x);};
299  return f;
300  }
302  template<typename outValue, typename Fct, typename TFct>
303  std::vector<Fct> contractParaLA(const std::vector<TFct> &TF, double t) const
304  {
305  std::vector<Fct> f;
306  for (size_t i = 0; i < TF.size(); i++){
307  f.emplace_back(contractPara<outValue>(TF[i], t));
308  }
309  return f;
310  }
311 
313  double computeResidual(const Eigen::VectorXd & x) const;
314 
316  Eigen::VectorXd computeConstraint(const Eigen::VectorXd & E, const Eigen::VectorXd & A, const double nonlinear_coeff);
317 
319  double computeConstraintNorm(const Eigen::VectorXd & Ch, const size_t itersolver) const;
320 
321  std::pair<double, double> computeConditionNum() const;
322 
323  private:
324 
325  Eigen::VectorXd _compute_local_vec(size_t iT,
326  const Eigen::VectorXd & interp_f,
327  const Eigen::VectorXd & interp_dE,
328  const Eigen::VectorXd & interp_A,
329  const Eigen::VectorXd & E_i,
330  const Eigen::VectorXd & A_i,
331  double dt,
332  double theta,
333  double nonlinear_coeff
334  );
335 
336  void _assemble_local_vec(size_t iT,
337  const Eigen::VectorXd & v,
338  std::vector<std::list<Eigen::Triplet<double>>> & triplets_sys,
339  std::vector<Eigen::VectorXd> & vecs //< List of vectors for the system
340  );
341 
344  _compute_local_newton(
345  size_t iT,
346  size_t option_1FEk_2DF,
347  const Eigen::VectorXd & E_i,
348  const Eigen::VectorXd & A_i,
349  const Eigen::VectorXd & Elambda_k,
350  double dt,
351  double theta,
352  double nonlinear_coeff
353  );
354 
356  LocalStaticCondensation _compute_static_condensation(const size_t & iT) const;
357 
358  void _assemble_local_newton(
359  size_t iT,
360  const SystemVectors<Eigen::MatrixXd> & lsT,
361  std::vector<std::list<Eigen::Triplet<double> > > & triplets_sys,
362  std::vector<Eigen::VectorXd> & vecs //< List of vectors for the system
363  );
364 
366  std::vector<Eigen::MatrixXd>
367  _compute_local_contribution(
368  size_t iT,
369  double dt
370  );
371 
373  void _assemble_local_contribution(
374  size_t iT,
375  const std::vector<Eigen::MatrixXd> & lsT,
376  std::vector<std::list<Eigen::Triplet<double>>> & triplets,
377  std::vector<Eigen::VectorXd> & vecs //< List of vectors for the system
378  );
379 
381  boost::multi_array<double, 3> _compute_detnij_Pot_PkF(size_t iF) const;
382 
384  boost::multi_array<double, 3> _compute_detnij_PkF(size_t iF) const;
385 
389  boost::multi_array<double, 3> _compute_crossij_Pot_T(size_t iT) const;
390 
392  boost::multi_array<double, 3> _compute_crossij_T(size_t iT) const;
393 
395  boost::multi_array<double, 3> _integral_Pot_ijk(size_t iT) const;
396 
397  const DDRCore & m_ddrcore;
398  const LADDRCore & m_laddrcore;
399  bool m_use_threads;
400  std::ostream & m_output;
401  SerendipityProblem m_ser_pro;
402  const XDiv m_xdiv;
403  const SXGrad m_sxgrad;
404  const SXCurl m_sxcurl;
405  const SXDiv m_sxdiv;
406  const LieAlgebra m_liealg;
407  const LASXGrad m_lasxgrad;
408  const LASXCurl m_lasxcurl;
409  const LASXDiv m_lasxdiv;
410  const Eigen::VectorXd m_nloc_sc_E; // Nb of statically condensed DOFs for electric field in each cell (cell unknowns)
411  const Eigen::VectorXd m_nloc_sc_lambda; // Nb of statically condensed DOFs for lambda in each cell (cell unknowns)
412  SystemMatrixType m_A; // Matrix and RHS system
413  Eigen::VectorXd m_b; // RHS to Newton problem after static condensation
414  Eigen::VectorXd m_b_i; // RHS to nonlinear problem (F(x)=b)
415  Eigen::VectorXd m_b_k; // RHS/Residual to Newton problem (DF(dx)=b-F(x_(k-1))
416  SystemMatrixType m_sc_A; // Static condensation operator and RHS (to recover statically condensed DOFs)
417  Eigen::VectorXd m_sc_b;
418  SystemMatrixType m_laxgrad_L2; // Global laxgrad L2Product matrix
419  SystemMatrixType m_laxcurl_L2; // Global laxcurl L2Product matrix
420  SystemMatrixType m_laxdiv_L2; // Global laxdiv L2Product matrix
421  int m_nl_par;
422  double m_stab_par;
423  };
424 
426  template<typename outValue, typename Fct>
427  std::vector<Fct> sumLA(const std::vector<Fct> & LAF, const std::vector<Fct> & LAG, double lambda)
428  {
429  std::vector<Fct> values;
430  for (size_t i = 0; i < LAF.size(); i++){
431  values.emplace_back([&, i, lambda](double t, const Eigen::Vector3d & x)-> outValue {return LAF[i](t, x) + lambda * LAG[i](t, x);});
432  }
433  return values;
434  }
435 
436  //------------------------------------------------------------------------------
437  // Exact solutions
438  //------------------------------------------------------------------------------
439 
440  static const double PI = boost::math::constants::pi<double>();
441  using std::sin;
442  using std::cos;
443  using std::pow;
444 
445  //------------------------------------------------------------------------------
446  // Trigonometric solution
447  //------------------------------------------------------------------------------
448 
449  double pressure_scaling = 1.;
451  trigonometric_E1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
452  return Eigen::Vector3d(
453  -0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
454  sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
455  -0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
456  );
457  };
459  trigonometric_E2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
460  return Eigen::Vector3d(
461  0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
462  -sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
463  0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
464  );
465  };
466 
468  trigonometric_E3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
469  return Eigen::Vector3d(
470  0.5*pow(sin(PI*x(1)), 2)*cos(t),
471  sin(t)*pow(cos(PI*x(2)), 2),
472  0.5*cos(t)*pow(cos(PI*x(0)), 2)
473  );
474  };
475 
478 
480  trigonometric_A1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
481  return Eigen::Vector3d(
482  -0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
483  sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
484  -0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
485  );
486  };
487 
489  trigonometric_A2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
490  return Eigen::Vector3d(
491  -0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
492  sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
493  -0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
494  );
495  };
496 
498  trigonometric_A3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
499  return Eigen::Vector3d(-0.5*sin(t)*pow(sin(PI*x(1)), 2),
500  cos(t)*pow(cos(PI*x(2)), 2),
501  -0.5*sin(t)*pow(cos(PI*x(0)), 2)
502  );
503  };
504 
507 
509  trigonometric_B1_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
510  return Eigen::Vector3d(
511  1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)),
512  0,
513  -1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2))
514  );
515  };
516 
518  trigonometric_B2_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
519  return Eigen::Vector3d(
520  1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)),
521  0,
522  -1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2))
523  );
524  };
525 
527  trigonometric_B3_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
528  return Eigen::Vector3d(
529  2*PI*sin(PI*x(2))*cos(t)*cos(PI*x(2)),
530  -1.0*PI*sin(t)*sin(PI*x(0))*cos(PI*x(0)),
531  1.0*PI*sin(t)*sin(PI*x(1))*cos(PI*x(1))
532  );
533  };
534 
537 
539  trigonometric_B1_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
540  return Eigen::Vector3d(
541  -0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2),
542  -0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)),
543  0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3)
544  );
545  };
546 
548  trigonometric_B2_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
549  return Eigen::Vector3d(
550  0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2),
551  0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)),
552  -0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3)
553  );
554  };
555 
557  trigonometric_B3_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
558  return Eigen::Vector3d(
559  0,
560  0,
561  0
562  );
563  };
564 
567 
569  trigonometric_dtE1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
570  return Eigen::Vector3d(
571  -0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
572  sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
573  -0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
574  );
575  };
577  trigonometric_dtE2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
578  return Eigen::Vector3d(
579  -0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
580  sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
581  -0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
582  );
583  };
584 
586  trigonometric_dtE3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
587  return Eigen::Vector3d(
588  -0.5*sin(t)*pow(sin(PI*x(1)), 2),
589  cos(t)*pow(cos(PI*x(2)), 2),
590  -0.5*sin(t)*pow(cos(PI*x(0)), 2)
591  );
592  };
593 
596 
598  trigonometric_f1_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
599  return Eigen::Vector3d(
600  -0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)) + 1.5*pow(PI, 2)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
601  -3.0*pow(PI, 2)*sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
602  -0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)) + 1.5*pow(PI, 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
603  );
604  };
605 
607  trigonometric_f2_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
608  return Eigen::Vector3d(
609  -0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 1.5*pow(PI, 2)*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
610  -3.0*pow(PI, 2)*sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)) + sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
611  -0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)) + 1.5*pow(PI, 2)*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
612  );
613  };
614 
616  trigonometric_f3_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
617  return Eigen::Vector3d(
618  -0.5*sin(t)*pow(sin(PI*x(1)), 2) + 1.0*pow(PI, 2)*sin(t)*pow(sin(PI*x(1)), 2) - 1.0*pow(PI, 2)*sin(t)*pow(cos(PI*x(1)), 2),
619  2*pow(PI, 2)*pow(sin(PI*x(2)), 2)*cos(t) - 2*pow(PI, 2)*cos(t)*pow(cos(PI*x(2)), 2) + cos(t)*pow(cos(PI*x(2)), 2),
620  -1.0*pow(PI, 2)*sin(t)*pow(sin(PI*x(0)), 2) - 0.5*sin(t)*pow(cos(PI*x(0)), 2) + 1.0*pow(PI, 2)*sin(t)*pow(cos(PI*x(0)), 2)
621  );
622  };
623 
626 
628  trigonometric_f1_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
629  return Eigen::Vector3d(
630  0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(cos(PI*x(0)), 2) + (-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*cos(t)*pow(cos(PI*x(2)), 2) + 0.75*PI*pow(sin(t), 2)*sin(PI*x(0))*sin(PI*x(2))*pow(cos(PI*x(0)), 2)*cos(PI*x(1)) - 2.25*PI*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 0.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(t)*pow(cos(PI*x(2)), 3),
631  0.5*(-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*sin(t)*pow(sin(PI*x(1)), 2) - 0.5*(1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*sin(t)*pow(cos(PI*x(0)), 2) - 0.5*PI*pow(sin(t), 2)*sin(PI*x(0))*pow(sin(PI*x(1)), 3)*cos(PI*x(2)) - 0.5*PI*pow(sin(t), 2)*sin(PI*x(0))*sin(PI*x(1))*pow(cos(PI*x(1)), 2)*cos(PI*x(2)) - 0.5*PI*pow(sin(t), 2)*sin(PI*x(1))*sin(PI*x(2))*pow(cos(PI*x(0)), 3) + 2.0*PI*sin(t)*pow(sin(PI*x(2)), 2)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 1.0*PI*sin(t)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 3),
632  -0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(sin(PI*x(1)), 2) - (1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*cos(t)*pow(cos(PI*x(2)), 2) - 1.0*PI*pow(sin(t), 2)*pow(sin(PI*x(0)), 2)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 0.25*PI*pow(sin(t), 2)*sin(PI*x(0))*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(1)) - 0.25*PI*pow(sin(t), 2)*pow(cos(PI*x(0)), 3)*cos(PI*x(1))*cos(PI*x(2)) + 1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0))*pow(cos(PI*x(2)), 2)
633  );
634  };
635 
637  trigonometric_f2_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
638  return Eigen::Vector3d(
639  -0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(cos(PI*x(0)), 2) - (0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*cos(t)*pow(cos(PI*x(2)), 2) - 0.75*PI*sin(t)*sin(PI*x(0))*sin(PI*x(2))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1)) + 2.25*PI*sin(t)*pow(sin(PI*x(1)), 2)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 0.5*PI*sin(PI*x(0))*sin(PI*x(1))*pow(cos(t), 2)*pow(cos(PI*x(2)), 3),
640  0.5*(-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*sin(t)*pow(cos(PI*x(0)), 2) - 0.5*(0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*sin(t)*pow(sin(PI*x(1)), 2) + 0.5*PI*sin(t)*sin(PI*x(0))*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(2)) + 0.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(t)*pow(cos(PI*x(1)), 2)*cos(PI*x(2)) + 0.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(t)*pow(cos(PI*x(0)), 3) - 2.0*PI*pow(sin(PI*x(2)), 2)*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 1.0*PI*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 3),
641  0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(sin(PI*x(1)), 2) + (-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*cos(t)*pow(cos(PI*x(2)), 2) + 1.0*PI*sin(t)*pow(sin(PI*x(0)), 2)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 0.25*PI*sin(t)*sin(PI*x(0))*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(1)) + 0.25*PI*sin(t)*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(1))*cos(PI*x(2)) - 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*pow(cos(PI*x(2)), 2)
642  );
643  };
644 
646  trigonometric_f3_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
647  return Eigen::Vector3d(
648  0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)) - 0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)) + (0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)) - (-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
649  -0.5*(-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)) + 0.5*(0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 0.5*(-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)) + 0.5*(1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)),
650  -0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)) - (-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)) + (1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)));
651  };
652 
655 
656  //------------------------------------------------------------------------------
657  // Linear solution
658  //------------------------------------------------------------------------------
659 
661  linear_E1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
662  return Eigen::Vector3d(-x(2), -x(1), -x(0));
663  };
664 
666  linear_E2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
667  return Eigen::Vector3d(-x(2), -x(0), -x(1));
668  };
669 
671  linear_E3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
672  return Eigen::Vector3d(-x(1) - x(2), -x(0) - x(2), -x(0) - x(1));
673  };
674 
677 
679  linear_A1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
680  return Eigen::Vector3d(t*x(2), t*x(1), t*x(0));
681  };
682 
684  linear_A2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
685  return Eigen::Vector3d(t*x(2), t*x(0), t*x(1));
686  };
687 
689  linear_A3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
690  return Eigen::Vector3d(t*x(1) + t*x(2), t*x(0) + t*x(2), t*x(0) + t*x(1));
691  };
692 
695 
697  linear_B1_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
698  return Eigen::Vector3d(
699  0,
700  0,
701  0
702  );
703  };
704 
706  linear_B2_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
707  return Eigen::Vector3d(
708  t,
709  t,
710  t
711  );
712  };
713 
715  linear_B3_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
716  return Eigen::Vector3d(
717  0,
718  0,
719  0
720  );
721  };
722 
725 
727  linear_B1_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
728  return Eigen::Vector3d(
729  1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2)),
730  1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1)),
731  -1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))
732  );
733  };
734 
736  linear_B2_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
737  return Eigen::Vector3d(
738  1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)),
739  -1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)),
740  1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2))
741  );
742  };
743 
745  linear_B3_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
746  return Eigen::Vector3d(
747  -1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2),
748  1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2),
749  1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)
750  );
751  };
752 
755 
757  linear_dtE1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
758  return Eigen::Vector3d(0, 0, 0);
759  };
761  linear_dtE2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
762  return Eigen::Vector3d(0, 0, 0);
763  };
764 
766  linear_dtE3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
767  return Eigen::Vector3d(0, 0, 0);
768  };
769 
772 
774  linear_f1_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
775  return Eigen::Vector3d(0, 0, 0);
776  };
777 
779  linear_f2_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
780  return Eigen::Vector3d(0, 0, 0);
781  };
782 
784  linear_f3_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
785  return Eigen::Vector3d(0, 0, 0);
786  };
787 
790 
792  linear_f1_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
793  return Eigen::Vector3d(1.0*pow(t, 2)*x(0) + 1.0*pow(t, 2)*x(1) - t*x(0)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) + t*x(1)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(0) + t*x(1)) - (t*x(0) + t*x(1))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t) + (t*x(0) + t*x(2))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 1.0*pow(t, 2)*x(1) + 1.0*pow(t, 2)*x(2) - t*x(1)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) + t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(1) + t*x(2)) + (t*x(0) + t*x(1))*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) - (t*x(1) + t*x(2))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 1.0*pow(t, 2)*x(0) + 1.0*pow(t, 2)*x(2) + t*x(0)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) - t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(0) + t*x(2)) - (t*x(0) + t*x(2))*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) + (t*x(1) + t*x(2))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t)
794  );
795  };
796 
798  linear_f2_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
799  return Eigen::Vector3d(-1.0*pow(t, 2)*x(0) - 1.0*pow(t, 2)*x(1) - t*x(0)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) + t*x(1)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) + 1.0*t*(t*x(0) + t*x(1)) - 1.0*t*(t*x(1) + t*x(2)) + (t*x(0) + t*x(1))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1))) - (t*x(0) + t*x(2))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))), -1.0*pow(t, 2)*x(0) - 1.0*pow(t, 2)*x(2) + t*x(0)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) - t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - (t*x(0) + t*x(1))*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) + (t*x(1) + t*x(2))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))), -1.0*pow(t, 2)*x(1) - 1.0*pow(t, 2)*x(2) - t*x(1)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) + t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(0) + t*x(1)) + 1.0*t*(t*x(1) + t*x(2)) + (t*x(0) + t*x(2))*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) - (t*x(1) + t*x(2))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1)))
800  );
801  };
802 
804  linear_f3_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
805  return Eigen::Vector3d(1.0*pow(t, 2)*x(0) - 1.0*pow(t, 2)*x(1) + 1.0*pow(t, 2)*x(2) + t*x(0)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))) + t*x(0)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t) - t*x(1)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1))) - t*x(1)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 1.0*pow(t, 2)*x(2) - t*x(0)*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) + t*x(1)*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) - t*x(2)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))) + t*x(2)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 2.0*pow(t, 2)*x(1) - 1.0*pow(t, 2)*x(2) - t*x(0)*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) + t*x(1)*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) + t*x(2)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1))) - t*x(2)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t)
806  );
807  };
808 
811 
812  //------------------------------------------------------------------------------
813  // Constant solution
814  //------------------------------------------------------------------------------
815 
817  const_E1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
818  return Eigen::Vector3d::Zero();
819  };
820 
822  const_E2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
823  return Eigen::Vector3d::Zero();
824  };
825 
827  const_E3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
828  return Eigen::Vector3d::Zero();
829  };
830 
833 
835  const_A1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
836  return Eigen::Vector3d(0.2, 0.4, 0.6);
837  };
838 
840  const_A2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
841  return Eigen::Vector3d(0.8, 1.0, 1.2);
842  };
843 
845  const_A3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
846  return Eigen::Vector3d(1.4, 1.6, 1.8);
847  };
848 
851 
853  const_B1_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
854  return Eigen::Vector3d(0, 0, 0);
855  };
856 
858  const_B2_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
859  return Eigen::Vector3d(0, 0, 0);
860  };
861 
863  const_B3_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
864  return Eigen::Vector3d(0, 0, 0);
865  };
866 
869 
871  const_B1_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
872  return Eigen::Vector3d(-0.12, 0.24, -0.12);
873  };
874 
876  const_B2_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
877  return Eigen::Vector3d(0.24, -0.48, 0.24);
878  };
879 
881  const_B3_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
882  return Eigen::Vector3d(-0.12, 0.24, -0.12);
883  };
884 
887 
889  const_dtE1 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
890  return Eigen::Vector3d(0, 0, 0);
891  };
893  const_dtE2 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
894  return Eigen::Vector3d(0, 0, 0);
895  };
896 
898  const_dtE3 = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
899  return Eigen::Vector3d(0, 0, 0);
900  };
901 
904 
906  const_f1_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
907  return Eigen::Vector3d(0, 0, 0);
908  };
909 
911  const_f2_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
912  return Eigen::Vector3d(0, 0, 0);
913  };
914 
916  const_f3_linear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
917  return Eigen::Vector3d(0, 0, 0);
918  };
919 
922 
924  const_f1_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
925  return Eigen::Vector3d(1.656, 0.144, -1.368);
926  };
927 
929  const_f2_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
930  return Eigen::Vector3d(0.432, 0, -0.432);
931  };
932 
934  const_f3_nonlinear = [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
935  return Eigen::Vector3d(-0.792, -0.144, 0.504);
936  };
937 
940 
941  //------------------------------------------------------------------------------
942 
943 } // end of namespace HArDCore3D
944 
945 #endif
Construct all polynomial spaces for the DDR sequence.
Definition: ddrcore.hpp:62
Construct the spaces for the LADDR sequence.
Definition: laddrcore.hpp:19
Discrete Lie algebra valued Serendipity Hcurl space: local operators, L2 product and global interpola...
Definition: lasxcurl.hpp:18
Discrete Serendipity Hdiv space: local operators, L2 product and global interpolator.
Definition: lasxdiv.hpp:20
Discrete Serendipity Hgrad space: local operators, L2 product and global interpolator.
Definition: lasxgrad.hpp:18
Lie algebra class: mass matrix, structure constants and Lie bracket.
Definition: liealgebra.hpp:17
Discrete Serendipity Hcurl space: local operators, L2 product and global interpolator.
Definition: sxcurl.hpp:21
Discrete Serendipity Hdiv space: local operators, L2 product and global interpolator.
Definition: sxdiv.hpp:18
Discrete Serendipity Hgrad space: local operators, L2 product and global interpolator.
Definition: sxgrad.hpp:20
Construct all polynomial spaces for the DDR sequence.
Definition: serendipity_problem.hpp:40
Discrete Hdiv space: local operators, L2 product and global interpolator.
Definition: xdiv.hpp:20
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition: variabledofspace.hpp:158
static const double PI
Definition: ddr-magnetostatics.hpp:186
double pressure_scaling
Definition: ddr-stokes.hpp:249
static YangMills::TMagneticFieldType linear_B2_nonlinear
Definition: laddr-yangmills-lm.hpp:639
static YangMills::TElectricFieldType trigonometric_dtE3
Definition: laddr-yangmills-lm.hpp:489
static YangMills::TLAForcingTermType const_f_linear
Definition: laddr-yangmills-lm.hpp:824
double lambda
Norm of Lagrange multiplier.
Definition: laddr-yangmills-lm.hpp:48
static YangMills::TMagneticFieldType trigonometric_B2_linear
Definition: laddr-yangmills-lm.hpp:421
std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TForcingTermType
Base functions with time component.
Definition: laddr-yangmills-lm.hpp:60
static YangMills::TMagneticFieldType const_B1_nonlinear
Definition: laddr-yangmills-lm.hpp:774
static YangMills::TMagneticFieldType trigonometric_B1_linear
Definition: laddr-yangmills-lm.hpp:412
static YangMills::TForcingTermType linear_f1_nonlinear
Definition: laddr-yangmills-lm.hpp:695
static YangMills::TMagneticFieldType trigonometric_B3_linear
Definition: laddr-yangmills-lm.hpp:430
static YangMills::TElectricFieldType trigonometric_E1
Definition: laddr-yangmills-lm.hpp:354
static YangMills::TForcingTermType linear_f3_linear
Definition: laddr-yangmills-lm.hpp:687
static YangMills::TLAElectricFieldType linear_E
Definition: laddr-yangmills-lm.hpp:579
static YangMills::TMagneticFieldType const_B2_linear
Definition: laddr-yangmills-lm.hpp:761
static YangMills::TMagneticFieldType linear_B2_linear
Definition: laddr-yangmills-lm.hpp:609
static YangMills::TLAElectricFieldType trigonometric_A
Definition: laddr-yangmills-lm.hpp:409
static YangMills::TElectricFieldType const_A1
Definition: laddr-yangmills-lm.hpp:738
static YangMills::TForcingTermType const_f2_nonlinear
Definition: laddr-yangmills-lm.hpp:832
static YangMills::TForcingTermType trigonometric_f2_nonlinear
Definition: laddr-yangmills-lm.hpp:540
double A
Norm of potential.
Definition: laddr-yangmills-lm.hpp:47
static YangMills::TForcingTermType linear_f2_linear
Definition: laddr-yangmills-lm.hpp:682
static YangMills::TElectricFieldType trigonometric_A1
Definition: laddr-yangmills-lm.hpp:383
static YangMills::TForcingTermType trigonometric_f3_linear
Definition: laddr-yangmills-lm.hpp:519
static YangMills::TElectricFieldType trigonometric_A2
Definition: laddr-yangmills-lm.hpp:392
static YangMills::TForcingTermType linear_f1_linear
Definition: laddr-yangmills-lm.hpp:677
Eigen::MatrixXd L2v_Bkt(size_t iT, boost::multi_array< double, 3 > &intPciPcjPgk, const Eigen::VectorXd &v, const size_t &entry) const
Wrapper to plug vector into L2 integral product: int (Pcurl 1,[Pcurl 2, Pgrad 3])
Definition: laddr-yangmills-lm.cpp:1433
static YangMills::TMagneticFieldType linear_B1_nonlinear
Definition: laddr-yangmills-lm.hpp:630
static YangMills::TElectricFieldType linear_dtE2
Definition: laddr-yangmills-lm.hpp:664
static YangMills::TElectricFieldType linear_E1
Definition: laddr-yangmills-lm.hpp:564
std::vector< TForcingTermType > TLAForcingTermType
Lie algebra valued functions with time component.
Definition: laddr-yangmills-lm.hpp:68
static YangMills::TMagneticFieldType const_B2_nonlinear
Definition: laddr-yangmills-lm.hpp:779
static YangMills::TElectricFieldType const_dtE1
Definition: laddr-yangmills-lm.hpp:792
double E
Norm of electric field.
Definition: laddr-yangmills-lm.hpp:44
static YangMills::TElectricFieldType linear_dtE1
Definition: laddr-yangmills-lm.hpp:660
static YangMills::TElectricFieldType trigonometric_dtE1
Definition: laddr-yangmills-lm.hpp:472
static YangMills::TMagneticFieldType trigonometric_B3_nonlinear
Definition: laddr-yangmills-lm.hpp:460
static YangMills::TMagneticFieldType trigonometric_B1_nonlinear
Definition: laddr-yangmills-lm.hpp:442
static YangMills::TForcingTermType trigonometric_f3_nonlinear
Definition: laddr-yangmills-lm.hpp:549
static YangMills::TLAElectricFieldType trigonometric_E
Definition: laddr-yangmills-lm.hpp:380
static YangMills::TForcingTermType const_f1_nonlinear
Definition: laddr-yangmills-lm.hpp:827
static YangMills::TElectricFieldType linear_A3
Definition: laddr-yangmills-lm.hpp:592
static YangMills::TLAForcingTermType trigonometric_f_linear
Definition: laddr-yangmills-lm.hpp:528
static YangMills::TForcingTermType const_f3_nonlinear
Definition: laddr-yangmills-lm.hpp:837
static YangMills::TElectricFieldType linear_E2
Definition: laddr-yangmills-lm.hpp:569
static YangMills::TMagneticFieldType const_B3_nonlinear
Definition: laddr-yangmills-lm.hpp:784
static YangMills::TLAElectricFieldType const_A
Definition: laddr-yangmills-lm.hpp:753
static YangMills::TElectricFieldType linear_A2
Definition: laddr-yangmills-lm.hpp:587
std::vector< TElectricFieldType > TLAElectricFieldType
Definition: laddr-yangmills-lm.hpp:69
static YangMills::TElectricFieldType const_E1
Definition: laddr-yangmills-lm.hpp:720
static YangMills::TElectricFieldType linear_A1
Definition: laddr-yangmills-lm.hpp:582
static YangMills::TForcingTermType trigonometric_f1_nonlinear
Definition: laddr-yangmills-lm.hpp:531
std::vector< MagneticFieldType > LAMagneticFieldType
Definition: laddr-yangmills-lm.hpp:66
std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TElectricFieldType
Definition: laddr-yangmills-lm.hpp:61
static YangMills::TForcingTermType trigonometric_f1_linear
Definition: laddr-yangmills-lm.hpp:501
static YangMills::TElectricFieldType trigonometric_dtE2
Definition: laddr-yangmills-lm.hpp:480
static YangMills::TLAMagneticFieldType linear_B_nonlinear
Definition: laddr-yangmills-lm.hpp:657
static YangMills::TMagneticFieldType linear_B1_linear
Definition: laddr-yangmills-lm.hpp:600
static YangMills::TMagneticFieldType linear_B3_linear
Definition: laddr-yangmills-lm.hpp:618
static YangMills::TForcingTermType linear_f2_nonlinear
Definition: laddr-yangmills-lm.hpp:701
static YangMills::TMagneticFieldType const_B1_linear
Definition: laddr-yangmills-lm.hpp:756
static YangMills::TLAForcingTermType linear_f_linear
Definition: laddr-yangmills-lm.hpp:692
static YangMills::TElectricFieldType trigonometric_E2
Definition: laddr-yangmills-lm.hpp:362
static YangMills::TLAMagneticFieldType trigonometric_B_linear
Definition: laddr-yangmills-lm.hpp:439
static YangMills::TLAElectricFieldType trigonometric_dtE
Definition: laddr-yangmills-lm.hpp:498
static YangMills::TMagneticFieldType const_B3_linear
Definition: laddr-yangmills-lm.hpp:766
std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TMagneticFieldType
Definition: laddr-yangmills-lm.hpp:62
static YangMills::TElectricFieldType const_dtE2
Definition: laddr-yangmills-lm.hpp:796
static YangMills::TLAMagneticFieldType trigonometric_B_nonlinear
Definition: laddr-yangmills-lm.hpp:469
static YangMills::TElectricFieldType const_A2
Definition: laddr-yangmills-lm.hpp:743
static YangMills::TElectricFieldType const_E3
Definition: laddr-yangmills-lm.hpp:730
static YangMills::TLAMagneticFieldType linear_B_linear
Definition: laddr-yangmills-lm.hpp:627
static YangMills::TLAElectricFieldType linear_dtE
Definition: laddr-yangmills-lm.hpp:674
static YangMills::TMagneticFieldType trigonometric_B2_nonlinear
Definition: laddr-yangmills-lm.hpp:451
Eigen::SparseMatrix< double > SystemMatrixType
Definition: laddr-yangmills-lm.hpp:54
static YangMills::TLAElectricFieldType const_dtE
Definition: laddr-yangmills-lm.hpp:806
static YangMills::TForcingTermType const_f1_linear
Definition: laddr-yangmills-lm.hpp:809
static YangMills::TLAElectricFieldType linear_A
Definition: laddr-yangmills-lm.hpp:597
static YangMills::TForcingTermType const_f2_linear
Definition: laddr-yangmills-lm.hpp:814
YangMills(const DDRCore &ddrcore, const LieAlgebra &liealgebra, bool use_threads, std::ostream &output=std::cout)
Constructor.
Definition: laddr-yangmills-lm.cpp:443
static YangMills::TElectricFieldType const_dtE3
Definition: laddr-yangmills-lm.hpp:801
static YangMills::TElectricFieldType linear_E3
Definition: laddr-yangmills-lm.hpp:574
static YangMills::TLAElectricFieldType const_E
Definition: laddr-yangmills-lm.hpp:735
static YangMills::TElectricFieldType trigonometric_E3
Definition: laddr-yangmills-lm.hpp:371
Eigen::MatrixXd epsBkt_v(size_t iT, boost::multi_array< double, 3 > &ebkt_T, const Eigen::VectorXd &v) const
Calculates the matrix of *[v,.]^div in element index iT.
Definition: laddr-yangmills-lm.cpp:1323
static YangMills::TLAMagneticFieldType const_B_linear
Definition: laddr-yangmills-lm.hpp:771
static YangMills::TForcingTermType linear_f3_nonlinear
Definition: laddr-yangmills-lm.hpp:707
Eigen::MatrixXd L2v_epsBkt(size_t iT, boost::multi_array< double, 3 > &ebkt_T, const Eigen::VectorXd &v_T, const Eigen::MatrixXd &L2prod) const
Calculates the bracket inside the L2prod with v.
Definition: laddr-yangmills-lm.cpp:1466
static YangMills::TElectricFieldType const_E2
Definition: laddr-yangmills-lm.hpp:725
static YangMills::TForcingTermType const_f3_linear
Definition: laddr-yangmills-lm.hpp:819
static YangMills::TElectricFieldType const_A3
Definition: laddr-yangmills-lm.hpp:748
static YangMills::TForcingTermType trigonometric_f2_linear
Definition: laddr-yangmills-lm.hpp:510
static YangMills::TLAMagneticFieldType const_B_nonlinear
Definition: laddr-yangmills-lm.hpp:789
void assembleSystemNewton(const Eigen::VectorXd &E_i, const Eigen::VectorXd &A_i, const Eigen::VectorXd &Elambdac_k, const Eigen::VectorXd &sysVec, double dt, double theta, double nonlinear_coeff)
Assembles the system for Newton iterations.
Definition: laddr-yangmills-lm.cpp:879
std::vector< TMagneticFieldType > TLAMagneticFieldType
Definition: laddr-yangmills-lm.hpp:70
static YangMills::TMagneticFieldType linear_B3_nonlinear
Definition: laddr-yangmills-lm.hpp:648
static YangMills::TLAForcingTermType linear_f_nonlinear
Definition: laddr-yangmills-lm.hpp:713
static YangMills::TLAForcingTermType trigonometric_f_nonlinear
Definition: laddr-yangmills-lm.hpp:557
static YangMills::TLAForcingTermType const_f_nonlinear
Definition: laddr-yangmills-lm.hpp:842
static YangMills::TElectricFieldType trigonometric_A3
Definition: laddr-yangmills-lm.hpp:401
std::vector< Fct > sumLA(const std::vector< Fct > &LAF, const std::vector< Fct > &LAG, double lambda)
Template to evaluate a vector of functions.
Definition: laddr-yangmills-lm.hpp:330
static YangMills::TElectricFieldType linear_dtE3
Definition: laddr-yangmills-lm.hpp:669
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ForcingTermType
Base functions.
Definition: lasddr-yangmills.hpp:57
double computeConstraintNorm(const Eigen::VectorXd &Ch, const size_t itersolver) const
Solves a system to calculate the norm of the constraint (in the dual space of LaXgrad)
std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TForcingTermType
Base functions with time component.
Definition: lasddr-yangmills.hpp:61
const LieAlgebra & lieAlg() const
Returns the Lie algebra.
Definition: lasddr-yangmills.hpp:188
std::pair< double, double > computeConditionNum() const
size_t nbSCDOFs() const
Returns the number of statically condensed DOFs (both velocity and pressure)
Definition: lasddr-yangmills.hpp:176
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> ElectricFieldType
Definition: lasddr-yangmills.hpp:58
const Eigen::VectorXd & systemVectorNewton() const
Returns the linear system right-hand side vector.
Definition: lasddr-yangmills.hpp:240
Eigen::VectorXd computeInitialConditions(const Eigen::MatrixXd &Eh, const Eigen::MatrixXd &Ah, const double nonlinear_coeff, const size_t solver)
Computes the projected initial conditions that satisfy the constraint.
YangMillsNorms(double norm_E, double norm_A, double norm_lambda)
Constructor.
Definition: lasddr-yangmills.hpp:39
std::vector< YangMillsNorms > computeYangMillsNorms(const std::vector< Eigen::VectorXd > &list_dofs) const
Compute the discrete L2 norms, for a family of Eigen::VectorXd representing the electric field and po...
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: lasddr-yangmills.hpp:270
const Eigen::VectorXd & nonlinearRes() const
Returns the right-hand side to the nonlinear problem.
Definition: lasddr-yangmills.hpp:260
std::vector< TForcingTermType > TLAForcingTermType
Lie algebra valued functions with time component.
Definition: lasddr-yangmills.hpp:69
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: lasddr-yangmills.hpp:235
const SXCurl & xSCurl() const
Returns the space XDiv.
Definition: lasddr-yangmills.hpp:218
Eigen::VectorXd & systemVectorNewton()
Returns the linear system right-hand side vector.
Definition: lasddr-yangmills.hpp:245
size_t dimensionSpace() const
Returns the global problem dimension.
Definition: lasddr-yangmills.hpp:158
const LASXDiv & laSXDiv() const
Returns the space LAXCurl.
Definition: lasddr-yangmills.hpp:206
Eigen::VectorXd computeConstraint(const Eigen::VectorXd &E, const Eigen::VectorXd &A, const double nonlinear_coeff)
Computes the constraint [E, grad P'] + int <E, [A, P']> for the DOFs.
Eigen::VectorXd & nonlinearRes()
Returns the right-hand side to the nonlinear problem.
Definition: lasddr-yangmills.hpp:265
const SXGrad & xSGrad() const
Returns the space XDiv.
Definition: lasddr-yangmills.hpp:212
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: lasddr-yangmills.hpp:275
size_t sizeSystem() const
Returns the size of the system.
Definition: lasddr-yangmills.hpp:182
double computeResidual(const Eigen::VectorXd &x) const
Calculates residual with current matrix and rhs from given x (m_A*x = m_b)
size_t nbSCDOFs_E() const
Returns the number of statically condensed DOFs (Electric field)
Definition: lasddr-yangmills.hpp:164
std::vector< TElectricFieldType > TLAElectricFieldType
Definition: lasddr-yangmills.hpp:70
Eigen::VectorXd & systemVector()
Returns the right-hand side to the nonlinear problem.
Definition: lasddr-yangmills.hpp:255
double stoppingCrit(const Eigen::VectorXd &v, const Eigen::VectorXd &u)
Stops once changes become small.
void setNonlinearRes(const Eigen::VectorXd &Elambda_k, const Eigen::VectorXd &E_i, const Eigen::VectorXd &A_i, double dt, double theta, double nonlinear_coeff)
Sets the nonlinear residual vector b-F(x_k) in m_b_k.
Definition: lasddr-yangmills.cpp:632
const LASXCurl & laSXCurl() const
Returns the space LAXCurl.
Definition: lasddr-yangmills.hpp:200
std::vector< MagneticFieldType > LAMagneticFieldType
Definition: lasddr-yangmills.hpp:67
std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TElectricFieldType
Definition: lasddr-yangmills.hpp:62
std::function< Eigen::Vector3d(const Eigen::Vector3d &)> MagneticFieldType
Definition: lasddr-yangmills.hpp:59
const double & stabilizationParameter() const
Returns the stabilization parameter.
Definition: lasddr-yangmills.hpp:280
void setSystemVector(const Eigen::VectorXd &interp_f, const Eigen::VectorXd &interp_dE, const Eigen::VectorXd &interp_A, const Eigen::VectorXd &E_i, const Eigen::VectorXd &A_i, double dt, double theta, double nonlinear_coeff)
Sets system vector for the nonlinear problem.
void addBoundaryConditions(const LAMagneticFieldType &curl_A, double dt)
Adds boundary condition for chosen solution to the system vector.
std::function< outValue(const Eigen::Vector3d &)> contractPara(const TFct &F, double t) const
Takes a two parameter function and a value and returns a function with the first parameter fixed to v...
Definition: lasddr-yangmills.hpp:296
double & stabilizationParameter()
Returns the stabilization parameter.
Definition: lasddr-yangmills.hpp:285
std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> TMagneticFieldType
Definition: lasddr-yangmills.hpp:63
std::vector< ElectricFieldType > LAElectricFieldType
Definition: lasddr-yangmills.hpp:66
Eigen::SparseMatrix< double > SystemMatrixType
Definition: lasddr-yangmills.hpp:55
size_t nbSCDOFs_lambda() const
Returns the number of statically condensed DOFs (lambda)
Definition: lasddr-yangmills.hpp:170
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: lasddr-yangmills.hpp:230
const LASXGrad & laSXGrad() const
Returns the space LAXGrad.
Definition: lasddr-yangmills.hpp:194
const Eigen::VectorXd & systemVector() const
Returns the right-hand side to the nonlinear problem.
Definition: lasddr-yangmills.hpp:250
void assembleLinearSystem(double dt)
Assemble the global system
const SXDiv & xSDiv() const
Returns the space XDiv.
Definition: lasddr-yangmills.hpp:224
std::vector< ForcingTermType > LAForcingTermType
Lie algebra valued functions (vector of dim Lie algebra)
Definition: lasddr-yangmills.hpp:65
std::vector< TMagneticFieldType > TLAMagneticFieldType
Definition: lasddr-yangmills.hpp:71
std::vector< Fct > contractParaLA(const std::vector< TFct > &TF, double t) const
Takes a vector of two parameter functions and a value and returns a function with the first parameter...
Definition: lasddr-yangmills.hpp:303
Definition: ddr-magnetostatics.hpp:40
Structure to store information for, and perform, local static condensation.
Definition: local_static_condensation.hpp:25
Struct to store the systems and vector.
Definition: parallel_for.hpp:18