#1112 Calculation of accumulated well flow is now implemented. Tracers are handled, as well as MSW's.

This commit is contained in:
Jacob Støren 2017-01-31 15:28:15 +01:00
parent b75fdb4a6d
commit d3f1c9d6fa
7 changed files with 450 additions and 191 deletions

View File

@ -39,136 +39,15 @@
#include "RiuMainPlotWindow.h"
#include "RiuWellAllocationPlot.h"
#include "RigAccWellFlowCalculator.h"
CAF_PDM_SOURCE_INIT(RimWellAllocationPlot, "WellAllocationPlot");
//==================================================================================================
//--------------------------------------------------------------------------------------------------
///
///
//==================================================================================================
class RigAccWellFlowCalculator
{
public:
RigAccWellFlowCalculator(const std::vector< std::vector <cvf::Vec3d> >& pipeBranchesCLCoords,
const std::vector< std::vector <RigWellResultPoint> >& pipeBranchesCellIds):
m_pipeBranchesCLCoords(pipeBranchesCLCoords), m_pipeBranchesCellIds(pipeBranchesCellIds)
{
m_accConnectionFlowPrBranch.resize(m_pipeBranchesCellIds.size());
calculateAccumulatedFlowPrConnection(0,1);
}
// Returned pair is connection number from top of well with accumulated flow
const std::vector<std::pair <size_t, double> >& accumulatedFlowPrConnection(size_t branchIdx)
{
return m_accConnectionFlowPrBranch[branchIdx];
}
private:
const std::vector<std::pair <size_t, double> >& calculateAccumulatedFlowPrConnection(size_t branchIdx, size_t startConnectionNumberFromTop)
{
std::vector<std::pair <size_t, double> >& accConnectionFlow = m_accConnectionFlowPrBranch[branchIdx];
const std::vector<RigWellResultPoint>& branchCells = m_pipeBranchesCellIds[branchIdx];
int clSegIdx = static_cast<int>( branchCells.size()) - 1;
double accFlow = 0.0;
size_t prevGridIdx = -1;
size_t prevGridCellIdx = -1;
std::vector<size_t> resPointToConnectionIndexFromBottom;
resPointToConnectionIndexFromBottom.resize(branchCells.size(), -1);
size_t connIdxFromBottom = 0;
while (clSegIdx >= 0)
{
resPointToConnectionIndexFromBottom[clSegIdx] = connIdxFromBottom;
if ( branchCells[clSegIdx].m_gridIndex != prevGridIdx
&& branchCells[clSegIdx].m_gridCellIndex != prevGridIdx )
{
++connIdxFromBottom;
}
prevGridIdx = branchCells[clSegIdx].m_gridIndex ;
prevGridCellIdx = branchCells[clSegIdx].m_gridCellIndex;
--clSegIdx;
}
size_t prevConnIndx = -1;
clSegIdx = static_cast<int>( branchCells.size()) - 1;
while (clSegIdx >= 0)
{
// Skip point if referring to the same cell as in the previous point did
{
if (resPointToConnectionIndexFromBottom[clSegIdx] == prevConnIndx)
{
--clSegIdx;
continue;
}
prevConnIndx = resPointToConnectionIndexFromBottom[clSegIdx];
}
size_t connNumFromTop = connecionIdxFromTop(resPointToConnectionIndexFromBottom, clSegIdx) + startConnectionNumberFromTop;
accFlow += branchCells[clSegIdx].m_flowRate;
std::vector<size_t> downstreamBranches = findDownstreamBranchIdxs(branchCells[clSegIdx]);
for (size_t dsBidx : downstreamBranches )
{
if ( dsBidx != branchIdx && m_accConnectionFlowPrBranch[dsBidx].size() == 0) // Not this branch or already calculated
{
const std::pair <size_t, double>& brancFlowPair = calculateAccumulatedFlowPrConnection(dsBidx, connNumFromTop).back();
accFlow += brancFlowPair.second;
}
}
accConnectionFlow.push_back(std::make_pair(connNumFromTop, accFlow));
--clSegIdx;
}
return m_accConnectionFlowPrBranch[branchIdx];
}
static size_t connecionIdxFromTop( std::vector<size_t> resPointToConnectionIndexFromBottom, size_t clSegIdx)
{
return resPointToConnectionIndexFromBottom.front() - resPointToConnectionIndexFromBottom[clSegIdx];
}
std::vector<size_t> findDownstreamBranchIdxs(const RigWellResultPoint& connectionPoint)
{
std::vector<size_t> downStreamBranchIdxs;
for (size_t bIdx = 0; bIdx < m_pipeBranchesCellIds.size(); ++bIdx)
{
if ( m_pipeBranchesCellIds[bIdx][0].m_gridIndex == connectionPoint.m_gridIndex
&& m_pipeBranchesCellIds[bIdx][0].m_gridCellIndex == connectionPoint.m_gridCellIndex)
{
downStreamBranchIdxs.push_back(bIdx);
}
}
return downStreamBranchIdxs;
}
const std::vector< std::vector <cvf::Vec3d> >& m_pipeBranchesCLCoords;
const std::vector< std::vector <RigWellResultPoint> >& m_pipeBranchesCellIds;
std::vector< std::vector<std::pair <size_t, double> > > m_accConnectionFlowPrBranch;
};
//--------------------------------------------------------------------------------------------------
//--------------------------------------------------------------------------------------------------
///
@ -250,30 +129,6 @@ void RimWellAllocationPlot::deleteViewWidget()
//--------------------------------------------------------------------------------------------------
void RimWellAllocationPlot::updateFromWell()
{
if (!m_case) return;
const RigSingleWellResultsData* wellResults = nullptr;
wellResults = m_case->reservoirData()->findWellResult(m_wellName);
if (!wellResults) return;
const RigWellResultFrame* wellResultFrame = nullptr;
std::vector< std::vector <cvf::Vec3d> > pipeBranchesCLCoords;
std::vector< std::vector <RigWellResultPoint> > pipeBranchesCellIds;
RigSimulationWellCenterLineCalculator::calculateWellPipeCenterlineFromWellFrame(m_case->reservoirData(),
wellResults,
m_timeStep,
true,
true,
pipeBranchesCLCoords,
pipeBranchesCellIds);
accumulatedWellFlowPlot()->setDescription("Accumulated Well Flow (" + m_wellName + ")");
RigAccWellFlowCalculator wfCalculator(pipeBranchesCLCoords, pipeBranchesCellIds);
// Delete existing tracks
{
std::vector<RimWellLogTrack*> tracks;
@ -288,9 +143,75 @@ void RimWellAllocationPlot::updateFromWell()
CVF_ASSERT(accumulatedWellFlowPlot()->trackCount() == 0);
accumulatedWellFlowPlot()->setDescription("Accumulated Well Flow (" + m_wellName + ")");
if (!m_case) return;
const RigSingleWellResultsData* wellResults = nullptr;
wellResults = m_case->reservoirData()->findWellResult(m_wellName);
if (!wellResults) return;
std::vector< std::vector <cvf::Vec3d> > pipeBranchesCLCoords;
std::vector< std::vector <RigWellResultPoint> > pipeBranchesCellIds;
RigSimulationWellCenterLineCalculator::calculateWellPipeCenterlineFromWellFrame(m_case->reservoirData(),
wellResults,
m_timeStep,
true,
true,
pipeBranchesCLCoords,
pipeBranchesCellIds);
std::unique_ptr< RigAccWellFlowCalculator > wfCalculator;
std::map<QString, const std::vector<double>* > tracerCellFractionValues;
if ( m_flowDiagSolution )
{
RimFlowDiagSolution::TracerStatusType requestedTracerType = RimFlowDiagSolution::UNDEFINED;
const RigWellResultFrame& wellResultFrame = wellResults->wellResultFrame(m_timeStep);
if (wellResultFrame.m_productionType == RigWellResultFrame::PRODUCER)
{
requestedTracerType = RimFlowDiagSolution::INJECTOR;
}
else
{
requestedTracerType = RimFlowDiagSolution::PRODUCER;
}
std::vector<QString> tracerNames = m_flowDiagSolution->tracerNames();
for ( const QString& tracerName : tracerNames )
{
if (m_flowDiagSolution->tracerStatusInTimeStep(tracerName, m_timeStep) == requestedTracerType)
{
RigFlowDiagResultAddress resAddr(RIG_FLD_CELL_FRACTION_RESNAME, tracerName.toStdString());
const std::vector<double>* tracerCellFractions = m_flowDiagSolution->flowDiagResults()->resultValues(resAddr, m_timeStep);
tracerCellFractionValues[tracerName] = tracerCellFractions;
}
}
RigEclCellIndexCalculator cellIdxCalc(m_case->reservoirData()->mainGrid(), m_case->reservoirData()->activeCellInfo(RifReaderInterface::MATRIX_RESULTS));
wfCalculator.reset(new RigAccWellFlowCalculator(pipeBranchesCLCoords,
pipeBranchesCellIds,
tracerCellFractionValues,
cellIdxCalc));
}
else
{
wfCalculator.reset(new RigAccWellFlowCalculator(pipeBranchesCLCoords,
pipeBranchesCellIds));
}
size_t branchCount = pipeBranchesCLCoords.size();
for (size_t brIdx = 0; brIdx < branchCount; ++brIdx)
{
// Skip Tiny dummy branches
if (pipeBranchesCellIds[brIdx].size() <= 3) continue;
RimWellLogTrack* plotTrack = new RimWellLogTrack();
// TODO: Description is overwritten by RimWellLogPlot::updateTrackNames()
@ -299,54 +220,48 @@ void RimWellAllocationPlot::updateFromWell()
accumulatedWellFlowPlot()->addTrack(plotTrack);
#if 0
std::vector<cvf::Vec3d> curveCoords;
std::vector<double> flowRate;
if ( m_flowDiagSolution )
{
std::vector<cvf::Vec3d> branchCoords = pipeBranchesCLCoords[brIdx];
std::vector<RigWellResultPoint> branchResultPoints = pipeBranchesCellIds[brIdx];
RigWellResultPoint prevResultPoint;
for (size_t i = 0; i < branchResultPoints.size(); i++)
std::vector<double> connNumbers;
{
const RigWellResultPoint& resultPoint = branchResultPoints[i];
const std::vector<size_t>& connNumbersSize_t = wfCalculator->connectionNumbersFromTop(brIdx);
for ( size_t conNumb : connNumbersSize_t ) connNumbers.push_back(static_cast<double>(conNumb));
}
if (i > 0 && prevResultPoint.m_gridCellIndex != resultPoint.m_gridCellIndex)
{
// Add an extra curve point when a cell transition is detected
curveCoords.push_back(branchCoords[i]);
flowRate.push_back(prevResultPoint.m_flowRate);
}
std::vector<QString> tracerNames = wfCalculator->tracerNames();
for (const QString& tracerName: tracerNames)
{
std::vector<double> accFlow;
curveCoords.push_back(branchCoords[i]);
flowRate.push_back(branchResultPoints[i].m_flowRate);
accFlow = wfCalculator->accumulatedTracerFlowPrConnection(tracerName, brIdx);
prevResultPoint = resultPoint;
RimWellFlowRateCurve* curve = new RimWellFlowRateCurve;
curve->setFlowValues(tracerName, connNumbers, accFlow);
plotTrack->addCurve(curve);
curve->loadDataAndUpdate();
}
}
RigSimulationWellCoordsAndMD helper(curveCoords);
#endif
const std::vector<std::pair <size_t, double> >& flowPrConnection = wfCalculator.accumulatedFlowPrConnection(brIdx);
std::vector<double> connNumbers;
std::vector<double> accFlow;
for (const std::pair<size_t, double> & flowPair: flowPrConnection)
else
{
connNumbers.push_back(flowPair.first);
accFlow.push_back(flowPair.second);
std::vector<double> connNumbers;
std::vector<double> accFlow;
accFlow = wfCalculator->accumulatedTotalFlowPrConnection(brIdx);
const std::vector<size_t>& connNumbersSize_t = wfCalculator->connectionNumbersFromTop(brIdx);
for ( size_t conNumb : connNumbersSize_t ) connNumbers.push_back(static_cast<double>(conNumb));
RimWellFlowRateCurve* curve = new RimWellFlowRateCurve;
curve->setFlowValues("Total", connNumbers, accFlow);
plotTrack->addCurve(curve);
curve->loadDataAndUpdate();
}
RimWellFlowRateCurve* curve = new RimWellFlowRateCurve;
curve->setFlowValues(connNumbers, accFlow);
plotTrack->addCurve(curve);
curve->loadDataAndUpdate();
}
setDescription("Well Allocation (" + m_wellName + ")");

View File

@ -76,7 +76,7 @@ QString RimWellFlowRateCurve::wellName() const
//--------------------------------------------------------------------------------------------------
QString RimWellFlowRateCurve::wellLogChannelName() const
{
return "Flow Rate";
return "AccumulatedFlowRate";
}
//--------------------------------------------------------------------------------------------------
@ -84,8 +84,7 @@ QString RimWellFlowRateCurve::wellLogChannelName() const
//--------------------------------------------------------------------------------------------------
QString RimWellFlowRateCurve::createCurveAutoName()
{
return QString("%1 (%2)").arg(wellLogChannelName()).arg(wellName());
return m_tracerName;
}
//--------------------------------------------------------------------------------------------------
@ -127,9 +126,11 @@ RimWellAllocationPlot* RimWellFlowRateCurve::wellAllocationPlot() const
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
void RimWellFlowRateCurve::setFlowValues(const std::vector<double>& measuredDepths, const std::vector<double>& flowRates)
void RimWellFlowRateCurve::setFlowValues(const QString& tracerName, const std::vector<double>& measuredDepths, const std::vector<double>& flowRates)
{
m_curveData = new RigWellLogCurveData;
m_curveData->setValuesAndMD(flowRates, measuredDepths, RimDefines::UNIT_METER, false);
m_tracerName = tracerName;
}

View File

@ -38,7 +38,7 @@ public:
RimWellFlowRateCurve();
virtual ~RimWellFlowRateCurve();
void setFlowValues(const std::vector<double>& measuredDepths, const std::vector<double>& flowRates);
void setFlowValues(const QString& tracerName , const std::vector<double>& measuredDepths, const std::vector<double>& flowRates);
virtual QString wellName() const override;
virtual QString wellLogChannelName() const override;
@ -49,5 +49,7 @@ protected:
private:
RimWellAllocationPlot* wellAllocationPlot() const;
QString m_tracerName;
};

View File

@ -27,6 +27,7 @@ ${CEE_CURRENT_LIST_DIR}RigFlowDiagSolverInterface.h
${CEE_CURRENT_LIST_DIR}RigFlowDiagInterfaceTools.h
${CEE_CURRENT_LIST_DIR}RigFlowDiagStatCalc.h
${CEE_CURRENT_LIST_DIR}RigFlowDiagVisibleCellsStatCalc.h
${CEE_CURRENT_LIST_DIR}RigAccWellFlowCalculator.h
${CEE_CURRENT_LIST_DIR}RigWellLogExtractor.h
${CEE_CURRENT_LIST_DIR}RigEclipseWellLogExtractor.h
${CEE_CURRENT_LIST_DIR}RigLocalGrid.h
@ -73,6 +74,7 @@ ${CEE_CURRENT_LIST_DIR}RigFlowDiagResultFrames.cpp
${CEE_CURRENT_LIST_DIR}RigFlowDiagSolverInterface.cpp
${CEE_CURRENT_LIST_DIR}RigFlowDiagStatCalc.cpp
${CEE_CURRENT_LIST_DIR}RigFlowDiagVisibleCellsStatCalc.cpp
${CEE_CURRENT_LIST_DIR}RigAccWellFlowCalculator.cpp
${CEE_CURRENT_LIST_DIR}RigWellLogExtractor.cpp
${CEE_CURRENT_LIST_DIR}RigEclipseWellLogExtractor.cpp
${CEE_CURRENT_LIST_DIR}RigLocalGrid.cpp

View File

@ -0,0 +1,231 @@
/////////////////////////////////////////////////////////////////////////////////
//
// Copyright (C) 2017 Statoil ASA
//
// ResInsight is free software: you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// ResInsight is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or
// FITNESS FOR A PARTICULAR PURPOSE.
//
// See the GNU General Public License at <http://www.gnu.org/licenses/gpl.html>
// for more details.
//
/////////////////////////////////////////////////////////////////////////////////
#include "RigAccWellFlowCalculator.h"
#include "RigSingleWellResultsData.h"
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
RigAccWellFlowCalculator::RigAccWellFlowCalculator(const std::vector< std::vector <cvf::Vec3d> >& pipeBranchesCLCoords, const std::vector< std::vector <RigWellResultPoint> >& pipeBranchesCellIds, const std::map<QString, const std::vector<double>* >& tracerCellFractionValues, const RigEclCellIndexCalculator cellIndexCalculator):
m_pipeBranchesCLCoords(pipeBranchesCLCoords),
m_pipeBranchesCellIds(pipeBranchesCellIds),
m_tracerCellFractionValues(&tracerCellFractionValues),
m_cellIndexCalculator(cellIndexCalculator)
{
m_accConnectionFlowPrBranch.resize(m_pipeBranchesCellIds.size());
for ( const auto& it: (*m_tracerCellFractionValues) ) m_tracerNames.push_back(it.first);
calculateAccumulatedFlowPrConnection(0, 1);
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
RigAccWellFlowCalculator::RigAccWellFlowCalculator(const std::vector< std::vector <cvf::Vec3d> >& pipeBranchesCLCoords, const std::vector< std::vector <RigWellResultPoint> >& pipeBranchesCellIds):
m_pipeBranchesCLCoords(pipeBranchesCLCoords),
m_pipeBranchesCellIds(pipeBranchesCellIds),
m_tracerCellFractionValues(nullptr),
m_cellIndexCalculator(RigEclCellIndexCalculator(nullptr, nullptr))
{
m_accConnectionFlowPrBranch.resize(m_pipeBranchesCellIds.size());
m_tracerNames.push_back("GrandTotalOnly");
calculateAccumulatedFlowPrConnection(0, 1);
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
const std::vector<double>& RigAccWellFlowCalculator::accumulatedTotalFlowPrConnection(size_t branchIdx)
{
CVF_ASSERT(m_accConnectionFlowPrBranch[branchIdx].accConnFlowFractionsPrTracer.find("GrandTotalOnly") != m_accConnectionFlowPrBranch[branchIdx].accConnFlowFractionsPrTracer.end());
return m_accConnectionFlowPrBranch[branchIdx].accConnFlowFractionsPrTracer["GrandTotalOnly"];
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
const std::vector<double>& RigAccWellFlowCalculator::accumulatedTracerFlowPrConnection(const QString& tracerName, size_t branchIdx)
{
CVF_ASSERT(m_accConnectionFlowPrBranch[branchIdx].accConnFlowFractionsPrTracer.find(tracerName) != m_accConnectionFlowPrBranch[branchIdx].accConnFlowFractionsPrTracer.end());
return m_accConnectionFlowPrBranch[branchIdx].accConnFlowFractionsPrTracer[tracerName];
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
const std::vector<size_t>& RigAccWellFlowCalculator::connectionNumbersFromTop(size_t branchIdx)
{
return m_accConnectionFlowPrBranch[branchIdx].connectionNumbersFromTop;
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
void RigAccWellFlowCalculator::calculateAccumulatedFlowPrConnection(size_t branchIdx, size_t startConnectionNumberFromTop)
{
const std::vector<RigWellResultPoint>& branchCells = m_pipeBranchesCellIds[branchIdx];
std::vector<size_t> resPointToConnectionIndexFromBottom = wrpToConnectionIndexFromBottom(branchCells);
size_t prevConnIndx = -1;
int clSegIdx = static_cast<int>(branchCells.size()) - 1;
std::map<QString, std::vector<double> >& accConnFlowFractionsPrTracer = m_accConnectionFlowPrBranch[branchIdx].accConnFlowFractionsPrTracer;
std::vector<size_t>& connNumbersFromTop = m_accConnectionFlowPrBranch[branchIdx].connectionNumbersFromTop;
std::vector<double> accFlow;
accFlow.resize(m_tracerNames.size(), 0.0);
while ( clSegIdx >= 0 )
{
// Skip point if referring to the same cell as the previous centerline segment did
{
if ( resPointToConnectionIndexFromBottom[clSegIdx] == prevConnIndx )
{
--clSegIdx;
continue;
}
prevConnIndx = resPointToConnectionIndexFromBottom[clSegIdx];
}
// Accumulate the connection-cell's fraction flows
if ( m_tracerCellFractionValues )
{
size_t resCellIndex = m_cellIndexCalculator.resultCellIndex(branchCells[clSegIdx].m_gridIndex,
branchCells[clSegIdx].m_gridCellIndex);
size_t tracerIdx = 0;
for ( const auto & tracerFractionIt: (*m_tracerCellFractionValues) )
{
accFlow[tracerIdx] += (*tracerFractionIt.second)[resCellIndex] * branchCells[clSegIdx].flowRate();
tracerIdx++;
}
}
else
{
accFlow[0] += branchCells[clSegIdx].flowRate();
}
// Add the total accumulated (fraction) flows from any branches connected to this cell
size_t connNumFromTop = connectionIndexFromTop(resPointToConnectionIndexFromBottom, clSegIdx) + startConnectionNumberFromTop;
std::vector<size_t> downstreamBranches = findDownstreamBranchIdxs(branchCells[clSegIdx]);
for ( size_t dsBidx : downstreamBranches )
{
if ( dsBidx != branchIdx && m_accConnectionFlowPrBranch[dsBidx].connectionNumbersFromTop.size() == 0 ) // Not this branch or already calculated
{
calculateAccumulatedFlowPrConnection(dsBidx, connNumFromTop);
BranchResult& accConnFlowFractionsDsBranch = m_accConnectionFlowPrBranch[dsBidx];
size_t tracerIdx = 0;
for ( const auto & tracerName: m_tracerNames )
{
accFlow[tracerIdx] += accConnFlowFractionsDsBranch.accConnFlowFractionsPrTracer[tracerName].back();
tracerIdx++;
}
}
}
// Push back the accumulated result into the storage
size_t tracerIdx = 0;
for ( const auto & tracerName: m_tracerNames )
{
accConnFlowFractionsPrTracer[tracerName].push_back(accFlow[tracerIdx]);
tracerIdx++;
}
connNumbersFromTop.push_back(connNumFromTop);
--clSegIdx;
}
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
std::vector<size_t> RigAccWellFlowCalculator::wrpToConnectionIndexFromBottom(const std::vector<RigWellResultPoint> &branchCells)
{
std::vector<size_t> resPointToConnectionIndexFromBottom;
resPointToConnectionIndexFromBottom.resize(branchCells.size(), -1);
size_t connIdxFromBottom = 0;
int clSegIdx = static_cast<int>(branchCells.size()) - 1;
if (clSegIdx < 0) return resPointToConnectionIndexFromBottom;
size_t prevGridIdx = branchCells[clSegIdx].m_gridIndex;
size_t prevGridCellIdx = branchCells[clSegIdx].m_gridCellIndex;
int prevErtSegId = branchCells[clSegIdx].m_ertSegmentId;
int prevErtBranchId = branchCells[clSegIdx].m_ertBranchId;
while ( clSegIdx >= 0 )
{
if ( branchCells[clSegIdx].m_gridIndex != prevGridIdx
|| branchCells[clSegIdx].m_gridCellIndex != prevGridCellIdx
|| branchCells[clSegIdx].m_ertSegmentId != prevErtSegId
|| branchCells[clSegIdx].m_ertBranchId != prevErtBranchId)
{
++connIdxFromBottom;
prevGridIdx = branchCells[clSegIdx].m_gridIndex ;
prevGridCellIdx = branchCells[clSegIdx].m_gridCellIndex;
}
resPointToConnectionIndexFromBottom[clSegIdx] = connIdxFromBottom;
--clSegIdx;
}
return resPointToConnectionIndexFromBottom;
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
size_t RigAccWellFlowCalculator::connectionIndexFromTop(const std::vector<size_t>& resPointToConnectionIndexFromBottom, size_t clSegIdx)
{
return resPointToConnectionIndexFromBottom.front() - resPointToConnectionIndexFromBottom[clSegIdx];
}
//--------------------------------------------------------------------------------------------------
///
//--------------------------------------------------------------------------------------------------
std::vector<size_t> RigAccWellFlowCalculator::findDownstreamBranchIdxs(const RigWellResultPoint& connectionPoint)
{
std::vector<size_t> downStreamBranchIdxs;
for ( size_t bIdx = 0; bIdx < m_pipeBranchesCellIds.size(); ++bIdx )
{
if ( m_pipeBranchesCellIds[bIdx][0].m_gridIndex == connectionPoint.m_gridIndex
&& m_pipeBranchesCellIds[bIdx][0].m_gridCellIndex == connectionPoint.m_gridCellIndex
&& m_pipeBranchesCellIds[bIdx][0].m_ertBranchId == connectionPoint.m_ertBranchId
&& m_pipeBranchesCellIds[bIdx][0].m_ertSegmentId == connectionPoint.m_ertSegmentId)
{
downStreamBranchIdxs.push_back(bIdx);
}
}
return downStreamBranchIdxs;
}

View File

@ -0,0 +1,96 @@
/////////////////////////////////////////////////////////////////////////////////
//
// Copyright (C) 2017 Statoil ASA
//
// ResInsight is free software: you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// ResInsight is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or
// FITNESS FOR A PARTICULAR PURPOSE.
//
// See the GNU General Public License at <http://www.gnu.org/licenses/gpl.html>
// for more details.
//
/////////////////////////////////////////////////////////////////////////////////
#pragma once
//==================================================================================================
///
///
//==================================================================================================
#include "RigActiveCellInfo.h"
#include "RigMainGrid.h"
#include "RigFlowDiagResults.h"
struct RigWellResultPoint;
class RigEclCellIndexCalculator
{
public:
RigEclCellIndexCalculator(const RigMainGrid* mainGrid, const RigActiveCellInfo* activeCellInfo)
: m_mainGrid(mainGrid), m_activeCellInfo(activeCellInfo)
{}
size_t resultCellIndex(size_t gridIndex, size_t gridCellIndex) const
{
const RigGridBase* grid = m_mainGrid->gridByIndex(gridIndex);
size_t reservoirCellIndex = grid->reservoirCellIndex(gridCellIndex);
return m_activeCellInfo->cellResultIndex(reservoirCellIndex);
}
private:
const RigMainGrid* m_mainGrid;
const RigActiveCellInfo* m_activeCellInfo;
};
//==================================================================================================
///
///
//==================================================================================================
class RigAccWellFlowCalculator
{
public:
RigAccWellFlowCalculator(const std::vector< std::vector <cvf::Vec3d> >& pipeBranchesCLCoords,
const std::vector< std::vector <RigWellResultPoint> >& pipeBranchesCellIds,
const std::map<QString, const std::vector<double>* >& tracerCellFractionValues,
const RigEclCellIndexCalculator cellIndexCalculator);
RigAccWellFlowCalculator(const std::vector< std::vector <cvf::Vec3d> >& pipeBranchesCLCoords,
const std::vector< std::vector <RigWellResultPoint> >& pipeBranchesCellIds);
const std::vector<double>& accumulatedTotalFlowPrConnection( size_t branchIdx);
const std::vector<double>& accumulatedTracerFlowPrConnection(const QString& tracerName, size_t branchIdx);
const std::vector<size_t>& connectionNumbersFromTop(size_t branchIdx);
const std::vector<QString>& tracerNames() { return m_tracerNames;}
private:
void calculateAccumulatedFlowPrConnection( size_t branchIdx, size_t startConnectionNumberFromTop);
std::vector<size_t> wrpToConnectionIndexFromBottom( const std::vector<RigWellResultPoint> &branchCells);
static size_t connectionIndexFromTop( const std::vector<size_t>& resPointToConnectionIndexFromBottom, size_t clSegIdx);
std::vector<size_t> findDownstreamBranchIdxs( const RigWellResultPoint& connectionPoint);
const std::vector< std::vector <cvf::Vec3d> >& m_pipeBranchesCLCoords;
const std::vector< std::vector <RigWellResultPoint> >& m_pipeBranchesCellIds;
const std::map<QString, const std::vector<double>* >* m_tracerCellFractionValues;
RigEclCellIndexCalculator m_cellIndexCalculator;
std::vector<QString> m_tracerNames;
struct BranchResult
{
std::vector<size_t> connectionNumbersFromTop;
std::map<QString, std::vector<double> > accConnFlowFractionsPrTracer;
};
std::vector< BranchResult > m_accConnectionFlowPrBranch;
};

View File

@ -58,6 +58,18 @@ struct RigWellResultPoint
return isCell() || isPointValid();
}
double flowRate() const
{
if ( isCell() && m_isOpen)
{
return m_flowRate;
}
else
{
return 0.0;
}
}
size_t m_gridIndex;
size_t m_gridCellIndex; //< Index to cell which is included in the well