queso-0.51.1
VectorRV.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008,2009,2010,2011,2012,2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 #include <queso/VectorRV.h>
26 #include <queso/GslVector.h>
27 #include <queso/GslMatrix.h>
28 
29 namespace QUESO {
30 
31 // Default constructor -----------------------------
32 template<class V, class M>
34  const char* prefix,
35  const VectorSet<V,M>& imageSet)
36  :
37  m_env (imageSet.env()),
38  m_prefix ((std::string)(prefix)+"rv_"),
39  m_imageSet (imageSet),
40  m_pdf (NULL),
41  m_realizer (NULL),
42  m_subCdf (NULL),
43  m_unifiedCdf(NULL),
44  m_mdf (NULL)
45 {
46  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 54)) {
47  *m_env.subDisplayFile() << "Entering BaseVectorRV<V,M>::constructor()"
48  << ": prefix = " << m_prefix
49  << std::endl;
50  }
51 
52  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 54)) {
53  *m_env.subDisplayFile() << "Leaving BaseVectorRV<V,M>::constructor()"
54  << ": prefix = " << m_prefix
55  << std::endl;
56  }
57 }
58 // Destructor ---------------------------------------
59 template<class V, class M>
61 {
62  //if (m_mdf ) delete m_mdf;
63  //if (m_subCdf ) delete m_subCdf;
64  //if (m_unifiedCdf) delete m_unifiedCdf;
65  //if (m_realizer ) delete m_realizer;
66  //if (m_pdf ) delete m_pdf;
67 }
68 // RV handling-methods ------------------------------
69 
70 template <class V, class M>
71 const BaseEnvironment&
73 {
74  return m_env;
75 }
76 //---------------------------------------------------
77 template<class V, class M>
78 const VectorSet<V,M>&
80 {
81  return m_imageSet;
82 }
83 //---------------------------------------------------
84 template<class V, class M>
85 const BaseJointPdf<V,M>&
87 {
88  UQ_FATAL_TEST_MACRO(m_pdf == NULL,
89  m_env.worldRank(),
90  "BaseVectorRV<V,M>::pdf()",
91  "m_pdf is NULL");
92 
93  return *m_pdf;
94 }
95 //---------------------------------------------------
96 template<class V, class M>
99 {
100  UQ_FATAL_TEST_MACRO(m_realizer == NULL,
101  m_env.worldRank(),
102  (std::string)("BaseVectorRV<V,M>::realizer(), prefix=")+m_prefix,
103  "m_realizer is NULL");
104 
105  return *m_realizer;
106 }
107 //---------------------------------------------------
108 template<class V, class M>
109 const BaseVectorCdf<V,M>&
111 {
112  UQ_FATAL_TEST_MACRO(m_subCdf == NULL,
113  m_env.worldRank(),
114  (std::string)("BaseVectorRV<V,M>::subCdf(), prefix=")+m_prefix,
115  "m_subCdf is NULL");
116 
117  return *m_subCdf;
118 }
119 //---------------------------------------------------
120 template<class V, class M>
121 const BaseVectorCdf<V,M>&
123 {
124  UQ_FATAL_TEST_MACRO(m_unifiedCdf == NULL,
125  m_env.worldRank(),
126  (std::string)("BaseVectorRV<V,M>::unifiedCdf(), prefix=")+m_prefix,
127  "m_unifiedCdf is NULL");
128 
129  return *m_unifiedCdf;
130 }
131 //---------------------------------------------------
132 template<class V, class M>
133 const BaseVectorMdf<V,M>&
135 {
136  UQ_FATAL_TEST_MACRO(m_mdf == NULL,
137  m_env.worldRank(),
138  (std::string)("BaseVectorRV<V,M>::mdf(), prefix=")+m_prefix,
139  "m_mdf is NULL");
140 
141  return *m_mdf;
142 }
143 
144 //---------------------------------------------------
145 #ifdef QUESO_HAS_ANN
146 template <class V, class M>
147 double
149 {
150  ANNpointArray data;
151  double* dists;
152  double ENT_est;
153 
154  // FIXME: these default values should be stored in the
155  // QUESO input file ( create a InfoTheoryOptions )
156  unsigned int k = UQ_INFTH_ANN_KNN;
157  double eps = UQ_INFTH_ANN_EPS;
158 
159  // here it is assumed that the entropy for the
160  // entire joint RV will be computed
161  unsigned int dim = this->imageSet().vectorSpace().dimGlobal();
162 
163  // FIXME: get the number already stored, otherwise
164  // use the default value
165  unsigned int N = this->realizer().subPeriod();
166  if( N == 0 ) {
167  N = UQ_INFTH_ANN_NO_SMP;
168  }
169 
170  // allocate memory
171  data = annAllocPts(N,dim);
172  dists = new double[N];
173 
174  // copy samples in the ANN data structure
175  V smpRV( this->imageSet().vectorSpace().zeroVector() );
176  for( unsigned int i = 0; i < N; i++ ) {
177  // get a sample from the distribution
178  this->realizer().realization( smpRV );
179  // copy the vector values in the ANN data structure
180  for( unsigned int j = 0; j < dim; j++ ) {
181  data[ i ][ j ] = smpRV[ j ];
182  }
183  }
184 
185  // get distance to knn for each point
186  // (k+1) because the 1st nn is itself
187  distANN_XY( data, data, dists, dim, dim, N, N, k+1, eps );
188 
189  // compute the entropy estimate using the L-infinity (Max) norm
190  // so no need for the adjustment of the mass of the hyperball
191  // this has to be enforced before compiling the ANN lib
192  double sum_log_dist = 0.0;
193  for( unsigned int i = 0; i < N; i++ ) {
194  if( dists[ i ] > 0 ) {
195  sum_log_dist += log( 2.0*dists[ i ] );
196  }
197  }
198  ENT_est = - gsl_sf_psi_int( k ) + gsl_sf_psi_int( N ) + (double)dim / (double)N * sum_log_dist; // todo: take specificity of gsl_, i.e., make it general (gsl or boost or etc)
199 
200  // deallocate memory
201  delete [] dists;
202  annDeallocPts( data );
203 
204  return ENT_est;
205 }
206 #endif // QUESO_HAS_ANN
207 
208 //---------------------------------------------------
209 // Method declared outside class definition ---------
210 //---------------------------------------------------
211 template <class P_V, class P_M, class Q_V, class Q_M>
212 void
214  const BaseVectorRV<P_V,P_M>& paramRv,
215  const BaseVectorRV<Q_V,Q_M>& qoiRv,
216  unsigned int localNumSamples,
217  P_M& pqCovMatrix,
218  P_M& pqCorrMatrix)
219 {
220  // Check input data consistency
221  const BaseEnvironment& env = paramRv.env();
222 
223  bool useOnlyInter0Comm = (paramRv.imageSet().vectorSpace().numOfProcsForStorage() == 1) &&
224  (qoiRv.imageSet().vectorSpace().numOfProcsForStorage() == 1);
225 
226  UQ_FATAL_TEST_MACRO((useOnlyInter0Comm == false),
227  env.worldRank(),
228  "ComputeCovCorrMatricesBetweenVectorRvs()",
229  "parallel vectors not supported yet");
230 
231  unsigned int numRows = paramRv.imageSet().vectorSpace().dim();
232  unsigned int numCols = qoiRv.imageSet().vectorSpace().dim();
233 
234  UQ_FATAL_TEST_MACRO((numRows != pqCovMatrix.numRows()) || (numCols != pqCovMatrix.numCols()),
235  env.worldRank(),
236  "ComputeCovCorrMatricesBetweenVectorRvs()",
237  "inconsistent dimensions for covariance matrix");
238 
239  UQ_FATAL_TEST_MACRO((numRows != pqCorrMatrix.numRows()) || (numCols != pqCorrMatrix.numCols()),
240  env.worldRank(),
241  "ComputeCorrelationBetweenVectorRvs()",
242  "inconsistent dimensions for correlation matrix");
243 
244  UQ_FATAL_TEST_MACRO((localNumSamples > paramRv.realizer().period()) || (localNumSamples > qoiRv.realizer().period()),
245  env.worldRank(),
246  "ComputeCovCorrMatricesBetweenVectorRvs()",
247  "localNumSamples is too large");
248 
249  // For both P and Q vector sequences: fill them
250  P_V tmpP(paramRv.imageSet().vectorSpace().zeroVector());
251  Q_V tmpQ(qoiRv.imageSet().vectorSpace().zeroVector());
252 
253  SequenceOfVectors<P_V,P_M> localWorkingPSeq(paramRv.imageSet().vectorSpace(),
254  localNumSamples,
255  "covTmpP");
256  SequenceOfVectors<Q_V,Q_M> localWorkingQSeq(qoiRv.imageSet().vectorSpace(),
257  localNumSamples,
258  "covTmpQ");
259  for (unsigned int k = 0; k < localNumSamples; ++k) {
260  paramRv.realizer().realization(tmpP);
261  localWorkingPSeq.setPositionValues(k,tmpP);
262 
263  qoiRv.realizer().realization(tmpQ);
264  localWorkingQSeq.setPositionValues(k,tmpQ);
265  }
266 
268  localWorkingQSeq,
269  localNumSamples,
270  pqCovMatrix,
271  pqCorrMatrix);
272 
273  return;
274 }
275 
276 } // End namespace QUESO
277 
A templated class for handling sets.
Definition: VectorSet.h:49
virtual ~BaseVectorRV()
Virtual destructor.
Definition: VectorRV.C:60
A templated (base) class for handling CDFs of vector functions.
Definition: VectorCdf.h:56
A templated (base) class for handling sampling from vector RVs.
const BaseVectorCdf< V, M > & unifiedCdf() const
Finds the Cumulative Distribution Function of this vector RV, considering the unified sequence of dat...
Definition: VectorRV.C:122
const BaseJointPdf< V, M > & pdf() const
Posterior Density Function of the vector RV; access to private attribute m_pdf.
Definition: VectorRV.C:86
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:305
const BaseVectorMdf< V, M > & mdf() const
Finds the Mass Density Function of this vector RV; access to private attribute m_mdf.
Definition: VectorRV.C:134
const BaseEnvironment & m_env
Definition: VectorRV.h:112
This (virtual) class sets up the environment underlying the use of the QUESO library by an executable...
Definition: Environment.h:187
void ComputeCovCorrMatricesBetweenVectorRvs(const BaseVectorRV< P_V, P_M > &paramRv, const BaseVectorRV< Q_V, Q_M > &qoiRv, unsigned int localNumSamples, P_M &pqCovMatrix, P_M &pqCorrMatrix)
Definition: VectorRV.C:213
void ComputeCovCorrMatricesBetweenVectorSequences(const BaseVectorSequence< P_V, P_M > &subPSeq, const BaseVectorSequence< Q_V, Q_M > &subQSeq, unsigned int subNumSamples, P_M &pqCovMatrix, P_M &pqCorrMatrix)
const BaseEnvironment & env() const
QUESO environment; access to private attribute m_env.
Definition: VectorRV.C:72
A templated (base) class for handling MDFs of vector functions.
Definition: VectorMdf.h:56
BaseVectorRV(const char *prefix, const VectorSet< V, M > &imageSet)
Constructor.
Definition: VectorRV.C:33
const VectorSet< V, M > & imageSet() const
Image set of the vector RV; access to private attribute m_imageSet.
Definition: VectorRV.C:79
const BaseVectorRealizer< V, M > & realizer() const
Finds a realization (sample) of the PDF of this vector RV; access to private attribute m_realizer...
Definition: VectorRV.C:98
Class for handling vector samples (sequence of vectors).
const V & zeroVector() const
Returns a vector filled with zeros.
Definition: VectorSpace.C:218
A templated base class for handling vector RV.
Definition: VectorRV.h:51
unsigned int displayVerbosity() const
Definition: Environment.C:436
virtual const VectorSpace< V, M > & vectorSpace() const =0
Vector space to which this set belongs to. See template specialization.
unsigned int numOfProcsForStorage() const
Returns total number of processes.
Definition: VectorSpace.C:193
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
std::string m_prefix
Definition: VectorRV.h:113
A templated (base) class for handling joint PDFs.
Definition: JointPdf.h:60
const BaseVectorCdf< V, M > & subCdf() const
Finds the Cumulative Distribution Function of this vector RV, considering only the sub-sequence of da...
Definition: VectorRV.C:110
virtual void realization(V &nextValues) const =0
Performs a realization (sample) from a probability density function. See template specialization...

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