![]() |
This module contains functions and subroutines to compute tsource functions's time history.
Public Member Functions | |
real function, public | gabor (t, ts, l, a) |
function \(N^{\circ}3\) to compute real part of Gabor wavelet : \( f(t) = a \times \exp \left[ \frac{t-ts}{2\sigma} \right]^{2} \times \cos \left[ \frac{2 \pi (t-ts)}{\lambda} \right] \) More... | |
real function, public | expcos (t, ts, fp) |
function \(N^{\circ}4\) to compute euroseistest project source function for case can1 : \( f(t) = \exp \left[ - \frac{\omega(t-ts)}{\gamma} \right]^{2} \times \cos[\omega(t-ts)+\theta] \) More... | |
subroutine, public | ispli3 (t, du, s) |
function \(N^{\circ}5\) to compute order 3 spline : see Wikipedia More... | |
real function, public | ricker (t, ts, tp, a) |
function \(N^{\circ}6\) to compute order 2 Ricker wavelet : \( f(t) = 2 a \left[ \pi^2 \frac{(t-ts)^2}{tp^2} -0.5 \right] \times \exp \left[ -\pi^2 \frac{(t-ts)^2}{tp^2} \right] \) More... | |
real function, public | spiexp (t, tp, a) |
function \(N^{\circ}7\) to compute spice project exponential source function : \( f(t) = a \left[1-\left(1+\frac{t}{tp}\right) \right] \times \exp \left[ -\frac{t}{tp} \right] \) More... | |
real function, public | fctanh (t, ts, tp, a) |
function \(N^{\circ}8\) to compute hyperbolic tangent function (also called 'Bouchon pulse') : \( f(t) = 0.5 a \left[ 1 + \tanh \left( \frac{t-ts}{2/5 tp} \right) \right] \) More... | |
real function, public | fctanh_dt (t, ts, tp, a) |
function \(N^{\circ}9\) to compute first derivative of hyperbolic tangent function : \( f(t) = 2 a \times f_c \left[ 1 - \tanh^{2} \left( 4 f_c (t-ts) \right) \right] \) with \(f_c = \frac{1}{8/5 tp}\) More... | |
subroutine, public | interp_linear (dim_num, data_num, t_data, p_data, interp_num, t_interp, p_interp) |
INTERP_LINEAR applies piecewise linear interpolation to data. More... | |
Private Member Functions | |
logical function, private | rvec_ascends_strictly (n, x) |
This function determines if an R8VEC is strictly ascending. More... | |
subroutine, private | rvec_bracket (n, x, xval, left, right) |
This subroutine searches a sorted R8VEC for successive brackets of a value. More... | |
Definition at line 126 of file module_source_function.f90.
real function, public mod_source_function::expcos | ( | real, intent(in) | t, |
real, intent(in) | ts, | ||
real, intent(in) | fp | ||
) |
t | : current time of computation \( t_{simu} \) |
ts | : time shift of the function |
fp | : pseudo-frequency of the function |
Definition at line 186 of file module_source_function.f90.
Referenced by mod_solver::compute_external_force().
real function, public mod_source_function::fctanh | ( | real, intent(in) | t, |
real, intent(in) | ts, | ||
real, intent(in) | tp, | ||
real, intent(in) | a | ||
) |
t | : current time of computation \( t_{simu} \) |
ts | : time shift of the function |
tp | : pseudo-period of the function |
a | : amplitude factor |
Definition at line 303 of file module_source_function.f90.
Referenced by mod_solver::compute_external_force().
real function, public mod_source_function::fctanh_dt | ( | real, intent(in) | t, |
real, intent(in) | ts, | ||
real, intent(in) | tp, | ||
real, intent(in) | a | ||
) |
t | : current time of computation \( t_{simu} \) |
ts | : time shift of the function |
tp | : pseudo-period of the function |
a | : amplitude factor |
Definition at line 333 of file module_source_function.f90.
Referenced by mod_solver::compute_external_force().
real function, public mod_source_function::gabor | ( | real, intent(in) | t, |
real, intent(in) | ts, | ||
real, intent(in) | l, | ||
real, intent(in) | a | ||
) |
t | : current time of computation \( t_{simu} \) |
ts | : time shift of the function |
l | : parameter \(\lambda\). By default: \(\sigma = \lambda/20\). |
a | : amplitude |
Definition at line 157 of file module_source_function.f90.
Referenced by mod_solver::compute_external_force().
subroutine, public mod_source_function::interp_linear | ( | integer, intent(in) | dim_num, |
integer, intent(in) | data_num, | ||
real, dimension(data_num), intent(in) | t_data, | ||
real, dimension(dim_num,data_num), intent(in) | p_data, | ||
integer, intent(in) | interp_num, | ||
real, dimension(interp_num), intent(in) | t_interp, | ||
real, dimension(dim_num,interp_num), intent(out) | p_interp | ||
) |
From a space of DIM_NUM dimensions, we are given a sequence of DATA_NUM points, which are presumed to be successive samples from a curve of points P.
We are also given a parameterization of this data, that is, an associated sequence of DATA_NUM values of a variable T. The values of T are assumed to be strictly increasing.
Thus, we have a sequence of values P(T), where T is a scalar, and each value of P is of dimension DIM_NUM.
We are then given INTERP_NUM values of T, for which values P are to be produced, by linear interpolation of the data we are given.
Note that the user may request extrapolation. This occurs whenever a T_INTERP value is less than the minimum T_DATA or greater than the maximum T_DATA. In that case, linear extrapolation is used.
Licensing: This code is distributed under the GNU LGPL license.
dim_num | : the spatial dimension |
data_num | : the number of data points. |
t_data(data_num) | : values of the independent variable at the sample points. Values of T_DATA must be strictly increasing. |
p_data(dim_num,data_num) | : values of the dependent variables at the sample points. |
interp_num | : number of points at which interpolation is to be done. |
t_interp(interp_num) | : values of the independent variable at the interpolation points. |
p_interp(dim_num,data_num) | : interpolated values of the dependent variables at the interpolation points. |
Definition at line 394 of file module_source_function.f90.
References rvec_ascends_strictly(), and rvec_bracket().
Referenced by mod_source::init_double_couple_source(), and mod_source::init_single_force_source().
subroutine, public mod_source_function::ispli3 | ( | real, intent(in) | t, |
real, intent(in) | du, | ||
real, intent(inout) | s | ||
) |
s | : spline value |
t | : current time of computation \( t_{simu} \) |
du | : source duration |
Definition at line 217 of file module_source_function.f90.
Referenced by mod_solver::compute_external_force().
real function, public mod_source_function::ricker | ( | real, intent(in) | t, |
real, intent(in) | ts, | ||
real, intent(in) | tp, | ||
real, intent(in) | a | ||
) |
t | : current time of computation \( t_{simu} \) |
ts | : time shift of the function |
tp | : pseudo-period of the function |
a | : amplitude factor |
Definition at line 255 of file module_source_function.f90.
Referenced by mod_solver::compute_external_force().
|
private |
n | : size of the array |
x | : array to be examined |
Definition at line 455 of file module_source_function.f90.
Referenced by interp_linear().
|
private |
n | : size of input array |
x | : an array sorted into ascending order |
xval | : a value to be bracketed |
left | : x(left) <= xval <= x(right). xval < x(1), when left = 1, right = 2. x(n) < xval, when left = n-1, right = n |
right | : x(left) <= xval <= x(right). xval < x(1), when left = 1, right = 2. x(n) < xval, when left = n-1, right = n |
Definition at line 493 of file module_source_function.f90.
Referenced by interp_linear().
real function, public mod_source_function::spiexp | ( | real, intent(in) | t, |
real, intent(in) | tp, | ||
real, intent(in) | a | ||
) |
t | : current time of computation \( t_{simu} \) |
tp | : pseudo-period of the function |
a | : amplitude factor |
Definition at line 280 of file module_source_function.f90.
Referenced by mod_solver::compute_external_force().