queso-0.53.0
Private Member Functions | Private Attributes | List of all members
QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M > Class Template Reference

This templated class represents a Statistical Forward Problem. More...

#include <StatisticalForwardProblem.h>

Collaboration diagram for QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >:
Collaboration graph
[legend]

Public Member Functions

Constructor/Destructor methods
 StatisticalForwardProblem (const char *prefix, const SfpOptionsValues *alternativeOptionsValues, const BaseVectorRV< P_V, P_M > &paramRv, const BaseVectorFunction< P_V, P_M, Q_V, Q_M > &qoiFunction, GenericVectorRV< Q_V, Q_M > &qoiRv)
 Constructor. More...
 
 ~StatisticalForwardProblem ()
 Destructor. More...
 
Statistical methods
bool computeSolutionFlag () const
 Whether or not compute the solution. More...
 
void solveWithMonteCarlo (const McOptionsValues *alternativeOptionsValues)
 Solves the problem through Monte Carlo algorithm. More...
 
const GenericVectorRV< Q_V, Q_M > & qoiRv () const
 Returns the QoI RV; access to private attribute m_qoiRv. More...
 
const BaseVectorSequence< Q_V,
Q_M > & 
getParamChain () const
 Returns the parameter chain; access to private attribute m_paramChain. More...
 

Private Member Functions

void commonConstructor ()
 TODO: Common constructor. More...
 

Private Attributes

const BaseEnvironmentm_env
 
const BaseVectorRV< P_V, P_M > & m_paramRv
 
const BaseVectorFunction< P_V,
P_M, Q_V, Q_M > & 
m_qoiFunction
 
GenericVectorRV< Q_V, Q_M > & m_qoiRv
 
BaseVectorSequence< Q_V, Q_M > * m_paramChain
 
BaseVectorSequence< Q_V, Q_M > * m_qoiChain
 
MonteCarloSG< P_V, P_M, Q_V,
Q_M > * 
m_mcSeqGenerator
 
BaseVectorRealizer< Q_V, Q_M > * m_solutionRealizer
 
BaseJointPdf< Q_V, Q_M > * m_solutionPdf
 
const SfpOptionsValuesm_optionsObj
 
bool m_userDidNotProvideOptions
 

I/O methods

void print (std::ostream &os) const
 TODO: Prints the sequence. More...
 
std::ostream & operator<< (std::ostream &os, const StatisticalForwardProblem< P_V, P_M, Q_V, Q_M > &obj)
 

Detailed Description

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
class QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >

This templated class represents a Statistical Forward Problem.

This templated class represents a statistical forward problem. It is templated on the types 'P_V' and 'Q_V' of vectors and types 'P_M' and 'Q_M' of matrices, where 'P_' stands for 'parameter' and 'Q_' stands for 'quantities of interest'. Conceptually, a statistical forward problem has two input entities and one output entity.
The input entities of a statistical forward problem are:

  1. the input (parameter) RV, an instance of class 'BaseVectorRV<P_V,P_M>', and
  2. the QoI function, an instance of class 'BaseVectorFunction<P_V,P_M,Q_V,Q_M>'.

Let $ q() $ denote the mathematical QoI function and $ x $ denote a vector of parameters. The QoI function object stores the routine that computes $ q(x) $ and whatever data necessary by such routine. See file 'libs/basic/inc/VectorFunction.h' for more details.
The output entity of a statistical forward problem is:

  1. the QoI RV, another instance of class 'BaseVectorRV<P_V,P_M>'.

The QoI RV stores the solution according to the Bayesian approach. The solution of a SPF is computed by calling 'solveWithMonteCarlo()'. More operations, with different methods, will be available in the future.
The solution process might demand extra objects to be passed through the chosen solution operation interface. This distinction is important: this class separates 'what the problem is' from 'how the problem is solved'. Upon return from a solution operation, the QoI RV is available through the operation 'qoiRv()'. Such QoI RV is able to provide:

  1. a vector realizer through the operation 'qoiRv().realizer()', which returns an instance of the class 'BaseVectorRealizer<Q_V,Q_M>'.

