GeosLogReaderUtils functions

This package define functions dedicated to the GeosLogReader.

geos.pv.geosLogReaderUtils.GeosLogReaderAquifers module

class geos.pv.geosLogReaderUtils.GeosLogReaderAquifers.GeosLogReaderAquifers(filepath, propertiesUnit)[source]

Bases: object

Reader for Aquifer.

Parameters:
  • filepath (str) – path to geos log file.

  • propertiesUnit (dict[str, Unit]) – unit preferences

calculateExtraValues()[source]

Add cumulated columns for each aquifer volume and aquifer rate.

createDataframe()[source]

Create and fill and return dataframeAquifers.

Returns:

dataframe with values from Geos log.

Return type:

pd.DataFrame

readAll(filepath)[source]

Initialize all the attributes of the class by reading a Geos log file.

Parameters:

filepath (str) – Geos log filepath.

readAquiferNames(file)[source]

Initialize the m_aquiferNames attribute by reading log file.

Parameters:

file (TextIOBase) – Geos Log file

Returns:

The last line with time info read. The id of the last line read that contained the tag “_pressureInfluence_table”., which will be the line containing the first positive timestep at 0s.

Return type:

tuple(str, int)

readPropertiesValues(file, line, id_line, total_lines)[source]

Read aquifer property values from geos log file.

Initialize the m_aquifersPropertiesValues and m_timesteps attributes by reading the Geos log. If a timestep contains the tag m_computeStatisticsName, the current timestep is added to m_timesteps and we recover the property values in m_regionsPropertiesValues.

Parameters:
  • file (TextIOBase) – Geos Log file

  • line (str) – last line read in the file.

  • id_line (int) – The id of the last line read in readPhaseNames.

  • total_lines (int) – The number of lines in the file.

geos.pv.geosLogReaderUtils.geosLogReaderConvergence module

class geos.pv.geosLogReaderUtils.GeosLogReaderConvergence.GeosLogReaderConvergence(filepath, propertiesUnit)[source]

Bases: object

Reader for Convergence information.

Parameters:
  • filepath (str) – path to geos log file.

  • propertiesUnit (dict[str, Unit]) – unit preferences

calculateExtraValues()[source]

Add cumulated columns for newtonIter and linearIter.

createDataframe()[source]

Create and fill and return dataframeSolversIterations.

Returns:

dataframe with values from Geos log.

Return type:

pd.DataFrame

readAll(filepath)[source]

Initialize all the attributes of the class by reading a Geos log file.

Parameters:

filepath (str) – Geos log filepath.

readIterationsValues(file, total_lines)[source]

Read iteration values from Geos log file.

Initialize the m_aquifersPropertiesValues and m_timesteps attributes by reading the Geos log. If a timestep contains the tag m_computeStatisticsName, the current timestep is added to m_timesteps and we recover the property values in m_regionsPropertiesValues.

Parameters:
  • file (TextIOBase) – Geos Log file

  • total_lines (int) – The number of lines in the file.

geos.pv.geosLogReaderUtils.GeosLogReaderFlow module

class geos.pv.geosLogReaderUtils.GeosLogReaderFlow.GeosLogReaderFlow(filepath, propertiesUnit, phaseNames=None)[source]

Bases: object

A reader that reads .txt and .out files containing Geos logs.

To do that, we use specific tags in the current version of this code. Supposed tags are:

  • for region names: “Adding Object CellElementRegion” Supposed

    line:”Adding Object CellElementRegion named Reservoir from ObjectManager::Catalog”.

  • for phase names: “phaseModel” Supposed line: “ TableFunction:

    fluid_phaseModel1_PhillipsBrineDensity_table”.

  • for timesteps: “Time:” Supposed line: “Time: 0s, dt:100s, Cycle: 0”.

  • for CFL properties: “CFL”. Supposed line: “compflowStatistics: Max

    phase CFL number: 0.00696878”

Another important tag in the log will be the name of the flow statistics model used to output 2D data in the Geos Log. This one will be found automatically. The flow statistics model that can output flow data are:

  • “SinglePhaseStatistics”.

  • “CompositionalMultiphaseStatistics”.

Parameters:
  • filepath (str) – path to Geos log file

  • propertiesUnit (dict[str, Unit]) – unit preferences

  • phaseNames (list[str], optional) –

    Name of the phases.

    Defaults to [].

createDataframe()[source]

Create and fill and return dataframeFlow.

Returns:

dataframe with values from Geos log.

