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

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