Added well controls in spu2p
This commit is contained in:
parent
738ec64ac8
commit
338f5907a6
@ -516,6 +516,8 @@ main(int argc, char** argv)
|
|||||||
Opm::WellReport wellreport;
|
Opm::WellReport wellreport;
|
||||||
std::vector<double> well_bhp;
|
std::vector<double> well_bhp;
|
||||||
std::vector<double> well_perfrates;
|
std::vector<double> well_perfrates;
|
||||||
|
std::vector<double> fractional_flows;
|
||||||
|
std::vector<double> well_resflows;
|
||||||
if (wells->c_wells()) {
|
if (wells->c_wells()) {
|
||||||
const int nw = wells->c_wells()->number_of_wells;
|
const int nw = wells->c_wells()->number_of_wells;
|
||||||
well_bhp.resize(nw, 0.0);
|
well_bhp.resize(nw, 0.0);
|
||||||
@ -539,38 +541,61 @@ main(int argc, char** argv)
|
|||||||
if (wells->c_wells()) {
|
if (wells->c_wells()) {
|
||||||
Opm::computeWDP(*wells->c_wells(), *grid->c_grid(), state.saturation(), props->density(), gravity[2], true, wdp);
|
Opm::computeWDP(*wells->c_wells(), *grid->c_grid(), state.saturation(), props->density(), gravity[2], true, wdp);
|
||||||
}
|
}
|
||||||
pressure_timer.start();
|
if (check_well_controls) {
|
||||||
if (rock_comp->isActive()) {
|
computeFractionalFlow(*props, allcells, state.saturation(), fractional_flows);
|
||||||
rc.resize(num_cells);
|
}
|
||||||
std::vector<double> initial_pressure = state.pressure();
|
bool well_control_passed = !check_well_controls;
|
||||||
std::vector<double> prev_pressure;
|
int well_control_iteration = 0;
|
||||||
for (int iter = 0; iter < nl_pressure_maxiter; ++iter) {
|
do {
|
||||||
prev_pressure = state.pressure();
|
pressure_timer.start();
|
||||||
for (int cell = 0; cell < num_cells; ++cell) {
|
if (rock_comp->isActive()) {
|
||||||
rc[cell] = rock_comp->rockComp(state.pressure()[cell]);
|
rc.resize(num_cells);
|
||||||
|
std::vector<double> initial_pressure = state.pressure();
|
||||||
|
std::vector<double> prev_pressure;
|
||||||
|
for (int iter = 0; iter < nl_pressure_maxiter; ++iter) {
|
||||||
|
prev_pressure = state.pressure();
|
||||||
|
for (int cell = 0; cell < num_cells; ++cell) {
|
||||||
|
rc[cell] = rock_comp->rockComp(state.pressure()[cell]);
|
||||||
|
}
|
||||||
|
state.pressure() = initial_pressure;
|
||||||
|
psolver.solve(totmob, omega, src, wdp, bcs.c_bcs(), porevol, rc, simtimer.currentStepLength(),
|
||||||
|
state.pressure(), state.faceflux(), well_bhp, well_perfrates);
|
||||||
|
double max_change = 0.0;
|
||||||
|
for (int cell = 0; cell < num_cells; ++cell) {
|
||||||
|
max_change = std::max(max_change, std::fabs(state.pressure()[cell] - prev_pressure[cell]));
|
||||||
|
}
|
||||||
|
std::cout << "Pressure iter " << iter << " max change = " << max_change << std::endl;
|
||||||
|
if (max_change < nl_pressure_tolerance) {
|
||||||
|
break;
|
||||||
|
}
|
||||||
}
|
}
|
||||||
state.pressure() = initial_pressure;
|
computePorevolume(*grid->c_grid(), *props, *rock_comp, state.pressure(), porevol);
|
||||||
psolver.solve(totmob, omega, src, wdp, bcs.c_bcs(), porevol, rc, simtimer.currentStepLength(),
|
} else {
|
||||||
state.pressure(), state.faceflux(), well_bhp, well_perfrates);
|
psolver.solve(totmob, omega, src, wdp, bcs.c_bcs(), state.pressure(), state.faceflux(),
|
||||||
double max_change = 0.0;
|
well_bhp, well_perfrates);
|
||||||
for (int cell = 0; cell < num_cells; ++cell) {
|
}
|
||||||
max_change = std::max(max_change, std::fabs(state.pressure()[cell] - prev_pressure[cell]));
|
pressure_timer.stop();
|
||||||
|
double pt = pressure_timer.secsSinceStart();
|
||||||
|
std::cout << "Pressure solver took: " << pt << " seconds." << std::endl;
|
||||||
|
ptime += pt;
|
||||||
|
|
||||||
|
if (check_well_controls) {
|
||||||
|
Opm::computePhaseFlowRatesPerWell(*wells->c_wells(),
|
||||||
|
fractional_flows,
|
||||||
|
well_perfrates,
|
||||||
|
well_resflows);
|
||||||
|
std::cout << "Checking well conditions." << std::endl;
|
||||||
|
// For testing we set surface := reservoir
|
||||||
|
well_control_passed = wells->conditionsMet(well_bhp, well_resflows, well_resflows);
|
||||||
|
++well_control_iteration;
|
||||||
|
if (!well_control_passed && well_control_iteration > max_well_control_iterations) {
|
||||||
|
THROW("Could not satisfy well conditions in " << max_well_control_iterations << " tries.");
|
||||||
}
|
}
|
||||||
std::cout << "Pressure iter " << iter << " max change = " << max_change << std::endl;
|
if (!well_control_passed) {
|
||||||
if (max_change < nl_pressure_tolerance) {
|
std::cout << "Well controls not passed, solving again." << std::endl;
|
||||||
break;
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
computePorevolume(*grid->c_grid(), *props, *rock_comp, state.pressure(), porevol);
|
} while (!well_control_passed);
|
||||||
} else {
|
|
||||||
psolver.solve(totmob, omega, src, wdp, bcs.c_bcs(), state.pressure(), state.faceflux(),
|
|
||||||
well_bhp, well_perfrates);
|
|
||||||
}
|
|
||||||
pressure_timer.stop();
|
|
||||||
double pt = pressure_timer.secsSinceStart();
|
|
||||||
std::cout << "Pressure solver took: " << pt << " seconds." << std::endl;
|
|
||||||
ptime += pt;
|
|
||||||
|
|
||||||
// Process transport sources (to include bdy terms and well flows).
|
// Process transport sources (to include bdy terms and well flows).
|
||||||
Opm::computeTransportSource(*grid->c_grid(), src, state.faceflux(), 1.0,
|
Opm::computeTransportSource(*grid->c_grid(), src, state.faceflux(), 1.0,
|
||||||
wells->c_wells(), well_perfrates, reorder_src);
|
wells->c_wells(), well_perfrates, reorder_src);
|
||||||
|
Loading…
Reference in New Issue
Block a user