Return type:

pd.DataFrame

readAll(filepath)[source]

Initialize all the attributes of the class by reading a Geos log file.

Parameters:

filepath (str) – Geos log filepath.

readComputeStatisticsName(file, id_line, total_lines)[source]

Read flow statistics from the Geos log file.

Parameters:
  • file (TextIOBase) – Geos Log file

  • id_line (int) – The id of the last line read in readPhaseNames.

  • total_lines (int) – total number of lines in the file.

Returns:

Tuple containingt the id of the last line read and the line.

Return type:

tuple[int, str]

readPropertiesValues(file, id_line, total_lines, lineTagStats)[source]

Read property values from Geos log file.

Initialize the m_regionsPropertiesValues and m_timesteps attributes by reading the Geos log. If a timestep contains the tag m_computeStatisticsName, the current timestep is added to m_timesteps and we recover the property values in m_regionsPropertiesValues.

Parameters:
  • file (TextIOBase) – Geos Log file

  • id_line (int) – The id of the last line read in readPhaseNames.

  • total_lines (int) – The number of lines in the file.

  • lineTagStats (str) – The first line containing the tag of the flow statistics model.

readRegionNames(file)[source]

Initialize the m_regionNames attribute by reading log file.

Parameters:

file (TextIOBase) – Geos Log file

Returns:

The id of the last line read that contained the tag “Adding Object CellElementRegion”

Return type:

int

geos.pv.geosLogReaderUtils.GeosLogReaderFunctions module

Functions to read and process Geos log.

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.buildPropertiesNameForComponents(phasesName)[source]

Builds the list of component property names from the list of phases name.

Parameters:

phasesName (list) – [“CO2”,”Water”]

Returns:

[‘Dissolved mass CO2 in CO2’,’Dissolved mass Water in CO2’, ‘Dissolved mass CO2 in Water’,’Dissolved mass Water in Water’]

Return type:

list

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.buildPropertiesNameForPhases(nameBlock, phasesName)[source]

Replace phase by phase names.

Parameters:
  • nameBlock (str) – “ Mobile phase mass”

  • phasesName (list[str]) – [“CO2”,”Water”]

Returns:

[‘Mobile CO2 mass’, ‘Mobile Water mass’]

Return type:

list[str]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.buildPropertiesNameFromGeosProperties(geosProperties, phasesName)[source]

Extracts the property name and its extensions like min, max, average.

Parameters:
  • geosProperties (str) – “ Delta pressure (min, max)”

  • phasesName (list[str]) – [“CO2”,”Water”]

Returns:

[” Delta pressure min”, “ Delta pressure max”]

Return type:

list[str]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.buildPropertiesNameNoPhases(nameBlock, extensions='')[source]

From a name block and extensions, builds a list of properties name.

Parameters:
  • nameBlock (str) – “ Delta pressure “

  • extensions (str) – “min, max)”

Returns:

[” Delta pressure min”, “ Delta pressure max”]

Return type:

list

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.convertValues(propertyNames, propertyValues, propertiesUnit)[source]

Convert properties to the desired units.

Knowing two lists : 1) float numbers that are supposed to be in SI units ; 2) properties name linked to the float numbers. And that these lists are of same dimension, creates a new list of float values converted to a specific unit linked to the property name.

Parameters:
  • propertyNames (list[str]) – list of property names

  • propertyValues (list[float]) – list of property values.

  • propertiesUnit (dict[str, Unit]) – dictionary of desired units for each property {“pressure”: UnitPressure, …, “propertyTypeN”: UnitPropertyTypeN}

Returns:

list of converted values.

Return type:

list[float]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.correctZeroValuesInListOfValues(values)[source]

Replace orhphelin 0 values of input list.

If 0 values are found in a list of values, either replace them with the value found before in the list or keep it 0. We suppose that 2 zeros in a row correspond to a continuity of the values, hence keeping it 0.

Parameters:

values (list[float]) – list of ints or floats

Returns:

list of ints or floats

Return type:

list[float]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.countNumberLines(filepath)[source]

Reads a file to find the number of lines within it.

Parameters:

filepath (str) – Path to the file.

Returns:

Number of lines in file.

Return type:

int

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.elementsAreInLog(filepath, elements)[source]

Indicates if input file contains element from input list of string.

To do so, this reads a file and checks at every line if an element was found within the line. If an element is found, it is not checked again. The function returns True only when there is no more element to check.

