User Manual

Introduction

This manual describes the main data-structures and function accessible at the user level to interact with the Chronos HPC library. The library is mainly designed for linear algebra problems, that is the solution of linear systems or eigenproblems, through iterative methods. This manual, in particular, describes:

  • the general settings of the library such as verbosity level, number of MPI processes and threads;
  • the DDMat object and its methods, which is used to store dense matrices;
  • the DSMat object and its methods, which is used to store sparse matrices;
  • the Preconditioner class, which implements various preconditioning techniques such as adaptive Factored Sparse Approximate Inverse (aFSAI) and adaptive Algebraic MultiGrid (aAMG);
  • the iterative methods for the solution of linear systems and eigenproblems available from Chronos;

Moreover, the two sample drivers that are available opensource are explained in detail.

Types, Classes and Methods in Chronos

Chronos general settings

Chronos general settings are stored in the global variable Chronos which is an instance of the class ChronosEnvironment made accessible from the header Chronos.h.

While using Chronos on distributed memory sysyems, first of all, it is necessary to initialize the MPI library. Since Chronos supports both threads and MPI processes, MPI initialization should be done with command:

  • MPI_Init_thread using option MPI_THREAD_FUNNELED.

After the MPI initialization, before using any other function provided by Chronos, it is MANDATORY to initialize the Chronos framework, otherwise execution is not possible and the program behavior is undefined. The initialization takes place by simply invoking the init method in Chronos:

  • Chronos.init(): this method initializes all the required internal variables of Chronos and returns true if the initialization is successful or false otherwise.

Then, the general environment should be set-up for the run in the following way:

  • Chronos.Set_inpNthreads( type_OMP_iReg nthreads ): this is used to set the number of threads used in the multi-threaded execution regions;
  • Chronos.Set_Verbosity( int level ): this is used to set the level of information dumped by the library. level may range from 0 to 3, with 0 meaning that no information is displayed. The higher the verbosity level, the higher is the overhead introduced in execution;

The types defined in Chronos to store integers and reals

Chronos defines its specific types to store integers and reals. The predefined C++ types associated to each Chronos type is set at the time of compilation and cannot be modified in the binaries.

The Chronos types, that the user need to access, are the following:

  • iReg: regular integer type, which typically corresponds to int;
  • iExt: extended integer type, which typically corresponds to int, but for specific installations can be extended to long int;
  • iGlo: global integer type, which is used to store global indices and typically corresponds to long int;
  • rExt: extended real type, which typically corresponds to double;

The class VEC

Chronos provides the VEC template class to manage one dimensional arrays of data. This class is a light-weight implementation of the standard-library vector and provides the users with all the most commonly used constructors and functions such as size, data, resize, assign, etc.

The class CSRMAT

The Compressed Sparse Row (CSR) format is the Chronos standard format for shared sparse matrices and their management is demanded to CSRMAT class.

All the CSRMAT members are {\tt private}. The most relevant for the CSR format definition are:

  • iRegExt nrows
    • Description: number of rows.
  • iRegExt ncols
    • Description: number of columns.
  • iExt nterm

    • Description: number of stored entries.
  • VEC vec_coef

    • Description: real array containing the values of matrix entries.
  • VEC vec_ja
    • Description: integer array containing the column indexes of the matrix entries.
  • VEC vec_iat
    • Description: integer array containing the pointers to the beginning of each row in the arrays vec_coef and vec_ja.

