queso-0.51.1
Private Member Functions | Private Attributes | List of all members
QUESO::ArrayOfSequences< V, M > Class Template Reference

Class for handling array samples (arrays of scalar sequences). More...

#include <ArrayOfSequences.h>

Inheritance diagram for QUESO::ArrayOfSequences< V, M >:
Inheritance graph
[legend]
Collaboration diagram for QUESO::ArrayOfSequences< V, M >:
Collaboration graph
[legend]

Public Member Functions

Constructor/Destructor methods
 ArrayOfSequences (const VectorSpace< V, M > &vectorSpace, unsigned int subSequenceSize, const std::string &name)
 Default constructor. More...
 
 ~ArrayOfSequences ()
 Destructor. More...
 
Sequence methods
unsigned int subSequenceSize () const
 Size of the sub-sequence of arrays. More...
 
void resizeSequence (unsigned int newSubSequenceSize)
 Resizes the sequence. More...
 
void resetValues (unsigned int initialPos, unsigned int numPos)
 Resets a total of numPos values of the sequence starting at position initialPos. More...
 
void erasePositions (unsigned int initialPos, unsigned int numPos)
 Erases numPos elements of the sequence starting at position initialPos. More...
 
void getPositionValues (unsigned int posId, V &vec) const
 Gets the values of the sequence at position posId and stores them at vec. More...
 
void setPositionValues (unsigned int posId, const V &vec)
 Set the values of vec at position posId of the sequence. More...
 
void setGaussian (const V &meanVec, const V &stdDevVec)
 Sets the values of the sequence as a Gaussian distribution of mean given by meanVec and standard deviation by stdDevVec. More...
 
void setUniform (const V &aVec, const V &bVec)
 Sets the values of the sequence as a uniform distribution between the values given by vectors aVec and bVec. More...
 
void mean (unsigned int initialPos, unsigned int numPos, V &meanVec) const
 Finds the mean value of the sub-sequence, considering numPos positions starting at position initialPos. More...
 
void unifiedMean (unsigned int initialPos, unsigned int numPos, V &unifiedMeanVec) const
 Finds the mean value of the unified sequence of numPos positions starting at position initialPos. See template specialization. More...
 
void sampleVariance (unsigned int initialPos, unsigned int numPos, const V &meanVec, V &samVec) const
 Finds the sample variance of the sub-sequence, considering numPos positions starting at position initialPos and of mean meanVec. More...
 
void unifiedSampleVariance (unsigned int initialPos, unsigned int numPos, const V &meanVec, V &samVec) const
 Finds the sample variance of the unified sequence, considering numPos positions starting at position initialPos and of mean meanVec. More...
 
void populationVariance (unsigned int initialPos, unsigned int numPos, const V &meanVec, V &popVec) const
 Finds the population variance of the sub-sequence, considering numPos positions starting at position initialPos and of mean meanVec. More...
 
void autoCovariance (unsigned int initialPos, unsigned int numPos, const V &meanVec, unsigned int lag, V &covVec) const
 Calculates the autocovariance. More...
 
void autoCorrViaDef (unsigned int initialPos, unsigned int numPos, unsigned int lag, V &corrVec) const
 Calculates autocorrelation via definition. More...
 
void autoCorrViaFft (unsigned int initialPos, unsigned int numPos, const std::vector< unsigned int > &lags, std::vector< V * > &corrVecs) const
 
void autoCorrViaFft (unsigned int initialPos, unsigned int numPos, unsigned int numSum, V &autoCorrsSumVec) const
 
void minMax (unsigned int initialPos, V &minVec, V &maxVec) const
 Given an initial position initialPos, finds the minimum and the maximum values of the sequence. More...
 
void histogram (unsigned int initialPos, const V &minVec, const V &maxVec, std::vector< V * > &centersForAllBins, std::vector< V * > &quanttsForAllBins) const
 Calculates the histogram of the sequence. More...
 
void interQuantileRange (unsigned int initialPos, V &iqrs) const
 Returns the interquartile range of the values in the sequence. More...
 
void scalesForKDE (unsigned int initialPos, const V &iqrs, unsigned int kdeDimension, V &scales) const
 Selects the scales (bandwidth, scaleVec) for the kernel density estimation, of the sequence. More...
 
void gaussianKDE (const V &evaluationParamVec, V &densityVec) const
 Gaussian kernel for the KDE estimate of the sequence. More...
 
void gaussianKDE (unsigned int initialPos, const V &scales, const std::vector< V * > &evaluationParamVecs, std::vector< V * > &densityVecs) const
 TODO: Gaussian kernel for the KDE of the sequence. More...
 
void writeContents (std::ofstream &ofsvar) const
 Write contents of the chain to a file. More...
 
void unifiedWriteContents (std::ofstream &ofsvar) const
 Writes info of the unified sequence to a file. More...
 
void unifiedWriteContents (const std::string &fileName, const std::string &fileType) const
 Writes info of the unified sequence to a file. More...
 
void unifiedReadContents (const std::string &fileName, const std::string &fileType, const unsigned int subSequenceSize)
 Reads the unified sequence from a file. More...
 
void select (const std::vector< unsigned int > &idsOfUniquePositions)
 Select positions in the sequence of vectors. More...
 
void filter (unsigned int initialPos, unsigned int spacing)
 Filters positions in the sequence of vectors, starting at initialPos, and with spacing given by spacing. More...
 
void extractScalarSeq (unsigned int initialPos, unsigned int spacing, unsigned int numPos, unsigned int paramId, ScalarSequence< double > &scalarSeq) const
 Extracts a sequence of scalars of size numPos, from position paramId of the array of sequences. More...
 
- Public Member Functions inherited from QUESO::BaseVectorSequence< V, M >
 BaseVectorSequence (const VectorSpace< V, M > &vectorSpace, unsigned int subSequenceSize, const std::string &name)
 Default constructor. More...
 
virtual ~BaseVectorSequence ()
 Destructor. More...
 