Definition at line 80 of file StatisticalForwardProblem.h.

Constructor & Destructor Documentation

template<class P_V , class P_M , class Q_V , class Q_M >
QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::StatisticalForwardProblem ( const char *  prefix,
const SfpOptionsValues alternativeOptionsValues,
const BaseVectorRV< P_V, P_M > &  paramRv,
const BaseVectorFunction< P_V, P_M, Q_V, Q_M > &  qoiFunction,
GenericVectorRV< Q_V, Q_M > &  qoiRv 
)

Constructor.

Requirements: 1) the image set of the vector random variable 'paramRv' and the domain set of the QoI function 'qoiFunction' should belong to vector spaces of equal dimensions and 2) the image set of the QoI function 'qoiFunction' and the image set of the vector random variable 'qoiRv' should belong to vector spaces of equal dimensions. If the requirements are satisfied, the constructor then reads input options that begin with the string '<prefix>fp_'. Options reading is handled by class 'StatisticalForwardProblemOptions'. If no options input file is provided, the construction assigns alternativeOptionsValues to the options of the SFP.

Parameters
prefixThe prefix
alternativeOptionsValuesOptions (if no input file)
paramRvThe input RV
qoiFunctionThe QoI function
qoiRvThe QoI RV

Definition at line 34 of file StatisticalForwardProblem.C.

References queso_require_equal_to_msg.

40  :
41  m_env (paramRv.env()),
42  m_paramRv (paramRv),
43  m_qoiFunction (qoiFunction),
44  m_qoiRv (qoiRv),
45  m_paramChain (NULL),
46  m_qoiChain (NULL),
47  m_mcSeqGenerator (NULL),
48  m_solutionRealizer (NULL),
49 #ifdef UQ_ALSO_COMPUTE_MDFS_WITHOUT_KDE
50  m_subMdfGrids (NULL),
51  m_subMdfValues (NULL),
52 #endif
53 #ifdef QUESO_COMPUTES_EXTRA_POST_PROCESSING_STATISTICS
54  m_subSolutionMdf (NULL),
55  m_subCdfGrids (NULL),
56  m_subCdfValues (NULL),
57  m_subSolutionCdf (NULL),
58  m_unifiedCdfGrids (NULL),
59  m_unifiedCdfValues (NULL),
60  m_unifiedSolutionCdf (NULL),
61 #endif
62  m_solutionPdf (NULL),
63  m_optionsObj (alternativeOptionsValues),
65 {
66  if (m_env.subDisplayFile()) {
67  *m_env.subDisplayFile() << "Entering StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::constructor()"
68  << ": prefix = " << prefix
69  << ", alternativeOptionsValues = " << alternativeOptionsValues
70  << ", m_env.optionsInputFileName() = " << m_env.optionsInputFileName()
71  << std::endl;
72  }
73 
74  // If NULL, we create one
75  if (m_optionsObj == NULL) {
76  SfpOptionsValues * tempOptions = new SfpOptionsValues(&m_env, prefix);
77 
78  // We did this dance because scanOptionsValues is not a const method, but
79  // m_optionsObj is a pointer to const
80  m_optionsObj = tempOptions;
81 
82  // We do this so we don't delete the user's object in the dtor
84  }
85 
86  if (m_optionsObj->m_help != "") {
87  if (m_env.subDisplayFile()) {
88  *m_env.subDisplayFile() << (*m_optionsObj) << std::endl;
89  }
90  }
91 
92  queso_require_equal_to_msg(paramRv.imageSet().vectorSpace().dimLocal(), qoiFunction.domainSet().vectorSpace().dimLocal(), "'paramRv' and 'qoiFunction' are related to vector spaces of different dimensions");
93 
94  queso_require_equal_to_msg(qoiFunction.imageSet().vectorSpace().dimLocal(), qoiRv.imageSet().vectorSpace().dimLocal(), "'qoiFunction' and 'qoiRv' are related to vector spaces of different dimensions");
95 
96  if (m_env.subDisplayFile()) {
97  *m_env.subDisplayFile() << "Leaving StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::constructor()"
98  << ": prefix = " << m_optionsObj->m_prefix
99  << std::endl;
100  }
101 }
BaseVectorRealizer< Q_V, Q_M > * m_solutionRealizer
unsigned int dimLocal() const
Definition: VectorSpace.C:170
const BaseEnvironment & env() const
QUESO environment; access to private attribute m_env.
Definition: VectorRV.C:72
BaseJointPdf< Q_V, Q_M > * m_solutionPdf
const BaseVectorFunction< P_V, P_M, Q_V, Q_M > & m_qoiFunction
BaseVectorSequence< Q_V, Q_M > * m_paramChain
#define queso_require_equal_to_msg(expr1, expr2, msg)
Definition: asserts.h:85
const BaseVectorRV< P_V, P_M > & m_paramRv
std::string optionsInputFileName() const
Access to the attribute m_optionsInputFileName, which stores the name of the input file passed by the...
Definition: Environment.C:307
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:274
BaseVectorSequence< Q_V, Q_M > * m_qoiChain
std::string m_help
If non-empty string, options and values are printed to the output file.
const VectorSet< V, M > & imageSet() const
Image set of the vector RV; access to private attribute m_imageSet.
Definition: VectorRV.C:79
GenericVectorRV< Q_V, Q_M > & m_qoiRv
virtual const VectorSpace< V, M > & vectorSpace() const =0
Vector space to which this set belongs to. See template specialization.
MonteCarloSG< P_V, P_M, Q_V, Q_M > * m_mcSeqGenerator
template<class P_V , class P_M , class Q_V , class Q_M >
QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::~StatisticalForwardProblem ( )