The public methods of the CSRMAT class are:

  • CSRMAT()
    • Description: default constructor for an empty object.
  • CSRMAT(const iRegExt input_nrows, const iRegExt input_ncols, const iExt input_nterm)
    • Description: constructor with storage allocation.
    • Parameters: [in] input_nrows: number of rows; [in] input_ncols: number of columns; [in] input_nterm: number of stored entries.
  • ~CSRMAT()
    • Description: default destructor.
  • void resize ( const iRegExt new_nrows, const iRegExt new_ncols, const iExt new_nterm )
    • Description: resizes the objects.
    • Parameters: [in] new_nrows: number of rows; [in] new_ncols: number of columns; [in] new_nterm: number of stored entries.
  • iRegExt get_nrows()
    • Description: retrieves the number of rows.
    • return: number of rows.
  • iRegExt get_ncols()
    • Description: retrieves the number of columns.
    • return: number of columns.
  • iRegExt get_nterm()
    • Description: retrieves the number of stored terms.
    • return: number of stored terms.
  • rExt* get_ptr_coef_data()
    • Description: retrieves the pointer to vec_coef data.
    • return: pointer to vec_coef data.
  • iRegExt* get_ptr_ja_data()
    • Description: retrieves the pointer to vec_ja data.
    • return: pointer to vec_ja data.
  • iExt* get_ptr_iat_data()
    • Description: retrieves the pointer to vec_iat data.
    • return: pointer to vec_iat data.

The template function allows to define: CSRMAT<iReg>, CSRMAT<iExt>, VEC<CSRMAT<iReg> >, VEC<CSRMAT<iReg> >, VEC<CSRMAT<iExt> >, VEC<CSRMAT<iExt> > and VEC<VEC<CSRMAT<iReg> > >.

The class CommGraph

In Chronos, communications between MPI ranks are handled efficiently by the CommGraph class. The most relevant private member are:

  • iReg nProc
    • Description: number of MPI processes involved in the communications.
  • iReg ncon_L
    • Description: number of MPI processes communicating with the calling process whose rank indexes are smaller than that of the calling process.
  • iReg ncon_R
    • Description: number of MPI processes communicating with the calling process whose rank indexes are larger than that of the calling process.
  • VEC<iReg> vec_comList
    • Description: integer array containing the rank indexes of the MPI processes involved in the communications.
  • iReg* proc_L
    • Description: integer array containing the rank indexes of the MPI processes whose rank indexes are smaller than that of the calling process, proc_L = vec_comList.data().
  • iReg* proc_R
    • Description: integer array containing the rank indexes of the MPI processes whose rank indexes are larger than that of the calling process, proc_R = vec_comList.data() + ncon_L + 1.

The only public member is:

  • iReg* myid
    • Description: rank index of the calling process, myid = vec_comList.data() + ncon_L.

The public methods of the CommGraph class are:

  • CommGraph()
    • Description: default constructor for an empty object.
  • CommGraph(const iReg input_ncon_L, const iReg input_ncon_R)
    • Description: constructor with storage allocation.
    • Parameters: [in] input_ncon_L: number of MPI processes whose rank < calling process; [in] input_ncon_R: number of MPI processes whose rank > calling process.
  • ~CommGraph()
    • Description: default destructor.
  • void resize ( const iReg new_ncon_L, const iReg new_ncon_R )
    • Description: resizes the objects.\
    • Parameters: [in] new_ncon_L: number of MPI processes whose rank < calling process; [in] new_ncon_R: number of MPI processes whose rank > calling process.
  • void mk_OutCommGraph ( CommGraph & __restrict__ OutCommGraph )
    • Description: creates the outgoing communication graph (considering ‘this’ as the incoming communication graph).
    • Parameters: [out] OutCommGraph: outgoing communication graph.
  • iReg get_nProc()
    • Description: retrieves the number of rows.
    • return: number of rows.
  • iReg get_ncon_L()
    • Description: retrieves the number of MPI processes whose rank < calling process.
    • return: ncon_L value.
  • iReg get_ncon_R()
    • Description: retrieves the number of MPI processes whose rank > calling process.
    • return: ncon_R value.
  • iReg* get_ptr_comList()
    • Description: retrieves the pointer to vec_comList data.
    • return: pointer to vec_comList data.
  • iReg* get_proc_L()
    • Description: retrieves proc_L pointer.
    • return: proc_L.
  • iReg* get_proc_R()
    • Description: retrieves proc_R pointer.
    • return: proc_R.

The class DDMat