Parameters:
  • filepath (str) – Path to the file.

  • elements (list[str]) – Every string that needs to be find inside the file.

Return type:

bool

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractAquifer(geosLogLine)[source]

Extracts the name of the aquifer from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “ TableFunction: aquifer1_pressureInfluence_table”

Raises:

ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

“aquifer1”

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractFirstIntFromString(string)[source]

Extracts the first int value from a string.

Parameters:

string (str) – A random string.

Returns:

int or None if no int was found.

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractFloatsFromString(line)[source]

Extracts a list of float numbers from a string.

Parameters:

line (str) – A random string.

Returns:

[float1, …, floatN]

Return type:

list[float]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractLinearIter(geosLogLine)[source]

From a Geos log line, extracts the int value of linear iterations.

Parameters:

geosLogLine (str) – #expected line : “ Linear Solver | Success | Iterations: 23 | Final Rel Res: 5.96636e-05 | Make Restrictor Time: 0 | Compute Auu Time: 0 | SC Filter Time: 0 | Setup Time: 1.5156 s | Solve Time: 0.041093 s”

Raises:
  • ValueError – “Not enough elements to unpack in time line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

23

Return type:

int

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractListIntsFromString(string)[source]

Builds a list of int numbers from a string.

Parameters:

string (str) – A random string.

Returns:

[int1, …, intN]

Return type:

list[int]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractNewtonIter(geosLogLine)[source]

From a Geos log line, extracts the int value of NewtonIter.

Parameters:

geosLogLine (str) – #expected line : “ Attempt: 0, ConfigurationIter: 0, NewtonIter: 0”

Raises:
  • ValueError – “Not enough elements to unpack in time line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

NewtonIter

Return type:

int

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractPhaseId(geosLogLine)[source]

Extracts the phase number id from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “wellControls1: Phase 0 surface volumetric rate: 30.023748128796043 sm3/s”

Raises:
  • ValueError – “Not enough elements to unpack in line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

0

Return type:

int

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractPhaseModel(geosLogLine)[source]

Extracts the name of a phase model from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “ TableFunction: fluid_phaseModel1_PhillipsBrineDensity_table”

Raises:
  • ValueError – “Not enough elements to unpack in line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

“PhillipsBrineDensity”

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractPropertiesFlow(geosLogLine, phasesName)[source]

Extracts flow property from a Geos log line.

Parameters:
  • geosLogLine (str) – #expected line : “compflowStatistics, Reservoir: Delta pressure (min, max): 0, 0 Pa”

  • phasesName (list[str]) – [“CO2”,”Water”]

Raises:
  • ValueError – “Not enough elements to unpack in line <<” + geosLogLine + “>>.”

  • ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

[“Reservoir__DeltaPressureMin”, “Reservoir__DeltaPressureMax”]

Return type:

list[str]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractPropertiesWell(geosLogLine, wellName, phaseName)[source]

Extracts the well property presented from a Geos log line.

Parameters:
  • geosLogLine (str) – “wellControls1: Phase 0 surface volumetric rate: 30.023748128796043 sm3/s”

  • wellName (str) – well1

  • phaseName (str) – phase1

Returns:

[“Well1_SurfaceVolumetricRatePhase1”]

Return type:

list[str]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractRegion(geosLogLine)[source]

Extracts the name of the region from a Geos log line.

Parameters:
  • geosLogLine (str) – #expected line : “Adding Object CellElementRegion

  • ObjectManager::Catalog." (named Reservoir from)

Raises:
  • ValueError – “Not enough elements to unpack in line.”

  • ValueError – “An error has occured while parsing a line.”

Returns:

“Reservoir”

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractStatsName(geosLogLine)[source]

Extracts the name of the computed statistics name from a Geos log line.

Parameters:

geosLogLine (str) – #expected line :”compflowStatistics, Reservoir: Pressure (min, average, max): 2.86419e+07, 2.93341e+07, 3.006e+07 Pa”

Returns:

“compflowStatistics”

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractTimeAndDt(geosLogLine)[source]

From a Geos log line, extracts the float values of Time and dt.

Parameters:

geosLogLine (str) – #expected lines : “Time: {} years, {} days, {} hrs, {} min, {} s, dt: {} s, Cycle: {}” “Time: {:.2f} years, dt: {} s, Cycle: {}” “Time: {:.2f} days, dt: {} s, Cycle: {}” “Time: {:.2f} hrs, dt: {} s, Cycle: {}” “Time: {:.2f} min, dt: {} s, Cycle: {}” “Time: {:4.2e} s, dt: {} s, Cycle: {}” “Time: {]s, dt:{}s, Cycle: {}”

