47 #ifndef SORPTION_IMPL_HH_ 48 #define SORPTION_IMPL_HH_ 52 #include <boost/math/tools/roots.hpp> 69 eps = T(ldexp(1.0F, 1-bits));
74 return fabs(a - b) <= (
eps * (std::max)(fabs(a), fabs(b)));
100 Linear(
double mult_coef) : mult_coef_(mult_coef) {}
103 return (mult_coef_*x);
118 Langmuir(
double mult_coef,
double alpha) : mult_coef_(mult_coef), alpha_(alpha) {}
121 return (mult_coef_*(alpha_ * x)/(alpha_ *x + 1));
138 Freundlich(
double mult_coef,
double exponent) : mult_coef_(mult_coef), exponent_(exponent){}
141 return (mult_coef_*pow(x, exponent_));
158 "The Boost solver of nonlinear equation did not converge in sorption model.\n" 159 <<
"Check input at the following address for possible user error: \t" << Input::EI_Address::val <<
"\n" 160 "Solver message: \n" << EI_BoostMessage::val);
164 "The total mass in sorption model became negative during the computation (value: " 165 << EI_TotalMass::val <<
").\n" 166 <<
"Check input at the following address for possible user error: \t" << Input::EI_Address::val <<
"\n");
179 ConcPair(
double x,
double y) : fluid(x), solid(y) {}
197 inline void reinit(
enum SorptionType sorption_type,
bool limited_solubility_on,
198 double aqua_density,
double scale_aqua,
double scale_sorbed,
199 double c_aqua_limit,
double mult_coef,
double second_coef);
206 void make_table(
int n_points);
212 inline void compute(
double &c_aqua,
double &c_sorbed);
219 inline void interpolate(
double &c_aqua,
double &c_sorbed);
225 return interpolation_table.size() != 0;
234 void make_table(
const Func &isotherm,
int n_points);
309 template <
class Func>
313 CrossFunction(
const Func &func_,
double total_mass,
double scale_aqua,
double scale_sorbed,
double rho_aqua)
314 : func(func_), total_mass_(total_mass),
315 scale_sorbed_(scale_sorbed), scale_aqua_(scale_aqua), rho_aqua_(rho_aqua)
321 return scale_sorbed_*func( conc_aqua/rho_aqua_) + (scale_aqua_) * conc_aqua - total_mass_;
335 double rho_aqua,
double scale_aqua,
double scale_sorbed,
336 double c_aqua_limit,
double mult_coef,
double second_coef)
338 adsorption_type_ = adsorption_type;
339 rho_aqua_ = rho_aqua;
340 scale_aqua_ = scale_aqua;
341 scale_sorbed_ = scale_sorbed;
342 limited_solubility_on_ = limited_solubility_on;
343 mult_coef_ = mult_coef*rho_aqua;
344 second_coef_ = second_coef;
346 table_limit_ = c_aqua_limit;
347 inv_scale_aqua_ = scale_aqua_/(scale_aqua_*scale_aqua_ + scale_sorbed_*scale_sorbed_);
348 inv_scale_sorbed_ = scale_sorbed_/(scale_aqua_*scale_aqua_ + scale_sorbed_*scale_sorbed_);
356 if (limited_solubility_on_ && (c_pair.
fluid > table_limit_)) {
357 result = precipitate( c_pair );
359 result = solve_conc( c_pair );
362 c_sorbed=result.
solid;
370 result = compute_projection( c_pair );
373 c_sorbed=result.
solid;
379 double total_mass = (scale_aqua_* c_pair.
fluid + scale_sorbed_ * c_pair.
solid);
381 THROW( Isotherm::ExcNegativeTotalMass()
382 << EI_TotalMass(total_mass)
385 double total_mass_steps = total_mass / total_mass_step_;
386 unsigned int total_mass_idx = static_cast <
unsigned int>(std::floor(total_mass_steps));
388 if (total_mass_idx < (interpolation_table.size() - 1) ) {
389 double rot_sorbed = interpolation_table[total_mass_idx] + (total_mass_steps - total_mass_idx)*(interpolation_table[total_mass_idx+1] - interpolation_table[total_mass_idx]);
390 return ConcPair( (total_mass * inv_scale_aqua_ - rot_sorbed * inv_scale_sorbed_),
391 (total_mass * inv_scale_sorbed_ + rot_sorbed * inv_scale_aqua_) );
393 if (limited_solubility_on_) {
394 return precipitate( c_pair );
396 return solve_conc( c_pair );
404 double total_mass = (scale_aqua_*c_pair.
fluid + scale_sorbed_ * c_pair.
solid);
406 (total_mass - scale_aqua_ * table_limit_)/scale_sorbed_ );
414 boost::uintmax_t max_iter = 20;
416 double total_mass = (scale_aqua_*c_pair.
fluid + scale_sorbed_ * c_pair.
solid);
417 double mass_limit = table_limit_*scale_aqua_ +
const_cast<Func &
>(isotherm)(table_limit_ / this->rho_aqua_)*scale_sorbed_;
418 double upper_solution_bound;
420 if(total_mass >= mass_limit)
422 mass_limit = total_mass;
424 upper_solution_bound = mass_limit / scale_aqua_;
425 CrossFunction<Func> eq_func(isotherm, total_mass, scale_aqua_, scale_sorbed_, this->rho_aqua_);
426 pair<double,double> solution;
430 solution = boost::math::tools::toms748_solve(eq_func, 0.0, upper_solution_bound, toler, max_iter);
432 catch(boost::exception
const & e)
434 THROW( Isotherm::ExcBoostSolver()
435 << EI_BoostMessage(boost::diagnostic_information(e))
440 difference = (solution.second - solution.first)/2;
441 double c_aqua = solution.first + difference;
443 (total_mass - scale_aqua_ * c_aqua)/scale_sorbed_);
449 switch(adsorption_type_)
454 return solve_conc( conc, obj_isotherm);
459 Linear obj_isotherm(mult_coef_);
460 return solve_conc( conc, obj_isotherm);
465 Freundlich obj_isotherm(mult_coef_, second_coef_);
466 return solve_conc( conc, obj_isotherm);
471 Langmuir obj_isotherm(mult_coef_, second_coef_);
472 return solve_conc( conc, obj_isotherm);
483 double mass_limit = scale_aqua_ * table_limit_ + scale_sorbed_ *
const_cast<Func &
>(isotherm)(table_limit_ / this->rho_aqua_);
486 THROW( Isotherm::ExcNegativeTotalMass()
487 << EI_TotalMass(mass_limit)
489 total_mass_step_ = mass_limit / n_steps;
491 for(
int i=0; i<= n_steps; i++) {
493 ConcPair c_pair( mass/scale_aqua_, 0.0 );
495 ConcPair result = solve_conc( c_pair, isotherm);
496 double c_sorbed_rot = ( result.
solid * scale_aqua_ - result.
fluid * scale_sorbed_);
497 interpolation_table.push_back(c_sorbed_rot);
498 mass = mass+total_mass_step_;
double operator()(double x)
Isotherm definition.
double mult_coef_
Multiplication parameter of the isotherm.
CrossFunction(const Func &func_, double total_mass, double scale_aqua, double scale_sorbed, double rho_aqua)
void compute(double &c_aqua, double &c_sorbed)
double second_coef_
Optional secod parameter of the isotherm.
IntersectionPoint * interpolate(const IntersectionPoint &A1, const IntersectionPoint &A2, double t)
ConcPair(double x, double y)
double operator()(double x)
Isotherm definition.
ConcPair solve_conc(ConcPair c_pair, const Func &isotherm)
Langmuir(double mult_coef, double alpha)
Constructor to set parameters.
#define DECLARE_EXCEPTION(ExcName, Format)
Macro for simple definition of exceptions.
double operator()(double x)
Isotherm definition.
double mult_coef_
Parameters of the isotherm.
void reinit(enum SorptionType sorption_type, bool limited_solubility_on, double aqua_density, double scale_aqua, double scale_sorbed, double c_aqua_limit, double mult_coef, double second_coef)
Freundlich(double mult_coef, double exponent)
Constructor to set parameters.
double mult_coef_
Parameters of the isotherm.
vector< double > interpolation_table
double scale_sorbed_
coefficient that convert adsorbed molar concentration to mass; molar_weight * rho_rock * (1 - porosit...
ConcPair precipitate(ConcPair conc)
Linear(double mult_coef)
Constructor to set parameters.
#define TYPEDEF_ERR_INFO(EI_Type, Type)
Macro to simplify declaration of error_info types.
bool limited_solubility_on_
Solubility limit flag.
Pair of soluted and adsorbed concentration.
double operator()(double conc_aqua)
double mult_coef_
Parameters of the isotherm.
SorptionType
Type of adsorption isotherm.
bool operator()(const T &a, const T &b)
void make_table(int n_points)
None()
Constructor to set parameters.
void interpolate(double &c_aqua, double &c_sorbed)
ConcPair compute_projection(ConcPair conc)
double scale_aqua_
coefficient that convert soluted concentration to mass; porosity = k_W, originally rho_aqua*porosity ...
bool is_precomputed(void)
double operator()(double x)
Isotherm definition.
double rho_aqua_
density of the solvent
#define THROW(whole_exception_expr)
Wrapper for throw. Saves the throwing point.