Adhesion.System package
Submodules
Adhesion.System.Factory module
- Adhesion.System.Factory.make_system(surface, interaction, substrate=None, communicator=<NuMPI.MPIStub.Intracomm object>, physical_sizes=None, system_class=None, **kwargs)[source]
Factory function for contact systems. The returned object is always of a subtype of SystemBase.
Parameters:
- substrateContactMechanics.Substrate
Defines the solid mechanics in the substrate
- interactionAdhesion.Interaction
Defines the contact formulation
- surfaceSurfaceTopography.Topography
Defines the profile.
- returns:
system
- rtype:
child class of SystemBase
Adhesion.System.PlasticSystemSpecialisations module
implements plastic mapping algorithms for contact systems
- class Adhesion.System.PlasticSystemSpecialisations.PlasticSmoothContactSystem(substrate, interaction, surface)[source]
Bases:
SmoothContactSystem
This system implements a simple penetration hardness model.
Adhesion.System.SmoothSystemSpecialisations module
implements the periodic and non-periodic smooth contact systems
- class Adhesion.System.SmoothSystemSpecialisations.BndSet(large, small)
Bases:
tuple
- large
Alias for field number 0
- small
Alias for field number 1
- class Adhesion.System.SmoothSystemSpecialisations.FastSmoothContactSystem(substrate, interaction, surface, margin=4)[source]
Bases:
SmoothContactSystem
This proxy class tries to take advantage of the system-physical_sizes independence of non-periodic FFT-solved systems by determining the required minimum system physical_sizes, encapsulating a SmoothContactSystem of that physical_sizes and run it almost transparently instead of the full system. It’s almost transparent, because by its nature, the small system does not compute displacements everywhere the large system exists. Therefore, for full-domain output, an extra evaluation step must be taken. But since using the small system significantly increases optimisation speed, this tradeoff seems to be well worth the trouble.
- exception BabushkaBoundaryError[source]
Bases:
Exception
Called when the choosen physical_sizes of the Babushka System had an Influence on the Solution. In the Future, one may catch this error, increase the margin and restart minimization
- exception FreeBoundaryError(message, disp)[source]
Bases:
Exception
called when the supplied system cannot be computed with the current constraints
- callback(force=False)[source]
Simple callback function that can be handed over to scipy’s minimize to get updates during minimisation Parameters: ———- force: bool, optional
whether to include the norm of the force vector in the update message (default False)
- check_margins()[source]
Checks whether the safety margin is sufficient (i.e. the outer ring of the force array equals zero
- check_margins_babushka()[source]
Inspired from check_margins but acts on the babushka System: When the pressures are non-zero at the boundary of the babushka-area, boundaries have an effect on the Solution this raises an error because the area has been chosen too small and the
- compute_babushka_bounds(babushka_nb_grid_pts)[source]
returns a list of tuples that contain the equivalent slices in the small and the large array. It differentiates between nb_grid_pts and nb_domain_grid_pts. Parameters: babushka_nb_grid_pts – nb_grid_pts of smaller scale
- create_babushka(offset, disp0=None)[source]
Create a (much smaller) system with just the contacting patch plus the margin
- evaluate(disp, offset, pot=True, forces=False)[source]
Compute the energies and forces in the system for a given displacement field
Parameters:
- disp: ndarray
displacement field, in the shape of system.substrate.nb_subdomain_grid_pts
- offset: float
determines indentation depth, constant value added to the heights (system.topography)
- pot: bool, optional
Wether to evaluate the energy, default True
- forces: bool, optional
Wether to evaluate the forces, default False
- logger: ContactMechanics.Tools.Logger
informations of current state of the system will be passed to logger at every evaluation
- static handles(substrate_type, interaction_type, surface_type, is_domain_decomposed)[source]
returns whether this class (in practice a subclass) handles this combination of types Keyword Arguments: substrate_type – self-explanatory surface_type – self-explanatory is_domain_decomposed: some systems cannot handle parallel computation
- minimize_proxy(offset, disp0=None, method='L-BFGS-B', options=None, gradient=True, lbounds=None, ubounds=None, tol=None, callback=None, deproxify_everytime=True)[source]
Convenience function. Eliminates boilerplate code for most minimisation problems by encapsulating the use of scipy.minimize for common default options. In the case of smart proxy systems, this may also encapsulate things like dynamics computation of safety margins, extrapolation of results onto the proxied system, etc.
Parameters: offset – determines indentation depth disp0 – (default zero) initial guess for displacement field. If
not chosen appropriately, results may be unreliable.
- method – (defaults to L-BFGS-B, see scipy documentation). Be sure
to choose method that can handle high-dimensional parameter spaces.
- options – (default None) options to be passed to the minimizer
method
gradient – (default True) whether to use the gradient or not tol – (default None) tolerance for termination. For detailed
control, use solver-specific options.
- callback – (default None) callback function to be at each iteration
as callback(disp_k) where disp_k is the current displacement vector. Instead of a callable, it can be set to ‘True’, in which case the system’s default callback function is called.
- objective(offset, disp0=None, gradient=False)[source]
See super().objective for general description this method’s purpose. Difference for this class wrt ‘dumb’ systems: Needs an initial guess for the displacement field in order to estimate the contact area. returns both the objective and the adapted ininial guess as a tuple Parameters: ———– offset: float
determines indentation depth
- gradient: bool, optional
whether the gradient is supposed to be used (default False)
- disp0: ndarray of float, optional
initial guess for displacement field. If not chosen appropriately, results may be unreliable. (default zero)
- shape_minimisation_input(in_array)[source]
For minimisation of smart systems, the initial guess array (e.g. displacement) may have a non-intuitive shape and physical_sizes (The problem physical_sizes may be decreased, as for free, non-periodic systems, or increased as with augmented-lagrangian-type issues). Use the output of this function as argument x0 for scipy minimisation functions. Also, if your initial guess has a shape that makes no sense, this will tell you before you get caught in debugging scipy-code
Arguments: in_array – array with the initial guess. has the intuitive shape you
think it has
Adhesion.System.Systems module
Defines the interface for Adhesion systems
- class Adhesion.System.Systems.BoundedSmoothContactSystem(substrate, interaction, surface)[source]
Bases:
SmoothContactSystem
- compute_attractive_coordinates()[source]
returns an array of all coordinates, where contact pressure is attractive. Useful for evaluating the number of contact islands etc.
- compute_nb_attractive_pts()[source]
compute and return the number of contact points under attractive pressure.
- compute_nb_repulsive_pts()[source]
compute and return the number of contact points under repulsive pressure.
- compute_repulsive_coordinates()[source]
returns an array of all coordinates, where contact pressure is repulsive. Useful for evaluating the number of contact islands etc.
- compute_repulsive_force()[source]
computes and returns the sum of all repulsive forces
Assumptions: there
- property contact_zone
returns an array of all coordinates, where contact pressure is repulsive. Useful for evaluating the number of contact islands etc.
- static handles(*args, **kwargs)[source]
returns whether this class (in practice a subclass) handles this combination of types Keyword Arguments: substrate_type – self-explanatory surface_type – self-explanatory is_domain_decomposed: some systems cannot handle parallel computation
- minimize_proxy(offset=0, initial_displacements=None, method='L-BFGS-B', gradient=True, lbounds=None, ubounds=None, callback=None, logger=None, **kwargs)[source]
Convenience function. Eliminates boilerplate code for most minimisation problems by encapsulating the use of scipy.minimize for common default options. In the case of smart proxy systems, this may also encapsulate things like dynamics computation of safety margins, extrapolation of results onto the proxied system, etc.
Parameters: offset : float
determines indentation depth
- initial_displacements(default zero)
initial guess for displacement field. If not chosen appropriately, results may be unreliable.
- methodstring or callable
(defaults to L-BFGS-B, see scipy documentation). Be sure to choose method that can handle high-dimensional parameter spaces.
- optionsdict
(default None) options to be passed to the minimizer method
- gradientbool
(default True) whether to use the gradient or not
- lboundsarray of shape substrate.subdomain_nb_grid_pts or
substrate.topography_subdomain_nb_grid_pts or “auto”
- (default None)
nodal ceiling/floor
- uboundsarray of shape substrate.subdomain_nb_grid_pts or
substrate.topography_subdomain_nb_grid_pts
(default None)
- tolfloat
(default None) tolerance for termination. For detailed control, use solver-specific options.
- callbackcallable
(default None) callback function to be at each iteration
as callback(disp_k) where disp_k is the current displacement vector. Instead of a callable, it can be set to ‘True’, in which case the system’s default callback function is called.
- logger :
(default None) log information at every objective evaluation.
- class Adhesion.System.Systems.SmoothContactSystem(substrate, interaction, surface)[source]
Bases:
SystemBase
For smooth contact mechanics (i.e. the ones for which optimization is only kinda-hell
- callback(force=False)[source]
Simple callback function that can be handed over to scipy’s minimize to get updates during minimisation Parameters: ———- force: bool, optional
whether to include the norm of the force vector in the update message (default False)
- compute_attractive_contact_area()[source]
computes and returns the are where contact pressure is attractive
- compute_attractive_coordinates()[source]
returns an array of all coordinates, where contact pressure is attractive. Useful for evaluating the number of contact islands etc.
- compute_mean_gap()[source]
mean of the gap in the the physical domain (means excluding padding region for the FreeFFTElasticHalfspace)
- compute_nb_attractive_pts()[source]
compute and return the number of contact points under attractive pressure. Note that this is of no physical interest, as it is a purely numerical artefact
- compute_nb_contact_pts()[source]
compute and return the number of contact points. Note that this is of no physical interest, as it is a purely numerical artefact
- compute_nb_repulsive_pts()[source]
compute and return the number of contact points under repulsive pressure. Note that this is of no physical interest, as it is a purely numerical artefact
- compute_repulsive_contact_area()[source]
computes and returns the area where contact pressure is repulsive
- compute_repulsive_coordinates()[source]
returns an array of all coordinates, where contact pressure is repulsive. Useful for evaluating the number of contact islands etc.
- evaluate(disp, offset, pot=True, forces=False, logger=None)[source]
Compute the energies and forces in the system for a given displacement field
Parameters:
- disp: ndarray
displacement field, in the shape of system.substrate.nb_subdomain_grid_pts
- offset: float
determines indentation depth, constant value added to the heights (system.topography)
- pot: bool, optional
Wether to evaluate the energy, default True
- forces: bool, optional
Wether to evaluate the forces, default False
- logger: ContactMechanics.Tools.Logger
informations of current state of the system will be passed to logger at every evaluation
- evaluate_k(disp_k, gap, offset, mw=False, pot=True, forces=False, logger=None)[source]
Compute the energies and forces in the system for a given displacement field in fourier space.
- Parameters:
disp_k (ndarray) – displacement field in fourier space.
gap (ndarray) – displacement field in real space, in the shape of system.substrate.nb_subdomain_grid_pts
offset (float) – determines indentation depth, constant value added to the heights (system.topography)
mw (bool, optional) – when mass weighting is required then set this to TRUE.
pot (bool, optional) – Wether to evaluate the energy, default True
forces (bool, optional) – Wether to evaluate the forces, default False
logger (ContactMechanics.Tools.Logger) – informations of current state of the system will be passed to logger at every evaluation.
- static handles(substrate_type, interaction_type, surface_type, comm)[source]
returns whether this class (in practice a subclass) handles this combination of types Keyword Arguments: substrate_type – self-explanatory surface_type – self-explanatory is_domain_decomposed: some systems cannot handle parallel computation
- property nb_grid_pts
For systems, nb_grid_pts can become non-trivial
- objective(offset, disp0=None, gradient=False, logger=None)[source]
This helper method exposes a scipy.optimize-friendly interface to the evaluate() method. Use this for optimization purposes, it makes sure that the shape of disp is maintained and lets you set the offset and ‘forces’ flag without using scipy’s cumbersome argument passing interface. Returns a function of only disp
Parameters:
- disp0: ndarray
unused variable, present only for interface compatibility with inheriting classes
- offset: float
determines indentation depth, constant value added to the heights (system.topography)
- gradient: bool, optional
Wether to evaluate the gradient, default False
- logger: ContactMechanics.Tools.Logger
informations of current state of the system will be passed to logger at every evaluation
- returns:
- function(disp)
Parameters: disp: an ndarray of shape
system.substrate.nb_subdomain_grid_pts displacements
- Returns:
energy or energy, gradient
- objective_k_float(offset, gradient=False, logger=None)[source]
Returns callable objective as needed by scipy minimizers.
The optimisation varialbe is the halfcomplex fourier transform of the gap.
This helper method interface to the evaluate_k() method without preconditioning active.
\[\begin{split}\frac{1}{2(n_x n_y)} \tilde{u}\tilde{K} \bar{\tilde{u}} + \phi(F^{-1}(\tilde{u} - \tilde{h})) \\\end{split}\]preconditioned problem:
\[ \begin{align}\begin{aligned}\begin{split} \tilde{v} = \tilde{k}^{\frac{1}{2}} \tilde{u} \\\end{split}\\ \frac{1}{2(n_x n_y)} \tilde{v} \bar{\tilde{v}} + \phi(F^{-1}(\frac{\tilde{v}}{\tilde{k}^{\frac{1}{2}}} - \tilde{h}))\end{aligned}\end{align} \]\ we solve for variable \(\tilde{v}\).
Parameters:
- offset: float
determines indentation depth, constant value added to the heights (system.topography)
- gradient: bool, optional
Whether to evaluate the gradient, default False
- logger: ContactMechanics.Tools.Logger
information of current state of the system will be passed to logger at every evaluation
- returns:
- param disp_k:
- type disp_k:
an ndarray in fourier space
- rtype:
energy, gradient_k_float
- rtype:
function(disp_k)
- preconditioned_objective(offset, gradient=False, logger=None)[source]
This helper method interface to the evaluate_k() method with preconditioning active. That is, it tries to solve a simpler problem formulated using,
original problem:
\[\frac{1}{2(n_x n_y)} \tilde{u}\tilde{K} \bar{\tilde{u}} + \phi(F^{-1}(\tilde{u} - \tilde{h}))\]preconditioned problem:
\[ \begin{align}\begin{aligned}\begin{split}\tilde{v} = \tilde{k}^{\frac{1}{2}} \tilde{u} \\\end{split}\\ \frac{1}{2(n_x n_y)} \tilde{v}\bar{\tilde{v}} + \phi(F^{-1}(\frac{\tilde{v}}{\tilde{k}^{\frac{1}{2}}} - \tilde{h}))\end{aligned}\end{align} \]we solve for variable \(\tilde{v}\).
Parameters:
- offset: float
determines indentation depth, constant value added to the heights (system.topography)
- gradient: bool, optional
Whether to evaluate the gradient, default False
- logger: ContactMechanics.Tools.Logger
information of current state of the system will be passed to logger at every evaluation
- returns:
- param disp_k:
- type disp_k:
an ndarray in fourier halfcomplex space
- returns:
energy (scalar) – energy of the system
force_h (an halfcomplex array of shape(disp_k)) – force of the system
- rtype:
function(disp_k)
- primal_evaluate(disp, gap, pot=True, forces=False, logger=None)[source]
Compute the energies and forces in the system for a given displacement and gap..
Parameters:
- disp: ndarray
displacement field, in the shape of system.substrate.nb_subdomain_grid_pts
- gap: ndarray
gap , in the shape of system.substrate.nb_subdomain_grid_pts
- pot: bool, optional
Whether to evaluate the energy, default True
- forces: bool, optional
Whether to evaluate the forces, default False
- logger: ContactMechanics.Tools.Logger
information of current state of the system will be passed to logger at every evaluation
- primal_hessian_product(gap, des_dir)[source]
Returns the hessian product of the primal_objective function.
- primal_objective(offset, disp0=None, gradient=False, logger=None)[source]
To solve the primal objective using gap as the variable. Can be fed directly to standard solvers ex: scipy solvers etc and returns the elastic energy and it’s gradient (negative of the forces) as a function of the gap.
- Parameters:
gap (float) – gap between the contact surfaces.
offset (float) – constant value to add to the surface heights
pot ((default False)) –
gradient ((default True)) –
- Returns:
energy (float) – value of energy(scalar value).
force (float,array) – value of force(array).
Notes
Objective:
\[\begin{split}\min_u f = \frac{1}{2} u_i K_{ij} u_j + \phi (u_{ij})\\ \\ \nabla f = K_{ij} u_j + \phi^{\prime} \text{ which is the force} \\\end{split}\]