elastohydrodynamic_solver module

This file is part of PyFrac.

Created by Haseeb Zia on Wed Dec 28 14:43:38 2016. Copyright (c) “ECOLE POLYTECHNIQUE FEDERALE DE LAUSANNE, Switzerland, Geo-Energy Laboratory”, 2016-2020. All rights reserved. See the LICENSE.TXT file for more details.

elastohydrodynamic_solver.Anderson(sys_fun, guess, interItr_init, sim_prop, *args, perf_node=None)[source]

Anderson solver for non linear system.

Parameters
  • sys_fun (function) – – The function giving the system A, b for the Anderson solver to solve the linear system of the form Ax=b.

  • guess (ndarray) – – The initial guess.

  • interItr_init (ndarray) – – Initial value of the variable(s) exchanged between the iterations (if any).

  • sim_prop (SimulationProperties) – – the SimulationProperties object giving simulation parameters.

  • relax (float) – – The relaxation factor.

  • args (tuple) – – arguments given to the residual and systems functions.

  • perf_node (IterationProperties) – – the IterationProperties object passed to be populated with data.

  • -- value of the recursive time steps to consider for the anderson iteration (m_Anderson) –

Returns

  • Xks[mk+1] (ndarray) – final solution at the end of the iterations.

  • data (tuple) – any data to be returned

elastohydrodynamic_solver.Elastohydrodynamic_ResidualFun(solk, system_func, interItr, *args)[source]

This function gives the residual of the solution for the system of equations formed using the given function.

elastohydrodynamic_solver.Elastohydrodynamic_ResidualFun_nd(solk, system_func, interItr, InterItr_o, indices_o, *args)[source]

This function gives the residual of the solution for the system of equations formed using the given function.

elastohydrodynamic_solver.FiniteDiff_operator_turbulent_implicit(w, pf, EltCrack, fluidProp, matProp, simProp, mesh, InCrack, vkm1, to_solve, active, to_impose)[source]

The function evaluate the finite difference matrix, i.e. the A matrix in the ElastoHydrodynamic equations ( see e.g. Dontsov and Peirce 2008). The matrix is evaluated by taking turbulence into account.

Parameters
  • w (ndarray) – – the width of the trial fracture.

  • EltCrack (ndarray) – – the list of elements inside the fracture

  • fluidProp (object) – – FluidProperties class object giving the fluid properties.

  • matProp (object) – – An instance of the MaterialProperties class giving the material properties.

  • simProp (object) – – An object of the SimulationProperties class.

  • mesh (CartesianMesh) – – the mesh.

  • InCrack (ndarray) – – an array specifying whether elements are inside the fracture or not with 1 or 0 respectively.

  • vkm1 (ndarray) – – the velocity at cell edges from the previous iteration (if necessary). Here, it is used as the starting guess for the implicit solver.

  • to_solve (ndarray) – – the channel elements to be solved.

  • active (ndarray) – – the channel elements where width constraint is active.

  • to_impose (ndarray) – – the tip elements to be imposed.

Returns

  • FinDiffOprtr (ndarray) – the finite difference matrix.

  • vk (ndarray) – the velocity evaluated for current iteration.

elastohydrodynamic_solver.Gravity_term(w, EltCrack, fluidProp, mesh, InCrack, simProp)[source]

This function returns the gravity term (G in Zia and Lecampion 2019).

Parameters
  • w (ndarray) – – the width of the trial fracture.

  • EltCrack (ndarray) – – the list of elements inside the fracture.

  • fluidProp (object) – – FluidProperties class object giving the fluid properties.

  • Mesh (CartesianMesh) – – the mesh.

  • InCrack (ndarray) – – An array specifying whether elements are inside the fracture or not with 1 or 0 respectively.

  • simProp (object) – – An object of the SimulationProperties class.

Returns

– the matrix with the gravity terms.

Return type

G (ndarray)

elastohydrodynamic_solver.Jacobian(Residual_function, sys_func, x, TypValue, interItr, *args)[source]

This function returns the Jacobian of the given function.