In Chronos, the Distributed Dense Matrices are managed by the DDMat class. The storage schemes require the matrix to be subdivided into \(n_p\) horizontal stripes of consecutive rows, where \(n_p\) is the number of active MPI processes. Each stripe is stored row-wisely among the process to guarantee better access in memory during multiplication operations. DDMat results useful for linear systems with multiple right-hand-sides and eigenproblems, where distributed vectors are stored as one-column DDMat.

All the DDMat members are private. The most relevant for the storage scheme definition are:

  • iReg nrows
    • Description: number of rows (in the stripe).
  • iReg ncols
    • Description: number of columns.
  • VEC<rExt> coef
    • Description: stripe entries, stored row-wisely.

The public methods of the DDMat class are:

  • DDMat()
    • Description: default constructor for an empty object.
  • DDMat(const iReg input_nrows, const iReg input_ncols)
    • Description: constructor with storage allocation.
    • Parameters: [in] input_nrows: number of rows (in the stripe); [in] input_ncols: number of columns.
  • ~DDMat()
    • Description: default destructor.
  • void resize(const iReg new_nrows, const iReg new_ncols)
    • Description: resizes the objects.
    • Parameters: [in] new_nrows: new number of rows (in the stripe); [in] new_ncols: new number of columns.
  • DDMat& operator= ( DDMat& other )
    • Description: copy operator.
    • return: copied DDMat.
  • void Redist(const type_MPI_iReg rank, const type_MPI_iReg size, const VEC<iGlo> &vec_part, const VEC<iGlo> &vec_iperm, DDMat& mat_AP)
    • Description: redistribute the DDMat instance according to the input partition vec_part.
    • Parameters: [in] rank: rank id of this process; [in] size: number of MPI processes; [in] vec_part: the array containing for each unknown of the matrix, the ID of the final partition; [in] vec_iperm: the array containing the global numbering of each unknown; [out] mat_AP: the partitioned matrix.
  • void LocalPermute(const bool Inverse_Flag, const VEC<iReg> &vec_perm_in, const VEC<iReg> &vec_iperm_in, DDMat &mat_AP)
    • Description: applies a local permutation to the unknowns of this stripe of the DDMat instance.
    • Parameters: [in] Inverse_Flag: flag to apply inverse local permutation; [in] vec_perm_in: local permutation; [in] vec_iperm_in: local inverse permutation; [out] mat_AP: locally reordered matrix.
  • iReg get_nrows()
    • Description: retrieves the number of rows (in the stripe).
    • return: number of rows.
  • iReg get_ncols()
    • Description: retrieves the number of columns.
    • return: number of columns.
  • rExt* get_ptr_coef_data()
    • Description: retrieves the pointer to coef data.
    • return: pointer to coef data.
  • void set_all_coef(const rExt input_value)
    • Description: sets all coef entries equal to a given input value.
    • Parameters: [in] input_value: given value for setting coef entries.
  • set_col_coef(const iReg input_icol, const rExt input_value)
    • Description: sets the entries of a {\tt DDMat} column equal to a given input value.\
    • Parameters: [in] input_value: given value for setting entries of a DDMat column.

The class DSMat

In Chronos, the Distributed Sparse Matrices are managed by the DSMat class. As for the DDMat class, the storage schemes require the matrix to be subdivided into \(n_p\) horizontal stripes of consecutive rows, where \(n_p\) is the number of active MPI processes. Each stripe is subdivided into Compressed Sparse Row (CSR) matrices. The DSMat storage scheme is very effective in both the preconditioner computation and the SpMV product because it allows a large superposition between communication and computation.
The most relevant public members of the DSMat class are:

  • CommGraph *MatGraph
    • Description: Incoming (for SpMV and SpMM products) Communication graph. This object is defined in the CommGraph class described above.
  • CommGraph *OutMatGraph
    • Description: Outgoing (for SpMV and SpMM products) Communication graph. If the matrix is symmetric OutMatGraph just points to the same MatGraph data. This object is defined in the CommGraph class described above.
  • VEC<iReg> *ptDom
    • Description: vector of pointers to the compacted column indices of each incoming block. That is, ptDom[ib+1] - ptDom[ib] gives the number of columns of the incoming block ib (ib block on this row of blocks).
  • VEC<iReg> *OutptDom
    • Description: vector of pointers to the compacted row indices of each outgoing block. That is, OutptDom[ib+1] - OutptDom[ib] gives the number of rows of the outgoing block ib (ib block on the column of blocks including this rank). If the matrix is symmetric OutptDom just points to the same ptDom data
  • VEC /> csr
    • Description: vector of CSR.
  • CSRMAT<iReg>* matR
    • Description: array to the Right CSR. Each of this will be a reference to the previous vector of CSR (no memory associated).
  • CSRMAT<iReg>* matL
    • Description: array to the left CSR: each of this will be a reference to a CSR mat (no memory associated).
  • CSRMAT<iReg>* matC
    • Description: pointer to the diagonal entry of the vector csr (no memory associated).

