lugre_model Derived Type

type, public, extends(friction_model) :: lugre_model

Defines the Lu-Gre friction model.

The Lu-Gre model is a bristle-type model that attempts to describe friction using a bristle interpretation of the frictional surfaces.
The bristle-type models assume that the frictional behavior is represented by an average deflection of elastic springs. These springs have their own stiffness and damping properties and act as a typical spring-damper pair under small velocities; however, once sufficient velocity occurs, the bristles slip resulting in Coulomb-like sliding behavior.

The Lu-Gre model is defined as follows.

where:

Friction Force

Normal Force

Velocity

Coulomb Friction Coefficient

Static Friction Coefficient

Bristle Stiffness

Bristle Damping Coefficient

Viscous Damping Coefficient

Stribeck Curve Shape Factor

Stribeck Velocity Coefficient


Components

Type Visibility Attributes Name Initial
real(kind=real64), public :: coulomb_coefficient

The Coulomb (dynamic) friction coefficient.

real(kind=real64), public :: damping

The bristle damping coefficient.

real(kind=real64), public :: shape_parameter

The Stribeck curve shape parameter.

real(kind=real64), public :: static_coefficient

The static friction coefficient.

real(kind=real64), public :: stiffness

The bristle stiffness.

real(kind=real64), public :: stribeck_velocity

The Stribeck velocity term.

real(kind=real64), public :: viscous_damping

The viscous damping coefficient.


Type-Bound Procedures

procedure, public :: constraint_equations => fmdl_constraints

  • private subroutine fmdl_constraints(this, t, x, dxdt, nrm, f, rst)

    Overload this routine to establish constraings for the model to be enforced as part of the fitting operation.

    Arguments

    Type IntentOptional Attributes Name
    class(friction_model), intent(in) :: this

    The friction_model object.

    real(kind=real64), intent(in), dimension(:) :: t

    An N-element array containing the time points at which the data to be fit was sampled.

    real(kind=real64), intent(in), dimension(:) :: x

    An N-element array containing the relative motion data.

    real(kind=real64), intent(in), dimension(:) :: dxdt

    An N-element array containing the relative velocity data.

    real(kind=real64), intent(in), dimension(:) :: nrm

    An N-element array containing the normal force data.

    real(kind=real64), intent(in), dimension(:) :: f

    An N-element array containing the friction force data.

    real(kind=real64), intent(out), dimension(:) :: rst

    An M-element array where the results of the constraint equations will be written. M must be equal to the number of constraint equations for the model.

procedure, public :: evaluate => lg_eval

  • private function lg_eval(this, t, x, dxdt, nrm, svars) result(rst)

    Evaluates the Lu-Gre friction model.

    Arguments

    Type IntentOptional Attributes Name
    class(lugre_model), intent(inout) :: this

    The lugre_model object.

    real(kind=real64), intent(in) :: t

    The current simulation time value.

    real(kind=real64), intent(in) :: x

    The current value of the relative position between the contacting bodies.

    real(kind=real64), intent(in) :: dxdt

    The current value of the relative velocity between the contacting bodies.

    real(kind=real64), intent(in) :: nrm

    The current normal force between the contacting bodies.

    real(kind=real64), intent(in), optional, dimension(:) :: svars

    An optional array containing any internal state variables the model may rely upon.

    Return Value real(kind=real64)

    The friction force.

procedure, public :: fit => fmdl_fit

  • private subroutine fmdl_fit(this, t, x, v, f, n, weights, maxp, minp, alpha, integrator, controls, settings, info, stats, fmod, resid, err)

    Attempts to fit a friction model to the supplied data using a Levenberg-Marquardt solver.

    Arguments

    Type IntentOptional Attributes Name
    class(friction_model), intent(inout), target :: this

    The friction model. On output, the model is updated with the final, fitted parameters.

    real(kind=real64), intent(in), target, dimension(:) :: t

    An N-element array containing the time points at which the friction data was sampled. This array must contain monotonically increasing data.

    real(kind=real64), intent(in), target, dimension(:) :: x

    An N-element array containing the relative position data.

    real(kind=real64), intent(in), target, dimension(:) :: v

    An N-element array containing the relative velocity data.

    real(kind=real64), intent(in), target, dimension(:) :: f

    An N-element array containing the friction force data.

    real(kind=real64), intent(in), target, dimension(:) :: n

    An N-element array containing the normal force data.

    real(kind=real64), intent(in), optional, dimension(:) :: weights

    An optional N-element array that can be used to weight specific data points. The default is an array of all ones such that all points are weighted equally.

    real(kind=real64), intent(in), optional, dimension(:) :: maxp

    An M-element array (M = the number of model parameters) containing a maximum limit for each model parameter.

    real(kind=real64), intent(in), optional, dimension(:) :: minp

    An M-element array containing the minimum limit for each model parameter.

    real(kind=real64), intent(in), optional :: alpha

    An optional input that defines the significance level at which to evaluate the confidence intervals. The default value is 0.05 such that a 95% confidence interval is calculated.

    class(ode_integrator), intent(inout), optional, target :: integrator

    An optional input, used in the event the model has internal state variables, that provides integration of the state equations. The defaults is a 4th order Rosenbrock method.

    type(iteration_controls), intent(in), optional :: controls

    An optional input providing custom iteration controls.

    type(lm_solver_options), intent(in), optional :: settings

    An optional input providing custom settings for the solver.

    type(convergence_info), intent(out), optional :: info

    An optional output that can be used to gain information about the iterative solution and the nature of the convergence.

    type(regression_statistics), intent(out), optional, dimension(:) :: stats

    An optional output array of M-elements that can be used to retrieve statistical information regarding the fit of each of the M model parameters.

    real(kind=real64), intent(out), optional, target, dimension(:) :: fmod

    An optional N-element array used to provide the fitted model results.

    real(kind=real64), intent(out), optional, target, dimension(:) :: resid

    An optional N-element array containing the fitted residuals.

    class(errors), intent(inout), optional, target :: err

    An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling.

