KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVIDZAFromZGrid.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Tue Mar 8 10:00:16 2016
2 //Author: Diego Gruyer
3 
4 #include "KVIDZAFromZGrid.h"
5 #include "TMultiGraph.h"
6 #include "KVIDZALine.h"
7 #include "TCanvas.h"
8 
12 
13 
14 
15 
20 {
21  // Default constructor
22  // Grid is declared as a 'ZOnlyGrid' by default (this is internal mechanics)
23 
24  init();
25  fTables.SetOwner(kTRUE);
26  fFits.SetOwner();
27  SetOnlyZId();
28  fIgnoreMassID = false;
29 }
30 
31 
32 
33 
41 
43 {
44  // This method copies the current state of 'this' object into 'obj'
45  // You should add here any member variables, for example:
46  // (supposing a member variable KVIDZAFromZGrid::fToto)
47  // CastedObj.fToto = fToto;
48  // or
49  // CastedObj.SetToto( GetToto() );
50 
51  KVIDZAGrid::Copy(obj);
53  i.fZmaxInt = fZmaxInt;
54  i.fPIDRange = fPIDRange;
55  i.fZminInt = fZminInt;
56  fTables.Copy(i.fTables);
59 }
60 
61 
62 
63 
72 
73 void KVIDZAFromZGrid::ReadFromAsciiFile(std::ifstream& gridfile)
74 {
75 // fPIDRange = kFALSE;
76 // KVIDGraph::ReadFromAsciiFile(gridfile);
77 // if (GetParameters()->HasParameter("PIDRANGE")) {
78 // fPIDRange = kTRUE;
79 // fZmaxInt = GetParameters()->GetIntValue("PIDRANGE");
80 // LoadPIDRanges();
81 // }
82 
83  fPIDRange = kFALSE;
84 // fHasMassCut = kFALSE;
86 
87 // if (GetIdentifier("MassID")) fHasMassCut = kTRUE;
88 
89  if (GetParameters()->HasParameter("PIDRANGE")) {
90  fPIDRange = kTRUE;
91  TString pidrange = GetParameters()->GetStringValue("PIDRANGE");
92  if (pidrange.Contains("-")) {
93  TString min = (pidrange(0, pidrange.Index("-")));
94  fZminInt = min.Atoi();
95  min = (pidrange(pidrange.Index("-") + 1, pidrange.Length()));
96  fZmaxInt = min.Atoi();
97  }
98  else {
99  fZminInt = 1;
100  fZmaxInt = pidrange.Atoi();
101  }
102  LoadPIDRanges();
103  }
104  // if <PARAMETER> IgnoreMassID=1 appears in file, we are only using the PID intervals to clean
105  // up a messy de-e plot, not to give mass identification. particles will only be identified in Z.
106  if (GetParameters()->HasParameter("IgnoreMassID") && GetParameters()->GetIntValue("IgnoreMassID") == 1)
107  fIgnoreMassID = true;
108  else
109  fIgnoreMassID = false;
110 }
111 
112 
113 
115 
116 void KVIDZAFromZGrid::WriteToAsciiFile(std::ofstream& gridfile)
117 {
118  ExportToGrid();
119  KVIDGraph::WriteToAsciiFile(gridfile);
120 }
121 
122 
123 
125 
127 {
128  fZminInt = 100000;
129  fZmaxInt = 0;
130 
131  KVIDentifier* id = 0;
132  TIter it(GetIdentifiers());
133  while ((id = (KVIDentifier*)it())) {
134  int zz = id->GetZ();
135  if (!GetParameters()->HasParameter(Form("PIDRANGE%d", zz))) continue;
136  KVString mes = GetParameters()->GetStringValue(Form("PIDRANGE%d", zz));
137  if (mes.IsWhitespace()) continue;
138  int type = (mes.Contains(",") ? 2 : 1);
139  interval_set* itv = new interval_set(zz, type);
140  itv->SetName(GetName());
141  mes.Begin("|");
142  while (!mes.End()) {
143  KVString tmp = mes.Next();
144  tmp.Begin(":");
145  int aa = tmp.Next().Atoi();
146  KVString val = tmp.Next();
147  double pidmin, pidmax, pid;
148  if (type == 1) itv->add(aa, val.Atof());
149  else if (type == 2) {
150  val.Begin(",");
151  pidmin = val.Next().Atof();
152  pid = val.Next().Atof();
153  pidmax = val.Next().Atof();
154  itv->add(aa, pid, pidmin, pidmax);
155 // itv->add(aa, pid, pid-0.02, pid+0.02);
156  }
157  }
158  if (zz < fZminInt) fZminInt = zz;
159  if (zz > fZmaxInt) fZmaxInt = zz;
160  fTables.Add(itv);
161  }
162  fPIDRange = kTRUE;
163  // PrintPIDLimits();
164 }
165 
166 
167 
169 
171 {
172  fTables.Clear("all");
173  fPIDRange = kFALSE;
174 }
175 
176 
177 
179 
181 {
182  fTables.Clear("all");
183  LoadPIDRanges();
184 }
185 
186 
187 
189 
191 {
192  interval_set* itv = 0;
193  TIter it(&fTables);
194  while ((itv = (interval_set*)it())) if (itv->GetZ() == zint) return itv;
195  return 0;
196 }
197 
198 
199 
206 
208 {
209 // ((interval_set*)fTables.At(12))->fIntervals.ls();
210 
211 // for (int zz = fZminInt; zz <= fZmaxInt; zz++) {
212 // Info("PrintPIDLimits", "Z=%2d [%.4lf %.4lf]", zz, ((interval_set*)fTables.At(zz - fZminInt))->fPIDmins.at(0),
213 // ((interval_set*)fTables.At(zz - fZminInt))->fPIDmaxs.at(((interval_set*)fTables.At(zz - fZminInt))->fNPIDs - 1));
214 // }
215 }
216 
217 
218 
220 
222 {
223  if (GetParameters()->HasParameter("PIDRANGE")) GetParameters()->RemoveParameter("PIDRANGE");
224  for (int ii = 1; ii < 30; ii++) {
225  if (GetParameters()->HasParameter(Form("PIDRANGE%d", ii))) GetParameters()->RemoveParameter(Form("PIDRANGE%d", ii));
226  }
227 }
228 
229 
230 
237 
238 int KVIDZAFromZGrid::is_inside(double pid) const
239 {
240  // Look for a set of mass-interval definitions in which the given PID
241  // falls (PID from linearisation of Z identification).
242  // In principle this should be the set corresponding to Z=nint(PID),
243  // but if not Z+/-1 are also tried.
244  // Returns the value of Z for the set found (or 0 if no set found)
245 
246  int zint = TMath::Nint(pid);
247  interval_set* it = GetIntervalSet(zint);
248  if (it) {
249  if (it->is_inside(pid)) return zint;
250  else if (it->is_above(pid)) {
251 
252  it = GetIntervalSet(zint + 1);
253  if (it && it->is_inside(pid)) return zint + 1;
254  else return 0;
255  }
256  else {
257  it = GetIntervalSet(zint - 1);
258  if (it && it->is_inside(pid)) return zint - 1;
259  else return 0;
260  }
261  }
262  else return 0;
263 }
264 
265 
266 
273 
275 {
276  // General initialisation method for identification grid.
277  // This method MUST be called once before using the grid for identifications.
278  // The ID lines are sorted.
279  // The natural line widths of all ID lines are calculated.
280  // The line with the largest Z (Zmax line) is found.
281 
282  SetOnlyZId();
284  // Zmax should be Z of last line in sorted list
285 
286  TIter it(GetIdentifiers());
287  KVIDentifier* id = 0;
288 
289  fZMax = 0;
290  while ((id = (KVIDentifier*)it())) {
291  if (!id->InheritsFrom("KVIDZALine")) continue;
292  int zz = ((KVIDZALine*)id)->GetZ();
293  if (zz > fZMax) fZMax = zz;
294  }
295 
296  // set to true if grid has a limited region for mass identification, indicated by an info "MassID"
297  fHasMassIDRegion = (GetInfos()->FindObject("MassID") != nullptr);
298 
299  // set up mass fits (if any)
300  fFits.Clear();
301  if (GetParameters()->HasStringParameter("MASSFITS")) {
302  KVNumberList zlist(GetParameters()->GetStringValue("MASSFITS"));
303  for (auto z : zlist) {
304  auto massfit = GetParameters()->GetTStringValue(Form("MASSFIT_%d", z));
305  massfit.ReplaceAll(":", "=");
306  KVNameValueList fitparams;
307  fitparams.Set(massfit);
308  fFits.Add(new KVMultiGaussIsotopeFit(z, fitparams));
309  }
310  }
311 }
312 
313 
314 
316 
318 {
319  double P;
320  auto A = fitfunc->GetA(idr->PID, P);
321  if (A) {
322  idr->A = A;
323  idr->PID = fitfunc->GetInterpolatedA(idr->PID);
324  if (P > 0.5) idr->IDquality = KVIDZAGrid::kICODE0; // probability of A is >50%
325  else idr->IDquality = KVIDZAGrid::kICODE3;// OK, slight ambiguity of A
326  }
327  else {
328  // returned A=0 => background noise
330  return false;
331  }
332  return true;
333 }
334 
335 
336 
348 
350 {
351  // Fill the KVIdentificationResult object with the results of identification for point (x,y)
352  // corresponding to some physically measured quantities related to a reconstructed nucleus.
353  //
354  // If identification is successful, idr->IDOK = true.
355  // In this case, idr->Aident and idr->Zident indicate whether isotopic or only Z identification
356  // was acheived.
357  //
358  // In case of unsuccessful identification, idr->IDOK = false,
359  // BUT idr->Zident and/or idr->Aident may be true: this is to indicate which kind of
360  // identification was attempted but failed (this changes the meaning of the quality code)
361 
362  idr->Aident = idr->Zident = kFALSE;
363 
364  KVIDZAGrid::Identify(x, y, idr);
365  idr->Zident = kTRUE; // meaning Z identification was attempted, even if it failed
366  if (!idr->IDOK) return;
367 
368  bool have_pid_range_for_Z = fPIDRange && (idr->Z <= fZmaxInt) && (idr->Z > fZminInt - 1);
369  auto mass_fit_for_Z = GetMultiGaussFit(idr->Z);
370  bool have_mass_fit_for_Z = (mass_fit_for_Z != nullptr);
371  bool mass_id_success = false;
372 
373  if ((have_mass_fit_for_Z || have_pid_range_for_Z)
374  && (!fHasMassIDRegion || idr->HasFlag(GetName(), "MassID"))) { // if a mass ID region is defined, we must be inside it
375  // try mass identification
376  if (have_mass_fit_for_Z)
377  mass_id_success = MassIdentificationFromMultiGaussFit(mass_fit_for_Z, idr);
378  else
379  mass_id_success = (DeduceAfromPID(idr) > 0);
380  if (mass_id_success) {
381  // mass identification was at least attempted
382  // make sure grid's quality code is consistent with KVIdentificationResult
383  const_cast<KVIDZAFromZGrid*>(this)->fICode = idr->IDquality;
384  idr->Aident = kTRUE; // meaning A identification was attempted, even if it failed
385  }
386  else {
387  // the pid falls outside of any mass ranges for a Z which has assigned isotopes
388  // therefore although the Z identification was good, we cannot consider this
389  // particle to be identified
390  const_cast<KVIDZAFromZGrid*>(this)->fICode = kICODE4;
391  idr->IDquality = fICode; // otherwise identfication result quality code is not coherent with comment (see below)
392  }
393  idr->IDOK = (fICode < kICODE4);
394  }
395 
396  // ignore isotopic successful isotopic identification if fIgnoreMassID=true
397  if (fIgnoreMassID && idr->IDOK && idr->Aident) idr->Aident = false;
398 
399  // set comments in identification result
400  switch (fICode) {
401  case kICODE0:
402  idr->SetComment("ok");
403  break;
404  case kICODE1:
405  if (mass_id_success) idr->SetComment("slight ambiguity of A, which could be larger");
406  else idr->SetComment("slight ambiguity of Z, which could be larger");
407  break;
408  case kICODE2:
409  if (mass_id_success) idr->SetComment("slight ambiguity of A, which could be smaller");
410  else idr->SetComment("slight ambiguity of Z, which could be smaller");
411  break;
412  case kICODE3:
413  if (mass_id_success) idr->SetComment("slight ambiguity of A, which could be larger or smaller");
414  else idr->SetComment("slight ambiguity of Z, which could be larger or smaller");
415  break;
416  case kICODE4:
417  if (mass_id_success) idr->SetComment("point is outside of mass identification range");
418  else idr->SetComment("point is in between two lines of different Z, too far from either to be considered well-identified");
419  break;
420  case kICODE5:
421  if (mass_id_success) idr->SetComment("point is in between two isotopes A & A+2 (e.g. 5He, 8Be, 9B)");
422  else idr->SetComment("point is in between two lines of different Z, too far from either to be considered well-identified");
423  break;
424  case kICODE6:
425  idr->SetComment("(x,y) is below first line in grid");
426  break;
427  case kICODE7:
428  idr->SetComment("(x,y) is above last line in grid");
429  break;
430  default:
431  idr->SetComment("no identification: (x,y) out of range covered by grid");
432  }
433 }
434 
435 
436 
437 
443 
445 {
446  // First look for a set of mass intervals in which the PID of the identification result falls,
447  // if there is one (see KVIDZAFromZGrid::is_inside).
448  // If an interval set is found for a Z different to the original identification, idr->Z is changed.
449  // Then call interval_set::eval for the mass interval for this Z.
450 
451  int zint = is_inside(idr->PID);
452  if (!zint) return -1;
453  if (zint != idr->Z) idr->Z = zint;
454 
455  double res = 0.;
456  interval_set* it = GetIntervalSet(zint);
457  if (it) res = it->eval(idr);
458  return res;
459 }
460 
461 
462 
463 
465 
467 {
469  KVNumberList pids;
470  interval_set* itvs = 0;
471  TIter npid(GetIntervalSets());
472  while ((itvs = (interval_set*)npid())) {
473  if (!itvs->GetNPID()) continue;
474  pids.Add(itvs->GetZ());
475  }
476  GetParameters()->SetValue("PIDRANGE", pids.AsString());
477 
478  itvs = 0;
479  TIter next(GetIntervalSets());
480  while ((itvs = (interval_set*)next())) {
481  if (!itvs->GetNPID()) continue;
482  KVString par = Form("PIDRANGE%d", itvs->GetZ());
483  KVString val = "";
484  interval* itv = 0;
485  TIter ni(itvs->GetIntervals());
486  while ((itv = (interval*)ni())) {
487  val += Form("%d:%lf,%lf,%lf|", itv->GetA(), itv->GetPIDmin(), itv->GetPID(), itv->GetPIDmax());
488  }
489  val.Remove(val.Length() - 1);
490  GetParameters()->SetValue(par.Data(), val.Data());
491  }
492 }
493 
494 
495 
496 
497 
499 
501 {
502  double pid = idr->PID;
503  if (pid < 0.5) return 0.;
504  // calculate interpolated mass from PID
505  double res = fPIDs.Eval(pid);
506  int ares = 0;
507 
509 
510  // look for mass interval PID is in
511  // in case it falls between two intervals remember also the interval
512  // immediately to the left & right of the PID
513  interval* left_int(nullptr), *right_int(nullptr);
514  interval* inter;
515  TIter it(&fIntervals);
516  while ((inter = (interval*)it())) {
517  if (inter->is_inside(pid)) {
518  ares = inter->GetA();
519  break;
520  }
521  else if (inter->is_left_of(pid)) {
522  left_int = inter;
523  }
524  else if (!right_int && inter->is_right_of(pid)) {
525  right_int = inter;
526  }
527  }
528  if (ares != 0) {
529  // the PID is inside a defined mass interval
530  idr->A = ares;
531  idr->PID = res;
533  }
534  else {
535  // the PID is not inside a defined mass interval
536  //
537  // * if it is in between two consecutive masses i.e. A and A+1 then it is
538  // Z- and A-identified with a slight ambiguity of A
539  // * if it is in between two non-consecutive masses i.e. A and A+2 then it
540  // is not identified (e.g. 5He, 8Be, 9B)
541  if (!right_int || !left_int) {
542  // case where no left or right interval were found
543  // to prevent from crashes but should not appen
544  idr->A = ares;
545  idr->PID = res;
547  }
548  else {
549  int dA = right_int->GetA() - left_int->GetA();
550  if (dA == 1) {
551  // OK, slight ambiguity of A
552  ares = TMath::Nint(res);
553  idr->A = ares;
554  idr->PID = res;
556  }
557  else {
558  // in a hole where no isotopes should be (e.g. 5He, 8Be, 9B)
559  idr->A = ares;
560  idr->PID = res;
562  }
563  }
564  }
565  }
566  else {
567  ares = TMath::Nint(res);
568  idr->A = ares;
569  idr->PID = res;
570  if (ares > fPIDs.GetX()[0] && ares < fPIDs.GetX()[fNPIDs - 1]) {
572  }
573  else {
575  }
576  }
577  return res;
578 }
579 
580 
581 
583 
584 bool interval_set::is_inside(double pid)
585 {
586  if (fType != KVIDZAFromZGrid::kIntType) return kTRUE;
587 
588 // Info("is_inside","min: %d max:%d npids:%d", ((interval*)fIntervals.At(0))->GetA(), ((interval*)fIntervals.At(fNPIDs-1))->GetA(), fNPIDs);
589 
590  if (pid > ((interval*)fIntervals.At(0))->GetPIDmin() && pid < ((interval*)fIntervals.At(fNPIDs - 1))->GetPIDmax()) return kTRUE;
591  else return kFALSE;
592 }
593 
594 
595 
597 
598 bool interval_set::is_above(double pid)
599 {
600  if (fType != KVIDZAFromZGrid::kIntType) return kTRUE;
601 
602  if (pid > ((interval*)fIntervals.At(fNPIDs - 1))->GetPIDmax()) return kTRUE;
603  else return kFALSE;
604 }
605 
606 
607 
608 
610 
612 {
613  if (!GetNPID()) return "-";
614  KVNumberList alist;
615  for (int ii = 0; ii < GetNPID(); ii++) alist.Add(((interval*)fIntervals.At(ii))->GetA());
616  return alist.AsString();
617 }
618 
619 
620 
622 
623 interval_set::interval_set(int zz, int type)
624 {
625  fType = type;
626  fZ = zz;
627  fNPIDs = 0;
629 }
630 
631 
632 
638 
639 void interval_set::add(int aa, double pid, double pidmin, double pidmax)
640 {
641 // if (fNPIDs && pid < fPIDs.GetX()[fNPIDs - 1]) {
642 // Error("add", "Please give me peaks in the right order for Z=%d and A=%d...", fZ, aa);
643 // return;
644 // }
645  if (fType == KVIDZAFromZGrid::kIntType && !(pid > pidmin && pid < pidmax)) {
646  Error("add", "Wrong interval for Z=%d and A=%d: [%.4lf %.4lf %.4lf] (%s)", fZ, aa, pidmin, pid, pidmax, GetName());
647  return;
648  }
649 
650  fPIDs.SetPoint(fNPIDs, pid, aa);
652  if (pid) fIntervals.AddLast(new interval(fZ, aa, pid, pidmin, pidmax));
653  }
654  fNPIDs++;
655 }
656 
657 
658 
659 
660 
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
int type
char * Form(const char *fmt,...)
virtual void WriteToAsciiFile(std::ofstream &gridfile)
Definition: KVIDGraph.cpp:442
const KVList * GetInfos() const
Definition: KVIDGraph.h:317
virtual void ReadFromAsciiFile(std::ifstream &gridfile)
Definition: KVIDGraph.cpp:601
const KVNameValueList * GetParameters() const
Definition: KVIDGraph.h:287
const Char_t * GetName() const
Definition: KVIDGraph.cpp:1332
const KVList * GetIdentifiers() const
Definition: KVIDGraph.h:297
void Initialize()
Definition: KVIDGrid.cpp:218
Hybrid charge & mass identification grid.
Bool_t fHasMassIDRegion
set to true if grid has a limited region for mass identification, indicated by an info "MassID"
bool MassIdentificationFromMultiGaussFit(KVMultiGaussIsotopeFit *, KVIdentificationResult *) const
KVList fTables
intervals for mass id
interval_set * GetIntervalSet(int zint) const
virtual void WriteToAsciiFile(std::ofstream &gridfile)
void SetOnlyZId(Bool_t=kTRUE)
void Copy(TObject &obj) const
virtual double DeduceAfromPID(KVIdentificationResult *idr) const
KVUniqueNameList fFits
multi-gaussian fits for mass id
KVMultiGaussIsotopeFit * GetMultiGaussFit(int z) const
virtual void ReadFromAsciiFile(std::ifstream &gridfile)
KVList * GetIntervalSets()
int is_inside(double pid) const
virtual void Identify(Double_t x, Double_t y, KVIdentificationResult *) const
virtual void Copy(TObject &) const
Copy this to 'obj'.
Definition: KVIDZAGrid.cpp:66
virtual void Identify(Double_t x, Double_t y, KVIdentificationResult *) const
UShort_t fZMax
largest Z of lines in grid
Definition: KVIDZAGrid.h:69
Int_t fICode
code de retour
Definition: KVIDZAGrid.h:84
Base class for identification ridge lines corresponding to different nuclear species.
Definition: KVIDZALine.h:32
Base class for graphical cuts used in particle identification.
Definition: KVIDentifier.h:27
Full result of one attempted particle identification.
Bool_t IDOK
general quality of identification, =kTRUE if acceptable identification made
void SetComment(const Char_t *c)
Bool_t Aident
= kTRUE if A of particle established
Double_t PID
= "real" Z if Zident==kTRUE and Aident==kFALSE, "real" A if Zident==Aident==kTRUE
Int_t A
A of particle found (if Aident==kTRUE)
Int_t Z
Z of particle found (if Zident==kTRUE)
Int_t IDquality
specific quality code returned by identification procedure
Bool_t HasFlag(std::string grid_name, TString flag)
Bool_t Zident
=kTRUE if Z of particle established
Function for fitting PID mass spectra.
double GetInterpolatedA(double PID) const
int GetA(double PID, double &P) const
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
void SetValue(const Char_t *name, value_type value)
void RemoveParameter(const Char_t *name)
Bool_t HasStringParameter(const Char_t *name) const
const Char_t * GetStringValue(const Char_t *name) const
bool Set(const KVString &)
TString GetTStringValue(const Char_t *name) const
Strings used to represent a set of ranges of values.
Definition: KVNumberList.h:83
const Char_t * AsString(Int_t maxchars=0) const
void Add(Int_t)
Add value 'n' to the list.
virtual void Copy(TObject &obj) const
virtual void AddLast(TObject *obj)
virtual void SetOwner(Bool_t enable=kTRUE)
virtual void Clear(Option_t *option="")
virtual TObject * At(Int_t idx) const
virtual void Add(TObject *obj)
virtual TObject * FindObject(const char *name) const
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
void Begin(TString delim) const
Definition: KVString.cpp:565
Bool_t End() const
Definition: KVString.cpp:634
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
virtual void Add(TObject *obj)
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Double_t * GetX() const
virtual Double_t Eval(Double_t x, TSpline *spline=nullptr, Option_t *option="") const
virtual const char * GetName() const
virtual void SetName(const char *name)
virtual Bool_t InheritsFrom(const char *classname) const
virtual void Error(const char *method, const char *msgfmt,...) const
Ssiz_t Length() const
Int_t Atoi() const
Double_t Atof() const
const char * Data() const
Bool_t IsWhitespace() const
TString & Remove(EStripType s, char c)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TString & ReplaceAll(const char *s1, const char *s2)
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
TString GetListOfMasses()
bool is_inside(double pid)
bool is_above(double pid)
void add(int aa, double pid, double pidmin=-1., double pidmax=-1.)
double eval(KVIdentificationResult *idr)
interval_set(int zz, int type)
KVList * GetIntervals()
bool is_right_of(double pid)
double GetPID()
bool is_left_of(double pid)
double GetPIDmin()
double GetPIDmax()
bool is_inside(double pid)
Double_t y[n]
Double_t x[n]
Int_t Nint(T x)