BlockIPPlatform
 All Classes Files Functions Variables
Classes | Public Member Functions | Public Attributes | Private Member Functions | List of all members
MatrixBlockIP Class Reference

Class for manipulating matrices, and interfacing SparseChol. More...

#include <MatrixBlockIP.h>

Collaboration diagram for MatrixBlockIP:
Collaboration graph
[legend]

Classes

struct  Order_ija
 Auxiliary struct for sorting matrices in ija format. More...
 
struct  Order_vector
 Auxiliary struct for sorting vectors. More...
 

Public Member Functions

 MatrixBlockIP (int blocks=1)
 Constructor.
 
 MatrixBlockIP (MatrixBlockIP *mbip)
 Copy constructor, does not copy SparseChol.
 
 ~MatrixBlockIP ()
 Destructor.
 
void copy (MatrixBlockIP *mbip)
 Copy, does not copy SparseChol.
 
void reset ()
 Comes back to the initial state.
 
void restore ()
 Comes back to the state after create the matrix.
 
void create_general_matrix_row_wise (int m, int n, int nz, int *&inirowa, int *&icola, double *&a)
 Native method to create a general row-wise packed matrix.
 
void create_general_matrix_column_wise (int m, int n, int nz, int *&inicola, int *&irowa, double *&a)
 Native method to create a general column-wise packed matrix.
 
void create_network_matrix (int num_arcs, int num_nodes, int *&src, int *&dst, bool oriented=true)
 Native method to create a network matrix.
 
void create_identity_matrix (int dim)
 Native method to create a identity matrix of dimension dim.
 
void create_idty_idty_matrix (int nrows)
 Native method to create a identity-identity matrix with two identities [I I].
 
void create_diagonal_matrix (int dim, double *&d)
 Native method to create a diagonal matrix D= diag(d) of dimension dim.
 
void create_diag_diag_matrix (int nrows, double *&d1, double *&d2)
 Native method to create a diagonal-diagonal matrix D= [D1 D2].
 
void create_general_matrix_format_ija (int m, int n, int nz, int *&row_index, int *&col_index, double *&values)
 Creates a general matrix in format ija.
 
void compute_full_matrix (int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[], int numActiveLnk=0, int *listActiveLnk=NULL)
 Builds a packed rowwise format structured matrix based on diagonal blocks and linking constraints blocks.
 
void analyze_D (int numBlocks, bool sameL, MatrixBlockIP L[])
 
void analyze_D_diagonal (int numBlocks, bool sameL, MatrixBlockIP L[])
 
void analyze_D_gen_sym_uptr (int numBlocks, bool sameL, MatrixBlockIP L[])
 
void compute_D (int numBlocks, bool sameL, MatrixBlockIP L[], bool isActive[], int iniTheta[], double Theta[])
 Compute D = Theta_(numBlocks+1) + sum{i..numBlocks} L_i*Theta_i*L_i'.
 
void compute_D_diagonal (int numBlocks, bool sameL, MatrixBlockIP L[], bool isActive[], int iniTheta[], double Theta[])
 Compute D = Theta_(numBlocks+1) + sum{i in 1..numBlocks} L_i*Theta_i*L_i' when D is diagonal.
 
void compute_D_gen_sym_uptr (int numBlocks, bool sameL, MatrixBlockIP L[], bool isActive[], int iniTheta[], double Theta[])
 Compute D = Theta_(numBlocks+1) + sum{i..numBlocks} L_i*Theta_i*L_i' when D is a symmetric general matrix.
 
void native_to_general (bool delete_native_format=false)
 Calculates and stores the matrix in packed rowwise format.
 
void ija_to_rowwise ()
 Calculates and stores the matrix in packed rowwise format.
 
void network_to_general ()
 Calculates and stores the network matrix (either oriented or nonoriented) in packed rowwise format.
 
void identity_to_general ()
 Calculates and stores the I matrix in packed rowwise format.
 
void diagonal_to_general ()
 Calculates and stores the diagonal D matrix in packed rowwise format.
 
void idty_idty_to_general ()
 Calculates and stores the matrix [I I] in packed rowwise format.
 
void diag_diag_to_general ()
 Calculates and stores the matrix [D1 D2] in packed rowwise format.
 
void network_to_ija_format ()
 Calculates and stores the network matrix in ija format. It considers both oriented and nonoriented cases.
 