unsigned int unifiedSequenceSize () const
 Calculates the size of the unified sequence of vectors. More...
 
unsigned int vectorSizeLocal () const
 Local dimension (size) of the vector space. More...
 
unsigned int vectorSizeGlobal () const
 Global dimension (size) of the vector space. More...
 
const VectorSpace< V, M > & vectorSpace () const
 Vector space; access to protected attribute VectorSpace<V,M>& m_vectorSpace. More...
 
const std::string & name () const
 Access to protected attribute m_name: name of the sequence of vectors. More...
 
void setName (const std::string &newName)
 Changes the name of the sequence of vectors. More...
 
void clear ()
 Reset the values and the size of the sequence of vectors. More...
 
const V & subMinPlain () const
 Finds the minimum value of the sub-sequence. More...
 
const V & unifiedMinPlain () const
 Finds the minimum value of the unified sequence. More...
 
const V & subMaxPlain () const
 Finds the maximum value of the sub-sequence. More...
 
const V & unifiedMaxPlain () const
 Finds the maximum value of the unified sequence. More...
 
const V & subMeanPlain () const
 Finds the mean value of the sub-sequence. More...
 
const V & unifiedMeanPlain () const
 Finds the mean value of the unified sequence. More...
 
const V & subMedianPlain () const
 Finds the median value of the sub-sequence. More...
 
const V & unifiedMedianPlain () const
 Finds the median value of the unified sequence. More...
 
const V & subSampleVariancePlain () const
 Finds the variance of a sample of the sub-sequence. More...
 
const V & unifiedSampleVariancePlain () const
 Finds the variance of a sample of the unified sequence. More...
 
const BoxSubset< V, M > & subBoxPlain () const
 Finds a box subset of the sub-sequence (given by its min and max values calculated via subMinPlain and subMaxPlain). More...
 
const BoxSubset< V, M > & unifiedBoxPlain () const
 Finds a box subset of the unified-sequence (given by the min and max values of the unified sequence calculated via unifiedMinPlain and unifiedMaxPlain). More...
 
void deleteStoredVectors ()
 Deletes all the stored vectors. More...
 
void append (const BaseVectorSequence< V, M > &src, unsigned int initialPos, unsigned int numPos)
 Appends the vector src to this vector. More...
 
double subPositionsOfMaximum (const ScalarSequence< double > &subCorrespondingScalarValues, BaseVectorSequence< V, M > &subPositionsOfMaximum)
 Finds the positions where the maximum element occurs in the sub-sequence. More...
 
double unifiedPositionsOfMaximum (const ScalarSequence< double > &subCorrespondingScalarValues, BaseVectorSequence< V, M > &unifiedPositionsOfMaximum)
 Finds the positions where the maximum element occurs in the unified sequence. More...
 
void setGaussian (const V &meanVec, const V &stdDevVec)
 Sets the values of the sequence as a Gaussian distribution of mean given by meanVec and standard deviation by stdDevVec. More...
 
void setUniform (const V &aVec, const V &bVec)
 Sets the values of the sequence as a uniform distribution between the values given by vectors aVec and bVec. More...
 
virtual void subMeanExtra (unsigned int initialPos, unsigned int numPos, V &meanVec) const =0
 Finds the mean value of the sub-sequence, considering numPos positions starting at position initialPos. See template specialization. More...
 
virtual void unifiedMeanExtra (unsigned int initialPos, unsigned int numPos, V &unifiedMeanVec) const =0
 Finds the mean value of the unified sequence of numPos positions starting at position initialPos. See template specialization. More...
 
virtual void subMedianExtra (unsigned int initialPos, unsigned int numPos, V &medianVec) const =0
 Finds the median value of the sub-sequence, considering numPos positions starting at position initialPos. See template specialization. More...
 
virtual void unifiedMedianExtra (unsigned int initialPos, unsigned int localNumPos, V &unifiedMedianVec) const =0
 Finds the median value of the unified sequence, considering numPos positions starting at position initialPos. See template specialization. More...
 
virtual void subSampleVarianceExtra (unsigned int initialPos, unsigned int numPos, const V &meanVec, V &samVec) const =0
 Finds the sample variance of the sub-sequence, considering numPos positions starting at position initialPos and of mean meanVec. See template specialization. More...
 
virtual void unifiedSampleVarianceExtra (unsigned int initialPos, unsigned int numPos, const V &unifiedMeanVec, V &unifiedSamVec) const =0
 Finds the sample variance of the unified sequence, considering numPos positions starting at position initialPos and of mean meanVec. See template specialization. More...
 
virtual void subPopulationVariance (unsigned int initialPos, unsigned int numPos, const V &meanVec, V &popVec) const =0
 Finds the population variance of the sub-sequence, considering numPos positions starting at position initialPos and of mean meanVec. See template specialization. More...
 
virtual void unifiedPopulationVariance (unsigned int initialPos, unsigned int numPos, const V &unifiedMeanVec, V &unifiedPopVec) const =0
 Finds the population variance of the unified-sequence, considering numPos positions starting at position initialPos and of mean meanVec. See template specialization. More...
 
virtual void subMinMaxExtra (unsigned int initialPos, unsigned int numPos, V &minVec, V &maxVec) const =0
 Finds the minimum and the maximum values of the sub-sequence, considering numPos positions starting at position initialPos. See template specialization. More...
 
virtual void unifiedMinMaxExtra (unsigned int initialPos, unsigned int numPos, V &unifiedMinVec, V &unifiedMaxVec) const =0
 Finds the minimum and the maximum values of the unified sequence, considering numPos positions starting at position initialPos. See template specialization. More...
 
virtual void subHistogram (unsigned int initialPos, const V &minVec, const V &maxVec, std::vector< V * > &centersForAllBins, std::vector< V * > &quanttsForAllBins) const =0
 Calculates the histogram of the sub-sequence. See template specialization. More...
 
