Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
tline.h File Reference

Go to the source code of this file.

Classes

class  EqualStepCoorMesh< T >
 
class  PointCoorMesh< T, D >
 
class  CopiedPointCoorMesh< T >
 
class  AbsCoorMesh< T >
 
class  VirtEqualStepCoorMesh< T >
 
class  VirtCopiedPointCoorMesh< T >
 

Functions

template<class T >
std::ostream & operator<< (std::ostream &file, const EqualStepCoorMesh< T > &f)
 
template<class T >
std::istream & operator>> (std::istream &file, EqualStepCoorMesh< T > &f)
 
template<class T >
int operator== (const EqualStepCoorMesh< T > &f1, const EqualStepCoorMesh< T > &f2)
 
template<class T >
int apeq_mant (const EqualStepCoorMesh< T > &f1, const EqualStepCoorMesh< T > &f2, T prec)
 
template<class T >
int operator!= (const EqualStepCoorMesh< T > &f1, const EqualStepCoorMesh< T > &f2)
 
template<class T , class D >
long t_find_interval (double x, long q, const D &coor)
 
template<class T , class D >
long t_find_interval_end (double x, long q, const D &coor, long n_start)
 
template<class T , class D >
std::ostream & operator<< (std::ostream &file, const PointCoorMesh< T, D > &f)
 
template<class T >
std::ostream & operator<< (std::ostream &file, CopiedPointCoorMesh< T > &f)
 
template<class T >
std::istream & operator>> (std::istream &file, CopiedPointCoorMesh< T > &f)
 
template<class T >
int operator== (const CopiedPointCoorMesh< T > &f1, const CopiedPointCoorMesh< T > &f2)
 
template<class T >
int apeq_mant (const CopiedPointCoorMesh< T > &f1, const CopiedPointCoorMesh< T > &f2, T prec)
 
template<class T >
int operator!= (const CopiedPointCoorMesh< T > &f1, const CopiedPointCoorMesh< T > &f2)
 
template<class T , class B >
int verify_types (const T *ths, const B *b2, const T **tb2)
 
template<class T >
std::ostream & operator<< (std::ostream &file, AbsCoorMesh< T > &f)
 
template<class T >
std::ostream & operator<< (std::ostream &file, VirtEqualStepCoorMesh< T > &f)
 
template<class T >
std::istream & operator>> (std::istream &file, VirtEqualStepCoorMesh< T > &f)
 
template<class T >
std::ostream & operator<< (std::ostream &file, VirtCopiedPointCoorMesh< T > &f)
 
template<class T >
std::istream & operator>> (std::istream &file, VirtCopiedPointCoorMesh< T > &f)
 
template<class T >
ActivePtr< AbsCoorMesh< T > > read_AbsCoorMesh (std::istream &file)
 
template<class T , class D , class M >
t_value_step_ar (const M &mesh, const D &y, T x, int s_include_last_point=0)
 
template<class T , class D , class M1 , class M2 >
t_value_step_ar (const M1 &mesh1, const M2 &mesh2, const D &y, T x1, T x2, int s_include_last_point=0)
 
template<class T , class D , class M >
void t_hfill_step_ar (const M &mesh, const D &y, T x, T val=1, int s_include_last_point=0)
 
template<class T , class D , class M >
void t_hfill_step_ar_ac (const M &mesh, const D &y, T x, T val=1, int s_include_last_point=0)
 
template<class T , class D , class M1 , class M2 >
void t_hfill_step_ar_ac (const M1 &mesh1, const M2 &mesh2, const D &y, T x1, T x2, T val=1, int s_include_last_point=0)
 
template<class T , class D , class M >
t_integ_step_ar (const M &mesh, const D &y, T x1, T x2, int xpower)
 
template<class T , class D , class M >
t_integ_generic_step_ar (const M &mesh, const D &y, T(*fun)(long np, T xp1, T xp2, T yp, T xmin, T xmax, T x1, T x2), T x1, T x2)
 
template<class T , class D , class M >
t_total_integ_step_ar (const M &mesh, const D &y)
 
template<class T , class D , class M1 , class M2 >
t_total_integ_step_ar (const M1 &mesh1, const M2 &mesh2, const D &y)
 
template<class T , class M1 , class M2 >
t_total_integ_step_ar (const M1 &mesh1, const M2 &mesh2, const DynArr< T > &y)
 
template<class T , class D , class M >
t_find_x_for_integ_step_ar (const M &mesh, const D &y, T integ, int *s_err)
 
template<class T , class D , class M >
t_find_x_for_already_integ_step_ar (const M &mesh, const D &y, T integ, int *s_err)
 
template<class T , class D , class M >
long t_find_entire_x_for_already_integ_step_ar (const M &mesh, const D &y, T integ, int *s_err)
 
template<class T , class D , class M >
t_hispre_step_ar (const M &mesh, const D &y, D &integ_y)
 
template<class T , class D , class M >
t_hisran_step_ar (const M &mesh, const D &integ_y, T rannum)
 
template<class T , class D , class M >
t_opposite_hisran_step_ar (const M &mesh, const D &integ_y, T x)
 
template<class T , class D , class M >
long t_entire_hisran_step_ar (const M &mesh, const D &integ_y, T rannum)
 
template<class T , class D , class M >
t_mean_step_ar (const M &mesh, const D &y, T x1, T x2, int &s_err)
 
template<class T >
t_value_straight_2point (T x1, T y1, T x2, T y2, T x, int s_ban_neg)
 
template<class T >
t_integ_straight_2point (T x1, T y1, T x2, T y2, T xl, T xr, int xpower, int s_ban_neg)
 
template<class T , class D , class M >
t_value_straight_point_ar (const M &mesh, const D &y, T x, int s_ban_neg, int s_extrap_left, T left_bond, int s_extrap_right, T right_bond)
 
template<class T , class D , class M >
t_value_generic_point_ar (const M &mesh, const D &y, T(*funval)(T xp1, T yp1, T xp2, T yp2, T xmin, T xmax, T x), T x, int s_extrap_left, T left_bond, int s_extrap_right, T right_bond)
 
template<class T >
t_value_power_2point (T x1, T y1, T x2, T y2, T x)
 
template<class T >
t_value_power_extended_2point (T x1, T y1, T x2, T y2, T x)
 
template<class T >
t_value_exp_2point (T x1, T y1, T x2, T y2, T x)
 
template<class T >
t_integ_power_2point (T x1, T y1, T x2, T y2, T xl, T xr)
 
template<class T , class D , class M >
t_integ_straight_point_ar (const M &mesh, const D &y, T x1, T x2, int xpower, int s_ban_neg, int s_extrap_left, T left_bond, int s_extrap_right, T right_bond)
 
template<class T , class D , class M >
t_mean_straight_point_ar (const M &mesh, const D &y, T x1, T x2, int s_extrap_left, T left_bond, int s_extrap_right, T right_bond, int &s_err)
 
template<class T , class D , class M >
t_integ_generic_point_ar (const M &mesh, const D &y, T(*fun)(T xp1, T yp1, T xp2, T yp2, T xmin, T xmax, T x1, T x2), T x1, T x2, int s_extrap_left, T left_bond, int s_extrap_right, T right_bond)
 
template<class T , class D , class M >
t_width_at_hheight_step_ar (const M &mesh, const D &y)
 

Function Documentation

◆ apeq_mant() [1/2]

template<class T >
int apeq_mant ( const CopiedPointCoorMesh< T > &  f1,
const CopiedPointCoorMesh< T > &  f2,
prec 
)

Definition at line 1150 of file tline.h.

1151 {
1152 if (f1.get_qi() != f2.get_qi() ||
1153 !apeq_mant(f1.get_xmin(), f2.get_xmin(), prec) ||
1154 !apeq_mant(f1.get_xmax(), f2.get_xmax(), prec))
1155 return 0;
1156 if (!apeq_mant(f1.get_mesh(), f2.get_mesh(), prec))
1157 return 0;
1158 else
1159 return 1;
1160}
const DynLinArr< T > & get_mesh(void) const
Definition: tline.h:766
long get_qi(void) const
Definition: tline.h:763
T get_xmax(void) const
Definition: tline.h:765
T get_xmin(void) const
Definition: tline.h:764
int apeq_mant(const EqualStepCoorMesh< T > &f1, const EqualStepCoorMesh< T > &f2, T prec)
Definition: tline.h:290

◆ apeq_mant() [2/2]

template<class T >
int apeq_mant ( const EqualStepCoorMesh< T > &  f1,
const EqualStepCoorMesh< T > &  f2,
prec 
)

Definition at line 290 of file tline.h.

291 {
292 if (f1.get_qi() != f2.get_qi() ||
293 !apeq_mant(f1.get_xmin(), f2.get_xmin(), prec) ||
294 !apeq_mant(f1.get_xmax(), f2.get_xmax(), prec)) {
295 Iprintn(mcout, !apeq_mant(f1.get_xmin(), f2.get_xmin(), prec));
296 Iprintn(mcout, !apeq_mant(f1.get_xmax(), f2.get_xmax(), prec));
297 return 0;
298 } else
299 return 1;
300}
long get_qi(void) const
Definition: tline.h:62
T get_xmin(void) const
Definition: tline.h:65
T get_xmax(void) const
Definition: tline.h:66
#define mcout
Definition: prstream.h:133
#define Iprintn(file, name)
Definition: prstream.h:216

