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

Generated on Thu Dec 15 2016 13:23:10 for queso-0.56.1 by  doxygen 1.8.5