24#ifndef DUMUX_INCOMPRESSIBLE_ONEP_CONVERGENCETEST_SOLVER_HH
25#define DUMUX_INCOMPRESSIBLE_ONEP_CONVERGENCETEST_SOLVER_HH
30#include <dune/common/timer.hh>
31#include <dune/grid/io/file/dgfparser/dgfexception.hh>
32#include <dune/grid/io/file/vtk.hh>
50template<
class TypeTag>
69template<
class TypeTag>
72 using namespace Dumux;
75 const auto c = std::to_string(numCells);
76 Parameters::init([&c] (
auto& tree) { tree[
"Grid.Cells"] = c +
" " + c; });
85 const auto& leafGridView = gridManager.grid().leafGridView();
92 auto gridGeometry = std::make_shared<GridGeometry>(leafGridView);
93 gridGeometry->update();
97 auto problem = std::make_shared<Problem>(gridGeometry);
101 auto x = std::make_shared<SolutionVector>(gridGeometry->numDofs());
105 auto gridVariables = std::make_shared<GridVariables>(problem, gridGeometry);
106 gridVariables->init(*x);
110 auto assembler = std::make_shared<Assembler>(problem, gridGeometry, gridVariables);
113 auto linearSolver = std::make_shared<LinearSolver>();
120 if (getParam<bool>(
"IO.WriteVTK",
false))
124 IOFields::initOutputModule(vtkWriter);
128 std::vector<Scalar> exact(gridGeometry->numDofs());
130 for (
const auto& element : elements(gridGeometry->gridView()))
132 auto fvGeometry =
localView(*gridGeometry);
133 fvGeometry.bindElement(element);
135 for (
const auto& scv : scvs(fvGeometry))
136 exact[scv.dofIndex()] = problem->exact(scv.dofPosition());
139 vtkWriter.
addField(exact,
"p_exact");
140 vtkWriter.
write(1.0);
144 storage.gridGeometry = gridGeometry;
145 storage.solution = x;
An enum class to define various differentiation methods available in order to compute the derivatives...
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Dumux sequential linear solver backends.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:52
The grid manager (this is the class used by the user) for all supported grid managers that constructs...
Definition: gridmanager_base.hh:312
A VTK output module to simplify writing dumux simulation data to VTK format.
Definition: io/vtkoutputmodule.hh:66
void write(double time, Dune::VTK::OutputType type=Dune::VTK::ascii)
Writing data.
Definition: io/vtkoutputmodule.hh:210
void addField(const Vector &v, const std::string &name, FieldType fieldType=FieldType::automatic)
Definition: io/vtkoutputmodule.hh:161
An implementation of a linear PDE solver.
Definition: linear/pdesolver.hh:62
void solve(SolutionVector &uCurrentIter) override
Solve a linear PDE system.
Definition: linear/pdesolver.hh:88
Sequential ILU(0)-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:580
Base class for linear solvers.
Definition: dumux/linear/solver.hh:37
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:52
std::shared_ptr< SolutionVector > solution
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:60
Dumux::GridManager< Grid > gridManager
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:58
Dumux::GetPropType< TypeTag, Dumux::Properties::GridGeometry > GridGeometry
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:54
std::shared_ptr< GridGeometry > gridGeometry
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:59
Dumux::GetPropType< TypeTag, Dumux::Properties::Grid > Grid
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:53
Dumux::GetPropType< TypeTag, Dumux::Properties::SolutionVector > SolutionVector
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:55
A linear system assembler (residual and Jacobian) for finite volume schemes.
A high-level solver interface for a linear PDE solver.
Declares all properties used in Dumux.
Definition of a problem for the MaxwellStefan problem: A rotating velocity field mixes a MaxwellStefa...
Convience header that includes all grid manager specializations.
A VTK output module to simplify writing dumux simulation data to VTK format.
SolutionStorage< TypeTag > solveRefinementLevel(int numCells)
Solves the problem for a given number of cells per direction.
Definition: test/porousmediumflow/1p/implicit/convergence/solver.hh:70