Referenced by apeq_mant(), t_entire_hisran_step_ar(), t_hisran_step_ar(), and t_opposite_hisran_step_ar().

◆ operator!=() [1/2]

template<class T >
int operator!= ( const CopiedPointCoorMesh< T > &  f1,
const CopiedPointCoorMesh< T > &  f2 
)

Definition at line 1163 of file tline.h.

1164 {
1165 if (f1.get_qi() != f2.get_qi() || f1.get_xmin() != f2.get_xmin() ||
1166 f1.get_xmax() != f2.get_xmax())
1167 return 1;
1168 if (f1.get_mesh() != f2.get_mesh())
1169 return 1;
1170 else
1171 return 0;
1172}

◆ operator!=() [2/2]

template<class T >
int operator!= ( const EqualStepCoorMesh< T > &  f1,
const EqualStepCoorMesh< T > &  f2 
)

Definition at line 303 of file tline.h.

304 {
305 if (f1.get_qi() != f2.get_qi() || f1.get_xmin() != f2.get_xmin() ||
306 f1.get_xmax() != f2.get_xmax())
307 return 1;
308 else
309 return 0;
310}

◆ operator<<() [1/6]

template<class T >
std::ostream & operator<< ( std::ostream &  file,
AbsCoorMesh< T > &  f 
)

Definition at line 1273 of file tline.h.

1273 {
1274 f.print(file);
1275 return file;
1276}
virtual void print(std::ostream &file) const
Definition: tline.h:1232

◆ operator<<() [2/6]

template<class T >
std::ostream & operator<< ( std::ostream &  file,
const EqualStepCoorMesh< T > &  f 
)

Definition at line 257 of file tline.h.

258 {
259 f.print(file);
260 return file;
261}
void print(std::ostream &file) const
Definition: tline.h:248

◆ operator<<() [3/6]

template<class T , class D >
std::ostream & operator<< ( std::ostream &  file,
const PointCoorMesh< T, D > &  f 
)

Definition at line 739 of file tline.h.

739 {
740 f.print(file);
741 return file;
742}
virtual void print(std::ostream &file) const
Definition: tline.h:715

◆ operator<<() [4/6]

template<class T >
std::ostream & operator<< ( std::ostream &  file,
CopiedPointCoorMesh< T > &  f 
)

Definition at line 1121 of file tline.h.

1121 {
1122 f.print(file);
1123 return file;
1124}
void print(std::ostream &file) const
Definition: tline.h:1096

◆ operator<<() [5/6]

template<class T >
std::ostream & operator<< ( std::ostream &  file,
VirtCopiedPointCoorMesh< T > &  f 
)

Definition at line 1501 of file tline.h.

1502 {
1503 f.print(file);
1504 return file;
1505}
virtual void print(std::ostream &file) const
Definition: tline.h:1460

◆ operator<<() [6/6]

template<class T >
std::ostream & operator<< ( std::ostream &  file,
VirtEqualStepCoorMesh< T > &  f 
)

Definition at line 1368 of file tline.h.

1368 {
1369 f.print(file);
1370 return file;
1371}
virtual void print(std::ostream &file) const
Definition: tline.h:1326

◆ operator==() [1/2]

template<class T >
int operator== ( const CopiedPointCoorMesh< T > &  f1,
const CopiedPointCoorMesh< T > &  f2 
)

Definition at line 1138 of file tline.h.

1139 {
1140 if (f1.get_qi() != f2.get_qi() || f1.get_xmin() != f2.get_xmin() ||
1141 f1.get_xmax() != f2.get_xmax())
1142 return 0;
1143 if (f1.get_mesh() != f2.get_mesh())
1144 return 0;
1145 else
1146 return 1;
1147}

◆ operator==() [2/2]

template<class T >
int operator== ( const EqualStepCoorMesh< T > &  f1,
const EqualStepCoorMesh< T > &  f2 
)

Definition at line 280 of file tline.h.

281 {
282 if (f1.get_qi() != f2.get_qi() || f1.get_xmin() != f2.get_xmin() ||
283 f1.get_xmax() != f2.get_xmax())
284 return 0;
285 else
286 return 1;
287}

◆ operator>>() [1/4]

template<class T >
std::istream & operator>> ( std::istream &  file,
CopiedPointCoorMesh< T > &  f 
)

Definition at line 1127 of file tline.h.

1127 {
1128 mfunname("istream& operator>>(istream& file, CopiedPointCoorMesh<T>& f)");
1129 definp_endpar dep(&file, 0, 1, 0);
1130 set_position("Type of T is (in internal notations)", *dep.istrm, dep.s_rewind,
1131 dep.s_req_sep);
1132 DynLinArr<T> mesh;
1133 f = CopiedPointCoorMesh<T>(mesh);
1134 return file;
1135}
#define mfunname(string)
Definition: FunNameStack.h:67
long set_position(const String &word, std::istream &istrm, int s_rewind, int s_req_sep)
Definition: definp.cpp:71

◆ operator>>() [2/4]

template<class T >
std::istream & operator>> ( std::istream &  file,
EqualStepCoorMesh< T > &  f 
)

Definition at line 264 of file tline.h.

264 {
265 mfunname("istream& operator>>(istream& file, EqualStepCoorMesh<T>& f)");
266 definp_endpar dep(&file, 0, 1, 0);
267 set_position("Type of T is (in internal notations)", *dep.istrm, dep.s_rewind,
268 dep.s_req_sep);
269 long q;
270 T xmin;
271 T xmax;
272 DEFINPAP(q);
273 DEFINPAP(xmin);
274 DEFINPAP(xmax);
275 f = EqualStepCoorMesh<T>(q, xmin, xmax);
276 return file;
277}
#define DEFINPAP(name)
Definition: definp.h:91

◆ operator>>() [3/4]

template<class T >
std::istream & operator>> ( std::istream &  file,
VirtCopiedPointCoorMesh< T > &  f 
)

Definition at line 1508 of file tline.h.

1509 {
1510 mfunname("istream& operator>>(istream& file, VirtCopiedPointCoorMesh<T>& f)");
1511 definp_endpar dep(&file, 0, 1, 0);
1512 set_position("CopiedPointCoorMesh<T>:", *dep.istrm, dep.s_rewind,
1513 dep.s_req_sep);
1514 // file>>f; // this would lead to loop
1515 //file>>(static_cast< CopiedPointCoorMesh<T> >(f));
1516 CopiedPointCoorMesh<T>& ft = f;
1517 file >> ft;
1518 return file;
1519}

◆ operator>>() [4/4]

template<class T >
std::istream & operator>> ( std::istream &  file,
VirtEqualStepCoorMesh< T > &  f 
)

Definition at line 1374 of file tline.h.

1374 {
1375 mfunname("istream& operator>>(istream& file, VirtEqualStepCoorMesh<T>& f)");
1376 definp_endpar dep(&file, 0, 1, 0);
1377 set_position("EqualStepCoorMesh<T>:", *dep.istrm, dep.s_rewind,
1378 dep.s_req_sep);
1379 //file>>f;
1380 //file>>(static_cast< EqualStepCoorMesh<T> >(f));
1381 EqualStepCoorMesh<T>& ft = f;
1382 file >> ft;
1383
1384 return file;
1385}

◆ read_AbsCoorMesh()

template<class T >
ActivePtr< AbsCoorMesh< T > > read_AbsCoorMesh ( std::istream &  file)

Definition at line 1604 of file tline.h.

1604 {
1605 mfunnamep("ActivePtr< AbsCoorMesh<T> > read_AbsCoorMesh(istream& file)");
1606 const long q = 28;
1607 long n;
1608 char keyline[q];
1609 for (n = 0; n < q - 1; n++) {
1610 //file.get(keyline[n]);
1611 keyline[n] = file.get();
1612 }
1613 keyline[n] = '\0';
1614 if (!strcmp(&(keyline[0]), "AbsCoorMesh: no content.")) {
1615 funnw.ehdr(mcerr);
1616 mcout << "attempt to read AbsCoorMesh\n";
1617 spexit(mcerr);
1618 }
1619 if (!strcmp(&(keyline[0]), "VirtEqualStepCoorMesh<T>: ")) {
1621 file >> (*amesh);
1622 return ActivePtr<AbsCoorMesh<T> >(amesh, dont_clone);
1623 }
1624 if (!strcmp(&(keyline[0]), "VirtCopiedPointCoorMesh<T>:")) {
1626 file >> (*amesh);
1627 return ActivePtr<AbsCoorMesh<T> >(amesh, dont_clone);
1628 }
1629 funnw.ehdr(mcerr);
1630 mcout << "cannot identifiy the type of derived AbsCoorMesh\n";
1631 Iprintn(mcout, keyline);
1632 spexit(mcerr);
1633#ifdef INS_CRETURN
1634 return ActivePtr<AbsCoorMesh<T> >();
1635#endif
1636}
@ dont_clone
Definition: AbsPtr.h:519
#define mfunnamep(string)
Definition: FunNameStack.h:77
#define spexit(stream)
Definition: FunNameStack.h:536
#define mcerr
Definition: prstream.h:135

◆ t_entire_hisran_step_ar()

template<class T , class D , class M >
long t_entire_hisran_step_ar ( const M &  mesh,
const D &  integ_y,
rannum 
)

Definition at line 2394 of file tline.h.

