24#ifndef OPM_COUNTGLOBALCELLS_HEADER_INCLUDED
25#define OPM_COUNTGLOBALCELLS_HEADER_INCLUDED
27#include <opm/core/props/BlackoilPhases.hpp>
29#include <dune/grid/common/gridview.hh>
38 std::vector<int> buildAllCells(
const int nc);
44 activePhases(
const PU& pu)
46 const int maxnp = BlackoilPhases::MaxNumPhases;
47 std::vector<bool> active(maxnp,
false);
49 for (
int p = 0; p < pu.MaxNumPhases; ++p) {
50 active[ p ] = pu.phase_used[ p ] != 0;
60 active2Canonical(
const PU& pu)
62 const int maxnp = BlackoilPhases::MaxNumPhases;
63 std::vector<int> act2can(maxnp, -1);
65 for (
int phase = 0; phase < maxnp; ++phase) {
66 if (pu.phase_used[ phase ]) {
67 act2can[ pu.phase_pos[ phase ] ] = phase;
76 double getGravity(
const double* g,
const int dim);
84 std::size_t countLocalInteriorCells(
const Grid& grid)
86 if ( grid.comm().size() == 1)
90 std::size_t count = 0;
91 const auto& gridView = grid.leafGridView();
92 for(
auto cell = gridView.template begin<0, Dune::Interior_Partition>(),
93 endCell = gridView.template end<0, Dune::Interior_Partition>();
94 cell != endCell; ++cell)
109 std::size_t countGlobalCells(
const Grid& grid)
111 if ( grid.comm().size() == 1)
115 std::size_t count = countLocalInteriorCells(grid);
116 return grid.comm().sum(count);
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition: BlackoilPhases.hpp:27