|
Cyclops Tensor Framework
parallel arithmetic on multidimensional arrays
|
Matrix class which encapsulates a 2D tensor. More...
#include <matrix.h>


Public Member Functions | |
| Matrix () | |
| default constructor for a matrix More... | |
| Matrix (Matrix< dtype > const &A) | |
| copy constructor for a matrix More... | |
| Matrix (Tensor< dtype > const &A) | |
| casts a tensor to a matrix More... | |
| Matrix (int nrow, int ncol, World &wrld, CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, int profile=0) | |
| constructor for a matrix More... | |
| Matrix (int nrow, int ncol, World &wrld, char const *name, int profile=0, CTF_int::algstrct const &sr=Ring< dtype >()) | |
| constructor for a matrix More... | |
| Matrix (int nrow, int ncol, int atr=0, World &wrld=get_universe(), CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, int profile=0) | |
| constructor for a matrix More... | |
| Matrix (int nrow, int ncol, int atr, World &wrld, char const *name, int profile=0, CTF_int::algstrct const &sr=Ring< dtype >()) | |
| constructor for a matrix More... | |
| Matrix (int nrow, int ncol, char const *idx, Idx_Partition const &prl, Idx_Partition const &blk=Idx_Partition(), int atr=0, World &wrld=get_universe(), CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, int profile=0) | |
| constructor for a matrix with a given initial cyclic distribution More... | |
| void | write_mat (int mb, int nb, int pr, int pc, int rsrc, int csrc, int lda, dtype const *data) |
| writes a nonsymmetric matrix from a block-cyclic initial distribution this is `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix) More... | |
| Matrix (int nrow, int ncol, int mb, int nb, int pr, int pc, int rsrc, int csrc, int lda, dtype const *data, World &wrld=get_universe(), CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, int profile=0) | |
| constructor for a nonsymmetric matrix with a block-cyclic initial distribution this is `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix) More... | |
| Matrix (int const *desc, dtype const *data, World &wrld=get_universe(), CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, int profile=0) | |
| construct Matrix from ScaLAPACK array descriptor `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix) More... | |
| void | read_mat (int mb, int nb, int pr, int pc, int rsrc, int csrc, int lda, dtype *data) |
| reads a nonsymmetric matrix into a block-cyclic initial distribution this is `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix) More... | |
| void | get_desc (int &ictxt, int *&desc) |
| get a ScaLAPACK descriptor for this Matrix, will always be in pure cyclic layout More... | |
| void | read_mat (int const *desc, dtype *data) |
| read Matrix into ScaLAPACK array descriptor `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix) More... | |
| void | print_matrix () |
| void | qr (Matrix< dtype > &Q, Matrix< dtype > &R) |
| void | svd (Matrix< dtype > &U, Vector< dtype > &S, Matrix< dtype > &VT, int rank=0) |
Public Member Functions inherited from CTF::Tensor< dtype > | |
| Tensor () | |
| default constructor More... | |
| Tensor (int order, int const *len, int const *sym, World &wrld=get_universe(), char const *name=NULL, bool profile=0, CTF_int::algstrct const &sr=Ring< dtype >()) | |
| defines tensor filled with zeros on the default algstrct More... | |
| Tensor (int order, int const *len, int const *sym, World &wrld, CTF_int::algstrct const &sr, char const *name=NULL, bool profile=0) | |
| defines a tensor filled with zeros on a specified algstrct More... | |
| Tensor (int order, int const *len, World &wrld=get_universe(), CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, bool profile=0) | |
| defines a nonsymmetric tensor filled with zeros on a specified algstrct More... | |
| Tensor (int order, bool is_sparse, int const *len, int const *sym, World &wrld=get_universe(), CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, bool profile=0) | |
| defines a (sparse) tensor on a specified algstrct More... | |
| Tensor (int order, bool is_sparse, int const *len, World &wrld=get_universe(), CTF_int::algstrct const &sr=Ring< dtype >(), char const *name=NULL, bool profile=0) | |
| defines a nonsymmetric tensor filled with zeros on a specified algstrct More... | |
| Tensor (Tensor< dtype > const &A) | |
| copies a tensor, copying the data of A More... | |
| Tensor (tensor const &A) | |
| copies a tensor, copying the data of A (same as above) More... | |
| Tensor (bool copy, tensor const &A) | |
| copies a tensor (setting data to zero or copying A) More... | |
| Tensor (tensor &A, int const *new_sym) | |
| repacks the tensor other to a different symmetry (assumes existing data contains the symmetry and keeps only values with indices in increasing order) WARN: LIMITATION: new_sym must cause unidirectional structural changes, i.e. {NS,NS}->{SY,NS} OK, {SY,NS}->{NS,NS} OK, {NS,NS,SY,NS}->{SY,NS,NS,NS} NOT OK! More... | |
| Tensor (tensor const &A, World &wrld) | |
| creates a zeroed out copy (data not copied) of a tensor in a different world More... | |
| Tensor (int order, int const *len, int const *sym, World &wrld, char const *idx, Idx_Partition const &prl, Idx_Partition const &blk=Idx_Partition(), char const *name=NULL, bool profile=0, CTF_int::algstrct const &sr=Ring< dtype >()) | |
| defines tensor filled with zeros on the default algstrct on a user-specified distributed layout More... | |
| Tensor (int order, bool is_sparse, int const *len, int const *sym, World &wrld, char const *idx, Idx_Partition const &prl, Idx_Partition const &blk=Idx_Partition(), char const *name=NULL, bool profile=0, CTF_int::algstrct const &sr=Ring< dtype >()) | |
| defines tensor filled with zeros on the default algstrct on a user-specified distributed layout More... | |
| Typ_Idx_Tensor< dtype > | operator[] (char const *idx_map) |
| associated an index map with the tensor for future operation More... | |
| Typ_Idx_Tensor< dtype > | i (char const *idx_map) |
| void | read (int64_t npair, Pair< dtype > *pairs) |
| Gives the values associated with any set of indices. More... | |
| void | read (int64_t npair, int64_t const *global_idx, dtype *data) |
| Gives the values associated with any set of indices The sparse data is defined in coordinate format. The tensor index (i,j,k,l) of a tensor with edge lengths {m,n,p,q} is associated with the global index g via the formula g=i+j*m+k*m*n+l*m*n*p. The row index is first and the column index is second for matrices, which means they are column major. More... | |
| void | read (int64_t npair, dtype alpha, dtype beta, Pair< dtype > *pairs) |
| sparse read with accumulation: pairs[i].d = alpha*A[pairs[i].k]+beta*pairs[i].d More... | |
| void | read (int64_t npair, dtype alpha, dtype beta, int64_t const *global_idx, dtype *data) |
| sparse read with accumulation: data[i] = alpha*A[global_idx[i]]+beta*data[i] More... | |
| void | get_local_data (int64_t *npair, int64_t **global_idx, dtype **data, bool nonzeros_only=false, bool unpack_sym=false) const |
| Gives the global indices and values associated with the local data. More... | |
| void | read_local (int64_t *npair, int64_t **global_idx, dtype **data, bool unpack_sym=false) const |
| Using get_local_data(), which returns an array that must be freed with delete [], is more efficient, this version exists for backwards compatibility. gives the global indices and values associated with the local data WARNING-1: for sparse tensors this includes the zeros to maintain consistency with the behavior for dense tensors, use get_local_pairs to get only nonzeros. More... | |
| void | get_local_pairs (int64_t *npair, Pair< dtype > **pairs, bool nonzeros_only=false, bool unpack_sym=false) const |
| gives the global indices and values associated with the local data More... | |
| void | read_local (int64_t *npair, Pair< dtype > **pairs, bool unpack_sym=false) const |
| gives the global indices and values associated with the local data WARNING: for sparse tensors this includes the zeros to maintain consistency with the behavior for dense tensors, use get_lcoal_pairs to get only nonzeros More... | |
| void | get_all_data (int64_t *npair, dtype **data, bool unpack=false) |
| collects the entire tensor data on each process (not memory scalable) More... | |
| void | read_all (int64_t *npair, dtype **data, bool unpack=false) |
| collects the entire tensor data on each process (not memory scalable) More... | |
| int64_t | read_all (dtype *data, bool unpack=false) |
| collects the entire tensor data on each process (not memory scalable) More... | |
| void | write (int64_t npair, int64_t const *global_idx, dtype const *data) |
| writes in values associated with any set of indices The sparse data is defined in coordinate format. The tensor index (i,j,k,l) of a tensor with edge lengths {m,n,p,q} is associated with the global index g via the formula g=i+j*m+k*m*n+l*m*n*p. The row index is first and the column index is second for matrices, which means they are column major. More... | |
| void | write (int64_t npair, Pair< dtype > const *pairs) |
| writes in values associated with any set of indices More... | |
| void | write (int64_t npair, dtype alpha, dtype beta, int64_t const *global_idx, dtype const *data) |
| sparse add: A[global_idx[i]] = beta*A[global_idx[i]]+alpha*data[i] More... | |
| void | write (int64_t npair, dtype alpha, dtype beta, Pair< dtype > const *pairs) |
| sparse add: A[pairs[i].k] = alpha*A[pairs[i].k]+beta*pairs[i].d More... | |
| void | contract (dtype alpha, CTF_int::tensor &A, char const *idx_A, CTF_int::tensor &B, char const *idx_B, dtype beta, char const *idx_C) |
| contracts C[idx_C] = beta*C[idx_C] + alpha*A[idx_A]*B[idx_B] More... | |
| void | contract (dtype alpha, CTF_int::tensor &A, char const *idx_A, CTF_int::tensor &B, char const *idx_B, dtype beta, char const *idx_C, Bivar_Function< dtype > fseq) |
| contracts computes C[idx_C] = beta*C[idx_C] fseq(alpha*A[idx_A],B[idx_B]) More... | |
| void | sum (dtype alpha, CTF_int::tensor &A, char const *idx_A, dtype beta, char const *idx_B) |
| sums B[idx_B] = beta*B[idx_B] + alpha*A[idx_A] More... | |
| void | sum (dtype alpha, CTF_int::tensor &A, char const *idx_A, dtype beta, char const *idx_B, Univar_Function< dtype > fseq) |
| sums B[idx_B] = beta*B[idx_B] + fseq(alpha*A[idx_A]) More... | |
| void | scale (dtype alpha, char const *idx_A) |
| scales A[idx_A] = alpha*A[idx_A] More... | |
| void | scale (dtype alpha, char const *idx_A, Endomorphism< dtype > fseq) |
| scales A[idx_A] = fseq(alpha*A[idx_A]) More... | |
| dtype * | get_mapped_data (char const *idx, Idx_Partition const &prl, Idx_Partition const &blk=Idx_Partition(), bool unpack=true) |
| returns local data of tensor with parallel distribution prl and local blocking blk More... | |
| void | write (char const *idx, dtype const *data, Idx_Partition const &prl, Idx_Partition const &blk=Idx_Partition(), bool unpack=true) |
| writes data to tensor from parallel distribution prl and local blocking blk More... | |
| double | estimate_time (CTF_int::tensor &A, char const *idx_A, CTF_int::tensor &B, char const *idx_B, char const *idx_C) |
| estimate the time of a contraction C[idx_C] = A[idx_A]*B[idx_B] More... | |
| double | estimate_time (CTF_int::tensor &A, char const *idx_A, char const *idx_B) |
| estimate the time of a sum B[idx_B] = A[idx_A] More... | |
| Tensor< dtype > | slice (int const *offsets, int const *ends) const |
| cuts out a slice (block) of this tensor A[offsets,ends) result will always be fully nonsymmetric More... | |
| Tensor< dtype > | slice (int64_t corner_off, int64_t corner_end) const |
| cuts out a slice (block) of this tensor with corners specified by global index result will always be fully nonsymmetric More... | |
| Tensor< dtype > | slice (int const *offsets, int const *ends, World *oworld) const |
| cuts out a slice (block) of this tensor A[offsets,ends) result will always be fully nonsymmetric More... | |
| Tensor< dtype > | slice (int64_t corner_off, int64_t corner_end, World *oworld) const |
| cuts out a slice (block) of this tensor with corners specified by global index result will always be fully nonsymmetric More... | |
| void | slice (int const *offsets, int const *ends, dtype beta, CTF_int::tensor const &A, int const *offsets_A, int const *ends_A, dtype alpha) |
| adds to a slice (block) of this tensor = B B[offsets,ends)=beta*B[offsets,ends) + alpha*A[offsets_A,ends_A) More... | |
| void | slice (int64_t corner_off, int64_t corner_end, dtype beta, CTF_int::tensor const &A, int64_t corner_off_A, int64_t corner_end_A, dtype alpha) |
| adds to a slice (block) of this tensor = B B[offsets,ends)=beta*B[offsets,ends) + alpha*A[offsets_A,ends_A) More... | |
| void | permute (dtype beta, CTF_int::tensor &A, int *const *perms_A, dtype alpha) |
| Apply permutation to matrix, potentially extracting a slice B[i,j,...] = beta*B[...] + alpha*A[perms_A[0][i],perms_A[1][j],...]. More... | |
| void | permute (int *const *perms_B, dtype beta, CTF_int::tensor &A, dtype alpha) |
| Apply permutation to matrix, potentially extracting a slice B[perms_B[0][i],perms_B[0][j],...] = beta*B[...] + alpha*A[i,j,...]. More... | |
| void | sparsify () |
| reduce tensor to sparse format, storing only nonzero data, or data above a specified threshold. makes dense tensors sparse. cleans sparse tensors of any 'computed' zeros. More... | |
| void | sparsify (dtype threshold, bool take_abs=true) |
| reduce tensor to sparse format, storing only nonzero data, or data above a specified threshold. makes dense tensors sparse. cleans sparse tensors of any small values More... | |
| void | sparsify (std::function< bool(dtype)> filter) |
| sparsifies tensor keeping only values v such that filter(v) = true More... | |
| void | read_sparse_from_file (const char *fpath, bool with_vals=true) |
| read sparse tensor from file, entries of tensor must be stored one per line, as i_1 ... i_order v, to create entry T[i_1, ..., i_order] = v or as i_1 ... i_order, to create entry T[i_1, ..., i_order] = mulid More... | |
| void | write_sparse_to_file (const char *fpath, bool with_vals=true) |
| write sparse tensor to file, entries of tensor will be stored one per line, as i_1 ... i_order v, corresponding to entry T[i_1, ..., i_order] = v or as i_1 ... i_order if with_vals =false More... | |
| void | add_to_subworld (Tensor< dtype > *tsr, dtype alpha, dtype beta) |
| accumulates this tensor to a tensor object defined on a different world More... | |
| void | add_to_subworld (Tensor< dtype > *tsr) |
| accumulates this tensor to a tensor object defined on a different world More... | |
| void | add_from_subworld (Tensor< dtype > *tsr, dtype alpha, dtype beta) |
| accumulates this tensor from a tensor object defined on a different world More... | |
| void | add_from_subworld (Tensor< dtype > *tsr) |
| accumulates this tensor from a tensor object defined on a different world More... | |
| void | align (CTF_int::tensor const &A) |
| aligns data mapping with tensor A More... | |
| dtype | reduce (OP op) |
| performs a reduction on the tensor More... | |
| dtype | norm1 () |
| computes the entrywise 1-norm of the tensor More... | |
| dtype | norm2 () |
| computes the frobenius norm of the tensor (needs sqrt()!) More... | |
| dtype | norm_infty () |
| finds the max absolute value element of the tensor More... | |
| void | norm1 (double &nrm) |
| computes the entrywise 1-norm of the tensor More... | |
| void | norm2 (double &nrm) |
| computes the frobenius norm of the tensor More... | |
| void | norm_infty (double &nrm) |
| finds the max absolute value element of the tensor More... | |
| dtype * | get_raw_data (int64_t *size) const |
| gives the raw current local data with padding included More... | |
| const dtype * | raw_data (int64_t *size) const |
| gives a read-only copy of the raw current local data with padding included More... | |
| void | get_max_abs (int n, dtype *data) const |
| obtains a small number of the biggest elements of the tensor in sorted order (e.g. eigenvalues) More... | |
| void | fill_random (dtype rmin, dtype rmax) |
| fills local unique tensor elements to random values in the range [min,max] works only for dtype in {float,double,int,int64_t}, for others you can use Transform() uses mersenne twister seeded every time a global world is created (reseeded if all worlds destroyed) More... | |
| void | fill_sp_random (dtype rmin, dtype rmax, double frac_sp) |
| generate roughly frac_sp*dense_tensor_size nonzeros between rmin and rmax, works only for dtype in {float,double,int,int64_t}, for others you can use Transform() uses mersenne twister seeded every time a global world is created (reseeded if all worlds destroyed) More... | |
| void | profile_on () |
| turns on profiling for tensor More... | |
| void | profile_off () |
| turns off profiling for tensor More... | |
| void | set_name (char const *name) |
| sets tensor name More... | |
| Tensor< dtype > & | operator= (dtype val) |
| sets all values in the tensor to val More... | |
| Tensor< dtype > & | operator= (const Tensor< dtype > A) |
| sets the tensor More... | |
| Sparse_Tensor< dtype > | operator[] (std::vector< int64_t > indices) |
| gives handle to sparse index subset of tensors More... | |
| void | print (FILE *fp, dtype cutoff) const |
| prints tensor data to file using process 0 (modify print(...) overload in set.h if you would like a different print format) More... | |
| void | print (FILE *fp=stdout) const |
| void | prnt () const |
| void | compare (const Tensor< dtype > &A, FILE *fp=stdout, double cutoff=-1.0) |
| prints two sets of tensor data side-by-side to file using process 0 More... | |
| ~Tensor () | |
| frees CTF tensor More... | |
| template<> | |
| void | read_sparse_from_file (const char *fpath, bool with_vals) |
| template<> | |
| void | read_sparse_from_file (const char *fpath, bool with_vals) |
| template<> | |
| void | read_sparse_from_file (const char *fpath, bool with_vals) |
| template<> | |
| void | read_sparse_from_file (const char *fpath, bool with_vals) |
| template<> | |
| void | write_sparse_to_file (const char *fpath, bool with_vals) |
| template<> | |
| void | write_sparse_to_file (const char *fpath, bool with_vals) |
| template<> | |
| void | write_sparse_to_file (const char *fpath, bool with_vals) |
| template<> | |
| void | write_sparse_to_file (const char *fpath, bool with_vals) |
| template<> | |
| void | fill_random (double rmin, double rmax) |
| template<> | |
| void | fill_random (float rmin, float rmax) |
| template<> | |
| void | fill_random (int64_t rmin, int64_t rmax) |
| template<> | |
| void | fill_random (int rmin, int rmax) |
| template<> | |
| void | fill_sp_random (double rmin, double rmax, double frac_sp) |
| template<> | |
| void | fill_sp_random (float rmin, float rmax, double frac_sp) |
| template<> | |
| void | fill_sp_random (int rmin, int rmax, double frac_sp) |
| template<> | |
| void | fill_sp_random (int64_t rmin, int64_t rmax, double frac_sp) |
Public Member Functions inherited from CTF_int::tensor | |
| CTF::Idx_Tensor | operator[] (char const *idx_map) |
| associated an index map with the tensor for future operation More... | |
| tensor () | |
| ~tensor () | |
| class free self More... | |
| void | free_self () |
| destructor More... | |
| tensor (algstrct const *sr, int order, int const *edge_len, int const *sym, CTF::World *wrld, bool alloc_data=true, char const *name=NULL, bool profile=1, bool is_sparse=0) | |
| defines a tensor object with some mapping (if alloc_data) More... | |
| tensor (algstrct const *sr, int order, bool is_sparse, int const *edge_len, int const *sym, CTF::World *wrld, char const *idx, CTF::Idx_Partition const &prl, CTF::Idx_Partition const &blk, char const *name=NULL, bool profile=1) | |
| defines a tensor object with some mapping (if alloc_data) More... | |
| tensor (tensor const *other, bool copy=1, bool alloc_data=1) | |
| creates tensor copy, unfolds other if other is folded More... | |
| tensor (tensor *other, int const *new_sym) | |
| repacks the tensor other to a different symmetry (assumes existing data contains the symmetry and keeps only values with indices in increasing order) WARN: LIMITATION: new_sym must cause unidirectional structural changes, i.e. {NS,NS}->{SY,NS} OK, {SY,NS}->{NS,NS} OK, {NS,NS,SY,NS}->{SY,NS,NS,NS} NOT OK! More... | |
| int * | calc_phase () const |
| compute the cyclic phase of each tensor dimension More... | |
| int | calc_tot_phase () const |
| calculate the total number of blocks of the tensor More... | |
| int64_t | calc_nvirt () const |
| calculate virtualization factor of tensor return virtualization factor More... | |
| int64_t | calc_npe () const |
| calculate the number of processes this tensor is distributed over return number of processes owning a block of the tensor More... | |
| void | set_padding () |
| sets padding and local size of a tensor given a mapping More... | |
| int | set_zero () |
| elects a mapping and sets tensor data to zero More... | |
| int | set (char const *val) |
| sets tensor data to val More... | |
| int | zero_out_padding () |
| sets padded portion of tensor to zero (this should be maintained internally) More... | |
| void | scale_diagonals (int const *sym_mask) |
| scales each element by 1/(number of entries equivalent to it after permutation of indices for which sym_mask is 1) More... | |
| void | addinv () |
| void | print_map (FILE *stream=stdout, bool allcall=1) const |
| displays mapping information More... | |
| void | set_name (char const *name) |
| set the tensor name More... | |
| char const * | get_name () const |
| get the tensor name More... | |
| void | profile_on () |
| turn on profiling More... | |
| void | profile_off () |
| turn off profiling More... | |
| void | get_raw_data (char **data, int64_t *size) const |
| get raw data pointer without copy WARNING: includes padding More... | |
| int | write (int64_t num_pair, char const *alpha, char const *beta, char *mapped_data, char const rw='w') |
| Add tensor data new=alpha*new+beta*old with <key, value> pairs where key is the global index for the value. More... | |
| void | write (int64_t num_pair, char const *alpha, char const *beta, int64_t const *inds, char const *data) |
| Add tensor data new=alpha*new+beta*old with <key, value> pairs where key is the global index for the value. More... | |
| void | read (int64_t num_pair, char const *alpha, char const *beta, int64_t const *inds, char *data) |
| read tensor data with <key, value> pairs where key is the global index for the value, which gets filled in with beta times the old values plus alpha times the values read from the tensor. More... | |
| int | read (int64_t num_pair, char const *alpha, char const *beta, char *mapped_data) |
| read tensor data with <key, value> pairs where key is the global index for the value, which gets filled in with beta times the old values plus alpha times the values read from the tensor. More... | |
| char * | read (char const *idx, CTF::Idx_Partition const &prl, CTF::Idx_Partition const &blk, bool unpack) |
| returns local data of tensor with parallel distribution prl and local blocking blk More... | |
| int | read (int64_t num_pair, char *mapped_data) |
| read tensor data with <key, value> pairs where key is the global index for the value, which gets filled in. More... | |
| int64_t | get_tot_size (bool packed) |
| get number of elements in whole tensor More... | |
| int | allread (int64_t *num_pair, char **all_data, bool unpack) |
| read entire tensor with each processor (in packed layout). WARNING: will use an 'unscalable' amount of memory. More... | |
| int | allread (int64_t *num_pair, char *all_data, bool unpack=true) |
| read entire tensor with each processor (in packed layout). WARNING: will use an 'unscalable' amount of memory. More... | |
| void | slice (int const *offsets_B, int const *ends_B, char const *beta, tensor *A, int const *offsets_A, int const *ends_A, char const *alpha) |
| accumulates out a slice (block) of this tensor = B B[offsets,ends)=beta*B[offsets,ends) + alpha*A[offsets_A,ends_A) More... | |
| int | permute (tensor *A, int *const *permutation_A, char const *alpha, int *const *permutation_B, char const *beta) |
| int | sparsify (char const *threshold=NULL, bool take_abs=true) |
| reduce tensor to sparse format, storing only nonzero data, or data above a specified threshold. makes dense tensors sparse. cleans sparse tensors of any 'computed' zeros. More... | |
| int | sparsify (std::function< bool(char const *)> f) |
| sparsifies tensor keeping only values v such that filter(v) = true More... | |
| int | read_local (int64_t *num_pair, char **mapped_data, bool unpack_sym=false) const |
| read tensor data pairs local to processor including those with zero values WARNING: for sparse tensors this includes the zeros to maintain consistency with the behavior for dense tensors, use read_local_nnz to get only nonzeros More... | |
| int | read_local_nnz (int64_t *num_pair, char **mapped_data, bool unpack_sym=false) const |
| read tensor data pairs local to processor that have nonzero values More... | |
| void | read_local (int64_t *num_pair, int64_t **inds, char **data, bool unpack_sym=false) const |
| read tensor data pairs local to processor including those with zero values WARNING: for sparse tensors this includes the zeros to maintain consistency with the behavior for dense tensors, use read_local_nnz to get only nonzeros More... | |
| void | read_local_nnz (int64_t *num_pair, int64_t **inds, char **data, bool unpack_sym=false) const |
| read tensor data pairs local to processor that have nonzero values More... | |
| int | align (tensor const *B) |
| align mapping of thisa tensor to that of B More... | |
| int | reduce_sum (char *result) |
| Performs an elementwise summation reduction on a tensor. More... | |
| int | reduce_sum (char *result, algstrct const *sr_other) |
| Performs an elementwise summation reduction on a tensor with summation defined by sr_other. More... | |
| int | reduce_sumabs (char *result) |
| Performs an elementwise absolute value summation reduction on a tensor. More... | |
| int | reduce_sumabs (char *result, algstrct const *sr_other) |
| Performs an elementwise absolute value summation reduction on a tensor. More... | |
| int | reduce_sumsq (char *result) |
| computes the sum of squares of the elements More... | |
| int | get_max_abs (int n, char *data) const |
| obtains the largest n elements (in absolute value) of the tensor More... | |
| void | print (FILE *fp=stdout, char const *cutoff=NULL) const |
| prints tensor data to file using process 0 More... | |
| void | prnt () const |
| void | compare (const tensor *A, FILE *fp, char const *cutoff) |
| prints two sets of tensor data side-by-side to file using process 0 More... | |
| void | orient_subworld (CTF::World *greater_world, int &bw_mirror_rank, int &fw_mirror_rank, distribution *&odst, char **sub_buffer_) |
| maps data from this world (subcomm) to the correct order of processors with respect to a parent (greater_world) comm More... | |
| void | add_to_subworld (tensor *tsr_sub, char const *alpha, char const *beta) |
| accumulates this tensor to a tensor object defined on a different world More... | |
| void | add_from_subworld (tensor *tsr_sub, char const *alpha, char const *beta) |
| accumulates this tensor from a tensor object defined on a different world More... | |
| void | unfold (bool was_mod=0) |
| undo the folding of a local tensor block unsets is_folded and deletes rec_tsr More... | |
| void | remove_fold () |
| removes folding without doing transpose unsets is_folded and deletes rec_tsr More... | |
| double | est_time_unfold () |
| estimate cost of potential transpose involved in undoing the folding of a local tensor block More... | |
| void | fold (int nfold, int const *fold_idx, int const *idx_map, int *all_fdim, int **all_flen) |
| fold a tensor by putting the symmetry-preserved portion in the leading dimensions of the tensor sets is_folded and creates rec_tsr with aliased data More... | |
| void | pull_alias (tensor const *other) |
| pulls data from an tensor with an aliased buffer More... | |
| void | clear_mapping () |
| zeros out mapping More... | |
| int | redistribute (distribution const &old_dist, int const *old_offsets=NULL, int *const *old_permutation=NULL, int const *new_offsets=NULL, int *const *new_permutation=NULL) |
| permutes the data of a tensor to its new layout More... | |
| double | est_redist_time (distribution const &old_dist, double nnz_frac) |
| int64_t | get_redist_mem (distribution const &old_dist, double nnz_frac) |
| int | map_tensor_rem (int num_phys_dims, CommData *phys_comm, int fill=0) |
| map the remainder of a tensor More... | |
| int | extract_diag (int const *idx_map, int rw, tensor *&new_tsr, int **idx_map_new) |
| extracts the diagonal of a tensor if the index map specifies to do so More... | |
| void | set_sym (int const *sym) |
| sets symmetry, WARNING: for internal use only !!!! More... | |
| void | set_new_nnz_glb (int64_t const *nnz_blk) |
| sets the number of nonzeros both locally (nnz_loc) and overall globally (nnz_tot) More... | |
| void | spmatricize (int m, int n, int nrow_idx, bool csr) |
| transposes local data in preparation for summation or contraction, transforms to COO or CSR format for sparse More... | |
| void | despmatricize (int nrow_idx, bool csr) |
| transposes back local data from sparse matrix format to key-value pair format More... | |
| void | leave_home_with_buffer () |
| degister home buffer More... | |
| void | register_size (int64_t size) |
| register buffer allocation for this tensor More... | |
| void | deregister_size () |
| deregister buffer allocation for this tensor More... | |
| void | write_dense_to_file (MPI_File &file, int64_t offset=0) |
| write all tensor data to binary file in element order, unpacking from sparse or symmetric formats More... | |
| void | write_dense_to_file (char const *filename) |
| write all tensor data to binary file in element order, unpacking from sparse or symmetric formats More... | |
| void | read_dense_from_file (MPI_File &file, int64_t offset=0) |
| read all tensor data from binary file in element order, which should be stored as nonsymmetric and dense as done in write_dense_to_file() More... | |
| void | read_dense_from_file (char const *filename) |
| read all tensor data from binary file in element order, which should be stored as nonsymmetric and dense as done in write_dense_to_file() More... | |
| template<typename dtype_A , typename dtype_B > | |
| void | conv_type (tensor *B) |
| convert this tensor from dtype_A to dtype_B and store the result in B (primarily needed for python interface) More... | |
| template<typename dtype_A , typename dtype_B > | |
| void | exp_helper (tensor *A) |
| template<typename dtype > | |
| void | compare_elementwise (tensor *A, tensor *B) |
| do an elementwise comparison(==) of two tensors with elements of type dtype (primarily needed for python interface), store result in this tensor (has to be boolean tensor) More... | |
| template<typename dtype > | |
| void | not_equals (tensor *A, tensor *B) |
| do an elementwise comparison(!=) of two tensors with elements of type dtype (primarily needed for python interface), store result in this tensor (has to be boolean tensor) More... | |
| template<typename dtype > | |
| void | smaller_than (tensor *A, tensor *B) |
| do an elementwise comparison(<) of two tensors with elements of type dtype (primarily needed for python interface), store result in this tensor (has to be boolean tensor) More... | |
| template<typename dtype > | |
| void | smaller_equal_than (tensor *A, tensor *B) |
| do an elementwise comparison(<=) of two tensors with elements of type dtype (primarily needed for python interface), store result in this tensor (has to be boolean tensor) More... | |
| template<typename dtype > | |
| void | larger_than (tensor *A, tensor *B) |
| do an elementwise comparison(>) of two tensors with elements of type dtype (primarily needed for python interface), store result in this tensor (has to be boolean tensor) More... | |
| template<typename dtype > | |
| void | larger_equal_than (tensor *A, tensor *B) |
| do an elementwise comparison(>=) of two tensors with elements of type dtype (primarily needed for python interface), store result in this tensor (has to be boolean tensor) More... | |
| template<typename dtype > | |
| void | true_divide (tensor *A) |
| tensor * | self_reduce (int const *idx_A, int **new_idx_A, int order_B, int const *idx_B, int **new_idx_B, int order_C=0, int const *idx_C=NULL, int **new_idx_C=NULL) |
| performs a partial reduction on the tensor (used in summation and contraction) More... | |
Data Fields | |
| int | nrow |
| int | ncol |
| int | symm |
Data Fields inherited from CTF_int::tensor | |
| CTF::World * | wrld |
| distributed processor context on which tensor is defined More... | |
| algstrct * | sr |
| algstrct on which tensor elements and operations are defined More... | |
| int * | sym |
| symmetries among tensor dimensions More... | |
| int | order |
| number of tensor dimensions More... | |
| int * | lens |
| unpadded tensor edge lengths More... | |
| int * | pad_edge_len |
| padded tensor edge lengths More... | |
| int * | padding |
| padding along each edge length (less than distribution phase) More... | |
| char * | name |
| name given to tensor More... | |
| int | is_scp_padded |
| whether tensor data has additional padding More... | |
| int * | scp_padding |
| additional padding, may be greater than ScaLAPACK phase More... | |
| int * | sym_table |
| order-by-order table of dimensional symmetry relations More... | |
| bool | is_mapped |
| whether a mapping has been selected More... | |
| topology * | topo |
| topology to which the tensor is mapped More... | |
| mapping * | edge_map |
| mappings of each tensor dimension onto topology dimensions More... | |
| int64_t | size |
| current size of local tensor data chunk (mapping-dependent) More... | |
| int64_t | registered_alloc_size |
| size CTF keeps track of for memory usage More... | |
| bool | is_folded |
| whether the data is folded/transposed into a (lower-order) tensor More... | |
| int * | inner_ordering |
| ordering of the dimensions according to which the tensori s folded More... | |
| tensor * | rec_tsr |
| representation of folded tensor (shares data pointer) More... | |
| bool | is_cyclic |
| whether the tensor data is cyclically distributed (blocked if false) More... | |
| bool | is_data_aliased |
| whether the tensor data is an alias of another tensor object's data More... | |
| tensor * | slay |
| tensor object associated with tensor object whose data pointer needs to be preserved, needed for ScaLAPACK wrapper FIXME: home buffer should presumably take care of this... More... | |
| bool | has_zero_edge_len |
| if true tensor has a zero edge length, so is zero, which short-cuts stuff More... | |
| char * | data |
| tensor data, either the data or the key-value pairs should exist at any given time More... | |
| bool | has_home |
| whether the tensor has a home mapping/buffer More... | |
| char * | home_buffer |
| buffer associated with home mapping of tensor, to which it is returned More... | |
| int64_t | home_size |
| size of home buffer More... | |
| bool | is_home |
| whether the latest tensor data is in the home buffer More... | |
| bool | left_home_transp |
| whether the tensor left home to transpose More... | |
| bool | profile |
| whether profiling should be done for contractions/sums involving this tensor More... | |
| bool | is_sparse |
| whether only the non-zero elements of the tensor are stored More... | |
| bool | is_csr |
| whether CSR or COO if folded More... | |
| int | nrow_idx |
| how many modes are folded into matricized row More... | |
| int64_t | nnz_loc |
| number of local nonzero elements More... | |
| int64_t | nnz_tot |
| maximum number of local nonzero elements over all procs More... | |
| int64_t * | nnz_blk |
| nonzero elements in each block owned locally More... | |
Additional Inherited Members | |
Protected Member Functions inherited from CTF_int::tensor | |
| void | init (algstrct const *sr, int order, int const *edge_len, int const *sym, CTF::World *wrld, bool alloc_data, char const *name, bool profile, bool is_sparse) |
| initializes tensor data More... | |
| PairIterator | read_all_pairs (int64_t *num_pair, bool unpack) |
| read all pairs with each processor (packed) More... | |
| void | copy_tensor_data (tensor const *other) |
| copies all tensor data from other More... | |
| void | set_distribution (char const *idx, CTF::Idx_Partition const &prl, CTF::Idx_Partition const &blk) |
| set edge mappings as specified More... | |
Matrix class which encapsulates a 2D tensor.
| [in] | dtype | specifies tensor element type |
| CTF::Matrix< dtype >::Matrix | ( | ) |
default constructor for a matrix
Definition at line 29 of file matrix.cxx.
References CTF::Matrix< dtype >::ncol, and CTF::Matrix< dtype >::nrow.
| CTF::Matrix< dtype >::Matrix | ( | Matrix< dtype > const & | A | ) |
copy constructor for a matrix
| [in] | A | matrix to copy along with its data |
Definition at line 35 of file matrix.cxx.
References CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, and CTF::Matrix< dtype >::symm.
| CTF::Matrix< dtype >::Matrix | ( | Tensor< dtype > const & | A | ) |
casts a tensor to a matrix
| [in] | A | tensor object of order 2 |
Definition at line 43 of file matrix.cxx.
References AS, IASSERT, CTF_int::tensor::lens, CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, NS, CTF_int::tensor::order, SY, CTF_int::tensor::sym, and CTF::Matrix< dtype >::symm.
| CTF::Matrix< dtype >::Matrix | ( | int | nrow, |
| int | ncol, | ||
| World & | wrld, | ||
| CTF_int::algstrct const & | sr = Ring<dtype>(), |
||
| char const * | name = NULL, |
||
| int | profile = 0 |
||
| ) |
constructor for a matrix
| [in] | nrow | number of matrix rows |
| [in] | ncol | number of matrix columns |
| [in] | wrld | CTF world where the tensor will live |
| [in] | sr | defines the tensor arithmetic for this tensor |
| [in] | name | an optionary name for the tensor |
| [in] | profile | set to 1 to profile contractions involving this tensor |
Definition at line 66 of file matrix.cxx.
References CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, NS, and CTF::Matrix< dtype >::symm.
| CTF::Matrix< dtype >::Matrix | ( | int | nrow, |
| int | ncol, | ||
| World & | wrld, | ||
| char const * | name, | ||
| int | profile = 0, |
||
| CTF_int::algstrct const & | sr = Ring<dtype>() |
||
| ) |
constructor for a matrix
| [in] | nrow | number of matrix rows |
| [in] | ncol | number of matrix columns |
| [in] | wrld | CTF world where the tensor will live |
| [in] | name | an optionary name for the tensor |
| [in] | profile | set to 1 to profile contractions involving this tensor |
| [in] | sr | defines the tensor arithmetic for this tensor |
Definition at line 95 of file matrix.cxx.
References CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, and CTF::Matrix< dtype >::symm.
| CTF::Matrix< dtype >::Matrix | ( | int | nrow, |
| int | ncol, | ||
| int | atr = 0, |
||
| World & | wrld = get_universe(), |
||
| CTF_int::algstrct const & | sr = Ring<dtype>(), |
||
| char const * | name = NULL, |
||
| int | profile = 0 |
||
| ) |
constructor for a matrix
| [in] | nrow | number of matrix rows |
| [in] | ncol | number of matrix columns |
| [in] | atr | quantifier for sparsity and symmetry of matrix |
| [in] | wrld | CTF world where the tensor will live |
| [in] | sr | defines the tensor arithmetic for this tensor |
| [in] | name | an optionary name for the tensor |
| [in] | profile | set to 1 to profile contractions involving this tensor |
Definition at line 80 of file matrix.cxx.
References CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, and CTF::Matrix< dtype >::symm.
| CTF::Matrix< dtype >::Matrix | ( | int | nrow, |
| int | ncol, | ||
| int | atr, | ||
| World & | wrld, | ||
| char const * | name, | ||
| int | profile = 0, |
||
| CTF_int::algstrct const & | sr = Ring<dtype>() |
||
| ) |
constructor for a matrix
| [in] | nrow | number of matrix rows |
| [in] | ncol | number of matrix columns |
| [in] | atr | quantifier for sparsity and symmetry of matrix |
| [in] | wrld | CTF world where the tensor will live |
| [in] | name | an optionary name for the tensor |
| [in] | profile | set to 1 to profile contractions involving this tensor |
| [in] | sr | defines the tensor arithmetic for this tensor |
Definition at line 110 of file matrix.cxx.
References CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, and CTF::Matrix< dtype >::symm.
| CTF::Matrix< dtype >::Matrix | ( | int | nrow, |
| int | ncol, | ||
| char const * | idx, | ||
| Idx_Partition const & | prl, | ||
| Idx_Partition const & | blk = Idx_Partition(), |
||
| int | atr = 0, |
||
| World & | wrld = get_universe(), |
||
| CTF_int::algstrct const & | sr = Ring<dtype>(), |
||
| char const * | name = NULL, |
||
| int | profile = 0 |
||
| ) |
constructor for a matrix with a given initial cyclic distribution
| [in] | nrow | number of matrix rows |
| [in] | ncol | number of matrix columns |
| [in] | idx | assignment of characters to each dim |
| [in] | prl | mesh processor topology with character labels |
| [in] | blk | local blocking with processor labels |
| [in] | atr | quantifier for sparsity and symmetry of matrix |
| [in] | wrld | CTF world where the tensor will live |
| [in] | sr | defines the tensor arithmetic for this tensor |
| [in] | name | an optionary name for the tensor |
| [in] | profile | set to 1 to profile contractions involving this tensor |
Definition at line 126 of file matrix.cxx.
References CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, and CTF::Matrix< dtype >::symm.
| CTF::Matrix< dtype >::Matrix | ( | int | nrow, |
| int | ncol, | ||
| int | mb, | ||
| int | nb, | ||
| int | pr, | ||
| int | pc, | ||
| int | rsrc, | ||
| int | csrc, | ||
| int | lda, | ||
| dtype const * | data, | ||
| World & | wrld = get_universe(), |
||
| CTF_int::algstrct const & | sr = Ring<dtype>(), |
||
| char const * | name = NULL, |
||
| int | profile = 0 |
||
| ) |
constructor for a nonsymmetric matrix with a block-cyclic initial distribution this is `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix)
| [in] | nrow | number of matrix rows |
| [in] | ncol | number of matrix columns |
| [in] | mb | row block dimension |
| [in] | nb | col block dimension |
| [in] | pr | number of rows in processor grid |
| [in] | pc | number of cols in processor grid |
| [in] | rsrc | processor row holding first block row (0-based unlike ScaLAPACK) |
| [in] | csrc | processor col holding first block row (0-based unlike ScaLAPACK) |
| [in] | lda | leading dimension (length of buffer corresponding to row) |
| [in] | data | locally stored values |
| [in] | wrld | CTF world where the tensor will live, must contain pr*pc processors |
| [in] | sr | defines the tensor arithmetic for this tensor |
| [in] | name | an optionary name for the tensor |
| [in] | profile | set to 1 to profile contractions involving this tensor |
Definition at line 344 of file matrix.cxx.
References CTF_SCALAPACK::cblacs_gridinfo(), CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, NS, CTF::Matrix< dtype >::symm, and CTF::Matrix< dtype >::write_mat().
| CTF::Matrix< dtype >::Matrix | ( | int const * | desc, |
| dtype const * | data, | ||
| World & | wrld = get_universe(), |
||
| CTF_int::algstrct const & | sr = Ring<dtype>(), |
||
| char const * | name = NULL, |
||
| int | profile = 0 |
||
| ) |
construct Matrix from ScaLAPACK array descriptor `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix)
| [in] | desc | ScaLAPACK descriptor array: see ScaLAPACK docs for "Array Descriptor for In-core Dense Matrices" |
| [in] | data | locally stored values |
| [in] | wrld | CTF world where the tensor will live, must contain pr*pc processors |
| [in] | sr | defines the tensor arithmetic for this tensor |
| [in] | name | an optionary name for the tensor |
| [in] | profile | set to 1 to profile contractions involving this tensor |
Definition at line 377 of file matrix.cxx.
References CTF_SCALAPACK::cblacs_gridinfo(), IASSERT, CTF::Matrix< dtype >::ncol, CTF::World::np, CTF::Matrix< dtype >::nrow, NS, CTF::World::rank, CTF::Matrix< dtype >::symm, and CTF::Matrix< dtype >::write_mat().
| void CTF::Matrix< dtype >::get_desc | ( | int & | ictxt, |
| int *& | desc | ||
| ) |
get a ScaLAPACK descriptor for this Matrix, will always be in pure cyclic layout
| [out] | ictxt | index of newly created context |
| [out] | desc | array of integers of size 9, which will be filled with attributes see ScaLAPACK docs for "Array Descriptor for In-core Dense Matrices" |
Definition at line 296 of file matrix.cxx.
References CTF_int::mapping::calc_phase(), CTF_SCALAPACK::cblacs_get(), CTF_SCALAPACK::cblacs_gridinit(), CTF::World::comm, CTF_int::grid_wrapper::ctxt, CTF_int::tensor::edge_map, IASSERT, CTF::Matrix< dtype >::ncol, CTF::World::np, CTF::Matrix< dtype >::nrow, CTF_int::tensor::pad_edge_len, CTF_int::grid_wrapper::pc, CTF_int::grid_wrapper::pr, CTF_int::scalapack_grids, and CTF_int::tensor::wrld.
Referenced by CTF::Matrix< dtype >::qr(), and CTF::Matrix< dtype >::svd().
| void CTF::Matrix< dtype >::print_matrix | ( | ) |
Definition at line 144 of file matrix.cxx.
References CTF_int::tensor::data, ctf.core::dtype, CTF::Tensor< dtype >::i(), CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, CTF_int::algstrct::print(), CTF::World::rank, CTF::Tensor< dtype >::read_all(), CTF_int::tensor::sr, and CTF_int::tensor::wrld.
Referenced by main().
| void CTF::Matrix< dtype >::qr | ( | Matrix< dtype > & | Q, |
| Matrix< dtype > & | R | ||
| ) |
Definition at line 425 of file matrix.cxx.
References ctf.core::dtype, CTF::Matrix< dtype >::get_desc(), CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, NS, CTF::Matrix< dtype >::read_mat(), CTF_int::tensor::size, CTF::Tensor< dtype >::slice(), CTF_int::tensor::sr, SY, and CTF_int::tensor::wrld.
Referenced by CTF_int::matrix_qr(), CTF_int::matrix_qr_cmplx(), and qr().
| void CTF::Matrix< dtype >::read_mat | ( | int | mb, |
| int | nb, | ||
| int | pr, | ||
| int | pc, | ||
| int | rsrc, | ||
| int | csrc, | ||
| int | lda, | ||
| dtype * | data | ||
| ) |
reads a nonsymmetric matrix into a block-cyclic initial distribution this is `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix)
| [in] | mb | row block dimension |
| [in] | nb | col block dimension |
| [in] | pr | number of rows in processor grid |
| [in] | pc | number of cols in processor grid |
| [in] | rsrc | processor row holding first block row (0-based unlike ScaLAPACK) |
| [in] | csrc | processor col holding first block row (0-based unlike ScaLAPACK) |
| [in] | lda | leading dimension (length of buffer corresponding to row) |
| [out] | data | locally stored values |
Definition at line 247 of file matrix.cxx.
References CTF::Pair< dtype >::d, CTF_int::tensor::data, ctf.core::dtype, CTF_int::tensor::edge_map, CTF::get_my_kv_pair(), CTF::Tensor< dtype >::i(), CTF_int::tensor::is_sparse, CTF::Matrix< dtype >::ncol, ctf.core::np(), CTF::Matrix< dtype >::nrow, CTF::World::rank, CTF::Tensor< dtype >::read(), CTF::Matrix< dtype >::read_mat(), CTF_int::tensor::size, CTF_int::tensor::sr, and CTF_int::tensor::wrld.
Referenced by CTF::Matrix< dtype >::qr(), CTF::Matrix< dtype >::read_mat(), and CTF::Matrix< dtype >::svd().
| void CTF::Matrix< dtype >::read_mat | ( | int const * | desc, |
| dtype * | data | ||
| ) |
read Matrix into ScaLAPACK array descriptor `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix)
| [in] | desc | ScaLAPACK descriptor array: see ScaLAPACK docs for "Array Descriptor for In-core Dense Matrices" |
| [in] | data | locally stored values |
Definition at line 332 of file matrix.cxx.
References CTF_SCALAPACK::cblacs_gridinfo(), IASSERT, CTF::World::rank, CTF::Matrix< dtype >::read_mat(), and CTF_int::tensor::wrld.
| void CTF::Matrix< dtype >::svd | ( | Matrix< dtype > & | U, |
| Vector< dtype > & | S, | ||
| Matrix< dtype > & | VT, | ||
| int | rank = 0 |
||
| ) |
Definition at line 481 of file matrix.cxx.
References CTF_int::alloc(), CTF_int::mapping::calc_phase(), CTF_int::cdealloc(), CTF_SCALAPACK::cdescinit(), ctf.core::dtype, CTF_int::tensor::edge_map, CTF::Matrix< dtype >::get_desc(), CTF::Tensor< dtype >::get_raw_data(), CTF::Tensor< dtype >::i(), CTF::Matrix< dtype >::ncol, CTF::Matrix< dtype >::nrow, CTF::World::rank, CTF::Matrix< dtype >::read_mat(), CTF_int::tensor::size, CTF::Tensor< dtype >::slice(), CTF_int::tensor::topo, and CTF_int::tensor::wrld.
Referenced by get_factor_matrices(), CTF::get_factor_matrices(), main(), CTF_int::matrix_svd(), CTF_int::matrix_svd_cmplx(), and svd().
| void CTF::Matrix< dtype >::write_mat | ( | int | mb, |
| int | nb, | ||
| int | pr, | ||
| int | pc, | ||
| int | rsrc, | ||
| int | csrc, | ||
| int | lda, | ||
| dtype const * | data | ||
| ) |
writes a nonsymmetric matrix from a block-cyclic initial distribution this is `cheap' if mb=nb=1, nrowpr=0, ncolpc=0, rsrc=0, csrc=0, but is done via sparse read/write otherwise assumes processor grid is row-major (otherwise transpose matrix)
| [in] | mb | row block dimension |
| [in] | nb | col block dimension |
| [in] | pr | number of rows in processor grid |
| [in] | pc | number of cols in processor grid |
| [in] | rsrc | processor row holding first block row (0-based unlike ScaLAPACK) |
| [in] | csrc | processor col holding first block row (0-based unlike ScaLAPACK) |
| [in] | lda | leading dimension (length of buffer corresponding to row) |
| [out] | data | locally stored values |
Definition at line 203 of file matrix.cxx.
References CTF::Pair< dtype >::d, CTF_int::tensor::data, ctf.core::dtype, CTF_int::tensor::edge_map, CTF::get_my_kv_pair(), CTF::Tensor< dtype >::i(), CTF::Matrix< dtype >::ncol, ctf.core::np(), CTF::Matrix< dtype >::nrow, CTF::World::rank, CTF_int::tensor::size, CTF::Tensor< dtype >::write(), and CTF_int::tensor::wrld.
Referenced by CTF::Matrix< dtype >::Matrix().
| int CTF::Matrix< dtype >::ncol |
Definition at line 20 of file matrix.h.
Referenced by ao_mo_transf_naive(), ao_mo_transf_slice(), CTF::Matrix< dtype >::get_desc(), CTF::Matrix< dtype >::Matrix(), CTF::Matrix< dtype >::print_matrix(), CTF::Matrix< dtype >::qr(), CTF::Matrix< dtype >::read_mat(), CTF::Matrix< dtype >::svd(), and CTF::Matrix< dtype >::write_mat().
| int CTF::Matrix< dtype >::nrow |
Definition at line 20 of file matrix.h.
Referenced by ao_mo_transf_naive(), ao_mo_transf_slice(), btwn_cnt_fast(), btwn_cnt_naive(), CTF::Matrix< dtype >::get_desc(), CTF::Matrix< dtype >::Matrix(), mis(), mis2(), CTF::Matrix< dtype >::print_matrix(), CTF::Matrix< dtype >::qr(), CTF::Matrix< dtype >::read_mat(), CTF::Matrix< dtype >::svd(), and CTF::Matrix< dtype >::write_mat().
| int CTF::Matrix< dtype >::symm |
Definition at line 20 of file matrix.h.
Referenced by CTF::Matrix< dtype >::Matrix(), mis(), and mis2().