properties module

This file is part of PyFrac.

Created by Haseeb Zia on 03.04.17. Copyright (c) ECOLE POLYTECHNIQUE FEDERALE DE LAUSANNE, Switzerland, Geo-Energy Laboratory, 2016-2020. All rights reserved. See the LICENSE.TXT file for more details.

class properties.FluidProperties(viscosity=None, density=1000.0, rheology='Newtonian', turbulence=False, compressibility=0, n=None, k=None, T0=None)[source]

Bases: object

Class defining the fluid properties.

Parameters
  • viscosity (ndarray) – – viscosity of the fluid .

  • density (float) – – density of the fluid.

  • rheology (string) –

    – string specifying rheology of the fluid. Possible options:

    • ”Newtonian”

    • ”non-Newtonian”

  • turbulence (bool) – – turbulence flag. If true, turbulence will be taken into account

  • compressibility (float) – – the compressibility of the fluid.

viscosity

– Viscosity of the fluid (note its different from local viscosity, see fracture class for local viscosity).

Type

ndarray

muPrime

– 12 * viscosity (parallel plates viscosity factor).

Type

float

rheology

– string specifying rheology of the fluid. Possible options:

  • “Newtonian”

  • “Herschel-Bulkley” or “HBF”

  • “power-law” or “PLF”

Type

string

density

– density of the fluid.

Type

float

turbulence

– turbulence flag. If true, turbulence will be taken into account.

Type

bool

compressibility

– the compressibility of the fluid.

Type

float

n

– flow index of the Herschel-Bulkey fluid.

Type

float

k

– consistency index of the Herschel-Bulkey fluid.

Type

float

T0

– yield stress of the Herschel-Bulkey fluid.

Type

float

Mprime -- 2**
Type

