queso-0.51.1
1D1DFunction.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/1D1DFunction.h>
26 #include <queso/1DQuadrature.h>
27 
28 namespace QUESO {
29 
30 //*****************************************************
31 // Base 1D->1D class
32 //*****************************************************
34  double minDomainValue,
35  double maxDomainValue)
36  :
37  m_minDomainValue(minDomainValue),
38  m_maxDomainValue(maxDomainValue)
39 {
42  "Base1D1DFunction::constructor()",
43  "min >= max");
44 }
45 
47 {
48 }
49 
50 double
52 {
53  return m_minDomainValue;
54 }
55 
56 double
58 {
59  return m_maxDomainValue;
60 }
61 
62 double
63 Base1D1DFunction::multiplyAndIntegrate(const Base1D1DFunction& func, unsigned int quadratureOrder, double* resultWithMultiplicationByTAsWell) const
64 {
65  double value = 0.;
66 
69  "Base1D1DFunction::multiplyAndIntegrate()",
70  "not implemented yet");
71 
72  if (resultWithMultiplicationByTAsWell) { // Just to eliminate INTEL compiler warnings
73  func.value((double) quadratureOrder);
74  }
75 
76  return value;
77 }
78 
79 //*****************************************************
80 // Generic 1D->1D class
81 //*****************************************************
83  double minDomainValue,
84  double maxDomainValue,
85  double (*valueRoutinePtr)(double domainValue, const void* routinesDataPtr),
86  double (*derivRoutinePtr)(double domainValue, const void* routinesDataPtr),
87  const void* routinesDataPtr)
88  :
89  Base1D1DFunction(minDomainValue,maxDomainValue),
90  m_valueRoutinePtr (valueRoutinePtr),
91  m_derivRoutinePtr (derivRoutinePtr),
92  m_routinesDataPtr (routinesDataPtr)
93 {
94 }
95 
97 {
98 }
99 
100 double
101 Generic1D1DFunction::value(double domainValue) const
102 {
103  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
104  std::cerr << "In Generic1D1DFunction::value()"
105  << ": requested x (" << domainValue
106  << ") is out of the interval (" << m_minDomainValue
107  << ", " << m_maxDomainValue
108  << ")"
109  << std::endl;
110  }
111 
112  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
114  "Generic1D1DFunction::value()",
115  "x out of range");
116 
117  return (*m_valueRoutinePtr)(domainValue,m_routinesDataPtr);
118 }
119 
120 double
121 Generic1D1DFunction::deriv(double domainValue) const
122 {
123  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
124  std::cerr << "In Generic1D1DFunction::deriv()"
125  << ": requested x (" << domainValue
126  << ") is out of the interval (" << m_minDomainValue
127  << ", " << m_maxDomainValue
128  << ")"
129  << std::endl;
130  }
131 
132  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
134  "Generic1D1DFunction::deriv()",
135  "x out of range");
136 
137  return (*m_derivRoutinePtr)(domainValue,m_routinesDataPtr);
138 }
139 
140 //*****************************************************
141 // Constant 1D->1D class
142 //*****************************************************
144  double minDomainValue,
145  double maxDomainValue,
146  double constantValue)
147  :
148  Base1D1DFunction(minDomainValue,maxDomainValue),
149  m_constantValue (constantValue)
150 {
151 }
152 
154 {
155 }
156 
157 double
158 Constant1D1DFunction::value(double domainValue) const
159 {
160  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
161  std::cerr << "In Constant1D1DFunction::value()"
162  << ": requested x (" << domainValue
163  << ") is out of the interval (" << m_minDomainValue
164  << ", " << m_maxDomainValue
165  << ")"
166  << std::endl;
167  }
168 
169  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
171  "Constant1D1DFunction::value()",
172  "x out of range");
173 
174  return m_constantValue;
175 }
176 
177 double
178 Constant1D1DFunction::deriv(double domainValue) const
179 {
180  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
181  std::cerr << "In Constant1D1DFunction::deriv()"
182  << ": requested x (" << domainValue
183  << ") is out of the interval (" << m_minDomainValue
184  << ", " << m_maxDomainValue
185  << ")"
186  << std::endl;
187  }
188 
189  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
191  "Constant1D1DFunction::deriv()",
192  "x out of range");
193 
194  return 0.;
195 }
196 
197 //*****************************************************
198 // Linear 1D->1D class
199 //*****************************************************
201  double minDomainValue,
202  double maxDomainValue,
203  double referenceDomainValue,
204  double referenceImageValue,
205  double rateValue)
206  :
207  Base1D1DFunction(minDomainValue,maxDomainValue),
208  m_referenceDomainValue (referenceDomainValue),
209  m_referenceImageValue (referenceImageValue),
210  m_rateValue (rateValue)
211 {
212 }
213 
215 {
216 }
217 
218 double
219 Linear1D1DFunction::value(double domainValue) const
220 {
221  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
222  std::cerr << "In Linear1D1DFunction::value()"
223  << ": requested x (" << domainValue
224  << ") is out of the interval (" << m_minDomainValue
225  << ", " << m_maxDomainValue
226  << ")"
227  << std::endl;
228  }
229 
230  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
232  "Linear1D1DFunction::value()",
233  "x out of range");
234 
235  double imageValue = m_referenceImageValue + m_rateValue*(domainValue - m_referenceDomainValue);
236 
237  return imageValue;
238 }
239 
240 double
241 Linear1D1DFunction::deriv(double domainValue) const
242 {
243  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
244  std::cerr << "In Linear1D1DFunction::deriv()"
245  << ": requested x (" << domainValue
246  << ") is out of the interval (" << m_minDomainValue
247  << ", " << m_maxDomainValue
248  << ")"
249  << std::endl;
250  }
251 
252  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
254  "Linear1D1DFunction::deriv()",
255  "x out of range");
256 
257  return m_rateValue;
258 }
259 
260 //*****************************************************
261 // PiecewiseLinear 1D->1D class
262 //*****************************************************
264  double minDomainValue,
265  double maxDomainValue,
266  const std::vector<double>& referenceDomainValues,
267  double referenceImageValue0,
268  const std::vector<double>& rateValues)
269  :
270  Base1D1DFunction(minDomainValue,maxDomainValue),
271  m_numRefValues (referenceDomainValues.size()),
272  m_referenceDomainValues(referenceDomainValues),
273  m_rateValues (rateValues)
274 {
277  "PiecewiseLinear1D1DFunction::constructor()",
278  "num ref values = 0");
279 
280  UQ_FATAL_TEST_MACRO(m_numRefValues != rateValues.size(),
282  "PiecewiseLinear1D1DFunction::constructor()",
283  "num rate values is inconsistent");
284 
285  for (unsigned int i = 1; i < m_numRefValues; ++i) { // Yes, from '1'
288  "PiecewiseLinear1D1DFunction::constructor()",
289  "reference domain values are inconsistent");
290  }
291 
292  m_referenceImageValues.clear();
293  m_referenceImageValues.resize(m_numRefValues,0.);
294  m_referenceImageValues[0] = referenceImageValue0;
295  for (unsigned int i = 1; i < m_numRefValues; ++i) { // Yes, from '1'
297  }
298 
299  if (false) { // For debug only
300  std::cout << "In PiecewiseLinear1D1DFunction::constructor():"
301  << std::endl;
302  for (unsigned int i = 0; i < m_numRefValues; ++i) {
303  std::cout << "i = " << i
304  << ", m_referenceDomainValues[i] = " << m_referenceDomainValues[i]
305  << ", m_referenceImageValues[i] = " << m_referenceImageValues[i]
306  << ", m_rateValues[i] = " << m_rateValues[i]
307  << std::endl;
308  }
309  }
310 }
311 
313 {
314  m_rateValues.clear();
315  m_referenceImageValues.clear();
316  m_referenceDomainValues.clear();
317 }
318 
319 double
320 PiecewiseLinear1D1DFunction::value(double domainValue) const
321 {
322  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
323  std::cerr << "In PiecewiseLinear1D1DFunction::value()"
324  << ": requested x (" << domainValue
325  << ") is out of the interval (" << m_minDomainValue
326  << ", " << m_maxDomainValue
327  << ")"
328  << std::endl;
329  }
330 
331  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
333  "PiecewiseLinear1D1DFunction::value()",
334  "x out of range");
335 
336  unsigned int i = 0;
337  if (m_numRefValues == 1) {
338  // Nothing else to do
339  }
340  else {
341  bool referenceDomainValueFound = false;
342  while (referenceDomainValueFound == false) {
343  if (domainValue < m_referenceDomainValues[i+1]) {
344  referenceDomainValueFound = true;
345  }
346  else {
347  ++i;
348  if (i == (m_numRefValues-1)) {
349  referenceDomainValueFound = true;
350  }
351  }
352  }
353  }
354  double imageValue = m_referenceImageValues[i] + m_rateValues[i]*(domainValue - m_referenceDomainValues[i]);
355  if (false) { // For debug only
356  std::cout << "In PiecewiseLinear1D1DFunction::value()"
357  << ": domainValue = " << domainValue
358  << ", i = " << i
359  << ", m_referenceDomainValues[i] = " << m_referenceDomainValues[i]
360  << ", m_referenceImageValues[i] = " << m_referenceImageValues[i]
361  << ", m_rateValues[i] = " << m_rateValues[i]
362  << ", imageValue = " << imageValue
363  << std::endl;
364  }
365 
366  return imageValue;
367 }
368 
369 double
370 PiecewiseLinear1D1DFunction::deriv(double domainValue) const
371 {
372  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
373  std::cerr << "In PiecewiseLinear1D1DFunction::deriv()"
374  << ": requested x (" << domainValue
375  << ") is out of the interval (" << m_minDomainValue
376  << ", " << m_maxDomainValue
377  << ")"
378  << std::endl;
379  }
380 
381  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
383  "PiecewiseLinear1D1DFunction::deriv()",
384  "x out of range");
385 
386  unsigned int i = 0;
387  if (m_numRefValues == 1) {
388  // Nothing else to do
389  }
390  else {
391  bool referenceDomainValueFound = false;
392  while (referenceDomainValueFound == false) {
393  if (domainValue < m_referenceDomainValues[i+1]) {
394  referenceDomainValueFound = true;
395  }
396  else {
397  ++i;
400  "PiecewiseLinear1D1DFunction::deriv()",
401  "too big 'i'");
402  }
403  }
404  }
405 
406  return m_rateValues[i];
407 }
408 
409 //*****************************************************
410 // Quadratic 1D->1D class
411 //*****************************************************
413  double minDomainValue,
414  double maxDomainValue,
415  double a,
416  double b,
417  double c)
418  :
419  Base1D1DFunction(minDomainValue,maxDomainValue),
420  m_a (a),
421  m_b (b),
422  m_c (c)
423 {
424 }
425 
427 {
428 }
429 
430 double
431 Quadratic1D1DFunction::value(double domainValue) const
432 {
433  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
434  std::cerr << "In Quadratic1D1DFunction::value()"
435  << ": requested x (" << domainValue
436  << ") is out of the interval (" << m_minDomainValue
437  << ", " << m_maxDomainValue
438  << ")"
439  << std::endl;
440  }
441 
442  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
444  "Quadratic1D1DFunction::value()",
445  "x out of range");
446 
447  double imageValue = m_a*domainValue*domainValue + m_b*domainValue + m_c;
448 
449  return imageValue;
450 }
451 
452 double
453 Quadratic1D1DFunction::deriv(double domainValue) const
454 {
455  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
456  std::cerr << "In Quadratic1D1DFunction::deriv()"
457  << ": requested x (" << domainValue
458  << ") is out of the interval (" << m_minDomainValue
459  << ", " << m_maxDomainValue
460  << ")"
461  << std::endl;
462  }
463 
464  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
466  "Quadratic1D1DFunction::deriv()",
467  "x out of range");
468 
469  return 2.*m_a*domainValue + m_b;
470 }
471 
472 //*****************************************************
473 // Sampled 1D->1D class
474 //*****************************************************
476  :
477  Base1D1DFunction(-INFINITY,INFINITY)
478 {
479  //UQ_FATAL_TEST_MACRO(true,
480  // UQ_UNAVAILABLE_RANK,
481  // "SampledD1DFunction::deriv()",
482  // "invalid constructor");
483 }
484 
486  const std::vector<double>& domainValues,
487  const std::vector<double>& imageValues)
488  :
489  Base1D1DFunction(domainValues[0],domainValues[domainValues.size()-1]),
490  m_domainValues (domainValues.size(),0.),
491  m_imageValues (imageValues.size(), 0.)
492 {
493  unsigned int tmpSize = m_domainValues.size();
494  for (unsigned int i = 0; i < tmpSize; ++i) {
495  m_domainValues[i] = domainValues[i];
496  m_imageValues [i] = imageValues [i];
497  }
498 }
499 
501 {
502 }
503 
504 double
505 Sampled1D1DFunction::value(double domainValue) const
506 {
507  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
508  std::cerr << "In Sampled1D1DFunction::value()"
509  << ": requested x (" << domainValue
510  << ") is out of the interval (" << m_minDomainValue
511  << ", " << m_maxDomainValue
512  << ")"
513  << std::endl;
514  }
515 
516  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
518  "Sampled1D1DFunction::value()",
519  "x out of range");
520 
521  double returnValue = 0.;
522 
523  unsigned int tmpSize = m_domainValues.size();
524  //std::cout << "In Sampled1D1DFunction::value()"
525  // << ": domainValue = " << domainValue
526  // << ", tmpSize = " << tmpSize
527  // << ", m_domainValues[0] = " << m_domainValues[0]
528  // << ", m_domainValues[max] = " << m_domainValues[tmpSize-1]
529  // << std::endl;
530 
531  UQ_FATAL_TEST_MACRO(tmpSize == 0,
533  "Sampled1D1DFunction::value()",
534  "m_domainValues.size() = 0");
535 
536  UQ_FATAL_TEST_MACRO(domainValue < m_domainValues[0],
538  "Sampled1D1DFunction::value()",
539  "domainValue < m_domainValues[0]");
540 
541  UQ_FATAL_TEST_MACRO(m_domainValues[tmpSize-1] < domainValue,
543  "Sampled1D1DFunction::value()",
544  "m_domainValues[max] < domainValue");
545 
546  unsigned int i = 0;
547  for (i = 0; i < tmpSize; ++i) {
548  if (domainValue <= m_domainValues[i]) break;
549  }
550 
551  if (domainValue == m_domainValues[i]) {
552  //if (domainValueWasMatchedExactly) *domainValueWasMatchedExactly = true;
553  returnValue = m_imageValues[i];
554  }
555  else {
556  //if (domainValueWasMatchedExactly) *domainValueWasMatchedExactly = false;
557  double ratio = (domainValue - m_domainValues[i-1])/(m_domainValues[i]-m_domainValues[i-1]);
558  returnValue = m_imageValues[i-1] + ratio * (m_imageValues[i]-m_imageValues[i-1]);
559  }
560 
561  return returnValue;
562 }
563 
564 double
565 Sampled1D1DFunction::deriv(double domainValue) const
566 {
567  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
568  std::cerr << "In Sampled1D1DFunction::deriv()"
569  << ": requested x (" << domainValue
570  << ") is out of the interval (" << m_minDomainValue
571  << ", " << m_maxDomainValue
572  << ")"
573  << std::endl;
574  }
575 
576  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
578  "Sampled1D1DFunction::deriv()",
579  "x out of range");
580 
581  UQ_FATAL_TEST_MACRO(true,
583  "Sampled1d1DFunction::deriv()",
584  "this function makes no sense for this class");
585  return 0.;
586 }
587 
588 const std::vector<double>&
590 {
591  return m_domainValues;
592 }
593 
594 const std::vector<double>&
596 {
597  return m_imageValues;
598 }
599 
600 bool
602 {
603  bool result = false;
604 
605  unsigned int tmpSize = m_domainValues.size();
606  for (unsigned int i = 0; i < tmpSize; ++i) {
607  if (domainValue <= m_domainValues[i]) {
608  result = (domainValue == m_domainValues[i]);
609  break;
610  }
611  }
612 
613  return result;
614 }
615 
616 void
618  const std::vector<double>& domainValues,
619  const std::vector<double>& imageValues)
620 {
621  m_domainValues.clear();
622  m_imageValues.clear();
623 
624  unsigned int tmpSize = domainValues.size();
625  m_minDomainValue = domainValues[0];
626  m_maxDomainValue = domainValues[tmpSize-1];
627 
628  m_domainValues.resize(tmpSize,0.);
629  m_imageValues.resize(tmpSize,0.);
630  for (unsigned int i = 0; i < tmpSize; ++i) {
631  m_domainValues[i] = domainValues[i];
632  m_imageValues [i] = imageValues [i];
633  }
634 
635  return;
636 }
637 
638 void
640  const BaseEnvironment& env,
641  std::ofstream& ofsvar,
642  const std::string& prefixName) const
643 {
644  unsigned int tmpSize = m_domainValues.size();
645  if (tmpSize == 0) {
646  tmpSize = 1;
647  ofsvar << "\n" << prefixName << "Time_sub" << env.subIdString() << " = zeros(" << tmpSize << ",1);"
648  << "\n" << prefixName << "Value_sub" << env.subIdString() << " = zeros(" << tmpSize << ",1);";
649  }
650  else {
651  ofsvar << "\n" << prefixName << "Time_sub" << env.subIdString() << " = zeros(" << tmpSize << ",1);"
652  << "\n" << prefixName << "Value_sub" << env.subIdString() << " = zeros(" << tmpSize << ",1);";
653  for (unsigned int i = 0; i < tmpSize; ++i) {
654  ofsvar << "\n" << prefixName << "Time_sub" << env.subIdString() << "(" << i+1 << ",1) = " << m_domainValues[i] << ";"
655  << "\n" << prefixName << "Value_sub" << env.subIdString() << "(" << i+1 << ",1) = " << m_imageValues[i] << ";";
656  }
657  }
658 
659  return;
660 }
661 
662 //*****************************************************
663 // 'ScalarTimesFunc' 1D->1D class
664 //*****************************************************
666  double scalar,
667  const Base1D1DFunction& func)
668  :
669  Base1D1DFunction(func.minDomainValue(),func.maxDomainValue()),
670  m_scalar(scalar),
671  m_func (func)
672 {
673 }
674 
676 {
677 }
678 
679 double
680 ScalarTimesFunc1D1DFunction::value(double domainValue) const
681 {
682  double value = 0.;
683 
684  value += m_scalar*m_func.value(domainValue);
685 
686  return value;
687 }
688 
689 double
690 ScalarTimesFunc1D1DFunction::deriv(double domainValue) const
691 {
692  double value = 0.;
693 
694  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
695  std::cerr << "In ScalarTimes1D1DFunction::deriv()"
696  << ": requested x (" << domainValue
697  << ") is out of the interval (" << m_minDomainValue
698  << ", " << m_maxDomainValue
699  << ")"
700  << std::endl;
701  }
702 
703  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
705  "ScalarTimes1D1DFunction::deriv()",
706  "x out of range");
707 
708  UQ_FATAL_TEST_MACRO(true,
710  "ScalarTimesFunc1D1DFunction::deriv()",
711  "not implemented yet");
712 
713  return value;
714 }
715 
716 //*****************************************************
717 // 'FuncTimesFunc' 1D->1D class
718 //*****************************************************
720  const Base1D1DFunction& func1,
721  const Base1D1DFunction& func2)
722  :
723  Base1D1DFunction(std::max(func1.minDomainValue(),func2.minDomainValue()),std::min(func1.maxDomainValue(),func2.maxDomainValue())),
724  m_func1(func1),
725  m_func2(func2)
726 {
727 }
728 
730 {
731 }
732 
733 double
734 FuncTimesFunc1D1DFunction::value(double domainValue) const
735 {
736  double value = 0.;
737 
738  value += m_func1.value(domainValue)*m_func2.value(domainValue);
739 
740  return value;
741 }
742 
743 double
744 FuncTimesFunc1D1DFunction::deriv(double domainValue) const
745 {
746  double value = 0.;
747 
748  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
749  std::cerr << "In FuncTimes1D1DFunction::deriv()"
750  << ": requested x (" << domainValue
751  << ") is out of the interval (" << m_minDomainValue
752  << ", " << m_maxDomainValue
753  << ")"
754  << std::endl;
755  }
756 
757  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
759  "FuncTimes1D1DFunction::deriv()",
760  "x out of range");
761 
762  UQ_FATAL_TEST_MACRO(true,
764  "FuncTimesFunc1D1DFunction::deriv()",
765  "not implemented yet");
766 
767  return value;
768 }
769 
770 //*****************************************************
771 // 'FuncPlusFunc' 1D->1D class
772 //*****************************************************
774  const Base1D1DFunction& func1,
775  const Base1D1DFunction& func2)
776  :
777  Base1D1DFunction(std::max(func1.minDomainValue(),func2.minDomainValue()),std::min(func1.maxDomainValue(),func2.maxDomainValue())),
778  m_func1(func1),
779  m_func2(func2)
780 {
781 }
782 
784 {
785 }
786 
787 double
788 FuncPlusFunc1D1DFunction::value(double domainValue) const
789 {
790  double value = 0.;
791 
792  value += m_func1.value(domainValue) + m_func2.value(domainValue);
793 
794  return value;
795 }
796 
797 double
798 FuncPlusFunc1D1DFunction::deriv(double domainValue) const
799 {
800  double value = 0.;
801 
802  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
803  std::cerr << "In FuncPlus1D1DFunction::deriv()"
804  << ": requested x (" << domainValue
805  << ") is out of the interval (" << m_minDomainValue
806  << ", " << m_maxDomainValue
807  << ")"
808  << std::endl;
809  }
810 
811  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
813  "FuncPlus1D1DFunction::deriv()",
814  "x out of range");
815 
816  UQ_FATAL_TEST_MACRO(true,
818  "FuncPlusFunc1D1DFunction::deriv()",
819  "not implemented yet");
820 
821  return value;
822 }
823 
824 //*****************************************************
825 // Lagrange Polynomial 1D->1D class
826 //*****************************************************
828  const std::vector<double>& positionValues,
829  const std::vector<double>* functionValues)
830  :
831  Base1D1DFunction(-INFINITY,INFINITY),
832  m_positionValues(positionValues),
833  m_functionValues(positionValues.size(),1.)
834 {
835  if (functionValues) {
836  UQ_FATAL_TEST_MACRO(m_positionValues.size() != functionValues->size(),
838  "LagrangePolynomial1D1DFunction::constructor()",
839  "invalid input");
840  m_functionValues = *functionValues;
841  }
842 }
843 
845 {
846 }
847 
848 double
849 LagrangePolynomial1D1DFunction::value(double domainValue) const
850 {
851  double value = 0.;
852 
853  for (unsigned int k = 0; k < m_positionValues.size(); ++k) {
854  double scaleFactor = 1.;
855  double posK = m_positionValues[k];
856  for (unsigned int j = 0; j < m_positionValues.size(); ++j) {
857  if (j != k) {
858  double posJ = m_positionValues[j];
859  scaleFactor *= (domainValue-posJ)/(posK-posJ);
860  }
861  }
862 
863  //if (m_env.subDisplayFile()) {
864  // *m_env.subDisplayFile() << "In sddTG<K_V,K_M>::lagrange()"
865  // << ": k = " << k
866  // << ", scaleFactor = " << scaleFactor
867  // << std::endl;
868  //}
869 
870  value += scaleFactor * m_functionValues[k];
871 
872  //if (m_env.subDisplayFile()) {
873  // *m_env.subDisplayFile() << "In sddTG<K_V,K_M>::lagrange()"
874  // << ": k = " << k
875  // << ", value = " << value
876  // << std::endl;
877  //}
878  }
879 
880  return value;
881 }
882 
883 double
884 LagrangePolynomial1D1DFunction::deriv(double domainValue) const
885 {
886  double value = 0.;
887 
888  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
889  std::cerr << "In LagrangePolynomial1D1DFunction::deriv()"
890  << ": requested x (" << domainValue
891  << ") is out of the interval (" << m_minDomainValue
892  << ", " << m_maxDomainValue
893  << ")"
894  << std::endl;
895  }
896 
897  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
899  "LagrangePolynomial1D1DFunction::deriv()",
900  "x out of range");
901 
902  UQ_FATAL_TEST_MACRO(true,
904  "LagrangePolynomial1D1DFunction::deriv()",
905  "not implemented yet");
906 
907  return value;
908 }
909 
910 //*****************************************************
911 // Lagrange Basis 1D->1D class
912 //*****************************************************
914  const std::vector<double>& positionValues,
915  unsigned int basisIndex)
916  :
917  Base1D1DFunction(-INFINITY,INFINITY),
918  m_positionValues(positionValues),
919  m_basisIndex (basisIndex)
920 {
923  "LagrangeBasis1D1DFunction::constructor()",
924  "invalid input");
925 }
926 
928 {
929 }
930 
931 double
932 LagrangeBasis1D1DFunction::value(double domainValue) const
933 {
934  double scaleFactor = 1.;
935 
936  unsigned int k = m_basisIndex;
937  double posK = m_positionValues[k];
938  for (unsigned int j = 0; j < m_positionValues.size(); ++j) {
939  if (j != k) {
940  double posJ = m_positionValues[j];
941  scaleFactor *= (domainValue-posJ)/(posK-posJ);
942  }
943  }
944 
945  return scaleFactor;
946 }
947 
948 double
949 LagrangeBasis1D1DFunction::deriv(double domainValue) const
950 {
951  double value = 0.;
952 
953  if ((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)) {
954  std::cerr << "In LagrangeBasis1D1DFunction::deriv()"
955  << ": requested x (" << domainValue
956  << ") is out of the interval (" << m_minDomainValue
957  << ", " << m_maxDomainValue
958  << ")"
959  << std::endl;
960  }
961 
962  UQ_FATAL_TEST_MACRO(((domainValue < m_minDomainValue) || (domainValue > m_maxDomainValue)),
964  "LagrangeBasis1D1DFunction::deriv()",
965  "x out of range");
966 
967  UQ_FATAL_TEST_MACRO(true,
969  "LagrangeBasis1D1DFunction::deriv()",
970  "not implemented yet");
971 
972  return value;
973 }
974 
975 template <class T>
976 double
978  const ScalarSequence<T>& scalarSeq2,
979  unsigned int initialPos,
980  double scaleValue1,
981  double scaleValue2,
982  const Base1D1DFunction& func1,
983  const Base1D1DFunction& func2,
984  unsigned int quadratureOrder)
985 {
986  double resultValue = 0.;
987 
988  UQ_FATAL_TEST_MACRO(initialPos != 0,
989  scalarSeq1.env().worldRank(),
990  "SubF1F2Gaussian2dKdeIntegral()",
991  "not implemented yet for initialPos != 0");
992  UQ_FATAL_TEST_MACRO(scalarSeq1.subSequenceSize() != scalarSeq2.subSequenceSize(),
993  scalarSeq1.env().worldRank(),
994  "SubF1F2Gaussian2dKdeIntegral()",
995  "different sizes");
996 
997  GaussianHermite1DQuadrature quadObj(0.,1.,quadratureOrder);
998  const std::vector<double>& quadPositions = quadObj.positions();
999  const std::vector<double>& quadWeights = quadObj.weights ();
1000  UQ_FATAL_TEST_MACRO(quadPositions.size() != quadWeights.size(),
1002  "SubF1F2Gaussian2dKdeIntegral()",
1003  "quadObj has invalid state");
1004 
1005  unsigned int numQuadraturePositions = quadPositions.size();
1006  unsigned int dataSize = scalarSeq1.subSequenceSize();
1007  for (unsigned int k = 0; k < dataSize; ++k) {
1008  double value1 = 0.;
1009  double value2 = 0.;
1010  double x1k = scalarSeq1[k];
1011  double x2k = scalarSeq2[k];
1012  for (unsigned int j = 0; j < numQuadraturePositions; ++j) {
1013  value1 += func1.value(scaleValue1*quadPositions[j]+x1k)*quadWeights[j];
1014  value2 += func2.value(scaleValue2*quadPositions[j]+x2k)*quadWeights[j];
1015  }
1016  resultValue += value1*value2;
1017  }
1018  resultValue *= 1./(2.*M_PI)/((double) dataSize);
1019 
1020  return resultValue;
1021 }
1022 
1023 } // End namespace QUESO
Quadratic1D1DFunction(double minDomainValue, double maxDomainValue, double a, double b, double c)
Default constructor.
Definition: 1D1DFunction.C:412
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
double value(double domainValue) const
Returns the value of the piecewise-linear function at point domainValue.
Definition: 1D1DFunction.C:320
bool domainValueMatchesExactly(double domainValue) const
Checks whether the domain value domainValue matches exactly one of the values in the function domain ...
Definition: 1D1DFunction.C:601
double value(double domainValue) const
Returns the value of the multiplication of a function func1 by another function func2 at the point do...
Definition: 1D1DFunction.C:734
virtual double value(double domainValue) const =0
Returns the value of the (one-dimensional) function. See template specialization. ...
virtual double multiplyAndIntegrate(const Base1D1DFunction &func, unsigned int quadratureOrder, double *resultWithMultiplicationByTAsWell) const
TODO: Multiplies this function with function, and integrates it numerically. See template specializat...
Definition: 1D1DFunction.C:63
Class for Hermite-Gauss quadrature rule for one-dimensional functions.
Definition: 1DQuadrature.h:219
~Quadratic1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:426
Generic1D1DFunction(double minDomainValue, double maxDomainValue, double(*valueRoutinePtr)(double domainValue, const void *routinesDataPtr), double(*derivRoutinePtr)(double domainValue, const void *routinesDataPtr), const void *routinesDataPtr)
Default constructor.
Definition: 1D1DFunction.C:82
double value(double domainValue) const
Returns the value of the Lagrange basis at point domainValue.
Definition: 1D1DFunction.C:932
virtual double deriv(double domainValue) const =0
Returns the value of the derivative of the function. See template specialization. ...
PiecewiseLinear1D1DFunction(double minDomainValue, double maxDomainValue, const std::vector< double > &referenceDomainValues, double referenceImageValue0, const std::vector< double > &rateValues)
Default constructor.
Definition: 1D1DFunction.C:263
~FuncTimesFunc1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:729
std::vector< double > m_rateValues
Rate value; in , for each =1 .. m_numRefValues.
Definition: 1D1DFunction.h:277
Class for handling scalar samples.
void set(const std::vector< double > &domainValues, const std::vector< double > &imageValues)
Sets the values of the independent (domainValues) and dependent (imageValues) variables of this sampl...
Definition: 1D1DFunction.C:617
double value(double domainValue) const
Returns the value of the (one-dimensional) function at point domainValue.
Definition: 1D1DFunction.C:101
double value(double domainValue) const
Returns the value of the linear function at point domainValue.
Definition: 1D1DFunction.C:219
~Generic1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:96
~LagrangeBasis1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:927
int worldRank() const
Returns the process world rank.
Definition: Environment.C:235
Base1D1DFunction(double minDomainValue, double maxDomainValue)
Default constructor.
Definition: 1D1DFunction.C:33
std::vector< double > m_positionValues
Definition: 1D1DFunction.h:584
double maxDomainValue() const
Returns the maximum value of the domain of the (one-dimensional) function.
Definition: 1D1DFunction.C:57
const int UQ_UNAVAILABLE_RANK
Definition: Defines.h:74
const Base1D1DFunction & m_func
Definition: 1D1DFunction.h:449
double(* m_derivRoutinePtr)(double domainValue, const void *routinesDataPtr)
Definition: 1D1DFunction.h:128
double deriv(double domainValue) const
TODO: Returns the value of the derivative of the function multiplied by the given scalar at point dom...
Definition: 1D1DFunction.C:690
Linear1D1DFunction(double minDomainValue, double maxDomainValue, double referenceDomainValue, double referenceImageValue, double rateValue)
Default constructor.
Definition: 1D1DFunction.C:200
This (virtual) class sets up the environment underlying the use of the QUESO library by an executable...
Definition: Environment.h:187
double m_referenceDomainValue
Reference value in the function domain; in .
Definition: 1D1DFunction.h:220
double value(double domainValue) const
Returns the value of the constant function at point domainValue.
Definition: 1D1DFunction.C:158
~FuncPlusFunc1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:783
~Constant1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:153
double value(double domainValue) const
Returns the value of the Lagrange polynomial at point domainValue.
Definition: 1D1DFunction.C:849
ScalarTimesFunc1D1DFunction(double scalar, const Base1D1DFunction &func)
Default constructor.
Definition: 1D1DFunction.C:665
unsigned int m_numRefValues
Number of points which will be evaluated.
Definition: 1D1DFunction.h:268
double m_rateValue
Rate value; in .
Definition: 1D1DFunction.h:226
double deriv(double domainValue) const
TODO: Returns the value of the derivative of the addition of two functions.
Definition: 1D1DFunction.C:798
double deriv(double domainValue) const
Returns the value of the derivative of the function at point domainValue.
Definition: 1D1DFunction.C:121
double deriv(double domainValue) const
TODO: Returns the value of the derivative of the Lagrange polynomial at point domainValue.
Definition: 1D1DFunction.C:884
Sampled1D1DFunction()
Default constructor. It should not be called by the user.
Definition: 1D1DFunction.C:475
double SubF1F2Gaussian2dKdeIntegral(const ScalarSequence< T > &scalarSeq1, const ScalarSequence< T > &scalarSeq2, unsigned int initialPos, double scaleValue1, double scaleValue2, const Base1D1DFunction &func1, const Base1D1DFunction &func2, unsigned int quadratureOrder)
Calculates the integral of a 2D Gaussian KDE.
Definition: 1D1DFunction.C:977
double(* m_valueRoutinePtr)(double domainValue, const void *routinesDataPtr)
Definition: 1D1DFunction.h:127
const std::vector< double > & domainValues() const
Array of the domain values (values of the independent variable).
Definition: 1D1DFunction.C:589
const BaseEnvironment & env() const
Access to QUESO environment.
const void * m_routinesDataPtr
Definition: 1D1DFunction.h:129
const Base1D1DFunction & m_func2
Definition: 1D1DFunction.h:492
double value(double domainValue) const
Returns the value of the addition of function func1 and func2 evaluated at the point domainValue...
Definition: 1D1DFunction.C:788
std::vector< double > m_positionValues
Definition: 1D1DFunction.h:635
double deriv(double domainValue) const
TODO: Returns the value of the derivative of the function func1 by another function func2 at the poin...
Definition: 1D1DFunction.C:744
virtual double value(double domainValue) const
Returns the value of the sampled function at point domainValue.
Definition: 1D1DFunction.C:505
FuncPlusFunc1D1DFunction(const Base1D1DFunction &func1, const Base1D1DFunction &func2)
Default constructor.
Definition: 1D1DFunction.C:773
const Base1D1DFunction & m_func1
Definition: 1D1DFunction.h:531
double deriv(double domainValue) const
Returns the value of the derivative of the constant function at point domainValue.
Definition: 1D1DFunction.C:178
double m_referenceImageValue
Reference value in the function image; in .
Definition: 1D1DFunction.h:223
LagrangePolynomial1D1DFunction(const std::vector< double > &positionValues, const std::vector< double > *functionValues)
Default constructor.
Definition: 1D1DFunction.C:827
virtual ~Base1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:46
double value(double domainValue) const
Returns the value of the (one-dimensional) function at point domainValue, already multiplied by the s...
Definition: 1D1DFunction.C:680
double deriv(double domainValue) const
TODO: Returns the value of the derivative of the Lagrange basis at point domainValue.
Definition: 1D1DFunction.C:949
~Linear1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:214
const std::vector< double > & imageValues() const
Array of the image values (values of the dependent variable).
Definition: 1D1DFunction.C:595
double deriv(double domainValue) const
Bogus: Derivative of the function.
Definition: 1D1DFunction.C:565
Constant1D1DFunction(double minDomainValue, double maxDomainValue, double constantValue)
Default constructor.
Definition: 1D1DFunction.C:143
unsigned int subSequenceSize() const
Size of the sub-sequence of scalars.
virtual void printForMatlab(const BaseEnvironment &env, std::ofstream &ofsvar, const std::string &prefixName) const
Prints the values of the function in Matlab/Octave format.
Definition: 1D1DFunction.C:639
std::vector< double > m_domainValues
Array of the values in the domain of the function (values of the independent variable).
Definition: 1D1DFunction.h:404
std::vector< double > m_imageValues
Array of the values in the image of the function (values of the dependent variable).
Definition: 1D1DFunction.h:407
const Base1D1DFunction & m_func2
Definition: 1D1DFunction.h:532
double deriv(double domainValue) const
Returns the value of the derivative of the linear function at point domainValue.
Definition: 1D1DFunction.C:241
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
std::vector< double > m_referenceDomainValues
Reference values in the piecewise-linear function domain; in , for each =1 .. m_numRefValues.
Definition: 1D1DFunction.h:271
double deriv(double domainValue) const
Returns the value of the derivative of the piecewise-linear function at point domainValue.
Definition: 1D1DFunction.C:370
LagrangeBasis1D1DFunction(const std::vector< double > &positionValues, unsigned int basisIndex)
Default constructor.
Definition: 1D1DFunction.C:913
virtual ~Sampled1D1DFunction()
Destructor.
Definition: 1D1DFunction.C:500
const Base1D1DFunction & m_func1
Definition: 1D1DFunction.h:491
Class for one-dimensional functions.
Definition: 1D1DFunction.h:53
double minDomainValue() const
Returns the minimum value of the domain of the (one-dimensional) function.
Definition: 1D1DFunction.C:51
std::vector< double > m_referenceImageValues
Reference values in the piecewise-linear function image; in , for each =1 .. m_numRefValues.
Definition: 1D1DFunction.h:274
FuncTimesFunc1D1DFunction(const Base1D1DFunction &func1, const Base1D1DFunction &func2)
Default constructor.
Definition: 1D1DFunction.C:719
std::vector< double > m_functionValues
Definition: 1D1DFunction.h:585

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