Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
MCGIDI_functions.hpp
Go to the documentation of this file.
1/*
2# <<BEGIN-copyright>>
3# Copyright 2019, Lawrence Livermore National Security, LLC.
4# This file is part of the gidiplus package (https://github.com/LLNL/gidiplus).
5# gidiplus is licensed under the MIT license (see https://opensource.org/licenses/MIT).
6# SPDX-License-Identifier: MIT
7# <<END-copyright>>
8*/
9
10#ifndef MCGIDI_functions_hpp_included
11#define MCGIDI_functions_hpp_included 1
12
13#include <nf_utilities.h>
14#include <ptwXY.h>
15#include <LUPI_dataBuffer.hpp>
16
17namespace MCGIDI {
18
21enum class Function2dType { none, XYs };
25
27
28namespace Functions {
29
30/*
31============================================================
32====================== FunctionBase ========================
33============================================================
34*/
36
37 private:
38 int m_dimension;
39 double m_domainMin;
40 double m_domainMax;
41 Interpolation m_interpolation;
42 double m_outerDomainValue;
43
44 public:
47 LUPI_HOST_DEVICE FunctionBase( int a_dimension, double a_domainMin, double a_domainMax, Interpolation a_interpolation, double a_outerDomainValue = 0 );
48 LUPI_HOST_DEVICE virtual ~FunctionBase( ) = 0;
49
50 LUPI_HOST_DEVICE Interpolation interpolation( ) const { return( m_interpolation ); }
51 LUPI_HOST_DEVICE double domainMin( ) const { return( m_domainMin ); }
52 LUPI_HOST_DEVICE double domainMax( ) const { return( m_domainMax ); }
53 LUPI_HOST_DEVICE double outerDomainValue( ) const { return( m_outerDomainValue ); }
55};
56
57/*
58============================================================
59======================== Function1d ========================
60============================================================
61*/
62class Function1d : public FunctionBase {
63
64 protected:
66
67 public:
69 LUPI_HOST_DEVICE Function1d( double a_domainMin, double a_domainMax, Interpolation a_interpolation, double a_outerDomainValue = 0 );
71
74
75 template <typename RNG>
76 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION int sampleBoundingInteger( double a_x1, RNG && a_rng ) const ;
79};
80
81/*
82============================================================
83====================== Function1d_d1 =======================
84============================================================
85*/
86class Function1d_d1 : public Function1d {
87
88 public:
91 LUPI_HOST_DEVICE Function1d_d1( double a_domainMin, double a_domainMax, Interpolation a_interpolation, double a_outerDomainValue = 0 ) :
92 Function1d( a_domainMin, a_domainMax, a_interpolation, a_outerDomainValue ) { }
93
94 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
95};
96
97/*
98============================================================
99====================== Function1d_d2 =======================
100============================================================
101*/
103
104 public:
107 LUPI_HOST_DEVICE Function1d_d2( double a_domainMin, double a_domainMax, Interpolation a_interpolation, double a_outerDomainValue = 0 ) :
108 Function1d_d1( a_domainMin, a_domainMax, a_interpolation, a_outerDomainValue ) { }
109
110 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
111};
112
113/*
114============================================================
115======================== Constant1d ========================
116============================================================
117*/
118class Constant1d : public Function1d_d2 {
119
120 private:
121 double m_value;
122
123 public:
125 LUPI_HOST_DEVICE Constant1d( double a_domainMin, double a_domainMax, double a_value, double a_outerDomainValue = 0 );
126 LUPI_HOST Constant1d( GIDI::Functions::Constant1d const &a_constant1d );
128
129 LUPI_HOST_DEVICE double evaluate( LUPI_maybeUnused double a_x1 ) const { return( m_value ); }
131};
132
133/*
134============================================================
135=========================== XYs1d ==========================
136============================================================
137*/
138class XYs1d : public Function1d_d2 {
139
140 private:
141 Vector<double> m_Xs;
142 Vector<double> m_Ys;
143
144 public:
146 LUPI_HOST XYs1d( Interpolation a_interpolation, Vector<double> a_Xs, Vector<double> a_Ys, double a_outerDomainValue = 0 );
147 LUPI_HOST XYs1d( GIDI::Functions::XYs1d const &a_XYs1d );
149
150 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
152};
153
154/*
155============================================================
156======================= Polynomial1d =======================
157============================================================
158*/
160
161 private:
162 Vector<double> m_coefficients;
163 Vector<double> m_coefficientsReversed;
164
165 public:
167 LUPI_HOST Polynomial1d( double a_domainMin, double a_domainMax, Vector<double> const &a_coefficients, double a_outerDomainValue = 0 );
170
171 LUPI_HOST_DEVICE Vector<double> const &coefficients( ) const { return( m_coefficients ); }
172 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
174};
175
176/*
177============================================================
178========================= Gridded1d ========================
179============================================================
180*/
181class Gridded1d : public Function1d_d2 {
182
183 private:
184 Vector<double> m_grid;
185 Vector<double> m_data;
186
187 public:
189 LUPI_HOST Gridded1d( GIDI::Functions::Gridded1d const &a_gridded1d );
191
192 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
194};
195
196/*
197============================================================
198========================= Regions1d ========================
199============================================================
200*/
201class Regions1d : public Function1d_d1 {
202
203 private:
204 Vector<double> m_Xs;
205 Vector<Function1d_d2 *> m_functions1d;
206
207 public:
209 LUPI_HOST Regions1d( GIDI::Functions::Regions1d const &a_regions1d );
211
212 LUPI_HOST_DEVICE void append( Function1d_d2 *a_function1d );
213 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
215};
216
217/*
218============================================================
219======================== Branching1d =======================
220============================================================
221*/
223
224 private:
225 int m_initialStateIndex;
226
227 public:
229 LUPI_HOST Branching1d( SetupInfo &a_setupInfo, GIDI::Functions::Branching1d const &a_branching1d );
231
232 LUPI_HOST_DEVICE int initialStateIndex( ) const { return( m_initialStateIndex ); }
233
234 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
236};
237
238/*
239============================================================
240========== TerrellFissionNeutronMultiplicityModel ==========
241============================================================
242*/
244
245 private:
246 double m_width;
247 Function1d_d1 *m_multiplicity;
248
249 public:
251 LUPI_HOST TerrellFissionNeutronMultiplicityModel( double a_width, Function1d_d1 *a_multiplicity );
253
254 template <typename RNG>
255 LUPI_HOST_DEVICE int sampleBoundingInteger( double a_energy, RNG && a_rng ) const ;
256 LUPI_HOST_DEVICE double evaluate( double a_energy ) const ;
258};
259
260/*
261============================================================
262======================== Function2d ========================
263============================================================
264*/
265class Function2d : public FunctionBase {
266
267 protected:
269
270 public:
272 LUPI_HOST Function2d( double a_domainMin, double a_domainMax, Interpolation a_interpolation, double a_outerDomainValue = 0 );
274
277
278 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate( double a_x2, double a_x1 ) const MCGIDI_TRUE_VIRTUAL;
280};
281
282/*
283============================================================
284=========================== XYs2d ==========================
285============================================================
286*/
287class XYs2d : public Function2d {
288
289 private:
290 Vector<double> m_Xs;
291 Vector<Function1d_d1 *> m_functions1d;
292
293 public:
295 LUPI_HOST XYs2d( GIDI::Functions::XYs2d const &a_XYs2d );
297
298 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
300};
301
302/*
303============================================================
304========================== others ==========================
305============================================================
306*/
307LUPI_HOST Function1d *parseMultiplicityFunction1d( SetupInfo &a_setupInfo, Transporting::MC const &a_settings, GIDI::Suite const &a_suite );
313
314} // End of namespace Functions.
315
316/*
317============================================================
318============================================================
319================== namespace Probabilities ==================
320============================================================
321============================================================
322*/
323namespace Probabilities {
324
325/*
326============================================================
327===================== ProbabilityBase ======================
328============================================================
329*/
343
344/*
345============================================================
346===================== ProbabilityBase1d ====================
347============================================================
348*/
350
351 protected:
353
354 public:
358
361
363 template <typename RNG>
364 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double sample( double a_rngValue, RNG && a_rng ) const MCGIDI_TRUE_VIRTUAL;
366};
367
368/*
369============================================================
370======================= Xs_pdf_cdf1d =======================
371============================================================
372*/
374
375 private:
376 Vector<double> m_pdf;
377 Vector<double> m_cdf;
378
379 public:
383
384 LUPI_HOST_DEVICE double evaluate( double a_x1 ) const ;
385 template <typename RNG>
386 LUPI_HOST_DEVICE double sample( double a_rngValue, RNG && a_rng ) const ;
388};
389
390/*
391============================================================
392===================== ProbabilityBase2d ====================
393============================================================
394*/
396
397 protected:
399
400 public:
405
408
409 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate( double a_x2, double a_x1 ) const MCGIDI_TRUE_VIRTUAL;
410 template <typename RNG>
411 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double sample( double a_x2, double a_rngValue, RNG && a_rng ) const MCGIDI_TRUE_VIRTUAL;
413};
414
415/*
416============================================================
417=================== ProbabilityBase2d_d1 ===================
418============================================================
419*/
421
422 public:
426 ProbabilityBase2d( a_probabilty ) { }
428 ProbabilityBase2d( a_probabilty, a_Xs ) { }
429
430 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
431 template <typename RNG>
432 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
433 template <typename RNG>
434 LUPI_HOST_DEVICE double sample2dOf3d( double a_x2, double a_rngValue, RNG && a_rng, double *a_x1_1, double *a_x1_2 ) const ;
435};
436
437/*
438============================================================
439=================== ProbabilityBase2d_d2 ===================
440============================================================
441*/
443
444 public:
450 ProbabilityBase2d_d1( a_probabilty, a_Xs ) { }
451
452 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
453 template <typename RNG>
454 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
455 template <typename RNG>
456 LUPI_HOST_DEVICE double sample2dOf3d( double a_x2, double a_rngValue, RNG && a_rng, double *a_x1_1, double *a_x1_2 ) const ;
457};
458
459/*
460============================================================
461========================== XYs2d ===========================
462============================================================
463*/
465
466 private:
467 Vector<ProbabilityBase1d *> m_probabilities;
468
469 public:
471 LUPI_HOST XYs2d( GIDI::Functions::XYs2d const &a_XYs2d );
473
474 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
475 template <typename RNG>
476 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
477 template <typename RNG>
478 LUPI_HOST_DEVICE double sample2dOf3d( double a_x2, double a_rngValue, RNG && a_rng, double *a_x1_1, double *a_x1_2 ) const ;
480};
481
482/*
483============================================================
484======================== Regions2d =========================
485============================================================
486*/
488
489 private:
490 Vector<ProbabilityBase2d_d2 *> m_probabilities;
491
492 public:
494 LUPI_HOST Regions2d( GIDI::Functions::Regions2d const &a_regions2d );
496
497 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
498 template <typename RNG>
499 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
501};
502
503/*
504============================================================
505======================== Isotropic2d =======================
506============================================================
507*/
509
510 public:
512 LUPI_HOST Isotropic2d( GIDI::Functions::Isotropic2d const &a_isotropic2d );
514
515 LUPI_HOST_DEVICE double evaluate( LUPI_maybeUnused double a_x2, LUPI_maybeUnused double a_x1 ) const { return( 0.5 ); }
516 template <typename RNG>
517 LUPI_HOST_DEVICE double sample( LUPI_maybeUnused double a_x2, double a_rngValue, LUPI_maybeUnused RNG && a_rng ) const { return( 1. - 2. * a_rngValue ); }
520};
521
522/*
523============================================================
524====================== DiscreteGamma2d =====================
525============================================================
526*/
528
529 private:
530 double m_value;
531
532 public:
536
537 LUPI_HOST_DEVICE double evaluate( LUPI_maybeUnused double a_x2, LUPI_maybeUnused double a_x1 ) const { return( m_value ); } // FIXME This is wrong, should be something like 1 when domainMin <= a_x1 <= domainMax ), I think. I.e., should be a probability.
538 template <typename RNG>
539 LUPI_HOST_DEVICE double sample( LUPI_maybeUnused double a_x2, LUPI_maybeUnused double a_rngValue, LUPI_maybeUnused RNG && a_rng ) const { return( m_value ); }
541};
542
543/*
544============================================================
545====================== PrimaryGamma2d =====================
546============================================================
547*/
549
550 private:
551 double m_primaryEnergy;
552 double m_massFactor;
553 String m_finalState;
554 int m_initialStateIndex;
555
556 public:
558 LUPI_HOST PrimaryGamma2d( GIDI::Functions::PrimaryGamma2d const &a_primaryGamma2d, SetupInfo *a_setupInfo );
560
561 double primaryEnergy( ) const { return( m_primaryEnergy ); } /**< Returns the value of the *m_primaryEnergy* member. */
562 double massFactor( ) const { return( m_massFactor ); } /**< Returns the value of the *m_massFactor* member. */
563 String const &finalState( ) const { return( m_finalState ); } /**< Returns a const reference to the *m_finalState* member. */
564 int initialStateIndex( ) const { return( m_initialStateIndex ); } /**< Returns the value of the *m_initialStateIndex* member. */
565
566 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
567 template <typename RNG>
568 LUPI_HOST_DEVICE double sample( double a_x2, LUPI_maybeUnused double a_rngValue, LUPI_maybeUnused RNG && a_rng ) const { return( m_primaryEnergy + a_x2 * m_massFactor ); }
570};
571
572/*
573============================================================
574========================= Recoil2d =========================
575============================================================
576*/
578
579 private:
580 String m_xlink;
581
582 public:
584 LUPI_HOST Recoil2d( GIDI::Functions::Recoil2d const &a_recoil2d );
586
587 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
588 template <typename RNG>
589 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
591};
592
593/*
594============================================================
595==================== NBodyPhaseSpace2d =====================
596============================================================
597*/
599
600 private:
601 int m_numberOfProducts;
602 double m_mass;
603 double m_energy_in_COMFactor;
604 double m_massFactor;
605 double m_Q;
606 ProbabilityBase1d *m_dist;
607
608 public:
610 LUPI_HOST NBodyPhaseSpace2d( GIDI::Functions::NBodyPhaseSpace2d const &a_NBodyPhaseSpace2d, SetupInfo *a_setupInfo );
612
613 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
614 template <typename RNG>
615 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
617};
618
619/*
620============================================================
621====================== Evaporation2d =======================
622============================================================
623*/
625
626 private:
627 double m_U;
629
630 public:
632 LUPI_HOST Evaporation2d( GIDI::Functions::Evaporation2d const &a_generalEvaporation2d );
634
635 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
636 template <typename RNG>
637 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
639};
640
641/*
642============================================================
643=================== GeneralEvaporation2d ===================
644============================================================
645*/
647
648 private:
651
652 public:
656
657 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
658 template <typename RNG>
659 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
661};
662
663/*
664============================================================
665================= SimpleMaxwellianFission2d ================
666============================================================
667*/
669
670 private:
671 double m_U;
673
674 public:
678
679 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
680 template <typename RNG>
681 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
683};
684
685/*
686============================================================
687========================== Watt2d ==========================
688============================================================
689*/
691
692 private:
693 double m_U;
696
697 public:
699 LUPI_HOST Watt2d( GIDI::Functions::Watt2d const &a_Watt2d );
701
702 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
703 template <typename RNG>
704 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
706};
707
708/*
709============================================================
710=================== WeightedFunctionals2d ==================
711============================================================
712*/
714
715 private:
718
719 public:
723
724 LUPI_HOST_DEVICE double evaluate( double a_x2, double a_x1 ) const ;
725 template <typename RNG>
726 LUPI_HOST_DEVICE double sample( double a_x2, double a_rngValue, RNG && a_rng ) const ;
728};
729
730/*
731============================================================
732===================== ProbabilityBase3d ====================
733============================================================
734*/
736
737 protected:
739
740 public:
744
747
748 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate( double a_x3, double a_x2, double a_x1 ) const MCGIDI_TRUE_VIRTUAL;
749 template <typename RNG>
750 LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double sample( double a_x3, double a_x2_1, double a_x2_2, double a_rngValue, RNG && a_rng ) const MCGIDI_TRUE_VIRTUAL;
752};
753
754/*
755============================================================
756========================== XYs3d ===========================
757============================================================
758*/
759class XYs3d : public ProbabilityBase3d {
760
761 private:
762 Vector<ProbabilityBase2d_d1 *> m_probabilities;
763
764 public:
766 LUPI_HOST XYs3d( GIDI::Functions::XYs3d const &a_XYs3d );
768
769 LUPI_HOST_DEVICE double evaluate( double a_x3, double a_x2, double a_x1 ) const ;
770 template <typename RNG>
771 LUPI_HOST_DEVICE double sample( double a_x3, double a_x2_1, double a_x2_2, double a_rngValue, RNG && a_rng ) const ;
773};
774
775/*
776============================================================
777========================== others ==========================
778============================================================
779*/
782LUPI_HOST ProbabilityBase2d *parseProbability2d( Transporting::MC const &a_settings, GIDI::Suite const &a_suite, SetupInfo *a_setupInfo );
788
789
790} // End of namespace Probabilities.
791
792/*
793============================================================
794========================== others ==========================
795============================================================
796*/
798
802 LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d1 *a_function1d );
804 LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d2 *a_function1d );
805
808
811
814 Probabilities::ProbabilityBase2d *a_probability2d );
819
822
823} // End of namespace MCGIDI.
824
825#endif // End of MCGIDI_functions_hpp_included
#define LUPI_HOST_DEVICE
#define LUPI_HOST
#define LUPI_maybeUnused
#define MCGIDI_VIRTUAL_FUNCTION
Definition MCGIDI.hpp:33
#define MCGIDI_TRUE_VIRTUAL
Definition MCGIDI.hpp:34
LUPI_HOST_DEVICE int initialStateIndex() const
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(LUPI_maybeUnused double a_x1) const
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST_DEVICE Function1d_d1(double a_domainMin, double a_domainMax, Interpolation a_interpolation, double a_outerDomainValue=0)
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST_DEVICE Function1d_d2(double a_domainMin, double a_domainMax, Interpolation a_interpolation, double a_outerDomainValue=0)
LUPI_HOST_DEVICE String typeString() const
LUPI_HOST_DEVICE Function1dType type() const
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION int sampleBoundingInteger(double a_x1, RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate(double a_x1) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE Function2dType type() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE String typeString() const
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate(double a_x2, double a_x1) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE double domainMin() const
LUPI_HOST_DEVICE Interpolation interpolation() const
virtual LUPI_HOST_DEVICE ~FunctionBase()=0
LUPI_HOST_DEVICE double domainMax() const
LUPI_HOST_DEVICE double outerDomainValue() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE Vector< double > const & coefficients() const
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void append(Function1d_d2 *a_function1d)
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST_DEVICE double evaluate(double a_energy) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE int sampleBoundingInteger(double a_energy, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double sample(LUPI_maybeUnused double a_x2, LUPI_maybeUnused double a_rngValue, LUPI_maybeUnused RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(LUPI_maybeUnused double a_x2, LUPI_maybeUnused double a_x1) const
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(LUPI_maybeUnused double a_x2, LUPI_maybeUnused double a_x1) const
LUPI_HOST_DEVICE double sample(LUPI_maybeUnused double a_x2, double a_rngValue, LUPI_maybeUnused RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double sample(double a_x2, LUPI_maybeUnused double a_rngValue, LUPI_maybeUnused RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double sample(double a_rngValue, RNG &&a_rng) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate(double a_x1) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE String typeString() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE ProbabilityBase1dType type() const
LUPI_HOST ProbabilityBase2d_d1(GIDI::Functions::FunctionForm const &a_probabilty, Vector< double > const &a_Xs)
LUPI_HOST ProbabilityBase2d_d1(GIDI::Functions::FunctionForm const &a_probabilty)
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE double sample2dOf3d(double a_x2, double a_rngValue, RNG &&a_rng, double *a_x1_1, double *a_x1_2) const
LUPI_HOST ProbabilityBase2d_d2(GIDI::Functions::FunctionForm const &a_probabilty)
LUPI_HOST_DEVICE double sample2dOf3d(double a_x2, double a_rngValue, RNG &&a_rng, double *a_x1_1, double *a_x1_2) const
LUPI_HOST ProbabilityBase2d_d2(GIDI::Functions::FunctionForm const &a_probabilty, Vector< double > const &a_Xs)
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate(double a_x2, double a_x1) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double sample(double a_x2, double a_rngValue, RNG &&a_rng) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE String typeString() const
LUPI_HOST_DEVICE ProbabilityBase2dType type() const
LUPI_HOST_DEVICE ProbabilityBase3dType type() const
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double sample(double a_x3, double a_x2_1, double a_x2_2, double a_rngValue, RNG &&a_rng) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE String typeString() const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE MCGIDI_VIRTUAL_FUNCTION double evaluate(double a_x3, double a_x2, double a_x1) const MCGIDI_TRUE_VIRTUAL
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double sample2dOf3d(double a_x2, double a_rngValue, RNG &&a_rng, double *a_x1_1, double *a_x1_2) const
LUPI_HOST_DEVICE double sample(double a_x2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x2, double a_x1) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double evaluate(double a_x3, double a_x2, double a_x1) const
LUPI_HOST_DEVICE double sample(double a_x3, double a_x2_1, double a_x2_2, double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE void serialize(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode)
LUPI_HOST_DEVICE double sample(double a_rngValue, RNG &&a_rng) const
LUPI_HOST_DEVICE double evaluate(double a_x1) const
LUPI_HOST Function1d_d1 * parseFunction1d_d1(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST Function2d * parseFunction2d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST Function1d * parseMultiplicityFunction1d(SetupInfo &a_setupInfo, Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST Function1d_d2 * parseFunction1d_d2(GIDI::Functions::Function1dForm const *form1d)
LUPI_HOST ProbabilityBase2d_d2 * parseProbability2d_d2(GIDI::Functions::Function2dForm const *form2d, SetupInfo *a_setupInfo)
LUPI_HOST ProbabilityBase3d * parseProbability3d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST ProbabilityBase2d_d1 * parseProbability2d_d1(GIDI::Functions::Function2dForm const *form2d, SetupInfo *a_setupInfo)
LUPI_HOST ProbabilityBase1d * parseProbability1d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite)
LUPI_HOST ProbabilityBase2d * parseProbability2d(Transporting::MC const &a_settings, GIDI::Suite const &a_suite, SetupInfo *a_setupInfo)
Simple C++ string class, useful as replacement for std::string if this cannot be used,...
Definition MCGIDI.hpp:43
LUPI_HOST_DEVICE Function2dType Function2dClass(Functions::Function2d *funct)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d * serializeProbability2d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d *a_probability2d)
LUPI_HOST_DEVICE Functions::Function2d * serializeFunction2d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function2d *a_function2d)
LUPI_HOST_DEVICE ProbabilityBase1dType ProbabilityBase1dClass(Probabilities::ProbabilityBase1d *funct)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d1 * serializeProbability2d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d_d1 *a_probability2d)
LUPI_HOST_DEVICE Function1dType Function1dClass(Functions::Function1d *funct)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase3d * serializeProbability3d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase3d *a_probability3d)
LUPI_HOST_DEVICE Interpolation GIDI2MCGIDI_interpolation(ptwXY_interpolation a_interpolation)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase2d_d2 * serializeProbability2d_d2(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase2d_d2 *a_probability2d)
LUPI_HOST_DEVICE Functions::Function1d_d2 * serializeFunction1d_d2(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d2 *a_function1d)
LUPI_HOST_DEVICE Functions::Function1d * serializeFunction1d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d *a_function1d)
LUPI_HOST_DEVICE ProbabilityBase3dType ProbabilityBase3dClass(Probabilities::ProbabilityBase3d *funct)
LUPI_HOST_DEVICE ProbabilityBase2dType ProbabilityBase2dClass(Probabilities::ProbabilityBase2d *funct)
LUPI_HOST_DEVICE Functions::Function1d_d1 * serializeFunction1d_d1(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Functions::Function1d_d1 *a_function1d)
LUPI_HOST_DEVICE Probabilities::ProbabilityBase1d * serializeProbability1d(LUPI::DataBuffer &a_buffer, LUPI::DataBuffer::Mode a_mode, Probabilities::ProbabilityBase1d *a_probability1d)
enum ptwXY_interpolation_e ptwXY_interpolation