2394 {
2395 mfunname("double t_entire_hisran_step_ar(...)");
2396
2397 //mcout<<"start t_step_integ_ar\n";
2398 //check_econd11(xpower , != 0 , mcerr);
2399 long qi = mesh.get_qi();
2400 long s_same = apeq_mant(integ_y[qi - 1], 1.0, 1.0e-12);
2401 check_econd11a(s_same, != 1.0, "integ_y[qi-1]=" << integ_y[qi - 1] << '\n',
2402 mcerr);
2403 //check_econd11(integ_y[qi-1] , != 1.0 , mcerr);
2404 int s_err;
2405
2406 long ret =
2408 integ_y, // dimension q-1
2409 rannum, &s_err);
2410 check_econd11a(s_err, != 0, "mesh=" << mesh << " integ_y=" << integ_y
2411 << " rannum=" << rannum << '\n',
2412 mcerr);
2413 return ret;
2414 //return t_find_entire_x_for_already_integ_step_ar
2415 // (mesh, // dimension q
2416 // integ_y, // dimension q-1
2417 // rannum,
2418 // &s_err);
2419}
#define check_econd11a(a, signb, add, stream)
Definition: FunNameStack.h:395
long t_find_entire_x_for_already_integ_step_ar(const M &mesh, const D &y, T integ, int *s_err)
Definition: tline.h:2241

◆ t_find_entire_x_for_already_integ_step_ar()

template<class T , class D , class M >
long t_find_entire_x_for_already_integ_step_ar ( const M &  mesh,
const D &  y,
integ,
int *  s_err 
)

Definition at line 2241 of file tline.h.

2244 {
2245 mfunname("double t_find_entire_x_for_already_integ_step_ar(...)");
2246 //Iprintn(mcout, mesh);
2247 //Iprintn(mcout, integ);
2248 *s_err = 0;
2249 //mcout<<"start t_step_integ_ar\n";
2250 //check_econd11(xpower , != 0 , mcerr);
2251 check_econd11(integ, < 0.0, mcerr);
2252 long qi = mesh.get_qi();
2253 check_econd12(qi, <, 1, mcerr);
2254 //if(x1 > x2) return 0.0;
2255 long xmin = mesh.get_xmin();
2256 long xmax = mesh.get_xmax();
2257 if (integ == 0) return xmin;
2258 if (integ > y[qi - 1]) {
2259 *s_err = 1;
2260 return xmax;
2261 }
2262 if (integ == y[qi - 1]) return xmax;
2263 if (integ < y[0]) { // answer in the first bin
2264 long xp1(0);
2265 mesh.get_scoor(0, xp1);
2266 return xp1;
2267 }
2268 // binary search
2269 long nl = 0;
2270 long nr = qi - 1;
2271 long nc;
2272 while (nr - nl > 1) {
2273 nc = (nr + nl) / 2;
2274 if (integ < y[nc])
2275 nr = nc;
2276 else
2277 nl = nc;
2278 }
2279 //Iprint2n(mcout, nl, nr);
2280 //Iprint2n(mcout, y[nl], y[nr]);
2281 //Iprint2n(mcout, nl, nr);
2282 long x(0);
2283 mesh.get_scoor(nr, x);
2284 //Iprintn(mcout, x);
2285
2286 return x;
2287}
#define check_econd11(a, signb, stream)
Definition: FunNameStack.h:366
#define check_econd12(a, sign, b, stream)
Definition: FunNameStack.h:380

Referenced by t_entire_hisran_step_ar().

◆ t_find_interval()

template<class T , class D >
long t_find_interval ( double  x,
long  q,
const D &  coor 
)

Definition at line 314 of file tline.h.

314 {
315 long n1, n2, n3;
316#ifndef TLINE_REDUCE_TO_RAW_ARR
317 if (q <= 1) return -1;
318 if (x < coor[0] || x > coor[q - 1]) return -1;
319 if (x < coor[1]) return 0;
320 if (x >= coor[q - 2]) return q - 2;
321 n1 = 0;
322 n2 = q - 1;
323 while (n2 - n1 > 1) {
324 n3 = n1 + (n2 - n1) / 2;
325 if (x < coor[n3])
326 n2 = n3;
327 else
328 n1 = n3;
329 }
330 return n1;
331#else
332 T* arr = &(coor[0]); // take the address of the first element
333 if (q <= 1) return -1;
334 if (x < arr[0] || x > arr[q - 1]) return -1;
335 if (x < arr[1]) return 0;
336 if (x >= arr[q - 2]) return q - 2;
337 n1 = 0;
338 n2 = q - 1;
339 while (n2 - n1 > 1) {
340 n3 = n1 + (n2 - n1) / 2;
341 if (x < arr[n3])
342 n2 = n3;
343 else
344 n1 = n3;
345 }
346 return n1;
347
348#endif
349}

◆ t_find_interval_end()

template<class T , class D >
long t_find_interval_end ( double  x,
long  q,
const D &  coor,
long  n_start 
)

Definition at line 358 of file tline.h.

358 {
359 long n1, n2, n3;
360 if (n_start < 0 || n_start > q - 1) {
361 mcerr << " ERROR in t_find_interval_end(...):\n";
362 mcerr << "n_start < 0 || n_start > q-1\n";
363 Iprint2n(mcout, n_start, q);
364 spexit(mcerr);
365 }
366#ifndef TLINE_REDUCE_TO_RAW_ARR
367 //if(q <= 1) return -1;
368 if (q - n_start <= 1) return -1;
369 if (x < coor[n_start] || x > coor[q - 1]) return -1;
370 if (x < coor[n_start + 1]) return n_start;
371 if (x >= coor[q - 2]) return q - 2;
372 n1 = n_start;
373 n2 = q - 1;
374 while (n2 - n1 > 1) {
375 n3 = n1 + (n2 - n1) / 2;
376 if (x < coor[n3])
377 n2 = n3;
378 else
379 n1 = n3;
380 }
381 return n1;
382#else
383 T* arr = &(coor[0]); // take the address of the first element
384 //if(q <= 1) return -1;
385 if (q - n_start <= 1) return -1;
386 if (x < arr[n_start] || x > arr[q - 1]) return -1;
387 if (x < arr[n_start + 1]) return n_start;
388 if (x >= arr[q - 2]) return q - 2;
389 n1 = n_start;
390 n2 = q - 1;
391 while (n2 - n1 > 1) {
392 n3 = n1 + (n2 - n1) / 2;
393 if (x < arr[n3])
394 n2 = n3;
395 else
396 n1 = n3;
397 }
398 return n1;
399
400#endif
401}
#define Iprint2n(file, name1, name2)
Definition: prstream.h:236

◆ t_find_x_for_already_integ_step_ar()

template<class T , class D , class M >
T t_find_x_for_already_integ_step_ar ( const M &  mesh,
const D &  y,
integ,
int *  s_err 
)

Definition at line 2176 of file tline.h.

2179 {
2180 mfunname("double t_find_x_for_already_integ_step_ar(...)");
2181
2182 *s_err = 0;
2183 //mcout<<"start t_step_integ_ar\n";
2184 //check_econd11(xpower , != 0 , mcerr);
2185 check_econd11(integ, < 0.0, mcerr);
2186 long qi = mesh.get_qi();
2187 check_econd12(qi, <, 1, mcerr);
2188 //if(x1 > x2) return 0.0;
2189 double xmin = mesh.get_xmin();
2190 double xmax = mesh.get_xmax();
2191 if (integ == 0.0) return xmin;
2192 if (integ > y[qi - 1]) {
2193 *s_err = 1;
2194 return xmax;
2195 }
2196 if (integ == y[qi - 1]) return xmax;
2197 if (integ < y[0]) { // answer in the first bin
2198 T xp1(0.0);
2199 T xp2(0.0);
2200 mesh.get_scoor(0, xp1);
2201 mesh.get_scoor(1, xp2);
2202 return xp1 + (xp2 - xp1) * integ / y[0];
2203 }
2204 // binary search
2205 long nl = 0;
2206 long nr = qi - 1;
2207 long nc;
2208 while (nr - nl > 1) {
2209 nc = (nr + nl) / 2;
2210 if (integ < y[nc])
2211 nr = nc;
2212 else
2213 nl = nc;
2214 }
2215 //Iprint2n(mcout, nl, nr);
2216 T xl(0.0);
2217 T xr(0.0);
2218 mesh.get_scoor(nl + 1, xl);
2219 mesh.get_scoor(nr + 1, xr);
2220 // Note "+1" in the previous two expressions.
2221 // This arises from the fact that the nl'th element of
2222 // y-array contains integral of nl'th bin plus all previous bins.
2223 // So y[nl] is related to x of nl+1.
2224 T a = (xr - xl) / (y[nr] - y[nl]);
2225 T ret = xl + a * (integ - y[nl]);
2226
2227 return ret;
2228}

Referenced by t_hisran_step_ar().

◆ t_find_x_for_integ_step_ar()

template<class T , class D , class M >
T t_find_x_for_integ_step_ar ( const M &  mesh,
const D &  y,
integ,
int *  s_err 
)

Definition at line 2127 of file tline.h.

