Thanks to visit codestin.com
Credit goes to github.com

Skip to content

asmusskar/ALVA

Folders and files

NameName
Last commit message
Last commit date

Latest commit

 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 

Repository files navigation

DOI DOI GitHub release (latest by date including pre-releases) GitHub View ALVA on File Exchange

ALVA

Adaptive Layered Viscoelastic Analysis (ALVA) is a MATLAB-package for pavement modeling. The aim with the software is to equip the civil engineering community with an advanced pavement modeling tool and computer package that is highly adaptive, transparent and open-access, capable of supporting current and future pavement evaluation needs.

Figure 1: ALVA roadmap.

The components of the ALVA software (see Figure 1) is briefly described below:

  • main.m - main script for defining vehicle loading conditions, pavement structure geometry and material properties, numerical parameters (for balancing accuracy and efficiency), analysis type and evaluation points (i.e., location of the output response), and post-processing of the results. The main.m scripts can be found in the ../ALVA/examples_ folder.

  • init_LET.m - initialize the analysis of a layered elastic half space selecting between: (a) let_response_full.m - evaluates the response (displacements, stresses and strains) of a layered elastic half-space model or (b) let_response_polfit.m - evaluates the displacements of a layered elastic half space model at the surface only.

  • init_LET.m also works as an interface application for simulating moving loads in the Viscoelastic module; based on the asphalt material properties from VE_moduli.m the module is engaged to calculate the elastic responses in VE_response.m, which is subsequently used to predict the linear viscoelastic response of a moving vehicle in VE_simulation.m.

  • The main support functions in init_LET.m are (a) arb_funct.m - evaluates the coefficients of integration that embody the layered system properties and connectivity and (b) numint_coeff.m – organizes the integration points and weights for numerical integration for the user defined number of Bessel roots from besselroots.m and gauss points and weights between Bessel roots from lookup_gauss.m.

  • Additional support functions for efficient surface displacement calculations in let_response_polfit.m is arb_func_polfit.m – interpolates intermediate coefficients of integration, polfit_int.m – evaluates coefficient proportional integrals, and polfit_abc.m – evaluates coefficients a, b and c for the polynomial. Examples on ‘main.m’ scripts can be found in the ../ALVA/examples folder and all underlying components (i.e., source files) described above in the ../ALVA/basic folder.

Main input parameters

Layered elastic analysis

The core algorithm behind this package is based on linear elastic theory (LET), i.e., the classic formulation for an N-layered half-space, shown in Figure 2

Figure 2: N-layered half-space model

In this model all layers are assumed linear elastic, isotropic, homogeneous and weightless. The model inputs include Young’s modulus En, Poisson’s ratio νn, and layer thickness tn (where n denotes the layer number). In addition ALVA allows users to define horizontally oriented springs, kh, operating at the top or bottom layer interfaces to model imperfect interface conditions. This model is engaged to calculate the response at any point, Aj, of interest and for a given set of uniformly distributed circular loadings with load radius, a, and pressure q). An overview of LET model assumptions and solution procedure is given in Khazanovich and Wang (2007).

  • Analysis type
alva.analysis = 'Full';
% 1) 'Full'     : Conventional full integration with one-step Richardson
%                 extrapolation (for improved convergence near surface) 
%                 applied for evaluation of displacements and stresses
% 2) 'PolFit'   : Use polynomial fit technique to reduce integral length
%                 according to Andersen et al. (2018) for evaluation of
%                 surface displacements
  • Interface

Figure 3: Close-up of interface showing a mechanical representation of the horizontal springs kh

alva.bond = 'Bonded';
% 'Bonded'       : Full bonding between layers
% 'Slip'         : Interface bonding factor
% 'Frictionless' : No bonding between layers

 % Interface bonding/horizontal spring [MPa/mm]: high -> 'Bonded' / low -> 'Frictionless'
alva.kh = [1e6 1e6];        
  • Numerical parameters

