HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
HHO_MHD.hpp
Go to the documentation of this file.
1
2#include <Eigen/Dense>
3#include <Eigen/Sparse>
4
5#include "MHDTests.hpp"
6#include "mesh_builder.hpp"
7#include "vtu_writer.hpp"
8#include <basis.hpp>
9#include <chrono>
10#include <elementquad.hpp>
11#include <hybridcore.hpp>
12#include <mesh.hpp>
13#include <parallel_for.hpp>
14
16#include <TestCase/TestCase.hpp>
17
18#include <boost/program_options.hpp>
19#include <boost/timer/timer.hpp>
20
21#include <linearsolver.hpp>
22
23
30using namespace HArDCore3D;
31
32#ifndef _HHO_MHD_HPP
33#define _HHO_MHD_HPP
34
42{
43public:
44 // SolutionVector() : m_cell_deg(0), m_face_deg(0)
45 // {
46 // }
47
48 SolutionVector(Eigen::VectorXd values,
49 const Mesh *mesh_ptr,
50 const size_t cell_deg,
51 const size_t face_deg
52 )
53 : m_values(values), m_mesh_ptr(mesh_ptr), m_cell_deg(cell_deg),
54 m_face_deg(face_deg)
55 {
56 assert((size_t)(m_values.size()) == n_total_dofs);
57 }
58
60 inline Eigen::VectorXd asVectorXd() const { return m_values; }
62 inline void set_values(Eigen::VectorXd values) const { m_values = values; }
63
65 inline const size_t get_cell_deg() const { return m_cell_deg; }
66
68 inline const size_t get_face_deg() const { return m_face_deg; }
69
70 Eigen::VectorXd velocity_values() const
71 {
72 return m_values.head(n_total_cell_velocity_dofs +
73 n_total_face_velocity_dofs);
74 }
75
76 Eigen::VectorXd pressure_values() const
77 {
78 return m_values.segment(n_total_cell_velocity_dofs + n_total_face_velocity_dofs, n_total_pressure_dofs);
79 }
80
81 Eigen::VectorXd magnetic_values() const
82 {
83 return m_values.segment( n_total_cell_velocity_dofs + n_total_face_velocity_dofs + n_total_pressure_dofs, n_total_cell_magnetic_dofs + n_total_face_magnetic_dofs);
84 }
85
86 Eigen::VectorXd lagrange_values() const
87 {
88 return m_values.tail(n_total_lagrange_dofs);
89 }
90
91 Eigen::VectorXd velocity_restr(size_t iT) const
92 {
93 Eigen::VectorXd XTF = Eigen::VectorXd::Zero(n_local_velocity_dofs(iT));
94
95 XTF.head(n_local_cell_velocity_dofs) = m_values.segment( iT * n_local_cell_velocity_dofs, n_local_cell_velocity_dofs);
96 for (size_t iTF = 0; iTF < m_mesh_ptr->cell(iT)->n_faces(); iTF++)
97 {
98 size_t global_offset = n_total_cell_velocity_dofs + m_mesh_ptr->cell(iT)->face(iTF)->global_index() * n_local_face_velocity_dofs;
99 size_t local_offset = n_local_cell_velocity_dofs + iTF * n_local_face_velocity_dofs;
100 XTF.segment(local_offset, n_local_face_velocity_dofs) = m_values.segment(global_offset, n_local_face_velocity_dofs);
101 }
102
103 return XTF;
104 }
105
106 Eigen::VectorXd pressure_restr(size_t iT) const
107 {
108 size_t offset = n_total_cell_velocity_dofs + n_total_face_velocity_dofs + iT * n_local_pressure_dofs;
109
110 return m_values.segment(offset, n_local_pressure_dofs);
111 }
112
113 Eigen::VectorXd magnetic_restr(size_t iT) const
114 {
115 Eigen::VectorXd XTF = Eigen::VectorXd::Zero(n_local_magnetic_dofs(iT));
116
117 size_t offset = n_total_cell_velocity_dofs + n_total_face_velocity_dofs + n_total_pressure_dofs;
118
119 XTF.head(n_local_cell_magnetic_dofs) = m_values.segment( offset + iT * n_local_cell_magnetic_dofs, n_local_cell_magnetic_dofs);
120 for (size_t iTF = 0; iTF < m_mesh_ptr->cell(iT)->n_faces(); iTF++)
121 {
122 size_t global_offset = offset + n_total_cell_magnetic_dofs + m_mesh_ptr->cell(iT)->face(iTF)->global_index() * n_local_face_magnetic_dofs;
123 size_t local_offset = n_local_cell_magnetic_dofs + iTF * n_local_face_magnetic_dofs;
124 XTF.segment(local_offset, n_local_face_magnetic_dofs) = m_values.segment(global_offset, n_local_face_magnetic_dofs);
125 }
126
127 return XTF;
128 }
129
130 Eigen::VectorXd lagrange_restr(size_t iT) const
131 {
132 size_t offset = n_total_cell_velocity_dofs + n_total_face_velocity_dofs + n_total_pressure_dofs + n_total_cell_magnetic_dofs + n_total_face_magnetic_dofs + iT * n_local_lagrange_dofs;
133
134 return m_values.segment(offset, n_local_lagrange_dofs);
135 }
136
137 Eigen::VectorXd restr(size_t iT) const
138 {
139 Eigen::VectorXd XTF = Eigen::VectorXd::Zero( n_local_velocity_dofs(iT) + n_local_pressure_dofs + n_local_magnetic_dofs(iT) + n_local_lagrange_dofs);
140
141 XTF.head(n_local_velocity_dofs(iT)) = this->velocity_restr(iT);
142 XTF.segment(n_local_velocity_dofs(iT), n_local_pressure_dofs) = this->pressure_restr(iT);
143 XTF.segment(n_local_velocity_dofs(iT) + n_local_pressure_dofs, n_local_magnetic_dofs(iT)) = this->magnetic_restr(iT);
144 XTF.tail(n_local_lagrange_dofs) = this->lagrange_restr(iT);
145
146 return XTF;
147 }
148
151 {
152 assert(m_cell_deg == b.get_cell_deg() || m_face_deg == b.get_face_deg());
153 return SolutionVector(m_values + b.asVectorXd(), m_mesh_ptr, m_cell_deg, m_face_deg);
154 }
155
158 {
159 assert(m_cell_deg == b.get_cell_deg() || m_face_deg == b.get_face_deg());
160 return SolutionVector(m_values - b.asVectorXd(), m_mesh_ptr, m_cell_deg, m_face_deg);
161 }
162
164 double operator()(size_t index) const { return m_values(index); }
165
166 // //Assignment operator
167 // SolutionVector &operator=(const SolutionVector &a);
168
169private:
170 mutable Eigen::VectorXd m_values;
171 const Mesh *m_mesh_ptr;
172
173 const size_t m_cell_deg;
174 const size_t m_face_deg;
175
176 static const size_t dim = 3;
177
178 const size_t n_cells = m_mesh_ptr->n_cells();
179 const size_t n_faces = m_mesh_ptr->n_faces();
180
181 const size_t n_local_cell_velocity_dofs = dim * DimPoly<Cell>(m_cell_deg);
182 const size_t n_local_cell_magnetic_dofs = dim * DimPoly<Cell>(m_cell_deg);
183 const size_t n_local_face_velocity_dofs = dim * DimPoly<Face>(m_face_deg);
184 const size_t n_local_face_magnetic_dofs = dim * DimPoly<Face>(m_face_deg);
185 const size_t n_local_pressure_dofs = DimPoly<Cell>(m_cell_deg);
186 const size_t n_local_lagrange_dofs = DimPoly<Cell>(m_cell_deg);
187
188 const size_t n_total_cell_velocity_dofs = n_cells * n_local_cell_velocity_dofs;
189 const size_t n_total_cell_magnetic_dofs = n_cells * n_local_cell_magnetic_dofs;
190 const size_t n_total_face_velocity_dofs = n_faces * n_local_face_velocity_dofs;
191 const size_t n_total_face_magnetic_dofs = n_faces * n_local_face_magnetic_dofs;
192 const size_t n_total_pressure_dofs = n_cells * n_local_pressure_dofs;
193 const size_t n_total_lagrange_dofs = n_cells * n_local_lagrange_dofs;
194
195 const size_t n_total_dofs = n_total_cell_velocity_dofs + n_total_cell_magnetic_dofs + n_total_face_velocity_dofs + n_total_face_magnetic_dofs + n_total_pressure_dofs + n_total_lagrange_dofs;
196
197 inline size_t n_local_velocity_dofs(const size_t iT) const
198 {
199 return n_local_cell_velocity_dofs + m_mesh_ptr->cell(iT)->n_faces() * n_local_face_velocity_dofs;
200 }
201
202 inline size_t n_local_magnetic_dofs(const size_t iT) const
203 {
204 return n_local_cell_magnetic_dofs + m_mesh_ptr->cell(iT)->n_faces() * n_local_face_magnetic_dofs;
205 }
206};
207
209{
210public:
211 MHDModel(HybridCore &, size_t, size_t, char, char);
212
213 SolutionVector solve_with_static_cond(FType<VectorRd> f_source, FType<VectorRd> g_source, double visc, double diff, double tol, LinearSolver<Eigen::SparseMatrix<double>> & solver, bool threading = true);
215 std::vector<double> compute_errors(SolutionVector interpolant, SolutionVector discrete, double visc, double diff);
216
217private:
218 Eigen::MatrixXd gradient_term(Cell *cell, ElementQuad &elquad);
219 Eigen::MatrixXd stabilisation_term(Cell *cell, ElementQuad &elquad);
220 Eigen::MatrixXd divergence_term(Cell *cell, ElementQuad &elquad);
221 Eigen::VectorXd source_term(FType<VectorRd> f_vec, Cell *cell, ElementQuad &elquad);
222 Eigen::MatrixXd convective_term(Cell *cell, std::vector<Eigen::MatrixXd> tplus_k, std::vector<Eigen::MatrixXd> tminus_k, Eigen::VectorXd u);
223 Eigen::MatrixXd th_ijk_plus_th_ikj(Cell *cell, ElementQuad &elquad, size_t k);
224 Eigen::MatrixXd th_ijk_minus_th_ikj(Cell *cell, ElementQuad &elquad, size_t k);
225
226 Eigen::VectorXd local_interpolant(FType<VectorRd> func, Cell *cell, ElementQuad &elquad);
227 Eigen::VectorXd scalar_l2_projection(FType<double> func, Cell *cell, ElementQuad &elquad);
228
229 HybridCore &m_hho;
230 size_t m_L;
231 size_t m_K;
232 char m_u_bc;
233 char m_b_bc;
234
235 static const size_t dim = 3;
236
237 const Mesh *mesh_ptr = m_hho.get_mesh();
238
239 const size_t n_cells = mesh_ptr->n_cells();
240 const size_t n_faces = mesh_ptr->n_faces();
241 const size_t n_bdry_faces = mesh_ptr->n_b_faces();
242
243 const size_t n_local_cell_scalar_dofs = DimPoly<Cell>(m_L);
244 const size_t n_local_face_scalar_dofs = DimPoly<Face>(m_K);
245 const size_t n_local_highorder_scalar_dofs = DimPoly<Cell>(m_K + 1);
246
247 const size_t n_local_cell_vector_dofs = dim * DimPoly<Cell>(m_L);
248 const size_t n_local_face_vector_dofs = dim * DimPoly<Face>(m_K);
249 const size_t n_local_highorder_vector_dofs = dim * DimPoly<Cell>(m_K + 1);
250
251 const size_t n_total_cell_scalar_dofs = n_local_cell_scalar_dofs * n_cells;
252 // const size_t n_total_face_scalar_dofs = n_local_face_scalar_dofs * n_faces;
253
254 const size_t n_total_cell_vector_dofs = n_local_cell_vector_dofs * n_cells;
255 const size_t n_total_face_vector_dofs = n_local_face_vector_dofs * n_faces;
256 // const size_t n_total_vector_dofs = n_total_cell_vector_dofs +
257 // n_total_face_vector_dofs;
258
259 Eigen::SparseMatrix<double> inv_LTll_LTlg;
260 Eigen::SparseMatrix<double> LTgl;
261 Eigen::SparseMatrix<double> LTgg;
262
263 Eigen::VectorXd inv_LTll_rTl;
264
265 std::vector<Eigen::MatrixXd> AT;
266 std::vector<Eigen::MatrixXd> RT;
267};
268
270{
271
272 double vel_energy_difference = 0.0;
273 double vel_energy_exact = 0.0;
274
275 double mag_energy_difference = 0.0;
276 double mag_energy_exact = 0.0;
277
278 double pressure_energy_difference = 0.0;
279 double pressure_energy_exact = 0.0;
280
281 double lagrange_energy_difference = 0.0;
282 double lagrange_energy_exact = 0.0;
283
284 double vel_l2_difference = 0.0;
285 double vel_l2_exact = 0.0;
286
287 double mag_l2_difference = 0.0;
288 double mag_l2_exact = 0.0;
289
290 // Construct the local matrices using multithreading if use_threads is true
291 // std::function<void(size_t, size_t)> compute_local_errors = [&](size_t
292 // start, size_t end) -> void
293 // {
294 // for (size_t iT = start; iT < end; ++iT)
295 // {
296 for (size_t iT = 0; iT < n_cells; ++iT)
297 {
298 Eigen::VectorXd local_vel_difference = (interpolant - discrete).velocity_restr(iT);
299 Eigen::VectorXd local_vel_exact = interpolant.velocity_restr(iT);
300 Eigen::VectorXd local_mag_difference = (interpolant - discrete).magnetic_restr(iT);
301 Eigen::VectorXd local_mag_exact = interpolant.magnetic_restr(iT);
302
305 vel_energy_exact += local_vel_exact.transpose() * AT[iT] * local_vel_exact;
306 mag_energy_exact += local_mag_exact.transpose() * AT[iT] * local_mag_exact;
307
308 ElementQuad elquad(m_hho, iT, m_L + m_L, m_K + m_K);
309
310 QuadratureRule quadT = elquad.get_quadT();
311
312 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
313 BasisQuad<VectorRd> dphiT_quadT = elquad.get_dphiT_quadT();
314
315 Eigen::MatrixXd MTT = compute_gram_matrix(phiT_quadT, phiT_quadT, quadT, n_local_cell_scalar_dofs, n_local_cell_scalar_dofs, "sym");
316
317 Eigen::MatrixXd vec_MTT = Eigen::MatrixXd::Zero(n_local_cell_vector_dofs, n_local_cell_vector_dofs);
318
319 for (size_t d = 0; d < dim; d++)
320 {
321 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
322 {
323 for (size_t j = 0; j < n_local_cell_scalar_dofs; j++)
324 {
325 vec_MTT(dim * i + d, dim * j + d) = MTT(i, j);
326 }
327 }
328 }
329
330 vel_l2_difference += local_vel_difference.head(n_local_cell_vector_dofs).transpose() * vec_MTT * local_vel_difference.head(n_local_cell_vector_dofs);
331 mag_l2_difference += local_mag_difference.head(n_local_cell_vector_dofs).transpose() * vec_MTT * local_mag_difference.head(n_local_cell_vector_dofs);
332 vel_l2_exact += local_vel_exact.head(n_local_cell_vector_dofs).transpose() * vec_MTT * local_vel_exact.head(n_local_cell_vector_dofs);
333 mag_l2_exact += local_mag_exact.head(n_local_cell_vector_dofs).transpose() * vec_MTT * local_mag_exact.head(n_local_cell_vector_dofs);
334
335 Cell *cell = mesh_ptr->cell(iT);
336
337 size_t n_local_bdry_vector_dofs = cell->n_faces() * n_local_face_vector_dofs;
338 Eigen::MatrixXd M_BDRY_BDRY = Eigen::MatrixXd::Zero(n_local_bdry_vector_dofs, n_local_bdry_vector_dofs);
339
340 for (size_t iTF = 0; iTF < cell->n_faces(); iTF++)
341 {
342 const size_t offset_F = iTF * n_local_face_vector_dofs;
343
344 QuadratureRule quadF = elquad.get_quadF(iTF);
345 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
346
348 boost::extents[n_local_face_vector_dofs][quadF.size()]);
349
350 for (size_t d = 0; d < dim; d++)
351 {
352 VectorRd eF = VectorRd::Zero();
353
354 // Sort face unknowns into perpendicular, then paralell components
355 if (d == 0)
356 {
357 eF = cell->face(iTF)->normal();
358 }
359 else if (d == 1)
360 {
361 eF = cell->face(iTF)->edge(0)->tangent();
362 }
363 else if (d == 2)
364 {
365 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
366 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
367 }
368 else
369 {
370 assert(false); // should never reach here
371 }
372
373 for (size_t iqn = 0; iqn < quadF.size(); iqn++)
374 {
375 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
376 {
377 vec_phiF_quadF[dim * i + d][iqn] = phiF_quadF[i][iqn] * eF;
378 }
379 }
380 }
381
382 M_BDRY_BDRY.block(offset_F, offset_F, n_local_face_vector_dofs, n_local_face_vector_dofs) = compute_gram_matrix(vec_phiF_quadF, vec_phiF_quadF, quadF, n_local_face_vector_dofs, n_local_face_vector_dofs, "sym");
383 }
384
385 double hT = cell->diam();
390
391 Eigen::VectorXd local_pressure_difference = (interpolant - discrete).pressure_restr(iT);
392 Eigen::VectorXd local_pressure_exact = interpolant.pressure_restr(iT);
393 Eigen::VectorXd local_lagrange_difference = (interpolant - discrete).lagrange_restr(iT);
394 Eigen::VectorXd local_lagrange_exact = interpolant.lagrange_restr(iT);
395
400 }
401 // };
402
403 // // Running the local constructions in parallel
404 // parallel_for(n_cells, compute_local_errors, threading);
405
406 // If exact is zero, give absolute rather than relative errors
407 if (vel_energy_exact < 1E-12)
408 {
409 vel_energy_exact = 1.0;
410 }
411 if (mag_energy_exact < 1E-12)
412 {
413 mag_energy_exact = 1.0;
414 }
415 if (pressure_energy_exact < 1E-12)
416 {
418 }
419 if (lagrange_energy_exact < 1E-12)
420 {
422 }
423 if (vel_l2_exact < 1E-12)
424 {
425 vel_l2_exact = 1.0;
426 }
427 if (mag_l2_exact < 1E-12)
428 {
429 mag_l2_exact = 1.0;
430 }
431
432 // double velocity_energy_error = std::pow(visc, 0.5) * std::sqrt(vel_energy_difference / vel_energy_exact);
433 // double magnetic_energy_error = std::pow(diff, 0.5) * std::sqrt(mag_energy_difference / mag_energy_exact);
434
441
443 // return std::vector<double>{velocity_energy_error, magnetic_energy_error,
444 // pressure_energy_error, lagrange_energy_error};
445}
446
448{
449 Eigen::VectorXd values = Eigen::VectorXd::Zero(2 * (n_total_cell_vector_dofs + n_total_face_vector_dofs + n_total_cell_scalar_dofs));
450 for (size_t iT = 0; iT < n_cells; ++iT)
451 {
452 ElementQuad elquad(m_hho, iT, 14, 14);
453 Cell *cell = mesh_ptr->cell(iT);
454
455 Eigen::VectorXd local_Ikv = local_interpolant(velocity, cell, elquad);
456 Eigen::VectorXd local_Ikm = local_interpolant(magnetic, cell, elquad);
457 Eigen::VectorXd local_pikp = scalar_l2_projection(pressure, cell, elquad);
458
459 values.segment(iT * n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_Ikv.head(n_local_cell_vector_dofs);
460 for (size_t iTF = 0; iTF < cell->n_faces(); ++iTF)
461 {
462 values.segment(n_total_cell_vector_dofs + cell->face(iTF)->global_index() * n_local_face_vector_dofs, n_local_face_vector_dofs) = local_Ikv.segment(n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs, n_local_face_vector_dofs); // inefficient -- rewriting on internal // faces
463 }
464 values.segment(n_total_cell_vector_dofs + n_total_face_vector_dofs + iT * n_local_cell_scalar_dofs, n_local_cell_scalar_dofs) = local_pikp;
465
466 size_t magnetic_offset = n_total_cell_vector_dofs + n_total_face_vector_dofs + n_total_cell_scalar_dofs;
467 values.segment(magnetic_offset + iT * n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_Ikm.head(n_local_cell_vector_dofs);
468 for (size_t iTF = 0; iTF < cell->n_faces(); ++iTF)
469 {
470 values.segment(magnetic_offset + n_total_cell_vector_dofs + cell->face(iTF)->global_index() * n_local_face_vector_dofs, n_local_face_vector_dofs) = local_Ikm.segment(n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs, n_local_face_vector_dofs); // inefficient -- rewriting on internal // faces
471 }
472 }
473 return SolutionVector(values, mesh_ptr, m_L, m_K);
474}
475
476MHDModel::MHDModel(HybridCore &hho, size_t L, size_t K, char u_bc, char b_bc) : m_hho(hho), m_L(L), m_K(K), m_u_bc(u_bc), m_b_bc(b_bc)
477{
478 AT.resize(n_cells);
479 RT.resize(n_cells);
480}
481
482Eigen::VectorXd MHDModel::scalar_l2_projection(FType<double> func, Cell *cell, ElementQuad &elquad)
483{
484 QuadratureRule quadT = elquad.get_quadT();
485
486 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
487
488 Eigen::MatrixXd MTT = compute_gram_matrix(phiT_quadT, phiT_quadT, quadT, n_local_cell_scalar_dofs, n_local_cell_scalar_dofs, "sym");
489
490 Eigen::VectorXd RHS = integrate(func, phiT_quadT, quadT, n_local_cell_scalar_dofs);
491
492 return (MTT.ldlt()).solve(RHS);
493}
494
495Eigen::VectorXd MHDModel::local_interpolant(FType<VectorRd> func, Cell *cell, ElementQuad &elquad)
496{
497 const size_t n_local_faces = cell->n_faces();
498 const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces * n_local_face_vector_dofs;
499
500 Eigen::VectorXd ITk = Eigen::VectorXd::Zero(n_local_vector_dofs);
501
502 FType<double> func1 = [&](const VectorRd x) -> double { return func(x)(0); };
503 FType<double> func2 = [&](const VectorRd x) -> double { return func(x)(1); };
504 FType<double> func3 = [&](const VectorRd x) -> double { return func(x)(2); };
505
506 Eigen::VectorXd piTk_func1 = scalar_l2_projection(func1, cell, elquad);
507 Eigen::VectorXd piTk_func2 = scalar_l2_projection(func2, cell, elquad);
508 Eigen::VectorXd piTk_func3 = scalar_l2_projection(func3, cell, elquad);
509
510 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
511 {
512 ITk(dim * i + 0) = piTk_func1(i);
513 ITk(dim * i + 1) = piTk_func2(i);
514 ITk(dim * i + 2) = piTk_func3(i);
515 }
516
517 for (size_t iTF = 0; iTF < n_local_faces; iTF++)
518 {
519 const size_t offset_F = n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs;
520
521 QuadratureRule quadF = elquad.get_quadF(iTF);
522 size_t n_quadF = quadF.size();
523
524 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
525
526 BasisQuad<VectorRd> vec_phiF_quadF(boost::extents[n_local_face_vector_dofs][n_quadF]);
527
528 for (size_t k = 0; k < dim; k++)
529 {
530 VectorRd eF = VectorRd::Zero();
531
532 // Sort face unknowns into perpendicular, then paralell components
533 if (k == 0)
534 {
535 eF = cell->face(iTF)->normal();
536 }
537 else if (k == 1)
538 {
539 eF = cell->face(iTF)->edge(0)->tangent();
540 }
541 else if (k == 2)
542 {
543 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
544 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
545 }
546 else
547 {
548 assert(false); // should never reach here
549 }
550
551 for (size_t iqn = 0; iqn < n_quadF; iqn++)
552 {
553 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
554 {
555 vec_phiF_quadF[dim * i + k][iqn] = phiF_quadF[i][iqn] * eF;
556 }
557 }
558 }
559
560 Eigen::MatrixXd MFF = compute_gram_matrix(vec_phiF_quadF, vec_phiF_quadF, quadF, n_local_face_vector_dofs, n_local_face_vector_dofs, "sym");
561 Eigen::VectorXd RHS = integrate(func, vec_phiF_quadF, quadF, n_local_face_vector_dofs);
562
563 ITk.segment(offset_F, n_local_face_vector_dofs) = (MFF.ldlt()).solve(RHS);
564 }
565
566 return ITk;
567}
568
569Eigen::MatrixXd MHDModel::divergence_term(Cell *cell, ElementQuad &elquad)
570{
571 const size_t n_local_faces = cell->n_faces();
572 const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces * n_local_face_vector_dofs;
573
574 QuadratureRule quadT = elquad.get_quadT();
575 size_t n_quadT = quadT.size();
576
577 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
578 BasisQuad<VectorRd> dphiT_quadT = elquad.get_dphiT_quadT();
579
580 BasisQuad<double> div_vec_phiT_quadT(boost::extents[n_local_cell_vector_dofs][n_quadT]);
581
582 for (size_t k = 0; k < dim; k++)
583 {
584 // VectorRd ek = VectorRd::Zero();
585 // ek(k) = 1;
586 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
587 {
588 for (size_t iqn = 0; iqn < n_quadT; iqn++)
589 {
590 div_vec_phiT_quadT[dim * i + k][iqn] = dphiT_quadT[i][iqn](k);
591 }
592 }
593 }
594
595 Eigen::MatrixXd DT_RHS = Eigen::MatrixXd::Zero(n_local_cell_scalar_dofs, n_local_vector_dofs);
596
597 DT_RHS.topLeftCorner(n_local_cell_scalar_dofs, n_local_cell_vector_dofs) = compute_gram_matrix(phiT_quadT, div_vec_phiT_quadT, quadT, n_local_cell_scalar_dofs, n_local_cell_vector_dofs);
598
599 for (size_t iTF = 0; iTF < n_local_faces; iTF++)
600 {
601 const size_t offset_F = n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs;
602 const VectorRd &nTF = cell->face_normal(iTF);
603
604 QuadratureRule quadF = elquad.get_quadF(iTF);
605 size_t n_quadF = quadF.size();
606
607 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
608 BasisQuad<double> phiT_quadF = elquad.get_phiT_quadF(iTF);
609
610 BasisQuad<double> vec_phiF_quadF_dot_nTF(boost::extents[n_local_face_vector_dofs][n_quadF]);
611 BasisQuad<double> vec_phiT_quadF_dot_nTF(boost::extents[n_local_cell_vector_dofs][n_quadF]);
612
613 for (size_t k = 0; k < dim; k++)
614 {
615 VectorRd eT = VectorRd::Zero();
616 eT(k) = 1;
617
618 VectorRd eF = VectorRd::Zero();
619
620 // Sort face unknowns into perpendicular, then paralell components
621 if (k == 0)
622 {
623 eF = cell->face(iTF)->normal();
624 }
625 else if (k == 1)
626 {
627 eF = cell->face(iTF)->edge(0)->tangent();
628 }
629 else if (k == 2)
630 {
631 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
632 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
633 }
634 else
635 {
636 assert(false); // should never reach here
637 }
638
639 for (size_t iqn = 0; iqn < n_quadF; iqn++)
640 {
641 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
642 {
643 vec_phiF_quadF_dot_nTF[dim * i + k][iqn] = phiF_quadF[i][iqn] * eF.dot(nTF);
644 }
645 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
646 {
647 vec_phiT_quadF_dot_nTF[dim * i + k][iqn] = phiT_quadF[i][iqn] * eT.dot(nTF);
648 }
649 }
650 }
651
652 DT_RHS.topLeftCorner(n_local_cell_scalar_dofs, n_local_cell_vector_dofs) -= compute_gram_matrix(phiT_quadF, vec_phiT_quadF_dot_nTF, quadF, n_local_cell_scalar_dofs, n_local_cell_vector_dofs);
653
654 DT_RHS.block(0, offset_F, n_local_cell_scalar_dofs, n_local_face_vector_dofs) = compute_gram_matrix(phiT_quadF, vec_phiF_quadF_dot_nTF, quadF, n_local_cell_scalar_dofs, n_local_face_vector_dofs);
655 }
656
657 return -DT_RHS;
658}
659
660Eigen::MatrixXd MHDModel::gradient_term(Cell *cell, ElementQuad &elquad)
661{
662 const size_t n_local_faces = cell->n_faces();
663 const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces * n_local_face_vector_dofs;
664
665 QuadratureRule quadT = elquad.get_quadT();
666
667 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
668 BasisQuad<VectorRd> dphiT_quadT = elquad.get_dphiT_quadT();
669
670 Eigen::MatrixXd RT_RHS = Eigen::MatrixXd::Zero(n_local_highorder_vector_dofs, n_local_vector_dofs);
671
672 Eigen::MatrixXd MTT = compute_gram_matrix(phiT_quadT, phiT_quadT, quadT, n_local_cell_scalar_dofs, n_local_highorder_scalar_dofs, "sym");
673 Eigen::MatrixXd StiffT = compute_gram_matrix(dphiT_quadT, dphiT_quadT, quadT, "sym");
674
675 Eigen::MatrixXd vec_MTT = Eigen::MatrixXd::Zero(n_local_cell_vector_dofs, n_local_highorder_vector_dofs);
676
677 Eigen::MatrixXd vec_StiffT = Eigen::MatrixXd::Zero(n_local_highorder_vector_dofs, n_local_highorder_vector_dofs);
678
679 for (size_t k = 0; k < dim; k++)
680 {
681 for (size_t i = 0; i < n_local_highorder_scalar_dofs; i++)
682 {
683 for (size_t j = 0; j < n_local_highorder_scalar_dofs; j++)
684 {
685 if (i < n_local_cell_scalar_dofs)
686 {
687 vec_MTT(dim * i + k, dim * j + k) = MTT(i, j);
688 }
689 vec_StiffT(dim * i + k, dim * j + k) = StiffT(i, j);
690 }
691 }
692 }
693
694 Eigen::MatrixXd LT = (vec_MTT.topRightCorner(dim, n_local_highorder_vector_dofs)).transpose(); // need both components of the constant term
695 Eigen::MatrixXd LTtLT = LT * (LT.transpose());
696
697 double scalT = vec_StiffT.norm() / LTtLT.norm();
698
699 RT_RHS.topLeftCorner(n_local_highorder_vector_dofs, n_local_cell_vector_dofs) = vec_StiffT.topLeftCorner(n_local_highorder_vector_dofs, n_local_cell_vector_dofs) + scalT * LTtLT.topLeftCorner(n_local_highorder_vector_dofs, n_local_cell_vector_dofs);
700
701 for (size_t iTF = 0; iTF < n_local_faces; iTF++)
702 {
703 const size_t offset_F = n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs;
704 const VectorRd &nTF = cell->face_normal(iTF);
705
706 QuadratureRule quadF = elquad.get_quadF(iTF);
707 size_t n_quadF = quadF.size();
708
709 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
710 BasisQuad<double> phiT_quadF = elquad.get_phiT_quadF(iTF);
711 BasisQuad<VectorRd> dphiT_quadF = elquad.get_dphiT_quadF(iTF);
712
713 BasisQuad<VectorRd> vec_phiF_quadF(boost::extents[n_local_face_vector_dofs][n_quadF]);
714 BasisQuad<VectorRd> vec_phiT_quadF(boost::extents[n_local_cell_vector_dofs][n_quadF]);
715
716 BasisQuad<VectorRd> nTF_doT_grad_vec_phiT_quadF(boost::extents[n_local_highorder_vector_dofs][n_quadF]);
717
718 for (size_t k = 0; k < dim; k++)
719 {
720 VectorRd eT = VectorRd::Zero();
721 eT(k) = 1;
722
723 VectorRd eF = VectorRd::Zero();
724
725 // Sort face unknowns into perpendicular, then paralell components
726 if (k == 0)
727 {
728 eF = cell->face(iTF)->normal();
729 }
730 else if (k == 1)
731 {
732 eF = cell->face(iTF)->edge(0)->tangent();
733 }
734 else if (k == 2)
735 {
736 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
737 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
738 }
739 else
740 {
741 assert(false); // should never reach here
742 }
743
744 for (size_t iqn = 0; iqn < n_quadF; iqn++)
745 {
746 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
747 {
748 vec_phiF_quadF[dim * i + k][iqn] = phiF_quadF[i][iqn] * eF;
749 }
750 for (size_t i = 0; i < n_local_highorder_scalar_dofs; i++)
751 {
752 if (i < n_local_cell_scalar_dofs)
753 {
754 vec_phiT_quadF[dim * i + k][iqn] = phiT_quadF[i][iqn] * eT;
755 }
756 nTF_doT_grad_vec_phiT_quadF[dim * i + k][iqn] = (dphiT_quadF[i][iqn].dot(nTF)) * eT;
757 }
758 }
759 }
760
761 RT_RHS.topLeftCorner(n_local_highorder_vector_dofs, n_local_cell_vector_dofs) -= compute_gram_matrix(nTF_doT_grad_vec_phiT_quadF, vec_phiT_quadF, quadF, n_local_highorder_vector_dofs, n_local_cell_vector_dofs);
762
763 RT_RHS.block(0, offset_F, n_local_highorder_vector_dofs, n_local_face_vector_dofs) = compute_gram_matrix(nTF_doT_grad_vec_phiT_quadF, vec_phiF_quadF, quadF, n_local_highorder_vector_dofs, n_local_face_vector_dofs);
764 }
765
766 RT[cell->global_index()] = (vec_StiffT + scalT * LTtLT).ldlt().solve(RT_RHS);
767
768 return RT[cell->global_index()].transpose() * vec_StiffT * RT[cell->global_index()];
769}
770
771Eigen::MatrixXd MHDModel::th_ijk_plus_th_ikj(Cell *cell, ElementQuad &elquad, size_t k)
772{
773 const size_t n_local_faces = cell->n_faces();
774 const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces * n_local_face_vector_dofs;
775
776 Eigen::MatrixXd T = Eigen::MatrixXd::Zero(n_local_vector_dofs, n_local_vector_dofs);
777
779 if (k >= n_local_vector_dofs)
780 {
782 }
783
784 if (k < n_local_cell_vector_dofs)
785 {
786 QuadratureRule quadT = elquad.get_quadT();
787
788 size_t k_row = (k % dim);
789 size_t k_scalar = (k / dim); // integer division
790
791 VectorRd eK = VectorRd::Zero();
792 eK(k_row) = 1.0;
793
794 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
795 BasisQuad<VectorRd> dphiT_quadT = elquad.get_dphiT_quadT();
796
797 BasisQuad<VectorRd> vec_phiT_quadT(boost::extents[n_local_cell_vector_dofs][quadT.size()]);
798 BasisQuad<VectorRd> phiT_quadT_dot_grad_basis_k(boost::extents[n_local_cell_vector_dofs][quadT.size()]);
799 BasisQuad<VectorRd> basis_k_dot_grad_phiT_quadT(boost::extents[n_local_cell_vector_dofs][quadT.size()]);
800 BasisQuad<MatrixRd> grad_phiT_quadT(boost::extents[n_local_cell_vector_dofs][quadT.size()]);
801 BasisQuad<MatrixRd> basis_k_tensor_phiT_quadT(boost::extents[n_local_cell_vector_dofs][quadT.size()]);
802
803 for (size_t d = 0; d < dim; ++d)
804 {
805 VectorRd eT = VectorRd::Zero();
806 eT(d) = 1.0;
807 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
808 {
809 for (size_t iqn = 0; iqn < quadT.size(); iqn++)
810 {
811 vec_phiT_quadT[dim * i + d][iqn] = phiT_quadT[i][iqn] * eT;
814 grad_phiT_quadT[dim * i + d][iqn] = MatrixRd::Zero();
815 grad_phiT_quadT[dim * i + d][iqn].row(d) = dphiT_quadT[i][iqn].transpose();
816 // grad_phiT_quadT[dim * i + d][iqn].col(d) = dphiT_quadT[i][iqn];
817 basis_k_tensor_phiT_quadT[dim * i + d][iqn] = MatrixRd::Zero();
819 // basis_k_tensor_phiT_quadT[dim * i + d][iqn](d, k_row) =
820 // phiT_quadT[k_scalar][iqn] * phiT_quadT[i][iqn];
821 }
822 }
823 }
824
825 Eigen::MatrixXd uT_dot_grad_wT_dot_vT = compute_gram_matrix(vec_phiT_quadT, basis_k_dot_grad_phiT_quadT, quadT, n_local_cell_vector_dofs, n_local_cell_vector_dofs);
826
827 T.topLeftCorner(n_local_cell_vector_dofs, n_local_cell_vector_dofs) = compute_gram_matrix(vec_phiT_quadT, phiT_quadT_dot_grad_basis_k, quadT, n_local_cell_vector_dofs, n_local_cell_vector_dofs) + uT_dot_grad_wT_dot_vT - compute_gram_matrix(grad_phiT_quadT, basis_k_tensor_phiT_quadT, quadT, n_local_cell_vector_dofs, n_local_cell_vector_dofs) - uT_dot_grad_wT_dot_vT.transpose();
828 // T.topLeftCorner(n_local_cell_vector_dofs, n_local_cell_vector_dofs) =
829 // compute_gram_matrix(vec_phiT_quadT, phiT_quadT_dot_grad_basis_k,
830 // quadT).transpose() + uT_dot_grad_wT_dot_vT.transpose() -
831 // compute_gram_matrix(basis_k_tensor_phiT_quadT, grad_phiT_quadT,
832 // quadT).transpose() - uT_dot_grad_wT_dot_vT;
833
834 for (size_t iTF = 0; iTF < n_local_faces; ++iTF)
835 {
836 const size_t offset_F = n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs;
837 VectorRd nTF = cell->face_normal(iTF);
838
839 QuadratureRule quadF = elquad.get_quadF(iTF);
840
841 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
842 BasisQuad<double> phiT_quadF = elquad.get_phiT_quadF(iTF);
843
844 BasisQuad<VectorRd> basis_k_dot_nTF_phiT_quadF(boost::extents[n_local_cell_vector_dofs][quadF.size()]);
845 BasisQuad<double> phiT_quadF_dot_nTF(boost::extents[n_local_cell_vector_dofs][quadF.size()]);
846
847 BasisQuad<double> basis_k_dot_phiF_quadF(boost::extents[n_local_face_vector_dofs][quadF.size()]);
848 BasisQuad<VectorRd> vec_phiF_quadF(boost::extents[n_local_face_vector_dofs][quadF.size()]);
849
850 for (size_t d = 0; d < dim; ++d)
851 {
852 // VectorRd eF = ((d == 0) ? cell->face(iTF)->normal() :
853 // cell->face(iTF)->tangent());
854
855 VectorRd eF = VectorRd::Zero();
856
857 // Sort face unknowns into perpendicular, then paralell components
858 if (d == 0)
859 {
860 eF = cell->face(iTF)->normal();
861 }
862 else if (d == 1)
863 {
864 eF = cell->face(iTF)->edge(0)->tangent();
865 }
866 else if (d == 2)
867 {
868 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
869 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
870 }
871 else
872 {
873 assert(false); // should never reach here
874 }
875
876 VectorRd eT = VectorRd::Zero();
877 eT(d) = 1.0;
878
879 for (size_t iqn = 0; iqn < quadF.size(); iqn++)
880 {
881 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
882 {
884 phiT_quadF_dot_nTF[dim * i + d][iqn] = phiT_quadF[i][iqn] * eT.dot(nTF);
885 }
886 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
887 {
889 vec_phiF_quadF[dim * i + d][iqn] = phiF_quadF[i][iqn] * eF;
890 }
891 }
892 }
893
895
896 T.block(0, offset_F, n_local_cell_vector_dofs, n_local_face_vector_dofs) = uT_dot_nTF_vT_dot_wF;
897 T.block(offset_F, 0, n_local_face_vector_dofs, n_local_cell_vector_dofs) = -(compute_gram_matrix(basis_k_dot_phiF_quadF, phiT_quadF_dot_nTF, quadF) + uT_dot_nTF_vT_dot_wF.transpose());
898 }
899
900 return 0.5 * T;
901 }
902
903 // const size_t iT = cell->global_index();
904 // const size_t n_local_faces = cell->n_faces();
905 // const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces
906 // * n_local_face_vector_dofs;
907
909
910 size_t iTF = (k - n_local_cell_vector_dofs) / n_local_face_vector_dofs; // integer division
911 size_t k_basis = (k - n_local_cell_vector_dofs) % n_local_face_vector_dofs;
912 size_t k_scalar = k_basis / dim; // integer division
913
914 VectorRd eF = VectorRd::Zero();
915
916 // Sort face unknowns into perpendicular, then paralell components
917 if (k_basis % dim == 0)
918 {
919 eF = cell->face(iTF)->normal();
920 }
921 else if (k_basis % dim == 1)
922 {
923 eF = cell->face(iTF)->edge(0)->tangent();
924 }
925 else if (k_basis % dim == 2)
926 {
927 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
928 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
929 }
930 else
931 {
932 assert(false); // should never reach here
933 }
934
935 VectorRd nTF = cell->face_normal(iTF);
936
937 QuadratureRule quadF = elquad.get_quadF(iTF);
938
939 BasisQuad<double> phiT_quadF = elquad.get_phiT_quadF(iTF);
940 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
941
942 BasisQuad<double> phiT_quadF_dot_nTF(boost::extents[n_local_cell_vector_dofs][quadF.size()]);
943 BasisQuad<double> phiT_quadF_dot_basis_k(boost::extents[n_local_cell_vector_dofs][quadF.size()]);
944
945 for (size_t d = 0; d < dim; ++d)
946 {
947 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
948 {
949 for (size_t iqn = 0; iqn < quadF.size(); iqn++)
950 {
951 phiT_quadF_dot_nTF[dim * i + d][iqn] = phiT_quadF[i][iqn] * nTF(d); // eT.dot(nTF)
952 phiT_quadF_dot_basis_k[dim * i + d][iqn] = phiT_quadF[i][iqn] * phiF_quadF[k_scalar][iqn] * eF(d); // eT.dot(eF)
953 }
954 }
955 }
956
957 // Eigen::MatrixXd T = Eigen::MatrixXd::Zero(n_local_vector_dofs,
958 // n_local_vector_dofs);
959 T.topLeftCorner(n_local_cell_vector_dofs, n_local_cell_vector_dofs) = 0.5 * compute_gram_matrix(phiT_quadF_dot_basis_k, phiT_quadF_dot_nTF, quadF, n_local_cell_vector_dofs, n_local_cell_vector_dofs);
960 return T;
961
962 // return 0.5 * compute_gram_matrix(phiT_quadF_dot_basis_k,
963 // phiT_quadF_dot_nTF, quadF, n_local_cell_vector_dofs,
964 // n_local_cell_vector_dofs);
965}
966
967Eigen::MatrixXd MHDModel::th_ijk_minus_th_ikj(Cell *cell, ElementQuad &elquad, size_t k)
968{
969 const size_t n_local_faces = cell->n_faces();
970 const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces * n_local_face_vector_dofs;
971
972 Eigen::MatrixXd T = Eigen::MatrixXd::Zero(n_local_vector_dofs, n_local_vector_dofs);
973
975 if (k >= n_local_vector_dofs)
976 {
978 }
979
980 if (k < n_local_cell_vector_dofs)
981 {
982 QuadratureRule quadT = elquad.get_quadT();
983
984 size_t k_row = (k % dim);
985 size_t k_scalar = (k / dim); // integer division
986
987 VectorRd eK = VectorRd::Zero();
988 eK(k_row) = 1.0;
989
990 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
991 BasisQuad<VectorRd> dphiT_quadT = elquad.get_dphiT_quadT();
992
993 BasisQuad<VectorRd> vec_phiT_quadT( boost::extents[n_local_cell_vector_dofs][quadT.size()]);
994 BasisQuad<VectorRd> phiT_quadT_dot_grad_basis_k( boost::extents[n_local_cell_vector_dofs][quadT.size()]);
995 BasisQuad<VectorRd> basis_k_dot_grad_phiT_quadT( boost::extents[n_local_cell_vector_dofs][quadT.size()]);
996 BasisQuad<MatrixRd> grad_phiT_quadT( boost::extents[n_local_cell_vector_dofs][quadT.size()]);
997 BasisQuad<MatrixRd> basis_k_tensor_phiT_quadT( boost::extents[n_local_cell_vector_dofs][quadT.size()]);
998
999 for (size_t d = 0; d < dim; ++d)
1000 {
1001 VectorRd eT = VectorRd::Zero();
1002 eT(d) = 1.0;
1003 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
1004 {
1005 for (size_t iqn = 0; iqn < quadT.size(); iqn++)
1006 {
1007 vec_phiT_quadT[dim * i + d][iqn] = phiT_quadT[i][iqn] * eT;
1010 grad_phiT_quadT[dim * i + d][iqn] = MatrixRd::Zero();
1011 grad_phiT_quadT[dim * i + d][iqn].row(d) = dphiT_quadT[i][iqn].transpose();
1012 basis_k_tensor_phiT_quadT[dim * i + d][iqn] = MatrixRd::Zero();
1014 }
1015 }
1016 }
1017
1019
1021
1022 for (size_t iTF = 0; iTF < n_local_faces; ++iTF)
1023 {
1024 const size_t offset_F = n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs;
1025 VectorRd nTF = cell->face_normal(iTF);
1026
1027 QuadratureRule quadF = elquad.get_quadF(iTF);
1028
1029 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
1030 BasisQuad<double> phiT_quadF = elquad.get_phiT_quadF(iTF);
1031
1032 BasisQuad<VectorRd> basis_k_dot_nTF_phiT_quadF( boost::extents[n_local_cell_vector_dofs][quadF.size()]);
1033 BasisQuad<double> phiT_quadF_dot_nTF( boost::extents[n_local_cell_vector_dofs][quadF.size()]);
1034
1035 BasisQuad<double> basis_k_dot_phiF_quadF( boost::extents[n_local_face_vector_dofs][quadF.size()]);
1036 BasisQuad<VectorRd> vec_phiF_quadF( boost::extents[n_local_face_vector_dofs][quadF.size()]);
1037
1038 for (size_t d = 0; d < dim; ++d)
1039 {
1040 VectorRd eF = VectorRd::Zero();
1041
1042 // Sort face unknowns into perpendicular, then paralell components
1043 if (d == 0)
1044 {
1045 eF = cell->face(iTF)->normal();
1046 }
1047 else if (d == 1)
1048 {
1049 eF = cell->face(iTF)->edge(0)->tangent();
1050 }
1051 else if (d == 2)
1052 {
1053 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
1054 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
1055 }
1056 else
1057 {
1058 assert(false); // should never reach here
1059 }
1060
1061 VectorRd eT = VectorRd::Zero();
1062 eT(d) = 1.0;
1063
1064 for (size_t iqn = 0; iqn < quadF.size(); iqn++)
1065 {
1066 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
1067 {
1069 phiT_quadF_dot_nTF[dim * i + d][iqn] = phiT_quadF[i][iqn] * eT.dot(nTF);
1070 }
1071 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
1072 {
1074 vec_phiF_quadF[dim * i + d][iqn] = phiF_quadF[i][iqn] * eF;
1075 }
1076 }
1077 }
1078
1080
1081 T.block(0, offset_F, n_local_cell_vector_dofs, n_local_face_vector_dofs) = uT_dot_nTF_vT_dot_wF;
1082 T.block(offset_F, 0, n_local_face_vector_dofs, n_local_cell_vector_dofs) = -uT_dot_nTF_vT_dot_wF.transpose() + compute_gram_matrix(basis_k_dot_phiF_quadF, phiT_quadF_dot_nTF, quadF);
1083 }
1084
1085 return 0.5 * T;
1086 }
1087
1089
1090 size_t iTF = (k - n_local_cell_vector_dofs) / n_local_face_vector_dofs; // integer division
1091 size_t k_basis = (k - n_local_cell_vector_dofs) % n_local_face_vector_dofs;
1092 size_t k_scalar = k_basis / dim; // integer division
1093
1094 VectorRd eF = VectorRd::Zero();
1095
1096 // Sort face unknowns into perpendicular, then paralell components
1097 if (k_basis % dim == 0)
1098 {
1099 eF = cell->face(iTF)->normal();
1100 }
1101 else if (k_basis % dim == 1)
1102 {
1103 eF = cell->face(iTF)->edge(0)->tangent();
1104 }
1105 else if (k_basis % dim == 2)
1106 {
1107 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
1108 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
1109 }
1110 else
1111 {
1112 assert(false); // should never reach here
1113 }
1114 VectorRd nTF = cell->face_normal(iTF);
1115
1116 QuadratureRule quadF = elquad.get_quadF(iTF);
1117
1118 BasisQuad<double> phiT_quadF = elquad.get_phiT_quadF(iTF);
1119 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
1120
1121 BasisQuad<double> phiT_quadF_dot_nTF( boost::extents[n_local_cell_vector_dofs][quadF.size()]);
1122 BasisQuad<double> phiT_quadF_dot_basis_k( boost::extents[n_local_cell_vector_dofs][quadF.size()]);
1123
1124 for (size_t d = 0; d < dim; ++d)
1125 {
1126 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
1127 {
1128 for (size_t iqn = 0; iqn < quadF.size(); iqn++)
1129 {
1130 phiT_quadF_dot_nTF[dim * i + d][iqn] = phiT_quadF[i][iqn] * nTF(d); // eT.dot(nTF)
1131 phiT_quadF_dot_basis_k[dim * i + d][iqn] = phiT_quadF[i][iqn] * phiF_quadF[k_scalar][iqn] * eF(d); // eT.dot(eF)
1132 }
1133 }
1134 }
1135
1136 // Eigen::MatrixXd T = Eigen::MatrixXd::Zero(n_local_vector_dofs,
1137 // n_local_vector_dofs);
1138 T.topLeftCorner(n_local_cell_vector_dofs, n_local_cell_vector_dofs) = -0.5 * compute_gram_matrix(phiT_quadF_dot_basis_k, phiT_quadF_dot_nTF, quadF, n_local_cell_vector_dofs, n_local_cell_vector_dofs);
1139 return T;
1140}
1141
1142Eigen::MatrixXd MHDModel::convective_term(Cell *cell, std::vector<Eigen::MatrixXd> tplus_k, std::vector<Eigen::MatrixXd> tminus_k, Eigen::VectorXd u)
1143{
1144 const size_t n_local_faces = cell->n_faces();
1145 const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces * n_local_face_vector_dofs;
1146 const size_t n_local_dofs = 2 * n_local_vector_dofs;
1147 assert((size_t)u.size() == n_local_dofs);
1148
1149 Eigen::MatrixXd T = Eigen::MatrixXd::Zero(n_local_dofs, n_local_dofs);
1150
1151 for (size_t k = 0; k < n_local_vector_dofs; ++k)
1152 {
1153 T.topLeftCorner(n_local_vector_dofs, n_local_vector_dofs) += (u(k) + u(k + n_local_vector_dofs)) * tplus_k[k];
1154 T.topRightCorner(n_local_vector_dofs, n_local_vector_dofs) += (u(k) + u(k + n_local_vector_dofs)) * -tplus_k[k];
1155 T.bottomRightCorner(n_local_vector_dofs, n_local_vector_dofs) += (u(k) + u(k + n_local_vector_dofs)) * tminus_k[k];
1156 T.bottomLeftCorner(n_local_vector_dofs, n_local_vector_dofs) += (u(k) + u(k + n_local_vector_dofs)) * -tminus_k[k];
1157 }
1158
1159 return T;
1160}
1161
1162Eigen::MatrixXd MHDModel::stabilisation_term(Cell *cell, ElementQuad &elquad)
1163{
1164 const size_t n_local_faces = cell->n_faces();
1165 const size_t n_local_bdry_vector_dofs = n_local_faces * n_local_face_vector_dofs;
1166
1167 QuadratureRule quadT = elquad.get_quadT();
1168
1169 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
1170 BasisQuad<VectorRd> dphiT_quadT = elquad.get_dphiT_quadT();
1171
1172 Eigen::MatrixXd MTT = compute_gram_matrix(phiT_quadT, phiT_quadT, quadT, n_local_cell_scalar_dofs, n_local_highorder_scalar_dofs, "sym");
1173 Eigen::MatrixXd StiffT = compute_gram_matrix(dphiT_quadT, dphiT_quadT, quadT, n_local_cell_scalar_dofs, n_local_cell_scalar_dofs, "sym");
1174
1175 Eigen::MatrixXd vec_MTT = Eigen::MatrixXd::Zero(n_local_cell_vector_dofs, n_local_highorder_vector_dofs);
1176 Eigen::MatrixXd vec_StiffT = Eigen::MatrixXd::Zero(n_local_cell_vector_dofs, n_local_cell_vector_dofs);
1177
1178 for (size_t k = 0; k < dim; k++)
1179 {
1180 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
1181 {
1182 for (size_t j = 0; j < n_local_highorder_scalar_dofs; j++)
1183 {
1184 if (j < n_local_cell_scalar_dofs)
1185 {
1186 vec_StiffT(dim * i + k, dim * j + k) = StiffT(i, j);
1187 }
1188 vec_MTT(dim * i + k, dim * j + k) = MTT(i, j);
1189 }
1190 }
1191 }
1192
1193 Eigen::MatrixXd diffT = (vec_MTT.topLeftCorner(n_local_cell_vector_dofs, n_local_cell_vector_dofs).ldlt()).solve(vec_MTT) * RT[cell->global_index()];
1194 diffT.topLeftCorner(n_local_cell_vector_dofs, n_local_cell_vector_dofs) -= Eigen::MatrixXd::Identity(n_local_cell_vector_dofs, n_local_cell_vector_dofs);
1195
1196 std::vector<Eigen::Triplet<double>> triplets_M_BDRY_BDRY;
1197 std::vector<Eigen::Triplet<double>> triplets_M_BDRY_T;
1198 std::vector<Eigen::Triplet<double>> triplets_SCALE_MAT;
1199
1201 Eigen::SparseMatrix<double> M_BDRY_T(n_local_bdry_vector_dofs, n_local_highorder_vector_dofs);
1203
1204 for (size_t iTF = 0; iTF < n_local_faces; iTF++)
1205 {
1206 const size_t offset_F = iTF * n_local_face_vector_dofs;
1207
1208 QuadratureRule quadF = elquad.get_quadF(iTF);
1209 size_t n_quadF = quadF.size();
1210
1211 BasisQuad<double> phiF_quadF = elquad.get_phiF_quadF(iTF);
1212 BasisQuad<double> phiT_quadF = elquad.get_phiT_quadF(iTF);
1213
1214 BasisQuad<VectorRd> vec_phiF_quadF( boost::extents[n_local_face_vector_dofs][n_quadF]);
1215 BasisQuad<VectorRd> vec_phiT_quadF( boost::extents[n_local_highorder_vector_dofs][n_quadF]);
1216
1217 for (size_t d = 0; d < dim; d++)
1218 {
1219 VectorRd eT = VectorRd::Zero();
1220 eT(d) = 1;
1221
1222 VectorRd eF = VectorRd::Zero();
1223
1224 // Sort face unknowns into perpendicular, then paralell components
1225 if (d == 0)
1226 {
1227 eF = cell->face(iTF)->normal();
1228 }
1229 else if (d == 1)
1230 {
1231 eF = cell->face(iTF)->edge(0)->tangent();
1232 }
1233 else if (d == 2)
1234 {
1235 // eF = mesh_ptr->cell(iT)->face(iTF)->edge(1)->tangent();
1236 eF = (cell->face(iTF)->normal()).cross(cell->face(iTF)->edge(0)->tangent());
1237 }
1238 else
1239 {
1240 assert(false); // should never reach here
1241 }
1242
1243 for (size_t iqn = 0; iqn < n_quadF; iqn++)
1244 {
1245 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
1246 {
1247 vec_phiF_quadF[dim * i + d][iqn] = phiF_quadF[i][iqn] * eF;
1248 }
1249 for (size_t i = 0; i < n_local_highorder_scalar_dofs; i++)
1250 {
1251 vec_phiT_quadF[dim * i + d][iqn] = phiT_quadF[i][iqn] * eT;
1252 }
1253 }
1254 }
1255
1257 Eigen::MatrixXd vec_MFT = compute_gram_matrix(vec_phiF_quadF, vec_phiT_quadF, quadF, n_local_face_vector_dofs, n_local_highorder_vector_dofs);
1258
1259 Eigen::MatrixXd LOCAL_SCALING = (1.0 / (cell->diam())) * Eigen::MatrixXd::Identity(n_local_face_vector_dofs, n_local_face_vector_dofs);
1260
1261 for (size_t i = 0; i < n_local_face_vector_dofs; ++i)
1262 {
1263 for (size_t j = 0; j < n_local_face_vector_dofs; ++j)
1264 {
1265 triplets_M_BDRY_BDRY.emplace_back(offset_F + i, offset_F + j, vec_MFF(i, j));
1266 triplets_SCALE_MAT.emplace_back(offset_F + i, offset_F + j, LOCAL_SCALING(i, j));
1267 }
1268 for (size_t k = 0; k < n_local_highorder_vector_dofs; ++k)
1269 {
1270 triplets_M_BDRY_T.emplace_back(offset_F + i, k, vec_MFT(i, k));
1271 }
1272 }
1273 }
1274
1275 M_BDRY_BDRY.setFromTriplets(std::begin(triplets_M_BDRY_BDRY), std::end(triplets_M_BDRY_BDRY));
1276 M_BDRY_T.setFromTriplets(std::begin(triplets_M_BDRY_T), std::end(triplets_M_BDRY_T));
1277 BDRY_SCALE_MAT.setFromTriplets(std::begin(triplets_SCALE_MAT), std::end(triplets_SCALE_MAT));
1278
1279 Eigen::BiCGSTAB<Eigen::SparseMatrix<double>> solver;
1280 solver.compute(M_BDRY_BDRY);
1281
1282 Eigen::MatrixXd inv_M_BDRY_BDRY_M_BDRY_T = solver.solve(M_BDRY_T);
1283
1284 // Eigen::MatrixXd diffBDRY = solver.solve(M_BDRY_T) * RT[cell->global_index()];
1285 Eigen::MatrixXd diffBDRY = inv_M_BDRY_BDRY_M_BDRY_T * RT[cell->global_index()];
1287
1288 return diffT.transpose() * vec_StiffT * diffT + diffBDRY.transpose() * BDRY_SCALE_MAT * M_BDRY_BDRY * diffBDRY;
1289 // return diffBDRY.transpose() * BDRY_SCALE_MAT * M_BDRY_BDRY * diffBDRY;
1290
1291 // Eigen::MatrixXd diffT_on_bdry = inv_M_BDRY_BDRY_M_BDRY_T.topLeftCorner(n_local_bdry_vector_dofs, n_local_cell_vector_dofs) * diffT;
1292
1293 // return (diffT_on_bdry - diffBDRY).transpose() * BDRY_SCALE_MAT * M_BDRY_BDRY * (diffT_on_bdry - diffBDRY);
1294}
1295
1296Eigen::VectorXd MHDModel::source_term(FType<VectorRd> f_vec, Cell *cell, ElementQuad &elquad)
1297{
1298 const size_t n_local_faces = cell->n_faces();
1299 const size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_faces * n_local_face_vector_dofs;
1300
1301 // get quadrature rule and basis function on quadrature nodes
1302 QuadratureRule quadT = elquad.get_quadT();
1303 BasisQuad<double> phiT_quadT = elquad.get_phiT_quadT();
1304 BasisQuad<VectorRd> vec_phiT_quadT(boost::extents[n_local_highorder_vector_dofs][quadT.size()]);
1305
1306 Eigen::VectorXd source_vec = Eigen::VectorXd::Zero(n_local_vector_dofs);
1307
1308 for (size_t k = 0; k < dim; k++)
1309 {
1310 VectorRd eT = VectorRd::Zero();
1311 eT(k) = 1;
1312
1313 for (size_t iqn = 0; iqn < quadT.size(); iqn++)
1314 {
1315 for (size_t i = 0; i < n_local_cell_scalar_dofs; i++)
1316 {
1317 vec_phiT_quadT[dim * i + k][iqn] = phiT_quadT[i][iqn] * eT;
1318 }
1319 }
1320 }
1321
1322 // Integrate source term with basis functions
1323 source_vec.head(n_local_cell_vector_dofs) = integrate(f_vec, vec_phiT_quadT, quadT, n_local_cell_vector_dofs);
1324 return source_vec;
1325}
1326
1328{
1329
1330 std::vector<Eigen::VectorXd> LTf(n_cells);
1331 std::vector<Eigen::VectorXd> LTg(n_cells);
1332 std::vector<Eigen::MatrixXd> DT(n_cells);
1333 // std::vector<std::vector<Eigen::MatrixXd>> TTk(n_cells);
1334 std::vector<std::vector<Eigen::MatrixXd>> Tplus_k(n_cells);
1335 std::vector<std::vector<Eigen::MatrixXd>> Tminus_k(n_cells);
1336
1337 std::cout << " Assembling local matrices\n";
1338 boost::timer::cpu_timer assembly_timer;
1339 assembly_timer.start();
1340
1341 // Construct the local matrices using multithreading if use_threads is true
1342 std::function<void(size_t, size_t)> construct_all_local_contributions = [&](size_t start, size_t end) -> void
1343 {
1344 for (size_t iT = start; iT < end; ++iT)
1345 {
1346 // ElementQuad elquad(m_hho, iT, m_K + 1 + m_K + 1 + 2, m_K + m_K + 1);
1347
1348 ElementQuad elquad(m_hho, iT, std::max(3 * m_L, m_L + m_K + 1) + 1, std::max(3 * m_K, m_K + m_K + 1) + 1);
1349 // ElementQuad elquad(m_hho, iT, 15, 15);
1350 Cell *cell = mesh_ptr->cell(iT);
1351
1352 AT[iT] = gradient_term(cell, elquad);
1353 AT[iT] = AT[iT] + stabilisation_term(cell, elquad);
1354 LTf[iT] = (source_term(f_source, cell, elquad).head(n_local_cell_vector_dofs));
1355 LTg[iT] = (source_term(g_source, cell, elquad).head(n_local_cell_vector_dofs));
1356 DT[iT] = (divergence_term(cell, elquad));
1357 size_t n_local_vector_dofs = n_local_cell_vector_dofs + cell->n_faces() * n_local_face_vector_dofs;
1358 // for (size_t k = 0; k < 2 * n_local_vector_dofs; ++k)
1359 for (size_t k = 0; k < n_local_vector_dofs; ++k)
1360 {
1361 // TTk[iT].push_back(T_ijk_plus_T_ikj(cell, elquad, k));
1362 Tplus_k[iT].push_back(th_ijk_plus_th_ikj(cell, elquad, k));
1363 Tminus_k[iT].push_back(th_ijk_minus_th_ikj(cell, elquad, k));
1364 }
1365 }
1366 };
1367
1368 // Running the local constructions in parallel
1370
1371 double exact_rhs_norm = 0.0;
1372 for (size_t iT = 0; iT < n_cells; ++iT)
1373 {
1374 exact_rhs_norm += LTf[iT].squaredNorm() + LTg[iT].squaredNorm();
1375 }
1376 exact_rhs_norm = std::sqrt(exact_rhs_norm);
1377
1378 assembly_timer.stop();
1379 std::cout << " Assembly time = " << double(assembly_timer.elapsed().wall) * 1E-9 << "\n";
1380
1381 size_t n_magnetic_dofs = n_total_cell_vector_dofs + n_total_face_vector_dofs + n_total_cell_scalar_dofs;
1382 size_t n_velocity_dofs = n_total_cell_vector_dofs + n_total_face_vector_dofs + n_total_cell_scalar_dofs;
1383
1384 size_t n_total_dofs = n_magnetic_dofs + n_velocity_dofs;
1385
1386 size_t n_velocity_unknowns = n_total_face_vector_dofs + n_cells; // face dofs plus one scalar dof per cell
1387 size_t n_magnetic_unknowns = n_total_face_vector_dofs + n_cells; // face dofs plus one scalar dof per cell
1388
1389 size_t n_implicit_velocity_dofs = n_total_cell_vector_dofs + n_cells * (n_local_cell_scalar_dofs - 1);
1390 size_t n_implicit_magnetic_dofs = n_total_cell_vector_dofs + n_cells * (n_local_cell_scalar_dofs - 1);
1391
1394
1395 SolutionVector SOL(Eigen::VectorXd::Zero(n_total_dofs), mesh_ptr, m_L, m_K);
1396
1397 std::cout << " Solving scheme using Newton method\n";
1398
1399 boost::timer::cpu_timer solving_timer;
1400 solving_timer.start();
1401
1402 double residual_error = 0.0;
1403 int i = 0;
1404
1405 // to not re-compute the solver if already done
1406 bool solver_analyzed = false;
1407
1408 do
1409 {
1410 std::vector<Eigen::VectorXd> local_RHS(n_cells);
1411
1412 std::vector<Eigen::MatrixXd> local_inv_LTll_LTlg(n_cells);
1413 std::vector<Eigen::MatrixXd> local_LTgl(n_cells);
1414 std::vector<Eigen::MatrixXd> local_LTgg(n_cells);
1415
1416 Eigen::VectorXd inv_LTll_rTl = Eigen::VectorXd::Zero(n_implicit_dofs);
1417
1418 // Construct the local matrices using multithreading if use_threads is true
1419 std::function<void(size_t, size_t)> build_local_mats = [&](size_t start, size_t end) -> void
1420 {
1421 for (size_t iT = start; iT < end; ++iT)
1422 {
1423 Cell *cell = mesh_ptr->cell(iT);
1424 size_t n_local_faces = cell->n_faces();
1425 size_t n_local_bdry_vector_dofs = n_local_faces * n_local_face_vector_dofs;
1426 size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_bdry_vector_dofs;
1427 size_t n_local_unknowns = 2 * (n_local_vector_dofs + n_local_cell_scalar_dofs);
1428 size_t magnetic_offset = n_local_vector_dofs + n_local_cell_scalar_dofs;
1429
1430 Eigen::VectorXd uT = SOL.restr(iT);
1431 Eigen::VectorXd vector_terms = Eigen::VectorXd::Zero(2 * n_local_vector_dofs);
1434
1435 Eigen::MatrixXd A = Eigen::MatrixXd::Zero(n_local_unknowns, n_local_unknowns);
1436
1437 A.block(0, 0, n_local_vector_dofs, n_local_vector_dofs) = visc * AT[iT];
1438 A.block(0, n_local_vector_dofs, n_local_vector_dofs, n_local_cell_scalar_dofs) = DT[iT].transpose();
1439 A.block(n_local_vector_dofs, 0, n_local_cell_scalar_dofs, n_local_vector_dofs) = -DT[iT];
1440
1442 A.block(magnetic_offset, magnetic_offset + n_local_vector_dofs, n_local_vector_dofs, n_local_cell_scalar_dofs) = DT[iT].transpose();
1443 A.block(magnetic_offset + n_local_vector_dofs, magnetic_offset, n_local_cell_scalar_dofs, n_local_vector_dofs) = -DT[iT];
1444
1445 // Eigen::MatrixXd T = convective_term(cell, TTk[iT], vector_terms);
1446 Eigen::MatrixXd T = convective_term(cell, Tplus_k[iT], Tminus_k[iT], vector_terms);
1447 Eigen::MatrixXd BMAT = Eigen::MatrixXd::Zero(n_local_unknowns, n_local_unknowns);
1452
1453 Eigen::MatrixXd local_DG = A + BMAT;
1454 Eigen::MatrixXd GMAT = A + 0.5 * BMAT;
1455
1456 Eigen::VectorXd G = GMAT * uT;
1457 G.segment(0, n_local_cell_vector_dofs) -= LTf[iT];
1458 G.segment(magnetic_offset, n_local_cell_vector_dofs) -= LTg[iT];
1459
1460 local_RHS[iT] = -G;
1461
1462 // local_DG.block(0, 0, n_local_cell_vector_dofs, n_local_cell_vector_dofs) += std::pow(deltaN, -1) * vec_MTT[iT];
1463 // local_DG.block(magnetic_offset, magnetic_offset, n_local_cell_vector_dofs, n_local_cell_vector_dofs) += std::pow(deltaN, -1) * vec_MTT[iT];
1464
1465 // local_DG.block(n_local_vector_dofs, n_local_vector_dofs, n_local_cell_scalar_dofs, n_local_cell_scalar_dofs) += std::pow(deltaN, -1) * MTT[iT];
1466 // local_DG.block(n_local_vector_dofs + magnetic_offset, n_local_vector_dofs + magnetic_offset, n_local_cell_scalar_dofs, n_local_cell_scalar_dofs) += std::pow(deltaN, -1) * MTT[iT];
1467
1470
1472
1473 size_t n_local_implicit_velocity_dofs = n_local_cell_vector_dofs + n_local_cell_scalar_dofs - 1;
1474 size_t n_local_implicit_magnetic_dofs = n_local_cell_vector_dofs + n_local_cell_scalar_dofs - 1;
1475
1477
1478 size_t local_magnetic_offset = n_local_vector_dofs + n_local_cell_scalar_dofs;
1480 size_t local_implicit_magnetic_offset = n_local_cell_vector_dofs + n_local_cell_scalar_dofs - 1;
1481
1482 Eigen::MatrixXd local_LTll = Eigen::MatrixXd::Zero(n_local_implicit_dofs, n_local_implicit_dofs);
1483
1484 // velocity-velocity dofs
1485 local_LTll.block(0, 0, n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(0, 0, n_local_cell_vector_dofs, n_local_cell_vector_dofs);
1486 local_LTll.block(0, n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1) = local_DG.block(0, n_local_vector_dofs + 1, n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1);
1487 local_LTll.block(n_local_cell_vector_dofs, 0, n_local_cell_scalar_dofs - 1, n_local_cell_vector_dofs) = local_DG.block(n_local_vector_dofs + 1, 0, n_local_cell_scalar_dofs - 1, n_local_cell_vector_dofs);
1488
1489 // velocity pressure-pressure dofs (relaxation)
1490 // local_LTll.block(n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1, n_local_cell_scalar_dofs - 1) = local_DG.block(n_local_vector_dofs + 1, n_local_vector_dofs + 1, n_local_cell_scalar_dofs - 1, n_local_cell_scalar_dofs - 1);
1491
1492 // velocity-magnetic dofs
1493 local_LTll.block(0, local_implicit_magnetic_offset, n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(0, local_magnetic_offset, n_local_cell_vector_dofs, n_local_cell_vector_dofs);
1494
1495 // magnetic-velocity dofs
1496 local_LTll.block(local_implicit_magnetic_offset, 0, n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(local_magnetic_offset, 0, n_local_cell_vector_dofs, n_local_cell_vector_dofs);
1497
1498 // magnetic-magnetic dofs
1499 local_LTll.block(local_implicit_magnetic_offset, local_implicit_magnetic_offset, n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(local_magnetic_offset, local_magnetic_offset, n_local_cell_vector_dofs, n_local_cell_vector_dofs);
1500 local_LTll.block(local_implicit_magnetic_offset, local_implicit_magnetic_offset + n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1) = local_DG.block(local_magnetic_offset, local_magnetic_offset + n_local_vector_dofs + 1, n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1);
1501 local_LTll.block(local_implicit_magnetic_offset + n_local_cell_vector_dofs, local_implicit_magnetic_offset, n_local_cell_scalar_dofs - 1, n_local_cell_vector_dofs) = local_DG.block(local_magnetic_offset + n_local_vector_dofs + 1, local_magnetic_offset, n_local_cell_scalar_dofs - 1, n_local_cell_vector_dofs);
1502
1503 // magnetic pressure-pressure dofs (relaxation)
1504 // local_LTll.block(local_implicit_magnetic_offset + n_local_cell_vector_dofs, local_implicit_magnetic_offset + n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1, n_local_cell_scalar_dofs - 1) = local_DG.block(local_magnetic_offset + n_local_vector_dofs + 1, local_magnetic_offset + n_local_vector_dofs + 1, n_local_cell_scalar_dofs - 1, n_local_cell_scalar_dofs - 1);
1505
1506 Eigen::MatrixXd local_LTlg = Eigen::MatrixXd::Zero(n_local_implicit_dofs, n_local_condensed_dofs);
1507
1508 // velocity-velocity dofs
1509 local_LTlg.block(0, 0, n_local_cell_vector_dofs, n_local_bdry_vector_dofs) = local_DG.block(0, n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_bdry_vector_dofs);
1510 local_LTlg.block(n_local_cell_vector_dofs, 0, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs) = local_DG.block(n_local_vector_dofs + 1, n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs);
1511
1512 // velocity-magnetic dofs
1513 local_LTlg.block(0, local_condensed_magnetic_offset, n_local_cell_vector_dofs, n_local_bdry_vector_dofs) = local_DG.block(0, local_magnetic_offset + n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_bdry_vector_dofs);
1514 // local_LTlg.block(0 + n_local_cell_vector_dofs, local_condensed_magnetic_offset, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs) = local_DG.block(0 + n_local_vector_dofs + 1, local_magnetic_offset + n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs);
1515
1516 // magnetic-velocity dofs
1517 local_LTlg.block(local_implicit_magnetic_offset, 0, n_local_cell_vector_dofs, n_local_bdry_vector_dofs) = local_DG.block(local_magnetic_offset, 0 + n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_bdry_vector_dofs);
1518 // local_LTlg.block(local_implicit_magnetic_offset + n_local_cell_vector_dofs, 0, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs) = local_DG.block(local_magnetic_offset + n_local_vector_dofs + 1, 0 + n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs);
1519
1520 // magnetic-magnetic dofs
1521 local_LTlg.block(local_implicit_magnetic_offset, local_condensed_magnetic_offset, n_local_cell_vector_dofs, n_local_bdry_vector_dofs) = local_DG.block(local_magnetic_offset, local_magnetic_offset + n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_bdry_vector_dofs);
1522 local_LTlg.block(local_implicit_magnetic_offset + n_local_cell_vector_dofs, local_condensed_magnetic_offset, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs) = local_DG.block(local_magnetic_offset + n_local_vector_dofs + 1, local_magnetic_offset + n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1, n_local_bdry_vector_dofs);
1523
1524 local_LTgl[iT] = Eigen::MatrixXd::Zero(n_local_condensed_dofs, n_local_implicit_dofs);
1525
1526 // velocity-velocity dofs
1527 local_LTgl[iT].block(0, 0, n_local_bdry_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(n_local_cell_vector_dofs, 0, n_local_bdry_vector_dofs, n_local_cell_vector_dofs);
1528 local_LTgl[iT].block(0, n_local_cell_vector_dofs, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1) = local_DG.block(n_local_cell_vector_dofs, n_local_vector_dofs + 1, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1);
1529
1530 // velocity-magnetic dofs
1531 local_LTgl[iT].block(0, local_implicit_magnetic_offset, n_local_bdry_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(0 + n_local_cell_vector_dofs, local_magnetic_offset, n_local_bdry_vector_dofs, n_local_cell_vector_dofs);
1532 // local_LTgl.block(0, local_implicit_magnetic_offset + n_local_cell_vector_dofs, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1) = local_DG.block(0 + n_local_cell_vector_dofs, local_magnetic_offset + n_local_vector_dofs + 1, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1);
1533
1534 // magnetic-velocity dofs
1535 local_LTgl[iT].block(local_condensed_magnetic_offset, 0, n_local_bdry_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(local_magnetic_offset + n_local_cell_vector_dofs, 0, n_local_bdry_vector_dofs, n_local_cell_vector_dofs);
1536 // local_LTgl.block(local_condensed_magnetic_offset, 0 + n_local_cell_vector_dofs, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1) = local_DG.block(local_magnetic_offset + n_local_cell_vector_dofs, 0 + n_local_vector_dofs + 1, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1);
1537
1538 // magnetic-magnetic dofs
1539 local_LTgl[iT].block(local_condensed_magnetic_offset, local_implicit_magnetic_offset, n_local_bdry_vector_dofs, n_local_cell_vector_dofs) = local_DG.block(local_magnetic_offset + n_local_cell_vector_dofs, local_magnetic_offset, n_local_bdry_vector_dofs, n_local_cell_vector_dofs);
1540 local_LTgl[iT].block(local_condensed_magnetic_offset, local_implicit_magnetic_offset + n_local_cell_vector_dofs, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1) = local_DG.block(local_magnetic_offset + n_local_cell_vector_dofs, local_magnetic_offset + n_local_vector_dofs + 1, n_local_bdry_vector_dofs, n_local_cell_scalar_dofs - 1);
1541
1543
1544 // velocity-velocity dofs
1545 local_LTgg[iT].block(0, 0, n_local_bdry_vector_dofs, n_local_bdry_vector_dofs) = local_DG.block(n_local_cell_vector_dofs, n_local_cell_vector_dofs, n_local_bdry_vector_dofs, n_local_bdry_vector_dofs);
1546 local_LTgg[iT].block(0, n_local_bdry_vector_dofs, n_local_bdry_vector_dofs, 1) = local_DG.block(n_local_cell_vector_dofs, n_local_cell_vector_dofs + n_local_bdry_vector_dofs, n_local_bdry_vector_dofs, 1);
1547 local_LTgg[iT].block(n_local_bdry_vector_dofs, 0, 1, n_local_bdry_vector_dofs) = local_DG.block(n_local_cell_vector_dofs + n_local_bdry_vector_dofs, n_local_cell_vector_dofs, 1, n_local_bdry_vector_dofs);
1548
1549 // velocity pressure-pressure dofs (relaxation)
1550 // local_LTgg[iT](n_local_bdry_vector_dofs, n_local_bdry_vector_dofs) = local_DG(n_local_cell_vector_dofs + n_local_bdry_vector_dofs, n_local_cell_vector_dofs + n_local_bdry_vector_dofs);
1551
1552 // velocity-magnetic dofs
1554 // local_LTgg[iT].block(0, local_condensed_magnetic_offset + n_local_bdry_vector_dofs, n_local_bdry_vector_dofs, 1) = local_DG.block(0 + n_local_cell_vector_dofs, local_magnetic_offset + n_local_cell_vector_dofs + n_local_bdry_vector_dofs, n_local_bdry_vector_dofs, 1);
1555 // local_LTgg[iT].block(0 + n_local_bdry_vector_dofs, local_condensed_magnetic_offset, 1, n_local_bdry_vector_dofs) = local_DG.block(0 + n_local_cell_vector_dofs + n_local_bdry_vector_dofs, local_magnetic_offset + n_local_cell_vector_dofs, 1, n_local_bdry_vector_dofs);
1556
1557 // magnetic-velocity dofs
1559 // local_LTgg[iT].block(local_condensed_magnetic_offset, 0 + n_local_bdry_vector_dofs, n_local_bdry_vector_dofs, 1) = local_DG.block(local_magnetic_offset + n_local_cell_vector_dofs, 0 + n_local_cell_vector_dofs + n_local_bdry_vector_dofs, n_local_bdry_vector_dofs, 1);
1560 // local_LTgg[iT].block(local_condensed_magnetic_offset + n_local_bdry_vector_dofs, 0, 1, n_local_bdry_vector_dofs) = local_DG.block(local_magnetic_offset + n_local_cell_vector_dofs + n_local_bdry_vector_dofs, 0 + n_local_cell_vector_dofs, 1, n_local_bdry_vector_dofs);
1561
1562 // magnetic-magnetic dofs
1566
1567 // magnetic pressure-pressure dofs (relaxation)
1568 // local_LTgg[iT](local_condensed_magnetic_offset + n_local_bdry_vector_dofs, local_condensed_magnetic_offset + n_local_bdry_vector_dofs) = local_DG(local_magnetic_offset + n_local_cell_vector_dofs + n_local_bdry_vector_dofs, local_magnetic_offset + n_local_cell_vector_dofs + n_local_bdry_vector_dofs);
1569
1570 Eigen::PartialPivLU<Eigen::MatrixXd> solver;
1571 solver.compute(local_LTll);
1572
1573 Eigen::VectorXd local_rTl = Eigen::VectorXd::Zero(n_local_implicit_dofs);
1574 local_rTl.head(n_local_cell_vector_dofs) = local_RHS[iT].head(n_local_cell_vector_dofs);
1575 local_rTl.segment(local_implicit_magnetic_offset, n_local_cell_vector_dofs) = local_RHS[iT].segment(local_magnetic_offset, n_local_cell_vector_dofs);
1576
1577 Eigen::VectorXd local_inv_LTll_rTl = solver.solve(local_rTl);
1579
1580 inv_LTll_rTl.segment(iT * n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_inv_LTll_rTl.segment(0, n_local_cell_vector_dofs);
1581 inv_LTll_rTl.segment(n_total_cell_vector_dofs + iT * (n_local_cell_scalar_dofs - 1), n_local_cell_scalar_dofs - 1) = local_inv_LTll_rTl.segment(n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1);
1582
1583 inv_LTll_rTl.segment(n_implicit_velocity_dofs + iT * n_local_cell_vector_dofs, n_local_cell_vector_dofs) = local_inv_LTll_rTl.segment(local_implicit_magnetic_offset, n_local_cell_vector_dofs);
1584 inv_LTll_rTl.segment(n_implicit_velocity_dofs + n_total_cell_vector_dofs + iT * (n_local_cell_scalar_dofs - 1), n_local_cell_scalar_dofs - 1) = local_inv_LTll_rTl.segment(local_implicit_magnetic_offset + n_local_cell_vector_dofs, n_local_cell_scalar_dofs - 1);
1585 }
1586 };
1587
1588 // Running the local constructions in parallel
1590
1591 Eigen::SparseMatrix<double> inv_LTll_LTlg(n_implicit_dofs, n_unknowns + 2);
1592 Eigen::SparseMatrix<double> LTgl(n_unknowns + 2, n_implicit_dofs);
1593 Eigen::SparseMatrix<double> LTgg(n_unknowns + 2, n_unknowns + 2);
1594
1595 std::vector<Eigen::Triplet<double>> triplets_inv_LTll_LTlg;
1596 std::vector<Eigen::Triplet<double>> triplets_LTgl;
1597 std::vector<Eigen::Triplet<double>> triplets_LTgg;
1598
1599 Eigen::VectorXd rTg = Eigen::VectorXd::Zero(n_unknowns + 2);
1600
1601 for (size_t iT = 0; iT < n_cells; ++iT)
1602 {
1603 Cell *cell = mesh_ptr->cell(iT);
1604 size_t n_local_faces = cell->n_faces();
1605 size_t n_local_bdry_vector_dofs = n_local_faces * n_local_face_vector_dofs;
1606
1607 size_t local_velocity_cell_offset = 0;
1608 size_t local_magnetic_cell_offset = n_local_cell_vector_dofs + n_local_cell_scalar_dofs - 1;
1609 size_t global_velocity_cell_offset = iT * n_local_cell_vector_dofs;
1610 size_t global_magnetic_cell_offset = n_implicit_velocity_dofs + iT * n_local_cell_vector_dofs;
1611
1612 size_t local_velocity_pressure_offset = n_local_cell_vector_dofs;
1613 size_t local_magnetic_pressure_offset = n_local_cell_vector_dofs + n_local_cell_scalar_dofs - 1 + n_local_cell_vector_dofs;
1614 size_t global_velocity_pressure_offset = n_total_cell_vector_dofs + iT * (n_local_cell_scalar_dofs - 1);
1615 size_t global_magnetic_pressure_offset = n_implicit_velocity_dofs + n_total_cell_vector_dofs + iT * (n_local_cell_scalar_dofs - 1);
1616
1617 for (size_t iTF = 0; iTF < n_local_faces; iTF++)
1618 {
1619 size_t iF = cell->face(iTF)->global_index();
1620
1621 size_t local_velocity_face_iTF_offset = iTF * n_local_face_vector_dofs;
1622 size_t global_velocity_face_iTF_offset = iF * n_local_face_vector_dofs;
1623 size_t local_magnetic_face_iTF_offset = n_local_bdry_vector_dofs + 1 + iTF * n_local_face_vector_dofs;
1624 size_t global_magnetic_face_iTF_offset = n_velocity_unknowns + iF * n_local_face_vector_dofs;
1625
1626 rTg.segment(global_velocity_face_iTF_offset, n_local_face_vector_dofs) += local_RHS[iT].segment(n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs, n_local_face_vector_dofs);
1627 rTg.segment(global_magnetic_face_iTF_offset, n_local_face_vector_dofs) += local_RHS[iT].segment(n_local_cell_vector_dofs + n_local_bdry_vector_dofs + n_local_cell_scalar_dofs + n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs, n_local_face_vector_dofs);
1628
1629 for (size_t i = 0; i < n_local_face_vector_dofs; i++)
1630 {
1631 for (size_t j = 0; j < n_local_cell_vector_dofs; j++)
1632 {
1637
1642 }
1643 for (size_t j = 0; j < n_local_cell_scalar_dofs - 1; j++)
1644 {
1647
1650 }
1651 for (size_t jTF = 0; jTF < n_local_faces; jTF++)
1652 {
1653 size_t jF = cell->face(jTF)->global_index();
1654
1655 size_t local_velocity_face_jTF_offset = jTF * n_local_face_vector_dofs;
1656 size_t global_velocity_face_jTF_offset = jF * n_local_face_vector_dofs;
1657 size_t local_magnetic_face_jTF_offset = n_local_bdry_vector_dofs + 1 + jTF * n_local_face_vector_dofs;
1658 size_t global_magnetic_face_jTF_offset = n_velocity_unknowns + jF * n_local_face_vector_dofs;
1659 for (size_t j = 0; j < n_local_face_vector_dofs; j++)
1660 {
1665 }
1666 }
1667
1670
1673 }
1674 }
1675
1676 triplets_LTgg.emplace_back(n_total_face_vector_dofs + iT, n_unknowns, std::sqrt(cell->measure()));
1677 triplets_LTgg.emplace_back(n_unknowns, n_total_face_vector_dofs + iT, -std::sqrt(cell->measure()));
1678 triplets_LTgg.emplace_back(n_velocity_unknowns + n_total_face_vector_dofs + iT, n_unknowns + 1, std::sqrt(cell->measure()));
1679 triplets_LTgg.emplace_back(n_unknowns + 1, n_velocity_unknowns + n_total_face_vector_dofs + iT, -std::sqrt(cell->measure()));
1680
1681
1682 // triplets_LTgg.emplace_back(n_unknowns, n_unknowns, std::sqrt(cell->measure()));
1683 // triplets_LTgg.emplace_back(n_unknowns + 1, n_unknowns + 1, std::sqrt(cell->measure()));
1684 // triplets_LTgg.emplace_back(n_unknowns, n_unknowns, 1.0);
1685 // triplets_LTgg.emplace_back(n_unknowns + 1, n_unknowns + 1, 1.0);
1686
1687 // triplets_LTgg.emplace_back(n_total_face_vector_dofs + iT, n_unknowns, 1.0);
1688 // triplets_LTgg.emplace_back(n_unknowns, n_total_face_vector_dofs + iT, -1.0);
1689 // triplets_LTgg.emplace_back(n_velocity_unknowns + n_total_face_vector_dofs + iT, n_unknowns + 1, 1.0);
1690 // triplets_LTgg.emplace_back(n_unknowns + 1, n_velocity_unknowns + n_total_face_vector_dofs + iT, -1.0);
1691
1692 // triplets_LTgg.emplace_back(n_total_face_vector_dofs + iT, n_unknowns, 1.0);
1693
1694 triplets_LTgg.emplace_back(n_total_face_vector_dofs + iT, n_total_face_vector_dofs + iT, std::pow((cell->diam() / double(m_K + 1)), m_K + 2));
1695 triplets_LTgg.emplace_back(n_velocity_unknowns + n_total_face_vector_dofs + iT, n_velocity_unknowns + n_total_face_vector_dofs + iT, std::pow((cell->diam() / double(m_K + 1)), m_K + 2));
1696
1697 // triplets_LTgg.emplace_back(n_total_face_vector_dofs + iT, n_total_face_vector_dofs + iT, local_LTgg[iT](n_local_bdry_vector_dofs, n_local_bdry_vector_dofs));
1698 // triplets_LTgg.emplace_back(n_velocity_unknowns + n_total_face_vector_dofs + iT, n_velocity_unknowns + n_total_face_vector_dofs + iT, local_LTgg[iT](n_local_bdry_vector_dofs + 1 + n_local_bdry_vector_dofs, n_local_bdry_vector_dofs + 1 + n_local_bdry_vector_dofs));
1699 }
1700
1701 inv_LTll_LTlg.setFromTriplets(std::begin(triplets_inv_LTll_LTlg), std::end(triplets_inv_LTll_LTlg));
1702 LTgl.setFromTriplets(std::begin(triplets_LTgl), std::end(triplets_LTgl));
1703 LTgg.setFromTriplets(std::begin(triplets_LTgg), std::end(triplets_LTgg));
1704
1705 for (size_t ibF = 0; ibF < mesh_ptr->n_b_faces(); ++ibF)
1706 {
1707 size_t iF = mesh_ptr->n_i_faces() + ibF;
1708 size_t velocity_offset = iF * n_local_face_vector_dofs;
1709 size_t magnetic_offset = n_velocity_unknowns + iF * n_local_face_vector_dofs;
1710 for (size_t i = 0; i < n_local_face_vector_dofs; ++i)
1711 {
1712 int vel_index = int(velocity_offset + i);
1713 int mag_index = int(magnetic_offset + i);
1714
1715 if ((m_u_bc == 'D') || ((m_u_bc == 'H') && i % dim == 0))
1716 {
1717 LTgg.row(vel_index) *= 0.0;
1718 LTgg.col(vel_index) *= 0.0;
1719 LTgg.coeffRef(vel_index, vel_index) = 1.0;
1720
1721 LTgl.row(vel_index) *= 0.0;
1722 inv_LTll_LTlg.col(vel_index) *= 0.0;
1723 }
1724 if ((m_b_bc == 'D') || ((m_b_bc == 'H') && i % dim == 0))
1725 {
1726 LTgg.row(mag_index) *= 0.0;
1727 LTgg.col(mag_index) *= 0.0;
1728 LTgg.coeffRef(mag_index, mag_index) = 1.0;
1729
1730 LTgl.row(mag_index) *= 0.0;
1731 inv_LTll_LTlg.col(mag_index) *= 0.0;
1732 }
1733 }
1734 }
1735 LTgg.prune(0.0);
1736 inv_LTll_LTlg.prune(0.0);
1737 LTgl.prune(0.0);
1738
1739 LTgg.row(n_unknowns) *= visc;
1740 LTgg.col(n_unknowns) *= visc;
1741 LTgg.row(n_unknowns + 1) *= diff;
1742 LTgg.col(n_unknowns + 1) *= diff;
1743
1744 Eigen::SparseMatrix<double> condensed_SysMat = LTgg - LTgl * inv_LTll_LTlg;
1745
1746 if (!solver_analyzed){
1747 solver.analyzePattern(condensed_SysMat);
1748 solver_analyzed = true;
1749 }
1750 solver.factorize(condensed_SysMat);
1753
1754 Eigen::VectorXd RHS = rTg - LTgl * inv_LTll_rTl;
1755
1757 Eigen::VectorXd condensed_terms = solver.solve(RHS);
1759 // std::cout << "Residual of linear system = " << (condensed_SysMat * condensed_terms - RHS).norm() << "\n";
1760
1761 Eigen::VectorXd implicit_terms = inv_LTll_rTl - inv_LTll_LTlg * condensed_terms;
1762
1763 Eigen::VectorXd deltaU = Eigen::VectorXd::Zero(n_total_dofs);
1764
1765 deltaU.segment(0, n_total_cell_vector_dofs) = implicit_terms.segment(0, n_total_cell_vector_dofs);
1766 deltaU.segment(n_total_cell_vector_dofs, n_total_face_vector_dofs) = condensed_terms.segment(0, n_total_face_vector_dofs);
1767
1768 deltaU.segment(n_velocity_dofs, n_total_cell_vector_dofs) = implicit_terms.segment(n_implicit_velocity_dofs, n_total_cell_vector_dofs);
1769 deltaU.segment(n_velocity_dofs + n_total_cell_vector_dofs, n_total_face_vector_dofs) = condensed_terms.segment(n_velocity_unknowns, n_total_face_vector_dofs);
1770
1771 for (size_t iT = 0; iT < n_cells; ++iT)
1772 {
1773 deltaU(n_total_cell_vector_dofs + n_total_face_vector_dofs + iT * n_local_cell_scalar_dofs) = condensed_terms(n_total_face_vector_dofs + iT);
1774 deltaU(n_velocity_dofs + n_total_cell_vector_dofs + n_total_face_vector_dofs + iT * n_local_cell_scalar_dofs) = condensed_terms(n_velocity_unknowns + n_total_face_vector_dofs + iT);
1775
1776 deltaU.segment(n_total_cell_vector_dofs + n_total_face_vector_dofs + iT * n_local_cell_scalar_dofs + 1, n_local_cell_scalar_dofs - 1) = implicit_terms.segment(n_total_cell_vector_dofs + iT * (n_local_cell_scalar_dofs - 1), n_local_cell_scalar_dofs - 1);
1777 deltaU.segment(n_velocity_dofs + n_total_cell_vector_dofs + n_total_face_vector_dofs + iT * n_local_cell_scalar_dofs + 1, n_local_cell_scalar_dofs - 1) = implicit_terms.segment(n_implicit_velocity_dofs + n_total_cell_vector_dofs + iT * (n_local_cell_scalar_dofs - 1), n_local_cell_scalar_dofs - 1);
1778 }
1779
1780 size_t internal_vector_dofs = n_total_cell_vector_dofs + n_total_face_vector_dofs - mesh_ptr->n_b_faces() * n_local_face_vector_dofs;
1781 size_t bdry_vector_dofs = mesh_ptr->n_b_faces() * n_local_face_vector_dofs;
1782
1783 // set boundary dofs
1784
1785 if (m_u_bc == 'D')
1786 {
1787 deltaU.segment(internal_vector_dofs, bdry_vector_dofs) = Eigen::VectorXd::Zero(bdry_vector_dofs); // Homogeneous Dirchlet BCs
1788 }
1789 else
1790 {
1791 assert(m_u_bc == 'H');
1792 for (size_t ibF = 0; ibF < mesh_ptr->n_b_faces(); ++ibF)
1793 {
1794 size_t offset = internal_vector_dofs + ibF * n_local_face_vector_dofs;
1795
1796 // set normal components to zero
1797 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
1798 {
1799 deltaU(offset + dim * i) = 0.0; // homogenous bcs
1800 }
1801 }
1802 }
1803
1804 if (m_b_bc == 'D')
1805 {
1806 deltaU.segment(n_total_cell_vector_dofs + n_total_face_vector_dofs + n_total_cell_scalar_dofs + internal_vector_dofs, bdry_vector_dofs) = Eigen::VectorXd::Zero(bdry_vector_dofs);
1807 }
1808 else
1809 {
1810 assert(m_b_bc == 'H');
1811 for (size_t ibF = 0; ibF < mesh_ptr->n_b_faces(); ++ibF)
1812 {
1813 size_t offset = n_total_cell_vector_dofs + n_total_face_vector_dofs + n_total_cell_scalar_dofs + internal_vector_dofs + ibF * n_local_face_vector_dofs;
1814
1815 // set normal components to zero
1816 for (size_t i = 0; i < n_local_face_scalar_dofs; i++)
1817 {
1818 deltaU(offset + dim * i) = 0.0; // homogenous bcs
1819 }
1820 }
1821 }
1822
1823 // SOL.set_values(SOL.asVectorXd() + relax * deltaU);
1824 SOL.set_values(SOL.asVectorXd() + deltaU);
1825 Eigen::VectorXd residual = Eigen::VectorXd::Zero(n_total_dofs);
1826
1827 std::function<void(size_t, size_t)> compute_residual = [&](size_t start, size_t end) -> void
1828 {
1829 for (size_t iT = start; iT < end; ++iT)
1830 {
1831 Cell *cell = mesh_ptr->cell(iT);
1832 size_t n_local_faces = cell->n_faces();
1833 size_t n_local_bdry_vector_dofs = n_local_faces * n_local_face_vector_dofs;
1834 size_t n_local_vector_dofs = n_local_cell_vector_dofs + n_local_bdry_vector_dofs;
1835 size_t n_local_unknowns = 2 * (n_local_vector_dofs + n_local_cell_scalar_dofs);
1836 size_t local_magnetic_offset = n_local_vector_dofs + n_local_cell_scalar_dofs;
1837
1838 Eigen::VectorXd uT = SOL.restr(iT);
1839 Eigen::VectorXd vector_terms = Eigen::VectorXd::Zero(2 * n_local_vector_dofs);
1842
1843 Eigen::MatrixXd A = Eigen::MatrixXd::Zero(n_local_unknowns, n_local_unknowns);
1844
1845 A.block(0, 0, n_local_vector_dofs, n_local_vector_dofs) = visc * AT[iT];
1846 A.block(0, n_local_vector_dofs, n_local_vector_dofs, n_local_cell_scalar_dofs) = DT[iT].transpose();
1847 A.block(n_local_vector_dofs, 0, n_local_cell_scalar_dofs, n_local_vector_dofs) = -DT[iT];
1848
1850 A.block(local_magnetic_offset, local_magnetic_offset + n_local_vector_dofs, n_local_vector_dofs, n_local_cell_scalar_dofs) = DT[iT].transpose();
1852
1853 // Eigen::MatrixXd T = convective_term(cell, TTk[iT], vector_terms);
1854 Eigen::MatrixXd T = convective_term(cell, Tplus_k[iT], Tminus_k[iT], vector_terms);
1855 Eigen::MatrixXd BMAT = Eigen::MatrixXd::Zero(n_local_unknowns, n_local_unknowns);
1860
1861 Eigen::MatrixXd GMAT = A + 0.5 * BMAT;
1862
1863 Eigen::VectorXd G = GMAT * uT;
1864 G.segment(0, n_local_cell_vector_dofs) -= LTf[iT];
1865 G.segment(local_magnetic_offset, n_local_cell_vector_dofs) -= LTg[iT];
1866
1867 residual.segment(iT * n_local_cell_vector_dofs, n_local_cell_vector_dofs) = G.head(n_local_cell_vector_dofs);
1868 // residual.segment(n_total_cell_vector_dofs + n_total_face_vector_dofs + iT * n_local_cell_scalar_dofs, n_local_cell_scalar_dofs) = G.segment(n_local_vector_dofs, n_local_cell_scalar_dofs);
1869
1870 size_t global_magnetic_offset = n_total_cell_vector_dofs + n_total_face_vector_dofs + n_total_cell_scalar_dofs;
1871
1872 residual.segment(global_magnetic_offset + iT * n_local_cell_vector_dofs, n_local_cell_vector_dofs) = G.segment(local_magnetic_offset, n_local_cell_vector_dofs);
1873 // residual.segment(global_magnetic_offset + n_total_cell_vector_dofs + n_total_face_vector_dofs + iT * n_local_cell_scalar_dofs, n_local_cell_scalar_dofs) = G.tail(n_local_cell_scalar_dofs);
1874
1875 // for (size_t iTF = 0; iTF < n_local_faces; ++iTF)
1876 // {
1877 // if (cell->face(iTF)->is_boundary())
1878 // {
1879 // continue;
1880 // }
1881 // size_t iF = cell->face(iTF)->global_index();
1882 // residual.segment(n_total_cell_vector_dofs + iF * n_local_face_vector_dofs, n_local_face_vector_dofs) += G.segment(n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs, n_local_face_vector_dofs);
1883 // residual.segment(global_magnetic_offset + n_total_cell_vector_dofs + iF * n_local_face_vector_dofs, n_local_face_vector_dofs) += G.segment(local_magnetic_offset + n_local_cell_vector_dofs + iTF * n_local_face_vector_dofs, n_local_face_vector_dofs);
1884 // }
1885 }
1886 };
1887
1888 // Running the local constructions in parallel
1890
1891 residual_error = residual.norm() / exact_rhs_norm;
1892
1893 ++i;
1894
1895 std::cout << " It. " << i << ":";
1896 std::cout << " Res = " << residual_error << ", Solver error = " << residual_of_linear_system << "\n";
1897
1898 if (i > 500)
1899 {
1900 std::cout << "Failed to converge\n";
1901 // std::cout << "Failed to converge, best residual = " << best_residual << "\n";
1902 // SOL.set_values(best_guess);
1903 break;
1904 }
1905 } while (residual_error > tol);
1906
1907 solving_timer.stop();
1908 std::cout << " Solving time = " << double(solving_timer.elapsed().wall) * 1E-9 << "\n";
1909
1910 return SOL;
1911}
1912
1913#endif
constexpr double tol
Definition BoundaryConditions.cpp:15
Definition elementquad.hpp:55
Definition hybridcore.hpp:174
This structure is a wrapper to allow a given code to select various linear solvers.
Definition linearsolver.hpp:106
Definition HHO_MHD.hpp:209
Class to describe a mesh.
Definition MeshND.hpp:17
boost::multi_array< T, 2 > BasisQuad
type for bases evaluated on quadrature nodes
Definition basis.hpp:61
std::function< T(const VectorRd &)> FType
type for function of point. T is the return type of the function
Definition basis.hpp:55
Eigen::VectorXd integrate(const FType< T > &f, const BasisQuad< T > &B, const QuadratureRule &qr, size_t n_rows=0)
Computes the vector of integrals (f, phi_i)
Definition basis.hpp:2626
Eigen::MatrixXd compute_gram_matrix(const boost::multi_array< VectorRd, 2 > &B1, const boost::multi_array< double, 2 > &B2, const QuadratureRule &qr)
Compute the Gram-like matrix given a family of vector-valued and one of scalar-valued functions by te...
Definition basis.cpp:353
@ Matrix
Definition basis.hpp:67
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:58
size_t L
Definition HHO_DiffAdvecReac.hpp:46
size_t K
Definition HHO_DiffAdvecReac.hpp:46
Eigen::VectorXd lagrange_restr(size_t iT) const
Definition HHO_MHD.hpp:130
SolutionVector solve_with_static_cond(FType< VectorRd > f_source, FType< VectorRd > g_source, double visc, double diff, double tol, LinearSolver< Eigen::SparseMatrix< double > > &solver, bool threading=true)
Definition HHO_MHD.hpp:1327
MHDModel(HybridCore &, size_t, size_t, char, char)
Definition HHO_MHD.hpp:476
void set_values(Eigen::VectorXd values) const
Return the values as an Eigen vector.
Definition HHO_MHD.hpp:62
const size_t get_face_deg() const
Return the face degree.
Definition HHO_MHD.hpp:68
std::vector< double > compute_errors(SolutionVector interpolant, SolutionVector discrete, double visc, double diff)
Definition HHO_MHD.hpp:269
Eigen::VectorXd pressure_restr(size_t iT) const
Definition HHO_MHD.hpp:106
SolutionVector operator-(const SolutionVector &b) const
Overloads the subtraction: subtracts the coefficients.
Definition HHO_MHD.hpp:157
Eigen::VectorXd magnetic_restr(size_t iT) const
Definition HHO_MHD.hpp:113
Eigen::VectorXd lagrange_values() const
Definition HHO_MHD.hpp:86
Eigen::VectorXd pressure_values() const
Definition HHO_MHD.hpp:76
SolutionVector(Eigen::VectorXd values, const Mesh *mesh_ptr, const size_t cell_deg, const size_t face_deg)
Definition HHO_MHD.hpp:48
Eigen::VectorXd velocity_values() const
Definition HHO_MHD.hpp:70
Eigen::VectorXd magnetic_values() const
Definition HHO_MHD.hpp:81
SolutionVector operator+(const SolutionVector &b) const
Overloads the addition: adds the coefficients.
Definition HHO_MHD.hpp:150
const size_t get_cell_deg() const
Return the cell degree.
Definition HHO_MHD.hpp:65
double operator()(size_t index) const
Overloads the (): returns the corresponding coefficient.
Definition HHO_MHD.hpp:164
Eigen::VectorXd velocity_restr(size_t iT) const
Definition HHO_MHD.hpp:91
Eigen::VectorXd asVectorXd() const
Return the values as an Eigen vector.
Definition HHO_MHD.hpp:60
SolutionVector global_interpolant(FType< VectorRd > velocity, FType< double > pressure, FType< VectorRd > magnetic)
Definition HHO_MHD.hpp:447
Eigen::VectorXd restr(size_t iT) const
Definition HHO_MHD.hpp:137
const size_t DimPoly< Face >(const int m)
Compute the size of the basis of 2-variate polynomials up to degree m.
Definition hybridcore.hpp:68
const Mesh * get_mesh() const
Returns a pointer to the mesh.
Definition hybridcore.hpp:198
const size_t DimPoly< Cell >(const int m)
Compute the size of the basis of 3-variate polynomials up to degree m.
Definition hybridcore.hpp:61
std::size_t n_faces() const
number of faces in the mesh.
Definition MeshND.hpp:59
std::size_t n_cells() const
number of cells in the mesh.
Definition MeshND.hpp:60
Cell< dimension > * cell(std::size_t index) const
get a constant pointer to a cell using its global index
Definition MeshND.hpp:197
std::size_t n_b_faces() const
number of boundary faces in the mesh.
Definition MeshND.hpp:77
std::size_t n_i_faces() const
number of internal faces in the mesh.
Definition MeshND.hpp:82
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:61
Definition ddr-magnetostatics.hpp:41
MeshND::VectorRd< 2 > VectorRd
Definition Mesh2D.hpp:14
Definition HHO_MHD.hpp:42