ARTS 2.5.11 (git: 6827797f)
|
#include "math_funcs.h"
#include <cmath>
#include <iostream>
#include <stdexcept>
#include "array.h"
#include "arts_constants.h"
#include "arts_conversions.h"
#include "logic.h"
#include "mystring.h"
Go to the source code of this file.
Functions | |
Numeric | fac (const Index n) |
fac | |
Index | integer_div (const Index &x, const Index &y) |
integer_div | |
Numeric | LagrangeInterpol4 (ConstVectorView x, ConstVectorView y, const Numeric a) |
Lagrange Interpolation (internal function). | |
Numeric | last (ConstVectorView x) |
last | |
Index | last (const ArrayOfIndex &x) |
last | |
void | linspace (Vector &x, const Numeric start, const Numeric stop, const Numeric step) |
linspace | |
void | nlinspace (Vector &x, const Numeric start, const Numeric stop, const Index n) |
nlinspace | |
void | nlinspace (VectorView x, const Numeric start, const Numeric stop, const Index n) |
void | nlogspace (Vector &x, const Numeric start, const Numeric stop, const Index n) |
nlogspace | |
Numeric | trapz (ConstVectorView x, ConstVectorView y) |
trapz | |
void | cumsum (VectorView csum, ConstVectorView x) |
cumsum | |
Numeric | AngIntegrate_trapezoid (ConstMatrixView Integrand, ConstVectorView za_grid, ConstVectorView aa_grid) |
AngIntegrate_trapezoid. | |
Numeric | AngIntegrate_trapezoid_opti (ConstMatrixView Integrand, ConstVectorView za_grid, ConstVectorView aa_grid, ConstVectorView grid_stepsize) |
AngIntegrate_trapezoid_opti. | |
Numeric | AngIntegrate_trapezoid (ConstVectorView Integrand, ConstVectorView za_grid) |
AngIntegrate_trapezoid. | |
Numeric | sign (const Numeric &x) |
sign | |
void | mgd (VectorView psd, const Vector &x, const Numeric &n0, const Numeric &mu, const Numeric &la, const Numeric &ga) |
void | mgd_with_derivatives (VectorView psd, MatrixView jac_data, const Vector &x, const Numeric &n0, const Numeric &mu, const Numeric &la, const Numeric &ga, const bool &do_n0_jac, const bool &do_mu_jac, const bool &do_la_jac, const bool &do_ga_jac) |
void | delanoe_shape_with_derivative (VectorView psd, MatrixView jac_data, const Vector &x, const Numeric &alpha, const Numeric &beta) |
! Shape functions for normalized PSD. | |
Numeric | mod_gamma_dist (Numeric x, Numeric N0, Numeric Lambda, Numeric mu, Numeric gamma) |
Generalized Modified Gamma Distribution. | |
void | unitl (Vector &x) |
unitl | |
void | flat (VectorView x, ConstMatrixView X) |
flat | |
void | flat (VectorView x, ConstTensor3View X) |
flat | |
void | reshape (Tensor3View X, ConstVectorView x) |
reshape | |
void | reshape (MatrixView X, ConstVectorView x) |
reshape | |
void | calculate_weights_linear (Vector &x, Vector &w, const Index nph) |
calculate_weights_linear | |
void | calculate_int_weights_arbitrary_grid (Vector &w, const Vector &x) |
Calculates trapezoidal integration weights for arbitray grid. | |
Variables | |
constexpr Numeric | DEG2RAD =Conversion::deg2rad(1) |
constexpr Numeric | PI =Constant::pi |
Numeric AngIntegrate_trapezoid | ( | ConstMatrixView | Integrand, |
ConstVectorView | za_grid, | ||
ConstVectorView | aa_grid | ||
) |
AngIntegrate_trapezoid.
Performs an integration of a matrix over all directions defined in angular grids using the trapezoidal integration method.
Integrand | The Matrix to be integrated |
za_grid | The zenith angle grid |
aa_grid | The azimuth angle grid |
Definition at line 318 of file math_funcs.cc.
References ARTS_ASSERT, and DEG2RAD.
Referenced by AngIntegrate_trapezoid_opti(), doit_scat_fieldCalc(), doit_scat_fieldCalcLimb(), doit_scat_fieldNormalize(), get_stepwise_scattersky_source(), scat_dataCheck(), and scat_dataReduceT().
Numeric AngIntegrate_trapezoid | ( | ConstVectorView | Integrand, |
ConstVectorView | za_grid | ||
) |
AngIntegrate_trapezoid.
Performs an integration of a matrix over all directions defined in angular grids using the trapezoidal integration method. The integrand is independant of the azimuth angle. The integration over the azimuth angle gives a 2*PI
Integrand | Input : The vector to be integrated |
za_grid | Input : The zenith angle grid |
Definition at line 415 of file math_funcs.cc.
References ARTS_ASSERT, DEG2RAD, and PI.
Numeric AngIntegrate_trapezoid_opti | ( | ConstMatrixView | Integrand, |
ConstVectorView | za_grid, | ||
ConstVectorView | aa_grid, | ||
ConstVectorView | grid_stepsize | ||
) |
AngIntegrate_trapezoid_opti.
Performs an integration of a matrix over all directions defined in angular grids using the trapezoidal integration method.
In addition to the "old fashined" integration method, it checks whether the stepsize is constant. If it is, it uses a faster method, if not, it uses the old one.
Integrand | Input : The Matrix to be integrated |
za_grid | Input : The zenith angle grid |
aa_grid | Input : The azimuth angle grid |
grid_stepsize | Input : stepsize of the grid |
Definition at line 362 of file math_funcs.cc.
References AngIntegrate_trapezoid(), ARTS_ASSERT, and DEG2RAD.
Referenced by doit_scat_fieldCalc(), and doit_scat_fieldCalcLimb().
void calculate_int_weights_arbitrary_grid | ( | Vector & | w, |
const Vector & | x | ||
) |
Calculates trapezoidal integration weights for arbitray grid.
w | Vector integration weights |
x | Vector evaluation points |
Definition at line 833 of file math_funcs.cc.
References ARTS_USER_ERROR_IF, and w.
Referenced by iySurfaceLambertian().
void calculate_weights_linear | ( | Vector & | x, |
Vector & | w, | ||
const Index | nph | ||
) |
calculate_weights_linear
Function to set the evaluation points and the corresponding weights for numerical integration on the domain from [-1,1]
[out] | x | evaluation points \paran[out] w integration weights |
[in] | nph | number of evaluation points per hemisphere |
Definition at line 815 of file math_funcs.cc.
References nlinspace(), and w.
Referenced by AngularGridsSetFluxCalc().
void cumsum | ( | VectorView | csum, |
ConstVectorView | x | ||
) |
cumsum
Cumulative sum of a vector
csum[0] becomes x[0] and last value in csum equals the full sum of x
csum | Vector with cumulative values |
x | Input vector |
Definition at line 297 of file math_funcs.cc.
References ARTS_ASSERT.
Referenced by dlosGauss().
void delanoe_shape_with_derivative | ( | VectorView | psd, |
MatrixView | jac_data, | ||
const Vector & | x, | ||
const Numeric & | alpha, | ||
const Numeric & | beta | ||
) |
! Shape functions for normalized PSD.
This function implements the shape function F(X, alpha, beta) from as proposed by Delanoe et al. in "Normalized particle size distribution for remote sensing application".
[OUT] | psd On return contains the values of F corresponding to the values in x. |
[OUT] | jac_data On return contains the first derivative of F w.r.t x evaluated at the values in x. |
[IN] | x The values at which to evaluate the shape functions and and derivatives. |
[IN] | alpha The alpha parameter of the shape function. |
[IN] | beta The beta parameter of the shape function. |
Definition at line 623 of file math_funcs.cc.
Referenced by psdDelanoeEtAl14().
Numeric fac | ( | const Index | n | ) |
fac
Calculates the factorial.
The function asserts that n must be >= 0
n | Nominator |
Definition at line 46 of file math_funcs.cc.
Referenced by LineShape::RosenkranzQuadratic::dNdf(), MatrixGaussian(), MCRadar(), Absorption::LineMixing::SpeciesErrorCorrectedSuddenData::Omega(), LineShape::RosenkranzQuadratic::operator()(), specular_losCalc(), VectorGaussian(), and ze_cfac().
void flat | ( | VectorView | x, |
ConstMatrixView | X | ||
) |
flat
Flattens a matrix to a vector
The matrix is read from front, i.e. rows are looped first. In Matlab this equals x=X(:).
[out] | x | The vector. Should already be sized |
[in] | X | The matrix. |
Definition at line 709 of file math_funcs.cc.
References ARTS_ASSERT, and c.
void flat | ( | VectorView | x, |
ConstTensor3View | X | ||
) |
flat
Converts Tensor3 to a vector
The matrix is read from front, i.e. pages are looped first, followed by rows. In Matlab this equals x=X(:).
[out] | x | The vector. Should already be sized |
[in] | X | The tensor. |
Definition at line 735 of file math_funcs.cc.
References ARTS_ASSERT, and c.
Index integer_div | ( | const Index & | x, |
const Index & | y | ||
) |
integer_div
Performs an integer division.
The function asserts that the reminder of the division x/y is 0.
x | Nominator |
y | Denominator |
Definition at line 70 of file math_funcs.cc.
References ARTS_ASSERT.
Numeric LagrangeInterpol4 | ( | ConstVectorView | x, |
ConstVectorView | y, | ||
const Numeric | a | ||
) |
Lagrange Interpolation (internal function).
This function calculates the Lagrange interpolation of four interpolation points as described in Lagrange Interpolating Polynomial.
The input are the four x-axis values [x0,x1,x2,x3] and their associated y-axis values [y0,y1,y2,y3]. Furthermore the x-axis point "a" at which the interpolation should be calculated must be given as input. NOTE that the relation x2 =< x < x3 MUST hold!
x | x-vector with four elements [x0,x1,x2,x3] |
y | y-vector with four elements: yj = y(xj), j=0,1,2,3 |
a | interpolation point on the x-axis with x1 =< a < x2 |
Definition at line 96 of file math_funcs.cc.
References a, ARTS_USER_ERROR_IF, and b.
Index last | ( | const ArrayOfIndex & | x | ) |
last
Returns the last value of an index array.
x | An index array. |
Definition at line 157 of file math_funcs.cc.
References ARTS_ASSERT, and Array< base >::nelem().
Numeric last | ( | ConstVectorView | x | ) |
last
Returns the last value of a vector.
x | A vector. |
Definition at line 142 of file math_funcs.cc.
References ARTS_ASSERT.
Referenced by chk_scat_data(), ext_abs_pfun_from_tro(), find_effective_channel_boundaries(), Raw::Reduce::find_first_and_last_1(), mixer_matrix(), particle_bulkpropRadarOnionPeeling(), pos2refell_r(), refr_index_airFreeElectrons(), scat_data_singleTmatrix(), ScatSpeciesExtendTemperature(), sensor_responseMixer(), summation_by_vecmult(), VectorNLinSpaceVector(), ybatchDoublingMeanFocus(), yDoublingMeanFocus(), and ySimpleSpectrometer().
void linspace | ( | Vector & | x, |
const Numeric | start, | ||
const Numeric | stop, | ||
const Numeric | step | ||
) |
linspace
Linearly spaced vector with specified spacing.
The first element of x is always start. The next value is start+step etc. Note that the last value can deviate from stop. The step can be both positive and negative. (in Matlab notation: start:step:stop)
Size of result is adjusted within this function!
x | Output: linearly spaced vector |
start | first value in x |
stop | last value of x <= stop |
step | distance between values in x |
Definition at line 181 of file math_funcs.cc.
void mgd | ( | VectorView | psd, |
const Vector & | x, | ||
const Numeric & | n0, | ||
const Numeric & | mu, | ||
const Numeric & | la, | ||
const Numeric & | ga | ||
) |
Modified gamma distribution
Uses all four free parameters (n0, mu, la, ga) to calculate psd(D) = n0 * D^mu * exp( -la * x^ga )
Reference: Eq 1 of Petty & Huang, JAS, (2011).
psd | Particle number density per x-interval. Sizing of vector should be done before calling the function. |
x | Mass or size. |
n0 | See above. |
mu | See above. |
la | See above. |
ga | See above. |
Definition at line 473 of file math_funcs.cc.
References ARTS_ASSERT, and ARTS_USER_ERROR_IF.
void mgd_with_derivatives | ( | VectorView | psd, |
MatrixView | jac_data, | ||
const Vector & | x, | ||
const Numeric & | n0, | ||
const Numeric & | mu, | ||
const Numeric & | la, | ||
const Numeric & | ga, | ||
const bool & | do_n0_jac, | ||
const bool & | do_mu_jac, | ||
const bool & | do_la_jac, | ||
const bool & | do_ga_jac | ||
) |
Modified gamma distribution, and derivatives
As mgd, but this version can also return the derivate of psd with respect to the four parameters.
psd | Particle number density per x-interval. Sizing of vector should be done before calling the function. |
jac_data | Container for returning jacobian data. Shall be a matrix with four rows, where the rows match n0, mu, la and ga. Number of columns same as length of psd. |
x | Mass or size. |
n0 | See above. |
mu | See above. |
la | See above. |
ga | See above. |
do_n0_jac | Flag to actually calculate d_psd/d_n0 |
do_mu_jac | Flag to actually calculate d_psd/d_mu |
do_la_jac | Flag to actually calculate d_psd/d_la |
do_ga_jac | Flag to actually calculate d_psd/d_ga |
Definition at line 543 of file math_funcs.cc.
References ARTS_ASSERT, and ARTS_USER_ERROR_IF.
Referenced by psd_mgd_mass_and_something(), psd_mgd_smm_common(), psdModifiedGamma(), and psdModifiedGammaMass().
Numeric mod_gamma_dist | ( | Numeric | x, |
Numeric | N0, | ||
Numeric | Lambda, | ||
Numeric | mu, | ||
Numeric | gamma | ||
) |
Generalized Modified Gamma Distribution.
Returns number density per unit of 'x' as function of 'x'.
x | Numeric |
N0 | Numeric, Scaling parameter |
Lambda | Numeric, Shape parameter |
mu | Numeric, Shape parameter |
gamma | Numeric, Shape parameter |
Definition at line 657 of file math_funcs.cc.
References ARTS_USER_ERROR.
Referenced by psd_MY05(), and psd_SB06().
void nlinspace | ( | Vector & | x, |
const Numeric | start, | ||
const Numeric | stop, | ||
const Index | n | ||
) |
nlinspace
Linearly spaced vector with specified length.
Returns a vector equally and linearly spaced between start and stop of length n. (equals the Matlab function linspace)
The length must be > 1.
x | Output: linearly spaced vector |
start | first value in x |
stop | last value of x <= stop |
n | length of x |
Definition at line 208 of file math_funcs.cc.
References ARTS_ASSERT.
Referenced by AngularGridsSetFluxCalc(), calc_ssp_fixed_test(), calc_ssp_random_test(), calculate_weights_linear(), dlosGauss(), DOAngularGridsSet(), doit_scat_fieldCalcLimb(), DoitScatteringDataPrepare(), get_angs(), HydrotableCalc(), iyEmissionHybrid(), line_irradianceCalcForSingleSpeciesNonOverlappingLinesPseudo2D(), pha_mat_sptFromDataDOITOpt(), pha_mat_sptFromMonoData(), pnd_fieldZero(), run_cdisort(), run_cdisort_flux(), sensor_responseFillFgrid(), VectorNLinSpace(), and WriteBuiltinPartitionFunctionsXML().
void nlinspace | ( | VectorView | x, |
const Numeric | start, | ||
const Numeric | stop, | ||
const Index | n | ||
) |
Definition at line 218 of file math_funcs.cc.
void nlogspace | ( | Vector & | x, |
const Numeric | start, | ||
const Numeric | stop, | ||
const Index | n | ||
) |
nlogspace
Logarithmically spaced vector with specified length.
Returns a vector logarithmically spaced vector between start and stop of length n (equals the Matlab function logspace)
The length must be > 1.
x | Output: logarithmically spaced vector |
start | first value in x |
stop | last value of x <= stop |
n | length of x |
Definition at line 244 of file math_funcs.cc.
References a, and ARTS_ASSERT.
void reshape | ( | MatrixView | X, |
ConstVectorView | x | ||
) |
reshape
Converts vector to Matrix
The matrix is filled from front, i.e. rows are looped first, followed by cols. In Matlab this equals X = reshape( x, [ X.nrows(), X.ncols() ]
[out] | X | The matrix. Should already be sized |
[in] | x | The vector. |
Definition at line 791 of file math_funcs.cc.
References ARTS_ASSERT, and c.
void reshape | ( | Tensor3View | X, |
ConstVectorView | x | ||
) |
reshape
Converts vector to Tensor3
The tensor is filled from front, i.e. pages are looped first, followed by rows. In Matlab this equals X = reshape( x, [ X.npages(), X.nrows(), X.ncols() ]
[out] | X | The tensor. Should already be sized |
[in] | x | The vector. |
Definition at line 763 of file math_funcs.cc.
References ARTS_ASSERT, and c.
Referenced by AtmFieldPerturb(), and x2artsAtmAndSurf().
Numeric sign | ( | const Numeric & | x | ) |
sign
Returns the sign of a numeric value.
The function returns 1 if the value is greater than zero, 0 if it equals zero and -1 if it is less than zero.
x | A Numeric. |
Definition at line 445 of file math_funcs.cc.
Referenced by error_if_limb_ppath(), geodeticposlos2cart(), geompath_from_r1_to_r2(), poslos2cart(), r_crossing_2d(), raytrace_2d_linear_basic(), specular_losCalc(), and specular_losCalcNoTopography().
Numeric trapz | ( | ConstVectorView | x, |
ConstVectorView | y | ||
) |
trapz
Integration by the basic trapezoidal rule
x | Grid values |
y | Integrand |
Definition at line 274 of file math_funcs.cc.
References ARTS_ASSERT.
Referenced by asymmetry_parameter().
void unitl | ( | Vector & | x | ) |
unitl
Normalises a vector to have unit length.
The standard Euclidean norm is used (2-norm).
param x In/Out: A vector.
Definition at line 689 of file math_funcs.cc.
References ARTS_ASSERT.
|
inlineconstexpr |
Definition at line 27 of file math_funcs.cc.
Referenced by AngIntegrate_trapezoid(), and AngIntegrate_trapezoid_opti().
|
inlineconstexpr |
Definition at line 28 of file math_funcs.cc.
Referenced by AngIntegrate_trapezoid().