24#include <opm/common/TimingMacros.hpp>
26#include <opm/common/utility/numeric/RootFinders.hpp>
28#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
29#include <opm/material/fluidstates/SimpleModularFluidState.hpp>
30#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
40#include <fmt/format.h>
61namespace Miscibility {
63template<
class Flu
idSystem>
65 const std::vector<Scalar>& depth,
66 const std::vector<Scalar>& rs)
67 : pvtRegionIdx_(pvtRegionIdx)
68 , rsVsDepth_(depth, rs)
72template<
class Flu
idSystem>
80 const auto sat_rs = satRs(press, temp);
81 if (satGas > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
85 if (rsVsDepth_.xMin() > depth)
86 return std::min(sat_rs, rsVsDepth_.valueAt(0));
87 else if (rsVsDepth_.xMax() < depth)
88 return std::min(sat_rs, rsVsDepth_.valueAt(rsVsDepth_.numSamples() - 1));
89 return std::min(sat_rs, rsVsDepth_.eval(depth,
false));
93template<
class Flu
idSystem>
96satRs(
const Scalar press,
const Scalar temp)
const
98 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
101template<
class Flu
idSystem>
103 const std::vector<Scalar>& depth,
104 const std::vector<Scalar>& pbub)
105 : pvtRegionIdx_(pvtRegionIdx)
106 , pbubVsDepth_(depth, pbub)
110template<
class Flu
idSystem>
116 const Scalar satGas)
const
120 if (pbubVsDepth_.xMin() > depth)
121 press = pbubVsDepth_.valueAt(0);
122 else if (pbubVsDepth_.xMax() < depth)
123 press = pbubVsDepth_.valueAt(pbubVsDepth_.numSamples() - 1);
125 press = pbubVsDepth_.eval(depth,
false);
127 return satRs(std::min(press, cellPress), temp);
130template<
class Flu
idSystem>
133satRs(
const Scalar press,
const Scalar temp)
const
135 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
138template<
class Flu
idSystem>
140 const std::vector<Scalar>& depth,
141 const std::vector<Scalar>& pdew)
142 : pvtRegionIdx_(pvtRegionIdx)
143 , pdewVsDepth_(depth, pdew)
147template<
class Flu
idSystem>
148typename PDVD<FluidSystem>::Scalar
151 const Scalar cellPress,
153 const Scalar satOil)
const
155 Scalar press = cellPress;
157 if (pdewVsDepth_.xMin() > depth)
158 press = pdewVsDepth_.valueAt(0);
159 else if (pdewVsDepth_.xMax() < depth)
160 press = pdewVsDepth_.valueAt(pdewVsDepth_.numSamples() - 1);
162 press = pdewVsDepth_.eval(depth,
false);
164 return satRv(std::min(press, cellPress), temp);
167template<
class Flu
idSystem>
168typename PDVD<FluidSystem>::Scalar
170satRv(
const Scalar press,
const Scalar temp)
const
172 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);
175template<
class Flu
idSystem>
177 const std::vector<Scalar>& depth,
178 const std::vector<Scalar>& rv)
179 : pvtRegionIdx_(pvtRegionIdx)
180 , rvVsDepth_(depth, rv)
184template<
class Flu
idSystem>
185typename RvVD<FluidSystem>::Scalar
190 const Scalar satOil)
const
192 if (satOil < - std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
193 throw std::logic_error {
194 "Must not pass negative oil saturation"
197 const auto sat_rv = satRv(press, temp);
198 if (satOil > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
202 if (rvVsDepth_.xMin() > depth)
203 return std::min(sat_rv, rvVsDepth_.valueAt(0));
204 else if (rvVsDepth_.xMax() < depth)
205 return std::min(sat_rv, rvVsDepth_.valueAt(rvVsDepth_.numSamples() - 1));
206 return std::min(sat_rv, rvVsDepth_.eval(depth,
false));
210template<
class Flu
idSystem>
211typename RvVD<FluidSystem>::Scalar
213satRv(
const Scalar press,
const Scalar temp)
const
215 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);
218template<
class Flu
idSystem>
220 const std::vector<Scalar>& depth,
221 const std::vector<Scalar>& rvw)
222 : pvtRegionIdx_(pvtRegionIdx)
223 , rvwVsDepth_(depth, rvw)
227template<
class Flu
idSystem>
228typename RvwVD<FluidSystem>::Scalar
233 const Scalar satWat)
const
235 if (satWat < - std::sqrt(std::numeric_limits<double>::epsilon())) {
236 throw std::logic_error {
237 "Must not pass negative water saturation"
241 const auto sat_rvw = satRvw(press, temp);
242 if (satWat > std::sqrt(std::numeric_limits<double>::epsilon())) {
246 if (rvwVsDepth_.xMin() > depth)
247 return std::min(sat_rvw,rvwVsDepth_.valueAt(0));
248 else if (rvwVsDepth_.xMax() < depth)
249 return std::min(sat_rvw, rvwVsDepth_.valueAt(rvwVsDepth_.numSamples() - 1));
250 return std::min(sat_rvw, rvwVsDepth_.eval(depth,
false));
254template<
class Flu
idSystem>
255typename RvwVD<FluidSystem>::Scalar
257satRvw(
const Scalar press,
const Scalar temp)
const
259 return FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press);
262template<
class Flu
idSystem>
264RsSatAtContact(
const int pvtRegionIdx,
const Scalar pContact,
const Scalar T_contact)
265 : pvtRegionIdx_(pvtRegionIdx)
267 rsSatContact_ = satRs(pContact, T_contact);
270template<
class Flu
idSystem>
271typename RsSatAtContact<FluidSystem>::Scalar
276 const Scalar satGas)
const
278 if (satGas > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
279 return satRs(press, temp);
282 return std::min(satRs(press, temp), rsSatContact_);
286template<
class Flu
idSystem>
287typename RsSatAtContact<FluidSystem>::Scalar
289satRs(
const Scalar press,
const Scalar temp)
const
291 return FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press);
294template<
class Flu
idSystem>
296RvSatAtContact(
const int pvtRegionIdx,
const Scalar pContact,
const Scalar T_contact)
297 : pvtRegionIdx_(pvtRegionIdx)
299 rvSatContact_ = satRv(pContact, T_contact);
302template<
class Flu
idSystem>
303typename RvSatAtContact<FluidSystem>::Scalar
308 const Scalar satOil)
const
310 if (satOil > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
311 return satRv(press, temp);
314 return std::min(satRv(press, temp), rvSatContact_);
318template<
class Flu
idSystem>
319typename RvSatAtContact<FluidSystem>::Scalar
321satRv(
const Scalar press,
const Scalar temp)
const
323 return FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press);;
326template<
class Flu
idSystem>
328RvwSatAtContact(
const int pvtRegionIdx,
const Scalar pContact,
const Scalar T_contact)
329 : pvtRegionIdx_(pvtRegionIdx)
331 rvwSatContact_ = satRvw(pContact, T_contact);
334template<
class Flu
idSystem>
335typename RvwSatAtContact<FluidSystem>::Scalar
340 const Scalar satWat)
const
342 if (satWat > std::sqrt(std::numeric_limits<Scalar>::epsilon())) {
343 return satRvw(press, temp);
346 return std::min(satRvw(press, temp), rvwSatContact_);
350template<
class Flu
idSystem>
351typename RvwSatAtContact<FluidSystem>::Scalar
353satRvw(
const Scalar press,
const Scalar temp)
const
355 return FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press);;
358template <
class Flu
idSystem>
361 : rs_constant_(rs_constant)
362 , pb_constant_(pb_constant)
365 if (this->rs_constant_ < 0.0) {
366 throw std::invalid_argument(
"RSCONST RS value cannot be negative");
369 if (this->pb_constant_ <= 0.0) {
370 throw std::invalid_argument(
"RSCONST bubble point pressure must be positive");
374template <
class Flu
idSystem>
375typename FluidSystem::Scalar
377 [[maybe_unused]]
const Scalar temp)
const
380 return this->rs_constant_;
383template <
class Flu
idSystem>
384typename FluidSystem::Scalar
386 [[maybe_unused]]
const Scalar press,
387 [[maybe_unused]]
const Scalar temp,
388 [[maybe_unused]]
const Scalar satGas)
const
395 if (press < this->pb_constant_) {
396 throw std::invalid_argument {
397 fmt::format(
"Pressure {} is below bubble point pressure {}",
398 press, this->pb_constant_)
402 return this->rs_constant_;
407template<
class Scalar>
412 const TabulatedFunction& tempVdTable,
413 const TabulatedFunction& saltVdTable,
419 , tempVdTable_ (tempVdTable)
420 , saltVdTable_ (saltVdTable)
425template<
class Scalar>
428 return this->rec_.datumDepth();
431template<
class Scalar>
434 return this->rec_.datumDepthPressure();
437template<
class Scalar>
440 return this->rec_.waterOilContactDepth();
443template<
class Scalar>
446 return this->rec_.waterOilContactCapillaryPressure();
449template<
class Scalar>
452 return this->rec_.gasOilContactDepth();
455template<
class Scalar>
458 return this->rec_.gasOilContactCapillaryPressure();
461template<
class Scalar>
464 return this->rec_.initializationTargetAccuracy();
467template<
class Scalar>
474template<
class Scalar>
481template<
class Scalar>
488template<
class Scalar>
489const typename EquilReg<Scalar>::TabulatedFunction&
495template<
class Scalar>
496const typename EquilReg<Scalar>::TabulatedFunction&
502template<
class Scalar>
505 return this->pvtIdx_;
508template<
class Flu
idSystem,
class MaterialLawManager>
510PcEq(
const MaterialLawManager& materialLawManager,
514 : materialLawManager_(materialLawManager)
517 , targetPc_(targetPc)
521template<
class Flu
idSystem,
class MaterialLawManager>
526 const auto& matParams = materialLawManager_.materialLawParams(cell_);
528 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 0.0);
529 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
530 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
531 fluidState.setSaturation(phase_, s);
534 assert(phase_ != FluidSystem::oilPhaseIdx);
536 if (phase_ == FluidSystem::gasPhaseIdx) {
537 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 1.0 - s);
538 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 1.0 - s);
540 std::array<Scalar, FluidSystem::numPhases> pc{0.0};
541 using MaterialLaw =
typename MaterialLawManager::MaterialLaw;
542 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
543 Scalar sign = (phase_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
544 Scalar pcPhase = pc[FluidSystem::oilPhaseIdx] + sign * pc[phase_];
545 return pcPhase - targetPc_;
548template<
class Flu
idSystem,
class MaterialLawManager>
550PcEqSum(
const MaterialLawManager& materialLawManager,
555 : materialLawManager_(materialLawManager)
559 , targetPc_(targetPc)
563template<
class Flu
idSystem,
class MaterialLawManager>
568 const auto& matParams = materialLawManager_.materialLawParams(cell_);
570 fluidState.setSaturation(FluidSystem::waterPhaseIdx, 0.0);
571 fluidState.setSaturation(FluidSystem::oilPhaseIdx, 0.0);
572 fluidState.setSaturation(FluidSystem::gasPhaseIdx, 0.0);
573 fluidState.setSaturation(phase1_, s);
574 fluidState.setSaturation(phase2_, 1.0 - s);
576 std::array<Scalar, FluidSystem::numPhases> pc {0.0};
578 using MaterialLaw =
typename MaterialLawManager::MaterialLaw;
579 MaterialLaw::capillaryPressures(pc, matParams, fluidState);
580 Scalar sign1 = (phase1_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
581 Scalar pc1 = pc[FluidSystem::oilPhaseIdx] + sign1 * pc[phase1_];
582 Scalar sign2 = (phase2_ == FluidSystem::waterPhaseIdx)? -1.0 : 1.0;
583 Scalar pc2 = pc[FluidSystem::oilPhaseIdx] + sign2 * pc[phase2_];
584 return pc1 + pc2 - targetPc_;
587template <
class Flu
idSystem,
class MaterialLawManager>
588typename FluidSystem::Scalar
590 const int phase,
const int cell)
592 const auto& scaledDrainageInfo =
593 materialLawManager.oilWaterScaledEpsInfoDrainage(cell);
597 case FluidSystem::waterPhaseIdx:
598 return scaledDrainageInfo.Swl;
600 case FluidSystem::gasPhaseIdx:
601 return scaledDrainageInfo.Sgl;
603 case FluidSystem::oilPhaseIdx:
604 throw std::runtime_error(
"Min saturation not implemented for oil phase.");
607 throw std::runtime_error(
"Unknown phaseIdx .");
612template <
class Flu
idSystem,
class MaterialLawManager>
613typename FluidSystem::Scalar
615 const int phase,
const int cell)
617 const auto& scaledDrainageInfo =
618 materialLawManager.oilWaterScaledEpsInfoDrainage(cell);
622 case FluidSystem::waterPhaseIdx:
623 return scaledDrainageInfo.Swu;
625 case FluidSystem::gasPhaseIdx:
626 return scaledDrainageInfo.Sgu;
628 case FluidSystem::oilPhaseIdx:
629 throw std::runtime_error(
"Max saturation not implemented for oil phase.");
632 throw std::runtime_error(
"Unknown phaseIdx .");
637template <
class Flu
idSystem,
class MaterialLawManager>
638typename FluidSystem::Scalar
642 const typename FluidSystem::Scalar targetPc,
643 const bool increasing)
645 using Scalar =
typename FluidSystem::Scalar;
647 Scalar s0 = increasing ? maxSaturations<FluidSystem>(materialLawManager, phase, cell) : minSaturations<FluidSystem>(materialLawManager, phase, cell);
648 Scalar s1 = increasing ? minSaturations<FluidSystem>(materialLawManager, phase, cell) : maxSaturations<FluidSystem>(materialLawManager, phase, cell);
654 if (!std::isfinite(f0 + f1))
655 throw std::logic_error(fmt::format(
"The capillary pressure values {} and {} are not finite", f0, f1));
662 const Scalar tol = 1e-10;
664 const int maxIter = -2*
static_cast<int>(std::log2(tol)) + 10;
665 int usedIterations = -1;
666 const Scalar root = RegulaFalsiBisection<ThrowOnError>::solve(f, s0, s1, maxIter, tol, usedIterations);
670template<
class Flu
idSystem,
class MaterialLawManager>
671typename FluidSystem::Scalar
676 const typename FluidSystem::Scalar targetPc)
678 using Scalar =
typename FluidSystem::Scalar;
680 Scalar s0 = minSaturations<FluidSystem>(materialLawManager, phase1, cell);
681 Scalar s1 = maxSaturations<FluidSystem>(materialLawManager, phase1, cell);
692 assert(f0 > 0.0 && f1 < 0.0);
693 const Scalar tol = 1e-10;
695 const int maxIter = -2*
static_cast<int>(std::log2(tol)) + 10;
696 int usedIterations = -1;
697 const Scalar root = RegulaFalsiBisection<ThrowOnError>::solve(f, s0, s1, maxIter, tol, usedIterations);
701template<
class Flu
idSystem,
class MaterialLawManager>
702typename FluidSystem::Scalar
704 const typename FluidSystem::Scalar cellDepth,
705 const typename FluidSystem::Scalar contactDepth,
708 const bool increasing)
710 using Scalar =
typename FluidSystem::Scalar;
711 const Scalar s0 = increasing ? maxSaturations<FluidSystem>(materialLawManager, phase, cell) : minSaturations<FluidSystem>(materialLawManager, phase, cell);
712 const Scalar s1 = increasing ? minSaturations<FluidSystem>(materialLawManager, phase, cell) : maxSaturations<FluidSystem>(materialLawManager, phase, cell);
714 if (cellDepth <= contactDepth) {
722template<
class Flu
idSystem,
class MaterialLawManager>
723bool isConstPc(
const MaterialLawManager& materialLawManager,
727 using Scalar =
typename FluidSystem::Scalar;
730 const Scalar f0 = f(minSaturations<FluidSystem>(materialLawManager, phase, cell));
731 const Scalar f1 = f(maxSaturations<FluidSystem>(materialLawManager, phase, cell));
732 return std::abs(f0 - f1) < std::numeric_limits<Scalar>::epsilon();
Auxiliary routines that to solve the ODEs that emerge from the hydrostatic equilibrium problem.
Scalar datum() const
Definition: EquilibrationHelpers_impl.hpp:426
int pvtIdx() const
Definition: EquilibrationHelpers_impl.hpp:503
EquilReg(const EquilRecord &rec, std::shared_ptr< Miscibility::RsFunction< Scalar > > rs, std::shared_ptr< Miscibility::RsFunction< Scalar > > rv, std::shared_ptr< Miscibility::RsFunction< Scalar > > rvw, const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtIdx)
Definition: EquilibrationHelpers_impl.hpp:408
Scalar pcgoGoc() const
Definition: EquilibrationHelpers_impl.hpp:456
int equilibrationAccuracy() const
Definition: EquilibrationHelpers_impl.hpp:462
const CalcEvaporation & evaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:476
Scalar zgoc() const
Definition: EquilibrationHelpers_impl.hpp:450
Scalar pressure() const
Definition: EquilibrationHelpers_impl.hpp:432
const TabulatedFunction & saltVdTable() const
Definition: EquilibrationHelpers_impl.hpp:490
Scalar pcowWoc() const
Definition: EquilibrationHelpers_impl.hpp:444
const TabulatedFunction & tempVdTable() const
Definition: EquilibrationHelpers_impl.hpp:497
Scalar zwoc() const
Definition: EquilibrationHelpers_impl.hpp:438
const CalcDissolution & dissolutionCalculator() const
Definition: EquilibrationHelpers_impl.hpp:469
const CalcWaterEvaporation & waterEvaporationCalculator() const
Definition: EquilibrationHelpers_impl.hpp:483
Definition: EquilibrationHelpers.hpp:216
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:113
PBVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pbub)
Definition: EquilibrationHelpers_impl.hpp:102
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:218
Definition: EquilibrationHelpers.hpp:269
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:150
PDVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pdew)
Definition: EquilibrationHelpers_impl.hpp:139
RsConst(const Scalar rs_constant, const Scalar pb_constant)
Definition: EquilibrationHelpers_impl.hpp:359
Scalar satRs(const Scalar press, const Scalar temp) const
Definition: EquilibrationHelpers_impl.hpp:376
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:385
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:614
Definition: EquilibrationHelpers.hpp:101
Definition: EquilibrationHelpers.hpp:162
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:75
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:164
RsVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rs)
Definition: EquilibrationHelpers_impl.hpp:64
Definition: EquilibrationHelpers.hpp:322
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satOil=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:187
RvVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rv)
Definition: EquilibrationHelpers_impl.hpp:176
Definition: EquilibrationHelpers.hpp:376
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satWat=0.0) const override
Definition: EquilibrationHelpers_impl.hpp:230
RvwVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rvw)
Definition: EquilibrationHelpers_impl.hpp:219
BlackOilFluidSystem< double > FluidSystemSimple
Definition: EquilibrationHelpers_impl.hpp:45
FluidSystem::Scalar satFromDepth(const MaterialLawManager &materialLawManager, const typename FluidSystem::Scalar cellDepth, const typename FluidSystem::Scalar contactDepth, const int phase, const int cell, const bool increasing=false)
Compute saturation from depth. Used for constant capillary pressure function.
Definition: EquilibrationHelpers_impl.hpp:703
FluidSystem::Scalar satFromPc(const MaterialLawManager &materialLawManager, const int phase, const int cell, const typename FluidSystem::Scalar targetPc, const bool increasing=false)
Definition: EquilibrationHelpers_impl.hpp:639
SimpleModularFluidState< double, 3, 3, FluidSystemSimple, false, false, false, false, true, false, false, false > SatOnlyFluidState
Definition: EquilibrationHelpers_impl.hpp:59
FluidSystem::Scalar minSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:589
FluidSystem::Scalar maxSaturations(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Definition: EquilibrationHelpers_impl.hpp:614
bool isConstPc(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Return true if capillary pressure function is constant.
Definition: EquilibrationHelpers_impl.hpp:723
FluidSystem::Scalar satFromSumOfPcs(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const typename FluidSystem::Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:672
Definition: blackoilbioeffectsmodules.hh:45
Definition: EquilibrationHelpers.hpp:801
Scalar operator()(Scalar s) const
Definition: EquilibrationHelpers_impl.hpp:524
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:802
PcEq(const MaterialLawManager &materialLawManager, const int phase, const int cell, const Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:510
Definition: EquilibrationHelpers.hpp:842
PcEqSum(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const Scalar targetPc)
Definition: EquilibrationHelpers_impl.hpp:550
Scalar operator()(Scalar s) const
Definition: EquilibrationHelpers_impl.hpp:566
typename FluidSystem::Scalar Scalar
Definition: EquilibrationHelpers.hpp:843