virtual void unifiedHistogram (unsigned int initialPos, const V &unifiedMinVec, const V &unifiedMaxVec, std::vector< V * > &unifiedCentersForAllBins, std::vector< V * > &unifiedQuanttsForAllBins) const =0
 Calculates the histogram of the unified sequence. See template specialization. More...
 
virtual void subInterQuantileRange (unsigned int initialPos, V &iqrVec) const =0
 Returns the interquartile range of the values in the sub-sequence. See template specialization. More...
 
virtual void unifiedInterQuantileRange (unsigned int initialPos, V &unifiedIqrVec) const =0
 Returns the interquartile range of the values in the unified sequence. See template specialization. More...
 
virtual void subScalesForKde (unsigned int initialPos, const V &iqrVec, unsigned int kdeDimension, V &scaleVec) const =0
 Selects the scales (bandwidth, scaleVec) for the kernel density estimation, considering only the sub-sequence. See template specialization. More...
 
virtual void unifiedScalesForKde (unsigned int initialPos, const V &unifiedIqrVec, unsigned int kdeDimension, V &unifiedScaleVec) const =0
 Selects the scales (bandwidth) for the kernel density estimation, considering the unified sequence. See template specialization. More...
 
virtual void subGaussian1dKde (unsigned int initialPos, const V &scaleVec, const std::vector< V * > &evaluationParamVecs, std::vector< V * > &densityVecs) const =0
 Gaussian kernel for the KDE estimate of the sub-sequence. See template specialization. More...
 
virtual void unifiedGaussian1dKde (unsigned int initialPos, const V &unifiedScaleVec, const std::vector< V * > &unifiedEvaluationParamVecs, std::vector< V * > &unifiedDensityVecs) const =0
 Gaussian kernel for the KDE estimate of the unified sequence. See template specialization. More...
 
virtual void subWriteContents (unsigned int initialPos, unsigned int numPos, const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds) const =0
 Writes info of the sub-sequence to a file. See template specialization. More...
 
virtual void subWriteContents (unsigned int initialPos, unsigned int numPos, std::ofstream &ofsvar, const std::string &fileType) const =0
 Writes info of the sub-sequence to a file. See template specialization. More...
 
void computeFilterParams (std::ofstream *passedOfs, unsigned int &initialPos, unsigned int &spacing)
 Computes the filtering parameters spacing for the sequence of vectors. More...
 
virtual double estimateConvBrooksGelman (unsigned int initialPos, unsigned int numPos) const =0
 Estimates convergence rate using Brooks & Gelman method. See template specialization. More...
 

Private Member Functions

void extractRawData (unsigned int initialPos, unsigned int spacing, unsigned int numPos, unsigned int paramId, std::vector< double > &rawData) const
 Extracts the raw data. More...
 

Private Attributes

DistArray< ScalarSequence
< double > * > 
m_scalarSequences
 Sequence of scalars. More...
 

Additional Inherited Members

- Protected Member Functions inherited from QUESO::BaseVectorSequence< V, M >
void copy (const BaseVectorSequence< V, M > &src)
 Copies vector sequence src to this. More...
 
- Protected Attributes inherited from QUESO::BaseVectorSequence< V, M >
const BaseEnvironmentm_env
 
const VectorSpace< V, M > & m_vectorSpace
 
std::string m_name
 
Fft< double > * m_fftObj
 
V * m_subMinPlain
 
V * m_unifiedMinPlain
 
V * m_subMaxPlain
 
V * m_unifiedMaxPlain
 
V * m_subMeanPlain
 
V * m_unifiedMeanPlain
 
V * m_subMedianPlain
 
V * m_unifiedMedianPlain
 
V * m_subSampleVariancePlain
 
V * m_unifiedSampleVariancePlain
 
BoxSubset< V, M > * m_subBoxPlain
 
BoxSubset< V, M > * m_unifiedBoxPlain
 

Detailed Description

template<class V, class M>
class QUESO::ArrayOfSequences< V, M >

Class for handling array samples (arrays of scalar sequences).

This class handles array of sequences (samples) generated by an algorithm, as well as operations that can be carried over them, e.g., calculation of means, correlation and covariance matrices. It is derived from and implements BaseVectorSequence<V,M>.

Definition at line 44 of file ArrayOfSequences.h.

Constructor & Destructor Documentation

template<class V , class M >
QUESO::ArrayOfSequences< V, M >::ArrayOfSequences ( const VectorSpace< V, M > &  vectorSpace,
unsigned int  subSequenceSize,
const std::string &  name 
)

Default constructor.

Definition at line 32 of file ArrayOfSequences.C.

References QUESO::BaseVectorSequence< V, M >::m_env, QUESO::ArrayOfSequences< V, M >::m_scalarSequences, and QUESO::ArrayOfSequences< V, M >::subSequenceSize().

35  : BaseVectorSequence<V,M>(vectorSpace,subSequenceSize,name),
37 {
38 
39  //if (m_env.subDisplayFile()) {
40  // *m_env.subDisplayFile() << "Entering ArrayOfSequences<V,M>::constructor()"
41  // << std::endl;
42  //}
43 
44  //if (m_env.subDisplayFile()) {
45  // *m_env.subDisplayFile() << "In ArrayOfSequences<V,M>::constructor()"
46  // << "\n subSequenceSize = " << subSequenceSize
47  // << "\n m_scalarSequences.MyLength() = " << m_scalarSequences.MyLength()
48  // << std::endl;
49  //}
50 
51  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
52  m_scalarSequences(i,0) = new ScalarSequence<double>(m_env,subSequenceSize);
53  }
54 
55  //if (m_env.subDisplayFile()) {
56  // *m_env.subDisplayFile() << "Leaving ArrayOfSequences<V,M>::constructor()"
57  // << std::endl;
58  //}
59 }
const VectorSpace< V, M > & m_vectorSpace
const BaseEnvironment & m_env
const std::string & name() const
Access to protected attribute m_name: name of the sequence of vectors.
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
const VectorSpace< V, M > & vectorSpace() const
Vector space; access to protected attribute VectorSpace&lt;V,M&gt;&amp; m_vectorSpace.
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
template<class V , class M >
QUESO::ArrayOfSequences< V, M >::~ArrayOfSequences ( )