Destructor.

Definition at line 105 of file StatisticalForwardProblem.C.

106 {
107  if (m_solutionPdf ) delete m_solutionPdf;
108 
109 #ifdef QUESO_COMPUTES_EXTRA_POST_PROCESSING_STATISTICS
110  if (m_unifiedSolutionCdf) delete m_unifiedSolutionCdf;
111  if (m_unifiedCdfValues ) delete m_unifiedCdfValues;
112  if (m_unifiedCdfGrids ) delete m_unifiedCdfGrids;
113 
114  if (m_subSolutionCdf ) delete m_subSolutionCdf;
115  if (m_subCdfValues ) delete m_subCdfValues;
116  if (m_subCdfGrids ) delete m_subCdfGrids;
117 
118  if (m_subSolutionMdf ) delete m_subSolutionMdf;
119 #endif
120 #ifdef UQ_ALSO_COMPUTE_MDFS_WITHOUT_KDE
121  if (m_subMdfValues ) delete m_subMdfValues;
122  if (m_subMdfGrids ) delete m_subMdfGrids;
123 #endif
125 
126  if (m_mcSeqGenerator ) delete m_mcSeqGenerator;
127 
128  if (m_qoiChain) {
129  m_qoiChain->clear();
130  delete m_qoiChain;
131  }
132 
133  if (m_paramChain) {
134  m_paramChain->clear();
135  delete m_paramChain;
136  }
137 
138  if (m_optionsObj ) delete m_optionsObj;
139 }
BaseVectorRealizer< Q_V, Q_M > * m_solutionRealizer
BaseJointPdf< Q_V, Q_M > * m_solutionPdf
BaseVectorSequence< Q_V, Q_M > * m_paramChain
void clear()
Reset the values and the size of the sequence of vectors.
BaseVectorSequence< Q_V, Q_M > * m_qoiChain
MonteCarloSG< P_V, P_M, Q_V, Q_M > * m_mcSeqGenerator

Member Function Documentation

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
void QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::commonConstructor ( )
private

TODO: Common constructor.

Todo:
: implement me!
template<class P_V , class P_M , class Q_V , class Q_M >
bool QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::computeSolutionFlag ( ) const