elastohydrodynamic_solver.MakeEquationSystem_ViscousFluid_pressure_substituted(solk, interItr, *args)[source]

This function makes the linearized system of equations to be solved by a linear system solver. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly; see description of the ILSA algorithm). The pressure in the tip cells and the cells where width constraint is active are solved separately. The pressure in the channel cells to be solved for change in width is substituted with width using the elasticity relation (see Zia and Lecampion 2019).

Parameters
  • solk (ndarray) – – the trial change in width and pressure for the current iteration of fracture front.

  • interItr (ndarray) – – the information from the last iteration.

  • args (tupple) –

    – arguments passed to the function. A tuple containing the following in order:

    • EltChannel (ndarray) – list of channel elements

    • to_solve (ndarray) – the cells where width is to be solved (channel cells).

    • to_impose (ndarray) – the cells where width is to be imposed (tip cells).

    • imposed_vel (ndarray) – the values to be imposed in the above list (tip volumes)

    • wc_to_impose (ndarray) – the values to be imposed in the cells where the width constraint is active. These can be different then the minimum width if the overall fracture width is small and it has not reached the minimum width yet.

    • frac (Fracture) – fracture from last time step to get the width and pressure.

    • fluidProp (object): – FluidProperties class object giving the fluid properties.

    • matProp (object): – an instance of the MaterialProperties class giving the material properties.

    • sim_prop (object): – An object of the SimulationProperties class.

    • dt (float) – the current time step.

    • Q (float) – fluid injection rate at the current time step.

    • C (ndarray) – the elasticity matrix.

    • InCrack (ndarray) – an array with one for all the elements in the fracture and zero for rest.

    • LeakOff (ndarray) – the leaked off fluid volume for each cell.

    • active (ndarray) – index of cells where the width constraint is active.

    • neiInCrack (ndarray) – an ndarray giving indices(in the EltCrack list) of the neighbours of all the cells in the crack.

    • edgeInCrk_lst (ndarray) – this list provides the indices of those cells in the EltCrack list whose neighbors are not outside the crack. It is used to evaluate the conductivity on edges of only these cells who are inside. It consists of four lists, one for each edge.

Returns

  • A (ndarray) – the A matrix (in the system Ax=b) to be solved by a linear system solver.

  • S (ndarray) – the b vector (in the system Ax=b) to be solved by a linear system solver.

  • interItr_kp1 (list) – the information transferred between iterations. It has three ndarrays
    • fluid velocity at edges

    • cells where width is closed

    • effective newtonian viscosity

  • indices (list) – the list containing 3 arrays giving indices of the cells where the solution is obtained for width, pressure and active width constraint cells.

elastohydrodynamic_solver.MakeEquationSystem_ViscousFluid_pressure_substituted_deltaP(solk, interItr, *args)[source]

This function makes the linearized system of equations to be solved by a linear system solver. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly; see description of the ILSA algorithm). The change is pressure in the tip cells and the cells where width constraint is active are solved separately. The pressure in the channel cells to be solved for change in width is substituted with width using the elasticity relation (see Zia and Lecamption 2019).

Parameters
  • solk (ndarray) – – the trial change in width and pressure for the current iteration of fracture front.

  • interItr (ndarray) – – the information from the last iteration.

  • args (tupple) –

    – arguments passed to the function. A tuple containing the following in order:

    • EltChannel (ndarray) – list of channel elements

    • to_solve (ndarray) – the cells where width is to be solved (channel cells).

    • to_impose (ndarray) – the cells where width is to be imposed (tip cells).

    • imposed_vel (ndarray) – the values to be imposed in the above list (tip volumes)

    • wc_to_impose (ndarray) – the values to be imposed in the cells where the width constraint is active. These can be different then the minimum width if the overall fracture width is small and it has not reached the minimum width yet.

    • frac (Fracture) – fracture from last time step to get the width and pressure.

    • fluidProp (object): – FluidProperties class object giving the fluid properties.

    • matProp (object): – an instance of the MaterialProperties class giving the material properties.

    • sim_prop (object): – An object of the SimulationProperties class.

    • dt (float) – the current time step.

    • Q (float) – fluid injection rate at the current time step.

    • C (ndarray) – the elasticity matrix.

    • InCrack (ndarray) – an array with one for all the elements in the fracture and zero for rest.

    • LeakOff (ndarray) – the leaked off fluid volume for each cell.

    • active (ndarray) – index of cells where the width constraint is active.

    • neiInCrack (ndarray) – an ndarray giving indices(in the EltCrack list) of the neighbours of all the cells in the crack.

    • edgeInCrk_lst (ndarray) – this list provides the indices of those cells in the EltCrack list whose neighbors are not outside the crack. It is used to evaluate the conductivity on edges of only these cells who are inside. It consists of four lists, one for each edge.