n + 1) * (2 * n + 1

class properties.InjectionProperties(rate, mesh, source_coordinates=None, source_loc_func=None, sink_loc_func=None, sink_vel_func=None, rate_delayed_second_injpoint=None, delayed_second_injpoint_loc=None, initial_rate_delayed_second_injpoint=None, rate_delayed_inj_pt_func=None, delayed_second_injpoint_loc_func=None)[source]

Bases: object

Class defining the injection parameters.

Parameters
  • rate (ndarray) –

    – array specifying the time series (row 0) and the corresponding injection rates (row 1). The times are instant where the injection rate changes.

    Attention

    The first time should be zero. The corresponding injection rate would be taken for initialization of the fracture with an analytical solution, if required.

  • mesh (CartesianMesh) – – the CartesianMesh object defining mesh.

  • source_coordinates (ndarray) – – list or ndarray with a length of 2, specifying the x and y coordinates of the injection point. Not used if source_loc_func is provided (See below).

  • source_loc_func – – the source location function is used to get the elements in which the fluid is injected. It should take the x and y coordinates and return True or False depending upon if the source is present on these coordinates. This function is evaluated at each of the cell centre coordinates to determine if the cell is a source element. It should have to arguments (x, y) and return True or False. It is also called upon re-meshing to get the source elements on the coarse mesh.

injectionRate

– array specifying the time series (row 0) and the corresponding injection rates (row 1). The time series provide the time when the injection rate changes.

Type

ndarray

sourceCoordinates

– array with a single row and two columns specifying the x and y coordinate of the injection point coordinates. If there are more than one source elements, the average is taken to get an estimate injection cell at the center.

Type

ndarray

sourceElem

– the element(s) where the fluid is injected in the cartesian mesh.

Type

ndarray

sourceLocFunc

– the source location function is used to get the elements in which the fluid is injected. It should take the x and y coordinates and return True or False depending upon if the source is present on these coordinates. This function is evaluated at each of the cell centre coordinates to determine if the cell is a source element. It should have to arguments (x, y) and return True or False. It is also called upon re-meshing to get the source elements on the coarse mesh.

Type

function

sinkLocFunc

– see description of arguments.

Type

function

sink_vel_func

– see description of arguments.

Type

function

get_injection_rate(tm, frac)[source]

This function gives the current injection rate at all of the cells in the domain.

Parameters
  • tm (float) – – the time at which the injection rate is required.

  • frac (CartesianMesh) – – the Fracture object containing the mesh and the current fracture elements.

Returns

– an numpy array of the size of the mesh with injection rates in each of the cell

Return type

Qin (ndarray)

remesh(new_mesh, old_mesh)[source]

This function is called every time the domian is remeshed.

Parameters
  • new_mesh (CartesianMesh) – – the CartesianMesh object describing the new coarse mesh.

  • old_mesh (CartesianMesh) – – the CartesianMesh object describing the old mesh.

class properties.IterationProperties(itr_type='not initialized')[source]

Bases: object

This class stores performance data in the form of a tree.

Parameters

itr_type (string) – – currently, the following iterations are profiled: - ‘time step’ - ‘time step attempt’ - ‘same front’ - ‘extended front’ - ‘tip inversion’ - ‘tip width’ - ‘nonlinear system solve’ - ‘width constraint iteration’ - ‘linear system solve’ - ‘Brent method’

class properties.LabelProperties(variable, data_subset='whole mesh', projection='2D', use_latex=True)[source]

Bases: object

This class stores the labels of a plot figure.

class properties.LoadingProperties(displ_rate=0.0, loaded_elts=None)[source]

Bases: object

Class defining the mechanical loading properties

EltLoaded

– array of elements that are loaded.

Type

ndarray

displ_rate

– the rate at which the elements in the EltLoaded list are displaced due to the applied mechanical loading

Type

float

class properties.MaterialProperties(Mesh, Eprime, toughness=0.0, Carters_coef=0.0, confining_stress=0.0, grain_size=0.0, K1c_func=None, anisotropic_K1c=False, confining_stress_func=None, Carters_coef_func=None, TI_elasticity=False, Cij=None, free_surf=False, free_surf_depth=1e+300, TI_plane_angle=0.0, minimum_width=1e-06, pore_pressure=- 1e+100)[source]

Bases: object

Class defining the Material properties of the solid.

Parameters
  • Mesh (CartesianMesh) – – the CartesianMesh object describing the mesh.

  • Eprime (float) – – plain strain modulus.

  • Toughness (float) – – Linear-Elastic Plane-Strain Fracture Toughness.

  • Carters_coef (float) – – Carter’s leak off coefficient.

  • confining_stress (ndarray) – – in-situ confining stress field normal to fracture surface.

  • grain_size (float) – – the grain size of the rock; used to calculate the relative roughness.

  • K1c_func (function) – – the function giving the toughness on the domain. It takes one argument (angle) in case of anisotropic toughness and two arguments (x, y) in case of heterogeneous toughness. The function is also used to get the toughness if the domain is re-meshed.

  • anisotropic_K1c (bool) – – flag to specify if the fracture toughness is anisotropic.

  • confining_stress_func (function) – – the function giving the in-situ stress on the domain. It should takes two arguments (x, y) to give the stress on these coordinates. It is also used to get the stress if the domain is re-meshed.

  • Carters_coef_func (function) – – the function giving the in Carter’s leak off coefficient on the domain. It should takes two arguments (x, y) to give the coefficient on these coordinates. It is also used to get the leak off coefficient if the domain is re-meshed.

  • TI_elasticity (bool) – – if True, the medium is elastic transverse isotropic.

  • Cij (ndarray) – – the transverse isotropic stiffness matrix (in the canonical basis); needs to be provided if TI_elasticity=True.

  • free_surf (bool) – – the free surface flag. True if the effect of free surface is to be taken into account.

  • free_surf_depth (float) – – the depth of the fracture from the free surface.

  • TI_plane_angle (float) – – the angle of the plane of the fracture with respect to the free surface.

  • minimum_width (float) – – minimum width corresponding to the asperity of the material.

Eprime

– plain strain modulus.

Type

float

K1c

– Linear-Elastic Plane-Strain Toughness for each cell.

Type

ndarray

Kprime

– 4*(2/pi)**0.5 * K1c.

Type

ndarray

Cl

– Carter’s leak off coefficient.

Type

float

Cprime

– 2 * Carter’s leak off coefficient.

Type

ndarray

SigmaO

– in-situ confining stress field normal to fracture surface.

Type

ndarray

grainSize

– the grain size of the rock; used to calculate the relative roughness.

Type

float

anisotropic_K1c

– if True, the toughness is considered anisotropic.

Type

bool

Kc1

– the fracture toughness along the x-axis, in case it is anisotropic.

Type

float

TI_elasticity

– the flag specifying transverse isotropic elasticity.

Type

bool

Cij

– the transverse isotropic stiffness matrix (in the canonical basis).

Type

ndarray

freeSurf

– if True, the effect of free surface is to be taken into account.

Type

bool

FreeSurfDepth

– the depth of the fracture from the free surface.

Type

float

TI_PlaneAngle

– the angle of the plane of the fracture with respect to the free surface.

Type

float

K1cFunc

– the function giving the toughness on the domain. It takes one argument (angle) in case of anisotropic toughness and two arguments (x, y) in case of heterogeneous toughness. The function is also used to get the toughness if the domain is re-meshed.

Type

function

SigmaOFunc

– the function giving the in-situ stress on the domain. It should takes two arguments(x, y) to give the stress on these coordinates. It is also used to get the confining stress if the domain is re-meshed.

Type

function

ClFunc

– the function giving the in Carter’s leak off coefficient on the domain. It should takes two arguments (x, y) to give the coefficient on these coordinates. It is also used to get the leak off coefficient if the domain is re-meshed.

Type

function

remesh(mesh)[source]

This function evaluates the toughness, confining stress and leak off coefficient on the given mesh using the functions provided in the MaterialProperties object. It should be evaluated each time re-meshing is done.

Parameters

mesh (CartesianMesh) – – the CartesianMesh object describing the new mesh.

class properties.PlotProperties(color_map=None, line_color=None, line_style='-', line_width=1.0, line_style_anal='--', line_color_anal='r', interpolation='none', alpha=0.8, line_width_anal=None, text_size=None, disp_precision=3, mesh_color='yellowgreen', mesh_edge_color='grey', mesh_label_color='black', graph_scaling='linear', color_maps=None, colors_list=None, plot_legend=True, plot_FP_time=True, use_tex=False)[source]

Bases: object

This class stores the parameters used for plotting of the post-processed results

class properties.SimulationProperties(address=None)[source]

Bases: object

Class defining the simulation properties.

Parameters

address (str) – named ‘simul_param’. For the description of the arguments and there default values, see the py:module::default_SimParam .

tolFractFront

– tolerance for the fracture front loop.

Type

float

toleranceEHL

– tolerance for the Elastohydrodynamic solver.

Type

float

toleranceVStagnant

– tolerance on the velocity to decide if a cell is stagnant.

Type

float

toleranceProjection

– tolerance for projection iteration for anisotropic case

Type

float

maxFrontItrs

– maximum iterations to for the fracture front loop.

Type

int

maxSolverItrs

– maximum iterations for the EHL iterative solver (Picard-Newton hybrid) in this case.

Type

int

maxProjItrs

– maximum iterations for the loop to find projection on the front from ribbon.

Type

int

tmStpPrefactor

– factor for time-step adaptivity.

Type

float

maxTimeSteps

– maximum number of time steps.

Type

integer

finalTime

– time where the simulation ends.

Type

float

timeStepLimit

– limit above which time step will not exceed.

Type

float

fixedTmStp

– a float or an array giving the fixed time step. The array should have two rows, with the first row giving the time at which the time step would change and the second row giving the corresponding time step. If None is given as time step, appropriate time step would be calculated.

Type

ndarray

maxReattempts

– maximum number of reattempts in case of failure of a time step. A smaller time step will be attempted the given number of times.

Type

int

reAttemptFactor

– the factor by which time step is reduced on reattempts.

Type

float

plotFigure

– flag specifying to plot fracture trace after the given time period.

Type

boolean

saveToDisk

– flag specifying to save fracture to dist after the given time period.

Type

boolean

plotAnalytical

– if true, analytical solution will also be plotted along with the computed solution.

Type

boolean

analyticalSol

– the analytical solution of the radial fracture to be plotted on the fracture. Possible options:

  • K (toughness dominated regime, without leak off)

  • Kt (toughness dominated regime , with leak off)

  • M (viscosity dominated regime, without leak off)

  • Mt (viscosity dominated regime , with leak off)

  • E (elliptical, toughness dominated without leak off)

Type

String

bckColor

– the string specifying the parameter according to which the background of the domain is color coded. Possible options:

  • sigma0 or confining stress

  • K1c or fracture toughness

  • Cl or leak-off coefficient

Type

String

plotTimePeriod

– the time period after which the figures are plotted during simulation.

Type

float

blockFigure

– if True, the plotted figure(s) will be blocked after every time they are plotted. The simulation will advance when any key will be pressed from keyboard.

Type

bool

plotTSJump

– the number of time steps after which the variables given in plotVar attribute are plotted. E.g. a value of 4 will result in plotting every four time steps.

Type

int

plotVar

– a list of variable(s) to be plotted during simulation. The time / time steps after which the output is done can be controlled with a number of parameters ( see above).

Type

list

saveTimePeriod

– the time period after which the results are saved to disk during simulation.

Type

float

saveTSJump

– the number of time steps after which the results are saved to disk, e.g. a value of 4 will result in plotting every four time steps.

Type

int

elastohydrSolver

– the type of solver to solve the elasto-hydrodynamic system. At the moment, two main solvers can be specified.

  • ‘implicit_Picard’

  • ‘implicit_Anderson’

  • ‘RKL2’

Type

string

substitutePressure

– a flag specifying the solver to be used. If True, the pressure will be substituted in the channel elements (see Zia and Lecampion, 2019).

Type

bool

solveDeltaP

– a flag specifying the solver to be used. If True, the change in pressure, instead of pressure will be solved in the tip cells and the cells where the width constraint is active (see Zia and Lecampion, 2019).

Type

bool

solveStagnantTip

– if True, the stagnant tip cells will also be solved for width. This may result in more stable pressure as the elasticity equation will also be solved in those cells.

Type

bool

solveTipCorrRib

– if True, the tip cells corresponding to the closed ribbon cells will also be considered as closed and the width will be imposed on them.

Type

bool

solveSparse

– if True, the fluid conductivity matrix will be made with sparse matrix.

Type

bool

saveRegime

– if True, the regime of the propagation as observed in the ribbon cell (see Zia and Lecampion 2018, IJF) will be saved.

Type

boolean

verbosity

– the level of details about the ongoing simulation to be written on the log file (currently the levels ‘debug’, ‘info’, ‘warning’ and ‘error’ are supported).

Type

string

log2file

– True if you want to log to a file, otherwise set it to false

Type

bool

enableRemeshing

– if True, the computational domain will be compressed by the factor given by by the variable remeshFactor after the fracture front reaches the end of the domain.

Type

bool

remeshFactor

– the factor by which the domain is compressed on re-meshing.

Type

float

meshExtension

– an array of booleans defining if the mesh should be extended in the given direction or if it should get compressed. The distribution is bottom, top, left, right

Type

bool array

meshExtensionFactor

– factor by which the current mesh is extended in the extension direction

Type

float

meshExtendAllDir

– allow the mesh to extend in all directions

Type

bool

frontAdvancing

– The type of front advancing to be done. Possible options are:

  • ‘explicit’

  • ‘predictor-corrector’

  • ‘implicit’

Type

string

gravity

– if True, the effect of gravity will be taken into account.

Type

bool

collectPerfData

– if True, the performance data will be collected in the form of a tree.

Type

bool

paramFromTip

– if True, the space dependant parameters such as toughness and leak-off coefficients will be taken from the tip by projections instead of taking them from the ribbon cell center. The numerical scheme as a result will become unstable due to the complexities in finding the projection

Type

bool

saveReynNumb

– if True, the Reynold’s number at each edge of the cells inside the fracture will be saved.

Type

boolean

saveFluidFlux

– if True, the fluid flux at each edge of the cells inside the fracture will be saved.

Type

boolean

saveFluidVel

– if True, the fluid velocity at each edge of the cells inside the fracture will be saved.

Type

boolean

saveEffVisc

be saved.

Type

boolean

TI_KernelExecPath

– the folder containing the executable to calculate transverse isotropic kernel or kernel with free surface.

Type

string

explicitProjection

– if True, direction from last time step will be used to evaluate TI parameters.

Type

bool

symmetric

– if True, the four quadrant of the domain will be considered symmetric and only one will be solved for. The rest will be replaced by its reflection along the x and y axes.

Attention

The symmetric fracture is only implemented in the toughness dominated case. Use full domain if viscous fluid is injected.

Type

bool

enableGPU -- if True, the dense matrix vector product for the RKL scheme would be done using

GPU. If False, multithreaded dot product implemented in the explicit_RKL module will be used to do it.

nThreads -- The number of threads to be used for the dense matrix dot product in the RKL

scheme. By default set to 4.

projMethod

– the method by which the angle prescribed by the projections onto the front are evaluated. Possible options are:

  • ‘ILSA_orig’ (the method described in the original ILSA scheme).

  • ‘LS_grad’ (using gradient of the level set).

Type

string

height

– this parameters is only used in the case of height contained hydraulic fracture plots, e.g. to plot analytical solutions of PKN and KGD fractures.

Type

float

aspectRatio

– this parameters is only used in the case of elliptical hydraulic fracture plots, e.g. to plot analytical solutions of anisotropic toughness or TI elasticity.

Type

float

maxReattemptsFracAdvMore2Cells -- number of time reduction that are made if the fracture is advancing more than two cells
Type

e.g. because of an heterogeneity

Attention

These attributes below are private:

__outputAddress

– disk address of the files to be saved. If not given, a new ./Data/”tim stamp” folder will be automatically created.

Type

string

__solTimeSeries

– time series where the solution is required. The time stepping would be adjusted to get solution exactly at the given times.

Type

ndarray

__dryCrack_mechLoading

– if True, the mechanical loading solver will be used.

Type

bool

__viscousInjection

– if True, the the solver will also take the fluid viscosity into account.

Type

bool

__volumeControl

– if True, the the volume control solver will be used.

Type

bool

__simName

– the name of the simulation.

Type

string

__timeStamp

– the time at which the simulation properties was created.

Type

string

__tipAsymptote

– the tip asymptote to be used. Can be modified with the provided function.

Type

string

get_dryCrack_mechLoading()[source]
get_mesh_extension_direction()[source]
get_mesh_extension_factor()[source]
get_outputFolder()[source]
get_simulation_name()[source]
get_solTimeSeries()[source]
get_timeStamp()[source]
get_time_step_prefactor(t)[source]
get_tipAsymptote()[source]
get_viscousInjection()[source]
get_volumeControl()[source]
set_dryCrack_mechLoading(mech_loading)[source]
set_logging_to_file(address, verbosity_level=None)[source]
This function sets up the log, both to the file and to the file

Note: from any module in the code you can use the logging capabilities. You just have to:

  1. import the module

import logging

  1. create a child of the logger named ‘PyFrac’ defined in this function. Use a pertinent name as ‘Pyfrac.frontrec’

logger1 = logging.getLogger(‘PyFrac.frontrec’)

  1. use the object to send messages in the module, such as

logger1.debug(‘debug message’) logger1.info(‘info message’) logger1.warning(‘warn message’) logger1.error(‘error message’) logger1.critical(‘critical message’)

  1. IMPORTANT TO KNOW: 1-If you want to log only to the file in the abobe example you have to use: logger1 = logging.getLogger(‘PyFrac_LF.frontrec’) 2-SystemExit and KeyboardInterrupt exceptions are never swallowed by the logging package .

Parameters
  • address – string that defines the location where to save the log file

  • verbosity_level

    string that defines the level of logging concerning the file: ‘debug’ - Detailed information, typically of interest only when diagnosing problems. ‘info’ - Confirmation that things are working as expected. ‘warning’ - An indication that something unexpected happened, or indicative of some

    problem in the near future (e.g. ‘disk space low’). The software is still working as expected.

    ’error’ - Due to a more serious problem, the software has not been able to perform

    some function.

    ’critical’ - A serious error, indicating that the program itself may be unable to

    continue running.

Returns

set_mesh_extension_direction(direction)[source]

The function to set up in which directions the mesh should be extended

Parameters

direction (string) –

– direction where the mesh should be extended:

  • top (mesh extension towards positive y)

  • bottom (mesh extension towards negative y)

  • Left (mesh extension towards negative x)

  • Right (mesh extension towards positive x)

  • vertical (mesh extension up and down)

  • horizontal (mesh extension left and right)

  • all (extend the mesh in all directions)

set_mesh_extension_factor(ext_factor)[source]

The function to set up the factor deciding on the number of elements to add in the corresponding direction

Parameters

ext_factor (list or float) –

– the factor either given: - a float: all directions extend by the same amount - a list with two entries: the first gives the factor in x the second in

y direction

  • a list with four entries: the entries respectively correspond to

’left’, ‘right’, ‘bottom’, ‘top’.

set_outputFolder(output_address)[source]
set_simulation_name(simul_name)[source]
set_solTimeSeries(sol_t_srs)[source]
set_tipAsymptote(tip_asymptote)[source]

The function to set up the tip asymptote.

Parameters

tip_asymptote (string) –

– propagation regime. possible options are:

  • K (toughness dominated regime, without leak off)

  • M (viscosity dominated regime, without leak off)

  • Mt (viscosity dominated regime , with leak off)

  • U (Universal regime accommodating viscosity, toughness and leak off (see Donstov and Pierce, 2017), 0 order)

  • U1 (Universal regime accommodating viscosity, toughness and leak off (see Donstov and Pierce, 2017), delta correction)

  • MK (viscosity to toughness transition regime)

  • MDR (Maximum drag reduction asymptote, see Lecampion & Zia 2019)

  • M_MDR (Maximum drag reduction asymptote in viscosity sotrage

    regime, see Lecampion & Zia 2019)

  • HBF or HBF_aprox (Herschel-Bulkley fluid, see Bessmertnykh and Dontsov 2019; the tip volume is evaluated with a fast aproximation)

  • HBF_num_quad (Herschel-Bulkley fluid, see Bessmertnykh and Dontsov 2019; the tip volume is evaluated with numerical quadrature of the

    approximate function, which makes it very slow)

  • PLF or PLF_aprox (power law fluid, see Dontsov and Kresse 2017; the tip volume is evaluated with a fast aproximation)

  • PLF_num_quad (power law fluid, see Dontsov and Kresse 2017; the tip volume is evaluated with numerical quadrature of the

    approximate function, which makes it very slow)

= PLF_M (power law fluid in viscosity storage regime; see Desroche et al.)

set_viscousInjection(visc_injection)[source]
set_volumeControl(vol_control)[source]
properties.instrument_close(perfNode, perfNode_subItr, norm, n_elem, status, fail_cause, simulated_time)[source]
properties.instrument_start(itr_type, perfNode)[source]