The public methods of the DSMat class are:

  • DSMat()
    • Description: default constructor for an empty object.
  • ~DSMat()
    • Description: default destructor.
  • DSMat& operator= ( DSMat& other )
    • Description: copy operator.
    • return: copied DSMat.
  • void PartKway( const type_MPI_iReg rank, const type_MPI_iReg size, VEC<iGlo> &vec_part, VEC<iGlo> &vec_iperm)
    • Description: Redistributes the entries currently stored locally into this process to all the other process according to a partition generated by metis.
    • Parameters: [in] rank: rank id of this process; [in] size: number of MPI processes; [out] vec_part: the array containing for each unknown of the matrix, the ID of the final partition; [out] vec_iperm: the array containing the global numbering of each unknown.
  • void Redist( const type_MPI_iReg rank, const type_MPI_iReg size, const VEC<iGlo> &vec_part, const VEC<iGlo> &vec_iperm, DSMat &mat_AP)
    • Description: Creates through metis a partition of the DSMat instance.
    • Parameters: [in] rank: rank id of this process; [in] size: number of MPI processes; [in] vec_part: the array containing for each unknown of the matrix, the ID of the final partition; [in] vec_iperm: the array containing the global numbering of each unknown; [out mat_AP]: the partitioned matrix.
  • void LocalRCM( VEC &vec_perm, VEC &vec_iperm)
    • Description: Computes a local permutation (i.e. each subdomain has a local independent ordering) on an input DSMat using the reverse Cuthill-McKee (RCM) algorithm.
    • Parameters: [out] vec_perm: array containing the local permutation; [out] vec_iperm: array containing the local inverse permutation.
  • void LocalPermute(const type_MPI_iReg rank, const type_MPI_iReg size, const bool Inverse_Flag, const VEC<iReg> vec_perm_in, const VEC<iReg> vec_iperm_in, DSMat &mat_AP)
    • Description: Performs a local permutation (i.e. each subdomain has a local independent ordering) on an input DSMat. An input flag is provided to decide whether to produce a direct or an inverse ordering.
    • Parameters: [in] rank: rank id of this process; [in] size: number of MPI processes; [in] vec_perm_in: array containing the local permutation; [in] vec_iperm_in: array containing the local inverse permutation. [out mat_AP]: the permuted matrix.
  • void mk_OutptDom()
    • Description: Creates OutptDom from the Outgoing communication graph.
  • void set_sym_flag(bool input_sym_flag)
    • Description: Set the symmetric flag.
  • bool get_sym_flag()
    • Description: Get the symmetric flag.
    • return: the symmetric flag.
  • iReg get_nrows()
    • Description: Get nrows.
    • return: the number of rows in the process stripe.
  • iReg get_nterm()
    • Description: Get local nterm.
    • return: the number of terms in the process stripe.
  • iReg* get_ptr_ptDom_data()
    • Description: Get pointer to ptDom.
    • return: the pointer to {\tt ptDom}.
  • iReg* get_ptr_OutptDom_data()
    • Description: Get OutptDom.
    • return: the pointer to OutptDom.
  • iGlo Get_Globalnterm()
    • Description: Evaluates the total number of nterm of the distributed csr matrices.
    • return: the global number of terms.
  • iGlo Get_Globalnrows()
    • Description: Evaluates the total number of rows of the distributed matrix.
    • return: the global number of rows.
  • void fPrintAll_PARDDMAT_ASCII(const string &filename)
    • Description: Prints the distributed dense matrix in ASCII filename (overwrite if it exists).
    • Parameters: [in] filename: file name in which the matrix is printed.
  • void fPrintAll_PARDDMAT_BIN(const string &filename)
    • Description: Prints the distributed dense matrix in BINARY filename (overwrite if it exists).
    • Parameters: [in] filename: file name in which the matrix is printed.
  • void MxV(DDMat& x, DDMat& b, bool Barrier_flag)
    • Description: Computes the Matrix by Vector Product.
    • Parameters: [in] x vector that is multiplied by the matrix; [out] b vector where result is stored;\ [in] Barrier_flag if true global Barrier is called before return.
  • void MxM(DSMat& B, DSMat& C)
    • Description: Computes the Matrix by Matrix Product \(A \times B = C\) (with A equal to this).
    • Parameters: [in] B matrix that multiplies the matrix. [out] C matrix where result is stored.
  • void TreatBC()
    • Description: This function modifies a DSMat instance by detecting rows/columns associated to Dirichlet boundary conditions and setting to those rows/columns a diagonal entry equal to the maximum diagonal entry of the whole matrix. In such a way the minimum eigenvalue of the matrix is not erroneously associated to the arbitrary value given to the diagonal entries of boundary dofs. Detection is carried out by performing one Jacobi smoothing step on the constant vector and setting as BC those dofs whose value is zero.
  • void Transpose (DSMat& MatT)
    • Description: Transpose the matrix.
    • Parameters: [out] MatT Tranposed Matrix.