void identity_to_ija_format ()
 Calculates and stores the identity I matrix in ija format.
 
void diagonal_to_ija_format ()
 Calculates and stores the diagonal D1 matrix in ija format.
 
void idty_idty_to_ija_format ()
 Calculates and stores the matrix [I I] in ija format.
 
void diag_diag_to_ija_format ()
 Calculates and stores the matrix [D1 D2] in ija format.
 
void order_matrix ()
 Order the matrix when has been created in row-wise or column-wise format.
 
void mul_Mv (double vout[], const double vin[])
 Matrix-vector product vout=M*vin (vout(m),vin(n)) (driver)
 
void mul_Mtv (double vout[], const double vin[])
 MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (driver)
 
void add_mul_Mv (double vout[], const double vin[])
 Add matrix-vector product vout += M*vin (vout(m),vin(n)) (driver)
 
void add_mul_Mtv (double vout[], const double vin[])
 Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (driver)
 
void exist_var_in_row (int column, bool appear[], double coefs[])
 Given a column determines if exists an element for each row of the matrix, in case of exist returns the value.
 
void add_new_column (int size, int irows[], double values[])
 Add a new column into the matrix.
 
void add_new_columns (int num_columns, int size[], int *irows[], double *values[])
 Add new columns into the matrix.
 
void change_columns_sign (int size, int columns[])
 Changes the sign of some columns.
 
void change_rows_sign (int size, int rows[])
 Change the sign of some rows.
 
void delete_rows (int size, int rows[])
 Delete some rows.
 
int get_column (int column, int *&irows, double *&values, bool invert_sign=false)
 Gets a column of the matrix.
 
void symbolic_fact_MMt (CHOL_SOLVER chslv=SPRSBLKLLT, int *prov_pfa=NULL, int *prov_ipfa=NULL, NUMBERING prov_pfa_numbering=NOT_COMPUTED)
 Interface routine to SparseChol symbolic factorization.
 
void numeric_fact_MMt (double *Theta, int i_k=0)
 Interface routine to SparseChol numeric factorization.
 
void numeric_solve_MMt (double *rhs, int i_k=0, WHO_PERMUTES whoperm=CHOLESKY)
 Interface routine to SparseChol numeric solve.
 
void symbolic_fact_M (CHOL_SOLVER chslv=SPRSBLKLLT, int *prov_pfa=NULL, int *prov_ipfa=NULL, NUMBERING prov_pfa_numbering=NOT_COMPUTED)
 Interface routine to SparseChol symbolic factorization.
 
void numeric_fact_M ()
 Interface routine to SparseChol numeric factorization.
 
void numeric_solve_M (double *rhs, WHO_PERMUTES whoperm=CHOLESKY)
 Interface routine to SparseChol numeric solve.
 
int get_pfa (int i)
 Interface routine to SparseChol get_pfa().
 
int get_ipfa (int i)
 Interface routine to SparseChol get_ipfa().
 
int get_maxlnz ()
 Interface routine to SparseChol get_maxlnz().
 
int get_maxfillin ()
 Interface routine to SparseChol get_maxfillin().
 
int get_njka ()
 Interface routine to SparseChol get_njka().
 
int get_num_zero_pivots ()
 Interface routine to SparseChol get_num_zero_pivots().
 
int get_num_semidef_matrix ()
 Interface routine to SparseChol get_semidef_matrix().
 
void print_matrix (bool ija=true, bool start_one=true)
 Print the matrix.
 
void print_matrix (ofstream &outfile, bool print_ija=true, bool start_one=true)
 Write the matrix into a file.
 
void print_vector (double *v, int size, string name)
 Print some positions and name of a vector.
 
void column_wise_to_row_wise_format ()
 Convert a matrix in column-wise format to row-wise format.
 
void row_wise_to_column_wise_format ()
 Convert a matrix in column-wise format to row-wise format.
 

Public Attributes

TYPE_MATRIX type
 Type of matrix.
 
int m
 Number of rows.
 
int n
 Number of columns.
 
bool ija
 True if irowa, icola and a are stored.
 
bool rowwise
 True if inirowa, icola and a are stored.
 
bool columnwise
 True if the matrix is stored in column-wise packed format, incompatible with rowwise.
 
int nz
 Number of non-zero elements.
 
