void init(const Grid& grid, const double* perm, const double* porosity, const typename Grid::Vector& gravity)
    {
        // Build C grid structure.
        grid_.init(grid);

        // Initialize data.
        int num_phases = 3;
        well_t* w = 0;
        if (wells_.number_of_wells != 0) {
            w = &wells_;
        }
        data_ = cfs_tpfa_construct(grid_.c_grid(), w, num_phases);
        if (!data_) {
            throw std::runtime_error("Failed to initialize cfs_tpfa solver.");
        }

        // Compute half-transmissibilities
        int num_cells = grid.numCells();
        int ngconn  = grid_.c_grid()->cell_facepos[num_cells];
        ncf_.resize(num_cells);
        for (int cell = 0; cell < num_cells; ++cell) {
            int num_local_faces = grid.numCellFaces(cell);
            ncf_[cell] = num_local_faces;
        }
        htrans_.resize(ngconn);
        tpfa_htrans_compute(grid_.c_grid(), perm, &htrans_[0]);

        // Compute transmissibilities.
        trans_.resize(grid_.numFaces());
        tpfa_trans_compute(grid_.c_grid(), &htrans_[0], &trans_[0]);

        // Compute pore volumes.
        porevol_.resize(num_cells);
        for (int i = 0; i < num_cells; ++i) {
            porevol_[i] = porosity[i]*grid.cellVolume(i);
        }

        // Set gravity.
        if (Grid::dimension != 3) {
            throw std::logic_error("Only 3 dimensions supported currently.");
        }
        std::copy(gravity.begin(), gravity.end(), gravity_);

        state_ = Initialized;
    }