The class MatrixProd

This class is the base class for every sparse linear operator, that is for any object that maps a vector (or a dense matrix) into another vector (or dense matrix with the same number of columns). In particular, the following classes, that will be described below are derived from the MatrixProd base class:

  • DSMat: the class used to store a sparse matrix;
  • Preconditioner: the class used to store the available preconditioners (see below);
  • MatList: a list of MatrixProd instances that is used to concatenate multiple linear operators.

The interested user can code its own preconditioner by deriving this specific base class and the use all the iterative methods and linear operators available in Chronos.

The members of this class are:

  • bool Prepared_flag: flag indicating if the linear operator is ready or not for the MxV product. In fact, some operators need to be “prepared” to carry out the MxV product, that is they may need to undergo a transformation or allocate some temporary buffer. The true flag means that the operator is ready;
  • iReg nRHS_prep: this member indicates the number of rhs the operator is prepared for. If a MxV product is needed with a larger number of rhs, then eventual buffers need to be reallocated.

The methods provided by the class are:

  • void copy_MatrixProd ( MatrixProd& other ): creates a copy of the instance. This method is superposed with the “=” operator;
  • bool get_Prepared_flag(): returns the value of Prepared_flag;
  • iReg get_nRHS_prep(): returns the value of nRHS_prep;
    • virtual iReg get_nrows()*: returns the local number or rows;
  • virtual iGlo Get_Globalnrows(): returns the global number of rows;
  • virtual void MxV(DDMat& x, DDMat& b, bool Barrier_flag): Computes the Matrix by vector Product “b = A*x”. If {\tt Barrier_flag} is true, a global Barrier is called before return;
  • virtual ~MatrixProd(): deletes the instance of this class.

The class Preconditioner

This class, which is derived from MatrixProd serves as a base class for deriving preconditioners. Beyond all the methods inherited from MatrixProd, Preconditioner provides the following ones:

  • virtual void fPrintAll_ASCII(const string &filename): prints the preconditioner in ASCII format;
  • virtual void MxM ( DSMat& Mat, DSMat& PrecMat ): Applies the preconditioner to a (DSMat) matrix;
  • virtual void Compute( DSMat& mat_A , DDMat *const V0 = nullptr ): Computes the preconditioner for the input matrix mat_A, optionally using an input test space {\tt V0};
  • virtual iExt Get_nnz (): Returns an estimate of the number of non-zeroes of the preconditioner stored on this process;
  • virtual iExt Get_FlopEst (): Returns an estimate of the number of FLOPs needed to perform the preconditioner application on this process.