Destructor.

Definition at line 63 of file ArrayOfSequences.C.

64 {
65  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
66  if (m_scalarSequences(i,0)) delete m_scalarSequences(i,0);
67  }
68 }
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.

Member Function Documentation

template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::autoCorrViaDef ( unsigned int  initialPos,
unsigned int  numPos,
unsigned int  lag,
V &  corrVec 
) const
virtual

Calculates autocorrelation via definition.

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 346 of file ArrayOfSequences.C.

350 {
351  V subChainMean (m_vectorSpace.zeroVector());
352  V subChainAutoCovarianceLag0(m_vectorSpace.zeroVector());
353 
354  this->mean(initialPos,
355  numPos,
356  subChainMean);
357  this->autoCovariance(initialPos,
358  numPos,
359  subChainMean,
360  0, // lag
361  subChainAutoCovarianceLag0);
362 
363  this->autoCovariance(initialPos,
364  numPos,
365  subChainMean,
366  lag,
367  corrVec);
368  corrVec /= subChainAutoCovarianceLag0;
369 
370  return;
371 }
const VectorSpace< V, M > & m_vectorSpace
void mean(unsigned int initialPos, unsigned int numPos, V &meanVec) const
Finds the mean value of the sub-sequence, considering numPos positions starting at position initialPo...
void autoCovariance(unsigned int initialPos, unsigned int numPos, const V &meanVec, unsigned int lag, V &covVec) const
Calculates the autocovariance.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::autoCorrViaFft ( unsigned int  initialPos,
unsigned int  numPos,
const std::vector< unsigned int > &  lags,
std::vector< V * > &  corrVecs 
) const
virtual

Calculates autocorrelation via Fast Fourier transforms (FFT). TODO: Implement me!

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 374 of file ArrayOfSequences.C.

378 {
379  return;
380 }
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::autoCorrViaFft ( unsigned int  initialPos,
unsigned int  numPos,
unsigned int  numSum,
V &  autoCorrsSumVec 
) const
virtual

Calculates autocorrelation via Fast Fourier transforms (FFT). TODO: Implement me!

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 383 of file ArrayOfSequences.C.

References UQ_FATAL_TEST_MACRO.