double * a
 Value of non-zero elements.
 
int * inirowa
 Index to the first element of each row.
 
int * icola
 Column position for each element.
 
int * inicola
 Index to the first element of each column.
 
int * irowa
 Row position for each element.
 
TYPE_ORIENTATION type_orientation
 Type of orientation, by default, oriented.
 
int num_arcs
 Number of arcs.
 
int num_nodes
 Number of nodes.
 
int * src
 Source of each arc.
 
int * dst
 Destination of each arc.
 
double * d1
 d1 for DIAGONAL and 1st submatrix of DIAG_DIAG
 
double * d2
 d2 for 2nd submatrix of DIAG_DIAG
 
SparseChol chol
 Sparse Cholesky class.
 
int sizeL
 
bool made_symbfct_MMt
 
bool made_analyze_D
 
int num_blocks
 

Private Member Functions

void free_memory ()
 Free all the memory allocated by the application - not the user.
 
void mul_Mv_row_wise (double vout[], const double vin[])
 Matrix-vector product vout=M*vin (vout(m),vin(n)) (general rowwise matrix)
 
void mul_Mtv_row_wise (double vout[], const double vin[])
 MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (general rowwise matrix)
 
void add_mul_Mv_row_wise (double vout[], const double vin[])
 Add matrix-vector product vout += M*vin (vout(m),vin(n)) (general rowwise matrix)
 
void add_mul_Mtv_row_wise (double vout[], const double vin[])
 Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (general rowwise matrix)
 
void mul_Mv_column_wise (double vout[], const double vin[])
 Matrix-vector product vout=M*vin (vout(m),vin(n)) (general columnwise matrix)
 
void mul_Mtv_column_wise (double vout[], const double vin[])
 MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (general columnwise matrix)
 
void add_mul_Mv_column_wise (double vout[], const double vin[])
 Add matrix-vector product vout += M*vin (vout(m),vin(n)) (general columnwise matrix)
 
void add_mul_Mtv_column_wise (double vout[], const double vin[])
 Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (general columnwise matrix)
 
void mul_Mv_network (double vout[], const double vin[])
 Matrix-vector product vout=M*vin (vout(m),vin(n)) (network matrix, either oriented or nonoriented)
 
void mul_Mtv_network (double vout[], const double vin[])
 MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (network matrix, either oriented or nonoriented)
 
void add_mul_Mv_network (double vout[], const double vin[])
 Add matrix-vector product vout += M*vin (vout(m),vin(n)) (network matrix, either oriented or nonoriented)
 
void add_mul_Mtv_network (double vout[], const double vin[])
 Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (network matrix, either oriented or nonoriented)
 
void mul_Mv_identity (double vout[], const double vin[])
 Identity-vector product vout=M*vin (vout(m),vin(n))
 
void add_mul_Mv_identity (double vout[], const double vin[])
 Add Identity-vector product vout += M*vin (vout(m),vin(n))
 
void mul_Mv_diagonal (double vout[], const double vin[])
 Diagonal-vector product vout=D*vin (vout(m),vin(n)) D is diagonal, m=n.
 
void add_mul_Mv_diagonal (double vout[], const double vin[])
 Add Diagonal-vector product vout += D*vin (vout(m),vin(n)) D is diagonal, m=n.
 
void mul_Mv_idty_idty (double vout[], const double vin[])
 Matrix-vector product vout= [I I]*vin (vout(m),vin(n)) (M= [I I] matrix, n= 2*m)
 
void mul_Mtv_idty_idty (double vout[], const double vin[])
 MatrixTranspose-vector product vout= [I I]'*vin (vout(n),vin(m) (M= [I I] matrix, n= 2*m)
 
void add_mul_Mv_idty_idty (double vout[], const double vin[])
 Add matrix-vector product vout += [I I]*vin (vout(m),vin(n)) (M= [I I] matrix, n= 2*m)
 
void add_mul_Mtv_idty_idty (double vout[], const double vin[])
 Add matrixTranspose-vector product vout += [I I]'*vin (vout(n),vin(m) (M= [I I] matrix, n= 2*m)
 
void mul_Mv_diag_diag (double vout[], const double vin[])
 Matrix-vector product vout= [D1 D2]*vin (vout(m),vin(n)) (M= [D1 D2] matrix, n= 2*m)
 
