inline void EulerUpstream<GI, RP, BC>::initObj(const GI& g, const RP& r, const BC& b) { residual_computer_.initObj(g, r, b); porevol_.resize(g.numberOfCells()); for (CIt c = g.cellbegin(); c != g.cellend(); ++c) { porevol_[c->index()] = c->volume()*r.porosity(c->index()); } }
inline void EulerUpstreamImplicit<GI, RP, BC>::initObj(const GI& g, const RP& r, const BC& b) { //residual_computer_.initObj(g, r, b); mygrid_.init(g.grid()); porevol_.resize(mygrid_.numCells()); for (int i = 0; i < mygrid_.numCells(); ++i){ porevol_[i]= mygrid_.cellVolume(i)*r.porosity(i); } // int numf=mygrid_.numFaces(); int num_cells = mygrid_.numCells(); int ngconn = mygrid_.c_grid()->cell_facepos[num_cells]; //std::vector<double> htrans_(ngconn); htrans_.resize(ngconn); const double* perm = &(r.permeability(0)(0,0)); tpfa_htrans_compute(mygrid_.c_grid(), perm, &htrans_[0]); // int count = 0; myrp_= r; typedef typename GI::CellIterator CIt; typedef typename CIt::FaceIterator FIt; std::vector<FIt> bid_to_face; int maxbid = 0; for (CIt c = g.cellbegin(); c != g.cellend(); ++c) { for (FIt f = c->facebegin(); f != c->faceend(); ++f) { int bid = f->boundaryId(); maxbid = std::max(maxbid, bid); } } bid_to_face.resize(maxbid + 1); std::vector<int> egf_cf(mygrid_.numFaces()); int cix=0; for (CIt c = g.cellbegin(); c != g.cellend(); ++c) { int loc_fix=0; for (FIt f = c->facebegin(); f != c->faceend(); ++f) { if (f->boundary() && b.satCond(*f).isPeriodic()) { bid_to_face[f->boundaryId()] = f; } int egf=f->index(); int cf=mygrid_.cellFace(cix,loc_fix); egf_cf[egf]=cf; loc_fix+=1; } cix+=1; } #ifndef NDEBUG const UnstructuredGrid& c_grid=*mygrid_.c_grid(); #endif int hf_ind=0; int bf_ind=0; periodic_cells_.resize(0); periodic_faces_.resize(0); periodic_hfaces_.resize(0); periodic_nbfaces_.resize(0); //cell1 = cell0; direclet_cells_.resize(0); direclet_sat_.resize(0); direclet_sat_.resize(0); direclet_hfaces_.resize(0); assert(periodic_cells_.size()==0); for (CIt c = g.cellbegin(); c != g.cellend(); ++c) { int cell0 = c->index(); for (FIt f = c->facebegin(); f != c->faceend(); ++f) { // Neighbour face, will be changed if on a periodic boundary. // Compute cell[1], cell_sat[1] FIt nbface = f; if (f->boundary()) { bf_ind+=1; if (b.satCond(*f).isPeriodic()) { nbface = bid_to_face[b.getPeriodicPartner(f->boundaryId())]; assert(nbface != f); int cell1 = nbface->cellIndex(); assert(cell0 != cell1); int f_ind=f->index(); int fn_ind=nbface->index(); // mapping face indices f_ind=egf_cf[f_ind]; fn_ind=egf_cf[fn_ind]; assert((c_grid.face_cells[2*f_ind]==-1) || (c_grid.face_cells[2*f_ind+1]==-1)); assert((c_grid.face_cells[2*fn_ind]==-1) || (c_grid.face_cells[2*fn_ind+1]==-1)); assert((c_grid.face_cells[2*f_ind]==cell0) || (c_grid.face_cells[2*f_ind+1]==cell0)); assert((c_grid.face_cells[2*fn_ind]==cell1) || (c_grid.face_cells[2*fn_ind+1]==cell1)); periodic_cells_.push_back(cell0); periodic_cells_.push_back(cell1); periodic_faces_.push_back(f_ind); periodic_hfaces_.push_back(hf_ind); periodic_nbfaces_.push_back(fn_ind); } else if (!( b.flowCond(*f).isNeumann() && b.flowCond(*f).outflux() == 0.0)) { //cell1 = cell0; direclet_cells_.push_back(cell0); direclet_sat_.push_back(b.satCond(*f).saturation()); direclet_sat_.push_back(1-b.satCond(*f).saturation());//only work for 2 phases direclet_hfaces_.push_back(hf_ind); } } hf_ind+=1; } } mygrid_.makeQPeriodic(periodic_hfaces_,periodic_cells_); // use fractional flow instead of saturation as src TwophaseFluid myfluid(myrp_); int num_b=direclet_cells_.size(); for(int i=0; i <num_b; ++i){ std::array<double,2> sat = {{direclet_sat_[2*i] ,direclet_sat_[2*i+1] }}; std::array<double,2> mob; std::array<double,2*2> dmob; myfluid.mobility(direclet_cells_[i], sat, mob, dmob); double fl = mob[0]/(mob[0]+mob[1]); direclet_sat_[2*i] = fl; direclet_sat_[2*i+1] = 1-fl; } }
void ImplicitCapillarity<GI, RP, BC, IP>::transportSolve(std::vector<double>& saturation, const double /*time*/, const typename GI::Vector& gravity, const PressureSolution& pressure_sol, const Opm::SparseVector<double>& injection_rates) const { // Start a timer. Opm::time::StopWatch clock; clock.start(); // Compute capillary mobilities. typedef typename RP::Mobility Mob; int num_cells = saturation.size(); std::vector<Mob> cap_mob(num_cells); for (int c = 0; c < num_cells; ++c) { Mob& m = cap_mob[c]; residual_.reservoirProperties().phaseMobility(0, c, saturation[c], m.mob); Mob mob2; residual_.reservoirProperties().phaseMobility(1, c, saturation[c], mob2.mob); Mob mob_tot; mob_tot.setToSum(m, mob2); Mob mob_totinv; mob_totinv.setToInverse(mob_tot); m *= mob_totinv; m *= mob2; ImplicitCapillarityDetails::thresholdMobility(m.mob, 1e-10); // @@TODO: User-set limit. // std::cout << m.mob(0,0) << '\n'; } ReservoirPropertyFixedMobility<Mob> capillary_mobilities(cap_mob); // Set up boundary conditions. BC cap_press_bcs(residual_.boundaryConditions()); for (int i = 0; i < cap_press_bcs.size(); ++i) { if (cap_press_bcs.flowCond(i).isPeriodic()) { cap_press_bcs.flowCond(i) = FlowBC(FlowBC::Periodic, 0.0); } } // Compute injection rates from residual. std::vector<double> injection_rates_residual(num_cells); residual_.computeResidual(saturation, gravity, pressure_sol, injection_rates, method_viscous_, method_gravity_, false, injection_rates_residual); for (int i = 0; i < num_cells; ++i) { injection_rates_residual[i] = -injection_rates_residual[i]; } // Compute capillary pressure. // Note that the saturation is just a dummy for this call, since the mobilities are fixed. psolver_.solve(capillary_mobilities, saturation, cap_press_bcs, injection_rates_residual, residual_tolerance_, linsolver_verbosity_, linsolver_type_); // Solve for constant to change capillary pressure solution by. std::vector<double> cap_press(num_cells); const PressureSolution& pcapsol = psolver_.getSolution(); for (CIt c = residual_.grid().cellbegin(); c != residual_.grid().cellend(); ++c) { cap_press[c->index()] = pcapsol.pressure(c); } MatchSaturatedVolumeFunctor<GI, RP> functor(residual_.grid(), residual_.reservoirProperties(), saturation, cap_press); double min_cap_press = *std::min_element(cap_press.begin(), cap_press.end()); double max_cap_press = *std::max_element(cap_press.begin(), cap_press.end()); double cap_press_range = max_cap_press - min_cap_press; double mod_low = 1e100; double mod_high = -1e100; Opm::bracketZero(functor, 0.0, cap_press_range, mod_low, mod_high); const int max_iter = 40; const double nonlinear_tolerance = 1e-12; int iterations_used = -1; typedef Opm::RegulaFalsi<Opm::ThrowOnError> RootFinder; double mod_correct = RootFinder::solve(functor, mod_low, mod_high, max_iter, nonlinear_tolerance, iterations_used); std::cout << "Moved capillary pressure solution by " << mod_correct << " after " << iterations_used << " iterations." << std::endl; // saturation = functor.lastSaturations(); const std::vector<double>& sat_new = functor.lastSaturations(); for (int i = 0; i < num_cells; ++i) { saturation[i] = (1.0 - update_relaxation_)*saturation[i] + update_relaxation_*sat_new[i]; } // Optionally check and/or clamp results. if (check_sat_ || clamp_sat_) { checkAndPossiblyClampSat(saturation); } // Stop timer and optionally print seconds taken. clock.stop(); #ifdef VERBOSE std::cout << "Seconds taken by transport solver: " << clock.secsSinceStart() << std::endl; #endif // VERBOSE }
void operator()(const CIt& c) const { // This is constant for the whole run. const double delta_rho = s.preservoir_properties_->densityDifference(); int cell[2]; double cell_sat[2]; cell[0] = c->index(); cell_sat[0] = saturation[cell[0]]; // Loop over all cell faces. for (FIt f = c->facebegin(); f != c->faceend(); ++f) { // Neighbour face, will be changed if on a periodic boundary. FIt nbface = f; double dS = 0.0; // Compute cell[1], cell_sat[1] if (f->boundary()) { if (s.pboundary_->satCond(*f).isPeriodic()) { nbface = s.bid_to_face_[s.pboundary_->getPeriodicPartner(f->boundaryId())]; assert(nbface != f); cell[1] = nbface->cellIndex(); assert(cell[0] != cell[1]); // Periodic faces will be visited twice, but only once // should they contribute. We make sure that we skip the // periodic faces half the time. if (cell[0] > cell[1]) { // We skip this face. continue; } cell_sat[1] = saturation[cell[1]]; } else { assert(s.pboundary_->satCond(*f).isDirichlet()); cell[1] = cell[0]; cell_sat[1] = s.pboundary_->satCond(*f).saturation(); } } else { cell[1] = f->neighbourCellIndex(); assert(cell[0] != cell[1]); if (cell[0] > cell[1]) { // We skip this face. continue; } cell_sat[1] = saturation[cell[1]]; } // Get some local properties. const double loc_area = f->area(); const double loc_flux = pressure_sol.outflux(f); const Vector loc_normal = f->normal(); // We will now try to establish the upstream directions for each // phase. They may be the same, or different (due to gravity). // Recall the equation for v_w (water phase velocity): // v_w = lambda_w * (lambda_o + lambda_w)^{-1} // * (v + lambda_o * K * grad p_{cow} + lambda_o * K * (rho_w - rho_o) * g) // ^ ^^^^^^^^^^^^^^^^^^^^^^^^^^^ ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ // viscous term capillary term gravity term // // For the purpose of upstream weighting, we only consider the viscous and gravity terms. // The question is, in which direction does v_w and v_o point? That is, what is the sign // of v_w*loc_normal and v_o*loc_normal? // // For the case when the mobilities are scalar, the following analysis applies: // The viscous contribution to v_w is loc_area*loc_normal*f_w*v == f_w*loc_flux. // Then the phase fluxes become // flux_w = f_w*(loc_flux + loc_area*loc_normal*lambda_o*K*(rho_w - rho_o)*g) // ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ // := lambda_o*G (only scalar case) // flux_o = f_o*(loc_flux - lambda_w*G) // In the above, we must decide where to evaluate K, and for this purpose (deciding // upstream directions) we use a K averaged between the two cells. // Since all mobilities and fractional flow functions are positive, the sign // of one of these cases is trivial. If G >= 0, flux_w is in the same direction as // loc_flux, if G <= 0, flux_o is in the same direction as loc_flux. // The phase k for which flux_k and loc_flux are of same sign, is called the trivial // phase in the code below. // // Assuming for the moment that G >=0, we know the direction of the water flux // (same as loc_flux) and evaluate lambda_w in the upstream cell. Then we may use // that lambda_w to evaluate flux_o using the above formula. Knowing flux_o, we know // the direction of the oil flux, and can evaluate lambda_o in the corresponding // upstream cell. Finally, we can use the equation for flux_w to compute that flux. // The opposite case is similar. // // What about tensorial mobilities? In the following code, we make the assumption // that the directions of vectors are not so changed by the multiplication with // mobility tensors that upstream directions change. In other words, we let all // the upstream logic stand as it is. This assumption may need to be revisited. // A worse problem is that // 1) we do not have v, just loc_area*loc_normal*v, // 2) we cannot define G, since the lambdas do not commute with the dot product. typedef typename UpstreamSolver::RP::Mobility Mob; using Opm::utils::arithmeticAverage; // Doing arithmetic averages. Should we consider harmonic or geometric instead? const MutablePermTensor aver_perm = arithAver(s.preservoir_properties_->permeability(cell[0]), s.preservoir_properties_->permeability(cell[1])); // Computing the raw gravity influence vector = (rho_w - rho_o)Kg Vector grav_influence = prod(aver_perm, gravity); grav_influence *= delta_rho; // Computing G. Note that we do not multiply with the mobility, // so this G is wrong in case of anisotropic relperm. const double G = s.method_gravity_ ? loc_area*inner(loc_normal, grav_influence) : 0.0; const int triv_phase = G >= 0.0 ? 0 : 1; const int ups_cell = loc_flux >= 0.0 ? 0 : 1; // Compute mobility of the trivial phase. Mob m_ups[2]; s.preservoir_properties_->phaseMobility(triv_phase, cell[ups_cell], cell_sat[ups_cell], m_ups[triv_phase].mob); // Compute gravity flow of the nontrivial phase. double sign_G[2] = { -1.0, 1.0 }; double grav_flux_nontriv = sign_G[triv_phase]*loc_area *inner(loc_normal, m_ups[triv_phase].multiply(grav_influence)); // Find flow direction of nontrivial phase. const int ups_cell_nontriv = (loc_flux + grav_flux_nontriv >= 0.0) ? 0 : 1; const int nontriv_phase = (triv_phase + 1) % 2; s.preservoir_properties_->phaseMobility(nontriv_phase, cell[ups_cell_nontriv], cell_sat[ups_cell_nontriv], m_ups[nontriv_phase].mob); // Now we have the upstream phase mobilities in m_ups[]. Mob m_tot; m_tot.setToSum(m_ups[0], m_ups[1]); Mob m_totinv; m_totinv.setToInverse(m_tot); const double aver_sat = Opm::utils::arithmeticAverage<double, double>(cell_sat[0], cell_sat[1]); Mob m1c0, m1c1, m2c0, m2c1; s.preservoir_properties_->phaseMobility(0, cell[0], aver_sat, m1c0.mob); s.preservoir_properties_->phaseMobility(0, cell[1], aver_sat, m1c1.mob); s.preservoir_properties_->phaseMobility(1, cell[0], aver_sat, m2c0.mob); s.preservoir_properties_->phaseMobility(1, cell[1], aver_sat, m2c1.mob); Mob m_aver[2]; m_aver[0].setToAverage(m1c0, m1c1); m_aver[1].setToAverage(m2c0, m2c1); Mob m_aver_tot; m_aver_tot.setToSum(m_aver[0], m_aver[1]); Mob m_aver_totinv; m_aver_totinv.setToInverse(m_aver_tot); // Viscous (pressure driven) term. if (s.method_viscous_) { // v is not correct for anisotropic relperm. Vector v(loc_normal); v *= loc_flux; const double visc_change = inner(loc_normal, m_ups[0].multiply(m_totinv.multiply(v))); // const double visc_change = (m_ups[0].mob/(m_ups[1].mob + m_ups[0].mob))*loc_flux; // std::cout << "New: " << visc_change_2 << " old: " << visc_change << '\n'; dS += visc_change; } // Gravity term. if (s.method_gravity_) { if (cell[0] != cell[1]) { // We only add gravity flux on internal or periodic faces. const double grav_change = loc_area *inner(loc_normal, m_ups[0].multiply(m_totinv.multiply(m_ups[1].multiply(grav_influence)))); // const double grav_change = (lambda_one*lambda_two/(lambda_two+lambda_one))*G; // const double grav_change = (lambda_one*lambda_two/(lambda_two+lambda_one))*loc_gravity_flux; dS += grav_change; } } // Capillary term. if (s.method_capillary_) { // J(s_w) = \frac{p_c(s_w)\sqrt{k/\phi}}{\sigma \cos\theta} // p_c = \frac{J \sigma \cos\theta}{\sqrt{k/\phi}} Vector cap_influence = prod(aver_perm, s.estimateCapPressureGradient(f, nbface, saturation)); const double cap_change = loc_area *inner(loc_normal, m_aver[0].multiply(m_aver_totinv.multiply(m_aver[1].multiply(cap_influence)))); // const double cap_vel = inner(loc_normal, prod(aver_perm, estimateCapPressureGradient(f, nbface, saturation))); // const double loc_cap_flux = cap_vel*loc_area; // // const double cap_change = loc_cap_flux*(m_aver[1].mob*m_aver[0].mob // // /(m_aver[0].mob + m_aver[1].mob)); // const double cap_change = loc_cap_flux*(aver_lambda_two*aver_lambda_one // /(aver_lambda_one + aver_lambda_two)); dS += cap_change; } // Modify saturation. if (cell[0] != cell[1]){ residual[cell[0]] -= dS; residual[cell[1]] += dS; } else { assert(cell[0] == cell[1]); residual[cell[0]] -= dS; } } // Source term. double rate = s.pinjection_rates_->element(cell[0]); if (rate < 0.0) { // For anisotropic relperm, fractionalFlow does not really make sense // as a scalar rate *= s.preservoir_properties_->fractionalFlow(cell[0], cell_sat[0]); } residual[cell[0]] += rate; }