Figure 4: Visualisation of numerical parameters, i.e., Bessel function zero points 'N' and Gauss points 'n' for solving the inverse Hankel transform. In the Figure a Bessel function of the first kind of order zero is plotted with 'N'='n'=5. Zero points are shown as red circles, and Gauss points as blue and green circles for the first and second interval, respectively.

alva.N  = 300;  % Number of Bessel zero points in numerical integration
alva.n  = 30;   % Number of Gauss points points between zero points.
  • Pavement material properties (minimum two layers required) - see Figure 1
alva.zi = [200 760];         % Depth of first n-1 layers from the 
                             % surface [mm]: last z = inf, and should not 
                             % be added NB: zi(i) > zi(i-1) > z(i-2)...
alva.E  = [5000 200 50];     % Layer Young's moduli [MPa]
alva.nu = [0.35 0.40 0.45];  % Layer Poisson's ratio [-]

  • Load configuration

Figure 5: Load configuration example: dual wheel load.

alva.q  = [0.7
           0.7];         % Load pressure [MPa] (uniform vertical pressure)
alva.a  = [106.6
           106.6];       % Load radii [mm] (circular load)
alva.Xl = [-170 0
           170  0];      % Load positions [mm]: [x1 y1; x2 y2;..xi yi];
  • Location of evaluation points: (x1 y1 z1; x2 y2 z2;..)

Response at the centre of one tire load, i.e., x=170 and y=0 with varying depth z=(0,259.99,260.1,760.1), and midpoint between tires, i.e., x=0 and y=0 with varying depth z=(0,259.99,260.1,760.1)

alva.Xd = [170 0 0; 170 0 259.99; 170 0 260.1; 170 0 760.1;
             0 0 0;   0 0 259.99;   0 0 260.1;   0 0 760.1];   
  • Initialize system and get response for the layered elastic model

The output response is given in cartesian coordinates, i.e., for all points Aj we have:

Displacements:

Stresses:

Strains:

alva = init_LET(alva);

% Displacements [mm]
ux    = alva.ux;   
uy    = alva.uy;  
uz    = alva.uz;        

% Stresses [MPa]
sigx  = alva.sigx;      
sigy  = alva.sigy;       
sigz  = alva.sigz;       
sigxy = alva.sigxy;      
sigyz = alva.sigyz;      
sigxz = alva.sigxz;   

% Strains [10e-6]
epsx  = alva.epsx.*1e6;  
epsy  = alva.epsy.*1e6;  
epsz  = alva.epsz.*1e6; 
epsxy = alva.epsxy.*1e6; 
epsyz = alva.epsyz.*1e6;
epsxz = alva.epsxz.*1e6;

Layered viscoelastic analysis

The viscoelastic response is approximated based on the LET calculations utilizing the methodology and load scheme suggested by Levenberg (2016) (see Figure 6).

  • Creep compliance curve for the asphalt layer

Viscoelastic layers are associated with a creep compliance of the form (Smith, 1971):

wherein D0 and D are the short and long time compliances (respectively), and shape parameters τD and nD, controlling the transition between D0 and D.

% Viscoelastic properties of asphalt layer E(1)
D0   = 2.5e-5;    % Instantaneous or glassy compliance, [1/MPa]
Dinf = 1.0e-2;    % Long time equilibrium or rubbery compliance [1/MPa]
nD   = 0.35;      % Slope of the creep curve in the transient region [-]
tauD = 1e3;       % Retardation time of the material response [s]

% Determine creep compliance curve and relaxation modulus
alva.ti = 12;     % Number increments on compliance curve

alva    = VE_moduli(D0,Dinf,tauD,nD,alva);
  • Location of evaluation point and mesh for simulating moving load

Figure 6: Load scheme for simulating a moving load

The load moves in a straight line from x=-x0 (Start) to x=x0 (End). The travel path is decomposed into N intervals (i=1,…,N), each Δx long. The point of response evaluation Aj is indicated in the Figure; this point is located near the middle of the travel path (i.e., x-coordinate of zero), at y-coordinate y0 and depth z0 below the surface.

