mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
[bugfix] handle exception when linearizing the wells.
As multisegment wells may throw in applyUMFPack this is now needed and the exception needs to communicated to all processes. We do this in the linearize method of the well model. Before this change this is what could happen: - The process with the exception would have chopped the time step - The others would have successfully setup the systems and entered the linear solve This poduced a deadlock. One processes was waiting in OPM_END_PARALLEL_TRY during the setup of the shorter time step and in collective communication during the setup of the linear solver for the unchopped time step.
This commit is contained in:
parent
d07aed204f
commit
1520f9c85f
@ -335,14 +335,15 @@ namespace Opm {
|
||||
const int nc = UgGridHelpers::numCells(grid_);
|
||||
BVector x(nc);
|
||||
|
||||
// apply the Schur compliment of the well model to the reservoir linearized
|
||||
// equations
|
||||
wellModel().linearize(ebosSimulator().model().linearizer().jacobian(),
|
||||
ebosSimulator().model().linearizer().residual());
|
||||
|
||||
// Solve the linear system.
|
||||
linear_solve_setup_time_ = 0.0;
|
||||
try {
|
||||
// apply the Schur compliment of the well model to the reservoir linearized
|
||||
// equations
|
||||
// Note that linearize may throw for MSwells.
|
||||
wellModel().linearize(ebosSimulator().model().linearizer().jacobian(),
|
||||
ebosSimulator().model().linearizer().residual());
|
||||
|
||||
solveJacobianSystem(x);
|
||||
report.linear_solve_setup_time += linear_solve_setup_time_;
|
||||
report.linear_solve_time += perfTimer.stop();
|
||||
|
@ -152,13 +152,19 @@ namespace Opm {
|
||||
if (!localWellsActive())
|
||||
return;
|
||||
|
||||
if (!param_.matrix_add_well_contributions_) {
|
||||
// if the well contributions are not supposed to be included explicitly in
|
||||
// the matrix, we only apply the vector part of the Schur complement here.
|
||||
for (const auto& well: well_container_) {
|
||||
// r = r - duneC_^T * invDuneD_ * resWell_
|
||||
well->apply(res);
|
||||
if (!param_.matrix_add_well_contributions_)
|
||||
{
|
||||
OPM_BEGIN_PARALLEL_TRY_CATCH();
|
||||
{
|
||||
// if the well contributions are not supposed to be included explicitly in
|
||||
// the matrix, we only apply the vector part of the Schur complement here.
|
||||
for (const auto& well: well_container_) {
|
||||
// r = r - duneC_^T * invDuneD_ * resWell_
|
||||
well->apply(res);
|
||||
}
|
||||
}
|
||||
OPM_END_PARALLEL_TRY_CATCH("BlackoilWellModel::linearize failed: ",
|
||||
ebosSimulator_.gridView().comm());
|
||||
return;
|
||||
}
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user