Raises:

KeyError – “Cannot add time values for tag=<<timeFactor>>”

Returns:

(time, dt)

Return type:

tuple[float]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractValueAndNameAquifer(geosLogLine)[source]

Extract value and name of the aquifer contained in a Geos log line.

Parameters:

geosLogLine (str) – “FlowSolverBase compositionalMultiphaseFlow (SimuDeck_aquifer_pression_meme.xml, l.28): at time 100s, the <Aquifer> boundary condition ‘aquifer1’ produces a flux of -0.6181975187076816 kg (or moles if useMass=0).”

Returns:

a tuple with the name and the float value. e.g. (“aquifer1”, -0.6181975187076816)

Return type:

tuple[str, float]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractValuesFlow(geosLogLine)[source]

Extract values from a Geos log line.

Parameters:

geosLogLine (str) – Geos log line “compflowStatistics, Reservoir: ” “Dissolved component mass: { { 0, 0 }, { 0, -6.38235e+10 } } kg”

Returns:

list of values in the line. [0.0, 0.0, 0.0, -6.38235e+10]

Return type:

list[float]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractValuesWell(geosLogLine, numberProperties)[source]

Extract values from Geos log line and returns them as a list of floats.

The idea here is first to extract all floats values from the line. Now all of them are useful so we need to keep some of them. Luckily, only the last one or two floats are useful. And to determine that, we use the number of well properties found in the line which ranges from one to two.

Parameters:
  • geosLogLine (str) – “Rank 129: well.CO2010: The density of phase 0 at surface conditions is 1.86 kg/sm3.”

  • numberProperties (int) – number of well properties found in the line.

Returns:

value of the property. e.g. [1.86]

Return type:

list[float]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractWell(geosLogLine)[source]

Extracts the name of the well from a Geos log line.

Parameters:

geosLogLine (str) – #expected line : “ TableFunction: wellControls_ConstantBHP_table” Or “ TableFunction: wellControls_ConstantPhaseRate_table”

Raises:

ValueError – “An error has occured while parsing <<” + geosLogLine + “>>”

Returns:

“wellControls”

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.extractWellTags(geosLogLine)[source]

Extracts the list of well property tags available from a Geos log line.

Parameters:

geosLogLine (str) – line from geos log file.

Returns:

list of tags.

Return type:

list[str]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.findNumberPhasesSimulation(filepath)[source]

Find the number of phases from Geos log file.

Geos logs do not have explicit message telling you how many phases were used to perform the simulation, unlike regions, wells etc … Therefore, we need at least to identify the exact number of phases that can be find in this Geos log file to extract correctly properties regarding phase related data.

Parameters:

filepath (str) – Filepath to a Geos log file.

Returns:

The number of phases found in the Geos log.

Return type:

int

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.formatPropertyName(propertyName)[source]

Clean the string by replacing special characters and removing spaces.

Parameters:

propertyName (str) – name;of:the property

Returns:

NameOfTheProperty

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.identifyCurrentWell(geosLogLine, lastWellName)[source]

Identify the current name of the well rom a Geos log line.

Because properties values of wells can be output without specifying the name of the well which they belong to, we have to assume that the name of well for the current properties observed is either :

  • the name of the well that can be found inside the line

  • the name of the last well name found in former lines.

Parameters:
  • geosLogLine (str) – line from Geos log file #expected lines with well name : “Rank 18: well.CO2001: BHP (at the specified reference elevation): 19318538.400682557 Pa” Or “wellControls1: BHP (at the specified reference elevation): 12337146.157562563 Pa” #line with no well name : “The total rate is 0 kg/s, which corresponds to a total surface volumetric rate of 0 sm3/s”

  • lastWellName (str) – name of the last well found

Raises:

ValueError – “An error has occured while parsing <<geosLogLine>>.”

Returns:

“wellControls”

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.identifyProperties(properties)[source]

Identify properties and add identifer.

From a list of properties name, identifies each of them with a certain integer, to link it to a meta property by adding an id in front of the property name.

Parameters:

properties (list[str]) – [“CaprockPressureMax”, “CaprockPressureMin”]

Returns:

[1:”CaprockPressureMax”, 1:”CaprockPressureMin”]

Return type:

list[tuple[str, int]]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.isFloat(element)[source]

Check whether an element is float or not.

