Various functions and classes
Classes that are related to algebra and contained in the directory Share are included when MontjoieCommon.hxx is included. A basic example is given as
#include "Share/MontjoieCommon.hxx" using namespace Montjoie; int main(int argc, char** argv) { InitMontjoie(argc, argv); FinalizeMontjoie(); }
All the classes and functions presented in this section should work with this example. In this page, we have documented types defined in Montjoie (Real_wp, R2, etc), global variables (glob_chrono, rank_processor, etc), methods associated with the classes VarRandomGenerator, BasicTimer and finally standalone functions.
Defined types
Below we have listed the types defined in Montjoie (with typedef command), they are present in the files MontjoieTypes.hxx and Precision.hxx.
Real_wp | real number in double precision (multiple precision if mpfr is enabled) |
Complex_wp | complex number in double precision (multiple precision if mpfr is enabled) |
VectReal_wp | vector of real numbers in double precision (multiple precision if mpfr is enabled) |
VectComplex_wp | vector of complex numbers in double precision (multiple precision if mpfr is enabled) |
VectBool | vector of booleans |
VectString | vector of strings |
R2 | a tiny vector with two real numbers |
R2_Complex_wp | a tiny vector with two complex numbers |
VectR2 | vector of tiny vectors with two real numbers |
VectR2_Complex_wp | vector of tiny vectors with two complex numbers |
Matrix2_2 | a tiny real matrix with two rows and columns |
Matrix2_2_Complex_wp | a tiny complex matrix with two rows and columns |
Matrix2_2sym | a tiny real symmetric matrix with two rows and columns |
Matrix2_2sym_Complex_wp | a tiny complex symmetric matrix with two rows and columns |
R3 | a tiny vector with three real numbers |
R3_Complex_wp | a tiny vector with three complex numbers |
VectR3 | vector of tiny vectors with three real numbers |
VectR3_Complex_wp | vector of tiny vectors with three complex numbers |
Matrix3_3 | a tiny real matrix with three rows and columns |
Matrix3_3_Complex_wp | a tiny complex matrix with three rows and columns |
Matrix3_3sym | a tiny real symmetric matrix with three rows and columns |
Matrix3_3sym_Complex_wp | a tiny complex symmetric matrix with three rows and columns |
The classes Dimension1, Dimension2 and Dimension3 have a single static attribute (dim_N which is the dimension) and defined types such that you can specify a dimension-dependent type if you are using Dimension as a template parameter. See the example below :
// definition of a class depending on the dimension template<class Dimension> class MyClass { public : void MyMethod() { // considering a point in R^n typename Dimension::R_N point; // a matrix with n rows and columns typename Dimension::MatrixN_N mat_dfj; // a vector of points typename Dimension::VectR_N list_points; // element of a mesh typename Dimension::Element_Rn elt; } } // then in the main you instantiate this class with Dimension1, Dimension2 or Dimension3 MyClass<Dimension3> var;
Global variables
Below we have listed the global variables of Montjoie.
pi_wp | Pi in machine precision (double or multiple) |
Iwp | pure imaginary number |
epsilon_machine | machine precision |
rank_global_processor | rank of the current processor (within COMM_WORLD communicator) |
size_all_communications | size of datas exchanged during MPI calls |
Seed | random seed |
is_seed_defined | true if the random seed is defined |
glob_chrono | global timer that you can use everywhere |
Methods for VarRandomGenerator
We begin with an example of use
VarRandomGenerator var; // default random generator is the rand() function // you can select a generator from gsl (if flag MONTJOIE_WITH_GSL is set) var.SetRandomGenerator(var.GSL_MRG); // then you can obtain a single random number double x = var.GetRand(); // or generate a complete vector of random numbers int n = 20; VectReal_wp X; var.GenerateRandomNumbers(n, X); // you can also shuffle elements of X var.ApplyRandomPermutation(X);
Clear | releases memory used by random generator |
RoundToSquare | finds closest integers nx, ny such that n = nx ny |
GetRand | returns a random number between 0 and 1 |
GenerateRandomNumbers | generates a vector of random numbers |
ApplyRandomPermutation | performs a random permutation between elements of a vector |
GenerateRandomPermutation | computes a random permutation (i.e. array containing values 0, 1, 2, ..., N-1 in a random order) |
SetRandomGenerator | selects a random generator |
SetInputRandomGenerator | selects a random generator |
GenerateBitReversedNumbers | computes low-discrepancy random numbers with reversed-bit method |
Methods for BasicTimer/AccurateTimer/RealTimer
The class BasicTimer is the default class used to compute times elapsed in Montjoie, it uses the function clock. If the flag MONTJOIE_USE_ACCURATE_TIMING is set, the class AccurateTimer will be used (this class is based on function getrusage). If the flag MONTJOIE_USE_REAL_TIMING is defined, the class RealTimer will be used (this class is based on the function clock_gettime). Below an example of the class BasicTimer is detailed (the classes AccurateTimer and RealTimer can be used in a similar fashion)
// variable storing different timers (=stopwatches) BasicTimer my_chrono; // there are some predefined stopwatches : // ALL, COMM, PROD, STIFFNESS, FLUX, MASS, EXTRAPOL, SCHEME // PML, OUTPUT, JACOBIAN, SOLVE, FACTO my_chrono.Start(my_chrono.ALL); // then you can interrupt the stopwatch my_chrono.Stop(my_chrono.ALL); // if you start again, the stopwatch will continue and time will be incremented my_chrono.Start(my_chrono.ALL); // then you can retrieve the number of seconds my_chrono.Stop(my_chrono.ALL); DISP(my_chrono.GetSeconds(my_chrono.ALL)); // if you wish to reset the stopwatch (set the elapsed time to 0) my_chrono.Reset(my_chrono.ALL); // you can also create a new stopwatch (each stopwatch is associated with a number) int chrono_mesh = my_chrono.GetNumber(); my_chrono.Start(chrono_mesh); my_chrono.Stop(chrono_mesh); DISP(my_chrono.GetSeconds(chrono_mesh)); // and you can destroy the stopwatch my_chrono.ReleaseNumber(chrono_mesh);
Below, you can see all the methods associated with those classes
GetNumber | creates a new timer, returns the number associated with this timer |
ReleaseNumber | destroys timer |
ReserveNumber | reserves a number so that this number won't be used when creating a new timer |
Reset | sets to 0 the elapsed time of a timer |
Start | starts a timer |
Stop | stops a timer, the elapsed time will be kept for the next start |
GetSeconds | returns the elapsed time in seconds |
Methods for MontjoieTimer
The class MontjoieTimer enables the user to select the timer he prefers. It also can be used in parallel computing in order to display the sum of elapsed times (for each processor) and the maximal time. The timers stored in MontjoieTimer can also be accessed with string instead of numbers.
MontjoieTimer my_chrono; // you can switch to another timer my_chrono.SetTimer(MontjoieTimer::ACCURATE_TIMER); // then you can create a stopwatch with a string instead of numbers // you can associate a message with this stopwatch my_chrono.SetMessage("Factorization", "Factorize the finite element matrix"); // then the name given (here Factorization) should be used // for starting, stopping, resetting the timer my_chrono.Start("Factorization"); // the computational task is performed Matrix<double> A; IVect pivot; GetLU(A, pivot); my_chrono.Stop("Factorization"); // you can display the time spent with DisplayTime // it will use the message specified when SetMessage has been called my_chrono.DisplayTime("Factorization"); // all the created timers can be displayed: my_chrono.DisplayAll();
Below, you can see all the methods associated with class MontjoieTimer
SetTimer | specifies with timing system has to be used |
SetCommunicator | specifies the MPI communicator |
GetNumber | creates a new timer, returns the number associated with this timer |
ReleaseNumber | destroys timer |
ReserveNumber | reserves a number so that this number won't be used when creating a new timer |
Reset | sets to 0 the elapsed time of a timer |
Start | starts a timer |
Stop | stops a timer, the elapsed time will be kept for the next start |
GetSeconds | returns the elapsed time in seconds |
SetMessage | creates a timer and associates a message with it |
ReserveName | creates a timer with a given name |
ReleaseName | destroys timer with a given name |
NameExists | returns true if a timer already exists |
GetGlobalSeconds | retrieves the elapsed time for each processor, and also the sum, maximum and minimum |
DisplayTime | displays the elapsed time for a given timer |
DisplayAll | displays elapsed times for all timers |
Functions
Here a list of standalone functions implemented in Montjoie
Sort | sorts several integers or arrays |
Check_And_ReallocateIfNecessary | resizes an array if needed |
SetIncidentAngle | sets the wave vector with incident angles |
PrintNoBrackets | displays a tiny vector without brackets |
GetMemorySize | returns the memory used by complex structures |
FillZero | fills any structure with 0 |
square | returns the square of a number |
MakeLink | makes symbolic link between two files |
CopyFile | copies a file to another destination |
StringTokenize | splits a string into several small strings by specifying a delimiter |
StringTrim | Removes space and tab characters of a string |
DeleteSpaceAtExtremityOfString | space characters placed at the end or the beginning of strings are removed |
NumberToString | converts an integer into a string of four characters (e.g. '0023') |
Linspace | generates a vector of regularly spaced points |
GetHumanReadableMemory | returns a human-readable memory |
InitMontjoie | initialization of Montjoie static variables |
FinalizeMontjoie | end of Montjoie's session |
toInteger | conversion of a real number to an integer |
toDouble | conversion to a double |
GetPrecisionMachine | computes machine precision |
to_complex | converts a number into another one |
SetImaginaryI | tries to set y=i |
realpart | returns the real part |
imagpart | returns the imaginary part |
sign | returns the sign of a number |
Sqrt | returns the square root of a number |
abs_real | returns the absolute value of the real part of x |
abs_cplx | returns the modulus of a number |
isnan | returns true if the real or imaginary part of x is nan |
InitPrecisionConstants | initializes pi_wp, Iwp and epsilon_machine |
FindPrimeNumbers | computes prime numbers comprised between 2 and N |
FindPrimeDecomposition | computes decomposition of a number in prime numbers |
FindPrimeFactorization | increments an integer until to obtain a product of selected prime numbers |
FindTwoFactors | finds nx, ny such that n = nx*ny |
FindClosestPow2 | finds the power of 2 closest to a given integer |
FindClosestPow23 | finds the power of 2 and 3 closest to a given integer |
FindClosestPow235 | finds the power of 2 , 3 and 5 closest to a given integer |
FindClosestPow2357 | finds the power of 2, 3, 5 and 7 closest to a given integer |
ExtendIntegerInterval | increases an integer interval to contain a given number of elements |
EvaluateContinuedFraction | evaluates a number from its decomposition in continued fractions |
DecomposeContinuedFraction | computes decomposition in continued fractions |
GetNumeratorDenominator | evaluates numerator and denominator of a number expressed by its decomposition in continued fractions |
PrintFactorizationFailed | displays the status of factorization (when a direct solver is called) |
square
Syntax :
T square(x);
This function returns the square of x.
Example :
double ma_variable = 2.5; // for long argument, instead of x*x, you can use square double y = square(ma_variable);
Location :
Sort
Syntax :
void Sort(int& a, int& b); void Sort(TinyVector& a, TinyVector& b); void Sort(int& a, int& b, int& c); void Sort(int& a, int& b, int& c, int& d);
This function sorts the arguments (a, b, c and d for the last syntax) in ascending order.
Example :
int a = 3, b = 5, c = -4; // we should obtain a = -4, b = 3, c = 5 Sort(a, b, c);
Location :
Check_And_ReallocateIfNecessary
Syntax :
void Check_And_ReallocateIfNecessary(Vector& X, int n);
This function checks if the size of x is greater or equal to n. If the size of x is lower than n, the vector is resized to contain 1.5*n + 2 elements.
Example :
// in order to avoid too many reallocations // which is the case for PushBack // you can use this function Vector<double> num; int i = 0, n = 0; while (i != 4) { n++; Check_And_ReallocateIfNecessary(num, n); i = rand(); num(n-1) = i; } // once num is filled, you can resize it to its exact size num.Resize(n);
Location :
SetIncidentAngle
Syntax :
void SetIncidentAngle(const T& omega, TinyVector& kwave, const T& teta, const T& phi);
This function fills the wave vector kwave with the amplitude omega and the two incident angles teta and phi. This function is available in 2-D and 3-D. In 2-D the wave vector is equal to (ω cos θ, ω sin θ). In 3-D the wave vector is equal to (ω sin θ cos φ, ω sin θ sin φ, ω cos θ). In 2-D, this function is the reciprocal function of CartesianToPolar. In 3-D, this is the reciprocal function of CartesianToSpherical.
Example :
Real_wp omega = pi_wp; Real_wp teta = 0.5*pi_wp, phi = pi_wp; TinyVector<Real_wp, 3> kwave; // we compute the wave vector kwave from the wave number omega // and the two incident angles teta and phi SetIncidentAngle(omega, kwave, teta, phi); // if you want to find original omega(=r in spherical coordinates), teta and phi from kwave you can call CartesianToSpherical Real_wp theta0, phi0, cos_teta0, sin_teta0; CartesianToSpherical(kwave(0), kwave(1), kwave(2), theta0, phi0, cos_teta0, sin_teta0); // teta0 and phi0 should be equal to teta and phi
Location :
PrintNoBrackets
Syntax :
void PrintNoBrackets(ostream& out, const Vector& u); void PrintNoBrackets(ostream& out, const TinyVector& u);
This function prints vectors or tiny vectors without brackets and commas.
Example :
TinyVector<double, 5> x; // cout << x will display x with parenthesis and commas // PrintNoBrackets prints the values of x separated by spaces PrintNoBrackets(cout, x); // you can also display a vector of tiny vectors Vector<TinyVector<double, 5> > y; PrintNoBrackets(cout, y);
Location :
GetMemorySize
Syntax :
size_t GetMemorySize(const Vector<Vector>&); size_t GetMemorySize(const Vector<Matrix>&); size_t GetMemorySize(const Matrix<Matrix>&); size_t GetMemorySize(const Vector<Vector<Vector> >&); size_t GetMemorySize(const Vector<Vector<Matrix> >&); size_t GetMemorySize(const Vector<TinyVector<Vector> >&); size_t GetMemorySize(const TinyVector<Vector>&);
This function returns the memory (in bytes) used by composed structures such as vector of vectors, vector of matrices, etc.
Example :
Vector<Vector<double> > x; // if you want to know the size used to store x (all the values) size_t taille = GetMemorySize(x);
Location :
FillZero
Syntax :
template<class T> void FillZero(T&);
This function sets to 0 all the elements of the structure by calling function FillZero on each element of the structure.
Example :
// it works for any type Matrix<TinyVector<double, 3> > B(3, 3); FillZero(B);
Location :
MakeLink
Syntax :
void MakeLink(string source, string destination);
This function creates a symbolic link. source is the original file. destination is the created symbolic link pointing to the original file.
Example :
// disque_bis.mesh will be a link pointing to disque.mesh MakeLink("disque.mesh", "disque_bis.mesh");
Location :
CopyFile
Syntax :
void CopyFile(string source, string destination);
This function copies the source file to a destination file.
Example :
// disque_bis.mesh will be a hard copy of disque.mesh CopyFile("disque.mesh", "disque_bis.mesh");
Location :
StringTokenize
Syntax :
void StringTokenize(string phrase, Vector<string>& mots, string delim);
This function splits the string phrase into words contained in the variable mots. The words are assumed to be separated by one of the character contained in delim.
Example :
string phrase = "coucou me voila."; Vector<string> mots; // if the words are separated with space or tab character StringTokenize(phrase, mots, " \t"); // mots should contain "coucou", "me" and "voila." DISP(mots);
Location :
StringTrim
Syntax :
void StringTrim(string& s);
This function removes spaces contained in the string s.
Example :
string s = " Coucou Me Voila. "; StringTrim(s); // s should be equal to "CoucouMeVoila."
Location :
DeleteSpaceAtExtremityOfString
Syntax :
void DeleteSpaceAtExtremityOfString(string& s);
This function removes spaces located at the extremities of the string s.
Example :
string s = " Coucou Me Voila. "; DeleteSpaceAtExtremityOfString(s); // s should be equal to "Coucou Me Voila."
Location :
NumberToString
Syntax :
string NumberToString(int n); string NumberToString(int n, int nb_char);
This function returns the string corresponding to n with additional zeros on the left. By default, the returned string should contain 4 characters. The second argument is optional and specifies the number of characters of the returned string.
Example :
// s should be equal to "0023" string s = NumberToString(23); // s should be equal to "000023" s = NumberToString(23, 6);
Location :
Linspace
Syntax :
void Linspace(T a, T b, int N, Vector<T>& x);
This function computes a regular subdivision of the interval [a, b] with N points. The result is placed in vector x.
Example :
Vector<double> x; double a = -2.0, b = 3.0; // to produce a vector of 100 points uniformly distributed // in the interval [a, b] : Linspace(a, b, 100, x);
Location :
GetHumanReadableMemory
Syntax :
string GetHumanReadableMemory(size_t taille);
This function returns a string describing taille in human-readable string. The string contains the quantity with the unit (bytes, KiB, MiB, etc).
Example :
Vector<double> x; // if you want to display the memory used to store x // with readable units : cout << "Memory used to store x = " << GetHumanReadableMemory(x.GetMemorySize()) << endl;
Location :
InitMontjoie
Syntax :
void InitMontjoie(int argc, char** argv); void InitMontjoie(int argc, char** argv, int print_level);
This function should be called prior to any Montjoie function. Usually it is called in the main function of all targets. It calls MPI_Init (if MPI is used in the Makefile), and completes other initializations. The third argument is optional and sets the verbosity level of this function. By default, it will display the MPI_Init_thread level (MPI_THREAD_SINGLE, MPI_THREAD_FUNNELED, MPI_THREAD_SERIALIZED or MPI_THREAD_MULTIPLE).
Example :
// for any new main function #include "Montjoie.hxx" using namespace Montjoie; int main(int argc, char** argv) { // prior to any call to Montjoie functions : InitMontjoie should be called // this function calls MPI_Init InitMontjoie(argc, argv); // here you can call any function // finally, MPI_Finalize is called through FinalizeMontjoie return FinalizeMontjoie(); }
Location :
FinalizeMontjoie
Syntax :
int FinalizeMontjoie();
This function should be called before ending the program for a clean termination of all processes. It calls MPI_Finalize.
Example :
// for any new main function #include "Montjoie.hxx" using namespace Montjoie; int main(int argc, char** argv) { // prior to any call to Montjoie functions : InitMontjoie should be called // this function calls MPI_Init InitMontjoie(argc, argv); // here you can call any function // finally, MPI_Finalize is called through FinalizeMontjoie return FinalizeMontjoie(); }
Location :
toInteger
Syntax :
int toInteger(Real_wp x);
This function converts a float into an integer
Example :
Real_wp x = 2.8; // closest integer with round int n = toInteger(round(x));
Location :
toDouble
Syntax :
double toDouble(Real_wp x);
This function converts a float into a double.
Example :
Real_wp x = 2.8; // conversion in double precision double t = toDouble(x);
Location :
GetPrecisionMachine
Syntax :
T GetPrecisionMachine(T x);
This function computes the machine precision (i.e. the difference between 1 and its closest number in the representation of type T) of the type T. This function is used in InitMontjoie to compute epsilon_machine which is the precision for Real_wp.
Example :
// machine precision for long double long double x; long double eps = GetPrecisionMachine(x);
Location :
to_complex
Syntax :
to_complex(T x, T& y);
This function converts x in y. If x is complex and y real, it takes the real part of x.
Example :
Real_wp x = 0.3; Complex_wp y; // we copy x into y to_complex(x, y); // for the reverse operation x = real(y) to_complex(y, x);
Location :
SetImaginaryI
Syntax :
void SetImaginaryI(Real_wp& x); void SetImaginaryI(Complex_wp& x);
If the argument is complex, it is set to the pure imaginary number i. If the argument is real, it is set to one.
Example :
Real_wp x = 0.3; Complex_wp y; // y = i SetImaginaryI(y); // x = 1 SetImaginaryI(x);
Location :
realpart
Syntax :
T realpart(T x); T realpart(complex<T> x);
If the argument is complex, it returns the real part. If the argument is real, it returns the argument.
Example :
Real_wp x = 0.3; Complex_wp y(0.5, -0.2); // returns x since x is real Real_wp xr = realpart(x); // returns real(y) = 0.5 xr = realpart(y);
Location :
imagpart
Syntax :
T imagpart(T x); T imagpart(complex<T> x);
If the argument is complex, it returns the imaginary part. If the argument is real, it returns zero.
Example :
Real_wp x = 0.3; Complex_wp y(0.5, -0.2); // returns 0 since x is real Real_wp xi = imagpart(x); // returns imag(y) = -0.2 xr = imagpart(y);
Location :
sign
Syntax :
T sign(T x);
It returns -1 if x is nonpositive, 1 otherwise.
Example :
Real_wp x = 0.3; // returns 1.0 Real_wp s = sign(x); // returns -1.0 s = sign(-x);
Location :
Sqrt
Syntax :
T Sqrt(T x);
It returns the square root of x.
Example :
Real_wp x = 0.3; // equivalent to sqrt Real_wp s = Sqrt(x);
Location :
abs_real
Syntax :
T abs_real(T x);
It returns abs(real(x)) if x is complex, otherwise it returns abs(x) if x is real.
Example :
Complex_wp z(-0.3, 0.8); // the result should be 0.3 Real_wp s = abs_real(z);
Location :
abs_cplx
Syntax :
T abs_cplx(T x);
It is equivalent to function abs.
Example :
Complex_wp z(-0.3, 0.8); // the result should be the modulus of z Real_wp s = abs(z);
Location :
abs_cplx
Syntax :
bool isnan(complex<T> x);
It returns true if the real part or imaginary part is nan (not a number).
Example :
Complex_wp z(-0.3, nan); // the result should be true bool s = isnan(z);
Location :
InitPrecisionConstants
Syntax :
void InitPrecisionConstants();
It is called by InitMontjoie and sets pi_wp, Iwp and epsilon_machine. Since InitMontjoie is usually called, this function should not be called.
Example :
InitPrecisionConstants();
Location :
Clear
Syntax :
void Clear();
This method releases the memory used by the object.
Example :
VarRandomGenerator var_random; // selecting a random generator var_random.SetRandomGenerator(var_random.GSL_RANLUX); // then you use the random generator Real_wp y = var_random.GetRand(); // you can erase the object var_random.Clear();
Location :
RoundToSquare
Syntax :
void RoundToSquare(int& n, Real_wp alpha, int& nx, int& ny); void RoundToSquare(int& n, Real_wp alpha, Real_wp beta, int& nx, int& ny, int& nz);
This static method modifies the integer n such that n = nx*ny where ny is close to alpha*nx. The second syntax is designed for 3-D, by searching n = nx*ny*nz where ny is close to alpha*nx and nz close to beta*nx.
Example :
// you have a target (number of points) int N = 40000; // ratio between the number of points in y and in x Real_wp alpha = 2.0; // then you want to construct a 2-D grid with approximatively this number of points // the 2-D grid contains nx points along x and ny points along y // and you want to have ny = alpha*nx (for a square, alpha is equal to 1) // RoundToSquare will give you N such that N = nx*ny and ny approximatively equal to alpha*nx int nx, ny; VarRandomGenerator::RoundToSquare(N, alpha, nx, ny); cout << "Number of points in x, y : " << nx << " " << ny << endl; // for 3-D grid, you need to know a ratio beta such that nz \approx beta*nx Real_wp beta = 1.5; int nz; VarRandomGenerator::RoundToSquare(N, alpha, beta, nx, ny, nz);
Location :
GetRand
Syntax :
Real_wp GetRand();
This method returns a random number between 0 and 1.
Example :
VarRandomGenerator var_random; // selecting a random generator (otherwise it will use rand() function) var_random.SetRandomGenerator(var_random.GSL_RANLUX); // then you use the random generator Real_wp y = var_random.GetRand(); Real_wp z = var_random.GetRand();
Location :
GenerateRandomNumbers
Syntax :
void GenerateRandomNumbers(int N, VectReal_wp& R);
This method fills the vector R with N random numbers between 0 and 1.
Example :
VarRandomGenerator var_random; // selecting a random generator (otherwise it will use rand() function) var_random.SetRandomGenerator(var_random.GSL_RANLUX); // if you want to generate N random values int N = 100; VectReal_wp val; var_random.GenerateRandomNumbers(N, val);
Location :
ApplyRandomPermutation
Syntax :
void ApplyRandomPermutation(VectReal_wp& R);
This method shuffles values contained in R. R is replaced by R(P) where P is a random permutation.
Example :
VarRandomGenerator var_random; // selecting a random generator (otherwise it will use rand() function) var_random.SetRandomGenerator(var_random.GSL_RANLUX); // a vector x is constructed VectReal_wp x; Linspace(-1.0, 1.0, 100, x); // then if you want to permute randomly values contained in x var_random.ApplyRandomPermutation(x);
Location :
GenerateRandomPermutation
Syntax :
void GenerateRandomPermutation(int N, IVect& permut);
This method generates a random permutation, and stores this permutation in the output argument permut.
Example :
VarRandomGenerator var_random; // selecting a random generator (otherwise it will use rand() function) var_random.SetRandomGenerator(var_random.GSL_RANLUX); int N = 100; // then if you want to obtain a random permutation // i.e. numbers 0, 1, 2, .., N-1 shuffled IVect permut; var_random.GenerateRandomPermutation(N, permut);
Location :
SetRandomGenerator
Syntax :
void SetRandomGenerator(int type);
This method selects the random generator to be used. It also initializes the random generator. For the default random generator, it will call srand(time(NULL)). The different types are RANDOM, BIT_REVERSED, GSL_MT19937, GSL_RANLXD1, GSL_RANLXD2, GSL_RANLUX, GSL_RANLUX389, GSL_CMRG, GSL_MRG, GSL_TAUS, GSL_TAUS2, GSL_GFSR4. RANDOM corresponds to the function rand().
Example :
VarRandomGenerator var_random; // selecting a random generator (otherwise it will use rand() function) var_random.SetRandomGenerator(var_random.GSL_RANLUX); // once the random generator is selected, you can generate // as many random numbers as you want Real_wp x = var_random.GetRand();
Location :
SetInputRandomGenerator
Syntax :
void SetInputRandomGenerator(string nom);
This method selects the random generator to be used. It also initializes the random generator. For the default random generator, it will call srand(time(NULL)). The different types are RANDOM, BIT_REVERSED, GSL_MT19937, GSL_RANLXD1, GSL_RANLXD2, GSL_RANLUX, GSL_RANLUX389, GSL_CMRG, GSL_MRG, GSL_TAUS, GSL_TAUS2, GSL_GFSR4. RANDOM corresponds to the function rand().
Example :
VarRandomGenerator var_random; // selecting a random generator (otherwise it will use rand() function) var_random.SetInputRandomGenerator("GSL_RANLUX"); // once the random generator is selected, you can generate // as many random numbers as you want Real_wp x = var_random.GetRand();
Location :
GenerateBitReversedNumbers
Syntax :
void GenerateBitReversedNumbers(int n, VectReal_wp& x);
This method generates random numbers with the bit-reversed method.
Example :
VarRandomGenerator var_random; // selecting a random generator (otherwise it will use rand() function) var_random.SetInputRandomGenerator("GSL_RANLUX"); // if GenerateBitReversedNumbers is called, the selected random generator // is not called, but the random numbers are generated with the bit-reversed method int N = 100; VectReal_wp x; var_random.GenerateBitReversedNumbers(N, x);
Location :
FindPrimeNumbers
Syntax :
void FindPrimeNumbers(IVect& prime_numbers, int N);
This function computes all the prime numbers less than 2*N. These prime numbers are stored in the array prime_numbers. They are computed with a Eratosthene's sieve.
Example :
IVect prime_numbers; // for example, if you want to know all prime numbers p <= 2*N int N = 100; FindPrimeNumbers(prime_numbers, N); // the prime numbers are stored in the array prime_numbers;
Location :
FindPrimeDecomposition
Syntax :
void FindPrimeDecomposition(int n, IVect& prime_decomp, const IVect& prime_num); void FindPrimeDecomposition(int n, IVect& prime_decomp);
This function finds the decomposition of an integer n into prime numbers whose list is provided in the array prime_num. The result is stored in prime_decomp. If n is not factorizable with the prime numbers contained in prime_num, prime_decomp will store in the last element the irreducible factor. The second syntax computes the decomposition without providing prime numbers. In that case, prime numbers are automatically computed such that there is no irreducible factor.
Example :
IVect prime_numbers; // for example, if you want to know all prime numbers p <= 2*N int N = 100; FindPrimeNumbers(prime_numbers, N); // the prime numbers are stored in the array prime_numbers int p = 2*2*2*3*3*11*17*17; // p is factorized with prime_numbers that have been computed IVect prime_decomp; FindPrimeDecomposition(p, prime_decomp, prime_numbers); // here prime_decomp should give (2, 2, 2, 3, 3, 11, 17, 17) // the prime number is repeated with its power // here prime_numbers will be equal to (2, 3, 5, 7) FindPrimeNumbers(prime_numbers, 5); // p is not decomposable with only these prime numbers // p = 2*2*2*3*3* (11*17*17) // the irreducible factor is 11*17*17 = 3179 FindPrimeDecomposition(p, prime_decomp, prime_numbers); // prime_decomp should be equal to (2, 2, 2, 3, 3, 3179) // if you do not want to provide prime numbers // and compute the decomposition will all prime numbers // that can appear in the decomposition of p FindPrimeDecomposition(p, prime_decomp); // here prime_decomp should give (2, 2, 2, 3, 3, 11, 17, 17)
Location :
FindPrimeFactorizaion
Syntax :
void FindPrimeFactorization(int& n, IVect& prime_decomp, const IVect& prime_num); void FindPrimeFactorization(int& n, IVect& prime_decomp, int pmax);
This function increases n until that n is an exact product of prime numbers stored in prime_num. In the second syntax, pmax is the maximal prime number allowed.
Example :
IVect prime_numbers; // prime_numbers = (2, 3, 5, 7) FindPrimeNumbers(prime_numbers, 5); // 37 is a prime number that cannot be decomposed with 2, 3, 5 and 7 int p = 37; // you can ask for the nearest integer greater than p, such that // p_new = pow(2, n1)*pow(3, n2)*pow(5, n3)*pow(6, n4) int p_new = p; IVect prime_decomp; FindPrimeFactorization(p_new, prime_decomp, prime_numbers); // here it should find p_new = 40 and prime_decomp = (2, 2, 2, 5) // the result should be the same, the third argument if the maximum prime number, you allow in the decomposition of p_new p_new = 37; FindPrimeFactorization(p_new, prime_decomp, 7);
Location :
FindTwoFactors
Syntax :
void FindTwoFactors(int n, int& nx, int& ny);
This function computes nx and ny such that n = nx*ny (exactly) with nx and ny as close as possible.
Example :
int N = 20, nx ny; // here it should find nx = 4 and ny = 5 FindTwoFactors(N, nx, ny);
Location :
FindClosestPow2
Syntax :
int FindClosestPow2(int n);
This function returns the power of 2 that is the closest integer greater than n.
Example :
// approximate number of points int N = 1000; int nb_points = FindClosestPow2(N); // here you should have nb_points = 1024 = pow(2, 10);
Location :
FindClosestPow23
Syntax :
int FindClosestPow23(int n);
This function returns an integer 2i 3j that is the closest integer greater than n.
Example :
// approximate number of points int N = 900; int nb_points = FindClosestPow23(N); // here you should have nb_points = 972 = 2*2*3*3*3*3*3
Location :
FindClosestPow235
Syntax :
int FindClosestPow235(int n);
This function returns an integer 2i 3j 5k that is the closest integer greater than n.
Example :
// approximate number of points int N = 950; int nb_points = FindClosestPow235(N); // here you should have nb_points = 960 = 2*2*2*2*2*2*3*5
Location :
FindClosestPow2357
Syntax :
int FindClosestPow2357(int n);
This function returns an integer 2i 3j 5k 7m that is the closest integer greater than n.
Example :
// approximate number of points int N = 975; int nb_points = FindClosestPow2357(N); // here you should have nb_points = 980 = 2*2*5*7*7
Location :
ExtendIntegerInterval
Syntax :
void ExtendIntegerInterval(int& i0, int& i1, int Nt, int n);
This function extends the interval [i0, i1] such that they stay within the interval [0, Nt-1] and i1-i0+1 = N.
Example :
int i0 = 2, i1 = 10; int Nt = 20; int N = 12; // i0 is decremented and i1 incremented such that // the interval [i0, i1] contains N elements // i0 and i1 must stay between 0 and Nt-1 ExtendIntegerInterval(i0, i1, Nt, N);
Location :
EvaluateContinuedFraction
Syntax :
void EvaluateContinuedFraction(const Vector<int>& num, Real_wp& x);
This function evaluates the continued fraction stored in num. The result is placed in x.
Example :
// continued fraction a0 + 1 / (a1 + 1/(a2 + 1/...)) // coefficients a0, a1, a2 are stored in num Vector<int> num(3); num(0) = 2; num(1) = 4; num(2) = 3; // evaluation of the continued fraction Real_wp x; EvaluateContinuedFraction(num, x);
Location :
DecomposeContinuedFraction
Syntax :
void DecomposeContinuedFraction(Real_wp x, Real_wp epsilon, IVect& decomp);
This function constructs the continued fraction that is an approximation of x. The continued fraction is stopped as soon as the result is close enough to x (with a threshold epsilon).
Example :
// continued fraction of pi with an accuracy of 1e-6 IVect decomp; DecomposeContinuedFraction(pi_wp, Real_wp(1e-6), decomp);
Location :
GetNumeratorDenominator
Syntax :
void GetNumeratorDenominator(const IVect& decomp, int& num, int& denom);
This function computes the numerator and denominator of a continued fraction.
Example :
// continued fraction of pi with an accuracy of 1e-6 IVect decomp; DecomposeContinuedFraction(pi_wp, Real_wp(1e-6), decomp); // then you can compute p and q, such that p/q is the fraction approximating pi int p, q; GetNumeratorDenominator(decomp, p, q);
Location :
GetNumber
Syntax :
int GetNumber();
This method creates a new timer and returns the number associated with this timer.
Example :
// object containing a set of timers BasicTimer timer; // if you want to create a new timer int jeton = timer.GetNumber(); // then you can start the timer with the retrieved number timer.Start(jeton); // perform operations ... // then stop the timer timer.Stop(jeton); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds(jeton) << endl;
Location :
ReleaseNumber
Syntax :
void ReleaseNumber(int n);
This method destroys the timer n.
Example :
// object containing a set of timers BasicTimer timer; // if you want to create a new timer int jeton = timer.GetNumber(); // then you can start the timer with the retrieved number timer.Start(jeton); // perform operations ... // then stop the timer timer.Stop(jeton); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds(jeton) << endl; // once you no longer need this timer, you can release the number timer.ReleaseNumber(jeton);
Location :
ReserveNumber
Syntax :
void ReserveNumber(int n);
This method reserves the timer n. It may destroy an already existing number or be in concurrence with predefined timers (ALL, PROD, etc). So this functionality should be used cautiously.
Example :
// object containing a set of timers BasicTimer timer; // if you want to reserve some numbers timer.ReserveNumber(0); timer.ReserveNumber(1); // then you can start the timer with the reserved numbers timer.Start(0); // perform operations ... // then stop the timer timer.Stop(0); timer.Start(1); // perform operations ... // then stop the timer timer.Stop(1); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds(0) << endl; cout << "Time spent = " << timer.GetSeconds(1) << endl;
Location :
Reset
Syntax :
void Reset(int n);
This method resets the timer n, i.e. the elapsed time is set to 0 for this timer.
Example :
// object containing a set of timers BasicTimer timer; // if you want to create a new timer int jeton = timer.GetNumber(); // then you can start the timer with the retrieved number timer.Start(jeton); // perform operations ... // then stop the timer timer.Stop(jeton); // if you restart the timer, the time is incremented timer.Start(jeton); // other operations ... // final stop timer.Stop(jeton); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds(jeton) << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset(jeton);
Location :
Start
Syntax :
void Start(int n);
This method starts the timer n. The elapsed time is incremented.
Example :
// object containing a set of timers BasicTimer timer; // if you want to create a new timer int jeton = timer.GetNumber(); // then you can start the timer with the retrieved number timer.Start(jeton); // perform operations ... // then stop the timer timer.Stop(jeton); // if you start again the timer, the time is incremented timer.Start(jeton); // other operations ... // final stop timer.Stop(jeton); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds(jeton) << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset(jeton);
Location :
Stop
Syntax :
void Stop(int n);
This method stops the timer n. The elapsed time is updated only when this method is called.
Example :
// object containing a set of timers BasicTimer timer; // if you want to create a new timer int jeton = timer.GetNumber(); // then you can start the timer with the retrieved number timer.Start(jeton); // perform operations ... // then stop the timer timer.Stop(jeton); // if you start again the timer, the time is incremented timer.Start(jeton); // other operations ... // final stop timer.Stop(jeton); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds(jeton) << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset(jeton);
Location :
GetSeconds
Syntax :
double GetSeconds(int n);
This method returns the elapsed time between a call to Start and a call to Stop. The elapsed time is updated only when the method Stop is called.
Example :
// object containing a set of timers BasicTimer timer; // if you want to create a new timer int jeton = timer.GetNumber(); // then you can start the timer with the retrieved number timer.Start(jeton); // perform operations ... // then stop the timer timer.Stop(jeton); // if you start again the timer, the time is incremented timer.Start(jeton); // other operations ... // final stop timer.Stop(jeton); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds(jeton) << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset(jeton);
Location :
SetTimer
Syntax :
void SetTimer(int type);
This method changes the timing system to use. The choices are BASIC_TIMER, ACCURATE_TIMER and REAL_TIMER (they use respectively the classes BasicTimer, AccurateTimer and RealTimer). When this method is called, all previous timers are erased.
Example :
// object containing a set of timers MontjoieTimer timer; // the choice of the timer can be changed before creating timers timer.SetTimer(MontjoieTimer::ACCURATE_TIMER); // once the choice is made, it will use the asked timing system for all functions // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with the retrieved number timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl;
Location :
SetCommunicator
Syntax :
void SetCommunicator(MPI::Comm);
This method sets the MPI communicator to use when GetGlobalSeconds is called.
Example :
// object containing a set of timers MontjoieTimer timer; // the choice of the timer can be changed before creating timers timer.SetTimer(MontjoieTimer::ACCURATE_TIMER); // the communicator can be set anywhere before a call to GetGlobalSeconds timer.SetCommunicator(MPI::COMM_WORLD); // once the choice is made, it will use the asked timing system for all functions // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with the retrieved number timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // if you want to know the sum of elapsed times for all involved processors // you retrieve also the minimum and maximum double dt_loc, dt_sum, dt_min, dt_max; timer.GetGlobalSeconds("Product", dt_loc, dt_sum, dt_min, dt_max);
Location :
Reset
Syntax :
void Reset(int n); void Reset(const string& name);
This method resets the timer (identified either with a number or a string), i.e. the elapsed time is set to 0 for this timer.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with its name for example timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // if you start again the timer, the time is incremented timer.Start("Product"); // other operations ... // final stop timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset("Product");
Location :
Start
Syntax :
void Start(int n); void Start(const string& name);
This method starts the timer (identified either with a number or a string). The elapsed time is incremented.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with its name timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // if you start again the timer, the time is incremented timer.Start("Product"); // other operations ... // final stop timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset("Product");
Location :
Stop
Syntax :
void Stop(int n); void Stop(const string& name);
This method stops the timer (identified either with a number or a string). The elapsed time is updated only when this method is called.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with its name timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // if you start again the timer, the time is incremented timer.Start("Product"); // other operations ... // final stop timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset("Product");
Location :
GetSeconds
Syntax :
double GetSeconds(int n); double GetSeconds(const string& name);
This method returns the elapsed time between a call to Stop and a call to Start. The elapsed time is updated only when the method Stop is called.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with its name timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // if you start again the timer, the time is incremented timer.Start("Product"); // other operations ... // final stop timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset("Product");
Location :
SetMessage
Syntax :
void SetMessage(const string& name, const string& message);
This method creates a new timer associated with the string name. If a timer with the same name already exists, this method will just modify the message associated with the timer. A message is also required, this message is used when elapsed times are displayed with the method DisplayTime.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with its name timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // if you start again the timer, the time is incremented timer.Start("Product"); // other operations ... // final stop timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // if you want to reuse the timer, but restart from 0 seconds : timer.Reset("Product");
Location :
ReleaseName
Syntax :
void ReleaseName(const string& name);
This method destroys the timer associated with the string name.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create a new timer timer.SetMessage("Product", "Perform this operation"); // then you can start the timer with its name timer.Start("Product"); // perform operations ... // then stop the timer timer.Stop("Product"); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds("Product") << endl; // once you no longer need this timer, you can release its name timer.ReleaseName("Product");
Location :
ReserveName
Syntax :
void ReserveName(const string& name);
This method creates a new timer associated with the string name. If a timer already exists with this name, the program is aborted.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create timers without specifying a message timer.ReserveName("Product"); timer.ReserveName("Factorize"); // the message can be set afterwards timer.SetMessage("Product", "perform the matrix vector product"); // then you can start the wished timer timer.Start("Product"); // perform operations ... // then stop the timer timer.Stop("Product"); timer.Start("Factorize"); // perform operations ... // then stop the timer timer.Stop("Factorize"); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds("Product") << endl; cout << "Time spent = " << timer.GetSeconds("Factorize") << endl;
Location :
NameExists
Syntax :
bool NameExists(const string& name);
This method returns true if a timer with a given name already exists.
Example :
// object containing a set of timers MontjoieTimer timer; // if you want to create timers without specifying a message timer.ReserveName("Product"); timer.ReserveName("Factorize"); // you can test the existence of a timer if (timer.NameExists("Product")) cout << "Timer Product exists" << endl; // the message can be set afterwards timer.SetMessage("Product", "perform the matrix vector product"); // then you can start the wished timer timer.Start("Product"); // perform operations ... // then stop the timer timer.Stop("Product"); timer.Start("Factorize"); // perform operations ... // then stop the timer timer.Stop("Factorize"); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds("Product") << endl; cout << "Time spent = " << timer.GetSeconds("Factorize") << endl;
Location :
GetGlobalSeconds
Syntax :
void GetGlobalSeconds(int n, double& dt_loc, double& dt_sum, double& dt_min, double& dt_max); void GetGlobalSeconds(const string& name, double& dt_loc, double& dt_sum, double& dt_min, double& dt_max);
This method retrieves the local time elapsed (i.e. time elapsed on the current processor), the sum of elapsed times on all processors, the maximal time among all processors, and minimal time. The timer is given through its name or its identification number.
Example :
// object containing a set of timers MontjoieTimer timer; // the choice of the timer can be changed before creating timers timer.SetTimer(MontjoieTimer::ACCURATE_TIMER); // the communicator can be set anywhere before a call to GetGlobalSeconds timer.SetCommunicator(MPI::COMM_WORLD); // once the choice is made, it will use the asked timing system for all functions // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with the retrieved number timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // if you want to know the sum of elapsed times for all involved processors // you retrieve also the minimum and maximum double dt_loc, dt_sum, dt_min, dt_max; timer.GetGlobalSeconds("Product", dt_loc, dt_sum, dt_min, dt_max);
Location :
DisplayTime
Syntax :
void DisplayTime(int n); void DisplayTime(const string& name); void DisplayTime(int n, const string& message); void DisplayTime(const string& name, const string& message);
This method displays the elapsed time for a given timer (either by giving its number n or its name). The message can be given as a second argument, otherwise the displayed message will be the message given when SetMessage has been called.
Example :
// object containing a set of timers MontjoieTimer timer; // the choice of the timer can be changed before creating timers timer.SetTimer(MontjoieTimer::ACCURATE_TIMER); // the communicator can be set anywhere before a call to GetGlobalSeconds/DisplayTime timer.SetCommunicator(MPI::COMM_WORLD); // once the choice is made, it will use the asked timing system for all functions // if you want to create a new timer with a string timer.SetMessage("Product", "Perform the matrix vector product"); // then you can start the timer with the retrieved number timer.Start("Product"); // perform operations Matrix<double> A; Vector<double> x, b; Mlt(A, x, b); // then stop the timer timer.Stop("Product"); // and display the elapsed time (it will be the sum of Stop - Start) cout << "Time spent = " << timer.GetSeconds("Product") << endl; // for a nice display (and with parallel information) timer.DisplayTime("Product"); // if you want to specify another message timer.DisplayTime("Product", "Perform my product");
Location :
DisplayAll
Syntax :
void DisplayAll();
This method displays the elapsed times for all the created timers (through the methods ReserveName or SetMessage). The timers created with GetNumber/ReserveNumber are not displayed since they do not have a name.
Example :
// object containing a set of timers MontjoieTimer timer; // the communicator can be set anywhere before a call to DisplayAll timer.SetCommunicator(MPI::COMM_WORLD); // if you want to create timers without specifying a message timer.ReserveName("Product"); timer.ReserveName("Factorize"); // the message can be set afterwards timer.SetMessage("Product", "perform the matrix vector product"); // then you can start the wished timer timer.Start("Product"); // perform operations ... // then stop the timer timer.Stop("Product"); timer.Start("Factorize"); // perform operations ... // then stop the timer timer.Stop("Factorize"); // and display the elapsed time cout << "Time spent = " << timer.GetSeconds("Product") << endl; cout << "Time spent = " << timer.GetSeconds("Factorize") << endl; // if you want to do that automatically (display elapsed times for all timers) timer.DisplayAll();
Location :
PrintFactorizationFailed
Syntax
void PrintFactorizationFailed(int type, int ierr);
This function prints errors potentially encountered in the factorisation/resolution of sparse linear systems. If an error occured during the factorisation, the program is aborted.
Example :
// factorisation of a sparse matrix Matrix<double, General, ArrayRowSparse> A; All_MatrixLU<double> mat_lu; mat_lu.Factorize(A); // then checking if the factorisation was successful int ierr; int type = mat_lu.GetInfoFactorization(ierr); PrintFactorizationFailed(type, ierr);