2130 {
2131 mfunname("double t_integ_step_ar(...)");
2132
2133 *s_err = 0;
2134 //mcout<<"start t_step_integ_ar\n";
2135 //check_econd11(xpower , != 0 , mcerr);
2136 check_econd11(integ, < 0.0, mcerr);
2137 long qi = mesh.get_qi();
2138 check_econd12(qi, <, 1, mcerr);
2139 //if(x1 > x2) return 0.0;
2140 double xmin = mesh.get_xmin();
2141 double xmax = mesh.get_xmax();
2142 if (integ == 0.0) return xmin;
2143 T s(0.0);
2144 long n = 0;
2145 T xp1(0.0);
2146 T xp2(0.0);
2147 mesh.get_scoor(n, xp2);
2148 for (n = 0; n < qi; n++) {
2149 xp1 = xp2;
2150 mesh.get_scoor(n + 1, xp2);
2151 T step = xp2 - xp1;
2152 T s1 = s + y[n] * step;
2153 //Iprint3n(mcout, n, s1, integ);
2154 if (s1 > integ) break;
2155 if (s1 == integ) return xp2;
2156 s = s1;
2157 }
2158
2159 if (n == qi) {
2160 *s_err = 1;
2161 return xmax;
2162 }
2163 double x = xp1;
2164 //Iprint3n(mcout, n, s, x);
2165 x += (integ - s) / y[n];
2166 //Iprintn(mcout, x);
2167 return x;
2168}

◆ t_hfill_step_ar()

template<class T , class D , class M >
void t_hfill_step_ar ( const M &  mesh,
const D &  y,
x,
val = 1,
int  s_include_last_point = 0 
)

Definition at line 1714 of file tline.h.

1717 {
1718 mfunname("double t_hfill_step_ar(...)");
1719 double xmin = mesh.get_xmin();
1720 double xmax = mesh.get_xmax();
1721 //Iprint3n(mcout, x, xmin, xmax);
1722 if (x < xmin) return;
1723 if (s_include_last_point == 0) {
1724 if (x >= xmax) return;
1725 } else {
1726 if (x > xmax) return;
1727 }
1728 long n1;
1729 int i_ret = 0;
1730 i_ret = mesh.get_interval(x, n1);
1731 check_econd11(i_ret, != 1, mcerr);
1732 y[n1] += val;
1733 return;
1734}

◆ t_hfill_step_ar_ac() [1/2]

template<class T , class D , class M >
void t_hfill_step_ar_ac ( const M &  mesh,
const D &  y,
x,
val = 1,
int  s_include_last_point = 0 
)

Definition at line 1740 of file tline.h.

1743 {
1744 mfunname("double t_hfill_step_ar(...)");
1745 double xmin = mesh.get_xmin();
1746 double xmax = mesh.get_xmax();
1747 //Iprint3n(mcout, x, xmin, xmax);
1748 if (x < xmin) return;
1749 if (s_include_last_point == 0) {
1750 if (x >= xmax) return;
1751 } else {
1752 if (x > xmax) return;
1753 }
1754 long n1;
1755 int i_ret = 0;
1756 i_ret = mesh.get_interval(x, n1);
1757 check_econd11(i_ret, != 1, mcerr);
1758 y.ac(n1) += val;
1759 return;
1760}

◆ t_hfill_step_ar_ac() [2/2]

template<class T , class D , class M1 , class M2 >
void t_hfill_step_ar_ac ( const M1 &  mesh1,
const M2 &  mesh2,
const D &  y,
x1,
x2,
val = 1,
int  s_include_last_point = 0 
)

Definition at line 1764 of file tline.h.

1768 {
1769 mfunname("double t_hfill_step_ar(...)");
1770 double x1min = mesh1.get_xmin();
1771 double x1max = mesh1.get_xmax();
1772 double x2min = mesh2.get_xmin();
1773 double x2max = mesh2.get_xmax();
1774 //Iprint3n(mcout, x, xmin, xmax);
1775 if (x1 < x1min) return;
1776 if (s_include_last_point == 0) {
1777 if (x1 >= x1max) return;
1778 } else {
1779 if (x1 > x1max) return;
1780 }
1781 if (x2 < x2min) return;
1782 if (s_include_last_point == 0) {
1783 if (x2 >= x2max) return;
1784 } else {
1785 if (x2 > x2max) return;
1786 }
1787 long n1;
1788 int i_ret1 = 0;
1789 i_ret1 = mesh1.get_interval(x1, n1);
1790 check_econd11(i_ret1, != 1, mcerr);
1791 long n2;
1792 int i_ret2 = 0;
1793 i_ret2 = mesh2.get_interval(x2, n2);
1794 check_econd11(i_ret2, != 1, mcerr);
1795 y.ac(n1, n2) += val;
1796 return;
1797}

◆ t_hispre_step_ar()

template<class T , class D , class M >
T t_hispre_step_ar ( const M &  mesh,
const D &  y,
D &  integ_y 
)

Definition at line 2295 of file tline.h.

2297 {
2298 mfunname("double t_hispre_step_ar(...)");
2299
2300 //mcout<<"start t_step_integ_ar\n";
2301 //check_econd11(xpower , != 0 , mcerr);
2302 long qi = mesh.get_qi();
2303 check_econd12(qi, <, 1, mcerr);
2304
2305 T s(0.0);
2306 long n = 0;
2307 T xp1(0.0);
2308 T xp2(0.0);
2309 mesh.get_scoor(n, xp2);
2310 for (n = 0; n < qi; n++) {
2311 xp1 = xp2;
2312 mesh.get_scoor(n + 1, xp2);
2313 T step = xp2 - xp1;
2314 check_econd11a(y[n], < 0.0,
2315 "n=" << n << " xp1=" << xp1 << " xp2=" << xp2 << '\n',
2316 mcerr);
2317 s = s + y[n] * step;
2318 integ_y[n] = s;
2319 //Iprint3n(mcout, n, s1, integ);
2320 }
2321 check_econd11a(s, <= 0.0, "y=" << y << " integ_y=" << integ_y << '\n', mcerr);
2322 for (n = 0; n < qi; n++) {
2323 integ_y[n] /= s;
2324 }
2325 return s;
2326}

◆ t_hisran_step_ar()

template<class T , class D , class M >
T t_hisran_step_ar ( const M &  mesh,
const D &  integ_y,
rannum 
)

Definition at line 2331 of file tline.h.

2331 {
2332 mfunname("double t_hisran_step_ar(...)");
2333 //mcout<<"start t_step_integ_ar\n";
2334 //check_econd11(xpower , != 0 , mcerr);
2335 long qi = mesh.get_qi();
2336 long s_same = apeq_mant(integ_y[qi - 1], 1.0, 1.0e-12);
2337 check_econd11a(s_same, != 1.0, "integ_y[qi-1]=" << integ_y[qi - 1] << '\n',
2338 mcerr);
2339
2340 //Iprintn(mcout, rannum);
2341 //check_econd11(integ_y[qi-1] , != 1.0 , mcerr);
2342 int s_err;
2343
2344 T ret = t_find_x_for_already_integ_step_ar(mesh, // dimension q
2345 integ_y, // dimension q-1
2346 rannum, &s_err);
2347 check_econd11a(s_err, != 0, "mesh=" << mesh << " integ_y=" << integ_y
2348 << " rannum=" << rannum << '\n',
2349 mcerr);
2350 return ret;
2351 //return t_find_x_for_already_integ_step_ar
2352 // (mesh, // dimension q
2353 // integ_y, // dimension q-1
2354 // rannum,
2355 // &s_err);
2356}
T t_find_x_for_already_integ_step_ar(const M &mesh, const D &y, T integ, int *s_err)
Definition: tline.h:2176

◆ t_integ_generic_point_ar()

template<class T , class D , class M >
T t_integ_generic_point_ar ( const M &  mesh,
const D &  y,
T(*)(T xp1, T yp1, T xp2, T yp2, T xmin, T xmax, T x1, T x2)  fun,
x1,
x2,
int  s_extrap_left,
left_bond,
int  s_extrap_right,
right_bond 
)

Definition at line 2786 of file tline.h.

2790 {
2791 mfunname("double t_integ_generic_point_ar(...)");
2792
2793 //mcout<<"Strart t_integ_straight_point_ar\n";
2794 //check_econd21(xpower , != 0 && , != 1 , mcerr);
2795 check_econd12(x1, >, x2, mcerr);
2796 long qi = mesh.get_qi();
2797 check_econd12(qi, <, 1, mcerr);
2798 //if(x1 > x2) return 0.0;
2799 double xmin = mesh.get_xmin();
2800 double xmax = mesh.get_xmax();
2801 if (x2 <= xmin && s_extrap_left == 0) return 0.0;
2802 if (x1 >= xmax && s_extrap_right == 0) return 0.0;
2803 if (x2 <= left_bond) return 0.0;
2804 if (x1 >= right_bond) return 0.0;
2805 // long istart, iafterend; // indexes to sum total intervals
2806 if (x1 < left_bond) x1 = left_bond;
2807 if (x2 > right_bond) x2 = right_bond;
2808 if (x1 <= xmin && s_extrap_left == 0) x1 = xmin;
2809 if (x2 > xmax && s_extrap_left == 0) x2 = xmax;
2810 long np1, np2;
2811 T bp1, bp2;
2812 int i_ret = 0;
2813 // restore the interval in which x1 reside
2814 i_ret = mesh.get_interval_extrap(x1, np1, bp1, np2, bp2);
2815 // restore the x-coordinates of given points
2816 T xp1;
2817 mesh.get_scoor(np1, xp1);
2818 T xp2;
2819 mesh.get_scoor(np2, xp2);
2820 T res;
2821 T yp1 = y[np1];
2822 T yp2 = y[np2];
2823 if (i_ret == 2 || x2 <= xp2) // then all in one interval
2824 {
2825 res = fun(xp1, yp1, xp2, yp2, xmin, xmax, x1, x2);
2826 } else {
2827 // integrate only till end of the current interval
2828 T x1i = x1;
2829 T x2i = xp2;
2830 res = fun(xp1, yp1, xp2, yp2, xmin, xmax, x1i, x2i);
2831 //x2i = x1; // prepere for loop
2832 int s_stop = 0;
2833 do {
2834 np1 = np2;
2835 np2++;
2836 xp1 = xp2;
2837 mesh.get_scoor(np2, xp2);
2838 x1i = x2i;
2839 if (xp2 >= x2) {
2840 x2i = x2; // till end of integral
2841 s_stop = 1;
2842 } else {
2843 if (np2 == qi) // end of the mesh, but x2 is farther
2844 {
2845 x2i = x2; // till end of integral
2846 s_stop = 1;
2847 } else {
2848 x2i = xp2; // till end of current interval
2849 s_stop = 0;
2850 }
2851 }
2852 yp1 = yp2;
2853 yp2 = y[np2];
2854 res += fun(xp1, yp1, xp2, yp2, xmin, xmax, x1i, x2i);
2855 //Iprint2n(mcout, xp1, xp2);
2856 //Iprint2n(mcout, x1i, x2i);
2857 //Iprint2n(mcout, yp1, yp2);
2858 //Iprint2n(mcout, res, s_stop);
2859
2860 } while (s_stop == 0);
2861 }
2862 return res;
2863}
void(* fun)(T &f))
Definition: AbsArr.h:591