Currently, Chronos provides the user with 3 preconditioning methods:

  • Jacobi: the simplest possible preconditioner, that is just a scaling by the diagonal elements of the system matrix;
  • aFSAI: the adaptive Factored Sparse Approximate Inverse preconditioner.
  • aAMG: the adaptive AMG preconditioner.

The pubblic functions of these 3 preconditioners are those inherited from the base class Preconditioner.

The set-up of aFSAI and aAMG can be customized by the user to increase performance for specific problems using the following methods:

  • iReg aFSAI::Set_default(): sets the aFSAI preconditioner to its default;
  • iReg aFSAI::Set_nstep(iReg nstep): sets to nstep the number of steps used to build aFSAI;
  • iReg aFSAI::Set_step_size(iReg step_size): sets to step_size the number of entries added for each row at each step of the adaptive procedure;
  • iReg aFSAI::Set_epsilon(rExt epsilon): sets to epsilon the exit tolerance used in the adaptive procedure;
  • iReg aAMG::Set_mech(): sets the aAMG preconditioner to its default for mechanical problems;
  • iReg aAMG::Set_mech_SRQCG(iReg niter_SRQCG): sets to niter_in the number of iterations of SRQGC for the construction of the test space;
  • iReg aAMG::Set_mech_smoother_type(iReg smoother_type): sets to smoother_type the smoother used in aAMG;
  • iReg aAMG::Set_mech_prol_smoothing(bool prol_smoothing): sets on or off the prolongation smoothing in aAMG;
  • iReg aAMG::Set_cfd(): sets the aAMG preconditioner to its default for CFD problems;
  • iReg aAMG::Set_cfd_smoother_type(iReg smoother_type): sets to smoother_type the smoother used in aAMG;
  • iReg aAMG::Set_cfd_prolongation_type(iReg prolongation_type): sets to prolongation_type the prolongation used in aAMG;

The above methods return 0 in case of success.

Iterative methods available in Chronos

Chronos provides iterative methods to solve both linear systems and eigenproblems. There are two base classes, one for the iterative methods for linear systems and the other for the iterative methods for eigenproblems.

Iterative methods for the solution to linear systems

The base class for these iterative methods is named LinSolver and has the following protected members:

  • iReg maxITER: maximum number of iterations;
  • rExt exitTOL: exit tolerance on the relative iterative residual;
  • iReg initial_sol: flag to set the initial solution. 0 \(\Rightarrow\) sol0 = Prec x rhs, 1 \(\Rightarrow\) sol0 = 0, 2 \(\Rightarrow\) initial solution is available from input;
  • bool printConvProfile: flag for Convergence Profile printing;
  • DDMat RES: residual array;
  • iReg ITER: total number of iterations;
  • VEC normRES0: euclidean norm of the initial residual(s);
  • VEC normRES: euclidean norm of the final iterative residual(s);
  • VEC normRESreal: euclidean norm of the final real residual(s);
  • VEC normRHS: euclidean norm of the right hand side(s);

It also has the following public functions:

  • virtual ~LinSolver(): Deletes the object;
  • virtual void Set_Solver(MatrixProd& MAT, DDMat& RHS, MatrixProd& PREC): sets the solver for monolithic preconditioning;
  • virtual void Set_Solver(MatrixProd& MAT, DDMat& RHS, MatrixProd& PREC_L, MatrixProd& PREC_R): sets the solver for split preconditioning;
  • virtual void Solve(MatrixProd& MAT, DDMat& RHS, DDMat& SOL, MatrixProd& PREC, const string filename): solves the linear system with monolithic preconditioning;
  • virtual void Solve(MatrixProd& MAT, DDMat& RHS, DDMat& SOL, MatrixProd& PREC_L, MatrixProd& PREC_R, const string filename): solves the linear system with split preconditioning;
  • void Set_maxITER(const iReg input_maxITER): sets the exit tolerance;
  • void Set_exitTOL(const rExt input_exitTOL): sets the exit tolerance;
  • void Set_initial_sol(const iReg input_initial_sol): sets the flag for initial solution;
  • Set_printConvProfile(const bool input_printConvProfile): Sets the flag for Convergence Profile printing;
  • iReg get_ITER(): Returns the number of iterations performed in previous solution;
  • const rExt* get_ptr_normRES(): Returns an handle to the Euclidean norm of the final iterative residual for each RHS.