387 {
388 #if 0
389  bool bRC = ((initialPos < this->subSequenceSize()) &&
390  (0 < numPos ) &&
391  ((initialPos+numPos) <= this->subSequenceSize()) &&
392  (autoCorrsSumVec.size() == this->vectorSize() ));
393  UQ_FATAL_TEST_MACRO(bRC == false,
394  m_env.worldRank(),
395  "ArrayOfSequences<V,M>::autoCorrViaFft(), for sum",
396  "invalid input data");
397 
398  ScalarSequence<double> data(m_env,0);
399 
400  unsigned int numParams = this->vectorSize();
401  for (unsigned int i = 0; i < numParams; ++i) {
402  this->extractScalarSeq(initialPos,
403  1, // spacing
404  numPos,
405  i,
406  data);
407 
408  data.autoCorrViaFft(0,
409  numPos,
410  autoCorrsSumVec[i]);
411  }
412 #endif
413  return;
414 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
void extractScalarSeq(unsigned int initialPos, unsigned int spacing, unsigned int numPos, unsigned int paramId, ScalarSequence< double > &scalarSeq) const
Extracts a sequence of scalars of size numPos, from position paramId of the array of sequences...
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::autoCovariance ( unsigned int  initialPos,
unsigned int  numPos,
const V &  meanVec,
unsigned int  lag,
V &  covVec 
) const
virtual

Calculates the autocovariance.

The autocovariance is the covariance of a variable with itself at some other time. It is calculated over a sequence of arrays with initial position initialPos, considering numPos positions, a lag of lag, with mean given by meanVec. The results are saved in the output vector covVec/

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 292 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences, and UQ_FATAL_TEST_MACRO.

297 {
298  bool bRC = ((0 <= initialPos ) &&
299  (0 < numPos ) &&
300  ((initialPos+numPos-1) <= (this->subSequenceSize()-1)));
301  UQ_FATAL_TEST_MACRO(bRC == false,
302  m_env.worldRank(),
303  "VectorSequenceAutoCovariance<V,M>()",
304  "invalid initial position or number of positions");
305 
306  bRC = (numPos > lag);
307  UQ_FATAL_TEST_MACRO(bRC == false,
308  m_env.worldRank(),
309  "VectorSequenceAutoCovariance<V,M>()",
310  "lag is too large");
311 
312  bRC = (this->vectorSize() == covVec.size());
313  UQ_FATAL_TEST_MACRO(bRC == false,
314  m_env.worldRank(),
315  "VectorSequenceAutoCovariance<V,M>()",
316  "incompatible sizes between covVec vector and vectors in sequence");
317 
318  bRC = (this->vectorSize() == meanVec.size());
319  UQ_FATAL_TEST_MACRO(bRC == false,
320  m_env.worldRank(),
321  "VectorSequenceAutoCovariance<V,M>()",
322  "incompatible sizes between meanVec vector and vectors in sequence");
323 
324  unsigned int loopSize = numPos - lag;
325  unsigned int finalPosPlus1 = initialPos + loopSize;
326  double doubleLoopSize = (double) loopSize;
327  covVec.cwSet(0.);
328 
329  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
330  for (unsigned int i = 0; i < covVec.size(); ++i) {
331  const ScalarSequence<double>& seq = *(tmp->m_scalarSequences(i,0));
332  double meanValue = meanVec[i];
333  double result = 0;
334  for (unsigned int j = initialPos; j < finalPosPlus1; ++j) {
335  double diff1 = seq[j] - meanValue;
336  double diff2 = seq[j+lag] - meanValue;
337  result += diff1*diff2;
338  }
339  covVec[i] = result/doubleLoopSize;
340  }
341 
342  return;
343 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::erasePositions ( unsigned int  initialPos,
unsigned int  numPos 
)
virtual

Erases numPos elements of the sequence starting at position initialPos.

This routine deletes all stored computed vectors

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 106 of file ArrayOfSequences.C.

References QUESO::BaseVectorSequence< V, M >::deleteStoredVectors(), and QUESO::ScalarSequence< T >::erasePositions().

108 {
109  if (initialPos < this->subSequenceSize()) {
110  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
111  ScalarSequence<double>& seq = *(m_scalarSequences(i,0));
112  seq.erasePositions(initialPos,numPos);
113  }
114  }
115 
117 
118  return;
119 }
void deleteStoredVectors()
Deletes all the stored vectors.
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::extractRawData ( unsigned int  initialPos,
unsigned int  spacing,
unsigned int  numPos,
unsigned int  paramId,
std::vector< double > &  rawData 
) const
privatevirtual

Extracts the raw data.

Extracts the data from sequence of array at position paramId, with spacing spacing, and saves is in rawData. The vector rawData will have size numPos. Note that, in fact, spacing = 1.

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 673 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences.

678 {
679  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
680  ScalarSequence<double>& seq = *(tmp->m_scalarSequences(paramId,0));
681 
682  rawData.resize(numPos);
683  if (spacing == 1) {
684  for (unsigned int j = 0; j < numPos; ++j) {
685  rawData[j] = seq[paramId];
686  }
687  }
688  else {
689  for (unsigned int j = 0; j < numPos; ++j) {
690  rawData[j] = seq[paramId];
691  }
692  }
693 
694  return;
695 }
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::extractScalarSeq ( unsigned int  initialPos,
unsigned int  spacing,
unsigned int  numPos,
unsigned int  paramId,
ScalarSequence< double > &  scalarSeq 
) const
virtual

Extracts a sequence of scalars of size numPos, from position paramId of the array of sequences.

In this method, the parameter initialPos is unused and the parameter spacing is always equal to one

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 647 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences, and QUESO::ScalarSequence< T >::resizeSequence().

652 {
653  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
654  ScalarSequence<double>& seq = *(tmp->m_scalarSequences(paramId,0));
655 
656  scalarSeq.resizeSequence(numPos);
657  if (spacing == 1) {
658  for (unsigned int j = 0; j < numPos; ++j) {
659  scalarSeq[j] = seq[paramId];
660  }
661  }
662  else {
663  for (unsigned int j = 0; j < numPos; ++j) {
664  scalarSeq[j] = seq[paramId];
665  }
666  }
667 
668  return;
669 }
void resizeSequence(unsigned int newSequenceSize)
Resizes the size of the sequence of scalars.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::filter ( unsigned int  initialPos,
unsigned int  spacing 
)
virtual

Filters positions in the sequence of vectors, starting at initialPos, and with spacing given by spacing.

TODO: implement me!

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 640 of file ArrayOfSequences.C.

642 {
643  return;
644 }
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::gaussianKDE ( const V &  evaluationParamVec,
V &  densityVec 
) const

Gaussian kernel for the KDE estimate of the sequence.

TODO: implement me!

Definition at line 539 of file ArrayOfSequences.C.

541 {
542  return;
543 }
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::gaussianKDE ( unsigned int  initialPos,
const V &  scales,
const std::vector< V * > &  evaluationParamVecs,
std::vector< V * > &  densityVecs 
) const

TODO: Gaussian kernel for the KDE of the sequence.

Todo:
: implement me! The density estimator will be:

\[ \hat{f}(x) = \frac{1}{nh} \sum_{j=1}^n K\Big(\frac{x-x_j}{h}\Big),\]

where $ K $ is a Gaussian kernel.

Definition at line 546 of file ArrayOfSequences.C.

References QUESO::MiscGaussianDensity().

550 {
551 #if 0
552  unsigned int dataSize = sequence.size() - initialPos;
553  unsigned int numEstimationsPerParam = evaluationParamVecs.size();
554 
555  for (unsigned int j = 0; j < numEstimationsPerParam; ++j) {
556  densityVecs[j] = new V(*(sequence[0]));
557  }
558 
559  unsigned int numParams = sequence[0]->size();
560  for (unsigned int i = 0; i < numParams; ++i) {
561  double scaleInv = 1./scales[i];
562  for (unsigned int j = 0; j < numEstimationsPerParam; ++j) {
563  double x = (*(evaluationParamVecs[j]))[i];
564  double value = 0.;
565  for (unsigned int k = 0; k < dataSize; ++k) {
566  double xk = (*(sequence[initialPos+k]))[i];
567  value += MiscGaussianDensity((x-xk)*scaleInv,0.,1.);
568  }
569  (*(densityVecs[j]))[i] = scaleInv * (value/(double) numEstimationsPerParam);
570  }
571  }
572 
573 #endif
574  return;
575 }
double MiscGaussianDensity(double x, double mu, double sigma)
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::getPositionValues ( unsigned int  posId,
V &  vec 
) const
virtual

Gets the values of the sequence at position posId and stores them at vec.

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 122 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences.

123 {
124  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
125  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
126  vec[i] = (*(tmp->m_scalarSequences(i,0)))[posId];
127  }
128 
129  return;
130 }
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::histogram ( unsigned int  initialPos,
const V &  minVec,
const V &  maxVec,
std::vector< V * > &  centersForAllBins,
std::vector< V * > &  quanttsForAllBins 
) const

Calculates the histogram of the sequence.

TODO: implement me!

Definition at line 432 of file ArrayOfSequences.C.

References UQ_FATAL_TEST_MACRO.