Referenced by Heed::glin_integ_ar().

◆ t_integ_generic_step_ar()

template<class T , class D , class M >
T t_integ_generic_step_ar ( const M &  mesh,
const D &  y,
T(*)(long np, T xp1, T xp2, T yp, T xmin, T xmax, T x1, T x2)  fun,
x1,
x2 
)

Definition at line 1907 of file tline.h.

1913 {
1914 mfunname("double t_integ_step_ar(...)");
1915
1916 //mcout<<"start t_step_integ_ar\n";
1917 check_econd12(x1, >, x2, mcerr);
1918 long qi = mesh.get_qi();
1919 check_econd12(qi, <, 1, mcerr);
1920 //if(x1 > x2) return 0;
1921 double xmin = mesh.get_xmin();
1922 double xmax = mesh.get_xmax();
1923 if (x2 <= xmin) return 0;
1924 if (x1 >= xmax) return 0;
1925 if (x1 == x2) return 0;
1926 long istart, iafterend; // indexes to sum total intervals
1927 T s(0);
1928 if (x1 <= xmin) {
1929 x1 = xmin;
1930 istart = 0;
1931 } else {
1932 long n1, n2;
1933 T b1, b2;
1934 int i_ret = 0;
1935 i_ret = mesh.get_interval(x1, n1, b1, n2, b2);
1936 //Iprint2n(mcout, x1, i_ret);
1937 //Iprint4n(mcout, n1, b1, n2, b2);
1938 check_econd11(i_ret, != 1, mcerr);
1939 if (b2 - x1 > 0) // otherwise it could be only equal to 0
1940 {
1941 if (x2 <= b2) // if x2 in the same interval
1942 {
1943 s = fun(n1, b1, b2, y[n1], xmin, xmax, x1, x2);
1944 return s;
1945 }
1946 s = fun(n1, b1, b2, y[n1], xmin, xmax, x1, x2);
1947 }
1948 istart = n2;
1949 }
1950 if (x2 >= xmax) {
1951 x2 = xmax;
1952 iafterend = qi;
1953 } else {
1954 long n1, n2;
1955 T b1, b2;
1956 int i_ret = 0;
1957 i_ret = mesh.get_interval(x2, n1, b1, n2, b2);
1958 //Iprint2n(mcout, x2, i_ret);
1959 //Iprint4n(mcout, n1, b1, n2, b2);
1960 check_econd11(i_ret, != 1, mcerr);
1961 if (x2 - b1 > 0) {
1962 s += fun(n1, b1, b2, y[n1], xmin, xmax, b1, x2);
1963 }
1964 iafterend = n1;
1965 }
1966 //Iprint2n(mcout, istart, iafterend);
1967 long i;
1968 double b;
1969 mesh.get_scoor(istart, b);
1970 for (i = istart; i < iafterend; i++) {
1971 double a = b;
1972 mesh.get_scoor(i + 1, b);
1973 s += fun(i, a, b, y[i], xmin, xmax, a, b);
1974 }
1975 //Iprintn(mcout, s);
1976
1977 //T t;
1978 return s;
1979
1980}

◆ t_integ_power_2point()

template<class T >
T t_integ_power_2point ( x1,
y1,
x2,
y2,
xl,
xr 
)

Definition at line 2661 of file tline.h.

2663 {
2664 mfunname("double t_integ_power_2point(...)");
2665
2666 check_econd11(y1, <= 0.0, mcerr);
2667 check_econd11(y2, <= 0.0, mcerr);
2668 check_econd12(y1, ==, y2, mcerr);
2669 check_econd12(x1, ==, x2, mcerr);
2670 T pw = log(y1 / y2) / log(x1 / x2);
2671 check_econd11(pw, == -1.0, mcerr);
2672 T k = y1 * pow(x1, -pw);
2673 T t = k / (1 + pw) * (pow(xr, (pw + 1)) - pow(xl, (pw + 1)));
2674 return t;
2675}
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:336

◆ t_integ_step_ar()

template<class T , class D , class M >
T t_integ_step_ar ( const M &  mesh,
const D &  y,
x1,
x2,
int  xpower 
)

Definition at line 1810 of file tline.h.

1812 {
1813 mfunname("double t_integ_step_ar(...)");
1814
1815 // mcout<<"start t_step_integ_ar\n";
1816 check_econd21(xpower, != 0 &&, != 1, mcerr);
1817 check_econd12(x1, >, x2, mcerr);
1818 long qi = mesh.get_qi();
1819 check_econd12(qi, <, 1, mcerr);
1820 //if(x1 > x2) return 0;
1821 double xmin = mesh.get_xmin();
1822 double xmax = mesh.get_xmax();
1823 if (x2 <= xmin) return 0;
1824 if (x1 >= xmax) return 0;
1825 if (x1 == x2) return 0;
1826 long istart, iafterend; // indexes to sum total intervals
1827 T s(0);
1828 if (x1 <= xmin) {
1829 x1 = xmin;
1830 istart = 0;
1831 } else {
1832 long n1, n2;
1833 T b1, b2;
1834 int i_ret = 0;
1835 i_ret = mesh.get_interval(x1, n1, b1, n2, b2);
1836 //Iprint2n(mcout, x1, i_ret);
1837 //Iprint4n(mcout, n1, b1, n2, b2);
1838 check_econd11(i_ret, != 1, mcerr);
1839 if (b2 - x1 > 0) { // otherwise it could be only equal to 0
1840 if (x2 <= b2) { // if x2 in the same interval
1841 if (xpower == 0) {
1842 s = (x2 - x1) * y[n1];
1843 } else {
1844 s = 0.5 * (x2 * x2 - x1 * x1) * y[n1];
1845 }
1846 return s;
1847 }
1848 if (xpower == 0) {
1849 s += (b2 - x1) * y[n1];
1850 } else {
1851 s += 0.5 * (b2 * b2 - x1 * x1) * y[n1];
1852 }
1853 }
1854 istart = n2;
1855 }
1856 if (x2 >= xmax) {
1857 x2 = xmax;
1858 iafterend = qi;
1859 } else {
1860 long n1, n2;
1861 T b1, b2;
1862 int i_ret = 0;
1863 i_ret = mesh.get_interval(x2, n1, b1, n2, b2);
1864 //Iprint2n(mcout, x2, i_ret);
1865 //Iprint4n(mcout, n1, b1, n2, b2);
1866 check_econd11(i_ret, != 1, mcerr);
1867 if (x2 - b1 > 0) {
1868 if (xpower == 0) {
1869 s += (x2 - b1) * y[n1];
1870 } else {
1871 s += 0.5 * (x2 * x2 - b1 * b1) * y[n1];
1872 }
1873 }
1874 iafterend = n1;
1875 }
1876 //Iprint2n(mcout, istart, iafterend);
1877 long i;
1878 double b;
1879 mesh.get_scoor(istart, b);
1880 if (xpower == 0) {
1881 for (i = istart; i < iafterend; i++) {
1882 double a = b;
1883 mesh.get_scoor(i + 1, b);
1884 s += (b - a) * y[i];
1885 }
1886 } else {
1887 for (i = istart; i < iafterend; i++) {
1888 double a = b;
1889 mesh.get_scoor(i + 1, b);
1890 s += 0.5 * (b * b - a * a) * y[i];
1891 }
1892 }
1893 //Iprintn(mcout, s);
1894
1895 //T t;
1896 return s;
1897
1898}
#define check_econd21(a, sign1_b1_sign0, sign2_b2, stream)
Definition: FunNameStack.h:428

Referenced by Heed::EnTransfCS::EnTransfCS(), and t_mean_step_ar().

◆ t_integ_straight_2point()

template<class T >
T t_integ_straight_2point ( x1,
y1,
x2,
y2,
xl,
xr,
int  xpower,
int  s_ban_neg 
)

Definition at line 2473 of file tline.h.