procedure, public :: from_array => lg_from_array

  • private subroutine lg_from_array(this, x, err)

    Converts an array into the parameters for the friction model.

    Arguments

    Type IntentOptional Attributes Name
    class(lugre_model), intent(inout) :: this

    The lugre_model object.

    real(kind=real64), intent(in), dimension(:) :: x

    The array of parameters. See parameter_count to determine the size of this array. The parameter order is as follows:

    1. static_coefficient

    2. coulomb_coefficient

    3. stribeck_velocity

    4. stiffness

    5. damping

    6. viscous_damping

    7. shape_parameter

    class(errors), intent(inout), optional, target :: err

    An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling.

procedure, public :: get_constraint_equation_count => fmdl_get_constraint_count

  • private pure function fmdl_get_constraint_count(this) result(rst)

    Gets the number of constraint equations the model requires to be satisfied when fitting to data.

    Arguments

    Type IntentOptional Attributes Name
    class(friction_model), intent(in) :: this

    The friction_model object.

    Return Value integer(kind=int32)

    The number of constraint equations.

procedure, public :: get_state_variable_count => lg_get_state_var_count

  • private pure function lg_get_state_var_count(this) result(rst)

    Gets the number of internal state variables used by the model.

    Arguments

    Type IntentOptional Attributes Name
    class(lugre_model), intent(in) :: this

    The lugre_model object.

    Return Value integer(kind=int32)

    The internal state variable count.

procedure, public :: has_internal_state => lg_has_state_vars

  • private pure function lg_has_state_vars(this) result(rst)

    Returns a value stating if the model relies upon internal state variables.

    Arguments

    Type IntentOptional Attributes Name
    class(lugre_model), intent(in) :: this

    The lugre_model object.

    Return Value logical

    Returns true if the model utilizes internal state variables; else, returns false.

procedure, public :: parameter_count => lg_parameter_count

  • private pure function lg_parameter_count(this) result(rst)

    Gets the number of model parameters.

    Arguments

    Type IntentOptional Attributes Name
    class(lugre_model), intent(in) :: this

    The lugre_model object.

    Return Value integer(kind=int32)

    The number of model parameters.

procedure, public :: reset => fmdl_reset

  • private subroutine fmdl_reset(this)

    Resets the friction model to it's original state.

    Arguments

    Type IntentOptional Attributes Name
    class(friction_model), intent(inout) :: this

    The friction_model object.

procedure, public :: state => lg_state_model

  • private subroutine lg_state_model(this, t, x, dxdt, nrm, svars, dsdt)

    Evaluates the time derivatives of the internal friction state model.

    Arguments

    Type IntentOptional Attributes Name
    class(lugre_model), intent(inout) :: this

    The lugre_model object.

    real(kind=real64), intent(in) :: t

    The current simulation time value.

    real(kind=real64), intent(in) :: x

    The current value of the relative position between the contacting bodies.

    real(kind=real64), intent(in) :: dxdt

    The current value of the relative velocity between the contacting bodies.

    real(kind=real64), intent(in) :: nrm

    The current normal force between the contacting bodies.

    real(kind=real64), intent(in), dimension(:) :: svars

    An N-element array containing any internal state variables the model may rely upon.

    real(kind=real64), intent(out), dimension(:) :: dsdt

    An N-element array where the state variable derivatives are to be written.

procedure, public :: to_array => lg_to_array

  • private subroutine lg_to_array(this, x, err)

    Converts the parameters of the friction model into an array.

    Arguments

    Type IntentOptional Attributes Name
    class(lugre_model), intent(in) :: this

    The lugre_model object.

    real(kind=real64), intent(out), dimension(:) :: x

    The array used to store the parameters. See @ref parameter_count to determine the size of this array. The parameter order is as follows:

    1. static_coefficient

    2. coulomb_coefficient

    3. stribeck_velocity

    4. stiffness

    5. damping

    6. viscous_damping

    7. shape_parameter

    class(errors), intent(inout), optional, target :: err

    An optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling.