Test for the one-phase facet coupling model. More...
#include <config.h>
#include <iostream>
#include <dune/common/parallel/mpihelper.hh>
#include <dune/geometry/quadraturerules.hh>
#include "problem_1p_bulk.hh"
#include "problem_1p_lowdim.hh"
#include "problem_tracer_bulk.hh"
#include "problem_tracer_lowdim.hh"
#include <dumux/common/properties.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/dumuxmessage.hh>
#include <dumux/assembly/diffmethod.hh>
#include <dumux/discretization/method.hh>
#include <dumux/discretization/elementsolution.hh>
#include <dumux/discretization/evalgradients.hh>
#include <dumux/linear/seqsolverbackend.hh>
#include <dumux/multidomain/newtonsolver.hh>
#include <dumux/multidomain/fvassembler.hh>
#include <dumux/multidomain/traits.hh>
#include <dumux/multidomain/facet/gridmanager.hh>
#include <dumux/multidomain/facet/couplingmapper.hh>
#include <dumux/multidomain/facet/couplingmanager.hh>
#include <dumux/multidomain/facet/codimonegridadapter.hh>
#include <dumux/io/vtkoutputmodule.hh>
Test for the one-phase facet coupling model.
Namespaces | |
namespace | Dumux |
make the local view function available whenever we use the grid geometry | |
namespace | Dumux::Properties |
Typedefs | |
using | Dumux::Properties::BoxTracerTraits = TestTraits< Properties::TTag::TracerBulkBox, Properties::TTag::TracerLowDimBox > |
using | Dumux::Properties::TpfaTracerTraits = TestTraits< Properties::TTag::TracerBulkTpfa, Properties::TTag::TracerLowDimTpfa > |
using | Dumux::Properties::MpfaTracerTraits = TestTraits< Properties::TTag::TracerBulkMpfa, Properties::TTag::TracerLowDimMpfa > |
Functions | |
template<class GridGeometry , class GridManager , class LowDimGridView , std::enable_if_t< GridGeometry::discMethod==Dumux::DiscretizationMethod::box, int > = 0> | |
void | updateBulkFVGridGeometry (GridGeometry &gridGeometry, const GridManager &gridManager, const LowDimGridView &lowDimGridView) |
Updates the finite volume grid geometry for the box scheme. More... | |
template<class FV , class Storage , class CM , class Assembler , class Prob , class GV , class Sol , std::size_t id> | |
void | computeVolumeFluxes (Storage &volumeFluxes, CM &couplingManager, const Assembler &assembler, const Prob &problem, const typename GV::GridGeometry &gridGeometry, const GV &gridVariables, const Sol &sol, Dune::index_constant< id > domainId) |
Computes the volume fluxes on all scvfs for a sub-domain. More... | |
int | main (int argc, char **argv) |
void computeVolumeFluxes | ( | Storage & | volumeFluxes, |
CM & | couplingManager, | ||
const Assembler & | assembler, | ||
const Prob & | problem, | ||
const typename GV::GridGeometry & | gridGeometry, | ||
const GV & | gridVariables, | ||
const Sol & | sol, | ||
Dune::index_constant< id > | domainId | ||
) |
Computes the volume fluxes on all scvfs for a sub-domain.
int main | ( | int | argc, |
char ** | argv | ||
) |
the problem (initial and boundary conditions)
the linear system
get some time loop parameters
instantiate time loop
set some check points for the vtk output
start the time loop
void updateBulkFVGridGeometry | ( | GridGeometry & | gridGeometry, |
const GridManager & | gridManager, | ||
const LowDimGridView & | lowDimGridView | ||
) |
Updates the finite volume grid geometry for the box scheme.
Updates the finite volume grid geometry for the cell centered schemes.
This is necessary as the finite volume grid geometry for the box scheme with facet coupling requires additional data for the update. The reason is that we have to create additional faces on interior boundaries, which wouldn't be created in the standard scheme.