2477 {
2478 mfunname("double t_integ_straight_2point(...)");
2479 check_econd12(x1, ==, x2, mcerr);
2480
2481 T a = (y2 - y1) / (x2 - x1);
2482 T b = y1;
2483 T yl = a * (xl - x1) + b;
2484 T yr = a * (xr - x1) + b;
2485 if (s_ban_neg == 1) {
2486 if (yl <= 0.0 && yr <= 0.0) return 0.0;
2487 if (yl < 0.0 || yr < 0.0) {
2488 T xz = x1 - b / a;
2489 if (yl < 0.0) {
2490 xl = xz;
2491 yl = 0.0;
2492 } else {
2493 xr = xz;
2494 yr = 0.0;
2495 }
2496 }
2497 }
2498 T res;
2499 if (xpower == 0)
2500 res = 0.5 * a * (xr * xr - xl * xl) + (b - a * x1) * (xr - xl);
2501 else
2502 res = a * (xr * xr * xr - xl * xl * xl) / 3.0 +
2503 0.5 * (b - a * x1) * (xr * xr - xl * xl);
2504
2505 return res;
2506}

◆ t_integ_straight_point_ar()

template<class T , class D , class M >
T t_integ_straight_point_ar ( const M &  mesh,
const D &  y,
x1,
x2,
int  xpower,
int  s_ban_neg,
int  s_extrap_left,
left_bond,
int  s_extrap_right,
right_bond 
)

Definition at line 2678 of file tline.h.

2683 {
2684 mfunname("double t_integ_straight_point_ar(...)");
2685
2686 //mcout<<"Strart t_integ_straight_point_ar\n";
2687 check_econd21(xpower, != 0 &&, != 1, mcerr);
2688 check_econd12(x1, >, x2, mcerr);
2689 long qi = mesh.get_qi();
2690 check_econd12(qi, <, 1, mcerr);
2691 //if(x1 > x2) return 0.0;
2692 double xmin = mesh.get_xmin();
2693 double xmax = mesh.get_xmax();
2694 if (x2 <= xmin && s_extrap_left == 0) return 0.0;
2695 if (x1 >= xmax && s_extrap_right == 0) return 0.0;
2696 if (x2 <= left_bond) return 0.0;
2697 if (x1 >= right_bond) return 0.0;
2698 T s(0.0);
2699 if (x1 < left_bond) x1 = left_bond;
2700 if (x2 > right_bond) x2 = right_bond;
2701 if (x1 <= xmin && s_extrap_left == 0) x1 = xmin;
2702 if (x2 > xmax && s_extrap_left == 0) x2 = xmax;
2703 long np1, np2;
2704 T bp1, bp2;
2705 int i_ret = 0;
2706 // restore the interval in which x1 reside
2707 i_ret = mesh.get_interval_extrap(x1, np1, bp1, np2, bp2);
2708 // restore the x-coordinates of given points
2709 T xp1;
2710 mesh.get_scoor(np1, xp1);
2711 T xp2;
2712 mesh.get_scoor(np2, xp2);
2713 T res;
2714 T yp1 = y[np1];
2715 T yp2 = y[np2];
2716 if (i_ret == 2 || x2 <= xp2) // then all in one interval
2717 {
2718 res = t_integ_straight_2point<T>(xp1, yp1, xp2, yp2, x1, x2, xpower,
2719 s_ban_neg);
2720 } else {
2721 // integrate only till end of the current interval
2722 T x1i = x1;
2723 T x2i = xp2;
2724 res = t_integ_straight_2point<T>(xp1, yp1, xp2, yp2, x1i, x2i, xpower,
2725 s_ban_neg);
2726 //x2i = x1; // prepere for loop
2727 int s_stop = 0;
2728 do {
2729 np1 = np2;
2730 np2++;
2731 xp1 = xp2;
2732 mesh.get_scoor(np2, xp2);
2733 x1i = x2i;
2734 if (xp2 >= x2) {
2735 x2i = x2; // till end of integral
2736 s_stop = 1;
2737 } else {
2738 if (np2 == qi) // end of the mesh, but x2 is farther
2739 {
2740 x2i = x2; // till end of integral
2741 s_stop = 1;
2742 } else {
2743 x2i = xp2; // till end of current interval
2744 s_stop = 0;
2745 }
2746 }
2747 yp1 = yp2;
2748 yp2 = y[np2];
2749 res += t_integ_straight_2point<T>(xp1, yp1, xp2, yp2, x1i, x2i, xpower,
2750 s_ban_neg);
2751 //Iprint2n(mcout, xp1, xp2);
2752 //Iprint2n(mcout, x1i, x2i);
2753 //Iprint2n(mcout, yp1, yp2);
2754 //Iprint2n(mcout, res, s_stop);
2755
2756 } while (s_stop == 0);
2757 }
2758 return res;
2759}

Referenced by t_mean_straight_point_ar().

◆ t_mean_step_ar()

template<class T , class D , class M >
T t_mean_step_ar ( const M &  mesh,
const D &  y,
x1,
x2,
int &  s_err 
)

Definition at line 2425 of file tline.h.

2426 {
2427 mfunname("double t_mean_step_ar(...)");
2428 s_err = 0;
2429 T integ = t_integ_step_ar(mesh, y, x1, x2, 0);
2430 if (integ == 0) {
2431 s_err = 1;
2432 return 0.0;
2433 }
2434 T xinteg = t_integ_step_ar(mesh, y, x1, x2, 1);
2435 return xinteg / integ;
2436}
T t_integ_step_ar(const M &mesh, const D &y, T x1, T x2, int xpower)
Definition: tline.h:1810

◆ t_mean_straight_point_ar()

template<class T , class D , class M >
T t_mean_straight_point_ar ( const M &  mesh,
const D &  y,
x1,
x2,
int  s_extrap_left,
left_bond,
int  s_extrap_right,
right_bond,
int &  s_err 
)

Definition at line 2762 of file tline.h.

2765 {
2766 mfunname("double t_mean_straight_point_ar(...)");
2767 s_err = 0;
2768 T integ = t_integ_straight_point_ar(mesh, y, x1, x2, 0, 1, s_extrap_left,
2769 left_bond, s_extrap_right, right_bond);
2770 if (integ == 0) {
2771 s_err = 1;
2772 return 0.0;
2773 }
2774 T xinteg = t_integ_straight_point_ar(mesh, y, x1, x2, 1, 1, s_extrap_left,
2775 left_bond, s_extrap_right, right_bond);
2776 return xinteg / integ;
2777}
T t_integ_straight_point_ar(const M &mesh, const D &y, T x1, T x2, int xpower, int s_ban_neg, int s_extrap_left, T left_bond, int s_extrap_right, T right_bond)
Definition: tline.h:2678

◆ t_opposite_hisran_step_ar()

template<class T , class D , class M >
T t_opposite_hisran_step_ar ( const M &  mesh,
const D &  integ_y,
x 
)

Definition at line 2362 of file tline.h.

2362 {
2363 mfunname("double t_opposite_hisran_step_ar(...)");
2364
2365 //mcout<<"start t_step_integ_ar\n";
2366 //check_econd11(xpower , != 0 , mcerr);
2367 long qi = mesh.get_qi();
2368 long s_same = apeq_mant(integ_y[qi - 1], 1.0, 1.0e-12);
2369 check_econd11a(s_same, != 1.0, "integ_y[qi-1]=" << integ_y[qi - 1] << '\n',
2370 mcerr);
2371
2372 long n1;
2373 T b1;
2374 long n2;
2375 T b2;
2376
2377 mesh.get_interval_extrap(x, n1, b1, n2, b2);
2378 T y1;
2379 T y2;
2380 if (n1 == 0) {
2381 y1 = 0.0;
2382 y2 = integ_y[0];
2383 } else {
2384 y1 = integ_y[n1 - 1];
2385 y2 = integ_y[n2 - 1];
2386 }
2387 T res = y1 + ((y2 - y1) / (b2 - b1)) * (x - b1);
2388 return res;
2389}

◆ t_total_integ_step_ar() [1/3]

template<class T , class D , class M >
T t_total_integ_step_ar ( const M &  mesh,
const D &  y 
)

Definition at line 1987 of file tline.h.

1988 {
1989 mfunname("double t_total_integ_step_ar(...)");
1990
1991 //mcout<<"start t_step_integ_ar\n";
1992 long qi = mesh.get_qi();
1993 check_econd12(qi, <, 1, mcerr);
1994 //if(x1 > x2) return 0;
1995 long istart, iafterend; // indexes to sum total intervals
1996 T s(0);
1997 istart = 0;
1998 iafterend = qi;
1999 //Iprint2n(mcout, istart, iafterend);
2000 long i;
2001 double b;
2002 mesh.get_scoor(istart, b);
2003 for (i = istart; i < iafterend; i++) {
2004 double a = b;
2005 mesh.get_scoor(i + 1, b);
2006 s += (b - a) * y[i];
2007 }
2008
2009 //T t;
2010 return s;
2011
2012}

◆ t_total_integ_step_ar() [2/3]

template<class T , class D , class M1 , class M2 >
T t_total_integ_step_ar ( const M1 &  mesh1,
const M2 &  mesh2,
const D &  y 
)

Definition at line 2017 of file tline.h.

