#include <DelayedStream.h>
|
| DelayedStream (const shared_ptr< tatami::NumericMatrix > &ptr, const vector< string > &rowNames, const int &chunkSize) |
| ~DelayedStream () |
void | setRegions (const vector< string > ®ions) override |
int | n_samples () override |
int | n_rows () |
vector< string > | getSampleNames () override |
string | getStreamType () override |
GenomicRanges | getChromRanges () override |
bool | getNextChunk (DataChunk< arma::mat > &chunk, const bool &useFilter=true) override |
bool | getNextChunk (DataChunk< arma::sp_mat > &chunk, const bool &useFilter=true) override |
bool | getNextChunk (DataChunk< Eigen::MatrixXd > &chunk, const bool &useFilter=true) override |
bool | getNextChunk (DataChunk< Eigen::SparseMatrix< double > > &chunk, const bool &useFilter=true) override |
bool | getNextChunk (DataChunk< Rcpp::NumericMatrix > &chunk, const bool &useFilter=true) override |
bool | getNextChunk (DataChunk< vector< double > > &chunk, const bool &useFilter=true) override |
void | getNextChunk (DataChunk< Eigen::MatrixXd > &chunk, int start, int len) |
| GenomicDataStream () |
| GenomicDataStream (const Param ¶m) |
virtual | ~GenomicDataStream () |
double | getMinVariance () |
double | getMAF () const |
void | setMinVariance (const double &value) |
void | setChunkSize (const int &chunkSize) |
Param | getParam () const |
◆ DelayedStream()
gds::DelayedStream::DelayedStream |
( |
const shared_ptr< tatami::NumericMatrix > & | ptr, |
|
|
const vector< string > & | rowNames, |
|
|
const int & | chunkSize ) |
|
inline |
◆ ~DelayedStream()
gds::DelayedStream::~DelayedStream |
( |
| ) |
|
|
inline |
◆ getChromRanges()
◆ getNextChunk() [1/7]
bool gds::DelayedStream::getNextChunk |
( |
DataChunk< arma::mat > & | chunk, |
|
|
const bool & | useFilter = true ) |
|
inlineoverridevirtual |
◆ getNextChunk() [2/7]
bool gds::DelayedStream::getNextChunk |
( |
DataChunk< arma::sp_mat > & | chunk, |
|
|
const bool & | useFilter = true ) |
|
inlineoverridevirtual |
◆ getNextChunk() [3/7]
bool gds::DelayedStream::getNextChunk |
( |
DataChunk< Eigen::MatrixXd > & | chunk, |
|
|
const bool & | useFilter = true ) |
|
inlineoverridevirtual |
◆ getNextChunk() [4/7]
void gds::DelayedStream::getNextChunk |
( |
DataChunk< Eigen::MatrixXd > & | chunk, |
|
|
int | start, |
|
|
int | len ) |
|
inline |
Get chunks based on start and end
◆ getNextChunk() [5/7]
bool gds::DelayedStream::getNextChunk |
( |
DataChunk< Eigen::SparseMatrix< double > > & | chunk, |
|
|
const bool & | useFilter = true ) |
|
inlineoverridevirtual |
◆ getNextChunk() [6/7]
bool gds::DelayedStream::getNextChunk |
( |
DataChunk< Rcpp::NumericMatrix > & | chunk, |
|
|
const bool & | useFilter = true ) |
|
inlineoverridevirtual |
◆ getNextChunk() [7/7]
bool gds::DelayedStream::getNextChunk |
( |
DataChunk< vector< double > > & | chunk, |
|
|
const bool & | useFilter = true ) |
|
inlineoverridevirtual |
◆ getSampleNames()
vector< string > gds::DelayedStream::getSampleNames |
( |
| ) |
|
|
inlineoverridevirtual |
◆ getStreamType()
string gds::DelayedStream::getStreamType |
( |
| ) |
|
|
inlineoverridevirtual |
◆ n_rows()
int gds::DelayedStream::n_rows |
( |
| ) |
|
|
inline |
Get number of rows in data matrix
◆ n_samples()
int gds::DelayedStream::n_samples |
( |
| ) |
|
|
inlineoverridevirtual |
◆ setRegions()
void gds::DelayedStream::setRegions |
( |
const vector< string > & | regions | ) |
|
|
inlineoverridevirtual |
The documentation for this class was generated from the following file: