Utility functions and classes. More...

Files

file  OFELI.h
 Header file that includes all kernel classes of the library.
 
file  OFELI_Config.h
 File that contains some macros.
 
file  constants.h
 File that contains some widely used constants.
 

Classes

class  Funct
 A simple class to parse real valued functions. More...
 
class  Tabulation
 To read and manipulate tabulated functions. More...
 
class  Point< T_ >
 Defines a point with arbitrary type coordinates. More...
 
class  Point2D< T_ >
 Defines a 2-D point with arbitrary type coordinates. More...
 
class  SpaceTime
 Defines a space-time point. More...
 
class  Gauss
 Calculate data for Gauss integration. More...
 
class  Timer
 To handle elapsed time counting. More...
 

Macros

#define OFELI_E   2.71828182845904523536028747135
 
#define OFELI_PI   3.14159265358979323846264338328
 
#define OFELI_THIRD   0.33333333333333333333333333333
 
#define OFELI_SIXTH   0.16666666666666666666666666667
 
#define OFELI_TWELVETH   0.08333333333333333333333333333
 
#define OFELI_SQRT2   1.41421356237309504880168872421
 
#define OFELI_SQRT3   1.73205080756887729352744634151
 
#define OFELI_ONEOVERPI   0.31830988618379067153776752675
 
#define OFELI_GAUSS2   0.57735026918962576450914878050196
 
#define OFELI_EPSMCH   DBL_EPSILON
 
#define OFELI_TOLERANCE   OFELI_EPSMCH*10000
 
#define VLG   1.e10
 
#define OFELI_IMAG   std::complex<double>(0.,1.);
 
#define LIGHT_SPEED   299792458.2;
 
#define CATCH_EXCEPTION
 

Functions

ostream & operator<< (ostream &s, const complex_t &x)
 Output a complex number.
 
ostream & operator<< (ostream &s, const std::string &c)
 Output a string.
 
template<class T_ >
ostream & operator<< (ostream &s, const vector< T_ > &v)
 Output a vector instance.
 
template<class T_ >
ostream & operator<< (ostream &s, const std::list< T_ > &l)
 Output a vector instance.
 
template<class T_ >
ostream & operator<< (ostream &s, const std::pair< T_, T_ > &a)
 Output a pair instance.
 
void saveField (Vect< real_t > &v, string output_file, int opt)
 Save a vector to an output file in a given file format.
 
void saveField (const Vect< real_t > &v, const Mesh &mesh, string output_file, int opt)
 Save a vector to an output file in a given file format.
 
void saveField (Vect< real_t > &v, const Grid &g, string output_file, int opt)
 Save a vector to an output file in a given file format, for a structured grid data.
 
void saveGnuplot (string input_file, string output_file, string mesh_file, int f=1)
 Save a vector to an input Gnuplot file.
 
void saveGnuplot (Mesh &mesh, string input_file, string output_file, int f=1)
 Save a vector to an input Gnuplot file.
 
void saveTecplot (string input_file, string output_file, string mesh_file, int f=1)
 Save a vector to an output file to an input Tecplot file.
 
void saveTecplot (Mesh &mesh, string input_file, string output_file, int f=1)
 Save a vector to an output file to an input Tecplot file.
 
void saveVTK (string input_file, string output_file, string mesh_file, int f=1)
 Save a vector to an output VTK file.
 
void saveVTK (Mesh &mesh, string input_file, string output_file, int f=1)
 Save a vector to an output VTK file.
 
void saveGmsh (string input_file, string output_file, string mesh_file, int f=1)
 Save a vector to an output Gmsh file.
 
void saveGmsh (Mesh &mesh, string input_file, string output_file, int f=1)
 Save a vector to an output Gmsh file.
 
ostream & operator<< (ostream &s, const Tabulation &t)
 Output Tabulated function data.
 
template<class T_ >
bool operator== (const Point< T_ > &a, const Point< T_ > &b)
 Operator ==
 
template<class T_ >
Point< T_ > operator+ (const Point< T_ > &a, const Point< T_ > &b)
 Operator +
 
template<class T_ >
Point< T_ > operator+ (const Point< T_ > &a, const T_ &x)
 Operator +
 
template<class T_ >
Point< T_ > operator- (const Point< T_ > &a)
 Unary Operator -
 
template<class T_ >
Point< T_ > operator- (const Point< T_ > &a, const Point< T_ > &b)
 Operator -
 
template<class T_ >
Point< T_ > operator- (const Point< T_ > &a, const T_ &x)
 Operator -
 
template<class T_ >
Point< T_ > operator* (const T_ &a, const Point< T_ > &b)
 Operator *
 
template<class T_ >
Point< T_ > operator* (const int &a, const Point< T_ > &b)
 Operator *.
 
template<class T_ >
Point< T_ > operator* (const Point< T_ > &b, const T_ &a)
 Operator /
 
template<class T_ >
Point< T_ > operator* (const Point< T_ > &b, const int &a)
 Operator *
 
template<class T_ >
T_ operator* (const Point< T_ > &a, const Point< T_ > &b)
 Operator *
 
template<class T_ >
Point< T_ > operator/ (const Point< T_ > &b, const T_ &a)
 Operator /
 
bool areClose (const Point< real_t > &a, const Point< real_t > &b, real_t toler=OFELI_TOLERANCE)
 Return true if both instances of class Point<double> are distant with less then toler
 
real_t SqrDistance (const Point< real_t > &a, const Point< real_t > &b)
 Return squared euclidean distance between points a and b
 
real_t Distance (const Point< real_t > &a, const Point< real_t > &b)
 Return euclidean distance between points a and b
 
bool operator< (const Point< size_t > &a, const Point< size_t > &b)
 Comparison operator. Returns true if all components of first vector are lower than those of second one.
 
template<class T_ >
std::ostream & operator<< (std::ostream &s, const Point< T_ > &a)
 Output point coordinates.
 
template<class T_ >
bool operator== (const Point2D< T_ > &a, const Point2D< T_ > &b)
 Operator ==.
 
