27#ifndef _HHO_DIFFUSION_HPP
28#define _HHO_DIFFUSION_HPP
34#include <boost/timer/timer.hpp>
37#include <Eigen/Sparse>
79 std::string solver_type,
81 std::ostream & output = std::cout
99 return double(_assembly_time) * pow(10, -9);
104 return double(_solving_time) * pow(10, -9);
109 return _solving_error;
114 return double(_itime[idx]) * pow(10, -9);
134 Eigen::MatrixXd diffusion_operator(
HybridCore& hho,
const size_t iT,
const ElementQuad& elquad)
const;
153 const std::string solver_type;
154 const bool m_use_threads;
155 std::ostream & m_output;
158 const size_t m_nlocal_cell_dofs;
159 const size_t m_nlocal_edge_dofs;
160 const size_t m_nhighorder_dofs;
161 const size_t m_ntotal_cell_dofs;
162 const size_t m_ntotal_edge_dofs;
163 const size_t m_ndir_edge_dofs;
164 const size_t m_nnondir_edge_dofs;
165 const size_t m_ntotal_dofs;
168 std::vector<Eigen::MatrixXd> aT;
170 Eigen::SparseMatrix<double> GlobMat;
171 Eigen::SparseMatrix<double> SysMat;
174 Eigen::SparseMatrix<double> ScBeMat;
176 Eigen::VectorXd GlobRHS;
177 Eigen::VectorXd ScRHS;
180 size_t _assembly_time;
181 size_t _solving_time;
182 double _solving_error;
183 mutable std::vector<size_t> _itime = std::vector<size_t>(10, 0);
187 HHO_Diffusion::HHO_Diffusion(
HybridCore& hho,
size_t K,
int L,
CellFType<MatrixRd> kappa,
CellFType<double> source,
BoundaryConditions BC,
FType<double> exact_solution,
CellFType<VectorRd> grad_exact_solution, std::string solver_type,
bool use_threads, std::ostream & output)
191 m_Ldeg(std::max(
L,0)),
195 exact_solution(exact_solution),
196 grad_exact_solution(grad_exact_solution),
197 solver_type(solver_type),
200 m_nlocal_cell_dofs(
DimPoly<Cell>(m_Ldeg)),
201 m_nlocal_edge_dofs(
DimPoly<Edge>(m_K)),
202 m_nhighorder_dofs(
DimPoly<Cell>(m_K+1)),
203 m_ntotal_cell_dofs(m_nlocal_cell_dofs * m_hho.get_mesh()->n_cells()),
204 m_ntotal_edge_dofs(m_nlocal_edge_dofs * m_hho.get_mesh()->n_edges()),
205 m_ndir_edge_dofs(m_nlocal_edge_dofs * m_BC.n_dir_edges()),
206 m_nnondir_edge_dofs(m_nlocal_edge_dofs * m_hho.get_mesh()->n_edges() - m_ndir_edge_dofs),
207 m_ntotal_dofs(m_ntotal_cell_dofs + m_ntotal_edge_dofs),
208 GlobRHS(
Eigen::VectorXd::Zero(m_ntotal_edge_dofs)),
209 ScRHS(
Eigen::VectorXd::Zero(m_ntotal_cell_dofs))
211 m_output <<
"[HHO_Diffusion] Initializing" << std::endl;
212 GlobMat.resize(m_ntotal_edge_dofs, m_ntotal_edge_dofs);
213 ScBeMat.resize(m_ntotal_cell_dofs, m_ntotal_edge_dofs);
219 boost::timer::cpu_timer timer;
226 std::vector<Eigen::Triplet<double>> triplets_GlobMat;
227 std::vector<Eigen::Triplet<double>> triplets_ScBe;
231 std::vector<std::vector<Eigen::Triplet<double>>> cell_triplets_GlobMat;
232 std::vector<std::vector<Eigen::Triplet<double>>> cell_triplets_ScBe;
233 std::vector<Eigen::VectorXd> cell_source(mesh->
n_cells());
234 cell_triplets_GlobMat.resize(mesh->
n_cells());
235 cell_triplets_ScBe.resize(mesh->
n_cells());
236 size_t size_triplets_GlobMat = 0;
237 size_t size_triplets_ScBe = 0;
242 std::function<void(
size_t,
size_t)> construct_all_local_contributions
243 = [&](
size_t start,
size_t end)->
void
245 for (
size_t iT = start; iT <
end; iT++) {
246 Cell* iCell = mesh->
cell(iT);
249 size_t nlocal_edges = iCell->
n_edges();
250 size_t edge_dofs = nlocal_edges * m_nlocal_edge_dofs;
253 size_t doeT = m_Ldeg + m_K + 1;
254 size_t doeF = 2*m_K + 1;
257 aT[iT] = diffusion_operator(hho, iT, elquad);
258 Eigen::VectorXd bT = load_operator(hho, iT, elquad);
261 Eigen::MatrixXd MatF = Eigen::MatrixXd::Zero(edge_dofs,edge_dofs);
262 cell_source[iT] = Eigen::VectorXd::Zero(edge_dofs);
268 Eigen::MatrixXd ATT = aT[iT].topLeftCorner(m_nlocal_cell_dofs, m_nlocal_cell_dofs);
269 Eigen::MatrixXd ATF = aT[iT].topRightCorner(m_nlocal_cell_dofs, edge_dofs);
270 Eigen::MatrixXd AFF = aT[iT].bottomRightCorner(edge_dofs, edge_dofs);
272 Eigen::PartialPivLU<Eigen::MatrixXd> invATT;
275 Eigen::MatrixXd invATT_ATF = invATT.solve(ATF);
276 Eigen::VectorXd invATT_bTcell = invATT.solve(bT.head(m_nlocal_cell_dofs));
277 MatF = AFF - ATF.transpose() * invATT_ATF;
279 cell_source[iT] = bT.tail(edge_dofs) - ATF.transpose() * invATT_bTcell;
282 ScRHS.segment(iT * m_nlocal_cell_dofs, m_nlocal_cell_dofs) = invATT_bTcell;
283 for (
size_t i = 0;
i < m_nlocal_cell_dofs;
i++) {
284 for (
size_t jlF = 0; jlF < nlocal_edges; jlF++) {
285 const size_t jF = iCell->edge(jlF)->global_index();
286 for (
size_t jk = 0; jk < m_nlocal_edge_dofs; jk++) {
287 const size_t jLocal = jlF * m_nlocal_edge_dofs + jk;
288 const size_t jGlobal = jF * m_nlocal_edge_dofs + jk;
289 cell_triplets_ScBe[iT].emplace_back(iT * m_nlocal_cell_dofs +
i, jGlobal, invATT_ATF(
i, jLocal));
293 size_triplets_ScBe += cell_triplets_ScBe[iT].size();
299 Eigen::MatrixXd red_matT = Eigen::MatrixXd::Zero(1+nlocal_edges,nlocal_edges);
303 for (
size_t ilF = 0; ilF < nlocal_edges; ilF++){
304 VectorRd xF = iCell->edge(ilF)->center_mass();
305 size_t iF = iCell->edge(ilF)->global_index();
307 red_matT(0,ilF) *= phiF_cst / phiT_cst;
309 red_matT.bottomRightCorner(nlocal_edges,nlocal_edges) = Eigen::MatrixXd::Identity(nlocal_edges,nlocal_edges);
311 cell_source[iT] = red_matT.transpose() * bT;
312 MatF = red_matT.transpose() * aT[iT] * red_matT;
315 for (
size_t jlF = 0; jlF < nlocal_edges; jlF++) {
316 const size_t jF = iCell->edge(jlF)->global_index();
317 const size_t jGlobal = jF * m_nlocal_edge_dofs;
318 cell_triplets_ScBe[iT].emplace_back(iT, jGlobal, red_matT(0,jlF));
320 size_triplets_ScBe += cell_triplets_ScBe[iT].size();
325 for (
size_t ilF = 0; ilF < nlocal_edges; ilF++) {
326 const size_t iF = iCell->edge(ilF)->global_index();
327 for (
size_t ik = 0; ik < m_nlocal_edge_dofs; ik++) {
328 const size_t iLocal = ilF * m_nlocal_edge_dofs + ik;
329 const size_t iGlobal = iF * m_nlocal_edge_dofs + ik;
330 for (
size_t jlF = 0; jlF < nlocal_edges; jlF++) {
331 const size_t jF = iCell->edge(jlF)->global_index();
332 for (
size_t jk = 0; jk < m_nlocal_edge_dofs; jk++) {
333 const size_t jLocal = jlF * m_nlocal_edge_dofs + jk;
334 const size_t jGlobal = jF * m_nlocal_edge_dofs + jk;
335 cell_triplets_GlobMat[iT].emplace_back(iGlobal, jGlobal, MatF(iLocal, jLocal));
340 size_triplets_GlobMat += cell_triplets_GlobMat[iT].size();
349 triplets_ScBe.reserve(size_triplets_ScBe);
350 triplets_GlobMat.reserve(size_triplets_GlobMat);
351 for (
size_t iT = 0; iT < mesh->
n_cells(); iT++){
352 for (
size_t i = 0;
i < cell_triplets_ScBe[iT].size();
i++){
353 triplets_ScBe.push_back(cell_triplets_ScBe[iT][
i]);
355 for (
size_t i = 0;
i < cell_triplets_GlobMat[iT].size();
i++){
356 triplets_GlobMat.push_back(cell_triplets_GlobMat[iT][
i]);
358 Cell&
T = *mesh->
cell(iT);
359 for (
size_t ilF = 0; ilF <
T.n_edges(); ilF++) {
360 const size_t iF =
T.
edge(ilF)->global_index();
361 for (
size_t ik = 0; ik < m_nlocal_edge_dofs; ik++) {
362 const size_t iLocal = ilF * m_nlocal_edge_dofs + ik;
363 const size_t iGlobal = iF * m_nlocal_edge_dofs + ik;
364 GlobRHS(iGlobal) += cell_source[iT](iLocal);
369 if (m_BC.
name()==
"Neumann"){
371 triplets_GlobMat.erase(std::remove_if(std::begin(triplets_GlobMat), std::end(triplets_GlobMat),
372 [](
const auto& x) {
return (x.row() == 0); }), std::end(triplets_GlobMat));
373 triplets_GlobMat.emplace_back(0, 0, 1);
378 GlobMat.setFromTriplets(std::begin(triplets_GlobMat), std::end(triplets_GlobMat));
379 ScBeMat.setFromTriplets(std::begin(triplets_ScBe), std::end(triplets_ScBe));
383 _assembly_time = timer.elapsed().wall;
390 boost::timer::cpu_timer timer;
399 size_t n_unknowns = 0;
400 size_t n_fixed_dofs = 0;
402 Eigen::VectorXd UDir;
404 if (m_BC.
name() !=
"Neumann"){
406 n_unknowns = m_nnondir_edge_dofs;
407 n_fixed_dofs = m_ndir_edge_dofs;
408 SysMat = GlobMat.topLeftCorner(n_unknowns, n_unknowns);
411 UDir = Eigen::VectorXd::Zero(n_fixed_dofs);
413 size_t n_nondir_edges = mesh->
n_edges() - n_dir_edges;
414 for (
size_t idF = 0; idF < n_dir_edges; idF++){
415 Edge* edge = mesh->
edge(n_nondir_edges + idF);
419 UDir.segment(idF * m_nlocal_edge_dofs, m_nlocal_edge_dofs) = l2_projection<HybridCore::PolyEdgeBasisType>(exact_solution, hho.
EdgeBasis(iF), quadF, phiF_quadF);
422 B = GlobRHS.segment(0, n_unknowns) - GlobMat.topRightCorner(n_unknowns, n_fixed_dofs) * UDir;
426 n_unknowns = m_ntotal_edge_dofs;
429 UDir = Eigen::VectorXd::Zero(n_fixed_dofs);
434 Eigen::VectorXd xF = Eigen::VectorXd::Zero(n_unknowns);
442 Eigen::BiCGSTAB<Eigen::SparseMatrix<double> > solver;
443 solver.compute(SysMat);
444 xF = solver.solve(B);
445 std::cout <<
" [solver] #iterations: " << solver.iterations() <<
", estimated error: " << solver.error() << std::endl;
447 _solving_error = (SysMat * xF - B).norm();
449 Eigen::VectorXd Xh = Eigen::VectorXd::Zero(m_ntotal_dofs);
450 Xh.tail(n_fixed_dofs) = UDir;
451 Xh.segment(m_ntotal_cell_dofs, n_unknowns) = xF;
453 Xh.head(m_ntotal_cell_dofs) = ScRHS - ScBeMat * Xh.tail(m_ntotal_edge_dofs);
455 Xh.head(m_ntotal_cell_dofs) = ScBeMat * Xh.tail(m_ntotal_edge_dofs);
459 if (m_BC.
name()==
"Neumann"){
462 double total_measure = 0;
463 for (
size_t iT = 0; iT < mesh->
n_cells(); iT++){
464 Cell&
T = *mesh->
cell(iT);
470 for (
size_t iqn = 0; iqn < quadT.size(); iqn++){
471 average += quadT[iqn].w * Xh(iT * m_nlocal_cell_dofs +
i) * basisT.
function(
i, quadT[iqn].vector());
475 double average_exact_sol = 0.0;
479 average_exact_sol += qT.w * exact_solution(qT.vector());
482 average_exact_sol /= total_measure;
487 std::function<double(
VectorRd)> AveDiff = [&average_exact_sol,&average](
VectorRd x)->
double
488 {
return average_exact_sol - average;};
494 _solving_time = timer.elapsed().user + timer.elapsed().system;
503 Eigen::MatrixXd HHO_Diffusion::diffusion_operator(
HybridCore &hho,
const size_t iT,
const ElementQuad &elquad)
const {
505 boost::timer::cpu_timer timeint;
508 Cell* cell = mesh->cell(iT);
509 const size_t nedgesT = cell->n_edges();
512 size_t local_dofs = m_nlocal_cell_dofs + nedgesT * m_nlocal_edge_dofs;
519 return kappa(cell->center_mass(), cell);
528 boost::multi_array<double, 2> phiT_quadT = elquad.
get_phiT_quadT();
529 boost::multi_array<VectorRd, 2> dphiT_quadT = elquad.
get_dphiT_quadT();
532 boost::multi_array<VectorRd, 2> kappaT_dphiT_quadT( boost::extents[dphiT_quadT.shape()[0]][quadT.size()] );
533 for (
size_t i = 0;
i < dphiT_quadT.shape()[0];
i++){
534 for (
size_t iqn = 0; iqn < quadT.size(); iqn++){
535 kappaT_dphiT_quadT[
i][iqn] = kappaT(quadT[iqn].vector()) * dphiT_quadT[
i][iqn];
542 std::vector<Eigen::MatrixXd> MFF(nedgesT, Eigen::MatrixXd::Zero(m_nlocal_edge_dofs, m_nlocal_edge_dofs));
543 std::vector<Eigen::MatrixXd> MFT(nedgesT, Eigen::MatrixXd::Zero(m_nlocal_edge_dofs, m_nhighorder_dofs));
548 Eigen::MatrixXd StiffT =
compute_gram_matrix(kappaT_dphiT_quadT, dphiT_quadT, quadT,
"sym");
550 _itime[1] += timeint.elapsed().user + timeint.elapsed().system;
554 Eigen::MatrixXd MTT =
compute_gram_matrix(phiT_quadT, phiT_quadT, quadT, m_nlocal_cell_dofs, m_nhighorder_dofs,
"sym");
556 _itime[2] += timeint.elapsed().user + timeint.elapsed().system;
560 Eigen::VectorXd LT = (MTT.row(0)).transpose();
561 Eigen::MatrixXd LTtLT = LT * (LT.transpose());
564 Eigen::MatrixXd BP = Eigen::MatrixXd::Zero(m_nhighorder_dofs, local_dofs);
565 BP.topLeftCorner(m_nhighorder_dofs, m_nlocal_cell_dofs) = StiffT.topLeftCorner(m_nhighorder_dofs, m_nlocal_cell_dofs);
568 for (
size_t ilF = 0; ilF < nedgesT; ilF++) {
570 const size_t offset_F = m_nlocal_cell_dofs + ilF * m_nlocal_edge_dofs;
571 const auto& nTF = cell->edge_normal(ilF);
575 size_t nbqF = quadF.size();
577 boost::multi_array<double, 2> phiT_quadF = elquad.
get_phiT_quadF(ilF);
578 boost::multi_array<double, 2> phiF_quadF = elquad.
get_phiF_quadF(ilF);
579 boost::multi_array<VectorRd, 2> dphiT_quadF = elquad.
get_dphiT_quadF(ilF);
584 std::vector<VectorRd> quadF_kappaT_nTF(nbqF, VectorRd::Zero());
585 for (
size_t iqn = 0; iqn < nbqF; iqn++){
586 quadF_kappaT_nTF[iqn] = quadF[iqn].w * kappaT(quadF[iqn].vector())*nTF;
588 for (
size_t i = 1;
i < m_nhighorder_dofs;
i++) {
590 for (
size_t iqn = 0; iqn < nbqF; iqn++){
592 double contrib_i = dphiT_quadF[
i][iqn].dot(quadF_kappaT_nTF[iqn]);
594 for (
size_t j = 0;
j < m_nlocal_cell_dofs;
j++) {
596 BP(
i,
j) -= contrib_i * phiT_quadF[
j][iqn];
598 for (
size_t j = 0;
j < m_nlocal_edge_dofs;
j++) {
600 BP(
i, offset_F +
j) += contrib_i * phiF_quadF[
j][iqn];
607 MFT[ilF] =
compute_gram_matrix(phiF_quadF, phiT_quadF, quadF, m_nlocal_edge_dofs, m_nhighorder_dofs,
"nonsym");
612 double scalT = StiffT.trace() / std::pow(LT.norm(), 2);
613 BP.topLeftCorner(m_nhighorder_dofs, m_nlocal_cell_dofs) += scalT * LTtLT.topLeftCorner(m_nhighorder_dofs, m_nlocal_cell_dofs);
614 Eigen::MatrixXd PT = ((StiffT+scalT*LTtLT).ldlt()).
solve(BP);
616 _itime[3] += timeint.elapsed().user + timeint.elapsed().system;
620 Eigen::MatrixXd ATF = PT.transpose() * StiffT * PT;
625 Eigen::MatrixXd STF = Eigen::MatrixXd::Zero(local_dofs, local_dofs);
628 Eigen::MatrixXd MTT_LL = MTT.topLeftCorner(m_nlocal_cell_dofs, m_nlocal_cell_dofs);
629 Eigen::MatrixXd deltaTL = MTT_LL.ldlt().solve( MTT * PT );
630 deltaTL.topLeftCorner(m_nlocal_cell_dofs, m_nlocal_cell_dofs) -= Eigen::MatrixXd::Identity(m_nlocal_cell_dofs, m_nlocal_cell_dofs);
632 for (
size_t ilF = 0; ilF < nedgesT; ilF++) {
634 double dTF = cell->edge(ilF)->diam();
637 VectorRd xF = cell->edge(ilF)->center_mass();
641 const VectorRd &nTF = cell->edge_normal(ilF);
642 const double kappa_TF = (kappa(xF, cell) * nTF).dot(nTF);
645 Eigen::MatrixXd MFFinv = MFF[ilF].inverse();
646 Eigen::MatrixXd deltaTFK = MFFinv * MFT[ilF] * PT;
647 deltaTFK.block(0, m_nlocal_cell_dofs + ilF * m_nlocal_edge_dofs, m_nlocal_edge_dofs, m_nlocal_edge_dofs) -=
648 Eigen::MatrixXd::Identity(m_nlocal_edge_dofs, m_nlocal_edge_dofs);
651 Eigen::MatrixXd deltaTFK_minus_deltaTL = deltaTFK - MFFinv * MFT[ilF].topLeftCorner(m_nlocal_edge_dofs, m_nlocal_cell_dofs) * deltaTL;
653 STF += (kappa_TF / dTF) * deltaTFK_minus_deltaTL.transpose() * MFF[ilF] * deltaTFK_minus_deltaTL;
656 _itime[5] += timeint.elapsed().user + timeint.elapsed().system;
659 ATF += mesh->dim() * STF;
670 Eigen::VectorXd HHO_Diffusion::load_operator(
HybridCore &hho,
const size_t iT,
const ElementQuad &elquad)
const {
673 Cell* cell = mesh->cell(iT);
674 size_t cell_edge_dofs = m_nlocal_cell_dofs + cell->n_edges()*m_nlocal_edge_dofs;
675 Eigen::VectorXd b = Eigen::VectorXd::Zero(cell_edge_dofs);
679 size_t nbq = quadT.size();
680 boost::multi_array<double, 2> phiT_quadT = elquad.
get_phiT_quadT();
683 Eigen::ArrayXd weight_source_quad = Eigen::ArrayXd::Zero(nbq);
684 for (
size_t iqn = 0; iqn < nbq; iqn++){
685 weight_source_quad(iqn) = quadT[iqn].w * source(quadT[iqn].vector(), cell);
688 for (
size_t i=0;
i < m_nlocal_cell_dofs;
i++){
689 for (
size_t iqn = 0; iqn < quadT.size(); iqn++){
690 b(
i) += weight_source_quad[iqn] * phiT_quadT[
i][iqn];
695 if (cell->is_boundary()){
697 for (
size_t ilF = 0; ilF < cell->n_edges(); ilF++) {
698 Edge* F = cell->edge(ilF);
699 if (m_BC.
type(*F)==
"neu"){
700 const size_t iF = F->global_index();
702 if (F->is_boundary()){
704 const size_t offset_F = m_nlocal_cell_dofs + ilF * m_nlocal_edge_dofs;
706 const auto& nTF = cell->edge_normal(ilF);
709 for (
size_t i = 0;
i < m_nlocal_edge_dofs;
i++){
712 return nTF.dot(kappa(p,cell) * grad_exact_solution(p,cell)) * basisF.function(
i,p);
715 b(offset_F +
i) += qF.w * Kgrad_n(qF.vector());
730 for (
size_t iT = 0; iT < mesh->n_cells(); iT++) {
731 Eigen::VectorXd XTF = Xh.
restr(iT);
732 value += XTF.transpose() * aT[iT] * XTF;
std::function< T(const VectorRd &, const Cell *)> CellFType
type for function of a point, that could be discontinuous between cells. T is the type of value of th...
Definition TestCase.hpp:42
The BoundaryConditions class provides definition of boundary conditions.
Definition BoundaryConditions.hpp:45
const std::string type(const Edge &edge) const
Test the boundary condition of an edge.
Definition BoundaryConditions.cpp:41
const std::string name() const
Returns the complete name of the boundary condition.
Definition BoundaryConditions.hpp:80
const size_t n_dir_edges() const
Returns the number of Dirichlet edges.
Definition BoundaryConditions.hpp:70
Definition elementquad.hpp:49
The HHO_Diffusion class provides tools to implement the HHO method for the diffusion problem.
Definition HHO_Diffusion.hpp:65
Definition hybridcore.hpp:179
Definition hybridcore.hpp:82
Compute max and min eigenvalues of all matrices for i
Definition compute_eigs.m:5
end
Definition compute_eigs.m:12
FunctionValue function(size_t i, const VectorRd &x) const
Evaluate the i-th function at point x.
Definition basis.hpp:351
Eigen::Vector2d VectorRd
Definition basis.hpp:55
size_t dimension() const
Dimension of the family. This is actually the number of functions in the family, not necessarily line...
Definition basis.hpp:345
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:225
Eigen::Matrix2d MatrixRd
Definition basis.hpp:54
std::function< T(const VectorRd &)> FType
type for function of point. T is the type of value of the function
Definition basis.hpp:62
static boost::multi_array< typename detail::basis_evaluation_traits< BasisType, BasisFunction >::ReturnValue, 2 > compute(const BasisType &basis, const QuadratureRule &quad)
Generic basis evaluation.
Definition basis.hpp:2288
static void parallel_for(unsigned nb_elements, std::function< void(size_t start, size_t end)> functor, bool use_threads=true, unsigned nb_threads_max=1e9)
Generic function to execute threaded processes.
Definition parallel_for.hpp:42
bool use_threads
Definition HHO_DiffAdvecReac.hpp:47
size_t L
Definition HHO_DiffAdvecReac.hpp:46
size_t K
Definition HHO_DiffAdvecReac.hpp:46
const size_t get_ntotal_edge_dofs()
Total number of edge DOFs over the entire mesh.
Definition HHO_Diffusion.hpp:126
double get_solving_time() const
cpu time to solve the scheme
Definition HHO_Diffusion.hpp:102
const size_t get_nlocal_edge_dofs()
Number of DOFs on each edge.
Definition HHO_Diffusion.hpp:120
const size_t get_ndir_edge_dofs()
Total number of edge DOFs for Dirichlet edges.
Definition HHO_Diffusion.hpp:128
const size_t get_nlocal_cell_dofs()
Number of DOFs in each cell.
Definition HHO_Diffusion.hpp:118
double get_itime(size_t idx) const
various intermediate assembly times
Definition HHO_Diffusion.hpp:112
double get_assembly_time() const
cpu time to assemble the scheme
Definition HHO_Diffusion.hpp:97
UVector solve(HybridCore &hho)
Definition HHO_Diffusion.hpp:387
Eigen::SparseMatrix< double > get_SysMat()
Return the (statically condensed) matrix system.
Definition HHO_Diffusion.hpp:92
const size_t get_ntotal_dofs()
Total number of degrees of freedom over the entire mesh.
Definition HHO_Diffusion.hpp:130
double EnergyNorm(HybridCore &hho, const UVector Xh)
Discrete energy norm (associated to the diffusion operator) of an hybrid function.
Definition HHO_Diffusion.hpp:726
void assemble(HybridCore &hho)
Assemble and solve the scheme.
Definition HHO_Diffusion.hpp:217
const size_t get_ntotal_cell_dofs()
Total number of cell DOFs over the entire mesh.
Definition HHO_Diffusion.hpp:124
double get_solving_error() const
residual after solving the scheme
Definition HHO_Diffusion.hpp:107
const size_t get_nhighorder_dofs()
Number of DOFs per cell for high-order (K+1) polynomials.
Definition HHO_Diffusion.hpp:122
HHO_Diffusion(HybridCore &hho, size_t K, int L, CellFType< MatrixRd > kappa, CellFType< double > source, BoundaryConditions BC, FType< double > exact_solution, CellFType< VectorRd > grad_exact_solution, std::string solver_type, bool use_threads, std::ostream &output=std::cout)
Constructor of the class.
Definition HHO_Diffusion.hpp:187
Eigen::VectorXd compute_weights(size_t iT) const
Computes the weights to get cell values from edge values when l=-1.
Definition hybridcore.cpp:138
Eigen::VectorXd restr(size_t iT) const
Extract the restriction of the unknowns corresponding to cell iT and its edges.
Definition hybridcore.cpp:29
const QuadratureRule & get_quadT() const
Returns quadrature rules in cell.
Definition elementquad.hpp:61
const boost::multi_array< VectorRd, 2 > & get_dphiT_quadT() const
Returns values of gradients of cell basis functions at cell quadrature nodes.
Definition elementquad.hpp:79
const boost::multi_array< VectorRd, 2 > & get_dphiT_quadF(size_t ilF) const
Returns values of gradients of cell basis functions at edge quadrature nodes, for edge with local num...
Definition elementquad.hpp:97
const QuadratureRule & get_quadF(size_t ilF) const
Returns quadrature rules on edge with local number ilF.
Definition elementquad.hpp:67
const PolyEdgeBasisType & EdgeBasis(size_t iE) const
Return edge basis for edge with global index iE.
Definition hybridcore.hpp:218
UVector interpolate(const ContinuousFunction &f, const int deg_cell, const size_t deg_edge, size_t doe) const
Compute the interpolant in the discrete space of a continuous function.
Definition hybridcore.hpp:290
Eigen::VectorXd & asVectorXd() const
Return the values as an Eigen vector.
Definition hybridcore.hpp:92
const boost::multi_array< double, 2 > & get_phiT_quadF(size_t ilF) const
Returns values of cell basis functions at edge quadrature nodes, for edge with local number ilF.
Definition elementquad.hpp:85
const boost::multi_array< double, 2 > get_phiF_quadF(size_t ilF) const
Returns values of edge basis functions at edge quadrature nodes, for edge with local number ilF.
Definition elementquad.hpp:91
size_t const DimPoly(int m)
const boost::multi_array< double, 2 > & get_phiT_quadT() const
Returns values of cell basis functions at cell quadrature nodes.
Definition elementquad.hpp:73
const Mesh * get_mesh() const
Returns a pointer to the mesh.
Definition hybridcore.hpp:201
const PolyCellBasisType & CellBasis(size_t iT) const
Return cell basis for element with global index iT.
Definition hybridcore.hpp:210
std::vector< Cell * > get_cells() const
lists the cells in the mesh.
Definition Mesh2D.hpp:78
double measure() const
Return the Lebesgue measure of the Polytope.
Definition Polytope2D.hpp:76
Polytope< DIMENSION > Cell
Definition Polytope2D.hpp:151
Cell * cell(std::size_t index) const
get a constant pointer to a cell using its global index
Definition Mesh2D.hpp:178
size_t n_edges() const
Return the number of edges of the Polytope.
Definition Polytope2D.hpp:88
Edge * edge(std::size_t index) const
get a constant pointer to a edge using its global index
Definition Mesh2D.hpp:168
Polytope< 1 > Edge
A Face is a Polytope with object_dim = DIMENSION - 1.
Definition Polytope2D.hpp:147
Polytope< 1 > * edge(const size_t i) const
Return the i-th edge of the Polytope.
Definition Polytope2D.hpp:303
std::size_t n_cells() const
number of cells in the mesh.
Definition Mesh2D.hpp:63
size_t global_index() const
Return the global index of the Polytope.
Definition Polytope2D.hpp:73
std::size_t n_edges() const
number of edges in the mesh.
Definition Mesh2D.hpp:61
std::vector< QuadratureNode > QuadratureRule
Definition quadraturerule.hpp:55
QuadratureRule generate_quadrature_rule(const Cell &T, const int doe, const bool force_split)
Generate quadrature rule on mesh element.
Definition quadraturerule.cpp:10
if(strcmp(field, 'real')) % real valued entries T
Definition mmread.m:93
for j
Definition mmread.m:174
Definition PastixInterface.hpp:16
Definition ddr-klplate.hpp:27
Description of one node and one weight from a quadrature rule.
Definition quadraturerule.hpp:41