Chronos presently offers two iterative linear solvers:

  • PCG: the preconditioned conjugate gradient, suitable only for Symmetric Positive Definite linear systems;
  • BiCGstab: the preconditioned biconjugate gradient stabilez, which can be used for general linear systems;

Iterative methods for the solution to eigenproblems

The base class for these iterative methods is named EigenSolver and has the following public functions:

  • virtual ~EigenSolver(): Deletes the object;
  • virtual void Set_Solver( MatrixProd& MAT , MatrixProd *const PREC = nullptr ): sets the solver allocating scratches for computation and preparing matrices;
  • virtual void Solve( MatrixProd& MAT, const DDMat& V0, VEC& EigValues, DDMat& EigVectors , MatrixProd *const PREC = nullptr ): solves the eigenvalue problem;

  • void set_maxITER(const iReg input_maxITER): Sets the maximum number of iterations;

  • void set_exitTOL(const rExt input_exitTOL): Sets the number of eigenpairs seeked in the solve stage;
  • bool get_solAVAIL(): returns true if a solution to the eigenproblem is available;
  • iReg get_maxITER(): returns the value of maxITER;
  • iReg get_exitTOL(): returns the value of exitTOL;
  • iReg get_nEigs(): returns the number of eigenpairs seeked in the solve stage;
  • iReg get_ITER(): returns the number of iterations performed in previous solution;
  • const VEC* get_ptr_normRES_eigvals(): returns an handle to the Euclidean norm of the final iterative residual for each eigenvalue;
  • const VEC* get_ptr_normRES_eigvect(): returns an handle to the Euclidean norm of the final iterative residual for each eigenvector;

Chronos presently offers two iterative eigensolvers:

  • PowMeth: the power method;
  • SRQCG: the Simultaneous Rayleigh Quotien minimization through Conjugate Gradient;

SRQCG offers the additional public methods:

  • void set_orth_freq (const iReg orth_freq_in): sets orthogonalization frequency;
  • void set_ritz_freq (const iReg ritz_freq_in): sets the Ritz projection frequency;
  • void get_orth_freq (): gets orthogonalization frequency;
  • void get_ritz_freq (): gets the Ritz projection frequency;

Utilities

Together with the Demo drivers, Chronos is equipped with some utilities to set-up the input files for the drivers. The following utilities can be used to convert sparse and dense matrices:

  • ascii2bin: a simple program to convert a matrix file from ascii to binary format;
  • bin2ascii: a simple program to convert a matrix file from binary to ascii format;
  • Tspace_ascii2bin: a simple program to convert a test space file from ascii to binary format;
  • Tspace_bin2ascii: a simple program to convert a test space file from binary to ascii format.

Users can call these utilities to convert their files into appropriate input files for Driver_Binary. The test space is stored in a DDMat as any other vector or set of vectors, so Tspace_ascii2bin and Tspace_bin2ascii can also be used to convert the input right-hand sides.

The utility:

  • MK_RBM

is used to transform 2D or 3D coordinates into a test space suitable for structural problems by computing the Rigid Body Modes (RBMs) of the model. This utility takes as argument an input file containing the following parameters:

  • nthreads: the number of threads used in creating RBM;
  • n_dim: the problem dimension, 2 or 3;
  • coord_fname: a string containing the name of the coordinate file (ascii);
  • RBM_fname: a string containing the name of the output RBMs file (binary);