Exemplo n.º 2
0
        void computeUpwindProperties(const Grid& grid,
                                     const BlackoilFluid& fluid,
                                     const typename Grid::Vector gravity,
                                     const std::vector<PhaseVec>& cell_pressure,
                                     const std::vector<PhaseVec>& face_pressure,
                                     const std::vector<CompVec>& cell_z,
                                     const CompVec& bdy_z)
        {
            int num_faces = face_pressure.size();
            ASSERT(num_faces == grid.numFaces());
            bool nonzero_gravity = gravity.two_norm() > 0.0;
            face_data.state_matrix.resize(num_faces);
            face_data.mobility.resize(num_faces);
            face_data.mobility_deriv.resize(num_faces);
            face_data.gravity_potential.resize(num_faces);
            face_data.surface_volume_density.resize(num_faces);
#pragma omp parallel for
            for (int face = 0; face < num_faces; ++face) {
            // Obtain properties from both sides of the face.
                typedef typename Grid::Vector Vec;
                Vec fc = grid.faceCentroid(face);
                int c[2] = { grid.faceCell(face, 0), grid.faceCell(face, 1) };

                // Get pressures and compute gravity contributions,
                // to decide upwind directions.
                PhaseVec phase_p[2];
                PhaseVec gravcontrib[2];
                for (int j = 0; j < 2; ++j) {
                    if (c[j] >= 0) {
                        // Pressures
                        phase_p[j] = cell_pressure[c[j]];
                        // Gravity contribution.
                        if (nonzero_gravity) {
                            Vec cdiff = fc;
                            cdiff -= grid.cellCentroid(c[j]);
                            gravcontrib[j] = fluid.phaseDensities(&cell_data.state_matrix[c[j]][0][0]);
                            gravcontrib[j] *= (cdiff*gravity);
                        } else {
                            gravcontrib[j] = 0.0;
                        }
                    } else {
                        // Pressures
                        phase_p[j] = face_pressure[face];
                        // Gravity contribution.
                        gravcontrib[j] = 0.0;
                    }
                }

                // Gravity contribution:
                //    gravcapf = rho_1*g*(z_12 - z_1) - rho_2*g*(z_12 - z_2)
                // where _1 and _2 refers to two neigbour cells, z is the
                // z coordinate of the centroid, and z_12 is the face centroid.
                // Also compute the potentials.
                PhaseVec pot[2];
                for (int phase = 0; phase < numPhases; ++phase) {
                    face_data.gravity_potential[face][phase] = gravcontrib[0][phase] - gravcontrib[1][phase];
                    pot[0][phase] = phase_p[0][phase] + face_data.gravity_potential[face][phase];
                    pot[1][phase] = phase_p[1][phase];
                }

                // Now we can easily find the upwind direction for every phase,
                // we can also tell which boundary faces are inflow bdys.

                // Compute face_z, which is averaged from the cells, unless on outflow or noflow bdy.
                // Get mobilities and derivatives.
                CompVec face_z(0.0);
                double face_z_factor = 0.5;
                PhaseVec phase_mob[2];
                PhaseJacobian phasemob_deriv[2];
                for (int j = 0; j < 2; ++j) {
                    if (c[j] >= 0) {
                        face_z += cell_z[c[j]];
                        phase_mob[j] = cell_data.mobility[c[j]];
                        phasemob_deriv[j] = cell_data.mobility_deriv[c[j]];
                    } else if (pot[j][Liquid] > pot[(j+1)%2][Liquid]) {
                        // Inflow boundary.
                        face_z += bdy_z;
                        FluidStateBlackoil bdy_state = fluid.computeState(face_pressure[face], bdy_z);
                        phase_mob[j] = bdy_state.mobility_;
                        phasemob_deriv[j] = bdy_state.dmobility_;
                    } else {
                        // For outflow or noflow boundaries, only cell z is used.
                        face_z_factor = 1.0;
                        // Also, make sure the boundary data are not used for mobilities.
                        pot[j] = -1e100;
                    }
                }
                face_z *= face_z_factor;
                face_data.surface_volume_density[face] = face_z;

                // Computing upwind mobilities and derivatives
                for (int phase = 0; phase < numPhases; ++phase) {
                    if (pot[0][phase] == pot[1][phase]) {
                        // Average.
                        double aver = 0.5*(phase_mob[0][phase] + phase_mob[1][phase]);
                        face_data.mobility[face][phase] = aver;
                        for (int p2 = 0; p2 < numPhases; ++p2) {
                            face_data.mobility_deriv[face][phase][p2] = phasemob_deriv[0][phase][p2]
                                + phasemob_deriv[1][phase][p2];
                        }
                    } else {
                        // Upwind.
                        int upwind = pot[0][phase] > pot[1][phase] ? 0 : 1;
                        face_data.mobility[face][phase] = phase_mob[upwind][phase];
                        for (int p2 = 0; p2 < numPhases; ++p2) {
                            face_data.mobility_deriv[face][phase][p2] = phasemob_deriv[upwind][phase][p2];
                        }
                    }
                }
            }
        }
        virtual void init(const Opm::parameter::ParameterGroup& param,
                          const Grid& grid,
                          const Fluid& fluid,
                          typename Grid::Vector gravity,
                          State& simstate)
        {
            typedef typename Fluid::CompVec CompVec;
            typedef typename Fluid::PhaseVec PhaseVec;

            if (param.getDefault("heterogenous_initial_mix", false)) {
                CompVec init_oil(0.0);
                init_oil[Fluid::Oil] = 1.0;
                CompVec init_water(0.0);
                init_water[Fluid::Water] = 1.0;
                simstate.cell_z_.resize(grid.numCells());
                std::fill(simstate.cell_z_.begin(), simstate.cell_z_.begin() + simstate.cell_z_.size()/2, init_oil);
                std::fill(simstate.cell_z_.begin() + simstate.cell_z_.size()/2, simstate.cell_z_.end(), init_water);
                OPM_MESSAGE("******* Assuming zero capillary pressures *******");
                PhaseVec init_p(100.0*Opm::unit::barsa);
                simstate.cell_pressure_.resize(grid.numCells(), init_p);
                //         if (gravity.two_norm() != 0.0) {
                //             double ref_gravpot = grid.cellCentroid(0)*gravity;
                //             double rho = init_z*fluid_.surfaceDensities();  // Assuming incompressible, and constant initial z.
                //             for (int cell = 1; cell < grid.numCells(); ++cell) {
                //                 double press = rho*(grid.cellCentroid(cell)*gravity - ref_gravpot) + simstate.cell_pressure_[0][0];
                //                 simstate.cell_pressure_[cell] = PhaseVec(press);
                //             }
                //         }
            } else if (param.getDefault("unstable_initial_mix", false)) {
                CompVec init_oil(0.0);
                init_oil[Fluid::Oil] = 1.0;
                init_oil[Fluid::Gas] = 0.0;
                CompVec init_water(0.0);
                init_water[Fluid::Water] = 1.0;
                CompVec init_gas(0.0);
                init_gas[Fluid::Gas] = 150.0;
                simstate.cell_z_.resize(grid.numCells());
                std::fill(simstate.cell_z_.begin(),
                          simstate.cell_z_.begin() + simstate.cell_z_.size()/3,
                          init_water);
                std::fill(simstate.cell_z_.begin() + simstate.cell_z_.size()/3,
                          simstate.cell_z_.begin() + 2*(simstate.cell_z_.size()/3),
                          init_oil);
                std::fill(simstate.cell_z_.begin() + 2*(simstate.cell_z_.size()/3),
                          simstate.cell_z_.end(),
                          init_gas);
                OPM_MESSAGE("******* Assuming zero capillary pressures *******");
                PhaseVec init_p(100.0*Opm::unit::barsa);
                simstate.cell_pressure_.resize(grid.numCells(), init_p);

                if (gravity.two_norm() != 0.0) {
            
                    typename Fluid::FluidState state = fluid.computeState(simstate.cell_pressure_[0], simstate.cell_z_[0]);
                    simstate.cell_z_[0] *= 1.0/state.total_phase_volume_density_;
                    for (int cell = 1; cell < grid.numCells(); ++cell) {
                        double fluid_vol_dens;
                        int cnt =0;    
                        do {
                            double rho = 0.5*((simstate.cell_z_[cell]+simstate.cell_z_[cell-1])*fluid.surfaceDensities());
                            double press = rho*((grid.cellCentroid(cell) - grid.cellCentroid(cell-1))*gravity) + simstate.cell_pressure_[cell-1][0];
                            simstate.cell_pressure_[cell] = PhaseVec(press);
                            state = fluid.computeState(simstate.cell_pressure_[cell], simstate.cell_z_[cell]);
                            fluid_vol_dens = state.total_phase_volume_density_;
                            simstate.cell_z_[cell] *= 1.0/fluid_vol_dens;
                            ++cnt;
                        } while (std::fabs((fluid_vol_dens-1.0)) > 1.0e-8 && cnt < 10);
                
                    }  
                } else {
                    std::cout << "---- Exit - BlackoilSimulator.hpp: No gravity, no fun ... ----" << std::endl;
                    exit(-1);
                } 
            } else if (param.getDefault("CO2-injection", false)) {
                CompVec init_water(0.0);
                // Initially water filled (use Oil-component for water in order
                // to utilise blackoil mechanisms for brine-co2 interaction)          
                init_water[Fluid::Oil] = 1.0;  
                simstate.cell_z_.resize(grid.numCells());
                std::fill(simstate.cell_z_.begin(),simstate.cell_z_.end(),init_water);

                double datum_pressure_barsa = param.getDefault<double>("datum_pressure", 200.0);
                double datum_pressure = Opm::unit::convert::from(datum_pressure_barsa, Opm::unit::barsa);
                PhaseVec init_p(datum_pressure);
                simstate.cell_pressure_.resize(grid.numCells(), init_p);

                // Simple initial condition based on "incompressibility"-assumption
                double zMin = grid.cellCentroid(0)[2];
                for (int cell = 1; cell < grid.numCells(); ++cell) {
                    if (grid.cellCentroid(cell)[2] < zMin)
                        zMin = grid.cellCentroid(cell)[2];
                }

                typename Fluid::FluidState state = fluid.computeState(init_p, init_water);
		simstate.cell_z_[0] *= 1.0/state.total_phase_volume_density_;
                double density = (init_water*fluid.surfaceDensities())/state.total_phase_volume_density_;

                for (int cell = 0; cell < grid.numCells(); ++cell) {
                    double pressure(datum_pressure + (grid.cellCentroid(cell)[2] - zMin)*gravity[2]*density);
                    simstate.cell_pressure_[cell] = PhaseVec(pressure);
                    state = fluid.computeState(simstate.cell_pressure_[cell], simstate.cell_z_[cell]);
                    simstate.cell_z_[cell] *= 1.0/state.total_phase_volume_density_;
                }       
            } else {
                CompVec init_z(0.0);
                double initial_mixture_gas = param.getDefault("initial_mixture_gas", 0.0);
                double initial_mixture_oil = param.getDefault("initial_mixture_oil", 1.0);
                double initial_mixture_water = param.getDefault("initial_mixture_water", 0.0);
                init_z[Fluid::Water] = initial_mixture_water;
                init_z[Fluid::Gas] = initial_mixture_gas;
                init_z[Fluid::Oil] = initial_mixture_oil;

                simstate.cell_z_.resize(grid.numCells(), init_z);
                OPM_MESSAGE("******* Assuming zero capillary pressures *******");
                PhaseVec init_p(param.getDefault("initial_pressure", 100.0*Opm::unit::barsa));
                simstate.cell_pressure_.resize(grid.numCells(), init_p);
                if (gravity.two_norm() != 0.0) {
                    double ref_gravpot = grid.cellCentroid(0)*gravity;
                    double rho = init_z*fluid.surfaceDensities();  // Assuming incompressible, and constant initial z.
                    for (int cell = 1; cell < grid.numCells(); ++cell) {
                        double press = rho*(grid.cellCentroid(cell)*gravity - ref_gravpot) + simstate.cell_pressure_[0][0];
                        simstate.cell_pressure_[cell] = PhaseVec(press);
                    }
                }
            }
        }