Returns

  • A (ndarray) – the A matrix (in the system Ax=b) to be solved by a linear system solver.

  • S (ndarray) – the b vector (in the system Ax=b) to be solved by a linear system solver.

  • interItr_kp1 (list) – the information transferred between iterations. It has three ndarrays
    • fluid velocity at edges

    • cells where width is closed

    • effective newtonian viscosity

  • indices (list) – the list containing 3 arrays giving indices of the cells where the solution is obtained for width, pressure and active width constraint cells.

elastohydrodynamic_solver.MakeEquationSystem_ViscousFluid_pressure_substituted_deltaP_sparse(solk, interItr, *args)[source]

This function makes the linearized system of equations to be solved by a linear system solver. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly; see description of the ILSA algorithm). The change is pressure in the tip cells and the cells where width constraint is active are solved separately. The pressure in the channel cells to be solved for change in width is substituted with width using the elasticity relation (see Zia and Lecamption 2019). The finite difference difference operator is saved as a sparse matrix.

Parameters
  • solk (ndarray) – – the trial change in width and pressure for the current iteration of fracture front.

  • interItr (ndarray) – – the information from the last iteration.

  • args (tupple) –

    – arguments passed to the function. A tuple containing the following in order:

    • EltChannel (ndarray) – list of channel elements

    • to_solve (ndarray) – the cells where width is to be solved (channel cells).

    • to_impose (ndarray) – the cells where width is to be imposed (tip cells).

    • imposed_vel (ndarray) – the values to be imposed in the above list (tip volumes)

    • wc_to_impose (ndarray) – the values to be imposed in the cells where the width constraint is active. These can be different then the minimum width if the overall fracture width is small and it has not reached the minimum width yet.

    • frac (Fracture) – fracture from last time step to get the width and pressure.

    • fluidProp (object): – FluidProperties class object giving the fluid properties.

    • matProp (object): – an instance of the MaterialProperties class giving the material properties.

    • sim_prop (object): – An object of the SimulationProperties class.

    • dt (float) – the current time step.

    • Q (float) – fluid injection rate at the current time step.

    • C (ndarray) – the elasticity matrix.

    • InCrack (ndarray) – an array with one for all the elements in the fracture and zero for rest.

    • LeakOff (ndarray) – the leaked off fluid volume for each cell.

    • active (ndarray) – index of cells where the width constraint is active.

    • neiInCrack (ndarray) – an ndarray giving indices(in the EltCrack list) of the neighbours of all the cells in the crack.

    • edgeInCrk_lst (ndarray) – this list provides the indices of those cells in the EltCrack list whose neighbors are not outside the crack. It is used to evaluate the conductivity on edges of only these cells who are inside. It consists of four lists, one for each edge.

Returns

  • A (ndarray) – the A matrix (in the system Ax=b) to be solved by a linear system solver.

  • S (ndarray) – the b vector (in the system Ax=b) to be solved by a linear system solver.

  • interItr_kp1 (list) – the information transferred between iterations. It has three ndarrays
    • fluid velocity at edges

    • cells where width is closed

    • effective newtonian viscosity

  • indices (list) – the list containing 3 arrays giving indices of the cells where the solution is obtained for width, pressure and active width constraint cells.