template<class T_ >
Point2D< T_ > operator+ (const Point2D< T_ > &a, const Point2D< T_ > &b)
 Operator +.
 
template<class T_ >
Point2D< T_ > operator+ (const Point2D< T_ > &a, const T_ &x)
 Operator +.
 
template<class T_ >
Point2D< T_ > operator- (const Point2D< T_ > &a)
 Unary Operator -
 
template<class T_ >
Point2D< T_ > operator- (const Point2D< T_ > &a, const Point2D< T_ > &b)
 Operator -
 
template<class T_ >
Point2D< T_ > operator- (const Point2D< T_ > &a, const T_ &x)
 Operator -
 
template<class T_ >
Point2D< T_ > operator* (const T_ &a, const Point2D< T_ > &b)
 Operator *.
 
template<class T_ >
Point2D< T_ > operator* (const int &a, const Point2D< T_ > &b)
 
template<class T_ >
Point2D< T_ > operator* (const Point2D< T_ > &b, const T_ &a)
 Operator /
 
template<class T_ >
Point2D< T_ > operator* (const Point2D< T_ > &b, const int &a)
 Operator *
 
template<class T_ >
T_ operator* (const Point2D< T_ > &b, const Point2D< T_ > &a)
 Operator *.
 
template<class T_ >
Point2D< T_ > operator/ (const Point2D< T_ > &b, const T_ &a)
 Operator /
 
bool areClose (const Point2D< real_t > &a, const Point2D< real_t > &b, real_t toler=OFELI_TOLERANCE)
 Return true if both instances of class Point2D<real_t> are distant with less then toler [Default: OFELI_EPSMCH].
 
real_t SqrDistance (const Point2D< real_t > &a, const Point2D< real_t > &b)
 Return squared euclidean distance between points a and b
 
real_t Distance (const Point2D< real_t > &a, const Point2D< real_t > &b)
 Return euclidean distance between points a and b
 
template<class T_ >
std::ostream & operator<< (std::ostream &s, const Point2D< T_ > &a)
 Output point coordinates.
 
bool operator== (const SpaceTime &a, const SpaceTime &b)
 Operator ==
 
SpaceTime operator+ (const SpaceTime &a, const SpaceTime &b)
 Operator +
 
SpaceTime operator+ (const SpaceTime &a, const real_t &x)
 Operator +
 
SpaceTime operator- (const SpaceTime &a)
 Unary Operator -
 
SpaceTime operator- (const SpaceTime &a, const SpaceTime &b)
 Operator -
 
SpaceTime operator- (const SpaceTime &a, const real_t &x)
 Operator -
 
SpaceTime operator* (const real_t &a, const SpaceTime &b)
 Operator *
 
SpaceTime operator* (const int &a, const SpaceTime &b)
 Operator *.
 
SpaceTime operator* (const SpaceTime &b, const real_t &a)
 Operator /
 
SpaceTime operator* (const SpaceTime &b, const int &a)
 Operator *
 
real_t operator* (const SpaceTime &a, const SpaceTime &b)
 Operator *
 
SpaceTime operator/ (const SpaceTime &b, const real_t &a)
 Operator /
 
bool areClose (const SpaceTime &a, const SpaceTime &b, real_t toler=OFELI_TOLERANCE)
 Return true if both instances of class Point<double> are distant with less then toler
 
real_t SqrDistance (const SpaceTime &a, const SpaceTime &b)
 Return squared euclidean distance between points a and b
 
real_t Distance (const SpaceTime &a, const SpaceTime &b)
 Return euclidean distance between points a and b
 
bool operator< (const SpaceTime &a, const SpaceTime &b)
 Comparison operator. Returns true if all components of first vector are lower than those of second one.
 
std::ostream & operator<< (std::ostream &s, const SpaceTime &a)
 Output space-time point coordinates.
 
real_t Discrepancy (Vect< real_t > &x, const Vect< real_t > &y, int n, int type=1)
 Return discrepancy between 2 vectors x and y
 
real_t Discrepancy (Vect< complex_t > &x, const Vect< complex_t > &y, int n, int type=1)
 Return discrepancy between 2 vectors x and y
 