void mul_Mtv_diag_diag (double vout[], const double vin[])
 MatrixTranspose-vector product vout= [D1 D2]'*vin (vout(n),vin(m) (M= [D1 D2] matrix, n= 2*m)
 
void add_mul_Mv_diag_diag (double vout[], const double vin[])
 Add matrix-vector product vout += [D1 D2]*vin (vout(m),vin(n)) (M= [D1 D2] matrix, n= 2*m)
 
void add_mul_Mtv_diag_diag (double vout[], const double vin[])
 Add matrixTranspose-vector product vout += [D1 D2]'*vin (vout(n),vin(m) (M= [D1 D2] matrix, n= 2*m)
 
void mul_Mv_gen_sym_uptr (double vout[], const double vin[])
 Matrix-vector product vout=M*vin (vout(m),vin(n)) (general symmetric upper triangular rowwise matrix)
 
void add_mul_Mv_gen_sym_uptr (double vout[], const double vin[])
 Add matrix-vector product vout += M*vin (vout(m),vin(n)) (general symmetric upper triangular rowwise matrix)
 
void order_packed_format (int begsize, int nz, int *&beg, int *&ind, double *&val)
 Order a matrix packed format (both column-wise and row-wise)
 

Detailed Description

Class for manipulating matrices, and interfacing SparseChol.

Constructor & Destructor Documentation

MatrixBlockIP::MatrixBlockIP ( int  blocks = 1)

Constructor.

Parameters
blocksnumber of blocks where this matrix will appear;
Note
if blocks unknown at construction time, do not pass this parameter to the constructor; it will be updated later, before symbolic factorization is performed, by the minimization algorithm
MatrixBlockIP::MatrixBlockIP ( MatrixBlockIP mbip)

Copy constructor, does not copy SparseChol.

Parameters
mbipMatrixBlockIP to copy

Member Function Documentation

