Redistancing¶
The redistancing algorithms to apply on a Level Set after a couple of time steps are implemented in _Solver/Redistancers
as a module.
The redistancing technique has to be implemented within a class, and the main method to redistance (given within the class)
should be named Redistance
in any case.
The various spatial application routine of the redistancing as well as the parsing of user wishes for the resitancing is performed in the wrapper class contained in the
_Solver/Redistancers/WrapRedistancer.py
file.
Page navigation
Redistancing choice and application subdomain selection¶
The module Redistancer
located in the file _Solver/SpatialModifiers/Redistancer.py`furnishes a single class
:code:`Redistancing
that interfaces the user instructions for redistancing the level set given a wished scheme and its parameters.
Upon instance creation, it converts the string-user input defining the redistancing method to use, its application range (everywhere, if __name__ == ‘__main__’:
an interface narrow band, etc) and defines an interface application routine that is to apply on any freshly computed Level set-residual.
The method Redistance
is therefore the generic interface to apply any LS residual to redistance.
Usage
Upon instance creation, the redistanciation routine is selected and constructed upon the user parameters that are specified in the following format
# Defining the parameters as RedistancingRoutineIndex_(ApplicationRangeIndex,ApplicationPeriod,ApplicationRangeParameter1,ApplicationRangeParameter2,...)<Parameters,of,selected,routine>
Parameters = 1_(1,1,0.5)<0.2,1e-8,1e-5>
# If no parameters is to be given or if you want to use the default ones (if specified in the selected RedistancingRoutine), then please use the syntax
Parameters = 1_(1,1,0.5)<>
You can then define the redistanciation routine by creating the wrapper instance;
Redistancer = Redistancing(Parameters, Problem, Mesh)
and finally call the redistanciation routine by using the field Redistance
:, by
# Within the lsrd code, always call with the "ConstructionPoints" argument, being here the Dofs where Phi is available
Phi = Solution.LSValues[0,:]
ConstructionPoints = Mesh.DofsXY[:,:]
NewValues = Redistancer.Redistance(Points, Phi, ConstructionPoints=ConstructionPoints)
Note
If the parameters of the selected redistancing routine requires a quadrature rule, please use the same quadrature rule format as for the properties of the spatial modifiers.
Note
The use of this class is triple: parse the user parameters to select the wished redistancer and its application range (spatially), allowing the use of a custom quadrature scheme in the redistancing function (if needed), and to selecting the points to redistance at.
-
class
3_LESARD.SourceCode._Solver.Redistancers.WrapRedistancer.
Redistancing
(Parameters, Problem, Mesh)¶ Class furnishing all the tools for a generic redistancing routine and registers the user’s wished redistancing parameters, both for its application spatial range andfor the method itself. The main iteration is accessible through the method
Redistance
, and is mapped in the methodSelectRedistancing
(see its description to known the parameter value to give for selecting the desired method).- Parameters
RedistancingParameters (string) – the redistancing parameters, containing the wished redistancing method, its spatial application range, frequency and the method’s parameters.
Problem (Problem) – the considered problem
Mesh (MeshStructure) – the considered mesh instance
Note
See
ParseRedistancingParams
for a description of the expected format.Methods
-
ParseRedistancingParams
()¶ Converting the string user input into sets of parameters: the index of the wished redistanciation routine, the spatial application properties and the properties of the redistancing method itself.
- Parameters
None –
the self.Params field of the class should be filled and containing a string of the format
RedistancingRoutineIndex_(ApplicationRangeIndex, ApplicationPeriod, ApplicationRangeParameter1,ApplicationRangeParameter2,...)<Parameters,of,selected,routine>
The minimal number of parameter to give is
RedistancingRoutineIndex_(ApplicationRangeIndex,ApplicationPeriod)<>
or0
if no redistancing is wished.
Note
See the module
Mappers
to get the Index corresponding to the desired redistancing routine methodSee the method
SelectRedistancing
to get the ApplicationRangeIndexRefer to the mapped routine’s documentation (
NarrowBandRedistancing
,GlobalRedistancing
, etc) for the knowing which ApplicationRangeParameters to useSee the selected RedistancingRoutine’s documentation for knowing its parameters
Warning
There should be NO space in the definition of the parameter string, even between the parameters specification
- Returns
Tuple contining in order:
the index of the wished redistancing routine (integer)
the spatial application range ([SpatialType (integer), ApplicationPeriod(integer), Paramter1 (string), Parameter2 (string), …])
the parameters of the redistancing routine (list containing strings and/or a quadrature scheme instance)
- Return type
Properties (tuple)
-
SelectRedistancing
(RedistanciationRange, *RangeParameters)¶ Mapping that maps to the wished redistancing method depending on the user wishes upon to the following indices to be given as a parameter.
- Parameters
RedistanciationRange (integer) –
The index corresponding to the wished spatial application range of the redistanciation routine, as follows
NoRedistancing
– Returns an identity functionNarrowBandRedistancing
– Performs a redistancing within a narrow band around the interfaceGlobalRedistancing
– Performs a redistancing over the full computational domainNarrowBandFixedRedistancing
– Performs a redistancing within a narrow band around the interface, keeping the level set values at the immediate vicinity of the interface untouched
RangeParameters (Optional) – If the selected redistanciation application metod requires arguments, its arguments (see the documentation of the selected routine to know them).
-
NoRedistancing
(Solution, *args)¶ Dummy routine returning the identity over the LevelSet values of the solution, only defined for compatibility isses of the general code.
- Parameters
Solution (Solution) – The solution structure of the investigated problem
- Returns
The level set values numpy array contained within the solution structures, without any change
- Return type
RedistancedSol (numpy array)
Note
*args is here for compatibility-on-call reasons
-
NarrowBandRedistancing
(Solution, BandWidth, *args)¶ Routines that applies the desired redistancing method over a narrow band around the interface.
- Parameters
Solution (Solution) – The solution structure of the investigated problem
BandWidth (single float-containing string) – The width of the band to consider for redistancing (= maximal value of the level set under which it will be redistanced)
- Returns
- The redistanced level set values to substitute to LSValues in the solution structure
(the substitution is not performed automatically)
- Return type
RedistancedSol (numpy array)
Note
*args is here for compatibility-on-call reasons
-
GlobalRedistancing
(Solution, *args)¶ Routines that applies the desired redistancing method over the full computational domain.
- Parameters
Solution (Solution) – The solution structure of the investigated problem
- Returns
- The redistanced level set values to substitute to LSValues in the solution structure
(the substitution is not performed automatically)
- Return type
RedistancedSol (numpy array)
Warning
Depending on the selected Redistanciation method, it may fail for points that are close to the boundary
Note
*args is here for compatibility-on-call reasons
Available redistancers and their parameter specifications¶
HopfLax Descent¶
Module that furnishes the tools to perform a redistancing through a Hopf-lax formula, solving the optimisation problem via a gradient-method approach. Following the work of Micheal Royston and Byungjoon Lee, “Parallel redistancing using the Hopf-Lax formula, 10.1016/j.jcp.2018.01.035”, but with a one-sided bisection method instead of the proposed secant method to allow saddle points safely
Note
Due to the used root minding method, the convergence is not strictly monotoneous in this implementation.
Due to the need of the gradient, when a bounded domain is given some points near the boundary may not have converged and are replaced by the initial value of Phi
The local extrema are avoided very dirtily and the convergence to the relevant local minima is not guaranteed
This toolbox is only valid for scalar-valued functions (representing a level set)
Warning
The code is not particularily robust
This method requires a gradient information, which if no gradient is furnished is done through a WENo scheme-like approximation.
- It is not robust for highly non-convex problem (it will be sensitive to the first guess in the Hopf-Lax solver via gradient descent).
Here the choice had been to take the initial guess at the evaluated point. However, for distorded initial level set functions this may yield to local-spurious values if one has been trapped in a wrond local minima in the optimisation process. For very distorded initial data, please use rather the reinitialisation HopfLax_GlobalBregman.
- When in the minimisation problem we arrive to a local maximum, the convergence may be slow due to the extra work to fall outside
the trap. From which side of the extrema we fall is not constrained and may yield to wrong local minima convergence (see the above note).
On coarse grid if the gradient has to be approximated we may be stuck in local maxima and redistancing may fail (esp. if the gradient make jumps to two local extrema)
-
class
3_LESARD.SourceCode._Solver.Redistancers.HopfLax_Descent.
Redistancing
(InitTime=0.2, eps=1e-08, gtol=1e-05, stol=1e-08, AvoidExtrema=0.01, smaxit=1000, gmaxit=50, verbose=False)¶ Class furnishing the routines to define a redistanciation routine based on a HopfLax formula.
Note
The routine is updating the values of the level set at the degrees of freedom of the given mesh from the values that are stored there. Therefore, it is only suitable to be used on solution computed from a point-values dofs and not mixed techniques as the reconstruction techniques back and forth from generic dofs to point-values are not considered.
- Parameters
InitTime (float) – Optional, default value = 0.2. Set the pseudo-time first iteration step in the Bregman split iteration process (w.r.t. hmin: t0=0, t1 = BregmanInitTime*hmax)
eps (float) – Optional, default value = 1e-8. Safety net in the division of the weno weight and gradient division
gtol (float) – Optional, default value = 1e-5. The tolerance to stop the gradient descent iteration
stol (float) – Optional, default value = 1e-8. The tolerance to stop the and the secant iteration
AvoidExtrema (float) – Optional, default value = 1e-2. Coefficient used in the local extrama get away (AvoidExtrema<1). If the curvature is very low, reduce the gtol. Has an impact on the desired gradient descent accuracy (refines it), usually not useful when no secant method is used.
smaxit (integer) – Optional, default value = 1000. Maximum number of secant method iteration
gmaxit (integer) – Optional, default value = 50. Maximum number of gradient descent iteration
verbose (boolean) – Optional, default value = False. True if the user wants a printout each point of the obtained values or False if no wished printout
Note
The arguments can be given as strings and will be converted accordingly if the string content is compatible with the exepected type.
-
W
(a, b, c, d, eps)¶ Function defining the weno weigths.
- Parameters
a (float) – The parameters to pass in order to compute the weno weights
b (float) – The parameters to pass in order to compute the weno weights
c (float) – The parameters to pass in order to compute the weno weights
d (float) – The parameters to pass in order to compute the weno weights
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
The Weno weight associated to the given coefficients
- Return type
Weight (float)
-
Deltas
(Phi, Points, eps)¶ Simple helper function that compute the differences involved in the Weno’s derivative approximation
- Parameters
Phi (Function callback) – An interpolation of the function to get the derivatives for
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns a tuple containing in order:
DeltaxMPhi (NbPoints float): Backward difference in x
DeltayMPhi (NbPoints float): Backward difference in y
DeltaxPPhi (NbPoints float): Forward difference in x
DeltayPPhi (NbPoints float): Forward difference in y
DDeltaxMPPhi (NbPoints float): Backward Forward difference in x
DDeltayMPPhi (NbPoints float): Backward Forward difference in y
- Return type
Differences (numpy arrays tuple)
-
GetGradient
(Phi, Points, eps)¶ Numerical gradient computation according to the WENO derivative computations
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback) – An interpolation of the function to get the derivatives for
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns np.array([Gradx, Grady]) where Gradx and Grady are as follows.
Gradx (NbPoints float): gradient value in the x direction at all the given points
Grady (NbPoints float): gradient value in the y direction at all the given points
- Return type
Grad (2d numpy array)
-
InterpolatePhi
(Points, Phi, *args)¶ Interpolator that considers either a function Phi or a set of Points and point-values Phi and returns an interpolator function. It also returns the Gradient interpolator of the function contructed from a gradient approximation at the given Points.
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback or NbPoints float) – Level set to reinitialise, either given as a function or a list of point values associated to points
args (Two function callback, optional) – The callback function defining the partial derivatives in x and y, in that order
- Returns
The interpolator of the function Phi InterpPhiDx (function callback): The interpolator of the derivative in x of Phi InterpPhiDy (function callback): The interpolator of the derivative in y of Phi
- Return type
InterpPhi (function callback)
Note
A fill value is set to max(abs(Phi))+0.1 for the function itself, but no fill value is set for the gradient
If function callbacks are given, they should be 2d->IR, taking the xy-values as a column list (NbPoints x 2) and return a 1d numpy vector
-
GetNextPhi
(Point, v, t, Interp, InterpDx, InterpDy, hmin)¶ Optimisation problem solver for the core Hopf-Lax formula, to be used within a zero finding method, using here the gradient descent method.
- Parameters
Point (float numpy array) – Point [x,y] at which we want to update the Level set value
v (float numpy array) – The point [x,y] to start the gradient descent algorithm at
t (float) – Distance from the Point point from which the value of the optimiser should be projected on the ball. It should be non-negative.
Interp (function callback) – The interpolator (or function itself) representing conitnously the level-set in the 2D plane (can be None or Nan outside the domain if there is no data there)
InterpDx (function callback) – The interpolator (or function itself) representing conitnously the derivative in x of level-set in the 2D plane (can be None or Nan outside the domain if there is no data there)
InterpDy (function callback) – The interpolator (or function itself) representing conitnously the derivative in y of level-set in the 2D plane (can be None or Nan outside the domain if there is no data there)
hmin (float) – The mininum distance between two points in the considered discretisation
- Returns
Value of the level set at the found optimizer v (float numpy array): The array [x,y] of the minimizer coordinates gfailed (float): Returns a failure flag: 0 if we could not converge to the wished tolerance within the given amount of steps, 1 otherwise
- Return type
Phi1 (float)
-
Redistance
(Points, Phi, *args, **kwargs)¶ Implementation of the main iteration strategy of the Bregman split for the redistancing problem.
- Parameters
Points (float 2D numpy array) – The points on which reinitialise the level-set on (NbPoints x 2)
Phi (function callback) – The (initial) level set function to redistance (should return a scalar value), given either as a float list associated to the points or as a function callback
args (Two function callback) – (Optional) The callback function defining the partial derivatives in x and y, in that order
kwargs (keyword argument) – The only keyword argument available is “ConstructionPoints”, only active when Phi is given as a float list. It corresponds to where the given point values Phi have been computed. If the argument is not passed when Phi is a float list, the function assumes that Phi has been constructed on Points. Note that the Points should be contained within the convex hull of ConstructionPoints.
- Returns
The value of the redistanced level set function at the dofs (NbDofs)
- Return type
RPhi (float numpy array)
HopfLax Bregman¶
Module that furnishes the tools to perform a redistancing through a Hopf-lax formula, solving the optimisation problem via a Bregman split approach. Following the work of Byungjoon Lee, “Revisiting the redistaning problem using the Hopf-Lax formula, 10.1016/j.jcp.2016.11.005” opf-Lax formula, 10.1016/j.jcp.2016.11.005”, with a one-sided bisection method instead of the proposed secant method to allow saddle points safely.
Note
The cheap way it is coded makes it slow
Requires a first coarse rough initial step, by default, this is done by a WENO procedure
Due to the root finding methd, the convergence is not strictly monotoneous in this implementation.
The local extrema are avoided only from a convexifying ratio, which slows down the convergence process for already convex functions
This toolbox is only valid for scalar-valued functions (representing a level set)
The optimisation subproblem is solved using the scipy library
Warning
The code is not particularily robust and may fail on saddle points if a too low convexifying ratio is set
When using a high convexity ratio, the implementation becomes considerabily slower and may result in unwanted results if the InitBregmanTime is too large.
-
class
3_LESARD.SourceCode._Solver.Redistancers.HopfLax_Bregman.
Redistancing
(InitTime=0.2, CoarseningRatio=1, ConvexifyingRatio=1, stol=1e-08, gtol=1e-05, smaxit=1000, gmaxit=50, eps=1e-06, cfl=0.2, verbose=False)¶ Class furnishing the routines to define a redistanciation routine based on a HopfLax formula.
Note
The routine is updating the values of the level set at the degrees of freedom of the given mesh from the values that are stored there. Therefore, it is only suitable to be used on solution computed from a point-values dofs and not mixed techniques as the reconstruction techniques back and forth from generic dofs to point-values are not considered.
- Parameters
InitTime (float) – Optional, default value = 0.2. Set the pseudo-time first iteration step in the Bregman split iteration process (w.r.t. hmin: t0=0, t1 = BregmanInitTime*hmax)
CoarseningRatio (float) – Optional, default value = 1. How much coarse the cartesian grid for the initialisation should be (from the unstructured hmax)
ConvexifyingRatio (float) – Optional, default value = 1. The convexification ratio to be used in the Bregman split iteration process
stol (float) – Optional, default value = 1e-8. The tolerance to stop the Secant iteration
gtol (float) – Optional, default value = 1e-5. The tolerance to stop the Bregman iteration
smaxit (float) – Optional, default value = 1000. The maximum iteration step for the secant method
gmaxit (float) – Optional, default value = 50. The maximum iteration step for the bregman method
eps (float) – Optional, default value = 1e-6. Safety net in the division of the weno weight and gradient division
cfl (float) – Optional, default value = 0.2. CFL for the intialisation of Bregman split
verbose (boolean) – Optional, default value = False. True if the user wants a printout each point of the obtained values or False if no wished printout
-
W
(a, b, c, d, eps)¶ Function defining the weno weigths.
- Parameters
a (float) – The parameters to pass in order to compute the weno weights
b (float) – The parameters to pass in order to compute the weno weights
c (float) – The parameters to pass in order to compute the weno weights
d (float) – The parameters to pass in order to compute the weno weights
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
The Weno weight associated to the given coefficients
- Return type
Weight (float)
-
Deltas
(Phi, Points, eps)¶ Simple helper function that compute the differences involved in the Weno’s derivative approximation
- Parameters
Phi (Function callback) – An interpolation of the function to get the derivatives for
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns a tuple containing in order:
DeltaxMPhi (NbPoints float): Backward difference in x
DeltayMPhi (NbPoints float): Backward difference in y
DeltaxPPhi (NbPoints float): Forward difference in x
DeltayPPhi (NbPoints float): Forward difference in y
DDeltaxMPPhi (NbPoints float): Backward Forward difference in x
DDeltayMPPhi (NbPoints float): Backward Forward difference in y
- Return type
Differences (numpy arrays tuple)
-
GetGradient
(Phi, Points, eps)¶ Numerical gradient computation according to the WENO derivative computations
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback) – An interpolation of the function to get the derivatives for
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns np.array([Gradx, Grady]) where Gradx and Grady are as follows.
Gradx (NbPoints float): gradient value in the x direction at all the given points
Grady (NbPoints float): gradient value in the y direction at all the given points
- Return type
Grad (2d numpy array)
-
InterpolatePhi
(Points, Phi, *args)¶ Interpolator that considers either a function Phi or a set of Points and point-values Phi and returns an interpolator function. It also returns the Gradient interpolator of the function contructed from a gradient approximation at the given Points.
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback or NbPoints float) – Level set to reinitialise, either given as a function or a list of point values associated to points
args (Two function callback, optional) – The callback function defining the partial derivatives in x and y, in that order
- Returns
The interpolator of the function Phi InterpPhiDx (function callback): The interpolator of the derivative in x of Phi InterpPhiDy (function callback): The interpolator of the derivative in y of Phi
- Return type
InterpPhi (function callback)
Note
If function callbacks are given, they should be 2d->IR, taking the xy-values as a column list (NbPoints x 2) and return a 1d numpy vector
-
RoughEstimationAssetEquation
(Points, Phi)¶ Getting a first rough approximation of the evaluation solution of Phi of the redistancing equation through a fifth order weno in space and euler in time scheme.
- Parameters
Points (Function callback) – The list of points where to consider the function
Phi (Function callback) – Level set to reinitialise, as a function or an interpolator
- Returns
The interpolator of the function Phi at the InitBregmanTime InterpPhiDx (function callback): The interpolator of the derivative in x of Phi at the InitBregmanTime InterpPhiDy (function callback): The interpolator of the derivative in y of Phi at the InitBregmanTime
- Return type
InterpPhi (function callback)
-
GetNextPhi
(Point, v, d, b, t, Interp)¶ Optimisation problem solver for the core Hopf-Lax formula, to be used within a zero finding method, using here the Bregman Split algorithm.
- Parameters
Point (float numpy array) – Point [x,y] at which we want to update the Level set value
v (float numpy array) – The initial value [x,y] of v to consider in the Bregman split
d (float numpy array) – The initial value [x,y] of d to consider in the Bregman split
b (float numpy array) – The initial value [x,y] of b to consider in the Bregman split
t (float) – Distance from the Point point from which the value of the optimiser should be projected on the ball. It should be non-negative.
Interp (function callback) – The interpolator (or function itself) representing conitnously the level-set in the 2D plane (can be None or Nan outside the domain if there is no data there)
- Returns
Value of the level set at the found optimizer d (float numpy array): The array [x,y] of the minimizer coordinates gfailed (float): Returns a failure flag: 0 if we could not converge to the wished tolerance within the given amount of steps, 1 otherwise
- Return type
Phi1 (float)
-
Redistance
(Points, Phi, *args, **kwargs)¶ Implementation of the main iteration strategy of the Bregman split for the redistancing problem.
- Parameters
Points (float 2D numpy array) – The points on which reinitialise the level-set on (NbPoints x 2)
Phi (function callback) – The (initial) level set function to redistance (should return a scalar value), given either as a float list associated to the points or as a function callback
args (Two function callback) – (Optional) The callback function defining the partial derivatives in x and y, in that order
kwargs (keyword argument) – The only keyword argument available is “ConstructionPoints”, only active when Phi is given as a float list. It corresponds to where the given point values Phi have been computed. If the argument is not passed when Phi is a float list, the function assumes that Phi has been constructed on Points. Note that the Points should be contained within the convex hull of ConstructionPoints.
- Returns
The value of the redistanced level set function at the dofs (NbDofs)
- Return type
RPhi (float numpy array)
HopfLax Global Bregman¶
Module that furnishes the tools to perform a redistancing through a Hopf-lax formula, solving the optimisation problem via a modified Bregman split approach that is supposed to deal with local minima (not always in practice though).
Following the work of Byungjoon Lee, “Revisiting the redistaning problem using the Hopf-Lax formula, 10.1016/j.jcp.2016.11.005” opf-Lax formula, 10.1016/j.jcp.2016.11.005”, with a one-sided bisection method instead of the proposed secant method to allow saddle points safely, and using an application of the work “Beyond Alternating Updates for Matrix Factorization with Inertial Bregman Proximal Gradient Algorithms” from Mahesh Chandra Mukkamala and Peter Ochs for the Bregman split instead.
Note
Due to the modification of the secant formula, the convergence is not strictly monotoneous in this implementation.
The local extrema are avoided but the saddle points stay problematic when trying to achieve global minima convergence.
This toolbox is only valid for scalar-valued functions (representing a level set)
Warning
The code is not particularily robust and may fail on saddle points
On coarse grid if the gradient has to be approximated we may be stuck in local maxima and redistancing may fail
-
class
3_LESARD.SourceCode._Solver.Redistancers.HopfLax_GlobalBregman.
Redistancing
(InitTime=0.2, stol=1e-08, gtol=1e-05, smaxit=1000, gmaxit=50, eps=1e-08, verbose=False)¶ Class furnishing the routines to define a redistanciation routine based on a HopfLax formula.
Note
The routine is updating the values of the level set at the degrees of freedom of the given mesh from the values that are stored there. Therefore, it is only suitable to be used on solution computed from a point-values dofs and not mixed techniques as the reconstruction techniques back and forth from generic dofs to point-values are not considered.
- Parameters
InitTime (float) – Optional, default value = 0.2. Set the pseudo-time first iteration step in the bisection iteration process
stol (float) – Optional, default value = 1e-8. The tolerance to stop the Secant iteration
gtol (float) – Optional, default value = 1e-5. The tolerance to stop the Bregman iteration
smaxit (float) – Optional, default value = 1000. The maximum iteration step for the secant method
gmaxit (float) – Optional, default value = 50. The maximum iteration step for the bregman method
eps (float) – Optional, default value = 1e-8. Safety net in the division of the weno weight and gradient division
verbose (boolean) – Optional, default value = False. True if the user wants a printout each point of the obtained values or False if no wished printout
Note
The arguments can be given as strings and will be converted accordingly if the string content is compatible with the exepected type.
-
W
(a, b, c, d, eps)¶ Function defining the weno weigths.
- Parameters
a (float) – The parameters to pass in order to compute the weno weights
b (float) – The parameters to pass in order to compute the weno weights
c (float) – The parameters to pass in order to compute the weno weights
d (float) – The parameters to pass in order to compute the weno weights
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
The Weno weight associated to the given coefficients
- Return type
Weight (float)
-
Deltas
(Phi, Points, eps)¶ Simple helper function that compute the differences involved in the Weno’s derivative approximation
- Parameters
Phi (Function callback) – An interpolation of the function to get the derivatives for
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns a tuple containing in order:
DeltaxMPhi (NbPoints float): Backward difference in x
DeltayMPhi (NbPoints float): Backward difference in y
DeltaxPPhi (NbPoints float): Forward difference in x
DeltayPPhi (NbPoints float): Forward difference in y
DDeltaxMPPhi (NbPoints float): Backward Forward difference in x
DDeltayMPPhi (NbPoints float): Backward Forward difference in y
- Return type
Differences (numpy arrays tuple)
-
GetGradient
(Phi, Points, eps)¶ Numerical gradient computation according to the WENO derivative computations
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback) – An interpolation of the function to get the derivatives for
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns np.array([Gradx, Grady]) where Gradx and Grady are as follows.
Gradx (NbPoints float): gradient value in the x direction at all the given points
Grady (NbPoints float): gradient value in the y direction at all the given points
- Return type
Grad (2d numpy array)
-
InterpolatePhi
(Points, Phi, *args)¶ Interpolator that considers either a function Phi or a set of Points and point-values Phi and returns an interpolator function. It also returns the Gradient interpolator of the function contructed from a gradient approximation at the given Points.
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback or NbPoints float) – Level set to reinitialise, either given as a function or a list of point values associated to points
args (Two function callback, optional) – The callback function defining the partial derivatives in x and y, in that order
- Returns
The interpolator of the function Phi InterpPhiDx (function callback): The interpolator of the derivative in x of Phi InterpPhiDy (function callback): The interpolator of the derivative in y of Phi
- Return type
InterpPhi (function callback)
Note
If function callbacks are given, they should be 2d->IR, taking the xy-values as a column list (NbPoints x 2) and return a 1d numpy vector
-
Redistance
(Points, Phi, *args, **kwargs)¶ Implementation of the main iteration strategy of the Bregman split for the redistancing problem.
- Parameters
Points (float 2D numpy array) – The points on which reinitialise the level-set on (NbPoints x 2)
Phi (function callback) – The (initial) level set function to redistance (should return a scalar value), given either as a float list associated to the points or as a function callback
args (Two function callback) – (optional) The callback function defining the partial derivatives in x and y, in that order
kwargs (keyword argument) – The only keyword argument available is “ConstructionPoints”, only active when Phi is given as a float list. It corresponds to where the given point values Phi have been computed. If the argument is not passed when Phi is a float list, the function assumes that Phi has been constructed on Points. Note that the Points should be contained within the convex hull of ConstructionPoints.
- Returns
The value of the redistanced level set function at the dofs (NbDofs)
- Return type
RPhi (float numpy array)
HopfLax Nested Bregman¶
Module that furnishes the tools to perform a redistancing through a Hopf-lax formula, solving the optimisation problem via a Bregman split approach. Following the work of Byungjoon Lee, “Revisiting the redistaning problem using the Hopf-Lax formula, 10.1016/j.jcp.2016.11.005” opf-Lax formula, 10.1016/j.jcp.2016.11.005”, with a one-sided bisection method instead of the proposed secant method to allow saddle points safely. In the minimasation of the sub-problem, it uses the global (also Bregman based) solver derived from the the work “Beyond Alternating Updates for Matrix Factorization with Inertial Bregman Proximal Gradient Algorithms” from Mahesh Chandra Mukkamala and Peter Ochs.
Note
The cheap way it is coded makes it slow
Requires a first coarse rough initial step, by default, this is done by a WENO procedure
Due to the root finding methd, the convergence is not strictly monotoneous in this implementation.
The local extrema and saddle points are avoided by the backtracking method of the nested solver
This toolbox is only valid for scalar-valued functions (representing a level set)
Warning
The code is not particularily robust and may fail on saddle points if a too low convexifying ratio is set
When using a high convexity ratio, the implementation becomes considerabily slower and may result in unwanted results if the InitBregmanTime is too large.
-
class
3_LESARD.SourceCode._Solver.Redistancers.HopfLax_NestedBregman.
Redistancing
(InitTime=0.2, CoarseningRatio=1, ConvexifyingRatio=1, stol=1e-08, gtol=1e-05, smaxit=1000, gmaxit=50, eps=1e-06, cfl=0.2, verbose=False)¶ Class furnishing the routines to define a redistanciation routine based on a HopfLax formula.
Note
The routine is updating the values of the level set at the degrees of freedom of the given mesh from the values that are stored there. Therefore, it is only suitable to be used on solution computed from a point-values dofs and not mixed techniques as the reconstruction techniques back and forth from generic dofs to point-values are not considered.
- Parameters
InitTime (float) – Optional, default value = 0.2. Set the pseudo-time first iteration step in the Bregman split iteration process (w.r.t. hmin: t0=0, t1 = BregmanInitTime*hmax)
CoarseningRatio (float) – Optional, default value = 1. How much coarse the cartesian grid for the initialisation should be (from the unstructured hmax)
ConvexifyingRatio (float) – Optional, default value = 1. The convexification ratio to be used in the Bregman split iteration process
stol (float) – Optional, default value = 1e-8. The tolerance to stop the Secant iteration
gtol (float) – Optional, default value = 1e-5. The tolerance to stop the Bregman iteration
smaxit (float) – Optional, default value = 1000. The maximum iteration step for the secant method
gmaxit (float) – Optional, default value = 50. The maximum iteration step for the bregman method
eps (float) – Optional, default value = 1e-6. Safety net in the division of the weno weight and gradient division
cfl (float) – Optional, default value = 0.2. CFL for the intialisation of Bregman split
verbose (boolean) – Optional, default value = False. True if the user wants a printout each point of the obtained values or False if no wished printout
Note
The arguments can be given as strings and will be converted accordingly if the string content is compatible with the exepected type.
-
W
(a, b, c, d, eps)¶ Function defining the weno weigths.
- Parameters
a (float) – The parameters to pass in order to compute the weno weights
b (float) – The parameters to pass in order to compute the weno weights
c (float) – The parameters to pass in order to compute the weno weights
d (float) – The parameters to pass in order to compute the weno weights
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
The Weno weight associated to the given coefficients
- Return type
Weight (float)
-
Deltas
(Phi, Points, eps)¶ Simple helper function that compute the differences involved in the Weno’s derivative approximation
- Parameters
Phi (Function callback) – An interpolation of the function to get the derivatives for
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns a tuple containing in order:
DeltaxMPhi (NbPoints float): Backward difference in x
DeltayMPhi (NbPoints float): Backward difference in y
DeltaxPPhi (NbPoints float): Forward difference in x
DeltayPPhi (NbPoints float): Forward difference in y
DDeltaxMPPhi (NbPoints float): Backward Forward difference in x
DDeltayMPPhi (NbPoints float): Backward Forward difference in y
- Return type
Differences (numpy arrays tuple)
-
GetGradient
(Phi, Points, eps)¶ Numerical gradient computation according to the WENO derivative computations
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback) – An interpolation of the function to get the derivatives for
eps (float) – The spacing to consider when selecting the neighborhing points in the gradient approximation (same in both x and y direction)
- Returns
- Returns np.array([Gradx, Grady]) where Gradx and Grady are as follows.
Gradx (NbPoints float): gradient value in the x direction at all the given points
Grady (NbPoints float): gradient value in the y direction at all the given points
- Return type
Grad (2d numpy array)
-
InterpolatePhi
(Points, Phi, *args)¶ Interpolator that considers either a function Phi or a set of Points and point-values Phi and returns an interpolator function. It also returns the Gradient interpolator of the function contructed from a gradient approximation at the given Points.
- Parameters
Points (float NbPoints x2) – The list of points where to consider the function for its interpolation
Phi (Function callback or NbPoints float) – Level set to reinitialise, either given as a function or a list of point values associated to points
args (Two function callback, optional) – The callback function defining the partial derivatives in x and y, in that order
- Returns
The interpolator of the function Phi InterpPhiDx (function callback): The interpolator of the derivative in x of Phi InterpPhiDy (function callback): The interpolator of the derivative in y of Phi
- Return type
InterpPhi (function callback)
Note
If function callbacks are given, they should be 2d->IR, taking the xy-values as a column list (NbPoints x 2) and return a 1d numpy vector
-
RoughEstimationAssetEquation
(Points, Phi)¶ Getting a first rough approximation of the evaluation solution of Phi of the redistancing equation through a fifth order weno in space and euler in time scheme.
- Parameters
Points (Function callback) – The list of points where to consider the function
Phi (Function callback) – Level set to reinitialise (either as a function or as an interpolator)
- Returns
The interpolator of the function Phi at the InitBregmanTime InterpPhiDx (function callback): The interpolator of the derivative in x of Phi at the InitBregmanTime InterpPhiDy (function callback): The interpolator of the derivative in y of Phi at the InitBregmanTime
- Return type
InterpPhi (function callback)
-
GetNextPhi
(Point, v, d, b, t, Interp, InterpDx, InterpDy)¶ Optimisation problem solver for the core Hopf-Lax formula, to be used within a zero finding method, using here the Bregman Split algorithm.
- Parameters
Point (float numpy array) – Point [x,y] at which we want to update the Level set value
v (float numpy array) – The initial value [x,y] of v to consider in the Bregman split
d (float numpy array) – The initial value [x,y] of d to consider in the Bregman split
b (float numpy array) – The initial value [x,y] of b to consider in the Bregman split
t (float) – Distance from the Point point from which the value of the optimiser should be projected on the ball. It should be non-negative.
Interp (function callback) – The interpolator (or function itself) representing conitnously the level-set in the 2D plane (can be None or Nan outside the domain if there is no data there)
InterpDx (function callback) – The interpolator (or function itself) representing conitnously the derivative in x of level-set in the 2D plane (can be None or Nan outside the domain if there is no data there)
InterpDy (function callback) – The interpolator (or function itself) representing conitnously the derivative in y of level-set in the 2D plane (can be None or Nan outside the domain if there is no data there)
- Returns
Value of the level set at the found optimizer d (float numpy array): The array [x,y] of the minimizer coordinates gfailed (float): Returns a failure flag: 0 if we could not converge to the wished tolerance within the given amount of steps, 1 otherwise
- Return type
Phi1 (float)
-
Redistance
(Points, Phi, *args, **kwargs)¶ Implementation of the main iteration strategy of the Bregman split for the redistancing problem.
- Parameters
Points (float 2D numpy array) – The points on which reinitialise the level-set on (NbPoints x 2)
Phi (function callback) – The (initial) level set function to redistance (should return a scalar value), given either as a float list associated to the points or as a function callback
args (Two function callback) – (Optional) The callback function defining the partial derivatives in x and y, in that order
kwargs (keyword argument) – The only keyword argument available is “ConstructionPoints”, only active when Phi is given as a float list. It corresponds to where the given point values Phi have been computed. If the argument is not passed when Phi is a float list, the function assumes that Phi has been constructed on Points. Note that the Points should be contained within the convex hull of ConstructionPoints.
- Returns
The value of the redistanced level set function at the dofs (NbDofs)
- Return type
RPhi (float numpy array)