Whether or not compute the solution.

Definition at line 144 of file StatisticalForwardProblem.C.

template<class P_V , class P_M , class Q_V , class Q_M >
const BaseVectorSequence< Q_V, Q_M > & QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::getParamChain ( ) const

Returns the parameter chain; access to private attribute m_paramChain.

Definition at line 402 of file StatisticalForwardProblem.C.

References queso_require_msg.

403 {
404 
405  // Make sure this runs after the forward propagation
406  // only then we obtain the actual realizations of the parameters
407  queso_require_msg(m_paramChain, "m_paramChain is NULL");
408 
409  return *m_paramChain;
410 
411 }
BaseVectorSequence< Q_V, Q_M > * m_paramChain
#define queso_require_msg(asserted, msg)
Definition: asserts.h:69
template<class P_V , class P_M , class Q_V , class Q_M >
void QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::print ( std::ostream &  os) const

TODO: Prints the sequence.

Todo:
: implement me!

Definition at line 415 of file StatisticalForwardProblem.C.

416 {
417  return;
418 }
template<class P_V , class P_M , class Q_V , class Q_M >
const GenericVectorRV< Q_V, Q_M > & QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::qoiRv ( ) const

Returns the QoI RV; access to private attribute m_qoiRv.

Definition at line 373 of file StatisticalForwardProblem.C.

374 {
375  return m_qoiRv;
376 }
GenericVectorRV< Q_V, Q_M > & m_qoiRv
template<class P_V , class P_M , class Q_V , class Q_M >
void QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::solveWithMonteCarlo ( const McOptionsValues alternativeOptionsValues)

Solves the problem through Monte Carlo algorithm.

  Requirements: none at this moment. If the requirements are satisfied, this operation checks
 the member flag 'm_computeSolution' (one of the options read from the input file during
 construction). If the flag is 'false', the operation returns immediately, computing nothing.
 If the flag is 'true', the operation sets the member variable 'm_qoiRv' accordingly.
 The operation:
  1. instantiates 'SequenceOfVectors<P_V,P_M>' (the input sequence of vectors),
  2. instantiates 'SequenceOfVectors<Q_V,Q_M>' (the output sequence of vectors),
  3. instantiates 'MonteCarloSG<P_V,P_M,Q_V,Q_M>' (the Monte Carlo algorithm),
  4. populates the output sequence with the Monte Carlo algorithm,
  5. sets the realizer of 'm_qoiRv' with the contents of the output sequence, and

Definition at line 151 of file StatisticalForwardProblem.C.

References QUESO::ComputeCovCorrMatricesBetweenVectorSequences(), QUESO::BaseVectorCdf< V, M >::m_prefix, QUESO::FilePtrSetStruct::ofsVar, queso_error_msg, QUESO::SequenceOfVectors< V, M >::subUniformlySampledCdf(), QUESO::SequenceOfVectors< V, M >::unifiedUniformlySampledCdf(), and UQ_FILE_EXTENSION_FOR_MATLAB_FORMAT.

153 {
154  m_env.fullComm().Barrier();
155  m_env.fullComm().syncPrintDebugMsg("Entering StatisticalForwardProblem<P_V,P_M>::solveWithMonteCarlo()",1,3000000);
156 
157  if (m_optionsObj->m_computeSolution == false) {
158  if ((m_env.subDisplayFile())) {
159  *m_env.subDisplayFile() << "In StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::solveWithMonteCarlo()"
160  << ": avoiding solution, as requested by user"
161  << std::endl;
162  }
163  return;
164  }
165  if ((m_env.subDisplayFile())) {
166  *m_env.subDisplayFile() << "In StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::solveWithMonteCarlo()"
167  << ": computing solution, as requested by user"
168  << std::endl;
169  }
170 
171  if (m_solutionPdf ) delete m_solutionPdf;
172 
173 #ifdef QUESO_COMPUTES_EXTRA_POST_PROCESSING_STATISTICS
174  if (m_unifiedSolutionCdf) delete m_unifiedSolutionCdf;
175  if (m_unifiedCdfValues ) delete m_unifiedCdfValues;
176  if (m_unifiedCdfGrids ) delete m_unifiedCdfGrids;
177 
178  if (m_subSolutionCdf ) delete m_subSolutionCdf;
179  if (m_subCdfValues ) delete m_subCdfValues;
180  if (m_subCdfGrids ) delete m_subCdfGrids;
181 
182  if (m_subSolutionMdf ) delete m_subSolutionMdf;
183 #endif
184 #ifdef UQ_ALSO_COMPUTE_MDFS_WITHOUT_KDE
185  if (m_subMdfValues ) delete m_subMdfValues;
186  if (m_subMdfGrids ) delete m_subMdfGrids;
187 #endif
189 
190  if (m_mcSeqGenerator ) delete m_mcSeqGenerator;
191 
192  if (m_qoiChain) {
193  m_qoiChain->clear();
194  delete m_qoiChain;
195  }
196 
197  if (m_paramChain) {
198  m_paramChain->clear();
199  delete m_paramChain;
200  }
201 
202  Q_V numEvaluationPointsVec(m_qoiRv.imageSet().vectorSpace().zeroVector());
203  numEvaluationPointsVec.cwSet(250.);
204 
205  // Compute output realizer: Monte Carlo approach
206  m_paramChain = new SequenceOfVectors<P_V,P_M>(m_paramRv.imageSet().vectorSpace(),0,m_optionsObj->m_prefix+"paramChain");
207  m_qoiChain = new SequenceOfVectors<Q_V,Q_M>(m_qoiRv.imageSet().vectorSpace(), 0,m_optionsObj->m_prefix+"qoiChain" );
208  m_mcSeqGenerator = new MonteCarloSG<P_V,P_M,Q_V,Q_M>(m_optionsObj->m_prefix.c_str(),
209  alternativeOptionsValues,
210  m_paramRv,
211  m_qoiFunction);
212  //m_qoiRv);
213  m_mcSeqGenerator->generateSequence(*m_paramChain,
214  *m_qoiChain);
215  m_solutionRealizer = new SequentialVectorRealizer<Q_V,Q_M>((m_optionsObj->m_prefix+"Qoi").c_str(),
216  *m_qoiChain);
218 
219  // Compute output mdf: uniform sampling approach
220 #ifdef UQ_ALSO_COMPUTE_MDFS_WITHOUT_KDE
221  m_subMdfGrids = new ArrayOfOneDGrids <Q_V,Q_M>((m_optionsObj->m_prefix+"QoiMdf_").c_str(),m_qoiRv.imageSet().vectorSpace());
222  m_subMdfValues = new ArrayOfOneDTables<Q_V,Q_M>((m_optionsObj->m_prefix+"QoiMdf_").c_str(),m_qoiRv.imageSet().vectorSpace());
223  m_qoiChain->subUniformlySampledMdf(numEvaluationPointsVec, // input
224  *m_subMdfGrids, // output
225  *m_subMdfValues); // output
226 
227  m_subSolutionMdf = new SampledVectorMdf<Q_V,Q_M>((m_optionsObj->m_prefix+"Qoi").c_str(),
228  *m_subMdfGrids,
229  *m_subMdfValues);
230  m_qoiRv.setMdf(*m_subSolutionMdf);
231 #endif
232 
233  // Compute output cdf: uniform sampling approach
234 #ifdef QUESO_COMPUTES_EXTRA_POST_PROCESSING_STATISTICS
235  std::string subCoreName_qoiCdf(m_optionsObj->m_prefix+ "QoiCdf_");
236  std::string uniCoreName_qoiCdf(m_optionsObj->m_prefix+"unifQoiCdf_");
237  if (m_env.numSubEnvironments() == 1) subCoreName_qoiCdf = uniCoreName_qoiCdf;
238 
239  std::string subCoreName_solutionCdf(m_optionsObj->m_prefix+ "Qoi");
240  std::string uniCoreName_solutionCdf(m_optionsObj->m_prefix+"unifQoi");
241  if (m_env.numSubEnvironments() == 1) subCoreName_solutionCdf = uniCoreName_solutionCdf;
242 
243  m_subCdfGrids = new ArrayOfOneDGrids <Q_V,Q_M>(subCoreName_qoiCdf.c_str(),m_qoiRv.imageSet().vectorSpace());
244  m_subCdfValues = new ArrayOfOneDTables<Q_V,Q_M>(subCoreName_qoiCdf.c_str(),m_qoiRv.imageSet().vectorSpace());
245  m_qoiChain->subUniformlySampledCdf(numEvaluationPointsVec, // input
246  *m_subCdfGrids, // output
247  *m_subCdfValues); // output
248 
249  m_subSolutionCdf = new SampledVectorCdf<Q_V,Q_M>(subCoreName_solutionCdf.c_str(),
250  *m_subCdfGrids,
251  *m_subCdfValues);
252  m_qoiRv.setSubCdf(*m_subSolutionCdf);
253 
254  // Compute unified cdf if necessary
255  if (m_env.numSubEnvironments() == 1) {
256  m_qoiRv.setUnifiedCdf(*m_subSolutionCdf);
257  }
258  else {
259  m_unifiedCdfGrids = new ArrayOfOneDGrids <Q_V,Q_M>(uniCoreName_qoiCdf.c_str(),m_qoiRv.imageSet().vectorSpace());
260  m_unifiedCdfValues = new ArrayOfOneDTables<Q_V,Q_M>(uniCoreName_qoiCdf.c_str(),m_qoiRv.imageSet().vectorSpace());
261  m_qoiChain->unifiedUniformlySampledCdf(numEvaluationPointsVec, // input
262  *m_unifiedCdfGrids, // output
263  *m_unifiedCdfValues); // output
264 
265  m_unifiedSolutionCdf = new SampledVectorCdf<Q_V,Q_M>(uniCoreName_solutionCdf.c_str(),
266  *m_unifiedCdfGrids,
267  *m_unifiedCdfValues);
268  m_qoiRv.setUnifiedCdf(*m_unifiedSolutionCdf);
269  }
270 #endif
271  // Compute (just unified one) covariance matrix, if requested
272  // Compute (just unified one) correlation matrix, if requested
273  P_M* pqCovarianceMatrix = NULL;
274  P_M* pqCorrelationMatrix = NULL;
276  if (m_env.subDisplayFile()) {
277  *m_env.subDisplayFile() << "In StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::solveWithMonteCarlo()"
278  << ", prefix = " << m_optionsObj->m_prefix
279  << ": instantiating cov and corr matrices"
280  << std::endl;
281  }
282  pqCovarianceMatrix = new P_M(m_env,
283  m_paramRv.imageSet().vectorSpace().map(), // number of rows
284  m_qoiRv.imageSet().vectorSpace().dimGlobal()); // number of cols
285  pqCorrelationMatrix = new P_M(m_env,
286  m_paramRv.imageSet().vectorSpace().map(), // number of rows
287  m_qoiRv.imageSet().vectorSpace().dimGlobal()); // number of cols
289  *m_qoiChain,
290  std::min(m_paramRv.realizer().subPeriod(),m_qoiRv.realizer().subPeriod()), // FIX ME: might be INFINITY
291  *pqCovarianceMatrix,
292  *pqCorrelationMatrix);
293  }
294 
295  // Write data out
296  if (m_env.subDisplayFile()) {
297  if (pqCovarianceMatrix ) {
298  *m_env.subDisplayFile() << "In StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::solveWithMonteCarlo()"
299  << ", prefix = " << m_optionsObj->m_prefix
300  << ": contents of covariance matrix are\n" << *pqCovarianceMatrix
301  << std::endl;
302  }
303  if (pqCorrelationMatrix) {
304  *m_env.subDisplayFile() << "In StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::solveWithMonteCarlo()"
305  << ", prefix = " << m_optionsObj->m_prefix
306  << ": contents of correlation matrix are\n" << *pqCorrelationMatrix
307  << std::endl;
308  }
309  }
310 
311  // Open data output file
312  if (m_env.subDisplayFile()) {
313  *m_env.subDisplayFile() << "In StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::solveWithMonteCarlo()"
314  << ", prefix = " << m_optionsObj->m_prefix
315  << ": checking necessity of opening data output file '" << m_optionsObj->m_dataOutputFileName
316  << "'"
317  << std::endl;
318  }
319  FilePtrSetStruct filePtrSet;
321  UQ_FILE_EXTENSION_FOR_MATLAB_FORMAT, // Yes, always ".m"
323  false,
324  filePtrSet)) {
325 #ifdef UQ_ALSO_COMPUTE_MDFS_WITHOUT_KDE
326  m_qoiRv.mdf().print(*filePtrSet.ofsVar);
327 #endif
328 #ifdef QUESO_COMPUTES_EXTRA_POST_PROCESSING_STATISTICS
329  *filePtrSet.ofsVar << m_qoiRv.subCdf();
330 #endif
331 
332  //if (pqCovarianceMatrix ) *filePtrSet.ofsVar << *pqCovarianceMatrix; // FIX ME: output matrix in matlab format
333  //if (pqCorrelationMatrix) *filePtrSet.ofsVar << *pqCorrelationMatrix; // FIX ME: output matrix in matlab format
334 
335  // Write unified cdf if necessary
336 #ifdef QUESO_COMPUTES_EXTRA_POST_PROCESSING_STATISTICS
337  if (m_env.numSubEnvironments() > 1) {
338  if (m_qoiRv.imageSet().vectorSpace().numOfProcsForStorage() == 1) {
339  if (m_env.inter0Rank() == 0) {
340  *filePtrSet.ofsVar << m_qoiRv.unifiedCdf(); //*m_unifiedSolutionCdf;
341  }
342  }
343  else {
344  queso_error_msg("unified cdf writing, parallel vectors not supported yet");
345  }
346  }
347 #endif
348  // Close data output file
350  if (m_env.subDisplayFile()) {
351  *m_env.subDisplayFile() << "In StatisticalForwardProblem<P_V,P_M,Q_V,Q_M>::solveWithMonteCarlo()"
352  << ", prefix = " << m_optionsObj->m_prefix
353  << ": closed data output file '" << m_optionsObj->m_dataOutputFileName
354  << "'"
355  << std::endl;
356  }
357  }
358  if (m_env.subDisplayFile()) {
359  *m_env.subDisplayFile() << std::endl;
360  }
361 
362  if (pqCovarianceMatrix ) delete pqCovarianceMatrix;
363  if (pqCorrelationMatrix) delete pqCorrelationMatrix;
364 
365  m_env.fullComm().syncPrintDebugMsg("Leaving StatisticalForwardProblem<P_V,P_M>::solveWithMonteCarlo()",1,3000000);
366  m_env.fullComm().Barrier();
367 
368  return;
369 }
BaseVectorRealizer< Q_V, Q_M > * m_solutionRealizer
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:95
void Barrier() const
Pause every process in *this communicator until all the processes reach this point.
Definition: MpiComm.C:122
BaseJointPdf< Q_V, Q_M > * m_solutionPdf
const BaseVectorFunction< P_V, P_M, Q_V, Q_M > & m_qoiFunction
BaseVectorSequence< Q_V, Q_M > * m_paramChain
#define queso_error_msg(msg)
Definition: asserts.h:47
const MpiComm & fullComm() const
Access function for MpiComm full communicator.
Definition: Environment.C:228
std::set< unsigned int > m_dataOutputAllowedSet
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)
unsigned int numSubEnvironments() const
Access function to the number of sub-environments.
Definition: Environment.C:288
const BaseVectorCdf< V, M > & unifiedCdf() const
Finds the Cumulative Distribution Function of this vector RV, considering the unified sequence of dat...
Definition: VectorRV.C:113
int inter0Rank() const
Returns the process inter0 rank.
Definition: Environment.C:261
void closeFile(FilePtrSetStruct &filePtrSet, const std::string &fileType) const
Closes the file.
Definition: Environment.C:1020
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:104
void clear()
Reset the values and the size of the sequence of vectors.
const BaseVectorRV< P_V, P_M > & m_paramRv
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:274
BaseVectorSequence< Q_V, Q_M > * m_qoiChain
#define UQ_FILE_EXTENSION_FOR_MATLAB_FORMAT
Definition: Defines.h:89
void setMdf(BaseVectorMdf< V, M > &mdf)
Sets the MDF of this vector RV to Mdf.
const VectorSet< V, M > & imageSet() const
Image set of the vector RV; access to private attribute m_imageSet.
Definition: VectorRV.C:79
void setRealizer(BaseVectorRealizer< V, M > &realizer)
Sets the realizer of this vector RV to realizer.
bool openOutputFile(const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds, bool writeOver, FilePtrSetStruct &filePtrSet) const
Opens an output file for each sub-environment that was chosen to send data to the file...
Definition: Environment.C:467
GenericVectorRV< Q_V, Q_M > & m_qoiRv
void setUnifiedCdf(BaseVectorCdf< V, M > &unifiedCdf)
Sets the CDF of the unified sequence of this vector RV to unifiedCdf.
virtual const VectorSpace< V, M > & vectorSpace() const =0
Vector space to which this set belongs to. See template specialization.
void syncPrintDebugMsg(const char *msg, unsigned int msgVerbosity, unsigned int numUSecs) const
Synchronizes all the processes and print debug message.
Definition: MpiComm.C:195
const BaseVectorMdf< V, M > & mdf() const
Finds the Mass Density Function of this vector RV; access to private attribute m_mdf.
Definition: VectorRV.C:122
MonteCarloSG< P_V, P_M, Q_V, Q_M > * m_mcSeqGenerator
void setSubCdf(BaseVectorCdf< V, M > &subCdf)
Sets the CDF of the sub-sequence of this vector RV to subCdf.

Friends And Related Function Documentation

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
std::ostream& operator<< ( std::ostream &  os,
const StatisticalForwardProblem< P_V, P_M, Q_V, Q_M > &  obj 
)
friend

Definition at line 144 of file StatisticalForwardProblem.h.

145  {
146  obj.print(os);
147  return os;
148  }

Member Data Documentation

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
const BaseEnvironment& QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_env
private

Definition at line 156 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
MonteCarloSG<P_V,P_M,Q_V,Q_M>* QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_mcSeqGenerator
private

Definition at line 164 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
const SfpOptionsValues* QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_optionsObj
private

Definition at line 170 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
BaseVectorSequence<Q_V,Q_M>* QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_paramChain
private

Definition at line 162 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
const BaseVectorRV<P_V,P_M>& QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_paramRv
private

Definition at line 158 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
BaseVectorSequence<Q_V,Q_M>* QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_qoiChain
private

Definition at line 163 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
const BaseVectorFunction<P_V,P_M,Q_V,Q_M>& QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_qoiFunction
private

Definition at line 159 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
GenericVectorRV<Q_V,Q_M>& QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_qoiRv
private

Definition at line 160 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
BaseJointPdf<Q_V,Q_M>* QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_solutionPdf
private

Definition at line 168 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
BaseVectorRealizer<Q_V,Q_M>* QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_solutionRealizer
private

Definition at line 166 of file StatisticalForwardProblem.h.

template<class P_V = GslVector, class P_M = GslMatrix, class Q_V = GslVector, class Q_M = GslMatrix>
bool QUESO::StatisticalForwardProblem< P_V, P_M, Q_V, Q_M >::m_userDidNotProvideOptions
private

Definition at line 187 of file StatisticalForwardProblem.h.


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

Generated on Thu Jun 11 2015 13:52:36 for queso-0.53.0 by  doxygen 1.8.5