rep  = [0 0 0];                    % [x y z]-coordinate of evaluation point
nels = 200;                        % Number of elements [-]
Vkh  = 60;                         % Vehicle speed [km/h]
V    = Vkh/3.6*1e3;                % Vehicle speed [mm/s]
x0   = 5000;                       % Start / end of mesh [mm]
dx0  = 2*x0/nels;                  % Mesh increment [mm]
dt   = dx0/V;  alva.dt = dt;       % Time increment [s]
tt   = 2*x0/V; alva.tt = tt;       % Total travel time [s]
xx   = (0:dx0:x0);                 % Mesh x-direction
yy   = rep(2)*ones(length(xx),1);  % Mesh y-direction
zz   = rep(3)*ones(length(xx),1);  % Mesh z-direction
Xr   = [xx' yy zz]; alva.Xr = Xr;  % Full mesh matrix
  • Calculate response
% Calculate linear elastic response at different times
alva = VE_response(alva.E,alva);

% Simulate moving load

% Select output response: displacements (dxm, dym or dzm), stresses (sigxm, 
% sigym, sigzm, sigxym, sigyzm or sigxzm), strains (epsxm, epsym, epszm, 
% epsxym, epsyzm or epsxzm)

letres = alva.dzm;                   % output response dz                                           
veres  = VE_simulation(letres,alva); % run simulation

Code validation

ALVA comes with six validation cases/examples (i.e., main.m scripts), comparing ALVA to existing pavement analysis sofware and analytical formulations. The files can be found in the in the ../ALVA/examples folder and the results obtained with independent codes can be found in the found in the ../ALVA/validation folder.

Case 1 - Basic validation at critical points in the pavemenet system (benchmark results)

In this section a basic validation of the ALVA package is presented, comparing ALVA repsonses to the Advanced Models for Analytical Design of European Pavement Structures (AMADEUS, 2000) report for benchmarking. As part of the AMADEUS study a number of popular response models were compared against the standard pavement structure shown in Table 1. Three packages, BISAR, KENLAYER and GAMES, that can facilitate layer interface slip situations, were further compared against ALVA for the proposed reference pavement system.

'ALVA_bonding_validation1.m' - tests the implementation of the interface spring model, as well as compares the ALVA model with a range of commonly used software at critical positions within the pavament system published by the European Commission.

Layer Thickness (mm) Young's moduli (MPa) Poisson's ratio
1 260 5000 0.35
2 500 200 0.40
3 50 0.45

Table 1: Reference pavement system used in basic validation of the code.

Load Radius (mm) Pressure, q (MPa) x-position (mm) y-position (mm)
Single 150.8 0.7 0 0
Dual (1) 106.6 0.7 -170 0
Dual (2) 106.6 0.7 170 0

Table 2: Load cases used for basic validation of the code.

Name Response Location Unit
R1 Vertical stress at the surface Load center MPa
R2 Horizontal strain at the bottom of layer 1 '' micron
R3 Vertical strain at the top of layer 2 '' micron
R4 Vertical strain at the top of layer 3 '' micron
------------- -------------------------------------------- ---------------------------- --------
R5 Vertical stress at the surface Edge of load MPa
R6 Horizontal strain at the bottom of layer 1 '' micron
R7 Vertical strain at the top of layer 2 '' micron
R8 Vertical strain at the top of layer 3 '' micron

Table 3: Description of key-point responses used for basic validation of the code.

Software R1 R2 R3 R4 R5 R6 R7 R8
BISAR 0.7 -100.5 251.7 185 0.4 -61.9 192.2 177.5
KENLAYER 0.8 -100.5 251.6 185.3 0.3 -62 192.2 177
GAMES 0.7 -100.5 251.6 185.1 0.3 -61.9 192.2 177.5
ALVA ('Slip') 0.7 -100.4 251.6 185.1 0.3 -61.8 192.3 177.5
ALVA ('Bonded') 0.7 -100.4 251.6 185.1 0.3 -61.8 192.3 177.5

Table 4: Reference pavement system subjeceted to a single wheel load - bonded interfaces.

Software R1 R2 R3 R4 R5 R6 R7 R8
BISAR 0.7 -120 1 217 0.4 -78 -10 205
KENLAYER 0.7 -120 1 216 0 -78 -10 205
GAMES 0.7 -119 11 217 0 -77 -1 205
ALVA ('Slip') 0.7 -119.5 0.7 216.5 0.3 -77.8 -9.7 204.7

Table 5: Reference pavement system subjeceted to a single wheel load - unbonded interface between layer 1 and 2, bonded interface between layer 2 and 3.

Software R1 R2 R3 R4 R5 R6 R7 R8
BISAR 0.7 N/A 186 170 0 N/A 182 177
KENLAYER 1.5 -85 186 170 0 -89 183 177
GAMES 0.7 -85 186 170 0 -89 183 177
ALVA ('Slip') 0.7 -84.9 186 169.7 0 -88.8 183.2 177.2
ALVA ('Bonded') 0.7 -84.9 186 169.7 0 -88.8 183.2 177.2

Table 6: Reference pavement system subjeceted to a dual wheel load - bonded interfaces.

Software R1 R2 R3 R4 R5 R6 R7 R8
BISAR 0.7 N/A 9 193 0 N/A -12 204
KENLAYER 0.7 -120 -1 216 0 -78 -10 205
GAMES 0.7 -101 -3 194 0 -106 1 205
ALVA ('Slip') 0.7 -102.6 -9.2 193.4 0 -107 -11.5 204.1

Table 7: Reference pavement system subjeceted to a dual wheel load - unbonded interface between layer 1 and 2, bonded interface between layer 2 and 3.

Case 2 - Vertical stresses and displacements for an elastic half-space with depth

Another basic validation step is carried by comparing vertical stresses and displacements in the line of loading with depth to the classical Boussinesq solution.

ALVA_let_validation1.m - tests the implementation of the ALVA LET model calculating vertical stresses and displacements with depth for a half-space subjected to a single circular load. The results are compared to the analytical Boussinesq solution and the computer programme ELLEA1.

Note: Minimum two layers is required for analysis of pavement systems in ALVA. For analysis of one-layer / half-space systems - select identical parameters for each layer, as well as "bonded" interface conditions

Figure 7: Validation: Half-space model.

Case 3 - Vertical and horizontal stresses and displacements with depth

ALVA_let_validation2.m - tests the implementation of the ALVA LET model, calculating stresses and displacements with depth for a multilayered pavement subjected to two circular loads utilizing the method proposed. The results are compared to the computer programme ELLEA1.

Figure 8: Validation: Multilayered model subjected to dual wheel load.

Case 4 - Surface displacements

In situ evaluation of mechanical pavement properties requires fitting measured surface displacements with model displacements. Such inverse analysis is guided by optimisation algorithms that entail re-execution of the underlying model many times over. For layered elasticity, which is the most commonly employed pavement model, this involves calculating computationally demanding semi-infinite integrals in every optimisation step. In this connection, a method was proposed (and implemented into ALVA) to improve the computational efficiency of surface displacement re-calculations in layered elastic systems.

ALVA_let_validation3.m - tests the implementation of the acclerated ALVA LET model, calculating the surface displacements with length for a multilayered pavement subjected to two circular loads utilizing the method proposed in Andersen et al. (2020). The results are compared to the computer programme ELLEA1.

Figure 9: Validation: Acceleration method for efficient calculation of surface displacements.

Case 5 - Shear stresses

ALVA_let_validation4.m - tests the implementation of the ALVA LET model calculating shear stresses with depth for a half-space subjected to a single circular load. The results are compared to the computer programme ELLEA1.

Figure 10: Validation: Shear stresses in multilayered model subjected to single wheel load.

Case 6 - Shear strains

ALVA_let_validation5.m - tests the implementation of the ALVA LET model calculating shear strain with depth for a multilayered pavement subjected to two circular loads. The results are compared to the computer programme ELLEA1.

Figure 11: Validation: Strains in multilayered model subjected to a dual wheel load.

Case 7 - Viscoelastic model (time dependent material behaviour)

  • ALVA_visco_validation1.m - tests the implementation of the ALVA VE model, calculating the displacements for a single evaluation point on the surface of a multilayered pavement considering a single circular load. The results are compared to the computer programme ELLVA1, see details in Levenberg (2016).

Figure 12: Validation: Surface displacements for a single evaluation point subjected to a moving load.

User examples

In this section we present relevant user examples with ALVA.

Example 1 - Backcalculation of layer moduli

ALVA_let_backcalculation.m - tests the implementation of the ALVA LET model for inferring layer moduli, or so called "backcalculation" of layer moduli, based on Falling Weight Deflectometer (FWD) measurements. The measurements was carried at DTU Lyngby Campus, Denmark Summer 2019.

The support script inv_loop.m was developed for this specific example. Moreover, for this optimisation problem, the fminsearch.m function in MATLAB was utilized, i.e., the Nelder–Mead simplex (NMS) multidimensional optimisation algorithm.

Input

Figure 13: Principle of operation for the Falling Weight Deflectometer (FWD) test: load and measurement configuration

Layer Thickness (mm) Young's moduli (MPa)* Poisson's ratio
1 161 200 0.35
2 260 200 0.35
3 200 0.35

Table 8: Measured pavement ticknesses. Note: * initial Young's moduli used in analysis.

Load Radius (mm) Pressure, q (MPa) x-position (mm) y-position (mm)
Single 150 0.7 0 0

Table 9: Falling Weight Deflectometer (FWD) load.

Geophone G1 G2 G3 G4 G5 G6 G7 G8 G9
Location along x-axis (mm): 0 200 300 450 600 900 1200 1500 1800
Displacement measurement (μm): 298.9 244.2 220.07 175.0 138.0 97.3 72.2 57.4 47.6

Table 10: Sensor location and displacement measurements.

ALVA 'Bonded' ALVA 'Frictionless'

Figure 14: Predicted versus measured pavement displacements.

Layer Thickness (mm) Young's moduli (MPa) - 'Bonded' Young's moduli (MPa) - 'Frictionless' Poisson's ratio
1 161 4885 2828 0.35
2 260 314 2392 0.35
3 179 162 0.35

Table 11: Backcalculated Young's moduli.

Installation

  • Download the package on your PC.
  • Open MATLAB
  • Go to the directory 'ALVA'
  • add the different directories of the ALVA on your MATLAB path — Now you are ready to run the validation examples provided and generate your own analysis.
  • ALVA is compatible with OCTAVE

For example validation of the library can be launched with

addpath("basic")
addpath("examples")
addpath("validation")
ALVA_let_validation1

How to contribute

  • To make changes or add a new function: (i) For the repository (make your own separate copy), (ii) make changes, and (iii) open a 'pull request'. Once approved, it can be merged into the master branch. If you wish to chat beforehand about your contribution, open an issue or email to [email protected].
  • If you find a bug in the code: open an 'issue' to notify contributors and create an official record.

Before contributing, please consider how your function fits into ALVA. At a minimum, functions must be well-documented and compatible with OCTAVE, not using any third party components.

References

Skar, Asmus, Andersen, Sebastian and Julius Nielsen (2020): Adaptive Layered Viscoelastic Analysis (ALVA). Technical University of Denmark. Software. https://doi.org/10.11583/DTU.12387305

Skar et al., (2020). ALVA: An adaptive MATLAB package for layered viscoelastic analysis. Journal of Open Source Software, 5(55), 2548. https://doi.org/10.21105/joss.02548

About

ALVA is a MATLAB-package for advanced pavement modeling that is highly adaptive

Resources

License

Stars

Watchers

Forks

Packages

No packages published

Contributors 2

  •  
  •