void getMesh (string file, ExternalFileFormat form, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in an external file format.
 
void getBamg (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in Bamg format.
 
void getEasymesh (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in Easymesh format.
 
void getGambit (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in Gambit neutral format.
 
void getGmsh (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in Gmsh format.
 
void getMatlab (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a Matlab mesh data.
 
void getNetgen (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in Netgen format.
 
void getTetgen (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in Tetgen format.
 
void getTriangle (string file, Mesh &mesh, size_t nb_dof=1)
 Construct an instance of class Mesh from a mesh file stored in Triangle format.
 
void saveMesh (const string &file, const Mesh &mesh, ExternalFileFormat form)
 This function saves mesh data a file for a given external format.
 
void saveGmsh (const string &file, const Mesh &mesh)
 This function outputs a Mesh instance in a file in Gmsh format.
 
void saveGnuplot (const string &file, const Mesh &mesh)
 This function outputs a Mesh instance in a file in Gmsh format.
 
void saveMatlab (const string &file, const Mesh &mesh)
 This function outputs a Mesh instance in a file in Matlab format.
 
void saveTecplot (const string &file, const Mesh &mesh)
 This function outputs a Mesh instance in a file in Tecplot format.
 
void saveVTK (const string &file, const Mesh &mesh)
 This function outputs a Mesh instance in a file in VTK format.
 
void BSpline (size_t n, size_t t, Vect< Point< real_t > > &control, Vect< Point< real_t > > &output, size_t num_output)
 Function to perform a B-spline interpolation.
 
void banner (const string &prog=" ")
 Outputs a banner as header of any developed program.
 
template<class T_ >
void QuickSort (std::vector< T_ > &a, int begin, int end)
 Function to sort a vector.
 
template<class T_ >
void qksort (std::vector< T_ > &a, int begin, int end)
 Function to sort a vector.
 
template<class T_ , class C_ >
void qksort (std::vector< T_ > &a, int begin, int end, C_ compare)
 Function to sort a vector according to a key function.
 
int Sgn (real_t a)
 Return sign of a: -1 or 1.
 
real_t Sgn (real_t a, real_t b)
 Return |a| if b>0, -|a| if not.
 
real_t Abs2 (complex_t a)
 Return square of modulus of complex number a
 
real_t Abs2 (real_t a)
 Return square of real number a
 
real_t Abs (real_t a)
 Return absolute value of a
 
real_t Abs (complex_t a)
 Return modulus of complex number a
 
real_t Abs (const Point< real_t > &p)
 Return norm of vector a
 
real_t Conjg (real_t a)
 Return complex conjugate of real number a
 
complex_t Conjg (complex_t a)
 Return complex conjugate of complex number a
 
string out_complex (complex_t z)
 Return string to conveniently display a complex number.
 
string out_complex (real_t x, real_t y)
 Return string to conviently display a complex number.
 
real_t Max (real_t a, real_t b, real_t c)
 Return maximum value of real numbers a, b and c
 
int Kronecker (int i, int j)
 Return Kronecker delta of i and j.
 
template<class T_ >
const T_ & Max (const T_ &a, const T_ &b)
 Return maximum value of elements a and b
 
template<class T_ >
const T_ & Min (const T_ &a, const T_ &b)
 Return minimum value of elements a and b
 
int Max (int a, int b, int c)
 Return maximum value of integer numbers a, b and c
 
real_t Min (real_t a, real_t b, real_t c)
 Return minimum value of real numbers a, b and c
 
int Min (int a, int b, int c)
 Return minimum value of integer numbers a, b and c
 
real_t Max (real_t a, real_t b, real_t c, real_t d)
 Return maximum value of integer numbers a, b, c and d
 
int Max (int a, int b, int c, int d)
 Return maximum value of integer numbers a, b, c and d
 
real_t Min (real_t a, real_t b, real_t c, real_t d)
 Return minimum value of real numbers a, b, c and d
 
int Min (int a, int b, int c, int d)
 Return minimum value of integer numbers a, b, c and d
 
real_t Arg (complex_t x)
 Return argument of complex number x
 
complex_t Log (complex_t x)
 Return principal determination of logarithm of complex number x
 
template<class T_ >
T_ Sqr (T_ x)
 Return square of value x
 
template<class T_ >
void Scale (T_ a, const vector< T_ > &x, vector< T_ > &y)
 Mutiply vector x by a and save result in vector y
 
template<class T_ >
void Scale (T_ a, const Vect< T_ > &x, Vect< T_ > &y)
 Mutiply vector x by a and save result in vector y
 
template<class T_ >
void Scale (T_ a, vector< T_ > &x)
 Mutiply vector x by a
 
template<class T_ >
void Xpy (size_t n, T_ *x, T_ *y)
 Add array x to y
 
template<class T_ >
void Xpy (const vector< T_ > &x, vector< T_ > &y)
 Add vector x to y
 
template<class T_ >
void Axpy (size_t n, T_ a, T_ *x, T_ *y)
 Multiply array x by a and add result to y
 
template<class T_ >
void Axpy (T_ a, const vector< T_ > &x, vector< T_ > &y)
 Multiply vector x by a and add result to y
 
template<class T_ >
void Copy (size_t n, T_ *x, T_ *y)
 Copy array x to y n is the arrays size.
 
real_t Error2 (const vector< real_t > &x, const vector< real_t > &y)
 Return absolute L2 error between vectors x and y
 
real_t RError2 (const vector< real_t > &x, const vector< real_t > &y)
 Return absolute L2 error between vectors x and y
 
real_t ErrorMax (const vector< real_t > &x, const vector< real_t > &y)
 Return absolute Max. error between vectors x and y
 
real_t RErrorMax (const vector< real_t > &x, const vector< real_t > &y)
 Return relative Max. error between vectors x and y
 
template<class T_ >
T_ Dot (size_t n, T_ *x, T_ *y)
 Return dot product of arrays x and y
 
real_t Dot (const vector< real_t > &x, const vector< real_t > &y)
 Return dot product of vectors x and y.
 
template<class T_ >
T_ Dot (const Point< T_ > &x, const Point< T_ > &y)
 Return dot product of x and y
 
real_t exprep (real_t x)
 Compute the exponential function with avoiding over and underflows.
 
template<class T_ >
void Assign (vector< T_ > &v, const T_ &a)
 Assign the value a to all entries of a vector v
 
template<class T_ >
void clear (vector< T_ > &v)
 Assign 0 to all entries of a vector.
 
template<class T_ >
void clear (Vect< T_ > &v)
 Assign 0 to all entries of a vector.
 
real_t Nrm2 (size_t n, real_t *x)
 Return 2-norm of array x
 
real_t Nrm2 (const vector< real_t > &x)
 Return 2-norm of vector x
 
template<class T_ >
real_t Nrm2 (const Point< T_ > &a)
 Return 2-norm of a
 
bool Equal (real_t x, real_t y, real_t toler=OFELI_EPSMCH)
 Function to return true if numbers x and y are close up to a given tolerance toler
 
char itoc (int i)
 Function to convert an integer to a character.
 
template<class T_ >
std::string toString (const T_ x)
 Function to convert any value to a string.
 
template<class T_ >
T_ stringTo (const std::string &s)
 Function to convert a string to a template type parameter.
 

Detailed Description

Utility functions and classes.

Macro Definition Documentation

◆ CATCH_EXCEPTION

#define CATCH_EXCEPTION
Value:
catch(OFELIException &e) { \
std::cout << "OFELI error: " << e.what() << endl; \
return 1; \
} \
catch(runtime_error &e) { \
std::cout << "OFELI Runtime error: " << e.what() << endl; \
return 1; \
} \
catch( ... ) { \
std::cout << "OFELI Unexpected error: " << endl; \
return 1; \
}

This macro can be inserted after a try loop to catch a thrown exception.

◆ LIGHT_SPEED

#define LIGHT_SPEED   299792458.2;

= Speed of light in SI units

◆ OFELI_E

#define OFELI_E   2.71828182845904523536028747135

Value of e or exp (with 28 digits)

◆ OFELI_EPSMCH

#define OFELI_EPSMCH   DBL_EPSILON

Value of Machine Epsilon

◆ OFELI_GAUSS2

#define OFELI_GAUSS2   0.57735026918962576450914878050196

Value of 1/sqrt(3) (with 32 digits)

◆ OFELI_IMAG

#define OFELI_IMAG   std::complex<double>(0.,1.);

= Unit imaginary number (i)

◆ OFELI_ONEOVERPI

#define OFELI_ONEOVERPI   0.31830988618379067153776752675

Value of 1/Pi (with 28 digits)

◆ OFELI_PI

#define OFELI_PI   3.14159265358979323846264338328

Value of Pi (with 28 digits)

◆ OFELI_SIXTH

#define OFELI_SIXTH   0.16666666666666666666666666667

Value of 1/6 (with 28 digits)

◆ OFELI_SQRT2

#define OFELI_SQRT2   1.41421356237309504880168872421

Value of sqrt(2) (with 28 digits)

◆ OFELI_SQRT3

#define OFELI_SQRT3   1.73205080756887729352744634151

Value of sqrt(3) (with 28 digits)

◆ OFELI_THIRD

#define OFELI_THIRD   0.33333333333333333333333333333

Value of 1/3 (with 28 digits)

◆ OFELI_TOLERANCE

#define OFELI_TOLERANCE   OFELI_EPSMCH*10000

Default tolerance for an iterative process = OFELI_EPSMCH * 10000

◆ OFELI_TWELVETH

#define OFELI_TWELVETH   0.08333333333333333333333333333

Value of 1/12 (with 28 digits)

◆ VLG

#define VLG   1.e10

Very large number: A real number for penalty

Function Documentation

◆ areClose() [1/2]

bool areClose ( const Point< real_t > &  a,
const Point< real_t > &  b,
real_t  toler = OFELI_TOLERANCE 
)

Return true if both instances of class Point<double> are distant with less then toler

Author
Rachid Touzani

◆ areClose() [2/2]

bool areClose ( const SpaceTime a,
const SpaceTime b,
real_t  toler = OFELI_TOLERANCE 
)

Return true if both instances of class Point<double> are distant with less then toler

Author
Rachid Touzani

◆ Axpy() [1/2]

template<class T_ >
void Axpy ( size_t  n,
T_  a,
T_ *  x,
T_ *  y 
)

Multiply array x by a and add result to y

n is the arrays size.

◆ Axpy() [2/2]

template<class T_ >
void Axpy ( T_  a,
const vector< T_ > &  x,
vector< T_ > &  y 
)

Multiply vector x by a and add result to y

x and y are instances of class vector<T_>

◆ banner()

void banner ( const string &  prog = " ")

Outputs a banner as header of any developed program.

Parameters
[in]progCalling program name. Enables writing a copyright notice accompanying the program.
Author
Rachid Touzani

◆ BSpline()

BSpline ( size_t  n,
size_t  t,
Vect< Point< real_t > > &  control,
Vect< Point< real_t > > &  output,
size_t  num_output 
)

Function to perform a B-spline interpolation.

This program is adapted from a free program ditributed by Keith Vertanen (verta.nosp@m.nkd@.nosp@m.cda.m.nosp@m.rs.u.nosp@m.mn.ed.nosp@m.u) in 1994.

Parameters
[in]nNumber of control points minus 1.
[in]tDegree of the polynomial plus 1.
[in]controlControl point array made up of Point stucture.
[out]outputVector in which the calculated spline points are to be put.
[in]num_outputHow many points on the spline are to be calculated.
Note
Condition: n+2>t (No curve results if n+2<=t) Control vector contains the number of points specified by n Output array is the proper size to hold num_output point structures
Author
Rachid Touzani

◆ clear() [1/2]

template<class T_ >
void clear ( Vect< T_ > &  v)

Assign 0 to all entries of a vector.

Parameters
[in]vVector to clear

◆ clear() [2/2]

template<class T_ >
void clear ( vector< T_ > &  v)

Assign 0 to all entries of a vector.

Parameters
[in]vVector to clear

◆ Discrepancy() [1/2]

real_t Discrepancy ( Vect< complex_t > &  x,
const Vect< complex_t > &  y,
int  n,
int  type = 1 
)

Return discrepancy between 2 vectors x and y

Parameters
[in,out]xFirst vector (Instance of class Vect). On output, x is assigned the vector y
[in]ySecond vector (Instance of class Vect)
[in]nType of norm
  • 1: Weighted 1-Norm
  • 2: Weighted 2-Norm
  • 0: Max-Norm
[in]typeDiscrepancy type (0: Absolute, 1: Relative [Default])
Returns
Computed discrepancy value

◆ Discrepancy() [2/2]

real_t Discrepancy ( Vect< real_t > &  x,
const Vect< real_t > &  y,
int  n,
int  type = 1 
)

Return discrepancy between 2 vectors x and y

Parameters
[in,out]xFirst vector (Instance of class Vect). On output, x is assigned the vector y
[in]ySecond vector (Instance of class Vect)
[in]nType of norm
  • 1: Weighted 1-Norm
  • 2: Weighted 2-Norm
  • 0: Max-Norm
[in]typeDiscrepancy type (0: Absolute, 1: Relative [Default])
Returns
Computed discrepancy value

◆ Distance() [1/3]

real_t Distance ( const Point2D< real_t > &  a,
const Point2D< real_t > &  b 
)

Return euclidean distance between points a and b

Author
Rachid Touzani

◆ Distance() [2/3]

double Distance ( const Point< real_t > &  a,
const Point< real_t > &  b 
)

Return euclidean distance between points a and b

Author
Rachid Touzani

◆ Distance() [3/3]

double Distance ( const SpaceTime a,
const SpaceTime b 
)

Return euclidean distance between points a and b

Author
Rachid Touzani

◆ Dot() [1/2]

double Dot ( const vector< real_t > &  x,
const vector< real_t > &  y 
)

Return dot product of vectors x and y.

x and y are instances of class vector<double>

◆ Dot() [2/2]

template<class T_ >
T_ Dot ( size_t  n,
T_ *  x,
T_ *  y 
)

Return dot product of arrays x and y

n is the arrays size.

◆ Equal()

bool Equal ( real_t  x,
real_t  y,
real_t  toler = OFELI_EPSMCH 
)

Function to return true if numbers x and y are close up to a given tolerance toler

Default value of tolerance is the constant OFELI_EPSMCH

◆ getBamg()

void getBamg ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in Bamg format.

Parameters
[in]fileName of a file written in the Bamg format.
Note
Bamg is a 2-D mesh generator. It allows to construct adapted meshes from a given metric. It was developed at INRIA, France. Available information can be found in the site:
http://raweb.inria.fr/rapportsactivite/RA2002/gamma/uid25.html
Parameters
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ getEasymesh()

void getEasymesh ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in Easymesh format.

Parameters
[in]fileName of a file (without extension) written in Easymesh format. Actually, the function Easymesh2MDF attempts to read mesh data from files file.e, file.n and file.s produced by Easymesh.
Note
Easymesh is a free program that generates 2-D, unstructured, Delaunay and constrained Delaunay triangulations in general domains. It can be downloaded from the site:
http://www-dinma.univ.trieste.it/nirftc/research/easymesh/Default.htm
Parameters
[in]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ getGambit()

void getGambit ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in Gambit neutral format.

Note
Gambit is a commercial mesh generator associated to the CFD code Fluent. Informations about Gambit can be found in the site:
http://www.fluent.com/software/gambit/
Parameters
[in]fileName of a file written in the Gambit neutral format.
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ getGmsh()

void getGmsh ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in Gmsh format.

Note
Gmsh is a free mesh generator that can be downloaded from the site:
http://www.geuz.org/gmsh/
Parameters
[in]fileName of a file written in the Gmsh format.
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ getMatlab()

void getMatlab ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a Matlab mesh data.

Note
Matlab is a language of scientific computing including visualization. It is developed by MathWorks. Available information can be found in the site:
http://www.mathworks.com/products/matlab/
Parameters
[in]fileName of a file created by Matlab by executing the script file Matlab2OFELI.m
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ getMesh()

void getMesh ( string  file,
ExternalFileFormat  form,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in an external file format.

Parameters
[in]fileInput mesh file name.
[in]formFormat of the mesh file. This one can be chosen among the enumerated values:
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ getNetgen()

void getNetgen ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in Netgen format.

Note
Netgen is a tetrahedral mesh generator that can be downloaded from the site:
http://www.hpfem.jku.at/netgen/
Parameters
[in]fileName of a file written in the Netgen format.
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. [ default = 1 ]
Author
Rachid Touzani

◆ getTetgen()

void getTetgen ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in Tetgen format.

Note
Tetgen is a free three-dimensional mesh generator that can be downloaded in the site:
http://tetgen.berlios.de/
Parameters
[in]fileName of a file written in the Tetgen format.
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ getTriangle()

void getTriangle ( string  file,
Mesh mesh,
size_t  nb_dof = 1 
)

Construct an instance of class Mesh from a mesh file stored in Triangle format.

Note
TRIANGLE is a C program that can generate meshes, Delaunay triangulations and Voronoi diagrams for 2D pointsets that can be downloaded in the site:
http://people.scs.fsu.edu/~burkardt/c_src/triangle/triangle.html/
Parameters
[in]fileName of a file written in the Tetgen format.
[out]meshMesh instance created by the function.
[in]nb_dofNumber of degrees of freedom for each node. This information is not provided, in general, by mesh generators. Its default value here is 1.
Author
Rachid Touzani

◆ Nrm2()

real_t Nrm2 ( size_t  n,
real_t *  x 
)

Return 2-norm of array x

Parameters
[in]nis Array length
[in]xArray to treat

◆ operator*() [1/15]

template<class T_ >
Point2D< T_ > operator* ( const int &  a,
const Point2D< T_ > &  b 
)

Operator *.

Return point b divided by integer constant a

Author
Rachid Touzani

◆ operator*() [2/15]

template<class T_ >
Point< T_ > operator* ( const int &  a,
const Point< T_ > &  b 
)

Operator *.

Return point b divided by integer constant a

Author
Rachid Touzani

◆ operator*() [3/15]

SpaceTime operator* ( const int &  a,
const SpaceTime b 
)

Operator *.

Return point b divided by integer constant a

Author
Rachid Touzani

◆ operator*() [4/15]

template<class T_ >
Point2D< T_ > operator* ( const Point2D< T_ > &  b,
const int &  a 
)

Operator *

Return point b postmultiplied by constant a

Author
Rachid Touzani

◆ operator*() [5/15]

template<class T_ >
T_ operator* ( const Point2D< T_ > &  b,
const Point2D< T_ > &  a 
)

Operator *.

Return point b postmultiplied by integer constant a.

Author
Rachid Touzani

◆ operator*() [6/15]

template<class T_ >
Point2D< T_ > operator* ( const Point2D< T_ > &  b,
const T_ &  a 
)

Operator /

Return point b postmultiplied by constant a

Author
Rachid Touzani

◆ operator*() [7/15]

template<class T_ >
T_ operator* ( const Point< T_ > &  b,
const Point< T_ > &  a 
)

Operator *

Return inner (scalar) product of points a and b

Author
Rachid Touzani

◆ operator*() [8/15]

template<class T_ >
Point< T_ > operator* ( const Point< T_ > &  b,
const int &  a 
)

Operator *

Return point b postmultiplied by constant a

Author
Rachid Touzani

◆ operator*() [9/15]

template<class T_ >
Point< T_ > operator* ( const Point< T_ > &  b,
const T_ &  a 
)

Operator /

Return point b multiplied by constant a

◆ operator*() [10/15]

SpaceTime operator* ( const real_t &  a,
const SpaceTime b 
)

Operator *

Return point b premultiplied by constant a

Author
Rachid Touzani

◆ operator*() [11/15]

real_t operator* ( const SpaceTime b,
const SpaceTime a 
)

Operator *

Return inner (scalar) product of points a and b

Author
Rachid Touzani

◆ operator*() [12/15]

SpaceTime operator* ( const SpaceTime b,
const int &  a 
)

Operator *

Return point b postmultiplied by constant a

Author
Rachid Touzani

◆ operator*() [13/15]

SpaceTime operator* ( const SpaceTime b,
const real_t &  a 
)

Operator /

Return point b multiplied by constant a

◆ operator*() [14/15]

template<class T_ >
Point2D< T_ > operator* ( const T_ &  a,
const Point2D< T_ > &  b 
)

Operator *.

Return point b premultiplied by constant a

Author
Rachid Touzani

◆ operator*() [15/15]

template<class T_ >
Point< T_ > operator* ( const T_ &  a,
const Point< T_ > &  b 
)

Operator *

Return point b premultiplied by constant a

Author
Rachid Touzani

◆ operator+() [1/6]

template<class T_ >
Point2D< T_ > operator+ ( const Point2D< T_ > &  a,
const Point2D< T_ > &  b 
)

Operator +.

Return sum of two points a and b

Author
Rachid Touzani

◆ operator+() [2/6]

template<class T_ >
Point2D< T_ > operator+ ( const Point2D< T_ > &  a,
const T_ &  x 
)

Operator +.

Translate a by x

Author
Rachid Touzani

◆ operator+() [3/6]

template<class T_ >
Point< T_ > operator+ ( const Point< T_ > &  a,
const Point< T_ > &  b 
)

Operator +

Return sum of two points a and b

Author
Rachid Touzani

◆ operator+() [4/6]

template<class T_ >
Point< T_ > operator+ ( const Point< T_ > &  a,
const T_ &  x 
)

Operator +

Translate a by x

Author
Rachid Touzani

◆ operator+() [5/6]

SpaceTime operator+ ( const SpaceTime a,
const real_t &  x 
)

Operator +

Translate a by x

Author
Rachid Touzani

◆ operator+() [6/6]

SpaceTime operator+ ( const SpaceTime a,
const SpaceTime b 
)

Operator +

Return sum of two points a and b

Author
Rachid Touzani

◆ operator-() [1/9]

template<class T_ >
Point2D< T_ > operator- ( const Point2D< T_ > &  a)

Unary Operator -

Return minus a

Author
Rachid Touzani

◆ operator-() [2/9]

template<class T_ >
Point2D< T_ > operator- ( const Point2D< T_ > &  a,
const Point2D< T_ > &  b 
)

Operator -

Return point a minus point b

Author
Rachid Touzani

◆ operator-() [3/9]

template<class T_ >
Point2D< T_ > operator- ( const Point2D< T_ > &  a,
const T_ &  x 
)

Operator -

Translate a by -x

Author
Rachid Touzani

◆ operator-() [4/9]

template<class T_ >
Point< T_ > operator- ( const Point< T_ > &  a)

Unary Operator -

Return minus a

Author
Rachid Touzani

◆ operator-() [5/9]

template<class T_ >
Point< T_ > operator- ( const Point< T_ > &  a,
const Point< T_ > &  b 
)

Operator -

Return point a minus point b

Author
Rachid Touzani

◆ operator-() [6/9]

template<class T_ >
Point< T_ > operator- ( const Point< T_ > &  a,
const T_ &  x 
)

Operator -

Translate a by -x

Author
Rachid Touzani

◆ operator-() [7/9]

SpaceTime operator- ( const SpaceTime a)

Unary Operator -

Return minus a

Author
Rachid Touzani

◆ operator-() [8/9]

SpaceTime operator- ( const SpaceTime a,
const real_t &  x 
)

Operator -

Translate a by -x

Author
Rachid Touzani

◆ operator-() [9/9]

SpaceTime operator- ( const SpaceTime a,
const SpaceTime b 
)

Operator -

Return point a minus point b

Author
Rachid Touzani

◆ operator/() [1/3]

template<class T_ >
Point2D< T_ > operator/ ( const Point2D< T_ > &  b,
const T_ &  a 
)

Operator /

Return point b divided by constant a

Author
Rachid Touzani

◆ operator/() [2/3]

template<class T_ >
Point< T_ > operator/ ( const Point< T_ > &  b,
const T_ &  a 
)

Operator /

Return point b divided by constant a

Author
Rachid Touzani

◆ operator/() [3/3]

SpaceTime operator/ ( const SpaceTime b,
const real_t &  a 
)

Operator /

Return point b divided by constant a

Author
Rachid Touzani

◆ operator<() [1/2]

bool operator< ( const Point< size_t > &  a,
const Point< size_t > &  b 
)

Comparison operator. Returns true if all components of first vector are lower than those of second one.

Return minus a

Author
Rachid Touzani

◆ operator<() [2/2]

bool operator< ( const SpaceTime a,
const SpaceTime b 
)

Comparison operator. Returns true if all components of first vector are lower than those of second one.

Return minus a

Author
Rachid Touzani

◆ operator<<() [1/8]

ostream & operator<< ( ostream &  s,
const complex_t &  x 
)

Output a complex number.

Author
Rachid Touzani

◆ operator<<() [2/8]

template<class T_ >
ostream & operator<< ( ostream &  s,
const std::list< T_ > &  l 
)

Output a vector instance.

Author
Rachid Touzani

◆ operator<<() [3/8]

template<class T_ >
ostream & operator<< ( ostream &  s,
const std::pair< T_, T_ > &  a 
)

Output a pair instance.

Author
Rachid Touzani

◆ operator<<() [4/8]

ostream & operator<< ( ostream &  s,
const std::string &  c 
)

Output a string.

Author
Rachid Touzani

◆ operator<<() [5/8]

template<class T_ >
ostream & operator<< ( ostream &  s,
const vector< T_ > &  v 
)

Output a vector instance.

Author
Rachid Touzani

◆ operator<<() [6/8]

template<class T_ >
ostream & operator<< ( std::ostream &  s,
const Point2D< T_ > &  a 
)

Output point coordinates.

Author
Rachid Touzani

◆ operator<<() [7/8]

template<class T_ >
ostream & operator<< ( std::ostream &  s,
const Point< T_ > &  a 
)

Output point coordinates.

Author
Rachid Touzani

◆ operator<<() [8/8]

ostream & operator<< ( std::ostream &  s,
const SpaceTime a 
)

Output space-time point coordinates.

Author
Rachid Touzani

◆ operator==() [1/3]

template<class T_ >
bool operator== ( const Point2D< T_ > &  a,
const Point2D< T_ > &  b 
)

Operator ==.

Return true if a=b, false if not.

Author
Rachid Touzani

◆ operator==() [2/3]

template<class T_ >
bool operator== ( const Point< T_ > &  a,
const Point< T_ > &  b 
)

Operator ==

Return true if a=b, false if not.

Author
Rachid Touzani

◆ operator==() [3/3]

bool operator== ( const SpaceTime a,
const SpaceTime b 
)

Operator ==

Return true if a=b, false if not.

Author
Rachid Touzani

◆ out_complex() [1/2]

string out_complex ( complex_t  z)

Return string to conveniently display a complex number.

Parameters
[in]zComplex number

◆ out_complex() [2/2]

string out_complex ( real_t  x,
real_t  y 
)

Return string to conviently display a complex number.

Parameters
[in]xReal part
[in]yImaginary part

◆ qksort() [1/2]

template<class T_ >
void qksort ( std::vector< T_ > &  a,
int  begin,
int  end 
)

Function to sort a vector.

qksort uses the famous quick sorting algorithm.

Parameters
[in,out]aVector to sort.
[in]beginindex of starting index (default value is 0)
[in]endindex of ending index (default value is the vector size - 1)
Author
Rachid Touzani

◆ qksort() [2/2]

template<class T_ , class C_ >
void qksort ( std::vector< T_ > &  a,
int  begin,
int  end,
C_  compare 
)

Function to sort a vector according to a key function.

qksort uses the famous quick sorting algorithm.

Parameters
[in,out]aVector to sort.
[in]beginindex of starting index (0 for the beginning of the vector)
[in]endindex of ending index
[in]compareA function object that implements the ordering. The user must provide this function that returns a boolean function that is true if the first argument is less than the second and false if not.
Author
Rachid Touzani

◆ QuickSort()

template<class T_ >
void QuickSort ( std::vector< T_ > &  a,
int  begin,
int  end 
)

Function to sort a vector.

qksort uses the famous quick sorting algorithm.

Parameters
[in,out]aVector to sort.
[in]beginindex of starting iterator
[in]endindex of ending iterator

The calling program must provide an overloading of the operator < for the type T_

Author
Rachid Touzani

◆ saveField() [1/3]

void saveField ( const Vect< real_t > &  v,
const Mesh mesh,
string  output_file,
int  opt 
)

Save a vector to an output file in a given file format.

Case where the vector does not contain mesh information

Parameters
[in]vVect instance to save
[in]meshMesh instance
[in]output_fileOutput file where to save the vector
[in]optOption to choose file format to save. This is to be chosen among enumerated values: GMSH, GNUPLOT, MATLAB, TECPLOT, VTK
Author
Rachid Touzani

◆ saveField() [2/3]

void saveField ( Vect< real_t > &  v,
const Grid g,
string  output_file,
int  opt = VTK 
)

Save a vector to an output file in a given file format, for a structured grid data.

Parameters
[in]vVect instance to save
[in]gGrid instance
[in]output_fileOutput file where to save the vector
[in]optOption to choose file format to save. This is to be chosen among enumerated values: GMSH, VTK
Author
Rachid Touzani

◆ saveField() [3/3]

void saveField ( Vect< real_t > &  v,
string  output_file,
int  opt 
)

Save a vector to an output file in a given file format.

Case where the vector contains mesh information

Parameters
[in]vVect instance to save
[in]output_fileOutput file where to save the vector
[in]optOption to choose file format to save. This is to be chosen among enumerated values: GMSH, GNUPLOT, MATLAB, TECPLOT, VTK
Author
Rachid Touzani

◆ saveGmsh() [1/3]

void saveGmsh ( const string &  file,
const Mesh mesh 
)

This function outputs a Mesh instance in a file in Gmsh format.

Note
Gmsh is a free mesh generator that can be downloaded from the site: http://www.geuz.org/gmsh/
Parameters
[in]fileOutput file in Gmsh format.
[in]meshMesh instance to save.
Author
Rachid Touzani

◆ saveGmsh() [2/3]

void saveGmsh ( Mesh mesh,
string  input_file,
string  output_file,
int  f = 1 
)

Save a vector to an output Gmsh file.

Gmsh is a free mesh generator and postprocessor that can be downloaded from the site:
http://www.geuz.org/gmsh/

Parameters
[in]meshReference to Mesh instance
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (Gmsh format file)
[in]fField is stored each f time step [Default: 1]
Author
Rachid Touzani

◆ saveGmsh() [3/3]

void saveGmsh ( string  input_file,
string  output_file,
string  mesh_file,
int  f = 1 
)

Save a vector to an output Gmsh file.

Gmsh is a free mesh generator and postprocessor that can be downloaded from the site:
http://www.geuz.org/gmsh/

Parameters
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (Gmsh format file)
[in]mesh_fileFile containing mesh data
[in]fField is stored each f time step [Default: 1]
Author
Rachid Touzani

◆ saveGnuplot() [1/3]

void saveGnuplot ( const string &  file,
const Mesh mesh 
)

This function outputs a Mesh instance in a file in Gmsh format.

Note
Gnuplot is a command-line driven program for producing 2D and 3D plots. It is under the GNU General Public License. Available information can be found in the site:
http://www.gnuplot.info/
Parameters
[out]fileOutput file in Gnuplot format.
[in]meshMesh instance to save.
Author
Rachid Touzani

◆ saveGnuplot() [2/3]

void saveGnuplot ( Mesh mesh,
string  input_file,
string  output_file,
int  f = 1 
)

Save a vector to an input Gnuplot file.

Gnuplot is a command-line driven program for producing 2D and 3D plots. It is under the GNU General Public License. Available information can be found in the site:
http://www.gnuplot.info/

Parameters
[in]meshReference to Mesh instance
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (gnuplot format file)
[in]fField is stored each f time step [Default: 1]
Author
Rachid Touzani

◆ saveGnuplot() [3/3]

void saveGnuplot ( string  input_file,
string  output_file,
string  mesh_file,
int  f = 1 
)

Save a vector to an input Gnuplot file.

Gnuplot is a command-line driven program for producing 2D and 3D plots. It is under the GNU General Public License. Available information can be found in the site:
http://www.gnuplot.info/

Parameters
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (gnuplot format file)
[in]mesh_fileFile containing mesh data
[in]fField is stored each f time step [Default: 1
Author
Rachid Touzani

◆ saveMatlab()

void saveMatlab ( const string &  file,
const Mesh mesh 
)

This function outputs a Mesh instance in a file in Matlab format.

Note
Matlab is a language of scientific computing including visualization. It is developed by MathWorks. Available information can be found in the site:
http://www.mathworks.com/products/matlab/
Parameters
[out]fileOutput file in Matlab format.
[in]meshMesh instance to save.
Author
Rachid Touzani

◆ saveMesh()

void saveMesh ( const string &  file,
const Mesh mesh,
ExternalFileFormat  form 
)

This function saves mesh data a file for a given external format.

Parameters
[in]fileFile where to store mesh
[in]meshMesh instance to save
[in]formFormat of the mesh file. This one can be chosen among the enumerated values:
Author
Rachid Touzani

◆ saveTecplot() [1/3]

void saveTecplot ( const string &  file,
const Mesh mesh 
)

This function outputs a Mesh instance in a file in Tecplot format.

Note
Tecplot is high quality post graphical commercial processing program developed by Amtec. Available information can be found in the site:
http://www.tecplot.com
Parameters
[out]fileOutput file in Tecplot format.
[in]meshMesh instance to save.
Author
Rachid Touzani

◆ saveTecplot() [2/3]

void saveTecplot ( Mesh mesh,
string  input_file,
string  output_file,
int  f = 1 
)

Save a vector to an output file to an input Tecplot file.

Tecplot is high quality post graphical commercial processing program developed by Amtec. Available information can be found in the site: http://www.tecplot.com

Parameters
[in]meshReference to Mesh instance
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (gnuplot format file)
[in]fField is stored each f time step [Default: 1]
Author
Rachid Touzani

◆ saveTecplot() [3/3]

void saveTecplot ( string  input_file,
string  output_file,
string  mesh_file,
int  f = 1 
)

Save a vector to an output file to an input Tecplot file.

Tecplot is high quality post graphical commercial processing program developed by Amtec. Available information can be found in the site: http://www.tecplot.com

Parameters
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (gnuplot format file)
[in]mesh_fileFile containing mesh data
[in]fField is stored each f time step [Default: 1]
Author
Rachid Touzani

◆ saveVTK() [1/3]

void saveVTK ( const string &  file,
const Mesh mesh 
)

This function outputs a Mesh instance in a file in VTK format.

Note
The Visualization ToolKit (VTK) is an open source, freely available software system for 3D computer graphics. Available information can be found in the site:
http://public.kitware.com/VTK/
Parameters
[out]fileOutput file in VTK format.
[in]meshMesh instance to save.
Author
Rachid Touzani

◆ saveVTK() [2/3]

saveVTK ( Mesh mesh,
string  input_file,
string  output_file,
int  f = 1 
)

Save a vector to an output VTK file.

The Visualization ToolKit (VTK) is an open source, freely available software system for 3D computer graphics. Available information can be found in the site:
http://public.kitware.com/VTK/

Parameters
[in]meshReference to Mesh instance
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (VTK format file)
[in]fField is stored each f time step [Default: 1]
Author
Rachid Touzani

◆ saveVTK() [3/3]

saveVTK ( string  input_file,
string  output_file,
string  mesh_file,
int  f = 1 
)

Save a vector to an output VTK file.

The Visualization ToolKit (VTK) is an open source, freely available software system for 3D computer graphics. Available information can be found in the site:
http://public.kitware.com/VTK/

Parameters
[in]input_fileInput file (OFELI XML file containing a field).
[in]output_fileOutput file (VTK format file)
[in]mesh_fileFile containing mesh data
[in]fField is stored each f time step [Default: 1]
Author
Rachid Touzani

◆ Scale() [1/3]

template<class T_ >
void Scale ( T_  a,
const Vect< T_ > &  x,
Vect< T_ > &  y 
)

Mutiply vector x by a and save result in vector y

x and y are instances of class Vect<T_>

◆ Scale() [2/3]

template<class T_ >
void Scale ( T_  a,
const vector< T_ > &  x,
vector< T_ > &  y 
)

Mutiply vector x by a and save result in vector y

x and y are instances of class vector<T_>

◆ Scale() [3/3]

template<class T_ >
void Scale ( T_  a,
vector< T_ > &  x 
)

Mutiply vector x by a

x is an instance of class vector<T_>

◆ SqrDistance() [1/3]

real_t SqrDistance ( const Point2D< real_t > &  a,
const Point2D< real_t > &  b 
)

Return squared euclidean distance between points a and b

Author
Rachid Touzani

◆ SqrDistance() [2/3]

double SqrDistance ( const Point< real_t > &  a,
const Point< real_t > &  b 
)

Return squared euclidean distance between points a and b

Author
Rachid Touzani

◆ SqrDistance() [3/3]

double SqrDistance ( const SpaceTime a,
const SpaceTime b 
)

Return squared euclidean distance between points a and b

Author
Rachid Touzani

◆ Xpy()

template<class T_ >
void Xpy ( const vector< T_ > &  x,
vector< T_ > &  y 
)

Add vector x to y

x and y are instances of class vector<T_>