22 #ifndef HHOBRINKMAN_HPP
23 #define HHOBRINKMAN_HPP
27 #include <boost/math/constants/constants.hpp>
29 #include <Eigen/Sparse>
30 #include <unsupported/Eigen/SparseExtra>
77 double muT =
mu.
value(T, T.center_mass());
83 return (muT > thres ? std::max(thres, kappainvT * std::pow(T.diam(), 2) / muT) : 1./thres);
97 energy( std::sqrt( std::pow(norm_u, 2) + std::pow(norm_p, 2) ) )
134 std::ostream & output = std::cout
143 Eigen::VectorXd & UDir
221 return m_vhhospace.
mesh();
268 const int doe_cell = -1,
269 const int doe_face = -1
274 const std::vector<Eigen::VectorXd> & list_dofs,
280 const Eigen::VectorXd & p
288 const Eigen::VectorXd & u,
289 const std::pair<std::vector<VectorRd>,
VectorRd> & surf
294 std::pair<Eigen::MatrixXd, Eigen::VectorXd>
295 _compute_local_contribution(
306 void _assemble_local_contribution(
308 const std::pair<Eigen::MatrixXd, Eigen::VectorXd> & lsT,
309 std::list<Eigen::Triplet<double> > & A1,
310 Eigen::VectorXd & b1,
311 std::list<Eigen::Triplet<double> > & A2,
321 std::ostream & m_output;
322 const size_t m_nloc_sc_u;
323 const size_t m_nloc_sc_p;
327 Eigen::VectorXd m_sc_b;
328 std::pair<double,double> m_stab_par;
336 static const double PI = boost::math::constants::pi<double>();
344 double constexpr
eps=1e-14;
406 M.row(0) << 2. * x(0),
421 return std::pow(x(0), 2) + x(1)*x(2) - 1./3. - 1./4.;
451 std::pow(x(0), 3) + std::pow(x(0), 2)*x(1),
452 std::pow(x(1), 2) + x(0) * x(1) * x(2),
460 M.row(0) << 3. * std::pow(x(0), 2) + 2 * x(0) * x(1),
463 M.row(1) << x(1) * x(2),
464 2 * x(1) + x(0) * x(2),
468 3 * std::pow(x(2), 2);
475 return std::pow(x(0), 3) + x(1)*x(2) - 1./4. - 1./4.;
480 return VectorRd( 3 * std::pow(x(0), 2),
507 0.5 * sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
508 0.5 * cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
509 -cos(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2))
516 M.row(0) <<
PI * cos(2 *
PI * x(0)) * cos(2 *
PI * x(1)) * cos(2 *
PI * x(2)),
517 -
PI * sin(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
518 -
PI * sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2));
519 M.row(1) << -
PI * sin(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
520 PI * cos(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
521 -
PI * cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * sin(2. *
PI * x(2));
522 M.row(2) << 2 *
PI * sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
523 2 *
PI * cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
524 - 2 *
PI * cos(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * cos(2. *
PI * x(2));
537 cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
538 sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
539 sin(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2))
546 sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
547 cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
548 -2. * cos(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2))
575 VectorRd uDvalue = VectorRd::Zero();
587 M.row(0) << sin(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
588 -cos(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
589 -cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2));
590 M.row(1) << -cos(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
591 sin(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * sin(2. *
PI * x(2)),
592 -sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * cos(2. *
PI * x(2));
593 M.row(2) << -cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
594 -sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
595 sin(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * sin(2. *
PI * x(2));
606 static std::function<double(
const VectorRd &)>
620 sin(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
621 cos(2. *
PI * x(0)) * sin(2. *
PI * x(1)) * cos(2. *
PI * x(2)),
622 -2. * cos(2. *
PI * x(0)) * cos(2. *
PI * x(1)) * sin(2. *
PI * x(2))
641 static std::function<double(
const VectorRd &)>
642 XiV = [](
const VectorRd & x)->
double {
return (x.z() > 2*std::abs(x.x()-1.)+0.5-
eps); };
660 }
else if ( (x.z() > 0.8) && (x.x() > 0.75-
eps) && (x.x() < 1.25+
eps) ){
671 return MatrixRd::Zero();
681 return VectorRd::Zero();
698 static std::function<double(
const VectorRd &)>
700 return ( (x.x()>0) && (x.x()<1) && (x.y()>0) && (x.y()<1) && (x.z()>-1) ? 1. : 0.);
702 static std::function<double(
const VectorRd &)>
704 return ( (x.x()>=1) && (x.x()<2) && (x.y()>0) && (x.y()<1) && (x.z()> -.75 + .25*(x.x()-1)) ? 1. : 0.);
706 static std::function<double(
const VectorRd &)>
728 return MatrixRd::Zero();
738 return VectorRd::Zero();
759 static std::function<double(
const VectorRd &)>
761 return ( (x.x()<0.5) ? 1. : 0.);
763 static std::function<double(
const VectorRd &)>
786 sin(
PI * x.y()) * sin(
PI * x.z()),
793 return cos(
PI * x.x()) * (x.x() - 0.5)
796 x.y() + cos(
PI * x.z()),
803 return cos(
PI * x.x()) * (x.x() - 0.5)
805 sin(
PI * x.y()) + cos(
PI * x.z()),
807 std::pow(x.y(), 2) * std::pow(x.z(), 2)
819 M.row(0) << 0, -exp(-x.y() - x.z()), -exp(-x.y() - x.z());
820 M.row(1) << 0,
PI*sin(
PI*x.z())*cos(
PI*x.y()),
PI*sin(
PI*x.y())*cos(
PI*x.z());
821 M.row(2) << 0, x.z(), x.y();
829 M.row(0) << -
PI*(x.x() - 0.5)*(x.y() + x.z())*sin(
PI*x.x()) + (x.y() + x.z())*cos(
PI*x.x()),
830 (x.x() - 0.5)*cos(
PI*x.x()),
831 (x.x() - 0.5)*cos(
PI*x.x());
832 M.row(1) << -
PI*(x.x() - 0.5)*(x.y() + cos(
PI*x.z()))*sin(
PI*x.x()) + (x.y() + cos(
PI*x.z()))*cos(
PI*x.x()),
833 (x.x() - 0.5)*cos(
PI*x.x()),
834 -
PI*(x.x() - 0.5)*sin(
PI*x.z())*cos(
PI*x.x());
835 M.row(2) << -
PI*(x.x() - 0.5)*sin(
PI*x.x())*sin(
PI*x.y()) + sin(
PI*x.y())*cos(
PI*x.x()),
836 PI*(x.x() - 0.5)*cos(
PI*x.x())*cos(
PI*x.y()),
845 M.row(0) << -
PI*(x.x() - 0.5)*(sin(
PI*x.y()) + cos(
PI*x.z()))*sin(
PI*x.x()) + (sin(
PI*x.y()) + cos(
PI*x.z()))*cos(
PI*x.x()),
846 PI*(x.x() - 0.5)*cos(
PI*x.x())*cos(
PI*x.y()),
847 -
PI*(x.x() - 0.5)*sin(
PI*x.z())*cos(
PI*x.x());
848 M.row(1) << -
PI*pow(x.z(), 3)*(x.x() - 0.5)*sin(
PI*x.x()) + pow(x.z(), 3)*cos(
PI*x.x()),
850 3*pow(x.z(), 2)*(x.x() - 0.5)*cos(
PI*x.x());
851 M.row(2) << -
PI*pow(x.y(), 2)*pow(x.z(), 2)*(x.x() - 0.5)*sin(
PI*x.x()) + pow(x.y(), 2)*pow(x.z(), 2)*cos(
PI*x.x()),
852 2*x.y()*pow(x.z(), 2)*(x.x() - 0.5)*cos(
PI*x.x()),
853 2*pow(x.y(), 2)*x.z()*(x.x() - 0.5)*cos(
PI*x.x());
882 2*exp(-x.y() - x.z()),
883 -2*pow(
PI, 2)*sin(
PI*x.y())*sin(
PI*x.z()),
891 -pow(
PI, 2)*(x.x() - 0.5)*(x.y() + x.z())*cos(
PI*x.x()) - 2*
PI*(x.y() + x.z())*sin(
PI*x.x()),
892 -pow(
PI, 2)*(x.x() - 0.5)*(x.y() + cos(
PI*x.z()))*cos(
PI*x.x()) - pow(
PI, 2)*(x.x() - 0.5)*cos(
PI*x.x())*cos(
PI*x.z()) - 2*
PI*(x.y() + cos(
PI*x.z()))*sin(
PI*x.x()),
893 -2*pow(
PI, 2)*(x.x() - 0.5)*sin(
PI*x.y())*cos(
PI*x.x()) - 2*
PI*sin(
PI*x.x())*sin(
PI*x.y())
901 -pow(
PI, 2)*(x.x() - 0.5)*(sin(
PI*x.y()) + cos(
PI*x.z()))*cos(
PI*x.x()) - pow(
PI, 2)*(x.x() - 0.5)*sin(
PI*x.y())*cos(
PI*x.x()) - pow(
PI, 2)*(x.x() - 0.5)*cos(
PI*x.x())*cos(
PI*x.z()) - 2*
PI*(sin(
PI*x.y()) + cos(
PI*x.z()))*sin(
PI*x.x()),
902 -pow(
PI, 2)*pow(x.z(), 3)*(x.x() - 0.5)*cos(
PI*x.x()) - 2*
PI*pow(x.z(), 3)*sin(
PI*x.x()) + 6*x.z()*(x.x() - 0.5)*cos(
PI*x.x()),
903 -pow(
PI, 2)*pow(x.y(), 2)*pow(x.z(), 2)*(x.x() - 0.5)*cos(
PI*x.x()) - 2*
PI*pow(x.y(), 2)*pow(x.z(), 2)*sin(
PI*x.x()) + 2*pow(x.y(), 2)*(x.x() - 0.5)*cos(
PI*x.x()) + 2*pow(x.z(), 2)*(x.x() - 0.5)*cos(
PI*x.x())
The BoundaryConditions class provides definition of boundary conditions.
Definition: BoundaryConditions.hpp:43
const size_t n_dir_faces() const
Returns the number of Dirichlet faces.
Definition: BoundaryConditions.hpp:63
Base class for global DOF spaces. Provides functions to manipulate global DOFs (the local version bei...
Definition: globaldofspace.hpp:16
Class definition: polynomial bases and operators.
Definition: vhhospace.hpp:52
Class to describe a mesh.
Definition: MeshND.hpp:17
Eigen::Matrix3d MatrixRd
Definition: basis.hpp:51
size_t dimension() const
Returns the dimension of the global space (all DOFs for all geometric entities)
Definition: localdofspace.hpp:80
size_t numLocalDofsFace() const
Returns the number of local face DOFs.
Definition: localdofspace.hpp:49
IntegralWeightValueType value
Definition: integralweight.hpp:69
static const double PI
Definition: ddr-magnetostatics.hpp:186
static Magnetostatics::SolutionPotentialType linear_u
Definition: ddr-magnetostatics.hpp:213
static Magnetostatics::PermeabilityType linear_mu
Definition: ddr-magnetostatics.hpp:232
static Magnetostatics::SolutionPotentialType trigonometric_u
Definition: ddr-magnetostatics.hpp:237
static Magnetostatics::PermeabilityType trigonometric_mu
Definition: ddr-magnetostatics.hpp:264
static Magnetostatics::ForcingTermType trigonometric_f
Definition: ddr-magnetostatics.hpp:255
static Magnetostatics::ForcingTermType linear_f
Definition: ddr-magnetostatics.hpp:227
static Stokes::PressureType trigonometric_p
Definition: ddr-stokes.hpp:269
static Stokes::PressureType linear_p
Definition: ddr-stokes.hpp:310
double pressure_scaling
Definition: ddr-stokes.hpp:249
const Mesh & mesh() const
Return a const reference to the mesh.
Definition: vhhospace.hpp:128
bool use_threads
Definition: HHO_DiffAdvecReac.hpp:47
static std::function< VectorRd(const VectorRd &)> beta_D
Definition: hho-brinkman.hpp:802
static Brinkman::CompressibilityForcingTermType vcracked_g
Definition: hho-brinkman.hpp:691
static Brinkman::PressureType regimes_p
Definition: hho-brinkman.hpp:567
constexpr double viscosity_in_cavity
Definition: hho-brinkman.hpp:710
size_t numDirDOFs() const
Returns the number of Dirichlet DOFs.
Definition: hho-brinkman.hpp:177
static BrinkmanParameters::PermeabilityInvType regimes_kappainv
Definition: hho-brinkman.hpp:564
static Brinkman::CompressibilityForcingTermType regimes_g
Definition: hho-brinkman.hpp:635
static Brinkman::MomentumForcingTermType cavity_f
Definition: hho-brinkman.hpp:743
Eigen::VectorXd interpolate(const VelocityType &u, const PressureType &p, const int doe_cell=-1, const int doe_face=-1) const
Interpolates velocity and pressure.
Definition: hho-brinkman.cpp:717
double energy
Global energy norm
Definition: hho-brinkman.hpp:104
static Brinkman::PressureType vcracked_p
Definition: hho-brinkman.hpp:675
static Brinkman::VelocityType cavity_u
Definition: hho-brinkman.hpp:722
std::function< VectorRd(const VectorRd &)> MomentumForcingTermType
Definition: hho-brinkman.hpp:121
static Brinkman::MomentumForcingTermType DivGrad_alpha
Definition: hho-brinkman.hpp:880
static Brinkman::PressureGradientType cubic_gradp
Definition: hho-brinkman.hpp:479
static Brinkman::VelocityType cubic_u
Definition: hho-brinkman.hpp:449
static Brinkman::VelocityGradientType regimes_graduD
Definition: hho-brinkman.hpp:583
double u
Energy norm of velocity.
Definition: hho-brinkman.hpp:100
const std::pair< double, double > & stabilizationParameter() const
Returns the stabilization parameter (scaling)
Definition: hho-brinkman.hpp:245
static const VectorRd vec_p
Definition: hho-brinkman.hpp:356
static Brinkman::VelocityType regimes_u
Definition: hho-brinkman.hpp:612
static std::function< double(const VectorRd &)> XiBox
Definition: hho-brinkman.hpp:707
IntegralWeight PermeabilityInvType
Definition: hho-brinkman.hpp:61
std::function< double(const VectorRd &)> PressureType
Definition: hho-brinkman.hpp:125
static Brinkman::PressureGradientType cavity_gradp
Definition: hho-brinkman.hpp:737
size_t sizeSystem() const
Returns the size of the final system with Lagrange multiplier, after application of SC and removal of...
Definition: hho-brinkman.hpp:201
static Brinkman::VelocityGradientType tworegions_gradu
Definition: hho-brinkman.hpp:859
constexpr double permeability_in_wedge
Definition: hho-brinkman.hpp:714
static std::function< double(const VectorRd &)> XiWedge
Definition: hho-brinkman.hpp:703
static BrinkmanParameters::ViscosityType cavity_mu
Definition: hho-brinkman.hpp:711
static Brinkman::VelocityGradientType quadratic_gradu
Definition: hho-brinkman.hpp:404
IntegralWeight ViscosityType
Definition: hho-brinkman.hpp:60
static BrinkmanParameters::PermeabilityInvType tworegions_kappainv
Definition: hho-brinkman.hpp:778
static BrinkmanParameters::PermeabilityInvType vcracked_kappainv
Definition: hho-brinkman.hpp:650
static Brinkman::VelocityGradientType regimes_gradu
Definition: hho-brinkman.hpp:615
static std::function< double(const VectorRd &)> Xi_D
Definition: hho-brinkman.hpp:764
constexpr double permeability_D
Definition: hho-brinkman.hpp:776
std::pair< double, double > computeFlux(const Eigen::VectorXd &u, const std::pair< std::vector< VectorRd >, VectorRd > &surf) const
Compute the velocity flux (integral of u.n) across a given surface, and the area of the surface,...
Definition: hho-brinkman.cpp:868
static BrinkmanParameters::PermeabilityInvType trigonometric_kappainv
Definition: hho-brinkman.hpp:500
size_t numNonSCDOFs() const
Returns the number of DOFs after SC and with Lagrange multiplier, but before eliminating Dirichlet DO...
Definition: hho-brinkman.hpp:195
static std::function< VectorRd(const VectorRd &)> beta_S
Definition: hho-brinkman.hpp:792
static std::function< double(const VectorRd &)> XiV
Definition: hho-brinkman.hpp:642
static std::function< VectorRd(const VectorRd &)> alpha
Definition: hho-brinkman.hpp:783
size_t dimVelocity() const
Returns the dimension of velocity space.
Definition: hho-brinkman.hpp:183
ViscosityType mu
Definition: hho-brinkman.hpp:86
PermeabilityInvType kappainv
Definition: hho-brinkman.hpp:87
static Brinkman::VelocityGradientType cavity_gradu
Definition: hho-brinkman.hpp:727
static Brinkman::VelocityType regimes_uD
Definition: hho-brinkman.hpp:574
constexpr double eps
Definition: hho-brinkman.hpp:344
static Brinkman::MomentumForcingTermType regimes_f
Definition: hho-brinkman.hpp:618
static BrinkmanParameters::ViscosityType tworegions_mu
Definition: hho-brinkman.hpp:772
static Brinkman::CompressibilityForcingTermType cubic_g
Definition: hho-brinkman.hpp:492
static BrinkmanParameters::PermeabilityInvType cubic_kappainv
Definition: hho-brinkman.hpp:446
constexpr double permeability_in_box
Definition: hho-brinkman.hpp:715
static Brinkman::PressureType cavity_p
Definition: hho-brinkman.hpp:732
static Brinkman::MomentumForcingTermType DivGrad_beta_S
Definition: hho-brinkman.hpp:889
static Brinkman::PressureType cubic_p
Definition: hho-brinkman.hpp:474
double Cf(const Cell &T) const
Returns the friction coefficient.
Definition: hho-brinkman.hpp:75
static BrinkmanParameters::ViscosityType cubic_mu
Definition: hho-brinkman.hpp:443
const Eigen::VectorXd & systemVector() const
Returns the linear system right-hand side vector.
Definition: hho-brinkman.hpp:235
std::function< MatrixRd(const VectorRd &)> VelocityGradientType
Definition: hho-brinkman.hpp:124
static Brinkman::PressureGradientType linear_gradp
Definition: hho-brinkman.hpp:374
static Brinkman::PressureType quadratic_p
Definition: hho-brinkman.hpp:420
std::pair< double, double > & stabilizationParameter()
Returns the stabilization parameter.
Definition: hho-brinkman.hpp:250
double scaling_mu
Definition: hho-brinkman.hpp:340
Brinkman(const VHHOSpace &vhho_space, const GlobalDOFSpace &p_space, const BoundaryConditions &BC, bool use_threads, std::ostream &output=std::cout)
Constructor.
Definition: hho-brinkman.cpp:459
std::vector< double > pressureVertexValues(const Eigen::VectorXd &p) const
Create vertex values for the pressure (from the element values), for plotting.
Definition: hho-brinkman.cpp:817
std::function< VectorRd(const VectorRd &)> VelocityType
Definition: hho-brinkman.hpp:123
static Brinkman::PressureType tworegions_p
Definition: hho-brinkman.hpp:864
static Brinkman::PressureGradientType vcracked_gradp
Definition: hho-brinkman.hpp:680
static std::function< double(const VectorRd &)> XiCavity
Definition: hho-brinkman.hpp:699
std::pair< double, double > computeConditionNum() const
Computes the condition number of the matrix.
Definition: hho-brinkman.cpp:843
static Brinkman::PressureGradientType regimes_gradp
Definition: hho-brinkman.hpp:570
static Brinkman::MomentumForcingTermType vcracked_f
Definition: hho-brinkman.hpp:686
static Brinkman::MomentumForcingTermType tworegions_f
Definition: hho-brinkman.hpp:908
static Brinkman::VelocityType regimes_uS
Definition: hho-brinkman.hpp:603
static Brinkman::CompressibilityForcingTermType cavity_g
Definition: hho-brinkman.hpp:748
size_t nloc_sc_u() const
Returns the local number of velocity statically condensed DOFs.
Definition: hho-brinkman.hpp:147
static Brinkman::MomentumForcingTermType quadratic_f
Definition: hho-brinkman.hpp:433
static BrinkmanParameters::ViscosityType quadratic_mu
Definition: hho-brinkman.hpp:389
static Brinkman::VelocityGradientType grad_beta_S
Definition: hho-brinkman.hpp:827
static BrinkmanParameters::PermeabilityInvType quadratic_kappainv
Definition: hho-brinkman.hpp:392
static Brinkman::PressureGradientType trigonometric_gradp
Definition: hho-brinkman.hpp:535
static BrinkmanParameters::ViscosityType regimes_mu
Definition: hho-brinkman.hpp:561
double scaling_kappainv
Definition: hho-brinkman.hpp:341
size_t numSCDOFs_p() const
Returns the number of pressure statically condensed DOFs.
Definition: hho-brinkman.hpp:165
std::vector< BrinkmanNorms > computeBrinkmanNorms(const std::vector< Eigen::VectorXd > &list_dofs, const BrinkmanParameters ¶) const
Compute the discrete energy norm of a family of vectors representing the dofs.
Definition: hho-brinkman.cpp:753
static Brinkman::VelocityGradientType linear_gradu
Definition: hho-brinkman.hpp:369
SystemMatrixType & systemMatrix()
Returns the linear system matrix.
Definition: hho-brinkman.hpp:230
static BrinkmanParameters::ViscosityType vcracked_mu
Definition: hho-brinkman.hpp:646
void assembleLinearSystem(const MomentumForcingTermType &f, const CompressibilityForcingTermType &g, const BrinkmanParameters ¶, const VelocityType &u, Eigen::VectorXd &UDir)
Assemble the global system
Definition: hho-brinkman.cpp:490
static Brinkman::MomentumForcingTermType cubic_f
Definition: hho-brinkman.hpp:487
static BrinkmanParameters::PermeabilityInvType cavity_kappainv
Definition: hho-brinkman.hpp:717
static Brinkman::VelocityType vcracked_u
Definition: hho-brinkman.hpp:655
static Brinkman::VelocityGradientType trigonometric_gradu
Definition: hho-brinkman.hpp:514
BrinkmanParameters(const ViscosityType &_mu, const PermeabilityInvType &_kappainv)
Constructor.
Definition: hho-brinkman.hpp:64
static BrinkmanParameters::PermeabilityInvType linear_kappainv
Definition: hho-brinkman.hpp:353
size_t dimPressure() const
Returns the dimension of pressure space.
Definition: hho-brinkman.hpp:189
static std::function< double(const VectorRd &)> Xi_S
Definition: hho-brinkman.hpp:760
static Brinkman::CompressibilityForcingTermType trigonometric_g
Definition: hho-brinkman.hpp:554
static const MatrixRd mat_u
Definition: hho-brinkman.hpp:355
std::function< VectorRd(const VectorRd &)> PressureGradientType
Definition: hho-brinkman.hpp:126
Eigen::VectorXd & scVector()
Returns the static condensation rhs.
Definition: hho-brinkman.hpp:260
static Brinkman::PressureGradientType tworegions_gradp
Definition: hho-brinkman.hpp:867
Eigen::SparseMatrix< double > SystemMatrixType
Definition: hho-brinkman.hpp:119
static Brinkman::MomentumForcingTermType DivGrad_beta_D
Definition: hho-brinkman.hpp:899
static Brinkman::PressureGradientType quadratic_gradp
Definition: hho-brinkman.hpp:425
static Brinkman::VelocityGradientType regimes_graduS
Definition: hho-brinkman.hpp:604
static Brinkman::VelocityType quadratic_u
Definition: hho-brinkman.hpp:395
std::function< double(const VectorRd &)> CompressibilityForcingTermType
Definition: hho-brinkman.hpp:122
constexpr double viscosity_S
Definition: hho-brinkman.hpp:769
static Brinkman::CompressibilityForcingTermType tworegions_g
Definition: hho-brinkman.hpp:917
static Brinkman::VelocityType tworegions_u
Definition: hho-brinkman.hpp:812
static Brinkman::VelocityGradientType grad_alpha
Definition: hho-brinkman.hpp:817
double p
L2 norm of p.
Definition: hho-brinkman.hpp:103
static Brinkman::VelocityGradientType cubic_gradu
Definition: hho-brinkman.hpp:458
Eigen::VectorXd & systemVector()
Returns the linear system right-hand side vector.
Definition: hho-brinkman.hpp:240
constexpr double viscosity_D
Definition: hho-brinkman.hpp:770
static Brinkman::VelocityGradientType grad_beta_D
Definition: hho-brinkman.hpp:843
const SystemMatrixType & scMatrix() const
Returns the static condensation recovery operator.
Definition: hho-brinkman.hpp:255
size_t numSCDOFs() const
Returns the number of statically condensed DOFs.
Definition: hho-brinkman.hpp:171
size_t numSCDOFs_u() const
Returns the number of velocity statically condensed DOFs.
Definition: hho-brinkman.hpp:159
constexpr double permeability_S
Definition: hho-brinkman.hpp:775
const VHHOSpace & vhhospace() const
Returns the velocity space.
Definition: hho-brinkman.hpp:207
BrinkmanNorms(double norm_u, double norm_p)
Constructor.
Definition: hho-brinkman.hpp:94
static Brinkman::CompressibilityForcingTermType linear_g
Definition: hho-brinkman.hpp:384
const GlobalDOFSpace & pspace() const
Returns the pressure space.
Definition: hho-brinkman.hpp:213
const SystemMatrixType & systemMatrix() const
Returns the linear system matrix.
Definition: hho-brinkman.hpp:225
size_t nloc_sc_p() const
Returns the local number of pressure statically condensed DOFs.
Definition: hho-brinkman.hpp:153
static Brinkman::VelocityGradientType vcracked_gradu
Definition: hho-brinkman.hpp:670
static Brinkman::CompressibilityForcingTermType quadratic_g
Definition: hho-brinkman.hpp:438
const Mesh & mesh() const
Returns the mesh.
Definition: hho-brinkman.hpp:219
static std::function< double(const VectorRd &)> XiS
Definition: hho-brinkman.hpp:607
std::size_t n_cells() const
number of cells in the mesh.
Definition: MeshND.hpp:60
Definition: ddr-magnetostatics.hpp:40
MeshND::VectorRd< 2 > VectorRd
Definition: Mesh2D.hpp:14
MeshND::Cell< 2 > Cell
Definition: Mesh2D.hpp:13
Structure to store norm components (energy for velocity, L2 pressure, and global)
Definition: hho-brinkman.hpp:92
Structure to store physical parameter (and compute derived parameters)
Definition: hho-brinkman.hpp:59
Structure for Brinkman model.
Definition: hho-brinkman.hpp:118
Structure for weights (scalar, at the moment) in integral.
Definition: integralweight.hpp:36
Structure to store information for, and perform, local static condensation.
Definition: local_static_condensation.hpp:25