elastohydrodynamic_solver.MakeEquationSystem_ViscousFluid_pressure_substituted_sparse(solk, interItr, *args)[source]

This function makes the linearized system of equations to be solved by a linear system solver. The finite difference difference opertator is saved as a sparse matrix. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly; see description of the ILSA algorithm). The pressure in the tip cells and the cells where width constraint is active are solved separately. The pressure in the channel cells to be solved for change in width is substituted with width using the elasticity relation (see Zia and Lecamption 2019).

Parameters
  • solk (ndarray) – – the trial change in width and pressure for the current iteration of fracture front.

  • interItr (ndarray) – – the information from the last iteration.

  • args (tupple) –

    – arguments passed to the function. A tuple containing the following in order:

    • EltChannel (ndarray) – list of channel elements

    • to_solve (ndarray) – the cells where width is to be solved (channel cells).

    • to_impose (ndarray) – the cells where width is to be imposed (tip cells).

    • imposed_vel (ndarray) – the values to be imposed in the above list (tip volumes)

    • wc_to_impose (ndarray) – the values to be imposed in the cells where the width constraint is active. These can be different then the minimum width if the overall fracture width is small and it has not reached the minimum width yet.

    • frac (Fracture) – fracture from last time step to get the width and pressure.

    • fluidProp (object): – FluidProperties class object giving the fluid properties.

    • matProp (object): – an instance of the MaterialProperties class giving the material properties.

    • sim_prop (object): – An object of the SimulationProperties class.

    • dt (float) – the current time step.

    • Q (float) – fluid injection rate at the current time step.

    • C (ndarray) – the elasticity matrix.

    • InCrack (ndarray) – an array with one for all the elements in the fracture and zero for rest.

    • LeakOff (ndarray) – the leaked off fluid volume for each cell.

    • active (ndarray) – index of cells where the width constraint is active.

    • neiInCrack (ndarray) – an ndarray giving indices(in the EltCrack list) of the neighbours of all the cells in the crack.

    • edgeInCrk_lst (ndarray) – this list provides the indices of those cells in the EltCrack list whose neighbors are not outside the crack. It is used to evaluate the conductivity on edges of only these cells who are inside. It consists of four lists, one for each edge.

Returns

  • A (ndarray) – the A matrix (in the system Ax=b) to be solved by a linear system solver.

  • S (ndarray) – the b vector (in the system Ax=b) to be solved by a linear system solver.

  • interItr_kp1 (list) – the information transferred between iterations. It has three ndarrays
    • fluid velocity at edges

    • cells where width is closed

    • effective newtonian viscosity

  • indices (list) – the list containing 3 arrays giving indices of the cells where the solution is obtained for width, pressure and active width constraint cells.

elastohydrodynamic_solver.MakeEquationSystem_mechLoading(wTip, EltChannel, EltTip, C, EltLoaded, w_loaded)[source]

This function makes the linear system of equations to be solved by a linear system solver. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly). The given width is imposed on the given loaded elements.

elastohydrodynamic_solver.MakeEquationSystem_volumeControl(w_lst_tmstp, wTip, EltChannel, EltTip, sigma_o, C, dt, Q, ElemArea, lkOff)[source]

This function makes the linear system of equations to be solved by a linear system solver. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly). The the volume of the fracture is imposed to be equal to the fluid injected into the fracture.

elastohydrodynamic_solver.MakeEquationSystem_volumeControl_double_fracture(w_lst_tmstp, wTipFR0, wTipFR1, EltChannel0, EltChannel1, EltTip0, EltTip1, sigma_o, C, dt, QFR0, QFR1, ElemArea, lkOff)[source]

This function makes the linear system of equations to be solved by a linear system solver. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly). The the volume of the fracture is imposed to be equal to the fluid injected into the fracture.

elastohydrodynamic_solver.MakeEquationSystem_volumeControl_symmetric(w_lst_tmstp, wTip_sym, EltChannel_sym, EltTip_sym, C_s, dt, Q, sigma_o, ElemArea, LkOff, vol_weights, sym_elements, dwTip)[source]

