12#ifndef DUMUX_IO_STRUCTURED_LATTICE_GRID_CREATOR_HH 
   13#define DUMUX_IO_STRUCTURED_LATTICE_GRID_CREATOR_HH 
   22#include <dune/common/concept.hh> 
   23#include <dune/common/exceptions.hh> 
   24#include <dune/grid/common/gridfactory.hh> 
   25#include <dune/grid/yaspgrid.hh> 
   26#include <dune/geometry/referenceelements.hh> 
   29#include <dune/foamgrid/foamgrid.hh> 
   30#include <dune/foamgrid/dgffoam.hh> 
   44template<
class GlobalPosition>
 
   45struct LowDimElementSelector
 
   48    auto require(F&& f) -> 
decltype(
 
   49        bool(f(std::declval<const GlobalPosition&>(), std::declval<const GlobalPosition&>()))
 
   59class StructuredLatticeGridCreator
 
   61    using GridType = Dune::FoamGrid<1, dimWorld>;
 
   62    using Element = 
typename GridType::template Codim<0>::Entity;
 
   63    using GlobalPosition = 
typename Element::Geometry::GlobalCoordinate;
 
   64    using CoordScalar = 
typename GridType::ctype;
 
   65    using HostGrid = Dune::YaspGrid<dimWorld, Dune::TensorProductCoordinates<CoordScalar, dimWorld>>;
 
   66    using HostGridView = 
typename HostGrid::LeafGridView;
 
   67    using ReferenceElements = 
typename Dune::ReferenceElements<CoordScalar, dimWorld>;
 
   69    using VertexIndexPair = std::pair<unsigned int, unsigned int>;
 
   73        VertexIndexPair localVertexIndices;
 
   74        unsigned int directionNumber;
 
   78    using Grid = GridType;
 
   80    void init(
const std::string& paramGroup = 
"")
 
   82        auto useAllElements = [](
const GlobalPosition& a, 
const GlobalPosition& b){ 
return true; };
 
   83        init(useAllElements, paramGroup);
 
   86    template<
class LowDimElementSelector, 
 
   87             typename std::enable_if_t<Dune::models<Concept::LowDimElementSelector<GlobalPosition>, LowDimElementSelector>(), 
int> = 0>
 
   88    void init(
const LowDimElementSelector& lowDimElementSelector,
 
   89              const std::string& paramGroup = 
"")
 
   91        paramGroup_ = paramGroup;
 
   92        removeThroatsOnBoundary_ = getParamFromGroup<std::vector<std::size_t>>(paramGroup, 
"Grid.RemoveThroatsOnBoundary",
 
   93                                                     std::vector<std::size_t>());
 
   95        setElementSelector_(lowDimElementSelector);
 
   96        initRandomNumberGenerator_();
 
   99        const auto hostGridInputData = getHostGridInputData_();
 
  100        using HostGrid = Dune::YaspGrid<dimWorld, Dune::TensorProductCoordinates<CoordScalar, dimWorld>>;
 
  101        using HastGridManager = GridManager<HostGrid>;
 
  102        HastGridManager hostGridManager;
 
  103        hostGridManager.init(hostGridInputData.positions, hostGridInputData.cells, hostGridInputData.grading, paramGroup_);
 
  104        hostGridView_ = std::make_unique<HostGridView>(hostGridManager.grid().leafGridView());
 
  106        hostGridLowerLeft_ = hostGridInputData.lowerLeft;
 
  107        hostGridUpperRight_ = hostGridInputData.upperRight;
 
  109        convertHostGridToLowDimGrid_();
 
  117        if (Dune::MPIHelper::getCommunication().size() > 1)
 
  118            gridPtr_->loadBalance();
 
  125    { 
return *gridPtr(); }
 
  130    std::shared_ptr<Grid>& gridPtr()
 
  135    template<
class LowDimElementSelector>
 
  136    void setElementSelector_(
const LowDimElementSelector& lowDimElementSelector)
 
  138        if (!removeThroatsOnBoundary_.empty())
 
  140            auto finalSelector = [&](
const GlobalPosition& a, 
const GlobalPosition& b)
 
  142                static const auto lambdaToRemoveThroatsOnBoundary = getLambdaToRemoveThroatsOnBoundary_();
 
  144                return min(lambdaToRemoveThroatsOnBoundary(a, b), lowDimElementSelector(a, b));
 
  146            considerLowDimElement_ = finalSelector;
 
  149            considerLowDimElement_ = lowDimElementSelector;
 
  152    void initRandomNumberGenerator_()
 
  154        directionProbability_ = getDirectionProbability_();
 
  158            const auto seed = getParamFromGroup<std::size_t>(paramGroup_, 
"Grid.DeletionRandomNumberSeed");
 
  159            generator_.seed(seed);
 
  163            std::random_device rd;
 
  164            generator_.seed(rd());
 
  168    void convertHostGridToLowDimGrid_()
 
  172        for (
const auto& element : elements(*hostGridView_))
 
  174            const auto geometry = 
element.geometry();
 
  175            const auto refElement = ReferenceElements::general(geometry.type());
 
  176            treatEdges_(element, refElement, geometry);
 
  177            treatDiagonals_(element, refElement, geometry);
 
  179        if (elementSet_.empty())
 
  180            DUNE_THROW(Dune::GridError, 
"Trying to create pore network with zero elements!");
 
  183        Dune::GridFactory<Grid> factory;
 
  184        for (
auto&& vertex : vertexSet_)
 
  185            factory.insertVertex(vertex);
 
  186        for (
auto&& element : elementSet_)
 
  187            factory.insertElement(Dune::GeometryTypes::cube(1), {
element.first, 
element.second});
 
  189        gridPtr_ = std::shared_ptr<Grid>(factory.createGrid());
 
  192    void resizeContainers_()
 
  194        vertexInserted_.resize(hostGridView_->size(dimWorld), 
false);
 
  195        hostVertexIdxToVertexIdx_.resize(hostGridView_->size(dimWorld));
 
  196        edgeInserted_.resize(hostGridView_->size(dimWorld-1), 
false);
 
  197        faceInserted_.resize(hostGridView_->size(dimWorld-2), 
false);
 
  200    template<
class HostGr
idElement>
 
  201    bool maybeInsertElementAndVertices_(
const HostGridElement& hostGridElement,
 
  202                                        const typename HostGridElement::Geometry& hostGridElementGeometry,
 
  203                                        std::size_t localVertexIdx0,
 
  204                                        std::size_t localVertexIdx1,
 
  205                                        double directionProbability)
 
  207        if (neglectElement_(hostGridElementGeometry, localVertexIdx0, localVertexIdx1, directionProbability))
 
  211            insertElementAndVertices_(hostGridElement, hostGridElementGeometry, localVertexIdx0, localVertexIdx1, directionProbability);
 
  216    template<
class HostGr
idElement>
 
  217    void insertElementAndVertices_(
const HostGridElement& hostGridElement,
 
  218                                   const typename HostGridElement::Geometry& hostGridElementGeometry,
 
  219                                   std::size_t localVertexIdx0,
 
  220                                   std::size_t localVertexIdx1,
 
  221                                   double directionProbability)
 
  224        const auto vIdxGlobal0 = hostGridView_->indexSet().subIndex(hostGridElement, localVertexIdx0, dimWorld);
 
  225        const auto vIdxGlobal1 = hostGridView_->indexSet().subIndex(hostGridElement, localVertexIdx1, dimWorld);
 
  227        auto getGobalVertexIdx = [&](
auto globalIdx, 
auto localIdx)
 
  229            if (!vertexInserted_[globalIdx])
 
  231                vertexInserted_[globalIdx] = 
true;
 
  232                hostVertexIdxToVertexIdx_[globalIdx] = vertexSet_.size();
 
  233                vertexSet_.push_back(hostGridElementGeometry.corner(localIdx));
 
  236            return hostVertexIdxToVertexIdx_[globalIdx];
 
  239        const auto newVertexIdx0 = getGobalVertexIdx(vIdxGlobal0, localVertexIdx0);
 
  240        const auto newVertexIdx1 = getGobalVertexIdx(vIdxGlobal1, localVertexIdx1);
 
  241        elementSet_.emplace_back(newVertexIdx0, newVertexIdx1);
 
  244    auto getDirectionProbability_()
 const 
  246        std::array<double, numDirections_()> directionProbability;
 
  247        directionProbability.fill(0.0); 
 
  250        if (getParamFromGroup<bool>(paramGroup_, 
"Grid.RegularLattice", 
false))
 
  252            directionProbability.fill(1.0); 
 
  253            for (
int i = 0; i < dimWorld; ++i)
 
  254                directionProbability[i] = 0.0; 
 
  256            return directionProbability;
 
  264                directionProbability = getParamFromGroup<decltype(directionProbability)>(paramGroup_, 
"Grid.DeletionProbability");
 
  266            catch(Dumux::ParameterException &e)
 
  268                throwDirectionError_();
 
  272        return directionProbability;
 
  275    void throwDirectionError_()
 const 
  278        Dune::FieldVector<std::string, numDirections_()> directions;
 
  282            directions[0] = 
"1: (1, 0)\n";
 
  283            directions[1] = 
"2: (0, 1)\n";
 
  285            directions[2] = 
"3: (1, 1)\n";
 
  286            directions[3] = 
"4: (1, -1)\n";
 
  291            directions[0] = 
" 1: (1, 0, 0)\n";
 
  292            directions[1] = 
"2: (0, 1, 0)\n";
 
  293            directions[2] = 
"3: (0, 0, 1)\n";
 
  295            directions[3] = 
"4: (1, 1, 0)\n";
 
  296            directions[4] = 
"5: (1, -1, 0)\n";
 
  297            directions[5] = 
"6: (1, 0, 1)\n";
 
  298            directions[6] = 
"7: (1, 0, -1)\n";
 
  299            directions[7] = 
"8: (0, 1, 1)\n";
 
  300            directions[8] = 
"9: (0, 1, -1)\n";
 
  302            directions[9] = 
"10: (1, 1, 1)\n";
 
  303            directions[10] = 
"11: (1, 1, -1)\n";
 
  304            directions[11] = 
"12: (-1, 1, 1)\n";
 
  305            directions[12] = 
"13: (-1, -1, 1)\n";
 
  307        DUNE_THROW(ParameterException, 
"You must specify probabilities for all directions (" << numDirections_() << 
") \n" << directions << 
"\nExample (3D):\n\n" 
  308        << 
"DeletionProbability = 0.5 0.5 0 0 0 0 0 0 0 0 0 0 0 \n\n" 
  309        << 
"deletes approximately 50% of all throats in x and y direction, while no deletion in any other direction takes place.\n" );
 
  312    static constexpr std::size_t numDirections_()
 
  313    { 
return (dimWorld < 3) ? 4 : 13; }
 
  315    template<
class Geometry>
 
  316    bool neglectElement_(Geometry& hostGridElementGeometry,
 
  317                         std::size_t localVertexIdx0,
 
  318                         std::size_t localVertexIdx1,
 
  319                         double directionProbability)
 
  321        if (randomNumer_() < directionProbability)
 
  326        if (!considerLowDimElement_(hostGridElementGeometry.corner(localVertexIdx0), hostGridElementGeometry.corner(localVertexIdx1)))
 
  333    { 
return uniformdist_(generator_); }
 
  335    auto getHostGridInputData_()
 const 
  339            std::array<std::vector<int>, dimWorld> cells;
 
  340            std::array<std::vector<CoordScalar>, dimWorld> positions;
 
  341            std::array<std::vector<CoordScalar>, dimWorld> grading;
 
  342            GlobalPosition lowerLeft;
 
  343            GlobalPosition upperRight;
 
  347        auto& cells = inputData.cells;
 
  348        auto& positions = inputData.positions;
 
  349        auto& grading = inputData.grading;
 
  350        auto& lowerLeft = inputData.lowerLeft;
 
  351        auto& upperRight = inputData.upperRight;
 
  354        for (
int i = 0; i < dimWorld; ++i)
 
  355            positions[i] = getParamFromGroup<std::vector<CoordScalar>>(paramGroup_, 
"Grid.Positions" + std::to_string(i), std::vector<CoordScalar>{});
 
  356        if (std::none_of(positions.begin(), positions.end(), [&](
auto& v){ return v.empty(); }))
 
  358            for (
int i = 0; i < dimWorld; ++i)
 
  360                cells[i].resize(positions[i].size()-1, 1.0);
 
  361                grading[i].resize(positions[i].size()-1, 1.0);
 
  366            const auto lowerLeft = getParamFromGroup<GlobalPosition>(paramGroup_, 
"Grid.LowerLeft", GlobalPosition(0.0));
 
  367            const auto upperRight = getParamFromGroup<GlobalPosition>(paramGroup_, 
"Grid.UpperRight");
 
  368            const auto numPores = getParamFromGroup<std::vector<unsigned int>>(paramGroup_, 
"Grid.NumPores");
 
  369            if (numPores.size() != dimWorld)
 
  370                DUNE_THROW(ParameterException, 
"Grid.NumPores has to be a space-separated list of " << dimWorld << 
" integers!");
 
  372            for (
int i = 0; i < dimWorld; ++i)
 
  374                positions[i].push_back(lowerLeft[i]);
 
  375                positions[i].push_back(upperRight[i]);
 
  376                cells[i].push_back(numPores[i] - 1);
 
  377                grading[i].resize(positions[i].size()-1, 1.0);
 
  378                grading[i] = getParamFromGroup<std::vector<CoordScalar>>(paramGroup_, 
"Grid.Grading" + std::to_string(i), grading[i]);
 
  385            GlobalPosition result;
 
  386            for (
int i = 0; i < dimWorld; ++i)
 
  387                result[i] = positions[i][0];
 
  394            GlobalPosition result;
 
  395            for (
int i = 0; i < dimWorld; ++i)
 
  396                result[i] = positions[i].back();
 
  403    template<
class HostGr
idElement, 
class ReferenceElement, 
class Geometry>
 
  404    void treatEdges_(
const HostGridElement& element, 
const ReferenceElement& refElement, 
const Geometry& geometry)
 
  407        for (
unsigned int edgeIdx = 0; edgeIdx < 
element.subEntities(dimWorld-1); ++edgeIdx)
 
  409            const auto vIdxLocal0 = refElement.subEntity(edgeIdx, dimWorld-1, 0, dimWorld);
 
  410            const auto vIdxLocal1 = refElement.subEntity(edgeIdx, dimWorld-1, 1, dimWorld);
 
  411            const auto edgeIdxGlobal = hostGridView_->indexSet().subIndex(element, edgeIdx, dimWorld-1);
 
  413            if (edgeInserted_[edgeIdxGlobal])
 
  416                edgeInserted_[edgeIdxGlobal] = 
true;
 
  418            std::size_t directionNumber = 0;
 
  431                else if(edgeIdx == 4 || edgeIdx == 5 || edgeIdx == 8 || edgeIdx == 9) 
 
  433                else if(edgeIdx == 6 || edgeIdx == 7 || edgeIdx == 10 || edgeIdx == 11) 
 
  436            maybeInsertElementAndVertices_(element, geometry, vIdxLocal0, vIdxLocal1, directionProbability_[directionNumber]);
 
  440    template<
class HostGr
idElement, 
class ReferenceElement, 
class Geometry>
 
  441    void treatDiagonals_(
const HostGridElement& element, 
const ReferenceElement& refElement, 
const Geometry& geometry)
 
  443        if constexpr (dimWorld < 3)
 
  444            treatDiagonalConnections2D_(element, geometry);
 
  446            treatDiagonalConnections3D_(element, refElement, geometry);
 
  449    template<
class HostGr
idElement, 
class Geometry>
 
  450    void treatDiagonalConnections2D_(
const HostGridElement& element,
 
  451                                     const Geometry& geometry)
 
  453         static constexpr std::array<Diagonal, 2> diagonals{ Diagonal{std::make_pair(0, 3), 2},
 
  454                                                             Diagonal{std::make_pair(1, 2), 3} };
 
  456        treatIntersectingDiagonals_(element, geometry, diagonals);
 
  459    template<
class HostGr
idElement, 
class ReferenceElement, 
class Geometry>
 
  460    void treatDiagonalConnections3D_(
const HostGridElement& element,
 
  461                                     const ReferenceElement& refElement,
 
  462                                     const Geometry& geometry)
 
  465        for (
auto faceIdx = 0; faceIdx < 
element.subEntities(dimWorld-2); ++faceIdx)
 
  467            const auto faceIdxGlobal = hostGridView_->indexSet().subIndex(element, faceIdx, dimWorld-2);
 
  469            if (faceInserted_[faceIdxGlobal])
 
  472                faceInserted_[faceIdxGlobal] = 
true;
 
  475            std::array<unsigned int, 4> vIdxLocal;
 
  476            for (
int i = 0; i < 4; i++)
 
  477                vIdxLocal[i]  = refElement.subEntity(faceIdx, dimWorld-2, i, dimWorld);
 
  479            const auto directionNumbers = [&]()
 
  482                    return std::make_pair<unsigned int, unsigned int>(8,7);
 
  483                else if (faceIdx < 4)
 
  484                    return std::make_pair<unsigned int, unsigned int>(6,5);
 
  486                    return std::make_pair<unsigned int, unsigned int>(4,3);
 
  489            const std::array<Diagonal, 2> diagonals{ Diagonal{std::make_pair(vIdxLocal[1], vIdxLocal[2]), directionNumbers.first},
 
  490                                                     Diagonal{std::make_pair(vIdxLocal[0], vIdxLocal[3]), directionNumbers.second} };
 
  492            treatIntersectingDiagonals_(element, geometry, diagonals);
 
  495        static constexpr std::array<Diagonal, 4> diagonals{ Diagonal{std::make_pair(0, 7), 9},
 
  496                                                            Diagonal{std::make_pair(3, 4), 10},
 
  497                                                            Diagonal{std::make_pair(1, 6), 11},
 
  498                                                            Diagonal{std::make_pair(2, 5), 12}, };
 
  500        treatIntersectingDiagonals_(element, geometry, diagonals);
 
  503    template<
class HostGr
idElement, 
class Geometry, 
class Diagonals>
 
  504    void treatIntersectingDiagonals_(
const HostGridElement& element,
 
  505                                     const Geometry& geometry,
 
  506                                     const Diagonals& diagonals)
 
  508        static const bool allowIntersectingDiagonals = getParamFromGroup<bool>(paramGroup_, 
"Grid.AllowIntersectingDiagonals", 
true);
 
  510        auto treat = [&](
const Diagonal& diagonal)
 
  512            return maybeInsertElementAndVertices_(element, geometry,
 
  513                                                  diagonal.localVertexIndices.first, diagonal.localVertexIndices.second,
 
  514                                                  directionProbability_[diagonal.directionNumber]);
 
  517        if (allowIntersectingDiagonals)
 
  520            for (
const auto& diagonal : diagonals)
 
  525            auto order = createOrderedList_(diagonals.size());
 
  526            std::shuffle(order.begin(), order.end(), generator_);
 
  529                const auto& diagonal = diagonals[i];
 
  530                const bool inserted = treat(diagonal);
 
  537    std::vector<std::size_t> createOrderedList_(
const std::size_t size)
 const 
  539        std::vector<std::size_t> result(size);
 
  540        std::iota(result.begin(), result.end(), 0);
 
  544    auto getLambdaToRemoveThroatsOnBoundary_()
 const 
  546        auto negletElementsOnBoundary = [&](
const GlobalPosition& a, 
const GlobalPosition& b)
 
  548            const auto center = 0.5 * (a + b);
 
  549            const auto eps = (a-b).two_norm() * 1e-5;
 
  551            bool neglectElement = 
false;
 
  552            for (
auto i : removeThroatsOnBoundary_)
 
  556                    case 0: neglectElement = 
center[0] < hostGridLowerLeft_[0] + eps; 
break;
 
  557                    case 1: neglectElement = 
center[0] > hostGridUpperRight_[0] - eps; 
break;
 
  558                    case 2: 
if constexpr (dimWorld > 1)
 
  560                                neglectElement = 
center[1] < hostGridLowerLeft_[1] + eps;
 
  563                    case 3: 
if constexpr (dimWorld > 1)
 
  565                                neglectElement = 
center[1] > hostGridUpperRight_[1] - eps;
 
  568                    case 4: 
if constexpr (dimWorld > 2)
 
  570                                neglectElement = 
center[2] < hostGridLowerLeft_[2] + eps;
 
  573                    case 5: 
if constexpr (dimWorld > 2)
 
  575                                neglectElement = 
center[2] > hostGridUpperRight_[2] - eps;
 
  587        return negletElementsOnBoundary;
 
  590    std::string paramGroup_;
 
  591    GlobalPosition hostGridLowerLeft_;
 
  592    GlobalPosition hostGridUpperRight_;
 
  593    std::vector<std::size_t> removeThroatsOnBoundary_;
 
  594    std::unique_ptr<const HostGridView> hostGridView_;
 
  595    std::function<bool(
const GlobalPosition&, 
const GlobalPosition&)> considerLowDimElement_;
 
  597    std::vector<GlobalPosition> vertexSet_;
 
  598    std::vector<VertexIndexPair> elementSet_;
 
  600    std::vector<bool> vertexInserted_;
 
  601    std::vector<std::size_t> hostVertexIdxToVertexIdx_;
 
  602    std::vector<std::size_t> edgeInserted_;
 
  603    std::vector<std::size_t> faceInserted_;
 
  605    mutable std::mt19937 generator_;
 
  606    std::uniform_real_distribution<> uniformdist_{0, 1};
 
  607    std::array<double, numDirections_()> directionProbability_;
 
  609    std::shared_ptr<Grid> gridPtr_;
 
Some exceptions thrown in DuMux
Grid manager specialization for YaspGrid.
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition center.hh:24
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition parameters.hh:165
Definition discretization/porenetwork/fvelementgeometry.hh:24
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.