void MatrixBlockIP::add_mul_Mtv ( double  vout[],
const double  vin[] 
)

Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (driver)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mtv_column_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (general columnwise matrix)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mtv_diag_diag ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrixTranspose-vector product vout += [D1 D2]'*vin (vout(n),vin(m) (M= [D1 D2] matrix, n= 2*m)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mtv_idty_idty ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrixTranspose-vector product vout += [I I]'*vin (vout(n),vin(m) (M= [I I] matrix, n= 2*m)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mtv_network ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (network matrix, either oriented or nonoriented)

Note
It assumes vin has dimension m+1; last equation, related to last node, is redundant and not used, but it is needed for an efficient computation.
Element vin[m] is used (but backed-up at beginning and restored at ending)
Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mtv_row_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrixTranspose-vector product vout += M(t)*vin (vout(n),vin(m) (general rowwise matrix)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mv ( double  vout[],
const double  vin[] 
)

Add matrix-vector product vout += M*vin (vout(m),vin(n)) (driver)

Parameters
voutResult of the product
vinVector to product

Here is the caller graph for this function:

void MatrixBlockIP::add_mul_Mv_column_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrix-vector product vout += M*vin (vout(m),vin(n)) (general columnwise matrix)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mv_diag_diag ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrix-vector product vout += [D1 D2]*vin (vout(m),vin(n)) (M= [D1 D2] matrix, n= 2*m)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mv_diagonal ( double  vout[],
const double  vin[] 
)
inlineprivate

Add Diagonal-vector product vout += D*vin (vout(m),vin(n)) D is diagonal, m=n.

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mv_gen_sym_uptr ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrix-vector product vout += M*vin (vout(m),vin(n)) (general symmetric upper triangular rowwise matrix)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::add_mul_Mv_identity ( double  vout[],
const double  vin[] 
)
inlineprivate

Add Identity-vector product vout += M*vin (vout(m),vin(n))

Identity matrix, m=n, just add input to output vector

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mv_idty_idty ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrix-vector product vout += [I I]*vin (vout(m),vin(n)) (M= [I I] matrix, n= 2*m)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mv_network ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrix-vector product vout += M*vin (vout(m),vin(n)) (network matrix, either oriented or nonoriented)

Note
It assumes vout has dimension m+1; last equation, related to last node, is redundant and not used, but it is needed for an efficient computation.
Element vout[m] is used (but backed-up at beginning and restored at ending)
Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_mul_Mv_row_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

Add matrix-vector product vout += M*vin (vout(m),vin(n)) (general rowwise matrix)

Parameters
voutResult of the product
vinVector to product
void MatrixBlockIP::add_new_column ( int  size,
int  irows[],
double  values[] 
)

Add a new column into the matrix.

Parameters
sizeNumber of non-zero elements of the new column
irowsRow position for each element, have to be ordered
valuesValue of non-zero elements
void MatrixBlockIP::add_new_columns ( int  num_columns,
int  size[],
int *  irows[],
double *  values[] 
)

Add new columns into the matrix.

Parameters
num_columnsNumber of columns to add
sizenumber of non-zero elements of the new column for each column
irowsRow position for each element and column, have to be ordered
valuesValue of non-zero elements for each column
void MatrixBlockIP::analyze_D ( int  numBlocks,
bool  sameL,
MatrixBlockIP  L[] 
)
inline

Create internal structure arrays to compute D = Theta_(numBlocks+1) + sum{i..numBlocks} L_i*Theta_i*L_i' It decides whether D is DIAGONAL or GEN_SYM_UPTR

Parameters
numBlocksNumber of diagonal blocks
sameLDefine if the same matrix is used for each L block
LLinking constraints blocks
void MatrixBlockIP::analyze_D_diagonal ( int  numBlocks,
bool  sameL,
MatrixBlockIP  L[] 
)

Create internal structure arrays to compute D = Theta_(numBlocks+1) + sum{1..numBlocks} L_i*Theta_i*L_i' when D is a diagonal

Parameters
numBlocksNumber of diagonal blocks
sameLDefine if the same matrix is used for each L block
LLinking constraints blocks
Note
valD (of analyze_D_gen_sym_uptr()) is reused, then it is allocated with ALLOC (see analyze_D_gen_sym_uptr for an explanation)
void MatrixBlockIP::analyze_D_gen_sym_uptr ( int  numBlocks,
bool  sameL,
MatrixBlockIP  L[] 
)

Create internal structure arrays to compute D = Theta_(numBlocks+1) + sum{i..numBlocks} L_i*Theta_i*L_i' when D is a general symmetric upper triangular matrix

Parameters
numBlocksNumber of diagonal blocks
sameLDefine if the same matrix is used for each L block
LLinking constraints blocks
Note
icolLLt, inivalD, blockD, valD and icolD are allocated with ALLOC for performance purposes, so they have to be freed with FREE
void MatrixBlockIP::change_columns_sign ( int  size,
int  columns[] 
)

Changes the sign of some columns.

Parameters
sizeNumber of columns to change the sign
columnsIndex to the columns to change the sign
void MatrixBlockIP::change_rows_sign ( int  size,
int  rows[] 
)

Change the sign of some rows.

Parameters
sizeNumber of rows to change the sign
rowsIndex to the rows to change the sign
void MatrixBlockIP::compute_D ( int  numBlocks,
bool  sameL,
MatrixBlockIP  L[],
bool  isActive[],
int  iniTheta[],
double  Theta[] 
)
inline

Compute D = Theta_(numBlocks+1) + sum{i..numBlocks} L_i*Theta_i*L_i'.

Parameters
numBlocksNumber of diagonal blocks
sameLDefine if the same matrix is used for each L block
isActiveSays what linking constrains are active
iniThetaIndice to the first element of the Theta matrix of each block
ThetanumBlocks diagonal matrices with dimension n_i x n_i
void MatrixBlockIP::compute_D_diagonal ( int  numBlocks,
bool  sameL,
MatrixBlockIP  L[],
bool  isActive[],
int  iniTheta[],
double  Theta[] 
)

Compute D = Theta_(numBlocks+1) + sum{i in 1..numBlocks} L_i*Theta_i*L_i' when D is diagonal.

Parameters
numBlocksNumber of diagonal blocks
sameLDefine if the same matrix is used for each L block
isActiveSays what linking constrains are active
iniThetaIndice to the first element of the Theta matrix of each block
ThetanumBlocks diagonal matrices with dimension n_i x n_i
void MatrixBlockIP::compute_D_gen_sym_uptr ( int  numBlocks,
bool  sameL,
MatrixBlockIP  L[],
bool  isActive[],
int  iniTheta[],
double  Theta[] 
)

Compute D = Theta_(numBlocks+1) + sum{i..numBlocks} L_i*Theta_i*L_i' when D is a symmetric general matrix.

Parameters
numBlocksNumber of diagonal blocks
sameLDefine if the same matrix is used for each L block
isActiveSays what linking constrains are active
iniThetaIndice to the first element of the Theta matrix of each block
ThetanumBlocks diagonal matrices with dimension n_i x n_i
void MatrixBlockIP::compute_full_matrix ( int  numBlocks,
bool  sameN,
MatrixBlockIP  N[],
bool  sameL,
MatrixBlockIP  L[],
int  numActiveLnk = 0,
int *  listActiveLnk = NULL 
)

Builds a packed rowwise format structured matrix based on diagonal blocks and linking constraints blocks.

Parameters
numBlocksNumber of diagonal blocks
sameNDefine if the same matrix is used for each N block
NDiagonal blocks
sameLDefine if the same matrix is used for each L block
LLinking constraints blocks
numActiveLnkNumber of active lnk
listActiveLnkList of active linking, j=listActiveLnk[i], i=1..numActiveLnk, j in {1,..,l_link}) If listActiveLnk is NULL all L rows are used

Here is the call graph for this function:

Here is the caller graph for this function:

void MatrixBlockIP::copy ( MatrixBlockIP mbip)

Copy, does not copy SparseChol.

Parameters
mbipMatrixBlockIP to copy
void MatrixBlockIP::create_diag_diag_matrix ( int  nrows,
double *&  d1,
double *&  d2 
)

Native method to create a diagonal-diagonal matrix D= [D1 D2].

where D1= diag(d1), D2= diag(d2) (dimension: nrows x (2nrows))

Parameters
nrowsNumber of rows
d1Values of the diagonal D1
d2Values of the diagonal D2
void MatrixBlockIP::create_diagonal_matrix ( int  dim,
double *&  d 
)

Native method to create a diagonal matrix D= diag(d) of dimension dim.

Parameters
dimDimension of the matrix
dValues of the diagonal elements

Here is the caller graph for this function:

void MatrixBlockIP::create_general_matrix_column_wise ( int  m,
int  n,
int  nz,
int *&  inicola,
int *&  irowa,
double *&  a 
)

Native method to create a general column-wise packed matrix.

Parameters
mNumber of rows
nNumber of columns
nzNumber of non-zero elements
inicolaIndex to the first element of each column
irowaRow position for each element
aValue of non-zero elements
Note
If icola is not ordered order_matrix function must be called when the arrays are filled

Here is the caller graph for this function:

void MatrixBlockIP::create_general_matrix_format_ija ( int  m,
int  n,
int  nz,
int *&  row_index,
int *&  col_index,
double *&  values 
)

Creates a general matrix in format ija.

Parameters
mNumber of rows
nNumber of columns
nzNumber of non-zero elements
row_indexRow position for each element
col_indexColumn position for each element
valuesValue of non-zero elements
void MatrixBlockIP::create_general_matrix_row_wise ( int  m,
int  n,
int  nz,
int *&  inirowa,
int *&  icola,
double *&  a 
)

Native method to create a general row-wise packed matrix.

Parameters
mNumber of rows
nNumber of columns
nzNumber of non-zero elements
inirowaIndex to the first element of each row
icolaColumn position for each element
aValue of non-zero elements
Note
If icola is not ordered order_matrix function must be called when the arrays are filled
void MatrixBlockIP::create_identity_matrix ( int  dim)

Native method to create a identity matrix of dimension dim.

Parameters
dimDimension of the matrix
void MatrixBlockIP::create_idty_idty_matrix ( int  nrows)

Native method to create a identity-identity matrix with two identities [I I].

Dimension: nrows x (2nrows)

Parameters
nrowsNumber of rows
void MatrixBlockIP::create_network_matrix ( int  num_arcs,
int  num_nodes,
int *&  src,
int *&  dst,
bool  oriented = true 
)

Native method to create a network matrix.

Parameters
num_arcsNumber of arcs
num_nodesNumber of nodes
srcSource of each arc. Dimension num_arcs
dstDestination of each arc. Dimension num_arcs
orientedFor oriented/nonoriented networks
void MatrixBlockIP::delete_rows ( int  size,
int  rows[] 
)

Delete some rows.

Parameters
sizeNumber of rows to be deleted
rowsIndex to the rows to be deleted, must be ordered
void MatrixBlockIP::exist_var_in_row ( int  column,
bool  appear[],
double  coefs[] 
)

Given a column determines if exists an element for each row of the matrix, in case of exist returns the value.

Parameters
columnThe column to examine
appearThe array that says if a element exists or not for each row. The user must allocate the space before call the function
coefsThe array that have the value of the element if exists, if not the content in that position is indeterminate. The user must allocate the space before call the function
int MatrixBlockIP::get_column ( int  column,
int *&  irows,
double *&  values,
bool  invert_sign = false 
)

Gets a column of the matrix.

Parameters
columnIndex to the column
irowsRow index for each ealement. Must be freed with delete[]
valuesValue for each non-zero element. Must be freed with delete[]
invert_signIf true change the sign of each element in the column
Returns
Number of non-zero elements in the column

Here is the caller graph for this function:

int MatrixBlockIP::get_ipfa ( int  i)
inline

Interface routine to SparseChol get_ipfa().

Returns
ipfa[i]
Note
The user must guarantee ipfa previously computed in call to symbolic_fact().
int MatrixBlockIP::get_maxfillin ( )
inline

Interface routine to SparseChol get_maxfillin().

Returns
maxfillin
Note
The user must guarantee maxfillin previously computed in call to symbolic_fact().
int MatrixBlockIP::get_maxlnz ( )
inline

Interface routine to SparseChol get_maxlnz().

Returns
maxlnz
Note
The user must guarantee maxlnz previously computed in call to symbolic_fact().
int MatrixBlockIP::get_njka ( )
inline

Interface routine to SparseChol get_njka().

Returns
njka
Note
The user must guarantee njka previously computed in call to symbolic_fact().
int MatrixBlockIP::get_num_semidef_matrix ( )
inline

Interface routine to SparseChol get_semidef_matrix().

Returns
num_semidef_matrix
Note
It will 0 if no numeric factorization made.
int MatrixBlockIP::get_num_zero_pivots ( )
inline

Interface routine to SparseChol get_num_zero_pivots().

Returns
num_zero_pivots
Note
It will 0 if no numeric factorization made.
int MatrixBlockIP::get_pfa ( int  i)
inline

Interface routine to SparseChol get_pfa().

Returns
pfa[i]
Note
The user must guarantee pfa previously computed in call to symbolic_fact().
void MatrixBlockIP::ija_to_rowwise ( )

Calculates and stores the matrix in packed rowwise format.

Note
ija format is loaded
void MatrixBlockIP::mul_Mtv ( double  vout[],
const double  vin[] 
)

MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (driver)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mtv_column_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (general columnwise matrix)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mtv_diag_diag ( double  vout[],
const double  vin[] 
)
inlineprivate

MatrixTranspose-vector product vout= [D1 D2]'*vin (vout(n),vin(m) (M= [D1 D2] matrix, n= 2*m)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mtv_idty_idty ( double  vout[],
const double  vin[] 
)
inlineprivate