Parameters:

element (Any) – input number to test.

Returns:

True if the number is a float.

Return type:

bool

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.phaseNamesBuilder(numberPhases, phasesFromUser)[source]

Build phase names.

When creating phase names, the user can or cannot have defined his own phase names when reading the log. Therefore, whether phase names were provided or not, if we have N phases found in the log, a list of N phase names will be created, starting from phase0 up to phaseN.

Parameters:
  • numberPhases (int) – Number of phases found in the log file.

  • phasesFromUser (list[str]) – names chosen by the user, can be more or less than numberPhases.

Returns:

[nameFromUser0, nameFromUser1, …, phaseN-1, phaseN]

Return type:

list[str]

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.replaceSpecialCharactersWithWhitespace(sentence)[source]

Replace every special characters in a string with whitespaces.

Parameters:

sentence (str) – Random string “hi ‘(_there(‘’&*$^,:;’”

Returns:

“hi there “

Return type:

str

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.timeInSecond(timeCounter)[source]

Calculates the time in s from a dict of different time quantities.

Parameters:

timeCounter (dict[str, float]) – timeCounter: {“years”: x0, “days”: x0, “hrs”: x0, “min”: x0, “s”: x0}

Returns:

Sum in seconds of all time quantities.

Return type:

float

geos.pv.geosLogReaderUtils.geosLogReaderFunctions.transformUserChoiceToListPhases(userChoice)[source]

Get a list of phase name from the input string.

When using GeosLogReader, the user can choose the names of the phases to use. The wished format is to specify each name in the good, separated by whitespaces or either commas.

Parameters:

userChoice (str | None) – Output from EnterPhaseNames string vector widget.

Returns:

[phase0, phase1, …, phaseN]

Return type:

list[str]

geos.pv.geosLogReaderUtils.GeosLogReaderWells module

class geos.pv.geosLogReaderUtils.GeosLogReaderWells.GeosLogReaderWells(filepath, propertiesUnit, phaseNames=None, numberWellsForMean=1)[source]

Bases: object

Read for Wells from Geos log file.

To do that, we use specific tags in the current version of this code. Supposed tags are :

  • for well names“Adding object WellElementRegion”

    and “_ConstantBHP_table” Supposed lines: “Adding Object WellElementRegion named wellRegion1 from ObjectManager::Catalog”. “ TableFunction: wellControls1_ConstantBHP_table”.

  • for phase names: “phaseModel”

    Supposed line: “ TableFunction: fluid_phaseModel1_PhillipsBrineDensity_table”.

  • for timesteps: “Time:”

    Supposed line : “Time: 0s, dt:100s, Cycle: 0” When it comes to well properties, special tags are used : “ BHP “ ; “ total rate” ; “ total surface volumetric rate” ; “phase surface volumetric rate” ; “well is shut” ; “density of phase” ; “total fluid density”.

Parameters:
  • filepath (str) – path of Geos log file

  • propertiesUnit (dict[str, Unit]) – unit preferences

  • phaseNames (list[str] | None, optional) –

    Name of the phases.

    Defaults to None.

  • numberWellsForMean (int, optional) – Number of wells. Defaults to 1.

calculateMeanValues()[source]

Calculate mean values of all wells.

createDataframe()[source]

Create and fill and return dataframeWells.

Returns:

dataframe with log values.

Return type:

pd.DataFrame

initWellPropertiesValues()[source]

Initialize the m_wellPropertiesValues.

readAll(filepath)[source]

Initialize all the attributes of the class by reading a Geos log file.

Parameters:
  • filepath (str) – Geos log filepath.

  • singlephase (bool) – True if its a singlephase simulation, False if multiphase.

readPropertiesValues(file, id_line, total_lines)[source]

Read property values from Geos log file.

Initialize the m_regionsPropertiesValues and m_timesteps attributes by reading the Geos log. If a timestep contains the tag m_computeStatisticsName, the current timestep is added to m_timesteps and we recover the property values in m_regionsPropertiesValues.

Parameters:
  • file (TextIOBase) – Geos Log file

  • id_line (int) – The id of the last line read in readPhaseNames.

  • total_lines (int) – The number of lines in the file.

readWellNames(file)[source]

Read well names from Geos log file.

Parameters:
  • file (TextIOBase) – Geos Log file

  • id_line (int) – The id of the last line read in readPhaseNames.

Returns:

The id of the last line read that contains the tag “Adding Object WellElementRegion”.

Return type:

int