2019 {
2020 mfunname("double t_total_integ_step_ar(...)");
2021
2022 //mcout<<"start t_step_integ_ar\n";
2023 long qi1 = mesh1.get_qi();
2024 check_econd12(qi1, <, 1, mcerr);
2025 long qi2 = mesh2.get_qi();
2026 check_econd12(qi2, <, 1, mcerr);
2027 //if(x1 > x2) return 0;
2028 long istart1, iafterend1; // indexes to sum total intervals
2029 T s1(0);
2030 istart1 = 0;
2031 iafterend1 = qi1;
2032 //Iprint2n(mcout, istart, iafterend);
2033 long i1;
2034 double b1;
2035 mesh1.get_scoor(istart1, b1);
2036 for (i1 = istart1; i1 < iafterend1; i1++) {
2037 double a1 = b1;
2038 mesh1.get_scoor(i1 + 1, b1);
2039 // time to obtain integral by the second dimension
2040 //if(x1 > x2) return 0;
2041 long istart2, iafterend2; // indexes to sum total intervals
2042 T s2(0);
2043 istart2 = 0;
2044 iafterend2 = qi2;
2045 //Iprint2n(mcout, istart, iafterend);
2046 long i2;
2047 double b2;
2048 mesh2.get_scoor(istart2, b2);
2049 for (i2 = istart2; i2 < iafterend2; i2++) {
2050 double a2 = b2;
2051 mesh2.get_scoor(i2 + 1, b2);
2052 s2 += (b2 - a2) * y[i1][i2];
2053 }
2054 // OK, integral = s2
2055 s1 += (b1 - a1) * s2;
2056 }
2057
2058 //T t;
2059 return s1;
2060
2061}

◆ t_total_integ_step_ar() [3/3]

template<class T , class M1 , class M2 >
T t_total_integ_step_ar ( const M1 &  mesh1,
const M2 &  mesh2,
const DynArr< T > &  y 
)

Definition at line 2066 of file tline.h.

2068 {
2069 mfunname("double t_total_integ_step_ar(...)");
2070
2071 //mcout<<"start t_step_integ_ar\n";
2072 long qi1 = mesh1.get_qi();
2073 check_econd12(qi1, <, 1, mcerr);
2074 check_econd12(qi1, !=, y.get_qel()[0], mcerr);
2075 long qi2 = mesh2.get_qi();
2076 check_econd12(qi2, <, 1, mcerr);
2077 check_econd12(qi2, !=, y.get_qel()[1], mcerr);
2078 //if(x1 > x2) return 0;
2079 long istart1, iafterend1; // indexes to sum total intervals
2080 T s1(0);
2081 istart1 = 0;
2082 iafterend1 = qi1;
2083 //Iprint2n(mcout, istart, iafterend);
2084 long i1;
2085 double b1;
2086 mesh1.get_scoor(istart1, b1);
2087 for (i1 = istart1; i1 < iafterend1; i1++) {
2088 double a1 = b1;
2089 mesh1.get_scoor(i1 + 1, b1);
2090
2091 // time to obtain integral by the second dimension
2092
2093 //if(x1 > x2) return 0.0;
2094 long istart2, iafterend2; // indexes to sum total intervals
2095 T s2(0.0);
2096 istart2 = 0;
2097 iafterend2 = qi2;
2098 //Iprint2n(mcout, istart, iafterend);
2099 long i2;
2100 double b2;
2101 mesh2.get_scoor(istart2, b2);
2102 for (i2 = istart2; i2 < iafterend2; i2++) {
2103 double a2 = b2;
2104 mesh2.get_scoor(i2 + 1, b2);
2105 s2 += (b2 - a2) * y.acu(i1, i2);
2106 }
2107
2108 // OK, integral = s2
2109
2110 s1 += (b1 - a1) * s2;
2111 }
2112
2113 //T t;
2114 return s1;
2115
2116}
T & acu(long i1)
Definition: AbsArr.h:2098
const DynLinArr< long > & get_qel(void) const
Definition: AbsArr.h:2548

◆ t_value_exp_2point()

template<class T >
T t_value_exp_2point ( x1,
y1,
x2,
y2,
x 
)

Definition at line 2643 of file tline.h.

2643 {
2644 mfunname("double t_value_exp_2point(...)");
2645
2646 check_econd11(y1, <= 0.0, mcerr);
2647 check_econd11(y2, <= 0.0, mcerr);
2648 check_econd12(y1, ==, y2, mcerr);
2649 check_econd12(x1, ==, x2, mcerr);
2650 T res;
2651
2652 T a = log(y1 / y2) / (x1 - x2);
2653 T c = y1 / exp(a * x1);
2654 //check_econd11(pw , == -1.0 , mcerr);
2655 res = c * exp(a * x);
2656 ;
2657 return res;
2658}
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:376

◆ t_value_generic_point_ar()

template<class T , class D , class M >
T t_value_generic_point_ar ( const M &  mesh,
const D &  y,
T(*)(T xp1, T yp1, T xp2, T yp2, T xmin, T xmax, T x)  funval,
x,
int  s_extrap_left,
left_bond,
int  s_extrap_right,
right_bond 
)

Definition at line 2572 of file tline.h.

2577 {
2578 // 0 - not include, 1 - include
2579 mfunname("double t_value_generic_point_ar(...)");
2580 double xmin = mesh.get_xmin();
2581 double xmax = mesh.get_xmax();
2582 //Iprint3n(mcout, x, xmin, xmax);
2583 if (x < left_bond) return 0.0;
2584 if (x > right_bond) return 0.0;
2585 if (x < xmin && s_extrap_left == 0) return 0.0;
2586 if (x > xmax && s_extrap_right == 0) return 0.0;
2587 long n1, n2;
2588 T b1, b2;
2589 mesh.get_interval_extrap(x, n1, b1, n2, b2);
2590 T x1;
2591 mesh.get_scoor(n1, x1);
2592 T x2;
2593 mesh.get_scoor(n2, x2);
2594 return funval(x1, y[n1], x2, y[n2], left_bond, right_bond, x);
2595}

Referenced by Heed::SimpleTablePhotoAbsCS::get_CS().

◆ t_value_power_2point()

template<class T >
T t_value_power_2point ( x1,
y1,
x2,
y2,
x 
)

Definition at line 2603 of file tline.h.

2603 {
2604 mfunname("double t_value_power_2point(...)");
2605
2606 check_econd11(y1, <= 0.0, mcerr);
2607 check_econd11(y2, <= 0.0, mcerr);
2608 check_econd12(y1, ==, y2, mcerr);
2609 check_econd12(x1, ==, x2, mcerr);
2610 T res = y1;
2611 if (x1 <= 0.0 && x2 >= 0.0) {
2612 mcerr << "T t_value_power_2point(...): \n";
2613 mcerr << "x's are of different sign, power cannot be drawn\n";
2614 spexit(mcerr);
2615 } else {
2616 T pw = log(y1 / y2) / log(x1 / x2);
2617 //check_econd11(pw , == -1.0 , mcerr);
2618 res = y1 * pow(x, pw) / pow(x1, pw);
2619 }
2620 return res;
2621}

◆ t_value_power_extended_2point()

template<class T >
T t_value_power_extended_2point ( x1,
y1,
x2,
y2,
x 
)

Definition at line 2625 of file tline.h.

2625 {
2626 mfunname("double t_value_power_2point(...)");
2627
2628 check_econd11(y1, <= 0.0, mcerr);
2629 check_econd11(y2, <= 0.0, mcerr);
2630 check_econd12(y1, ==, y2, mcerr);
2631 check_econd12(x1, ==, x2, mcerr);
2632 T res;
2633 if (x1 <= 0.0 && x2 >= 0.0) {
2634 res = y1 + (x - x1) * (y2 - y1) / (x2 - x1);
2635 } else {
2636 T pw = log(y1 / y2) / log(x1 / x2);
2637 //check_econd11(pw , == -1.0 , mcerr);
2638 res = y1 * pow(x, pw) / pow(x1, pw);
2639 }
2640 return res;
2641}

◆ t_value_step_ar() [1/2]

template<class T , class D , class M >
T t_value_step_ar ( const M &  mesh,
const D &  y,
x,
int  s_include_last_point = 0 
)

Definition at line 1649 of file tline.h.

1652 {
1653 mfunname("double t_value_step_ar(...)");
1654 double xmin = mesh.get_xmin();
1655 double xmax = mesh.get_xmax();
1656 //Iprint3n(mcout, x, xmin, xmax);
1657 if (x < xmin) return 0;
1658 if (s_include_last_point == 0) {
1659 if (x >= xmax) return 0;
1660 } else {
1661 if (x > xmax) return 0;
1662 }
1663 long n1, n2;
1664 T b1, b2;
1665 int i_ret = 0;
1666 i_ret = mesh.get_interval(x, n1, b1, n2, b2);
1667 check_econd11(i_ret, != 1, mcerr);
1668 return y[n1];
1669}

◆ t_value_step_ar() [2/2]

template<class T , class D , class M1 , class M2 >
T t_value_step_ar ( const M1 &  mesh1,
const M2 &  mesh2,
const D &  y,
x1,
x2,
int  s_include_last_point = 0 
)

Definition at line 1673 of file tline.h.

