mirror of
https://github.com/OPM/opm-simulators.git
synced 2024-12-18 21:43:27 -06:00
Sending *dynamic* pressures to assembly function...
This commit is contained in:
parent
9aa0c89603
commit
5abf755e2d
@ -108,7 +108,7 @@ namespace Opm
|
||||
computePerIterationDynamicData(dt, state, well_state);
|
||||
|
||||
// Assemble J and F.
|
||||
assemble();
|
||||
assemble(dt, state, well_state);
|
||||
|
||||
bool residual_ok = false; // Replace with tolerance check.
|
||||
while (!residual_ok) {
|
||||
@ -125,7 +125,7 @@ namespace Opm
|
||||
computePerIterationDynamicData(dt, state, well_state);
|
||||
|
||||
// Assemble J and F.
|
||||
assemble();
|
||||
assemble(dt, state, well_state);
|
||||
|
||||
// Check for convergence.
|
||||
// Include both tolerance check for residual
|
||||
@ -182,14 +182,10 @@ namespace Opm
|
||||
|
||||
|
||||
/// Compute per-solve dynamic properties.
|
||||
void CompressibleTpfa::computePerSolveDynamicData(const double dt,
|
||||
void CompressibleTpfa::computePerSolveDynamicData(const double /*dt*/,
|
||||
const BlackoilState& state,
|
||||
const WellState& well_state)
|
||||
const WellState& /*well_state*/)
|
||||
{
|
||||
dt_ = dt;
|
||||
initial_press_ = state.pressure();
|
||||
initial_bhp_ = well_state.bhp();
|
||||
cell_z_ = &state.surfacevol()[0];
|
||||
computeWellPotentials(state);
|
||||
}
|
||||
|
||||
@ -407,10 +403,13 @@ namespace Opm
|
||||
|
||||
|
||||
/// Compute the residual and Jacobian.
|
||||
void CompressibleTpfa::assemble()
|
||||
void CompressibleTpfa::assemble(const double dt,
|
||||
const BlackoilState& state,
|
||||
const WellState& well_state)
|
||||
{
|
||||
const double* cell_press = &initial_press_[0];
|
||||
const double* well_bhp = initial_bhp_.empty() ? NULL : &initial_bhp_[0];
|
||||
const double* cell_press = &state.pressure()[0];
|
||||
const double* well_bhp = well_state.bhp().empty() ? NULL : &well_state.bhp()[0];
|
||||
const double* z = &state.surfacevol()[0];
|
||||
UnstructuredGrid* gg = const_cast<UnstructuredGrid*>(&grid_);
|
||||
CompletionData completion_data;
|
||||
completion_data.gpot = &wellperf_gpot_[0];
|
||||
@ -428,7 +427,7 @@ namespace Opm
|
||||
cq.Af = &face_A_[0];
|
||||
cq.phasemobf = &face_phasemob_[0];
|
||||
cq.voldiscr = &cell_voldisc_[0];
|
||||
cfs_tpfa_res_assemble(gg, dt_, &forces, cell_z_, &cq, &trans_[0],
|
||||
cfs_tpfa_res_assemble(gg, dt, &forces, z, &cq, &trans_[0],
|
||||
&face_gravcap_[0], cell_press, well_bhp,
|
||||
&porevol_[0], h_);
|
||||
}
|
||||
|
@ -86,7 +86,9 @@ namespace Opm
|
||||
void computeWellDynamicData(const double dt,
|
||||
const BlackoilState& state,
|
||||
const WellState& well_state);
|
||||
void assemble();
|
||||
void assemble(const double dt,
|
||||
const BlackoilState& state,
|
||||
const WellState& well_state);
|
||||
void solveIncrement();
|
||||
|
||||
void computeResults(std::vector<double>& pressure,
|
||||
@ -109,10 +111,6 @@ namespace Opm
|
||||
struct cfs_tpfa_res_data* h_;
|
||||
|
||||
// ------ Data that will be modified for every solve. ------
|
||||
double dt_;
|
||||
std::vector<double> initial_press_;
|
||||
std::vector<double> initial_bhp_;
|
||||
const double* cell_z_;
|
||||
std::vector<double> wellperf_gpot_;
|
||||
|
||||
// ------ Data that will be modified for every solver iteration. ------
|
||||
|
Loading…
Reference in New Issue
Block a user