MatrixTranspose-vector product vout= [I I]'*vin (vout(n),vin(m) (M= [I I] matrix, n= 2*m)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mtv_network ( double  vout[],
const double  vin[] 
)
inlineprivate

MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (network matrix, either oriented or nonoriented)

Note
It assumes vin has dimension m+1; last equation, related to last node, is redundant and not used, but it is needed for an efficient computation.
Element vin[m] is used (but backed-up at beginning and restored at ending)
Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mtv_row_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

MatrixTranspose-vector product vout=M(t)*vin (vout(n),vin(m) (general rowwise matrix)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv ( double  vout[],
const double  vin[] 
)

Matrix-vector product vout=M*vin (vout(m),vin(n)) (driver)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_column_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

Matrix-vector product vout=M*vin (vout(m),vin(n)) (general columnwise matrix)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_diag_diag ( double  vout[],
const double  vin[] 
)
inlineprivate

Matrix-vector product vout= [D1 D2]*vin (vout(m),vin(n)) (M= [D1 D2] matrix, n= 2*m)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_diagonal ( double  vout[],
const double  vin[] 
)
inlineprivate

Diagonal-vector product vout=D*vin (vout(m),vin(n)) D is diagonal, m=n.

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_gen_sym_uptr ( double  vout[],
const double  vin[] 
)
inlineprivate

Matrix-vector product vout=M*vin (vout(m),vin(n)) (general symmetric upper triangular rowwise matrix)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_identity ( double  vout[],
const double  vin[] 
)
inlineprivate

Identity-vector product vout=M*vin (vout(m),vin(n))

Identity matrix, m=n, just copy input to output vector

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_idty_idty ( double  vout[],
const double  vin[] 
)
inlineprivate

Matrix-vector product vout= [I I]*vin (vout(m),vin(n)) (M= [I I] matrix, n= 2*m)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_network ( double  vout[],
const double  vin[] 
)
inlineprivate

Matrix-vector product vout=M*vin (vout(m),vin(n)) (network matrix, either oriented or nonoriented)

Note
It assumes vout has dimension m+1; last equation, related to last node, is redundant and not used, but it is needed for an efficient computation.
Element vout[m] is used (but backed-up at beginning and restored at ending)
Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::mul_Mv_row_wise ( double  vout[],
const double  vin[] 
)
inlineprivate

Matrix-vector product vout=M*vin (vout(m),vin(n)) (general rowwise matrix)

Parameters
voutResult of the product. Must be allocated
vinVector to product
void MatrixBlockIP::native_to_general ( bool  delete_native_format = false)
inline

Calculates and stores the matrix in packed rowwise format.

Parameters
delete_native_formatIf true only packed rowwise format will be stored

Here is the caller graph for this function:

void MatrixBlockIP::order_packed_format ( int  begsize,
int  nz,
int *&  beg,
int *&  ind,
double *&  val 
)
private

Order a matrix packed format (both column-wise and row-wise)

Parameters
begsizeNumber of rows (if row-wise) or columns (if column-wise) + 1
nzNumber of nonzeros in the matrix.
begRows (if row-wise) or columns (if column-wise) starts of the matrix.
indColumns (if row-wise) or rows (if column-wise) indices of nonzeros entries.
valValues of the nonzero entries.
Note
Parameters beg, ind and val assumes to be of appropriate dimensions before the method is called, namely beg[begsize], ind[nz], val[nz]
void MatrixBlockIP::print_matrix ( bool  print_ija = true,
bool  start_one = true 
)

Print the matrix.

Parameters
print_ijaIf print_ija is true then write triple (i,j,a) whenever possible (if matrix is in ija or packed rowwise format); if print_ija=false then writes (inirow, j, a) whenever possible (if matrix is in packed rowwise format).
start_oneStart vectors at position 1 (true) or zero (false)
void MatrixBlockIP::print_matrix ( ofstream &  outfile,
bool  print_ija = true,
bool  start_one = true 
)

Write the matrix into a file.

Parameters
outfileOutput file stream where the matrix will be printed
void MatrixBlockIP::print_vector ( double *  v,
int  size,
string  name 
)

Print some positions and name of a vector.

Parameters
vVector to print
sizeNumber of positions to print
nameVector name

Member Data Documentation

SparseChol MatrixBlockIP::chol

Sparse Cholesky class.

For sparse Cholesky

double* MatrixBlockIP::d1

d1 for DIAGONAL and 1st submatrix of DIAG_DIAG

Diagonal format attributes

int* MatrixBlockIP::inicola

Index to the first element of each column.

General column-wise packed format attributes

int* MatrixBlockIP::inirowa

Index to the first element of each row.

General row-wise packed format attributes

bool MatrixBlockIP::made_analyze_D

boolean to check whether analyze_D already made, needed to compute_D()

bool MatrixBlockIP::made_symbfct_MMt

boolean to check whether symbolic factorization already made

int MatrixBlockIP::num_blocks

number of replications of this matrix; this is needed to store space for num_blocks numerical factorization of M*Theta[i]*M', for different Theta[i], i=0,...,num_blocks-1 This value is set by the minimization algorithm, which needs the factorizations.

int MatrixBlockIP::nz

Number of non-zero elements.

Common packed format attributes

int MatrixBlockIP::sizeL

For D Matrix

TYPE_MATRIX MatrixBlockIP::type

Type of matrix.

Common in all types matrix

TYPE_ORIENTATION MatrixBlockIP::type_orientation

Type of orientation, by default, oriented.

Network format attributes


The documentation for this class was generated from the following files: