Skip to content
Snippets Groups Projects
Commit 22f58dd4 authored by Timo Koch's avatar Timo Koch
Browse files

[test][richards] Add two benchmark cases

parent b1f87621
No related branches found
No related tags found
1 merge request!2432[test][richards] Add two benchmark cases Vanderborght 2005
Showing
with 1670 additions and 0 deletions
add_subdirectory(benchmarks)
add_subdirectory(implicit)
dune_symlink_to_source_files(FILES
params_evaporation_sand.input
params_evaporation_loam1.input
params_evaporation_loam2.input
params_evaporation_clay.input
params_infiltration_sand.input
params_infiltration_loam.input
params_infiltration_clay.input
run_and_plot_m21.py
run_and_plot_m22.py
)
add_executable(test_richards_benchmark_tpfa EXCLUDE_FROM_ALL main.cc)
dumux_add_test(NAME test_richards_benchmark_infiltration_tpfa
TARGET test_richards_benchmark_tpfa
LABELS porousmediumflow richards
COMMAND ${CMAKE_SOURCE_DIR}/bin/testing/runtest.py
CMD_ARGS --script fuzzyData
--relative 0.05
--files ${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_ana_loam_0.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_ana_loam_0.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_loam_1.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_loam_1.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_loam_2.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_loam_2.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_loam_3.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_loam_3.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_ana_clay_0.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_ana_clay_0.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_clay_1.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_clay_1.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_clay_2.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_clay_2.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_clay_3.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_clay_3.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_ana_sand_0.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_ana_sand_0.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_sand_1.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_sand_1.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_sand_2.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_sand_2.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_theta_deltaeta_num_sand_3.dat
${CMAKE_CURRENT_BINARY_DIR}/theta_deltaeta_num_sand_3.dat
--command "${CMAKE_CURRENT_BINARY_DIR}/run_and_plot_m21.py")
dumux_add_test(NAME test_richards_benchmark_evaporation_tpfa
TARGET test_richards_benchmark_tpfa
LABELS porousmediumflow richards
COMMAND ${CMAKE_SOURCE_DIR}/bin/testing/runtest.py
CMD_ARGS --script fuzzyData
--files ${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_analytical_clay.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_analytical_clay.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_analytical_loam1.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_analytical_loam1.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_analytical_loam2.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_analytical_loam2.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_analytical_sand.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_analytical_sand.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_actual_clay.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_actual_clay.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_actual_loam1.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_actual_loam1.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_actual_loam2.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_actual_loam2.dat
${CMAKE_SOURCE_DIR}/test/references/test_richards_benchmark_tpfa_rate_actual_sand.dat
${CMAKE_CURRENT_BINARY_DIR}/rate_actual_sand.dat
--command "${CMAKE_CURRENT_BINARY_DIR}/run_and_plot_m22.py")
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup RichardsTests
* \brief Helper functions for analytical solutions
*/
#ifndef DUMUX_RICHARDS_BENCHMARKS_ANALYTICAL_HH
#define DUMUX_RICHARDS_BENCHMARKS_ANALYTICAL_HH
#include <cmath>
#include <algorithm>
#include <dumux/io/format.hh>
#include <dumux/io/gnuplotinterface.hh>
#include <dumux/common/math.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/integrate.hh>
#include <dumux/material/fluidmatrixinteractions/2p/vangenuchten.hh>
#include <dumux/material/components/simpleh2o.hh>
namespace Dumux {
/*!
* \file
* \ingroup RichardsTests
* \brief Analytical solution for infiltration scenario
* \note Root-soil benchmark paper Schnepf et al. (case M2.1, Eq. 4) https://doi.org/10.3389/fpls.2020.00316
* \note based on Vanderborght 2005 (see Fig. 4abc and Eq. 56-60) https://doi.org/10.2113/4.1.206
*
* Assumption is that the infiltration front is sufficiently far
* away from the boundary. The boundary conditions for the analytical
* solution are prescribed at z = ±∞.
*/
class AnalyticalSolutionM21
{
using PcKrSwCurve = FluidMatrix::VanGenuchtenNoReg<double>;
public:
AnalyticalSolutionM21()
: pcKrSwCurve_("SpatialParams")
{
gravity_ = 9.81;
density_ = Components::SimpleH2O<double>::liquidDensity(0,0);
viscosity_ = Components::SimpleH2O<double>::liquidViscosity(0,0);
porosity_ = getParam<double>("SpatialParams.Porosity");
permeability_ = getParam<double>("SpatialParams.Permeability");
const auto initialHead = getParam<double>("Problem.InitialHeadInCm", -400.0)*0.01; // cm to m
const auto pwInitial = initialHead*gravity_*density_ + 1.0e5;
const auto pcInitial = 1.0e5 - pwInitial;
const auto swIntial = pcKrSwCurve_.sw(pcInitial);
// compute θ_i, K_i, θ_sur, K_sur, θ_a
waterContentInitial_ = swIntial*porosity_;
conductivityInitial_ = conductivity(waterContentInitial_);
waterContentSurface_ = getParam<double>("Analytical.WaterContentSurface");
conductivitySurface_ = conductivity(waterContentSurface_);
waterContentRef_ = 0.5*(waterContentSurface_ + waterContentInitial_);
std::cout << std::endl << "Computed analytical solution (Vanderborght 2005 4abc, Schnepf 2020 M2.1) with:\n";
const auto soilType = getParam<std::string>("SpatialParams.Name");
std::cout << "Material: " << soilType << "\n";
std::cout << Fmt::format("-> θ_i: {}\n-> θ_sur: {}\n", waterContentInitial_, waterContentSurface_);
std::cout << Fmt::format("-> θ_a: {}\n", waterContentRef_);
std::cout << Fmt::format("-> K_i: {}\n-> K_sur: {}\n", conductivityInitial_, conductivitySurface_);
// create a table for fast evaluation of the inverse
// for the inverse we need to go from high water content to low water content so
// that deltaEta is monotonically increasing which is a requirement for the value key
// for the table interpolation
waterContents_ = linspace(waterContentSurface_ - 1e-8, waterContentInitial_ + 1e-8, 10000);
deltaEtas_.resize(waterContents_.size());
std::transform(waterContents_.begin(), waterContents_.end(), deltaEtas_.begin(),
[this](auto theta){ return deltaEtaExact(theta); });
// plot delta eta over water content
GnuplotInterface<double> gnuplot(false);
gnuplot.setOpenPlotWindow(false);
gnuplot.resetPlot();
gnuplot.setXlabel("water content");
gnuplot.setYlabel("delta eta");
gnuplot.addDataSetToPlot(waterContents_, deltaEtas_, "theta_vs_deltaeta_" + soilType + ".dat", "w l t 'analytical'");
gnuplot.plot("infiltration_theta_vs_deltaeta");
// TODO how to choose these? Is there an exact value?
waterContentRefDepth_ = getParam<double>("Analytical.WaterContentRefDepth");
waterContentRefTime_ = getParam<double>("Analytical.WaterContentRefTime");
std::cout << Fmt::format("-> η_a: {}\n", eta(waterContentRefDepth_, waterContentRefTime_));
}
// analytical solution for given x and t
double waterContent(double x /*in cm*/, double t /*in days*/) const
{
const auto deltaEta = eta(x, t) - eta(waterContentRefDepth_, waterContentRefTime_);
return interpolate<InterpolationPolicy::LinearTable>(deltaEta, deltaEtas_, waterContents_);
}
double deltaEta(double x /*in cm*/, double t /*in days*/) const
{ return eta(x, t) - eta(waterContentRefDepth_, waterContentRefTime_); }
double initialWaterContent() const
{ return waterContentInitial_; }
double surfaceWaterContent() const
{ return waterContentSurface_; }
double refWaterContent() const
{ return waterContentRef_; }
auto waterContentAndDeltaEtaPlot() const
{ return std::tie(waterContents_, deltaEtas_); }
private:
// hydraulic conductivity Kf (including relative permeability) in cm/day
double conductivity(double waterContent) const
{
const auto sw = waterContent/porosity_;
const auto krw = pcKrSwCurve_.krw(sw);
return krw*permeability_*gravity_*density_/viscosity_*100*86400;
}
// D = Kf*dh/dtheta Van Genuchten 1980 Eq. 10 in cm^2/day
double diffusivity(double waterContent) const
{
// pc = -h*gravity_*density_
// h = -pc/(gravity_*density_)
const auto sw = waterContent/porosity_;
const auto dh_dpc = -1.0/(0.01*gravity_*density_);
const auto dpc_dsw = pcKrSwCurve_.dpc_dsw(sw);
const auto dsw_dtheta = 1.0/porosity_;
const auto dh_dtheta = dh_dpc*dpc_dsw*dsw_dtheta;
return conductivity(waterContent)*dh_dtheta;
}
// integral travelling wave (Eq. 4 Schnepf 2020)
double integral(double waterContent) const
{
const auto integrand = [this](auto theta){
return diffusivity(theta)/(
(conductivitySurface_ - conductivityInitial_)*(theta - waterContentInitial_)
- (conductivity(theta) - conductivityInitial_)*(waterContentSurface_ - waterContentInitial_)
);
};
return integrateScalarFunction(integrand, waterContent, waterContentRef_);
}
// relative transformed coordinate position (Eq. 4 Schnepf 2020)
double deltaEtaExact(double waterContent) const
{
return (waterContentSurface_ - waterContentInitial_)*integral(waterContent);
}
// tranformed coordinate
double eta(double x /*in cm*/, double t /*in days*/) const
{
return std::fabs(x) - (conductivitySurface_ - conductivityInitial_)/(waterContentSurface_ - waterContentInitial_)*t;
}
PcKrSwCurve pcKrSwCurve_;
double waterContentSurface_, waterContentInitial_, waterContentRef_;
double waterContentRefDepth_, waterContentRefTime_;
double conductivitySurface_, conductivityInitial_;
double porosity_, permeability_, viscosity_, density_, gravity_;
std::vector<double> waterContents_;
std::vector<double> deltaEtas_;
};
/*!
* \file
* \ingroup RichardsTests
* \brief Analytical solution for evaporation scenario
* \note Root-soil benchmark paper Schnepf et al. (case M2.2) https://doi.org/10.3389/fpls.2020.00316
* \note based on Vanderborght 2005 (see Fig. 5abcd and Eq. 39-47) https://doi.org/10.2113/4.1.206
* \note Derivation in Lockington 1994 https://doi.org/10.1029/93WR03411
*
* The analytical solution is an approximative solution for an inifite domain
* and is neglects graviational effects. Lockington 1994 estimates less than 1% error
* against a given exact solution (without gravity) for a simplified pc-sw relationship.
*/
class AnalyticalSolutionM22
{
using PcKrSwCurve = FluidMatrix::VanGenuchtenNoReg<double>;
public:
AnalyticalSolutionM22()
: pcKrSwCurve_("SpatialParams")
{
gravity_ = 9.81;
density_ = Components::SimpleH2O<double>::liquidDensity(0,0);
viscosity_ = Components::SimpleH2O<double>::liquidViscosity(0,0);
porosity_ = getParam<double>("SpatialParams.Porosity");
permeability_ = getParam<double>("SpatialParams.Permeability");
const auto initialHead = getParam<double>("Problem.InitialHeadInCm", -40.0)*0.01; // cm to m
potentialEvaporationRate_ = getParam<double>("Problem.SurfaceFluxMilliMeterPerDay");
const auto pwInitial = initialHead*gravity_*density_ + 1.0e5;
const auto pcInitial = 1.0e5 - pwInitial;
const auto swIntial = pcKrSwCurve_.sw(pcInitial);
waterContentInitial_ = swIntial*porosity_;
const auto pwSurface = -10000*0.01*gravity_*density_ + 1.0e5;
const auto pcSurface = 1.0e5 - pwSurface;
const auto swSurface = pcKrSwCurve_.sw(pcSurface);
waterContentSurface_ = swSurface*porosity_;
std::cout << std::endl << "Computed analytical solution (Vanderborght 2005 5abcd, Schnepf 2020 M2.2) with:\n";
const auto soilType = getParam<std::string>("SpatialParams.Name");
std::cout << "Material: " << soilType << "\n";
std::cout << Fmt::format("-> θ_i: {}\n-> θ_sur: {}\n", waterContentInitial_, waterContentSurface_);
// water content from "effective" water content, Vanderborght 2005 Eq. 40
const auto theta = [this](auto thetaEff){
return thetaEff*(waterContentInitial_ - waterContentSurface_) + waterContentSurface_;
};
const auto dIntegral = integrateScalarFunction(
[&,this](auto thetaEff){ return diffusivity(theta(thetaEff)); }, 0.0, 1.0
);
std::cout << Fmt::format("D(0), D(1): {}, {}\n", diffusivity(waterContentSurface_), diffusivity(waterContentInitial_));
std::cout << Fmt::format("Kf(0), Kf(1): {}, {}\n", conductivity(waterContentSurface_), conductivity(waterContentInitial_));
std::cout << Fmt::format("-> Integral D(Θ): {}\n", dIntegral);
const auto dThetaIntegral = integrateScalarFunction(
[&,this](auto thetaEff){ return thetaEff*diffusivity(theta(thetaEff)); }, 0.0, 1.0
);
std::cout << Fmt::format("0D(0), 0.5D(0.5) 1D(1): {}, {}, {}\n",
0.0*diffusivity(theta(0.0)), 0.5*diffusivity(theta(0.5)), 1.0*diffusivity(theta(1.0)));
std::cout << Fmt::format("-> Integral Θ D(Θ): {}\n", dThetaIntegral);
// Eq. 43 Vanderborght 2005
const auto beta = dThetaIntegral*dThetaIntegral/(dIntegral*dIntegral);
const auto dThetaFactorIntegral = integrateScalarFunction(
[&,this](auto thetaEff){
const auto weight = (1.0-beta*thetaEff);
return weight*weight*diffusivity(theta(thetaEff));
}, 0.0, 1.0
);
std::cout << Fmt::format("-> Integral (1.0 - βΘ)² D(θ): {}\n", dThetaFactorIntegral);
// Eq. 42 Vanderborght 2005
const auto alpha = dThetaFactorIntegral/dIntegral;
// Eq. 41 Vanderborght 2005
const auto bracketFactor = 1.0 - alpha/((1.0 - beta)*(1.0 - beta));
const auto mu = 3*beta*(1.0 + std::sqrt(1.0 - 14.0/9.0*bracketFactor)) / ( 2*(beta - 1.0)*bracketFactor );
std::cout << Fmt::format("-> µ: {}, α: {}, β: {}\n", mu, alpha, beta);
// Eq. 39 Vanderborght 2005
desorptivity_ = (waterContentInitial_ - waterContentSurface_)*std::sqrt(4.0/mu*dIntegral);
std::cout << Fmt::format("-> Desorptivity, Sw(θ_i, θ_sur): {}\n", desorptivity_);
// Eq. 44 & 45 Vanderborght 2005
tPotential_ = desorptivity_*desorptivity_/(2.0*potentialEvaporationRate_*0.1*potentialEvaporationRate_*0.1);
tPrime_ = 0.5*tPotential_;
std::cout << Fmt::format("-> Critical time (t_pot): {} days\n", tPotential_);
std::cout << std::endl;
}
// Evaporationrate in mm/day (Eq. 46 & 47 Vanderborght 2005)
double evaporationRate(double t /*in days*/) const
{
if (t < tPotential_)
return potentialEvaporationRate_;
else
return 0.5*desorptivity_/std::sqrt(tPrime_ + t - tPotential_)*10;
}
private:
// hydraulic conductivity Kf (including relative permeability) in cm/day
double conductivity(double waterContent) const
{
const auto sw = waterContent/porosity_;
const auto krw = pcKrSwCurve_.krw(sw);
return krw*permeability_*gravity_*density_/viscosity_*100*86400;
}
// D = Kf*dh/dtheta Van Genuchten 1980 Eq. 10 in cm^2/day
double diffusivity(double waterContent) const
{
// pc = -h*gravity_*density_
// h = -pc/(gravity_*density_)
const auto sw = waterContent/porosity_;
const auto dh_dpc = -1.0/(0.01*gravity_*density_);
const auto dpc_dsw = pcKrSwCurve_.dpc_dsw(sw);
const auto dsw_dtheta = 1.0/porosity_;
const auto dh_dtheta = dh_dpc*dpc_dsw*dsw_dtheta;
return conductivity(waterContent)*dh_dtheta;
}
PcKrSwCurve pcKrSwCurve_;
double waterContentSurface_, waterContentInitial_;
double desorptivity_;
double tPrime_, tPotential_;
double potentialEvaporationRate_;
double porosity_, permeability_, viscosity_, density_, gravity_;
};
} // end namespace Dumux
#endif
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup RichardsTests
* \brief Richards benchmark test
*
* Infiltration benchmark:
* Root-soil benchmark paper Schnepf et al. (case M2.1, Eq. 4) https://doi.org/10.3389/fpls.2020.00316
* based on Vanderborght 2005 (see Fig. 4abc and Eq. 56-60) https://doi.org/10.2113/4.1.206
*
* Evaporation benchmark:
* Root-soil benchmark paper Schnepf et al. (case M2.2) https://doi.org/10.3389/fpls.2020.00316
* based on Vanderborght 2005 (see Fig. 5abcd and Eq. 39-47) https://doi.org/10.2113/4.1.206
*/
#include <config.h>
#include <iostream>
#include <dune/common/parallel/mpihelper.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/math.hh>
#include <dumux/linear/linearsolvertraits.hh>
#include <dumux/linear/amgbackend.hh>
#include <dumux/io/format.hh>
#include <dumux/io/gnuplotinterface.hh>
#include <dumux/io/vtkoutputmodule.hh>
#include <dumux/io/grid/gridmanager_yasp.hh>
#include <dumux/porousmediumflow/richards/newtonsolver.hh>
#include <dumux/assembly/fvassembler.hh>
#include <test/porousmediumflow/richards/benchmarks/analytical.hh>
#include <test/porousmediumflow/richards/benchmarks/properties.hh>
int main(int argc, char** argv)
{
using namespace Dumux;
// initialize MPI, finalize is done automatically on exit
const auto& mpiHelper = Dune::MPIHelper::instance(argc, argv);
// parse command line arguments and input file
Parameters::init(argc, argv);
// Choose tpfa scheme
using TypeTag = Properties::TTag::RichardsBenchmarkCCTpfa;
// try to create a grid (from the given grid file or the input file)
GridManager<GetPropType<TypeTag, Properties::Grid>> gridManager;
gridManager.init();
// we compute on the leaf grid view
const auto& leafGridView = gridManager.grid().leafGridView();
using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
auto gridGeometry = std::make_shared<GridGeometry>(leafGridView);
gridGeometry->update();
// the problem (initial and boundary conditions)
using Problem = GetPropType<TypeTag, Properties::Problem>;
auto problem = std::make_shared<Problem>(gridGeometry);
using SolutionVector = GetPropType<TypeTag, Properties::SolutionVector>;
SolutionVector x(gridGeometry->numDofs());
problem->applyInitialSolution(x);
auto xOld = x;
using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
auto gridVariables = std::make_shared<GridVariables>(problem, gridGeometry);
gridVariables->init(x);
// intialize the vtk output module
auto vtkWriter = std::make_shared<VtkOutputModule<GridVariables, SolutionVector>>(*gridVariables, x, problem->name());
vtkWriter->addVolumeVariable([](const auto& volVars){ return volVars.saturation(0)*volVars.porosity(); }, "water content");
vtkWriter->addVolumeVariable([](const auto& volVars){ return volVars.saturation(0); }, "saturation");
vtkWriter->addVolumeVariable([](const auto& volVars){ return volVars.pressure(0); }, "pressure");
using LinearSolver = AMGBiCGSTABBackend<LinearSolverTraits<GridGeometry>>;
auto linearSolver = std::make_shared<LinearSolver>(leafGridView, gridGeometry->dofMapper());
// using Assembler = FVAssembler<TypeTag, DiffMethod::analytic>; // TODO: something is wrong with analytic derivatives!
using Assembler = FVAssembler<TypeTag, DiffMethod::numeric>;
const auto dt = getParam<double>("TimeLoop.DtInitial");
const auto checkPoints = getParam<std::vector<double>>("TimeLoop.TEnd");
const auto tEnd = checkPoints.back();
// instantiate time loop
auto timeLoop = std::make_shared<CheckPointTimeLoop<double>>(0.0, dt, tEnd);
timeLoop->setMaxTimeStepSize(getParam<double>("TimeLoop.MaxTimeStepSize"));
timeLoop->setCheckPoint(checkPoints);
int checkPointCounter = 0;
auto assembler = std::make_shared<Assembler>(problem, gridGeometry, gridVariables, timeLoop, xOld);
using Newton = RichardsNewtonSolver<Assembler, LinearSolver>;
Newton nonLinearSolver(assembler, linearSolver);
// find out which scenario we are running (evaporation/infiltration)
const auto soilType = getParam<std::string>("SpatialParams.Name");
const double rate = getParam<double>("Problem.SurfaceFluxMilliMeterPerDay");
const BenchmarkScenario scenario = (rate > 0) ? BenchmarkScenario::evaporation : BenchmarkScenario::infiltration;
// compute analytical solution for given scenario
std::unique_ptr<AnalyticalSolutionM21> solutionInfiltration;
std::unique_ptr<AnalyticalSolutionM22> solutionEvaporation;
if (scenario == BenchmarkScenario::infiltration)
solutionInfiltration = std::make_unique<AnalyticalSolutionM21>();
else if (scenario == BenchmarkScenario::evaporation)
solutionEvaporation = std::make_unique<AnalyticalSolutionM22>();
///////////////////////////////////////
// for evaporation scenario output
///////////////////////////////////////
const bool openGnuplotWindow = getParam<bool>("Problem.OpenGnuplotWindow", false);
GnuplotInterface<double> gnuplotEvap(openGnuplotWindow);
gnuplotEvap.setOpenPlotWindow(openGnuplotWindow);
std::vector<double> outTime, outActual, outAnalytic;
///////////////////////////////////////
///////////////////////////////////////
///////////////////////////////////////
// for infiltration scenario output
///////////////////////////////////////
const auto updateAnalyticalWaterContent = [&](auto& theta)
{
const auto t = timeLoop->time()/86400.0;
for (const auto& element : elements(leafGridView))
{
const auto eIdx = leafGridView.indexSet().index(element);
const auto pos = element.geometry().center()[GridGeometry::GridView::dimensionworld -1];
theta[eIdx] = solutionInfiltration->waterContent(pos*100, t);
}
};
const auto updateDeltaEtaOutput = [&](auto& deltaEtaNum, auto& thetaNum)
{
const auto t = timeLoop->time()/86400.0;
for (const auto& element : elements(leafGridView))
{
const auto eIdx = leafGridView.indexSet().index(element);
const auto pos = element.geometry().center()[GridGeometry::GridView::dimensionworld -1];
deltaEtaNum[eIdx] = solutionInfiltration->deltaEta(pos*100, t);
auto fvGeometry = localView(*gridGeometry);
auto elemVolVars = localView(gridVariables->curGridVolVars());
fvGeometry.bindElement(element);
elemVolVars.bindElement(element, fvGeometry, x);
for (const auto& scv : scvs(fvGeometry))
thetaNum[eIdx] = elemVolVars[scv].saturation() * elemVolVars[scv].porosity();
}
};
std::vector<double> pos(leafGridView.size(0));
for (const auto& element : elements(leafGridView))
pos[leafGridView.indexSet().index(element)] = element.geometry().center()[GridGeometry::GridView::dimensionworld -1];
GnuplotInterface<double> gnuplotInfilt(openGnuplotWindow);
gnuplotInfilt.setOpenPlotWindow(openGnuplotWindow);
std::vector<double> deltaEtaNum(leafGridView.size(0)), thetaNum(leafGridView.size(0));
std::vector<double> analyticalWaterContent;
if (scenario == BenchmarkScenario::infiltration)
{
analyticalWaterContent.resize(leafGridView.size(0), solutionInfiltration->initialWaterContent());
vtkWriter->addField(analyticalWaterContent, "analytic water content");
// analytical solution
const auto& [theta, deltaEta] = solutionInfiltration->waterContentAndDeltaEtaPlot();
gnuplotInfilt.addDataSetToPlot(theta, deltaEta,
Fmt::format("theta_deltaeta_ana_{}_{:d}.dat", soilType, checkPointCounter), "w l t 'analytic'"
);
const auto yLimits = getParam<std::vector<double>>("Analytical.PlotLimits");
gnuplotInfilt.setYRange(yLimits[0], yLimits[1]);
gnuplotInfilt.setXRange(0.0, 0.5);
}
///////////////////////////////////////
///////////////////////////////////////
// initial solution
vtkWriter->write(0.0);
// time loop
timeLoop->start();
while (!timeLoop->finished())
{
// solve the non-linear system with time step control
nonLinearSolver.solve(x, *timeLoop);
// make the new solution the old solution
xOld = x;
gridVariables->advanceTimeStep();
// advance to the time loop to the next step
timeLoop->advanceTimeStep();
///////////////////////////////////////
// output for infiltration scenario
///////////////////////////////////////
if (scenario == BenchmarkScenario::infiltration && (timeLoop->finished() || timeLoop->isCheckPoint()))
{
++checkPointCounter;
// write analytical solution to ParaView
updateAnalyticalWaterContent(analyticalWaterContent);
vtkWriter->write(timeLoop->time());
updateDeltaEtaOutput(deltaEtaNum, thetaNum);
const auto posRef = interpolate<InterpolationPolicy::LinearTable>(solutionInfiltration->refWaterContent(), thetaNum, pos);
std::cout << Fmt::format("Simulated x_a: {}cm\n", posRef);
gnuplotInfilt.addDataSetToPlot(thetaNum, deltaEtaNum,
Fmt::format("theta_deltaeta_num_{}_{:d}.dat", soilType, checkPointCounter),
Fmt::format("w p t 'DuMux {:.1f} days'", timeLoop->time()/86400.0)
);
gnuplotInfilt.plot(Fmt::format("infiltration_{}_{:d}", soilType, checkPointCounter));
}
///////////////////////////////////////
// output for evaporation scenario
///////////////////////////////////////
if (scenario == BenchmarkScenario::evaporation)
{
if (timeLoop->finished())
vtkWriter->write(timeLoop->time());
const double rate = problem->computeActualRate(x, *gridVariables, false);
const double analyticRate = solutionEvaporation->evaporationRate(timeLoop->time()/86400.0);
std::cout << Fmt::format("Actual rate: {:.5g} (mm/day), analytic rate: {:.5g} (mm/day)\n", rate, analyticRate);
outTime.push_back(timeLoop->time()/86400.0);
outActual.push_back(rate);
outAnalytic.push_back(analyticRate);
gnuplotEvap.resetPlot();
gnuplotEvap.setXRange(0, tEnd/86400.0);
gnuplotEvap.setXlabel("time in days");
gnuplotEvap.setYlabel("evaporation rate in mm/days");
gnuplotEvap.addDataSetToPlot(outTime, outActual, Fmt::format("rate_actual_{}.dat", soilType), "w p t 'DuMux'");
gnuplotEvap.addDataSetToPlot(outTime, outAnalytic, Fmt::format("rate_analytical_{}.dat", soilType), "w l t 'reference'");
gnuplotEvap.plot(Fmt::format("evaporation_{}", soilType));
}
// report statistics of this time step
timeLoop->reportTimeStep();
// set new dt as suggested by the newton solver
timeLoop->setTimeStepSize(nonLinearSolver.suggestTimeStepSize(timeLoop->timeStepSize()));
}
timeLoop->finalize(leafGridView.comm());
nonLinearSolver.report();
// print dumux end message
if (mpiHelper.rank() == 0)
Parameters::print();
return 0;
}
[TimeLoop]
DtInitial = 10 # [s]
TEnd = 518400 # [s]
MaxTimeStepSize = 1200 # [s]
[Grid]
Cells0 = 1000
Grading0 = -0.99
Positions0 = -1.0 0
[Flux]
UpwindWeight = 0.5
[Problem]
Name = test_richards_benchmark_evaporation_clay
EnableGravity = false
InitialHeadInCm = -200.0
CriticalSurfaceHeadInCm = -10000.0
SurfaceFluxMilliMeterPerDay = 3.0
[Newton]
MaxRelativeShift = 1e-6
[Assembly.NumericDifference]
BaseEpsilon = 1e-7
[SpatialParams]
Name = clay
Swr = 0.25
Snr = 0.0
VanGenuchtenPcLowSweThreshold = 1e-6
VanGenuchtenPcHighSweThreshold = 1.0
VanGenuchtenAlpha = 1.01936799e-4
VanGenuchtenN = 1.1
VanGenuchtenL = 0.5
Porosity = 0.40
Permeability = 1.1798241e-13
[TimeLoop]
DtInitial = 10 # [s]
TEnd = 864000 # [s]
MaxTimeStepSize = 7200 # [s]
[Grid]
Cells0 = 1000
Grading0 = -0.99
Positions0 = -1.0 0
[Flux]
UpwindWeight = 0.5
[Problem]
Name = test_richards_benchmark_evaporation_loam1
EnableGravity = false
InitialHeadInCm = -200.0
CriticalSurfaceHeadInCm = -10000.0
SurfaceFluxMilliMeterPerDay = 1.0
[Newton]
MaxRelativeShift = 1e-6
[Assembly.NumericDifference]
BaseEpsilon = 1e-7
[SpatialParams]
Name = loam1
Swr = 0.18604651162
Snr = 0.0
VanGenuchtenPcLowSweThreshold = 1e-5
VanGenuchtenPcHighSweThreshold = 1.0
VanGenuchtenAlpha = 4.07747197e-4
VanGenuchtenN = 1.6
VanGenuchtenL = 0.5
Porosity = 0.43
Permeability = 5.8991203e-13
[TimeLoop]
DtInitial = 10 # [s]
TEnd = 172800 # [s]
MaxTimeStepSize = 7200 # [s]
[Grid]
Cells0 = 1000
Grading0 = -0.99
Positions0 = -1.0 0
[Flux]
UpwindWeight = 0.5
[Problem]
Name = test_richards_benchmark_evaporation_loam2
EnableGravity = false
InitialHeadInCm = -200.0
CriticalSurfaceHeadInCm = -10000.0
SurfaceFluxMilliMeterPerDay = 3.0
[Newton]
MaxRelativeShift = 1e-6
[Assembly.NumericDifference]
BaseEpsilon = 1e-7
[SpatialParams]
Name = loam2
Swr = 0.18604651162
Snr = 0.0
VanGenuchtenPcLowSweThreshold = 1e-5
VanGenuchtenPcHighSweThreshold = 1.0
VanGenuchtenAlpha = 4.07747197e-4
VanGenuchtenN = 1.6
VanGenuchtenL = 0.5
Porosity = 0.43
Permeability = 5.8991203e-13
[TimeLoop]
DtInitial = 1 # [s]
TEnd = 86400 # [s]
MaxTimeStepSize = 360 # [s]
[Grid]
Cells0 = 1000
Grading0 = -0.99
Positions0 = -1.0 0
[Flux]
UpwindWeight = 0.5
[Problem]
Name = test_richards_benchmark_evaporation_sand
EnableGravity = false
InitialHeadInCm = -40.0
CriticalSurfaceHeadInCm = -10000.0
SurfaceFluxMilliMeterPerDay = 1.0
[Newton]
MaxRelativeShift = 1e-6
[Assembly.NumericDifference]
BaseEpsilon = 1e-7
[SpatialParams]
Name = sand
Swr = 0.10465116279
Snr = 0.0
VanGenuchtenPcLowSweThreshold = 1e-5
VanGenuchtenPcHighSweThreshold = 1.0
VanGenuchtenAlpha = 1.52905199e-3
VanGenuchtenN = 3.0
VanGenuchtenL = 0.5
Porosity = 0.43
Permeability = 1.1798241e-11
[TimeLoop]
DtInitial = 1 # [s]
TEnd = 8640 17280 43200 # [s]
MaxTimeStepSize = 40 # [s]
[Grid]
Cells0 = 300
Grading0 = 1.0
Positions0 = -2.0 0
[Analytical]
WaterContentRefDepth = -24.41 # cm determined with 1500 cells
WaterContentRefTime = 0.1 # days
WaterContentSurface = 0.40
PlotLimits = 2.0 -5.0
[Flux]
UpwindWeight = 0.5
[Problem]
Name = test_richards_benchmark_infiltration_clay
EnableGravity = true
InitialHeadInCm = -400.0
CriticalSurfaceHeadInCm = 0
SurfaceFluxMilliMeterPerDay = -1000.0
[Newton]
MaxRelativeShift = 1e-6
[SpatialParams]
Name = clay
Swr = 0.25
Snr = 0.0
VanGenuchtenPcLowSweThreshold = 0.0
VanGenuchtenPcHighSweThreshold = 1.0
VanGenuchtenKrwHighSweThreshold = 0.9999
VanGenuchtenAlpha = 1.01936799e-4
VanGenuchtenN = 1.1
VanGenuchtenL = 0.5
Porosity = 0.40
Permeability = 1.1798241e-13
[TimeLoop]
DtInitial = 1 # [s]
TEnd = 17280 43200 86400 # [s]
MaxTimeStepSize = 200 # [s]
[Grid]
Cells0 = 300
Grading0 = 1.0
Positions0 = -2.0 0
[Analytical]
WaterContentRefDepth = -40.03 # cm determined with 1500 cells
WaterContentRefTime = 0.2 # days
WaterContentSurface = 0.43
PlotLimits = 5.0 -10.0
[Flux]
UpwindWeight = 0.5
[Problem]
Name = test_richards_benchmark_infiltration_loam
EnableGravity = true
InitialHeadInCm = -400.0
CriticalSurfaceHeadInCm = 0
SurfaceFluxMilliMeterPerDay = -1000.0
[Newton]
MaxRelativeShift = 1e-6
[SpatialParams]
Name = loam
Swr = 0.18604651162
Snr = 0.0
VanGenuchtenPcLowSweThreshold = 0.0
VanGenuchtenPcHighSweThreshold = 1.0
VanGenuchtenKrwHighSweThreshold = 1.0
VanGenuchtenAlpha = 4.07747197e-4
VanGenuchtenN = 1.6
VanGenuchtenL = 0.5
Porosity = 0.43
Permeability = 5.8991203e-13
[TimeLoop]
DtInitial = 1 # [s]
TEnd = 8640 17280 25920 # [s]
MaxTimeStepSize = 25 # [s]
[Grid]
Cells0 = 300
Grading0 = 1.0
Positions0 = -2.0 0
[Analytical]
WaterContentRefDepth = -42.78 # cm determined with 1500 cells
WaterContentRefTime = 0.1 # days
WaterContentSurface = 0.2824
PlotLimits = 5.0 -10.0
[Flux]
UpwindWeight = 0.5
[Problem]
Name = test_richards_benchmark_infiltration_sand
EnableGravity = true
InitialHeadInCm = -400.0
CriticalSurfaceHeadInCm = 0
SurfaceFluxMilliMeterPerDay = -1000.0
[Newton]
MaxRelativeShift = 1e-6
[SpatialParams]
Name = sand
Swr = 0.10465116279
Snr = 0.0
VanGenuchtenPcLowSweThreshold = 0.0
VanGenuchtenPcHighSweThreshold = 1.0
VanGenuchtenKrwHighSweThreshold = 1.0
VanGenuchtenAlpha = 1.52905199e-3
VanGenuchtenN = 3.0
VanGenuchtenL = 0.5
Porosity = 0.43
Permeability = 1.1798241e-11
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup RichardsTests
* \brief Richards benchmarks base problem
*
* Infiltration benchmark:
* Root-soil benchmark paper Schnepf et al. (case M2.1, Eq. 4) https://doi.org/10.3389/fpls.2020.00316
* based on Vanderborght 2005 (see Fig. 4abc and Eq. 56-60) https://doi.org/10.2113/4.1.206
*
* Evaporation benchmark:
* Root-soil benchmark paper Schnepf et al. (case M2.2) https://doi.org/10.3389/fpls.2020.00316
* based on Vanderborght 2005 (see Fig. 5abcd and Eq. 39-47) https://doi.org/10.2113/4.1.206
*/
#ifndef DUMUX_RICHARDS_BECHMARK_PROBLEM_HH
#define DUMUX_RICHARDS_BECHMARK_PROBLEM_HH
#include <dumux/common/properties.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/boundarytypes.hh>
#include <dumux/porousmediumflow/problem.hh>
#include <dumux/material/components/simpleh2o.hh>
namespace Dumux {
enum class BenchmarkScenario {
evaporation, infiltration
};
/*!
* \ingroup RichardsTests
* \brief Richards benchmarks base problem
*/
template <class TypeTag>
class RichardsBenchmarkProblem : public PorousMediumFlowProblem<TypeTag>
{
using ParentType = PorousMediumFlowProblem<TypeTag>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
using BoundaryTypes = Dumux::BoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>;
using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
using FVElementGeometry = typename GridGeometry::LocalView;
using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
using GridView = typename GridGeometry::GridView;
using Element = typename GridView::template Codim<0>::Entity;
using GlobalPosition = typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
static constexpr int dimWorld = GridView::dimensionworld;
public:
RichardsBenchmarkProblem(std::shared_ptr<const GridGeometry> gridGeometry)
: ParentType(gridGeometry)
{
name_ = getParam<std::string>("Problem.Name");
const auto density = Components::SimpleH2O<double>::liquidDensity(0,0);
const auto initialHead = getParam<Scalar>("Problem.InitialHeadInCm")*0.01;
const auto criticalHead = getParam<Scalar>("Problem.CriticalSurfaceHeadInCm")*0.01;
initialPressure_ = 1.0e5 + initialHead*9.81*density;
criticalSurfacePressure_ = 1.0e5 + criticalHead*9.81*density;
const auto& fm = this->spatialParams().fluidMatrixInteractionAtPos(0.0);
const auto criticalSaturation = fm.sw(1.0e5 - criticalSurfacePressure_);
criticalSurfaceKrw_ = fm.krw(criticalSaturation);
enableGravity_ = getParam<bool>("Problem.EnableGravity", true);
const auto potentialRate = getParam<Scalar>("Problem.SurfaceFluxMilliMeterPerDay");
potentialRate_ = density*potentialRate/(1000*86400.0); // mm/day -> kg/s/m^2
useKrwAverage_ = getParam<bool>("Problem.UseKrwAverage", false);
bottomDirichlet_ = getParam<bool>("Problem.BottomDirichlet", false);
scenario_ = (potentialRate > 0) ? BenchmarkScenario::evaporation : BenchmarkScenario::infiltration;
surfaceArea_ = (scenario_ == BenchmarkScenario::evaporation) ? 0.1*0.1 : 0.05*0.05;
}
// output name
const std::string& name() const
{ return name_; }
// reference temperature (unused but required)
Scalar temperature() const
{ return 273.15 + 10; };
// reference pressure
Scalar nonwettingReferencePressure() const
{ return 1.0e5; };
// column cross-section area
Scalar extrusionFactorAtPos(const GlobalPosition &globalPos) const
{ return surfaceArea_; }
/*!
* \brief Specifies which kind of boundary condition should be
* used for which equation on a given boundary segment.
*/
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
{
BoundaryTypes bcTypes;
if (onLowerBoundary(globalPos) && !bottomDirichlet_)
bcTypes.setAllNeumann();
else if (onLowerBoundary(globalPos) && bottomDirichlet_)
bcTypes.setAllDirichlet();
else if (onUpperBoundary(globalPos))
bcTypes.setAllNeumann();
else
DUNE_THROW(Dune::InvalidStateException, "Wrong boundary?");
return bcTypes;
}
/*!
* \brief Evaluates the boundary conditions for a Dirichlet boundary segment.
*/
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
{ return initialAtPos(globalPos); }
/*!
* \brief Evaluates the initial values for a control volume
*/
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
{
PrimaryVariables values(0.0);
values[Indices::pressureIdx] = initialPressure_;
values.setState(Indices::bothPhases);
return values;
}
/*!
* \brief Evaluates the boundary conditions for a Neumann boundary segment.
* Negative values mean influx.
*/
template<class ElementVolumeVariables, class ElementFluxVariablesCache>
NumEqVector neumann(const Element& element,
const FVElementGeometry& fvGeometry,
const ElementVolumeVariables& elemVolVars,
const ElementFluxVariablesCache& elemFluxVarsCache,
const SubControlVolumeFace& scvf) const
{
NumEqVector values(0.0);
const auto& globalPos = scvf.ipGlobal();
if (onUpperBoundary(globalPos))
{
const auto& volVars = elemVolVars[scvf.insideScvIdx()];
const auto dist = (fvGeometry.scv(scvf.insideScvIdx()).center() - globalPos).two_norm();
const auto cellPressure = volVars.pressure(0);
const auto density = volVars.density(0);
const auto viscosity = volVars.viscosity(0);
const auto relPerm = volVars.relativePermeability(0);
const auto K = volVars.permeability();
const auto gravity = enableGravity_ ? 9.81 : 0.0;
const auto avgRelPerm = 0.5*(relPerm + criticalSurfaceKrw_);
// kg/m^3 * m^2 * Pa / m / Pa / s = kg/s/m^2
auto criticalRate = density*K/viscosity*((cellPressure - criticalSurfacePressure_)/dist - density*gravity);
if (!std::signbit(criticalRate))
criticalRate *= useKrwAverage_ ? avgRelPerm : relPerm;
if (scenario_ == BenchmarkScenario::evaporation)
values[Indices::conti0EqIdx] = std::min(potentialRate_, criticalRate);
else
values[Indices::conti0EqIdx] = std::max(potentialRate_, criticalRate);
}
// free drainage (gravitational flux) for infiltration scenario
else if (onLowerBoundary(globalPos) && (scenario_ == BenchmarkScenario::infiltration))
{
const auto& volVars = elemVolVars[scvf.insideScvIdx()];
const auto gravity = enableGravity_ ? 9.81 : 0.0;
const auto density = volVars.density(0);
const auto viscosity = volVars.viscosity(0);
const auto relPerm = volVars.relativePermeability(0);
const auto K = volVars.permeability();
values[Indices::conti0EqIdx] = density*K*relPerm/viscosity*(density*gravity);
}
return values;
}
bool onLowerBoundary(const GlobalPosition &globalPos) const
{ return globalPos[dimWorld-1] < this->gridGeometry().bBoxMin()[dimWorld-1] + eps_; }
bool onUpperBoundary(const GlobalPosition &globalPos) const
{ return globalPos[dimWorld-1] > this->gridGeometry().bBoxMax()[dimWorld-1] - eps_; }
//! compute the actual evaporation/infiltration rate
template<class SolutionVector, class GridVariables>
Scalar computeActualRate(const SolutionVector& sol, const GridVariables& gridVars, bool verbose = true) const
{
Scalar rate = 0.0;
auto fvGeometry = localView(this->gridGeometry());
auto elemVolVars = localView(gridVars.curGridVolVars());
for (const auto& element : elements(this->gridGeometry().gridView()))
{
fvGeometry.bindElement(element);
elemVolVars.bindElement(element, fvGeometry, sol);
for (const auto& scvf : scvfs(fvGeometry))
if (scvf.boundary())
rate += this->neumann(element, fvGeometry, elemVolVars, 0.0, scvf)[0];
}
if (verbose)
std::cout << Fmt::format("Actual rate: {:.5g} (mm/day)\n", rate*86400*1000/1000);
return rate*86400*1000/1000;
}
private:
Scalar initialPressure_, criticalSurfacePressure_, potentialRate_;
Scalar criticalSurfaceKrw_;
static constexpr Scalar eps_ = 1.5e-7;
std::string name_;
bool enableGravity_;
bool bottomDirichlet_;
bool useKrwAverage_;
BenchmarkScenario scenario_;
Scalar surfaceArea_;
};
} // end namespace Dumux
#endif
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup RichardsTests
* \brief Common properties for the Richards benchmarks
*/
#ifndef DUMUX_RICHARDS_BENCHMARKS_PROPERTIES_HH
#define DUMUX_RICHARDS_BENCHMARKS_PROPERTIES_HH
#include <dune/grid/yaspgrid.hh>
#include <dumux/common/properties.hh>
#include <dumux/discretization/cctpfa.hh>
#include <dumux/porousmediumflow/richards/model.hh>
#include "spatialparams.hh"
#include "problem.hh"
namespace Dumux::Properties {
// Create new type tags
namespace TTag {
struct RichardsBenchmark { using InheritsFrom = std::tuple<Richards>; };
struct RichardsBenchmarkCCTpfa { using InheritsFrom = std::tuple<RichardsBenchmark, CCTpfaModel>; };
} // end namespace TTag
template<class TypeTag>
struct Grid<TypeTag, TTag::RichardsBenchmark>
{ using type = Dune::YaspGrid<1, Dune::TensorProductCoordinates<double, 1>>; };
// Set the physical problem to be solved
template<class TypeTag>
struct Problem<TypeTag, TTag::RichardsBenchmark>
{ using type = RichardsBenchmarkProblem<TypeTag>; };
// Set the spatial parameters
template<class TypeTag>
struct SpatialParams<TypeTag, TTag::RichardsBenchmark>
{
using type = RichardsBenchmarkSpatialParams<
GetPropType<TypeTag, Properties::GridGeometry>, GetPropType<TypeTag, Properties::Scalar>
>;
};
} // end namespace Dumux::Properties
#endif
#!/usr/bin/env python3
import sys
import numpy as np
import subprocess
# Run benchmark M2.1 infiltration scenarios
subprocess.run(["./test_richards_benchmark_tpfa", "params_infiltration_sand.input"])
subprocess.run(["./test_richards_benchmark_tpfa", "params_infiltration_loam.input"])
subprocess.run(["./test_richards_benchmark_tpfa", "params_infiltration_clay.input"])
try:
import matplotlib.pyplot as plt
except ImportError:
print("Matplotlib not found so no result plots will be created.")
sys.exit(0)
# Create Fig. 4def (right column) of Vanderborght 2005
fig, axes = plt.subplots(nrows=3, ncols=1, figsize=(5, 10), dpi=72)
# set labels
axes[0].set_title(r"$\theta$")
for ax in axes:
ax.set_ylabel(r"$\Delta\eta$ (cm)")
# set limits
axes[0].set_ylim([5.0, -10.0])
axes[0].set_xlim([0.0, 0.5])
axes[1].set_ylim([5.0, -10.0])
axes[1].set_xlim([0.0, 0.5])
axes[2].set_ylim([2.0, -5.0])
axes[2].set_xlim([0.0, 0.5])
# first plot
for i, (t, marker) in enumerate(zip([0.1, 0.2, 0.3], ["x", "D", "o"])):
theta, eta = np.genfromtxt("theta_deltaeta_num_sand_{}.dat".format(i+1)).T
axes[0].plot(theta, eta, linestyle='none', fillstyle='none', marker=marker, label="{:.1f} days".format(t))
theta, eta = np.genfromtxt("theta_deltaeta_ana_sand_0.dat").T
axes[0].plot(theta, eta, "k", label="Benchmark")
axes[0].legend()
# second plot
for i, (t, marker) in enumerate(zip([0.2, 0.5, 1.0], ["x", "D", "o"])):
theta, eta = np.genfromtxt("theta_deltaeta_num_loam_{}.dat".format(i+1)).T
axes[1].plot(theta, eta, linestyle='none', fillstyle='none', marker=marker, label="{:.1f} days".format(t))
theta, eta = np.genfromtxt("theta_deltaeta_ana_loam_0.dat").T
axes[1].plot(theta, eta, "k", label="Benchmark")
axes[1].legend()
# third plot
for i, (t, marker) in enumerate(zip([0.1, 0.2, 0.5], ["x", "D", "o"])):
theta, eta = np.genfromtxt("theta_deltaeta_num_clay_{}.dat".format(i+1)).T
axes[2].plot(theta, eta, linestyle='none', fillstyle='none', marker=marker, label="{:.1f} days".format(t))
theta, eta = np.genfromtxt("theta_deltaeta_ana_clay_0.dat").T
axes[2].plot(theta, eta, "k", label="Benchmark")
axes[2].legend()
fig.tight_layout()
plt.savefig("benchmark_infiltration.png")
#!/usr/bin/env python3
import sys
import numpy as np
import subprocess
# Run benchmark M2.2 evaporation scenarios
subprocess.run(["./test_richards_benchmark_tpfa", "params_evaporation_sand.input"])
subprocess.run(["./test_richards_benchmark_tpfa", "params_evaporation_loam1.input"])
subprocess.run(["./test_richards_benchmark_tpfa", "params_evaporation_loam2.input"])
subprocess.run(["./test_richards_benchmark_tpfa", "params_evaporation_clay.input"])
# Create equally spaced plots to improve test robustness.
# In the case of slightly different time steps this makes sure
# we always compare the same number of data points
for refDataFile in ["rate_analytical_sand.dat", "rate_analytical_loam1.dat", "rate_analytical_loam2.dat", "rate_analytical_clay.dat"]:
tRef, rateRef = np.genfromtxt(refDataFile).T
# sample output data at equally spaced times using linear interpolation
tEquallySpaced = np.linspace(np.min(tRef), np.max(tRef), 50, endpoint=True)
sampledAnaRate = np.interp(tEquallySpaced, tRef, rateRef)
sanitizedOutput = np.vstack((tEquallySpaced, sampledAnaRate)).T
np.savetxt(refDataFile, sanitizedOutput)
numericDataFile = refDataFile.replace("analytical", "actual")
t, rate = np.genfromtxt(numericDataFile).T
sampledRate = np.interp(tEquallySpaced, t, rate)
sanitizedOutput = np.vstack((tEquallySpaced, sampledRate)).T
np.savetxt(numericDataFile, sanitizedOutput)
try:
import matplotlib.pyplot as plt
except ImportError:
print("Matplotlib not found so no result plots will be created.")
sys.exit(0)
# Create Fig. 5abcd of Vanderborght 2005
fig, axes = plt.subplots(nrows=2, ncols=2, figsize=(8, 8))
axes[0,0].set_ylabel(r"$E_\mathrm{act}$ in (cm day${}^{-1})$")
axes[1,0].set_ylabel(r"$E_\mathrm{act}$ in (cm day${}^{-1})$")
axes[1,0].set_xlabel(r"t (days)")
axes[1,1].set_xlabel(r"t (days)")
axes[0,0].set_ylim([0.0, 0.1])
axes[1,0].set_ylim([0.0, 0.3])
plots = {"sand": (0,0), "loam1": (0,1), "loam2": (1,0), "clay": (1,1)}
for soil, index in plots.items():
t, e = np.genfromtxt("rate_actual_{}.dat".format(soil)).T
t, e = t, e*0.1
axes[index].plot(t, e, "k", label="DuMux")
t, e = np.genfromtxt("rate_analytical_{}.dat".format(soil)).T
t, e = t, e*0.1
axes[index].plot(t, e, marker="x", linestyle='none', label="Benchmark")
axes[index].legend()
plt.savefig("benchmark_evaporation.png")
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 3 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup RichardsTests
* \brief Spatial parameters for the Richards benchmarks
*/
#ifndef DUMUX_RICHARDS_BENCHMARKS_SPATIAL_PARAMETERS_HH
#define DUMUX_RICHARDS_BENCHMARKS_SPATIAL_PARAMETERS_HH
#include <dumux/common/parameters.hh>
#include <dumux/material/spatialparams/fv.hh>
#include <dumux/material/fluidmatrixinteractions/2p/vangenuchten.hh>
namespace Dumux {
/*!
* \ingroup RichardsTests
* \brief Spatial parameters for the Richards benchmarks
*/
template<class GridGeometry, class Scalar>
class RichardsBenchmarkSpatialParams
: public FVSpatialParams<GridGeometry, Scalar, RichardsBenchmarkSpatialParams<GridGeometry, Scalar>>
{
using ThisType = RichardsBenchmarkSpatialParams<GridGeometry, Scalar>;
using ParentType = FVSpatialParams<GridGeometry, Scalar, ThisType>;
using GridView = typename GridGeometry::GridView;
using GlobalPosition = typename GridView::template Codim<0>::Geometry::GlobalCoordinate;
using PcKrSwCurve = FluidMatrix::VanGenuchtenDefault<Scalar>;
public:
// export permeability type
using PermeabilityType = Scalar;
RichardsBenchmarkSpatialParams(std::shared_ptr<const GridGeometry> gridGeometry)
: ParentType(gridGeometry)
, pcKrSwCurve_("SpatialParams")
, permeability_(getParam<Scalar>("SpatialParams.Permeability"))
, porosity_(getParam<Scalar>("SpatialParams.Porosity"))
{}
/*!
* \brief Returns the intrinsic permeability tensor [m^2] at a given location
*/
PermeabilityType permeabilityAtPos(const GlobalPosition& globalPos) const
{ return permeability_; }
/*!
* \brief Returns the porosity [-] at a given location
*/
Scalar porosityAtPos(const GlobalPosition& globalPos) const
{ return porosity_; }
/*!
* \brief Returns the fluid-matrix interaction law for the sub-control volume
*/
auto fluidMatrixInteractionAtPos(const GlobalPosition& globalPos) const
{ return makeFluidMatrixInteraction(pcKrSwCurve_); }
private:
const PcKrSwCurve pcKrSwCurve_;
const Scalar permeability_, porosity_;
};
} // end namespace Dumux
#endif
1.157410000000000020e-04 3.000000000000000000e+00
1.225623585306122448e-01 3.000000000000000000e+00
2.450089760612244871e-01 3.000000000000000000e+00
3.674555935918367156e-01 3.000000000000000000e+00
4.899022111224489717e-01 3.000000000000000000e+00
6.123488286530611724e-01 3.000000000000000000e+00
7.347954461836734286e-01 3.000000000000000000e+00
8.572420637142856847e-01 2.356547345873095889e+00
9.796886812448979409e-01 2.009320758820943187e+00
1.102135298775509975e+00 1.813572569630184761e+00
1.224581916306122231e+00 1.674257710275444344e+00
1.347028533836734487e+00 1.566283578193238490e+00
1.469475151367346744e+00 1.478537553819808670e+00
1.591921768897959000e+00 1.405009279854189774e+00
1.714368386428571256e+00 1.342038502374267361e+00
1.836815003959183512e+00 1.287205304322299293e+00
1.959261621489795768e+00 1.238836683591719279e+00
2.081708239020408246e+00 1.195719185054402001e+00
2.204154856551020281e+00 1.156936192438849043e+00
2.326601474081632759e+00 1.121785504213962170e+00
2.449048091612244793e+00 1.089717977096938739e+00
2.571494709142857271e+00 1.060307494591189670e+00
2.693941326673469305e+00 1.033201638573059356e+00
2.816387944204081784e+00 1.008105194632021195e+00
2.938834561734693818e+00 9.847791777775084521e-01
3.061281179265306296e+00 9.630198236853880767e-01
3.183727796795918330e+00 9.426579647906333514e-01
3.306174414326530808e+00 9.235508530607631217e-01
3.428621031857142842e+00 9.055723904843258065e-01
3.551067649387755321e+00 8.886140918390854626e-01
3.673514266918367355e+00 8.725820234764813943e-01
3.795960884448979833e+00 8.573946521632799245e-01
3.918407501979591867e+00 8.429789740355959626e-01
4.040854119510203901e+00 8.292733565488841219e-01
4.163300737040816379e+00 8.162193004681020936e-01
4.285747354571428858e+00 8.037664485738968967e-01
4.408193972102040448e+00 7.918695582987409498e-01
4.530640589632652926e+00 7.804893870751060669e-01
4.653087207163265404e+00 7.695880387120920263e-01
4.775533824693877882e+00 7.591345019352786494e-01
4.897980442224489472e+00 7.490973992825950001e-01
5.020427059755101951e+00 7.394496605387167465e-01
5.142873677285714429e+00 7.301672105861153339e-01
5.265320294816326907e+00 7.212267210992168165e-01
5.387766912346938497e+00 7.126083635941140582e-01
5.510213529877550975e+00 7.042931380708070677e-01
5.632660147408163454e+00 6.962637162037216276e-01
5.755106764938775932e+00 6.885030113987010703e-01
5.877553382469387522e+00 6.809972017462628813e-01
6.000000000000000000e+00 6.737319999999999975e-01
1.157410000000000020e-04 1.000000000000000000e+00
2.041950115918367281e-01 1.000000000000000000e+00
4.082742821836734537e-01 1.000000000000000000e+00
6.123535527755101793e-01 1.000000000000000000e+00
8.164328233673469049e-01 1.000000000000000000e+00
1.020512093959183408e+00 1.000000000000000000e+00
1.224591364551020245e+00 1.000000000000000000e+00
1.428670635142857082e+00 7.778737875101386923e-01
1.632749905734693696e+00 6.728506609907902725e-01
1.836829176326530311e+00 6.082763028448093801e-01
2.040908446918367147e+00 5.621196125119106712e-01
2.244987717510204206e+00 5.257719936418311368e-01
2.449066988102040821e+00 4.962841255973509114e-01
2.653146258693877435e+00 4.714697127218852035e-01
2.857225529285714494e+00 4.502224935268295813e-01
3.061304799877551108e+00 4.317214223025953035e-01
3.265384070469387723e+00 4.153780819726516360e-01
3.469463341061224337e+00 4.008384118433797050e-01
3.673542611653060952e+00 3.877243180427907832e-01
3.877621882244898011e+00 3.758831851860980344e-01
4.081701152836734181e+00 3.650521808992008843e-01
4.285780423428571240e+00 3.551461893886338483e-01
4.489859694020408298e+00 3.460000093493642570e-01
4.693938964612244469e+00 3.375420541391472296e-01
4.898018235204081527e+00 3.296805789347614923e-01
5.102097505795918586e+00 3.223473319738066056e-01
5.306176776387754757e+00 3.154925970025837523e-01
5.510256046979591815e+00 3.090526926870446611e-01
5.714335317571428874e+00 3.030052356350920184e-01
5.918414588163265044e+00 2.972907827085818355e-01
6.122493858755102103e+00 2.919015811367194901e-01
6.326573129346938273e+00 2.867895227916334444e-01
6.530652399938775332e+00 2.819426499138086295e-01
6.734731670530612391e+00 2.773322473746159256e-01
6.938810941122448561e+00 2.729425771831885950e-01
7.142890211714285620e+00 2.687571204268193736e-01
7.346969482306121790e+00 2.647572017408195144e-01
7.551048752897958849e+00 2.609344305390497998e-01
7.755128023489795908e+00 2.572691248037059508e-01
7.959207294081632078e+00 2.537594278750826549e-01
8.163286564673468249e+00 2.503861091876726230e-01
8.367365835265305307e+00 2.471481366054413820e-01
8.571445105857142366e+00 2.440301818330647687e-01
8.775524376448979424e+00 2.410285468260425223e-01
8.979603647040816483e+00 2.381357633170685106e-01
9.183682917632653542e+00 2.353439878780564654e-01
9.387762188224488824e+00 2.326496473077268867e-01
9.591841458816325883e+00 2.300438760865215204e-01
9.795920729408162941e+00 2.275270441838346258e-01
1.000000000000000000e+01 2.250880000000000103e-01
1.157410000000000020e-04 3.000000000000000000e+00
4.092970546938776155e-02 3.000000000000000000e+00
8.174366993877552057e-02 3.000000000000000000e+00
1.225576344081632796e-01 3.000000000000000000e+00
1.633715988775510386e-01 2.271078521268931283e+00
2.041855633469387976e-01 1.845019560127752589e+00
2.449995278163265566e-01 1.620535773311463812e+00
2.858134922857142879e-01 1.469120155903942182e+00
3.266274567551020747e-01 1.364575394818577925e+00
3.674414212244898614e-01 1.269769535359806634e+00
4.082553856938775927e-01 1.199668822142912372e+00
4.490693501632653239e-01 1.132794855464952288e+00
4.898833146326531107e-01 1.081247230721113395e+00
5.306972791020408975e-01 1.031245389008651481e+00
5.715112435714285732e-01 9.913977398273021713e-01
6.123252080408163600e-01 9.522653823693737474e-01
6.531391725102041468e-01 9.203353474538268353e-01
6.939531369795919336e-01 8.886788876467560661e-01
7.347671014489797203e-01 8.623820076907375842e-01
7.755810659183673961e-01 8.361133308158121835e-01
8.163950303877551828e-01 8.139860330974717506e-01
8.572089948571429696e-01 7.918587353791312067e-01
8.980229593265306454e-01 7.728048946406237407e-01
9.388369237959184321e-01 7.538625681935013922e-01
9.796508882653062189e-01 7.372228182778443051e-01
1.020464852734693784e+00 7.207745883121543518e-01
1.061278817204081681e+00 7.060874955694643740e-01
1.102092781673469357e+00 6.916370798075467397e-01
1.142906746142857033e+00 6.785508611004817947e-01
1.183720710612244931e+00 6.657252943946687651e-01
1.224534675081632606e+00 6.539724509949772502e-01
1.265348639551020504e+00 6.424913500574166436e-01
1.306162604020408180e+00 6.318630305784083134e-01
1.346976568489795856e+00 6.215096825720809459e-01
1.387790532959183754e+00 6.118375969276407256e-01
1.428604497428571429e+00 6.024380914418371313e-01
1.469418461897959327e+00 5.935883457413171360e-01
1.510232426367347003e+00 5.850062950340975121e-01
1.551046390836734679e+00 5.768703738342368936e-01
1.591860355306122576e+00 5.689955412368316034e-01
1.632674319775510252e+00 5.614821234768757607e-01
1.673488284244897928e+00 5.542210287539828251e-01
1.714302248714285826e+00 5.472549995037219173e-01
1.755116213183673501e+00 5.405331604496251829e-01
1.795930177653061177e+00 5.340523402275998421e-01
1.836744142122449075e+00 5.278073237299559795e-01
1.877558106591836751e+00 5.217576565777397590e-01
1.918372071061224649e+00 5.159340871844718679e-01
1.959186035530612324e+00 5.102303896977580244e-01
2.000000000000000000e+00 5.046910000000000007e-01
1.157409999999999986e-05 1.000000000000000000e+00
2.041950115918367628e-02 3.295934362560769548e-01
4.082742821836735092e-02 2.252023654273773157e-01
6.123535527755102209e-02 1.812868377206766635e-01
8.164328233673469326e-02 1.558101310769424153e-01
1.020512093959183714e-01 1.387101880406292065e-01
1.224591364551020356e-01 1.262144359367440016e-01
1.428670635142857137e-01 1.165750720153176379e-01
1.632749905734693918e-01 1.088486914930787558e-01
1.836829176326530699e-01 1.024775306697169458e-01
2.040908446918367480e-01 9.710652280823442450e-02
2.244987717510204261e-01 9.249928496689792390e-02
2.449066988102040765e-01 8.849074952788624215e-02
2.653146258693878101e-01 8.496135964228500315e-02
2.857225529285714605e-01 8.182353441516870829e-02
3.061304799877551663e-01 7.900995070571752155e-02
3.265384070469388167e-01 7.646753624184916831e-02
3.469463341061224670e-01 7.415564181085756990e-02
3.673542611653061729e-01 7.204129322196131668e-02
3.877621882244898233e-01 7.009792845706301623e-02
4.081701152836735291e-01 6.830373612044857157e-02
4.285780423428571795e-01 6.664047525707428310e-02
4.489859694020408853e-01 6.509300225184083688e-02
4.693938964612245357e-01 6.364861068951874201e-02
4.898018235204081861e-01 6.229645695649491111e-02
5.102097505795918364e-01 6.102700401554703863e-02
5.306176776387755423e-01 5.983209111130089602e-02
5.510256046979591371e-01 5.870460149931776977e-02
5.714335317571428430e-01 5.763860743491050487e-02
5.918414588163265488e-01 5.662850661421763715e-02
6.122493858755102547e-01 5.566977970465630809e-02
6.326573129346938495e-01 5.475811460094261868e-02
6.530652399938775554e-01 5.388975863332265415e-02
6.734731670530612613e-01 5.306141255974963777e-02
6.938810941122448561e-01 5.227036887046056929e-02
7.142890211714285620e-01 5.151356661687750416e-02
7.346969482306122678e-01 5.078872530286561437e-02
7.551048752897959737e-01 5.009359084705073545e-02
7.755128023489795686e-01 4.942622865044665748e-02
7.959207294081632744e-01 4.878492419993582518e-02
8.163286564673469803e-01 4.816786308901593178e-02
8.367365835265305751e-01 4.757363260788292914e-02
8.571445105857142810e-01 4.700086177098793117e-02
8.775524376448979869e-01 4.644830267191920464e-02
8.979603647040816927e-01 4.591480203246568054e-02
9.183682917632652876e-01 4.539928969142054127e-02
9.387762188224489934e-01 4.490075197359985615e-02
9.591841458816326993e-01 4.441829415844209678e-02
9.795920729408162941e-01 4.395102290637325354e-02
1.000000000000000000e+00 4.349809999999999788e-02
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment