queso-0.57.0
GslVector.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-2017 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/Defines.h>
26 #include <queso/GslVector.h>
27 #include <queso/RngBase.h>
28 #include <algorithm>
29 #include <gsl/gsl_sort_vector.h>
30 #include <cmath>
31 
32 namespace QUESO {
33 
34 GslVector::GslVector(const BaseEnvironment& env, const Map& map)
35  :
36  Vector(env,map),
37  m_vec (gsl_vector_calloc(map.NumGlobalElements()))
38 {
39  //std::cout << "Entering GslVector::constructor(1)" << std::endl;
40 
41  queso_require_msg(m_vec, "null vector generated");
42 
43  queso_require_equal_to_msg(m_vec->size, (unsigned int) map.NumMyElements(), "incompatible local vec size");
44 
45  queso_require_equal_to_msg(m_vec->size, (unsigned int) map.NumGlobalElements(), "incompatible global vec size");
46 
47  queso_require_equal_to_msg(m_vec->size, (unsigned int) m_map.NumMyElements(), "incompatible own vec size");
48 
49  //std::cout << "In GslVector::constructor(env,map)"
50  // << "\n m_vec->size = " << m_vec->size
51  // << "\n map.NumGlobalElements() = " << map.NumGlobalElements()
52  // << "\n map.NumMyElements() = " << map.NumMyElements()
53  // << std::endl;
54 
55 }
56 
57 GslVector::GslVector(const BaseEnvironment& env, const Map& map, double value)
58  :
59  Vector(env,map),
60  m_vec (gsl_vector_calloc(map.NumGlobalElements()))
61 {
62  //std::cout << "Entering GslVector::constructor(2)" << std::endl;
63 
64  queso_require_msg(m_vec, "null vector generated");
65 
66  queso_require_equal_to_msg(m_vec->size, (unsigned int) map.NumMyElements(), "incompatible local vec size");
67 
68  queso_require_equal_to_msg(m_vec->size, (unsigned int) map.NumGlobalElements(), "incompatible global vec size");
69 
70  this->cwSet(value);
71 
72  queso_require_equal_to_msg(m_vec->size, (unsigned int) m_map.NumMyElements(), "incompatible own vec size");
73 
74  //std::cout << "Leaving GslVector::constructor(2)" << std::endl;
75 }
76 
77 GslVector::GslVector(const BaseEnvironment& env, double d1, double d2, const Map& map)
78  :
79  Vector(env,map),
80  m_vec (gsl_vector_calloc(map.NumGlobalElements()))
81 {
82  //std::cout << "Entering GslVector::constructor(3)" << std::endl;
83 
84  queso_require_msg(m_vec, "null vector generated");
85 
86  queso_require_equal_to_msg(m_vec->size, (unsigned int) map.NumMyElements(), "incompatible local vec size");
87 
88  queso_require_equal_to_msg(m_vec->size, (unsigned int) map.NumGlobalElements(), "incompatible global vec size");
89 
90  for (unsigned int i = 0; i < m_vec->size; ++i) {
91  double alpha = (double) i / ((double) m_vec->size - 1.);
92  (*this)[i] = (1.-alpha)*d1 + alpha*d2;
93  }
94 
95  queso_require_equal_to_msg(m_vec->size, (unsigned int) m_map.NumMyElements(), "incompatible own vec size");
96 
97  //std::cout << "Leaving GslVector::constructor(3)" << std::endl;
98 }
99 
100 GslVector::GslVector(const GslVector& v, double start, double end)
101  :
102  Vector(v.env(),v.map()),
103  m_vec (gsl_vector_calloc(v.sizeLocal()))
104 {
105  //std::cout << "Entering GslVector::constructor(4)" << std::endl;
106 
107  queso_require_msg(m_vec, "null vector generated");
108 
109  queso_require_equal_to_msg(m_vec->size, (unsigned int) v.map().NumMyElements(), "incompatible local vec size");
110 
111  queso_require_equal_to_msg(m_vec->size, (unsigned int) v.map().NumGlobalElements(), "incompatible global vec size");
112 
113  for (unsigned int i = 0; i < m_vec->size; ++i) {
114  double alpha = (double) i / ((double) m_vec->size - 1.);
115  (*this)[i] = (1. - alpha) * start + alpha * end;
116  }
117 
118  queso_require_equal_to_msg(m_vec->size, (unsigned int) m_map.NumMyElements(), "incompatible own vec size");
119 
120  //std::cout << "Leaving GslVector::constructor(4)" << std::endl;
121 }
122 
124  :
125  Vector(v.env(),v.map()),
126  m_vec (gsl_vector_calloc(v.sizeLocal()))
127 {
128  //std::cout << "Entering GslVector::constructor(5)" << std::endl;
129 
130  // prudenci 2010-06-17 mox
131  queso_require_msg(m_vec, "null vector generated");
132 
133  queso_require_equal_to_msg(m_vec->size, (unsigned int) v.map().NumMyElements(), "incompatible local vec size");
134 
135  queso_require_equal_to_msg(m_vec->size, (unsigned int) v.map().NumGlobalElements(), "incompatible global vec size");
136 
137  this->copy(v);
138 
139  queso_require_equal_to_msg(m_vec->size, (unsigned int) m_map.NumMyElements(), "incompatible own vec size");
140 
141  //std::cout << "Leaving GslVector::constructor(5)" << std::endl;
142 }
143 
145 {
146  if (m_vec) gsl_vector_free(m_vec);
147 }
148 
149 GslVector&
151 {
152  //std::cout << "In GslVector::operator=()" // mox
153  // << ": setting size1"
154  // << std::endl;
155  unsigned int size1 = this->sizeLocal();
156  //std::cout << "In GslVector::operator=()" // mox
157  // << ": setting size2"
158  // << std::endl;
159  unsigned int size2 = rhs.sizeLocal();
160  queso_require_equal_to_msg(size1, size2, "sizes are not compatible");
161  this->copy(rhs);
162  return *this;
163 }
164 
165 GslVector&
167 {
168  int iRC;
169  iRC = gsl_vector_scale(m_vec,a);
170  queso_require_msg(!(iRC), "failed");
171  return *this;
172 }
173 
174 GslVector&
176 {
177  *this *= (1./a);
178 
179  return *this;
180 }
181 
182 GslVector&
184 {
185  unsigned int size1 = this->sizeLocal();
186  unsigned int size2 = rhs.sizeLocal();
187  queso_require_equal_to_msg(size1, size2, "different sizes of this and rhs");
188 
189  for (unsigned int i = 0; i < size1; ++i) {
190  (*this)[i] *= rhs[i];
191  }
192 
193  return *this;
194 }
195 
196 GslVector&
198 {
199  unsigned int size1 = this->sizeLocal();
200  unsigned int size2 = rhs.sizeLocal();
201  queso_require_equal_to_msg(size1, size2, "different sizes of this and rhs");
202 
203  for (unsigned int i = 0; i < size1; ++i) {
204  (*this)[i] /= rhs[i];
205  }
206 
207  return *this;
208 }
209 
210 GslVector&
212 {
213  int iRC;
214  iRC = gsl_vector_add(m_vec,rhs.m_vec);
215  queso_require_msg(!(iRC), "failed");
216  return *this;
217 }
218 
219 GslVector&
221 {
222  int iRC;
223  iRC = gsl_vector_sub(m_vec,rhs.m_vec);
224  queso_require_msg(!(iRC), "failed");
225 
226  return *this;
227 }
228 
229 void
231 {
232  this->Vector::base_copy(src);
233  int iRC;
234  iRC = gsl_vector_memcpy(this->m_vec, src.m_vec);
235  queso_require_msg(!(iRC), "failed");
236 
237  return;
238 }
239 
240 unsigned int
242 {
243  // mox
244  //std::cout << "Entering GslVector::sizeLocal()"
245  // << ": &m_map = " << &m_map
246  // << std::endl;
247  //std::cout << ", m_map.NumMyElements() = " << m_map.NumMyElements()
248  // << std::endl;
249  //std::cout << ", m_vec = " << m_vec
250  // << std::endl;
251  //std::cout << ", m_vec->size = " << m_vec->size
252  // << std::endl;
253 
254  queso_require_equal_to_msg(m_vec->size, (unsigned int) m_map.NumMyElements(), "incompatible vec size");
255 
256  //std::cout << "Leaving GslVector::sizeLocal()"
257  // << ": m_vec = " << m_vec
258  // << ", m_vec->size = " << m_vec->size
259  // << ", &m_map = " << &m_map
260  // << ", m_map.NumMyElements() = " << m_map.NumMyElements()
261  // << std::endl;
262 
263  return m_vec->size;
264 }
265 
266 unsigned int
268 {
269  queso_require_equal_to_msg(m_vec->size, (unsigned int) m_map.NumGlobalElements(), "incompatible vec size");
270 
271  return m_vec->size;
272 }
273 
274 double
276 {
277  return scalarProduct(*this,*this);
278 }
279 
280 double
282 {
283  return std::sqrt(this->norm2Sq());
284 }
285 
286 double
288 {
289  double result = 0.;
290 
291  unsigned int size = this->sizeLocal();
292  for (unsigned int i = 0; i < size; ++i) {
293  result += fabs((*this)[i]);
294  }
295 
296  return result;
297 }
298 
299 double
301 {
302  double result = 0.;
303 
304  unsigned int size = this->sizeLocal();
305  double aux = 0.;
306  for (unsigned int i = 0; i < size; ++i) {
307  aux = fabs((*this)[i]);
308  if (aux > result) result = aux;
309  }
310 
311  return result;
312 }
313 
314 double
316 {
317  double result = 0.;
318  unsigned int size = this->sizeLocal();
319  for (unsigned int i = 0; i < size; ++i) {
320  result += (*this)[i];
321  }
322 
323  return result;
324 }
325 
326 void
327 GslVector::cwSet(double value)
328 {
329  unsigned int size = this->sizeLocal();
330  for (unsigned int i = 0; i < size; ++i) {
331  (*this)[i] = value;
332  }
333 
334  return;
335 }
336 
337 void
338 GslVector::cwSetGaussian(double mean, double stdDev)
339 {
340  for (unsigned int i = 0; i < this->sizeLocal(); ++i) {
341  (*this)[i] = mean + m_env.rngObject()->gaussianSample(stdDev);
342  }
343 
344  return;
345 }
346 
347 void
348 GslVector::cwSetGaussian(const GslVector& meanVec, const GslVector& stdDevVec)
349 {
350  for (unsigned int i = 0; i < this->sizeLocal(); ++i) {
351  (*this)[i] = meanVec[i] + m_env.rngObject()->gaussianSample(stdDevVec[i]);
352  }
353  return;
354 }
355 
356 void
357 GslVector::cwSetUniform(const GslVector& aVec, const GslVector& bVec)
358 {
359  for (unsigned int i = 0; i < this->sizeLocal(); ++i) {
360  (*this)[i] = aVec[i] + (bVec[i]-aVec[i])*m_env.rngObject()->uniformSample();
361  }
362  return;
363 }
364 
365 void
366 GslVector::cwSetBeta(const GslVector& alpha, const GslVector& beta)
367 {
368  queso_require_equal_to_msg(this->sizeLocal(), alpha.sizeLocal(), "incompatible alpha size");
369 
370  queso_require_equal_to_msg(this->sizeLocal(), beta.sizeLocal(), "incompatible beta size");
371 
372  double tmpSample = 0.;
373  for (unsigned int i = 0; i < this->sizeLocal(); ++i) {
374  tmpSample = m_env.rngObject()->betaSample(alpha[i],beta[i]);
375  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
376  *m_env.subDisplayFile() << "In GslVector::cwSetBeta()"
377  << ": fullRank " << m_env.fullRank()
378  << ", i = " << i
379  << ", alpha[i] = " << alpha[i]
380  << ", beta[i] = " << beta[i]
381  << ", sample = " << tmpSample
382  << std::endl;
383  }
384  if ((alpha[i] == 1. ) &&
385  (beta [i] == 0.1)) {
386  if (tmpSample == 1.) {
387  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
388  *m_env.subDisplayFile() << "Hitting 'sampe = 1' in GslVector::cwSetBeta()"
389  << ": fullRank " << m_env.fullRank()
390  << ", i = " << i
391  << ", alpha[i] = " << alpha[i]
392  << ", beta[i] = " << beta[i]
393  << ", sample = " << tmpSample
394  << std::endl;
395  }
396 #if 1
397  std::cerr << "Hitting 'sample = 1' in GslVector::cwSetBeta()"
398  << ": fullRank " << m_env.fullRank()
399  << ", i = " << i
400  << ", alpha[i] = " << alpha[i]
401  << ", beta[i] = " << beta[i]
402  << ", sample = " << tmpSample
403  << std::endl;
404  do {
405  tmpSample = m_env.rngObject()->betaSample(alpha[i],beta[i]);
406  } while (tmpSample == 1.);
407  std::cerr << "Code was able to get 'sample != 1' in GslVector::cwSetBeta()"
408  << ": fullRank " << m_env.fullRank()
409  << ", i = " << i
410  << ", alpha[i] = " << alpha[i]
411  << ", beta[i] = " << beta[i]
412  << ", sample = " << tmpSample
413  << std::endl;
414  }
415 #endif
416  }
417  (*this)[i] = tmpSample;
418  }
419  return;
420 }
421 
422 void
424 {
425  queso_require_equal_to_msg(this->sizeLocal(), a.sizeLocal(), "incompatible a size");
426 
427  queso_require_equal_to_msg(this->sizeLocal(), b.sizeLocal(), "incompatible b size");
428 
429  for (unsigned int i = 0; i < this->sizeLocal(); ++i) {
430  (*this)[i] = m_env.rngObject()->gammaSample(a[i],b[i]);
431  }
432  return;
433 }
434 
435 void
437 {
438  queso_require_equal_to_msg(this->sizeLocal(), alpha.sizeLocal(), "incompatible alpha size");
439 
440  queso_require_equal_to_msg(this->sizeLocal(), beta.sizeLocal(), "incompatible beta size");
441 
442  for (unsigned int i = 0; i < this->sizeLocal(); ++i) {
443  (*this)[i] = 1./m_env.rngObject()->gammaSample(alpha[i],1./beta[i]);
444  }
445  return;
446 }
447 
448 void
450 {
451  queso_require_equal_to_msg(this->sizeLocal(), v1.sizeLocal() + v2.sizeLocal(), "incompatible vector sizes");
452 
453  for (unsigned int i = 0; i < v1.sizeLocal(); ++i) {
454  (*this)[i] = v1[i];
455  }
456 
457  for (unsigned int i = 0; i < v2.sizeLocal(); ++i) {
458  (*this)[v1.sizeLocal()+i] = v2[i];
459  }
460 
461  return;
462 }
463 
464 void
465 GslVector::cwSetConcatenated(const std::vector<const GslVector* >& vecs)
466 {
467  unsigned int cummulativeSize = 0;
468  for (unsigned int i = 0; i < vecs.size(); ++i) {
469  GslVector tmpVec(*(vecs[i]));
470  for (unsigned int j = 0; j < vecs[i]->sizeLocal(); ++j) {
471  (*this)[cummulativeSize+j] = tmpVec[j];
472  }
473  cummulativeSize += vecs[i]->sizeLocal();
474  }
475 
476  queso_require_equal_to_msg(this->sizeLocal(), cummulativeSize, "incompatible vector sizes");
477  return;
478 }
479 
480 
481 void
482 GslVector::cwSet(unsigned int initialPos, const GslVector& vec)
483 {
484  queso_require_less_msg(initialPos, this->sizeLocal(), "invalid initialPos");
485 
486  queso_require_less_equal_msg((initialPos +vec.sizeLocal()), this->sizeLocal(), "invalid vec.sizeLocal()");
487 
488  for (unsigned int i = 0; i < vec.sizeLocal(); ++i) {
489  (*this)[initialPos+i] = vec[i];
490  }
491 
492  return;
493 }
494 
495 void
496 GslVector::cwExtract(unsigned int initialPos, GslVector& vec) const
497 {
498  queso_require_less_msg(initialPos, this->sizeLocal(), "invalid initialPos");
499 
500  queso_require_less_equal_msg((initialPos +vec.sizeLocal()), this->sizeLocal(), "invalid vec.sizeLocal()");
501 
502  for (unsigned int i = 0; i < vec.sizeLocal(); ++i) {
503  vec[i] = (*this)[initialPos+i];
504  }
505 
506  return;
507 }
508 
509 void
511 {
512  unsigned int size = this->sizeLocal();
513  for (unsigned int i = 0; i < size; ++i) {
514  (*this)[i] = 1./(*this)[i];
515  }
516 
517  return;
518 }
519 
520 void
522 {
523  unsigned int size = this->sizeLocal();
524  for (unsigned int i = 0; i < size; ++i) {
525  (*this)[i] = sqrt((*this)[i]);
526  }
527 
528  return;
529 }
530 
531 void
533  unsigned int firstPositionToStoreDiff,
534  double valueForRemainderPosition,
535  GslVector& outputVec) const
536 {
537  unsigned int size = this->sizeLocal();
538 
539  queso_require_less_equal_msg(firstPositionToStoreDiff, 1, "invalid firstPositionToStoreDiff");
540 
541  queso_require_equal_to_msg(size, outputVec.sizeLocal(), "invalid size of outputVecs");
542 
543  for (unsigned int i = 0; i < (size-1); ++i) {
544  outputVec[firstPositionToStoreDiff+i] = (*this)[i+1]-(*this)[i];
545  }
546  if (firstPositionToStoreDiff == 0) {
547  outputVec[size-1] = valueForRemainderPosition;
548  }
549  else {
550  outputVec[0] = valueForRemainderPosition;
551  }
552 
553  return;
554 }
555 
556 void
558  const GslVector& x1Vec,
559  const GslVector& y1Vec,
560  const GslVector& x2Vec)
561 {
562  queso_require_greater_msg(x1Vec.sizeLocal(), 1, "invalid 'x1' size");
563 
564  queso_require_equal_to_msg(x1Vec.sizeLocal(), y1Vec.sizeLocal(), "invalid 'x1' and 'y1' sizes");
565 
566  queso_require_equal_to_msg(x2Vec.sizeLocal(), this->sizeLocal(), "invalid 'x2' and 'this' sizes");
567 
568  for (unsigned int i = 1; i < x1Vec.sizeLocal(); ++i) { // Yes, '1'
569  queso_require_greater_msg(x1Vec[i], x1Vec[i-1], "invalid 'x1' values");
570  }
571 
572  for (unsigned int id2 = 0; id2 < x2Vec.sizeLocal(); ++id2) {
573  double x2 = x2Vec[id2];
574  unsigned int id1 = 0;
575  bool found1 = false;
576  for (id1 = 0; id1 < x1Vec.sizeLocal(); ++id1) {
577  if (x2 <= x1Vec[id1]) {
578  found1 = true;
579  break;
580  }
581  }
582  bool makeLinearModel = false;
583  double xa = 0.;
584  double xb = 0.;
585  double ya = 0.;
586  double yb = 0.;
587  if (x2 == x1Vec[id1]) {
588  (*this)[id2] = y1Vec[id1];
589  }
590  else if (x2 < x1Vec[0]) {
591  // Extrapolation case
592  makeLinearModel = true;
593  xa = x1Vec[0];
594  xb = x1Vec[1];
595  ya = y1Vec[0];
596  yb = y1Vec[1];
597  }
598  else if (found1 == true) {
599  // Interpolation case
600  makeLinearModel = true;
601  xa = x1Vec[id1-1];
602  xb = x1Vec[id1];
603  ya = y1Vec[id1-1];
604  yb = y1Vec[id1];
605  }
606  else {
607  // Extrapolation case
608  makeLinearModel = true;
609  xa = x1Vec[x1Vec.sizeLocal()-2];
610  xb = x1Vec[x1Vec.sizeLocal()-1];
611  ya = y1Vec[x1Vec.sizeLocal()-2];
612  yb = y1Vec[x1Vec.sizeLocal()-1];
613  }
614 
615  if (makeLinearModel) {
616  double rate = (yb-ya)/(xb-xa);
617  (*this)[id2] = ya + (x2-xa)*rate;
618  }
619  }
620 
621 
622 
623  return;
624 }
625 
626 void
628 {
629  gsl_sort_vector(m_vec);
630 
631  return;
632 }
633 
634 void
635 GslVector::mpiBcast(int srcRank, const MpiComm& bcastComm)
636 {
637  // Filter out those nodes that should not participate
638  if (bcastComm.MyPID() < 0) return;
639 
640  // Check 'srcRank'
641  queso_require_msg(!((srcRank < 0) || (srcRank >= bcastComm.NumProc())), "invalud srcRank");
642 
643  // Check number of participant nodes
644  double localNumNodes = 1.;
645  double totalNumNodes = 0.;
646  bcastComm.Allreduce<double>(&localNumNodes, &totalNumNodes, (int) 1, RawValue_MPI_SUM,
647  "GslVector::mpiBcast()",
648  "failed MPI.Allreduce() for numNodes");
649  queso_require_equal_to_msg(((int) totalNumNodes), bcastComm.NumProc(), "inconsistent numNodes");
650 
651  // Check that all participant nodes have the same vector size
652  double localVectorSize = this->sizeLocal();
653  double sumOfVectorSizes = 0.;
654  bcastComm.Allreduce<double>(&localVectorSize, &sumOfVectorSizes, (int) 1, RawValue_MPI_SUM,
655  "GslVector::mpiBcast()",
656  "failed MPI.Allreduce() for vectorSize");
657 
658  if ( ((unsigned int) sumOfVectorSizes) != ((unsigned int)(totalNumNodes*localVectorSize)) ) {
659  std::cerr << "rank " << bcastComm.MyPID()
660  << ": sumOfVectorSizes = " << sumOfVectorSizes
661  << ", totalNumNodes = " << totalNumNodes
662  << ", localVectorSize = " << localVectorSize
663  << std::endl;
664  }
665  bcastComm.Barrier();
666  queso_require_equal_to_msg(((unsigned int) sumOfVectorSizes), ((unsigned int)(totalNumNodes*localVectorSize)), "inconsistent vectorSize");
667 
668  // Ok, bcast data
669  std::vector<double> dataBuffer((unsigned int) localVectorSize, 0.);
670  if (bcastComm.MyPID() == srcRank) {
671  for (unsigned int i = 0; i < dataBuffer.size(); ++i) {
672  dataBuffer[i] = (*this)[i];
673  }
674  }
675 
676  bcastComm.Bcast((void *) &dataBuffer[0], (int) localVectorSize, RawValue_MPI_DOUBLE, srcRank,
677  "GslVector::mpiBcast()",
678  "failed MPI.Bcast()");
679 
680  if (bcastComm.MyPID() != srcRank) {
681  for (unsigned int i = 0; i < dataBuffer.size(); ++i) {
682  (*this)[i] = dataBuffer[i];
683  }
684  }
685 
686  return;
687 }
688 
689 void
690 GslVector::mpiAllReduce(RawType_MPI_Op mpiOperation, const MpiComm& opComm, GslVector& resultVec) const
691 {
692  // Filter out those nodes that should not participate
693  if (opComm.MyPID() < 0) return;
694 
695  unsigned int size = this->sizeLocal();
696  queso_require_equal_to_msg(size, resultVec.sizeLocal(), "different vector sizes");
697 
698  for (unsigned int i = 0; i < size; ++i) {
699  double srcValue = (*this)[i];
700  double resultValue = 0.;
701  opComm.Allreduce<double>(&srcValue, &resultValue, (int) 1, mpiOperation,
702  "GslVector::mpiAllReduce()",
703  "failed MPI.Allreduce()");
704  resultVec[i] = resultValue;
705  }
706 
707  return;
708 }
709 
710 void
711 GslVector::mpiAllQuantile(double probability, const MpiComm& opComm, GslVector& resultVec) const
712 {
713  // Filter out those nodes that should not participate
714  if (opComm.MyPID() < 0) return;
715 
716  queso_require_msg(!((probability < 0.) || (1. < probability)), "invalid input");
717 
718  unsigned int size = this->sizeLocal();
719  queso_require_equal_to_msg(size, resultVec.sizeLocal(), "different vector sizes");
720 
721  for (unsigned int i = 0; i < size; ++i) {
722  double auxDouble = (int) (*this)[i];
723  std::vector<double> vecOfDoubles(opComm.NumProc(),0.);
724  opComm.Gather<double>(&auxDouble, 1, &vecOfDoubles[0], (int) 1, 0,
725  "GslVector::mpiAllQuantile()",
726  "failed MPI.Gather()");
727 
728  std::sort(vecOfDoubles.begin(), vecOfDoubles.end());
729 
730  double result = vecOfDoubles[(unsigned int)( probability*((double)(vecOfDoubles.size()-1)) )];
731 
732  opComm.Bcast((void *) &result, (int) 1, RawValue_MPI_DOUBLE, 0,
733  "GslVector::mpiAllQuantile()",
734  "failed MPI.Bcast()");
735 
736  resultVec[i] = result;
737  }
738 
739  return;
740 }
741 
742 void
743 GslVector::print(std::ostream& os) const
744 {
745  //std::cout << "In GslVector::print(): before sizelocal()"
746  // << std::endl;
747  unsigned int size = this->sizeLocal();
748  //std::cout << "In GslVector::print(): after sizelocal()"
749  // << std::endl;
750 
751  //std::cout << "In GslVector::print(): before os.flags()"
752  // << std::endl;
753  std::ostream::fmtflags curr_fmt = os.flags();
754  //std::cout << "In GslVector::print(): after os.flags()"
755  // << std::endl;
756 
757  if (m_printScientific) {
758  unsigned int savedPrecision = os.precision();
759  os.precision(16);
760 
761  if (m_printHorizontally) {
762  for (unsigned int i = 0; i < size; ++i) {
763  os << std::scientific << (*this)[i]
764  << " ";
765  }
766  }
767  else {
768  for (unsigned int i = 0; i < size; ++i) {
769  os << std::scientific << (*this)[i]
770  << std::endl;
771  }
772  }
773 
774  os.precision(savedPrecision);
775  }
776  else {
777  if (m_printHorizontally) {
778  //std::cout << "In GslVector::print(): where expected"
779  // << std::endl;
780  for (unsigned int i = 0; i < size; ++i) {
781  os << std::dec << (*this)[i]
782  << " ";
783  }
784  }
785  else {
786  for (unsigned int i = 0; i < size; ++i) {
787  os << std::dec << (*this)[i]
788  << std::endl;
789  }
790  }
791  }
792 
793  //std::cout << "In GslVector::print(): before os.flags(curr_fmt)"
794  // << std::endl;
795  os.flags(curr_fmt);
796  //std::cout << "In GslVector::print(): after os.flags(curr_fmt)"
797  // << std::endl;
798 
799  return;
800 }
801 
802 void
804  const std::string& varNamePrefix,
805  const std::string& fileName,
806  const std::string& fileType,
807  const std::set<unsigned int>& allowedSubEnvIds) const
808 {
809  queso_require_greater_equal_msg(m_env.subRank(), 0, "unexpected subRank");
810 
811  queso_require_less_equal_msg(this->numOfProcsForStorage(), 1, "implemented just for sequential vectors for now");
812 
813  FilePtrSetStruct filePtrSet;
814  if (m_env.openOutputFile(fileName,
815  fileType, // "m or hdf"
816  allowedSubEnvIds,
817  false,
818  filePtrSet)) {
819  *filePtrSet.ofsVar << varNamePrefix << "_sub" << m_env.subIdString() << " = zeros(" << this->sizeLocal()
820  << "," << 1
821  << ");"
822  << std::endl;
823  *filePtrSet.ofsVar << varNamePrefix << "_sub" << m_env.subIdString() << " = [";
824 
825  bool savedVectorPrintScientific = this->getPrintScientific();
826  bool savedVectorPrintHorizontally = this->getPrintHorizontally();
827  this->setPrintScientific (true);
828  this->setPrintHorizontally(false);
829  *filePtrSet.ofsVar << *this;
830  //<< std::endl; // No need for 'endl' because horizontally = 'false'
831  this->setPrintHorizontally(savedVectorPrintHorizontally);
832  this->setPrintScientific (savedVectorPrintScientific);
833 
834  *filePtrSet.ofsVar << "];\n";
835 
836  m_env.closeFile(filePtrSet,fileType);
837  }
838 
839  return;
840 }
841 
842 void
844  const std::string& fileName,
845  const std::string& fileType,
846  const std::set<unsigned int>& allowedSubEnvIds)
847 {
848  queso_require_greater_equal_msg(m_env.subRank(), 0, "unexpected subRank");
849 
850  queso_require_less_equal_msg(this->numOfProcsForStorage(), 1, "implemented just for sequential vectors for now");
851 
852  FilePtrSetStruct filePtrSet;
853  if (m_env.openInputFile(fileName,
854  fileType, // "m or hdf"
855  allowedSubEnvIds,
856  filePtrSet)) {
857  double subReadSize = this->sizeLocal();
858 
859  // In the logic below, the id of a line' begins with value 0 (zero)
860  unsigned int idOfMyFirstLine = 1;
861  unsigned int idOfMyLastLine = this->sizeLocal();
862  unsigned int numParams = 1; // Yes, just '1'
863 
864  // Read number of chain positions in the file by taking care of the first line,
865  // which resembles something like 'variable_name = zeros(n_positions,m_params);'
866  std::string tmpString;
867 
868  // Read 'variable name' string
869  *filePtrSet.ifsVar >> tmpString;
870  //std::cout << "Just read '" << tmpString << "'" << std::endl;
871 
872  // Read '=' sign
873  *filePtrSet.ifsVar >> tmpString;
874  //std::cout << "Just read '" << tmpString << "'" << std::endl;
875  queso_require_equal_to_msg(tmpString, std::string("="), std::string("string should be the '=' sign"));
876 
877  // Read 'zeros(n_positions,n_params)' string
878  *filePtrSet.ifsVar >> tmpString;
879  //std::cout << "Just read '" << tmpString << "'" << std::endl;
880  unsigned int posInTmpString = 6;
881 
882  // Isolate 'n_positions' in a string
883  char nPositionsString[tmpString.size()-posInTmpString+1];
884  unsigned int posInPositionsString = 0;
885  do {
886  queso_require_less_msg(posInTmpString, tmpString.size(), "symbol ',' not found in first line of file");
887  nPositionsString[posInPositionsString++] = tmpString[posInTmpString++];
888  } while (tmpString[posInTmpString] != ',');
889  nPositionsString[posInPositionsString] = '\0';
890 
891  // Isolate 'n_params' in a string
892  posInTmpString++; // Avoid reading ',' char
893  char nParamsString[tmpString.size()-posInTmpString+1];
894  unsigned int posInParamsString = 0;
895  do {
896  queso_require_less_msg(posInTmpString, tmpString.size(), "symbol ')' not found in first line of file");
897  nParamsString[posInParamsString++] = tmpString[posInTmpString++];
898  } while (tmpString[posInTmpString] != ')');
899  nParamsString[posInParamsString] = '\0';
900 
901  // Convert 'n_positions' and 'n_params' strings to numbers
902  unsigned int sizeOfVecInFile = (unsigned int) strtod(nPositionsString,NULL);
903  unsigned int numParamsInFile = (unsigned int) strtod(nParamsString, NULL);
904  if (m_env.subDisplayFile()) {
905  *m_env.subDisplayFile() << "In GslVector::subReadContents()"
906  << ": fullRank " << m_env.fullRank()
907  << ", sizeOfVecInFile = " << sizeOfVecInFile
908  << ", numParamsInFile = " << numParamsInFile
909  << ", this->sizeLocal() = " << this->sizeLocal()
910  << std::endl;
911  }
912 
913  // Check if [size of vec in file] >= [requested sub vec size]
914  queso_require_greater_equal_msg(sizeOfVecInFile, subReadSize, "size of vec in file is not big enough");
915 
916  // Check if [num params in file] == [num params in current vec]
917  queso_require_equal_to_msg(numParamsInFile, numParams, "number of parameters of vec in file is different than number of parameters in this vec object");
918 
919  // Code common to any core in a communicator
920  unsigned int maxCharsPerLine = 64*numParams; // Up to about 60 characters to represent each parameter value
921 
922  unsigned int lineId = 0;
923  while (lineId < idOfMyFirstLine) {
924  filePtrSet.ifsVar->ignore(maxCharsPerLine,'\n');
925  lineId++;
926  };
927 
928  if (m_env.subDisplayFile()) {
929  *m_env.subDisplayFile() << "In GslVector::subReadContents()"
930  << ": beginning to read input actual data"
931  << std::endl;
932  }
933 
934  // Take care of initial part of the first data line,
935  // which resembles something like 'variable_name = [value1 value2 ...'
936  // Read 'variable name' string
937  *filePtrSet.ifsVar >> tmpString;
938  //std::cout << "Core 0 just read '" << tmpString << "'" << std::endl;
939 
940  // Read '=' sign
941  *filePtrSet.ifsVar >> tmpString;
942  //std::cout << "Core 0 just read '" << tmpString << "'" << std::endl;
943  queso_require_equal_to_msg(tmpString, std::string("="), std::string("in core 0, string should be the '=' sign"));
944 
945  // Take into account the ' [' portion
946  std::streampos tmpPos = filePtrSet.ifsVar->tellg();
947  filePtrSet.ifsVar->seekg(tmpPos+(std::streampos)2);
948 
949  if (m_env.subDisplayFile()) {
950  *m_env.subDisplayFile() << "In GslVector::subReadContents()"
951  << ": beginning to read lines with numbers only"
952  << ", lineId = " << lineId
953  << ", idOfMyFirstLine = " << idOfMyFirstLine
954  << ", idOfMyLastLine = " << idOfMyLastLine
955  << std::endl;
956  }
957 
958  while (lineId <= idOfMyLastLine) {
959  *filePtrSet.ifsVar >> (*this)[lineId - idOfMyFirstLine];
960  lineId++;
961  };
962 
963  m_env.closeFile(filePtrSet,fileType);
964  }
965 
966  return;
967 }
968 
969 gsl_vector*
971 {
972  return m_vec;
973 }
974 
975 bool
977 {
978  queso_require_equal_to_msg(this->sizeLocal(), rhs.sizeLocal(), "vectors have different sizes");
979 
980  bool result = false;
981  unsigned int i = 0;
982  unsigned int size = this->sizeLocal();
983  while ((i < size) && (result == false)) {
984  result = ( (*this)[i] < rhs[i] );
985  i++;
986  };
987 
988  return result;
989 }
990 
991 bool
993 {
994  queso_require_equal_to_msg(this->sizeLocal(), rhs.sizeLocal(), "vectors have different sizes");
995 
996  bool result = false;
997  unsigned int i = 0;
998  unsigned int size = this->sizeLocal();
999  while ((i < size) && (result == false)) {
1000  result = ( (*this)[i] > rhs[i] );
1001  i++;
1002  };
1003 
1004  return result;
1005 }
1006 
1007 bool
1009 {
1010  queso_require_equal_to_msg(this->sizeLocal(), rhs.sizeLocal(), "vectors have different sizes");
1011 
1012  bool result = false;
1013  unsigned int i = 0;
1014  unsigned int size = this->sizeLocal();
1015  while ((i < size) && (result == false)) {
1016  result = ( (*this)[i] <= rhs[i] ); // prudencio 2012-02-06
1017  i++;
1018  };
1019 
1020  return result;
1021 }
1022 
1023 bool
1025 {
1026  queso_require_equal_to_msg(this->sizeLocal(), rhs.sizeLocal(), "vectors have different sizes");
1027 
1028  bool result = false;
1029  unsigned int i = 0;
1030  unsigned int size = this->sizeLocal();
1031  while ((i < size) && (result == false)) {
1032  result = ( (*this)[i] >= rhs[i] ); // prudencio 2012-02-06
1033  i++;
1034  };
1035 
1036  return result;
1037 }
1038 
1039 double
1041 {
1042  return gsl_vector_max( m_vec );
1043 }
1044 
1045 double
1047 {
1048  return gsl_vector_min( m_vec );
1049 }
1050 
1051 int
1053 {
1054  return gsl_vector_max_index( m_vec );
1055 }
1056 
1057 int
1059 {
1060  return gsl_vector_min_index( m_vec );
1061 }
1062 
1063 void
1064 GslVector::getMaxValueAndIndex( double& max_value, int& max_value_index )
1065 {
1066  max_value = this->getMaxValue();
1067  max_value_index = this->getMaxValueIndex();
1068 
1069  return;
1070 }
1071 
1072 void
1073 GslVector::getMinValueAndIndex( double& min_value, int& min_value_index )
1074 {
1075  min_value = this->getMinValue();
1076  min_value_index = this->getMinValueIndex();
1077 
1078  return;
1079 }
1080 
1081 GslVector
1083 {
1084  GslVector abs_of_this_vec( *this );
1085 
1086  unsigned int size = abs_of_this_vec.sizeLocal();
1087 
1088  for( unsigned int i = 0; i < size; ++i )
1089  {
1090  abs_of_this_vec[i] = std::fabs( (*this)[i] );
1091  }
1092 
1093  return abs_of_this_vec;
1094 }
1095 
1096 std::ostream&
1097 operator<<(std::ostream& os, const GslVector& obj)
1098 {
1099  obj.print(os);
1100 
1101  return os;
1102 }
1103 
1104 GslVector operator/(double a, const GslVector& x)
1105 {
1106  GslVector answer(x);
1107  answer.cwInvert();
1108  answer *= a;
1109 
1110  return answer;
1111 }
1112 
1114 {
1115  GslVector answer(x);
1116  answer /= y;
1117 
1118  return answer;
1119 }
1120 
1121 GslVector operator*(double a, const GslVector& x)
1122 {
1123  GslVector answer(x);
1124  answer *= a;
1125 
1126  return answer;
1127 }
1128 
1130 {
1131  GslVector answer(x);
1132  answer *= y;
1133 
1134  return answer;
1135 }
1136 
1137 double scalarProduct(const GslVector& x, const GslVector& y)
1138 {
1139  unsigned int size1 = x.sizeLocal();
1140  unsigned int size2 = y.sizeLocal();
1141  queso_require_equal_to_msg(size1, size2, "different sizes of x and y");
1142 
1143  double result = 0.;
1144  for (unsigned int i = 0; i < size1; ++i) {
1145  result += x[i]*y[i];
1146  }
1147 
1148  return result;
1149 }
1150 
1152 {
1153  GslVector answer(x);
1154  answer += y;
1155 
1156  return answer;
1157 }
1158 
1160 {
1161  GslVector answer(x);
1162  answer -= y;
1163 
1164  return answer;
1165 }
1166 
1167 bool
1168 operator== (const GslVector& lhs, const GslVector& rhs)
1169 {
1170  bool answer = true;
1171 
1172  unsigned int size1 = lhs.sizeLocal();
1173  unsigned int size2 = rhs.sizeLocal();
1174  queso_require_equal_to_msg(size1, size2, "different sizes of lhs and rhs");
1175 
1176  for (unsigned int i = 0; i < size1; ++i) {
1177  if (lhs[i] != rhs[i]) {
1178  answer = false;
1179  break;
1180  }
1181  }
1182 
1183  return answer;
1184 }
1185 
1186 } // End namespace QUESO
GslMatrix operator+(const GslMatrix &m1, const GslMatrix &m2)
Definition: GslMatrix.C:2022
virtual double gammaSample(double a, double b) const =0
Samples a value from a Gamma distribution.
void matlabLinearInterpExtrap(const GslVector &x1Vec, const GslVector &y1Vec, const GslVector &x2Vec)
Definition: GslVector.C:557
void cwSetBeta(const GslVector &alpha, const GslVector &beta)
This function returns a random variate from the beta distribution, with vector parameters alpha and b...
Definition: GslVector.C:366
int NumGlobalElements() const
Returns the total number of elements across all processors.
Definition: Map.C:85
int getMinValueIndex() const
This function returns the index of the minimum value in the vector this.
Definition: GslVector.C:1058
double norm2() const
Returns the 2-norm (Euclidean norm) of the vector.
Definition: GslVector.C:281
void cwInvert()
This function inverts component-wise the element values of this.
Definition: GslVector.C:510
std::ostream & operator<<(std::ostream &os, const SequenceStatisticalOptions &obj)
std::ifstream * ifsVar
Provides a stream interface to read data from files.
Definition: Environment.h:89
const BaseEnvironment & m_env
Environment variable.
Definition: Vector.h:117
virtual void base_copy(const Vector &src)
Copies base data from vector src to this vector.
Definition: Vector.C:101
bool atLeastOneComponentSmallerThan(const GslVector &rhs) const
This function returns true if at least one component of this is smaller than the respective component...
Definition: GslVector.C:976
void mpiBcast(int srcRank, const MpiComm &bcastComm)
Definition: GslVector.C:635
void cwSet(double value)
Component-wise sets all values to this with value.
Definition: GslVector.C:327
virtual double gaussianSample(double stdDev) const =0
Samples a value from a Gaussian distribution with standard deviation given by stdDev.
GslVector abs() const
This function returns absolute value of elements in this.
Definition: GslVector.C:1082
bool m_printHorizontally
Flag for either or not print this matrix horizontally.
Definition: Vector.h:130
void cwSetUniform(const GslVector &aVec, const GslVector &bVec)
This function sets component-wise a number uniformly distributed in the range of elements of [aVec...
Definition: GslVector.C:357
A class for partitioning vectors and matrices.
Definition: Map.h:49
double sumOfComponents() const
Returns the sum of the components of the vector.
Definition: GslVector.C:315
The QUESO MPI Communicator Class.
Definition: MpiComm.h:203
void Gather(void *sendbuf, int sendcnt, RawType_MPI_Datatype sendtype, void *recvbuf, int recvcount, RawType_MPI_Datatype recvtype, int root, const char *whereMsg, const char *whatMsg) const
Gather values from each process to collect on all processes.
Definition: MpiComm.C:202
double getMaxValue() const
Returns the maximum value in the vector this.
Definition: GslVector.C:1040
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:348
void setPrintScientific(bool value) const
Determines whether vector should be printed in Scientific Notation.
Definition: Vector.C:88
double getMinValue() const
Returns minimum value in the vector this.
Definition: GslVector.C:1046
bool m_printScientific
Flag for either or not print this matrix in scientific notation.
Definition: Vector.h:133
GslVector & operator+=(const GslVector &rhs)
Stores in this the coordinate-wise addition of this and rhs.
Definition: GslVector.C:211
double norm2Sq() const
Returns the 2-norm squared of this vector.
Definition: GslVector.C:275
void cwSetGaussian(double mean, double stdDev)
This function sets component-wise Gaussian random variates, with mean mean and standard deviation std...
Definition: GslVector.C:338
int fullRank() const
Returns the rank of the MPI process in QUESO&#39;s full communicator.
Definition: Environment.C:268
void cwExtract(unsigned int initialPos, GslVector &vec) const
This function sets the values of this starting at position initialPos ans saves them in vector vec...
Definition: GslVector.C:496
bool getPrintScientific() const
Checks if the vector should be printed in Scientific Notation.
Definition: Vector.C:95
void cwSetConcatenated(const GslVector &v1, const GslVector &v2)
This function concatenates GslVector v1 and GslVector v2 into this.
Definition: GslVector.C:449
double scalarProduct(const GslVector &x, const GslVector &y)
Definition: GslVector.C:1137
void getMinValueAndIndex(double &value, int &index)
This function returns minimum value in the vector this and its the index.
Definition: GslVector.C:1073
void setPrintHorizontally(bool value) const
Determines whether vector should be printed horizontally.
Definition: Vector.C:75
double norm1() const
Returns the 1-norm of the vector.
Definition: GslVector.C:287
void subWriteContents(const std::string &varNamePrefix, const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds) const
Definition: GslVector.C:803
std::ofstream * ofsVar
Provides a stream interface to write data to files.
Definition: Environment.h:86
void copy(const GslVector &src)
This function copies the elements of the vector src into this.
Definition: GslVector.C:230
virtual double betaSample(double alpha, double beta) const =0
Samples a value from a Beta distribution.
const RngBase * rngObject() const
Access to the RNG object.
Definition: Environment.C:471
gsl_vector * m_vec
GSL vector.
Definition: GslVector.h:253
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:521
GslVector & operator*=(double a)
Stores in this the coordinate-wise multiplication of this and a.
Definition: GslVector.C:166
GslMatrix operator-(const GslMatrix &m1, const GslMatrix &m2)
Definition: GslMatrix.C:2029
bool operator==(const GslVector &lhs, const GslVector &rhs)
Definition: GslVector.C:1168
Class for vector operations using GSL library.
Definition: GslVector.h:49
void mpiAllReduce(RawType_MPI_Op mpiOperation, const MpiComm &opComm, GslVector &resultVec) const
Definition: GslVector.C:690
double normInf() const
Returns the infinity-norm (maximum norm) of the vector.
Definition: GslVector.C:300
const Map & map() const
Definition: Vector.C:61
Class for vector operations (virtual).
Definition: Vector.h:47
GslVector()
Default Constructor.
void Allreduce(void *sendbuf, void *recvbuf, int count, RawType_MPI_Datatype datatype, RawType_MPI_Op op, const char *whereMsg, const char *whatMsg) const
Combines values from all processes and distributes the result back to all processes.
Definition: MpiComm.C:143
void Barrier() const
Pause every process in *this communicator until all the processes reach this point.
Definition: MpiComm.C:174
void subReadContents(const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds)
Definition: GslVector.C:843
bool atLeastOneComponentBiggerOrEqualThan(const GslVector &rhs) const
This function returns true if at least one component of this is bigger than or equal to the respectiv...
Definition: GslVector.C:1024
void matlabDiff(unsigned int firstPositionToStoreDiff, double valueForRemainderPosition, GslVector &outputVec) const
Definition: GslVector.C:532
int getMaxValueIndex() const
This function returns the index of the maximum value in the vector this.
Definition: GslVector.C:1052
void cwSqrt()
This function returns component-wise the square-root of this.
Definition: GslVector.C:521
unsigned int sizeGlobal() const
Returns the global length of this vector.
Definition: GslVector.C:267
MPI_Op RawType_MPI_Op
Definition: MpiComm.h:45
gsl_vector * data() const
Definition: GslVector.C:970
GslMatrix operator*(double a, const GslMatrix &mat)
Definition: GslMatrix.C:1982
unsigned int numOfProcsForStorage() const
Definition: Vector.C:68
Struct for handling data input and output from files.
Definition: Environment.h:77
int subRank() const
Returns the rank of the MPI process in the sub-communicator subComm()
Definition: Environment.C:287
void getMaxValueAndIndex(double &value, int &index)
This function returns maximum value in the vector this and its the index.
Definition: GslVector.C:1064
GslVector & operator=(const GslVector &rhs)
Copies values from vector rhs to this.
Definition: GslVector.C:150
bool atLeastOneComponentBiggerThan(const GslVector &rhs) const
This function returns true if at least one component of this is bigger than the respective component ...
Definition: GslVector.C:992
GslVector & operator-=(const GslVector &rhs)
Stores in this the coordinate-wise subtraction of this by rhs.
Definition: GslVector.C:220
void closeFile(FilePtrSetStruct &filePtrSet, const std::string &fileType) const
Closes the file.
Definition: Environment.C:1084
GslVector operator/(double a, const GslVector &x)
Definition: GslVector.C:1104
int MyPID() const
Return my process ID.
Definition: MpiComm.C:124
virtual double uniformSample() const =0
Samples a value from a uniform distribution.
void sort()
This function sorts the elements of the vector this in ascending numerical order. ...
Definition: GslVector.C:627
int NumMyElements() const
Returns the number of elements owned by the calling processor.
Definition: Map.C:107
void mpiAllQuantile(double probability, const MpiComm &opComm, GslVector &resultVec) const
Definition: GslVector.C:711
int NumProc() const
Returns total number of processes.
Definition: MpiComm.C:133
MonteCarloSGOptions::MonteCarloSGOptions(const BaseEnvironment &env, const char *prefix, const McOptionsValues &alternativeOptionsValues queso_require_equal_to_msg)(m_env.optionsInputFileName(), std::string(""), std::string("this constructor is incompatible with the existence of an options input file"))
bool atLeastOneComponentSmallerOrEqualThan(const GslVector &rhs) const
This function returns true if at least one component of this is smaller than or equal to the respecti...
Definition: GslVector.C:1008
void print(std::ostream &os) const
Print method. Defines the behavior of the std::ostream &lt;&lt; operator inherited from the Object class...
Definition: GslVector.C:743
unsigned int sizeLocal() const
Returns the length of this vector.
Definition: GslVector.C:241
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:320
bool openInputFile(const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds, FilePtrSetStruct &filePtrSet) const
Opens an input file.
Definition: Environment.C:896
const Map m_map
Mapping variable.
Definition: Vector.h:122
bool getPrintHorizontally() const
Checks if vector is printed horizontally.
Definition: Vector.C:82
~GslVector()
Destructor.
Definition: GslVector.C:144
void cwSetInverseGamma(const GslVector &alpha, const GslVector &beta)
This function returns a random variate from the inverse gamma distribution with vector parameters alp...
Definition: GslVector.C:436
unsigned int displayVerbosity() const
Definition: Environment.C:450
void Bcast(void *buffer, int count, RawType_MPI_Datatype datatype, int root, const char *whereMsg, const char *whatMsg) const
Broadcast values from the root process to the slave processes.
Definition: MpiComm.C:191
GslVector & operator/=(double a)
Stores in this the coordinate-wise division of this by a.
Definition: GslVector.C:175
void cwSetGamma(const GslVector &a, const GslVector &b)
This function returns a random variate from the gamma distribution with vector parameters a and b...
Definition: GslVector.C:423
This (virtual) class sets up the environment underlying the use of the QUESO library by an executable...
Definition: Environment.h:198

Generated on Sat Apr 22 2017 14:04:35 for queso-0.57.0 by  doxygen 1.8.5