437 {
438 #if 0
439  UQ_FATAL_TEST_MACRO(centersForAllBins.size() != quanttsForAllBins.size(),
440  sequence[0]->env().worldRank(),
441  "VectorSequenceHistogram<V,M>()",
442  "vectors 'centers' and 'quantts' have different sizes");
443 
444  for (unsigned int j = 0; j < quanttsForAllBins.size(); ++j) {
445  centersForAllBins[j] = new V(*(sequence[0]));
446  quanttsForAllBins [j] = new V(*(sequence[0]));
447  }
448 
449  unsigned int dataSize = sequence.size() - initialPos;
450  unsigned int numParams = sequence[0]->size();
451  for (unsigned int i = 0; i < numParams; ++i) {
452  ScalarSequence<double> data(dataSize,0.);
453  for (unsigned int j = 0; j < dataSize; ++j) {
454  data[j] = (*(sequence[initialPos+j]))[i];
455  }
456 
457  std::vector<double> centers(centersForAllBins.size(),0.);
458  std::vector<double> quantts(quanttsForAllBins.size(),0.);
459  data.histogram(minVec[i],
460  maxVec[i],
461  centers,
462  quantts);
463 
464  for (unsigned int j = 0; j < quantts.size(); ++j) {
465  (*(centersForAllBins[j]))[i] = centers[j];
466  (*(quanttsForAllBins[j]))[i] = quantts[j];
467  }
468  }
469 
470 #endif
471  return;
472 }
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::interQuantileRange ( unsigned int  initialPos,
V &  iqrs 
) const

Returns the interquartile range of the values in the sequence.

TODO: implement me!

Definition at line 475 of file ArrayOfSequences.C.

477 {
478 #if 0
479  unsigned int dataSize = sequence.size() - initialPos;
480 
481  ArrayOfSequences sortedSequence(dataSize,m_vectorSpace.zeroVector());
482  this->sort(initialPos,
483  sortedSequence);
484 
485  unsigned int pos1 = (unsigned int) ( (((double) dataSize) + 1.)*1./4. - 1. );
486  unsigned int pos3 = (unsigned int) ( (((double) dataSize) + 1.)*3./4. - 1. );
487 
488  double fraction1 = (((double) dataSize) + 1.)*1./4. - 1. - ((double) pos1);
489  double fraction3 = (((double) dataSize) + 1.)*3./4. - 1. - ((double) pos3);
490 
491  unsigned int numParams = sequence[0]->size();
492  for (unsigned int i = 0; i < numParams; ++i) {
493  double value1 = (1.-fraction1) * (*sortedSequence[pos1])[i] + fraction1 * (*sortedSequence[pos1+1])[i];
494  double value3 = (1.-fraction3) * (*sortedSequence[pos3])[i] + fraction3 * (*sortedSequence[pos3+1])[i];
495  iqrs[i] = value3 - value1;
496  }
497 
498 #endif
499  return;
500 }
const VectorSpace< V, M > & m_vectorSpace
ArrayOfSequences(const VectorSpace< V, M > &vectorSpace, unsigned int subSequenceSize, const std::string &name)
Default constructor.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::mean ( unsigned int  initialPos,
unsigned int  numPos,
V &  meanVec 
) const

Finds the mean value of the sub-sequence, considering numPos positions starting at position initialPos.

Definition at line 172 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences, QUESO::ScalarSequence< T >::subMeanExtra(), and UQ_FATAL_TEST_MACRO.

174 {
175  bool bRC = ((0 <= initialPos ) &&
176  (0 < numPos ) &&
177  ((initialPos+numPos-1) <= (this->subSequenceSize()-1)));
178  UQ_FATAL_TEST_MACRO(bRC == false,
179  m_env.worldRank(),
180  "ArrayOfSequences<V,M>::mean()",
181  "invalid initial position or number of positions");
182 
183  bRC = (this->vectorSize() == meanVec.size());
184  UQ_FATAL_TEST_MACRO(bRC == false,
185  m_env.worldRank(),
186  "ArrayOfSequences<V,M>::mean()",
187  "incompatible sizes between meanVec vector and vectors in sequence");
188 
189  meanVec.cwSet(0.);
190  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
191  for (unsigned int i = 0; i < meanVec.size(); ++i) {
192  const ScalarSequence<double>& seq = *(tmp->m_scalarSequences(i,0));
193  meanVec[i] = seq.subMeanExtra(initialPos, numPos);
194  }
195 
196  return;
197 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::minMax ( unsigned int  initialPos,
V &  minVec,
V &  maxVec 
) const

Given an initial position initialPos, finds the minimum and the maximum values of the sequence.

Definition at line 417 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences, and QUESO::ScalarSequence< T >::subMinMaxExtra().

419 {
420  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
421 
422  unsigned int numParams = this->vectorSize();
423  for (unsigned int i = 0; i < numParams; ++i) {
424  ScalarSequence<double>& seq = *(tmp->m_scalarSequences(i,0));
425  seq.subMinMaxExtra(initialPos,minVec[i],maxVec[i]);
426  }
427 
428  return;
429 }
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::populationVariance ( unsigned int  initialPos,
unsigned int  numPos,
const V &  meanVec,
V &  popVec 
) const

Finds the population variance of the sub-sequence, considering numPos positions starting at position initialPos and of mean meanVec.

Definition at line 246 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences, and UQ_FATAL_TEST_MACRO.

250 {
251  bool bRC = ((0 <= initialPos ) &&
252  (0 < numPos ) &&
253  ((initialPos+numPos-1) <= (this->subSequenceSize()-1)));
254  UQ_FATAL_TEST_MACRO(bRC == false,
255  m_env.worldRank(),
256  "ArrayOfSequences<V,M>::populationVariance()",
257  "invalid initial position or number of positions");
258 
259  bRC = (this->vectorSize() == popVec.size());
260  UQ_FATAL_TEST_MACRO(bRC == false,
261  m_env.worldRank(),
262  "ArrayOfSequences<V,M>::populationVariance()",
263  "incompatible sizes between popVec vector and vectors in sequence");
264 
265  bRC = (this->vectorSize() == meanVec.size());
266  UQ_FATAL_TEST_MACRO(bRC == false,
267  m_env.worldRank(),
268  "ArrayOfSequences<V,M>::populationVariance()",
269  "incompatible sizes between meanVec vector and vectors in sequence");
270 
271  unsigned int loopSize = numPos;
272  unsigned int finalPosPlus1 = initialPos + loopSize;
273  double doubleLoopSize = (double) loopSize;
274  popVec.cwSet(0.);
275 
276  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
277  for (unsigned int i = 0; i < popVec.size(); ++i) {
278  const ScalarSequence<double>& seq = *(tmp->m_scalarSequences(i,0));
279  double tmpMean = meanVec[i];
280  double result = 0.;
281  for (unsigned int j = initialPos; j < finalPosPlus1; ++j) {
282  double diff = seq[j] - tmpMean;
283  result += diff*diff;
284  }
285  popVec[i] = result/doubleLoopSize;
286  }
287 
288  return;
289 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::resetValues ( unsigned int  initialPos,
unsigned int  numPos 
)
virtual

Resets a total of numPos values of the sequence starting at position initialPos.

This routine deletes all stored computed vectors

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 93 of file ArrayOfSequences.C.

References QUESO::BaseVectorSequence< V, M >::deleteStoredVectors().

95 {
96  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
97  m_scalarSequences(i,0)->resetValues(initialPos,numPos);
98  }
99 
101 
102  return;
103 }
void deleteStoredVectors()
Deletes all the stored vectors.
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::resizeSequence ( unsigned int  newSubSequenceSize)
virtual

Resizes the sequence.

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 80 of file ArrayOfSequences.C.

References QUESO::BaseVectorSequence< V, M >::deleteStoredVectors().

81 {
82  if (newSubSequenceSize != this->subSequenceSize()) {
83  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
84  m_scalarSequences(i,0)->resizeSequence(newSubSequenceSize);
85  }
86  }
87 
89  return;
90 }
void deleteStoredVectors()
Deletes all the stored vectors.
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::sampleVariance ( unsigned int  initialPos,
unsigned int  numPos,
const V &  meanVec,
V &  samVec 
) const