This function makes the linear system of equations to be solved by a linear system solver. The system is assembled with the extended footprint (treating the channel and the extended tip elements distinctly). The the volume of the fracture is imposed to be equal to the fluid injected into the fracture (see Zia and Lecampion 2018).

elastohydrodynamic_solver.Picard_Newton(Res_fun, sys_fun, guess, TypValue, interItr_init, sim_prop, *args, PicardPerNewton=1000, perf_node=None)[source]

Mixed Picard Newton solver for nonlinear systems.

Parameters
  • Res_fun (function) – – The function calculating the residual.

  • sys_fun (function) – – The function giving the system A, b for the Picard solver to solve the linear system of the form Ax=b.

  • guess (ndarray) – – The initial guess.

  • TypValue (ndarray) – – Typical value of the variable to estimate the Epsilon to calculate Jacobian.

  • interItr_init (ndarray) – – Initial value of the variable(s) exchanged between the iterations (if any).

  • sim_prop (SimulationProperties) – – the SimulationProperties object giving simulation parameters.

  • relax (float) – – The relaxation factor.

  • args (tuple) – – arguments given to the residual and systems functions.

  • PicardPerNewton (int) – – For hybrid Picard/Newton solution. Number of picard iterations for every Newton iteration.

  • perf_node (IterationProperties) – – the IterationProperties object passed to be populated with data.

Returns

  • solk (ndarray) – solution at the end of iteration.

  • data (tuple) – any data to be returned

elastohydrodynamic_solver.Velocity_Residual(v, *args)[source]

This function gives the residual of the velocity equation. It is used by the root finder.

Parameters
  • v (float) – – current velocity guess

  • args (tuple) –

    – a tuple consisting of the following:

    • w (float) width at the given cell edge

    • mu (float) viscosity at the given cell edge

    • rho (float) density of the injected fluid

    • dp (float) pressure gradient at the given cell edge

    • rough (float) roughness (width / grain size) at the cell center

Returns

– residual of the velocity equation

Return type

float

elastohydrodynamic_solver.calculate_fluid_flow_characteristics_laminar(w, pf, sigma0, Mesh, EltCrack, InCrack, muPrime, density)[source]

This function calculate fluid flux and velocity at the cell edges evaluated with the pressure calculated from the elasticity relation for the given fracture width and the poisoille’s Law.

elastohydrodynamic_solver.check_covergance(solk, solkm1, indices, tol)[source]

This function checks for convergence of the solution

Parameters
  • solk (ndarray) –

  • solkm1 (ndarray) –

  • indices (list) – for channel, tip and active width constraint cells.

  • tol (float) –

Returns

  • converged (bool) – True if converged

  • norm (float) – the evaluated norm which is checked against tolerance

elastohydrodynamic_solver.findBracket(func, guess, *args)[source]

This function can be used to find bracket for a root finding algorithm.

Parameters
  • func (callable function) – – the function giving the residual for which zero is to be found

  • guess (float) – – starting guess

  • args (tupple) – – arguments passed to the function

Returns

  • a (float) – the lower bracket

  • b (float) – the higher bracket

elastohydrodynamic_solver.finiteDiff_operator_Herschel_Bulkley(w, pf, EltCrack, fluidProp, Mesh, InCrack, neiInCrack, edgeInCrk_lst, simProp)[source]

The function evaluate the finite difference 5 point stencil matrix, i.e. the A matrix in the ElastoHydrodynamic equations in e.g. Dontsov and Peirce 2008. The matrix is evaluated for Herschel-Bulkley fluid rheology.

Parameters
  • w (ndarray) – – the width of the trial fracture.

  • pf (ndarray) – – the fluid pressure.

  • EltCrack (ndarray) – – the list of elements inside the fracture.

  • fluidProp (object) – – FluidProperties class object giving the fluid properties.

  • Mesh (CartesianMesh) – – the mesh.

  • InCrack (ndarray) – – an array specifying whether elements are inside the fracture or not with 1 or 0 respectively.

  • neiInCrack (ndarray) – – an ndarray giving indices of the neighbours of all the cells in the crack, in the EltCrack list.

  • edgeInCrk_lst (ndarray) – – this list provides the indices of those cells in the EltCrack list whose neighbors are not outside the crack. It is used to evaluate the conductivity on edges of only these cells who are inside. It consists of four lists, one for each edge.

  • simProp (object) – – An object of the SimulationProperties class.