1677 {
1678 mfunname("double t_value_step_ar(...)");
1679 double x1min = mesh1.get_xmin();
1680 double x1max = mesh1.get_xmax();
1681 //Iprint3n(mcout, x, xmin, xmax);
1682 if (x1 < x1min) return 0;
1683 if (s_include_last_point == 0) {
1684 if (x1 >= x1max) return 0;
1685 } else {
1686 if (x1 > x1max) return 0;
1687 }
1688 double x2min = mesh2.get_xmin();
1689 double x2max = mesh2.get_xmax();
1690 //Iprint3n(mcout, x, xmin, xmax);
1691 if (x2 < x2min) return 0;
1692 if (s_include_last_point == 0) {
1693 if (x2 >= x2max) return 0;
1694 } else {
1695 if (x2 > x2max) return 0;
1696 }
1697 long n11, n12;
1698 long n21, n22;
1699 T b1, b2;
1700 int i_ret = 0;
1701
1702 i_ret = mesh1.get_interval(x1, n11, b1, n12, b2);
1703 check_econd11(i_ret, != 1, mcerr);
1704
1705 i_ret = mesh2.get_interval(x2, n21, b1, n22, b2);
1706
1707 check_econd11(i_ret, != 1, mcerr);
1708 return y[n11][n21];
1709}

◆ t_value_straight_2point()

template<class T >
T t_value_straight_2point ( x1,
y1,
x2,
y2,
x,
int  s_ban_neg 
)

Definition at line 2439 of file tline.h.

2439 {
2440 mfunname("double t_value_straight_2point(...)");
2441 check_econd12(x1, ==, x2, mcerr);
2442
2443 T a = (y2 - y1) / (x2 - x1);
2444 // Less numerical precision
2445 //T b = y[n1];
2446 //T res = a * ( x - x1) + b;
2447 // More numerical precision (although it is not proved),
2448 // starting from what is closer to x
2449 T res;
2450 T dx1 = x - x1;
2451 T adx1 = (dx1 > 0) ? dx1 : -dx1; // absolute value
2452 //if(dx1 > 0)
2453 // adx1 = dx1;
2454 //else
2455 // adx1 = -dx1;
2456 T dx2 = x - x2;
2457 T adx2 = (dx2 > 0) ? dx2 : -dx2; // absolute value
2458 //if(dx2 > 0)
2459 // adx2 = dx2;
2460 //else
2461 // adx2 = -dx2;
2462 if (adx1 < adx2) // x is closer to x1
2463 {
2464 res = a * dx1 + y1;
2465 } else {
2466 res = a * dx2 + y2;
2467 }
2468 if (s_ban_neg == 1 && res < 0.0) res = 0.0;
2469 return res;
2470}

Referenced by t_value_straight_point_ar(), and t_width_at_hheight_step_ar().

◆ t_value_straight_point_ar()

template<class T , class D , class M >
T t_value_straight_point_ar ( const M &  mesh,
const D &  y,
x,
int  s_ban_neg,
int  s_extrap_left,
left_bond,
int  s_extrap_right,
right_bond 
)

Definition at line 2511 of file tline.h.

2514 {
2515 // 0 - not include, 1 - include
2516 mfunname("double t_value_straight_point_ar(...)");
2517 double xmin = mesh.get_xmin();
2518 double xmax = mesh.get_xmax();
2519 //Iprint3n(mcout, x, xmin, xmax);
2520 if (x < left_bond) return 0.0;
2521 if (x > right_bond) return 0.0;
2522 if (x < xmin && s_extrap_left == 0) return 0.0;
2523 if (x > xmax && s_extrap_right == 0) return 0.0;
2524 long n1, n2;
2525 T b1, b2;
2526 mesh.get_interval_extrap(x, n1, b1, n2, b2);
2527 T x1;
2528 mesh.get_scoor(n1, x1);
2529 T x2;
2530 mesh.get_scoor(n2, x2);
2531 return t_value_straight_2point(x1, y[n1], x2, y[n2], x, s_ban_neg);
2532 /*
2533 T a = (y[n2] - y[n1])/(x2 - x1);
2534 // Less numerical precision
2535 //T b = y[n1];
2536 //T res = a * ( x - x1) + b;
2537 // More numerical precision (although it is not proved),
2538 // starting from what is closer to x
2539 T res;
2540 T dx1 = x - x1;
2541 T adx1 = (dx1 > 0) ? dx1 : -dx1; // absolute value
2542 //if(dx1 > 0)
2543 // adx1 = dx1;
2544 //else
2545 // adx1 = -dx1;
2546 T dx2 = x - x2;
2547 T adx2 = (dx2 > 0) ? dx2 : -dx2; // absolute value
2548 //if(dx2 > 0)
2549 // adx2 = dx2;
2550 //else
2551 // adx2 = -dx2;
2552 if(adx1 < adx2) // x is closer to x1
2553 {
2554 res = a * dx1 + y[n1];
2555 }
2556 else
2557 {
2558 res = a * dx2 + y[n2];
2559 }
2560 //T res = a * ( x - x1) + b;
2561 if(s_ban_neg == 1 && res < 0.0) res = 0.0;
2562 //Iprint2n(mcout, x, res);
2563 //Iprint2n(mcout, x1, x2);
2564 //Iprint2n(mcout, y[n1], y[n2]);
2565 //Iprintn(mcout, a);
2566 return res;
2567 */
2568}
T t_value_straight_2point(T x1, T y1, T x2, T y2, T x, int s_ban_neg)
Definition: tline.h:2439

Referenced by Heed::HeedMatterDef::replace_epsi12().

◆ t_width_at_hheight_step_ar()

template<class T , class D , class M >
T t_width_at_hheight_step_ar ( const M &  mesh,
const D &  y 
)

Definition at line 2874 of file tline.h.

2874 {
2875 // 0 - not include, 1 - include
2876 mfunname("double t_width_at_hheight_step_ar(...)");
2877 //mcout<<"t_width_at_hheight_step_ar is started\n";
2878 long qi = mesh.get_qi();
2879 long n;
2880 T ymax = 0;
2881 long nmax;
2882 for (n = 0; n < qi; ++n) {
2883 if (y[n] > ymax) {
2884 check_econd11(y[n], < 0.0, mcerr);
2885 ymax = y[n];
2886 nmax = n;
2887 }
2888 }
2889 //Iprint2n(mcout, ymax, nmax);
2890 if (ymax == 0) return 0;
2891 T ylev = ymax / 2.0;
2892 T s2 = 0;
2893 long q = 0;
2894 for (n = nmax; n < qi; n++) {
2895
2896 if (y[n] > ylev && y[n + 1] <= ylev) {
2897 T x1, x2;
2898 mesh.get_interval(n, x1, x2);
2899 T step1, step2;
2900 mesh.get_step(n, step1);
2901 mesh.get_step(n + 1, step2);
2902 step1 = step1 / 2.0;
2903 step2 = step2 / 2.0;
2904 s2 += t_value_straight_2point(y[n], x1 + step1, y[n + 1], x2 + step2,
2905 ylev, 0);
2906 //Iprint2n(mcout, x1, x2);
2907 //Iprint2n(mcout, x1+step1, x2+step2);
2908 //Iprint2n(mcout, y[n], y[n+1]);
2909 //Iprint2n(mcout, n, t_value_straight_2point(y[n], x1+step1, y[n+1],
2910 //x2+step2, ylev, 0));
2911 q++;
2912 }
2913 }
2914 check_econd11(q, <= 0, mcerr);
2915 s2 = s2 / q;
2916 T s1 = 0;
2917 q = 0;
2918 for (n = nmax; n >= 0; n--) {
2919 if (y[n] > ylev && y[n - 1] <= ylev) {
2920 T x1, x2;
2921 mesh.get_interval(n - 1, x1, x2);
2922 T step1, step2;
2923 mesh.get_step(n - 1, step1);
2924 mesh.get_step(n, step2);
2925 step1 = step1 / 2.0;
2926 step2 = step2 / 2.0;
2927 s1 += t_value_straight_2point(y[n - 1], x1 + step1, y[n], x2 + step2,
2928 ylev, 0);
2929 //Iprint2n(mcout, x1, x2);
2930 //Iprint2n(mcout, x1+step1, x2+step2);
2931 //Iprint2n(mcout, y[n-1], y[n]);
2932 //Iprint2n(mcout, n, t_value_straight_2point(y[n-1], x1+step1, y[n],
2933 //x2+step2, ylev, 0));
2934 q++;
2935 }
2936 }
2937 check_econd11(q, <= 0, mcerr);
2938 s1 = s1 / q;
2939 //Iprint3n(mcout, s1, s2, s2 - s1);
2940 return s2 - s1;
2941}

◆ verify_types()

template<class T , class B >
int verify_types ( const T *  ths,
const B *  b2,
const T **  tb2 
)

Definition at line 1179 of file tline.h.

1179 {
1180 mfunnamep("template<class T> verify_types(const T* ths, const B* b2, const "
1181 "T* tb2)");
1182 //mcout<<"verify_types is working\n";
1183 //Iprintn(mcerr, typeid(T).name());
1184 //Iprintn(mcerr, typeid(B).name());
1185 if (typeid(*b2) != typeid(*ths)) return 0;
1186 if (typeid(*ths) != typeid(T)) {
1187 funnw.ehdr(cerr);
1188 mcerr << "Operator == for class T is called for first object of different "
1189 "class\n";
1190 mcerr << "This may mean that this operator is missed in one of derived "
1191 "classes\n";
1192 Iprintn(mcerr, typeid(*ths).name());
1193 spexit(mcerr);
1194 }
1195 *tb2 = dynamic_cast<const T*>(b2);
1196 check_econd11(*tb2, == NULL, mcerr);
1197 return 1;
1198}

Referenced by VirtEqualStepCoorMesh< T >::apeq_mant(), VirtCopiedPointCoorMesh< T >::apeq_mant(), VirtEqualStepCoorMesh< T >::operator==(), and VirtCopiedPointCoorMesh< T >::operator==().