Finds the sample variance of the sub-sequence, considering numPos positions starting at position initialPos and of mean meanVec.

Definition at line 200 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences, and UQ_FATAL_TEST_MACRO.

204 {
205  bool bRC = ((0 <= initialPos ) &&
206  (0 < numPos ) &&
207  ((initialPos+numPos-1) <= (this->subSequenceSize()-1)));
208  UQ_FATAL_TEST_MACRO(bRC == false,
209  m_env.worldRank(),
210  "ArrayOfSequences<V,M>::sampleVariance()",
211  "invalid initial position or number of positions");
212 
213  bRC = (this->vectorSize() == samVec.size());
214  UQ_FATAL_TEST_MACRO(bRC == false,
215  m_env.worldRank(),
216  "ArrayOfSequences<V,M>::sampleVariance()",
217  "incompatible sizes between samVec vector and vectors in sequence");
218 
219  bRC = (this->vectorSize() == meanVec.size());
220  UQ_FATAL_TEST_MACRO(bRC == false,
221  m_env.worldRank(),
222  "ArrayOfSequences<V,M>::sampleVariance()",
223  "incompatible sizes between meanVec vector and vectors in sequence");
224 
225  unsigned int loopSize = numPos;
226  unsigned int finalPosPlus1 = initialPos + loopSize;
227  double doubleLoopSize = (double) loopSize;
228  samVec.cwSet(0.);
229 
230  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
231  for (unsigned int i = 0; i < samVec.size(); ++i) {
232  const ScalarSequence<double>& seq = *(tmp->m_scalarSequences(i,0));
233  double tmpMean = meanVec[i];
234  double result = 0.;
235  for (unsigned int j = initialPos; j < finalPosPlus1; ++j) {
236  double diff = seq[j] - tmpMean;
237  result += diff*diff;
238  }
239  samVec[i] = result/(doubleLoopSize - 1.);
240  }
241 
242  return;
243 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::scalesForKDE ( unsigned int  initialPos,
const V &  iqrs,
unsigned int  kdeDimension,
V &  scales 
) const

Selects the scales (bandwidth, scaleVec) for the kernel density estimation, of the sequence.

TODO: implement me!

Definition at line 503 of file ArrayOfSequences.C.

507 {
508 #if 0
509  unsigned int dataSize = sequence.size() - initialPos;
510 
511  V mean(*(sequence[0]));
512  VectorSequenceMean(sequence,
513  initialPos,
514  dataSize,
515  mean);
516 
517  V samVec(*(sequence[0]));
518  VectorSequenceSampleVariance(sequence,
519  initialPos,
520  dataSize,
521  mean,
522  samVec);
523 
524  unsigned int numParams = sequence[0]->size();
525  for (unsigned int i = 0; i < numParams; ++i) {
526  if (iqrs[i] <= 0.) {
527  scales[i] = 1.06*std::sqrt(samVec[i])/std::pow(dataSize,1./5.);
528  }
529  else {
530  scales[i] = 1.06*std::min(std::sqrt(samVec[i]),iqrs[i]/1.34)/std::pow(dataSize,1./5.);
531  }
532  }
533 
534 #endif
535  return;
536 }
void mean(unsigned int initialPos, unsigned int numPos, V &meanVec) const
Finds the mean value of the sub-sequence, considering numPos positions starting at position initialPo...
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::select ( const std::vector< unsigned int > &  idsOfUniquePositions)
virtual

Select positions in the sequence of vectors.

TODO: implement me!

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 633 of file ArrayOfSequences.C.

635 {
636  return;
637 }
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::setGaussian ( const V &  meanVec,
const V &  stdDevVec 
)

Sets the values of the sequence as a Gaussian distribution of mean given by meanVec and standard deviation by stdDevVec.

This routine deletes all stored computed vectors

Definition at line 146 of file ArrayOfSequences.C.

References QUESO::BaseVectorSequence< V, M >::deleteStoredVectors(), and QUESO::ScalarSequence< T >::setGaussian().

147 {
148  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
149  ScalarSequence<double>& seq = *(m_scalarSequences(i,0));
150  seq.setGaussian(meanVec[i],stdDevVec[i]);
151  }
152 
154 
155  return;
156 }
void deleteStoredVectors()
Deletes all the stored vectors.
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::setPositionValues ( unsigned int  posId,
const V &  vec 
)
virtual