Returns

– the finite difference matrix.

Return type

FinDiffOprtr (ndarray)

elastohydrodynamic_solver.finiteDiff_operator_laminar(w, EltCrack, muPrime, Mesh, InCrack, neiInCrack, simProp)[source]

The function evaluate the finite difference 5 point stencil matrix, i.e. the A matrix in the ElastoHydrodynamic equations in e.g. Dontsov and Peirce 2008. The matrix is evaluated with the laminar flow assumption.

Parameters
  • w (ndarray) – – the width of the trial fracture.

  • EltCrack (ndarray) – – the list of elements inside the fracture.

  • muPrime (ndarray) – – the scaled local viscosity of the injected fluid (12 * viscosity).

  • Mesh (CartesianMesh) – – the mesh.

  • InCrack (ndarray) – – an array specifying whether elements are inside the fracture or not with 1 or 0 respectively.

  • neiInCrack (ndarray) – – an ndarray giving indices of the neighbours of all the cells in the crack, in the EltCrack list.

  • simProp (object) – – An object of the SimulationProperties class.

Returns

– the finite difference matrix.

Return type

FinDiffOprtr (ndarray)

elastohydrodynamic_solver.finiteDiff_operator_power_law(w, pf, EltCrack, fluidProp, Mesh, InCrack, neiInCrack, edgeInCrk_lst, simProp)[source]

The function evaluate the finite difference 5 point stencil matrix, i.e. the A matrix in the ElastoHydrodynamic equations in e.g. Dontsov and Peirce 2008. The matrix is evaluated for Herschel-Bulkley fluid rheology.

Parameters
  • w (ndarray) – – the width of the trial fracture.

  • pf (ndarray) – – the fluid pressure.

  • EltCrack (ndarray) – – the list of elements inside the fracture.

  • fluidProp (object) – – FluidProperties class object giving the fluid properties.

  • Mesh (CartesianMesh) – – the mesh.

  • InCrack (ndarray) – – an array specifying whether elements are inside the fracture or not with 1 or 0 respectively.

  • neiInCrack (ndarray) – – an ndarray giving indices of the neighbours of all the cells in the crack, in the EltCrack list.

  • edgeInCrk_lst (ndarray) – – this list provides the indices of those cells in the EltCrack list whose neighbors are not outside the crack. It is used to evaluate the conductivity on edges of only these cells who are inside. It consists of four lists, one for each edge.

  • simProp (object) – – An object of the SimulationProperties class.

Returns

– the finite difference matrix.

Return type

FinDiffOprtr (ndarray)

elastohydrodynamic_solver.get_complete_solution(sol, indices, *args)[source]
elastohydrodynamic_solver.get_finite_difference_matrix(wNplusOne, sol, frac_n, EltCrack, neiInCrack, fluid_prop, mat_prop, sim_prop, mesh, InCrack, C, interItr, to_solve, to_impose, active, interItr_kp1, list_edgeInCrack)[source]
elastohydrodynamic_solver.populate_full(indices, sol, values=None)[source]
elastohydrodynamic_solver.pressure_gradient(w, C, sigma0, Mesh, EltCrack, InCrack)[source]

This function gives the pressure gradient at the cell edges evaluated with the pressure calculated from the elasticity relation for the given fracture width.

elastohydrodynamic_solver.pressure_gradient_form_pressure(pf, Mesh, EltCrack, InCrack)[source]

This function gives the pressure gradient at the cell edges evaluated with the pressure

elastohydrodynamic_solver.velocity(w, EltCrack, Mesh, InCrack, muPrime, C, sigma0)[source]

This function gives the velocity at the cell edges evaluated using the Poiseuille flow assumption.