Set the values of vec at position posId of the sequence.

This routine deletes all stored computed vectors

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 133 of file ArrayOfSequences.C.

References QUESO::BaseVectorSequence< V, M >::deleteStoredVectors().

134 {
135  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
136  ScalarSequence<double>& seq = *(m_scalarSequences(i,0));
137  seq[posId] = vec[i];
138  }
139 
141 
142  return;
143 }
void deleteStoredVectors()
Deletes all the stored vectors.
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::setUniform ( const V &  aVec,
const V &  bVec 
)

Sets the values of the sequence as a uniform distribution between the values given by vectors aVec and bVec.

This routine deletes all stored computed vectors

Definition at line 159 of file ArrayOfSequences.C.

References QUESO::BaseVectorSequence< V, M >::deleteStoredVectors(), and QUESO::ScalarSequence< T >::setUniform().

160 {
161  for (unsigned int i = 0; i < (unsigned int) m_scalarSequences.MyLength(); ++i) {
162  ScalarSequence<double>& seq = *(m_scalarSequences(i,0));
163  seq.setUniform(aVec[i],bVec[i]);
164  }
165 
167 
168  return;
169 }
void deleteStoredVectors()
Deletes all the stored vectors.
DistArray< ScalarSequence< double > * > m_scalarSequences
Sequence of scalars.
template<class V , class M >
unsigned int QUESO::ArrayOfSequences< V, M >::subSequenceSize ( ) const
virtual

Size of the sub-sequence of arrays.

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 72 of file ArrayOfSequences.C.

References QUESO::ArrayOfSequences< V, M >::m_scalarSequences.

Referenced by QUESO::ArrayOfSequences< V, M >::ArrayOfSequences().

73 {
74  ArrayOfSequences<V,M>* tmp = const_cast<ArrayOfSequences<V,M>*>(this);
75 
76  return tmp->m_scalarSequences(0,0)->subSequenceSize();
77 }
template<class V, class M>
void QUESO::ArrayOfSequences< V, M >::unifiedMean ( unsigned int  initialPos,
unsigned int  numPos,
V &  unifiedMeanVec 
) const

Finds the mean value of the unified sequence of numPos positions starting at position initialPos. See template specialization.

template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::unifiedReadContents ( const std::string &  fileName,
const std::string &  fileType,
const unsigned int  subSequenceSize 
)
virtual

Reads the unified sequence from a file.

TODO: implement me!

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 621 of file ArrayOfSequences.C.

References UQ_FATAL_TEST_MACRO.

624 {
625  UQ_FATAL_TEST_MACRO(true,
626  m_env.worldRank(),
627  "ArrayOfSequences<V,M>::unifiedReadContents()",
628  "not implemented yet");
629  return;
630 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V, class M>
void QUESO::ArrayOfSequences< V, M >::unifiedSampleVariance ( unsigned int  initialPos,
unsigned int  numPos,
const V &  meanVec,
V &  samVec 
) const

Finds the sample variance of the unified sequence, considering numPos positions starting at position initialPos and of mean meanVec.

template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::unifiedWriteContents ( std::ofstream &  ofsvar) const

Writes info of the unified sequence to a file.

TODO: implement me!

Definition at line 600 of file ArrayOfSequences.C.

References UQ_FATAL_TEST_MACRO.

601 {
602  UQ_FATAL_TEST_MACRO(true,
603  m_env.worldRank(),
604  "ArrayOfSequences<V,M>::unifiedWriteContents(1)",
605  "not implemented yet");
606  return;
607 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::unifiedWriteContents ( const std::string &  fileName,
const std::string &  fileType 
) const
virtual

Writes info of the unified sequence to a file.

TODO: implement me!

Implements QUESO::BaseVectorSequence< V, M >.

Definition at line 610 of file ArrayOfSequences.C.

References UQ_FATAL_TEST_MACRO.

612 {
613  UQ_FATAL_TEST_MACRO(true,
614  m_env.worldRank(),
615  "ArrayOfSequences<V,M>::unifiedWriteContents(2)",
616  "not implemented yet");
617  return;
618 }
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
const BaseEnvironment & m_env
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
template<class V , class M >
void QUESO::ArrayOfSequences< V, M >::writeContents ( std::ofstream &  ofsvar) const

Write contents of the chain to a file.

Definition at line 578 of file ArrayOfSequences.C.

579 {
580  // Write chain
581  ofsvar << m_name << "_sub" << m_env.subIdString() << " = zeros(" << this->subSequenceSize()
582  << "," << this->vectorSize()
583  << ");"
584  << std::endl;
585  ofsvar << m_name << "_sub" << m_env.subIdString() << " = [";
586 
587  V tmpVec(m_vectorSpace.zeroVector());
588  unsigned int chainSize = this->subSequenceSize();
589  for (unsigned int j = 0; j < chainSize; ++j) {
590  this->getPositionValues(j,tmpVec);
591  ofsvar << tmpVec
592  << std::endl;
593  }
594  ofsvar << "];\n";
595 
596  return;
597 }
const std::string & subIdString() const
Access to the attribute m_subIdString; which stores the string for the sub-environment, and it will be used, for instance, to create the output files for each sub-environment.
Definition: Environment.C:335
const VectorSpace< V, M > & m_vectorSpace
const BaseEnvironment & m_env
unsigned int subSequenceSize() const
Size of the sub-sequence of arrays.
void getPositionValues(unsigned int posId, V &vec) const
Gets the values of the sequence at position posId and stores them at vec.

Member Data Documentation

template<class V, class M>
DistArray<ScalarSequence<double>*> QUESO::ArrayOfSequences< V, M >::m_scalarSequences
private

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

Generated on Thu Apr 23 2015 19:26:17 for queso-0.51.1 by  doxygen 1.8.5