CbmRoot
Loading...
Searching...
No Matches
CbmTrdFASP.cxx
Go to the documentation of this file.
1/* Copyright (C) 2018-2020 Horia Hulubei National Institute of Physics and Nuclear Engineering, Bucharest
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Alexandru Bercuci [committer] */
4
5#include "CbmTrdFASP.h"
6
7#include "CbmMatch.h" // for CbmMatch
8#include "CbmTrdDigi.h" // for CbmTrdDigi, CbmTrdDigi::eCbmTrdAsicType::kFASP
9
10#include <Logger.h> // for Logger, LOG
11
12#include <TAxis.h> // for TAxis
13#include <TCanvas.h> // for TCanvas
14#include <TGraph.h> // for TGraph
15#include <TH1.h> // for TH1, TH1F
16#include <TLine.h> // for TLine
17#include <TMath.h> // for Ceil, Floor
18#include <TMathBase.h> // for Min
19#include <TString.h> // for Form, TString
20#include <TVirtualPad.h> // for TVirtualPad
21
22#include <cstdio> // for printf, size_t
23#include <cstring> // for memcpy, memset, strcmp
24#include <iostream> // for operator<<, basic_ostream, cout, ostream
25#include <utility> // for pair, make_pair
26
27#define VERBOSE 0
28#define DRAW 0
29
30using namespace std;
31
32// Bool_t CbmTrdFASP::fgNeighbour = kTRUE; // by default enable neighbour trigger
33Bool_t CbmTrdFASP::fgNeighbour = kFALSE; // AB : match HW on mCBM22
34// Float_t CbmTrdFASP::fgShaperThr = 0.2; // [V]
35// Float_t CbmTrdFASP::fgNeighbourThr = 0.07; // [V]
36//Float_t CbmTrdFASP::fgShaperThr = 0.05; // [V]
37Float_t CbmTrdFASP::fgShaperThr = 0.08; //[V] AB : match HW on mCBM22
39Float_t CbmTrdFASP::fgBaseline = 0.25; //[V] AB : match HW on mCBM22
40Float_t CbmTrdFASP::fgOutGain = 2.025; // [V/ 4095 ADC chs]
41const Int_t CbmTrdFASP::fgkNclkFT = 14; // [clk]
42Int_t CbmTrdFASP::fgNclkLG = 31; // [clk]
43const Int_t CbmTrdFASP::fgkBufferKeep = 400; // [5*ns]
44//___________________________________________________________________
46 : TObject()
47{
50
51 if (uslice < 2 * fgkBufferKeep) {
52 LOG(warn) << "CbmTrdFASP::CbmTrdFASP() : Buffer should be at least " << 2 * fgkBufferKeep
53 << " [5*ns] long. Expand.";
54 uslice = 2 * fgkBufferKeep;
55 }
56 fNphys[0] = 0;
57 fNphys[1] = 0;
58 fHitThPrev.assign(uslice, 0);
59 fShaper.assign(uslice, 0.);
60 fShaperNext.assign(uslice, 0.);
61 if (DRAW) fOut.assign(uslice + 2, 0.);
62 memset(fGraph, 0, NGRAPH * sizeof(TGraph*));
63 memset(fGraphShp, 0, NGRAPH * sizeof(TGraph*));
64 memset(fGraphPhys, 0, NGRAPH * sizeof(TGraph*));
66}
67
68//___________________________________________________________________
70{
71 Int_t nalloc(0);
72 for (Int_t ig(0); ig < NGRAPH; ig++) {
73 if (fGraph[ig]) delete fGraph[ig];
74 if (fGraphShp[ig]) delete fGraphShp[ig];
75 if (fGraphPhys[ig]) delete fGraphPhys[ig];
76 nalloc++;
77 }
78 if (VERBOSE && DRAW) printf("CbmTrdFASP::~CbmTrdFASP() : allocated graphs.");
79 if (fMonitor) delete fMonitor;
80}
81
82//___________________________________________________________________
83void CbmTrdFASP::Clear(Option_t* opt)
84{
85 // printf("CbmTrdFASP::Clear : Nphys[0] = %d\n", fNphys[0]);
86
87 if (fNphys[1]) ProcessShaper('R'); // process last rect channel without interference from current tilt
88 fHitThPrev.assign(fHitThPrev.size(), 0); //*sizeof(Bool_t));
89 WriteDigi(); // finalize fDigi list
90 if (strcmp(opt, "draw") == 0) Draw();
91}
92
93//___________________________________________________________________
94void CbmTrdFASP::InitChannel(int id, const CbmTrdParFaspChannel* par, int asicId, int chId)
95{
96 if (id < 0 || id > 1) {
97 LOG(error) << "FASP simulator supports only two channels at time. Request for channel " << id << " was skipped.";
98 return;
99 }
100 fPar[id] = par;
101 fAsicId[id] = asicId;
102 if (asicId < 0) {
103 fChId[id] = -1;
104 return;
105 }
106
107 if (chId < 0 || chId >= NFASPCH) {
108 LOG(warn) << "FASP ASIC implements " << NFASPCH << " channels. Channel identifier " << chId << " not considered.";
109 fChId[id] = -1;
110 return;
111 }
112 else
113 fChId[id] = chId;
114
115 if (DRAW) {
116 if (fGraphMap.find(fAsicId[id]) == fGraphMap.end()) fGraphMap[fAsicId[id]] = {0};
117 else {
118 int jg = fGraphMap[fAsicId[id]][fChId[id]];
119 if (jg > 0) {
120 memset(fGraph[jg]->GetY(), 0, fGraph[jg]->GetN() * sizeof(double));
121 memset(fGraphShp[jg]->GetY(), 0, fGraphShp[jg]->GetN() * sizeof(double));
122 memset(fGraphPhys[jg]->GetY(), 0, fGraphPhys[jg]->GetN() * sizeof(double));
123 }
124 }
125 }
126}
127
128//___________________________________________________________________
130{
131 int gid(0);
132 if (fNgraph >= NGRAPH) {
133 LOG(warn) << "CbmTrdFASP::AddGraph : Draw buffer size " << NGRAPH << " exhausted. Expert setting.";
134 return gid;
135 }
136
137 int id(typ == 'T' ? 0 : 1);
138 gid = fNgraph;
139 fGraph[gid] = new TGraph(/*fOut.size()*/);
140 fGraph[gid]->SetNameTitle(Form("g%03d", gid),
141 Form("FASP id/ch [%3d / %2d] pad[%3d][%c]", fAsicId[id], fChId[id], fPad, typ));
142 fGraph[gid]->SetFillStyle(0);
143 fGraph[gid]->SetMarkerStyle(7);
144 fGraph[gid]->SetLineColor(kRed);
145 fGraph[gid]->SetMarkerColor(kRed);
146
147 fGraphShp[gid] = (TGraph*) fGraph[gid]->Clone();
148 fGraphShp[gid]->SetName(Form("gs%03d", gid));
149 fGraphShp[gid]->SetLineColor(kBlack);
150 fGraphShp[gid]->SetMarkerColor(kBlack);
151
152 fGraphPhys[gid] = (TGraph*) fGraph[gid]->Clone();
153 fGraphPhys[gid]->SetName(Form("gp%03d", gid));
154 fGraphPhys[gid]->SetMarkerStyle(20);
155 fGraphPhys[gid]->SetLineWidth(2);
156 for (UInt_t ip(0), tm(fStartTime / 5); ip < fOut.size(); ip++, tm += 5)
157 fGraphPhys[gid]->SetPoint(ip, tm, 0.);
158
159 fNgraph++;
160
161 return gid;
162}
163
164//___________________________________________________________________
165void CbmTrdFASP::Draw(Option_t* /*opt*/)
166{
167 if (!DRAW) return;
168 if (!fGraphMap.size()) return;
169
170 TH1* h(nullptr);
171 TVirtualPad* c1(nullptr);
172 if (!fMonitor) {
173 fMonitor = new TCanvas("c", "FASP Simulator", 10, 10, 1850, 2500);
174 fMonitor->Divide(4, 4, 1.e-5, 1.e-5);
175 for (Int_t ic(1); ic <= NFASPCH; ic++) {
176 c1 = fMonitor->cd(ic);
177 c1->SetLogy();
178 c1->SetLeftMargin(0.03580097);
179 c1->SetRightMargin(0.006067961);
180 c1->SetTopMargin(0.01476793);
181 }
182 fGthr = new TLine();
183 fGthr->SetLineStyle(2);
184 }
185
186 for (auto iasic : fGraphMap) {
187 for (Int_t ic(0); ic < NFASPCH; ic++)
188 fMonitor->cd(ic + 1)->Clear();
189 for (int ich(0), jg(0); ich < NFASPCH; ich++) {
190 if (!(jg = iasic.second[ich])) continue;
191
192 int rch = ich / 4, cch = ich % 4;
193 c1 = fMonitor->cd(4 * cch + rch + 1);
194
195 if (!fGraph[jg]) LOG(error) << "Missing graph for " << jg;
196 else
197 fGraph[jg]->Draw("al");
198 if (!fGraphShp[jg]) LOG(error) << "Missing shaper graph for " << jg;
199 else
200 fGraphShp[jg]->Draw("l");
201 if (!fGraphPhys[jg]) LOG(error) << "Missing physics graph for " << jg;
202 else
203 fGraphPhys[jg]->Draw("lp");
204 fGthr->SetLineColor(kGreen);
205 fGthr->DrawLine(0, fgShaperThr + 0.2, 5 * fOut.size(), fgShaperThr + 0.2);
206 fGthr->SetLineColor(kRed);
207 fGthr->DrawLine(0, fgNeighbourThr + 0.2, 5 * fOut.size(), fgNeighbourThr + 0.2);
208 h = fGraph[jg]->GetHistogram();
209 h->GetYaxis()->SetRangeUser(0.1, 6);
210 }
211 fMonitor->Modified();
212 fMonitor->Update();
213 fMonitor->SaveAs(Form("FASP%03d_sim.png", iasic.first));
214 }
215}
216
217//___________________________________________________________________
218Bool_t CbmTrdFASP::Go(ULong64_t time)
219{
223 if (fStartTime == 0) {
224 fStartTime = time;
225 return kFALSE;
226 }
227
228 if (time < fStartTime) {
229 LOG(warn) << "FASP simulator start time " << fStartTime << "[ns] is larger than event time " << time
230 << "[ns]. Skip";
231 return kFALSE;
232 }
233
234 // check to see if there are at most 2us left
235 if (time - fStartTime < fProcTime) return kFALSE;
236
237 return kTRUE; // go
238}
239
240//___________________________________________________________________
241void CbmTrdFASP::GetShaperSignal(Double_t charge)
242{
246
247 // printf("CbmTrdFASP::GetShaperSignal(%5.1f)\n", charge);
248 Int_t idx0(-1), idx1(0);
249 for (Int_t is0(1); is0 < fgkNDB; is0++) {
250 if (charge > fgkCharge[is0]) continue;
251 idx0 = is0 - 1;
252 break;
253 }
254
255 if (idx0 < 0) { // above calibrated region
256 if (VERBOSE) printf(" refMax[%2d]=%5.1f", fgkNDB - 1, fgkCharge[fgkNDB - 1]);
257 memcpy(fSignal, fgkShaper[fgkNDB - 1], FASP_WINDOW * sizeof(Float_t));
258 }
259 else if (idx0 == 0 && charge < fgkCharge[0]) { // below calibrated region
260 if (VERBOSE) printf(" refMin[0]={%5.1f}", fgkCharge[0]);
261 for (Int_t it(0); it < FASP_WINDOW; it++)
262 fSignal[it] = charge * fgkShaper[0][it] / fgkCharge[0];
263 }
264 else {
265 idx1 = idx0 + 1;
266 if (VERBOSE) printf(" ref={%5.1f[%2d] %5.1f[%2d]}", fgkCharge[idx0], idx0, fgkCharge[idx1], idx1);
267 // linear interpolation
268 Double_t dq(fgkCharge[idx0] - fgkCharge[idx1]);
269 for (Int_t it(0); it < FASP_WINDOW; it++)
270 fSignal[it] = (charge * (fgkShaper[idx0][it] - fgkShaper[idx1][it]) + fgkShaper[idx1][it] * fgkCharge[idx0]
271 - fgkShaper[idx0][it] * fgkCharge[idx1])
272 / dq;
273 }
274}
275
276//___________________________________________________________________
278{
279 if (time <= fTimeFT) return fFT;
280
281 Double_t out(0);
282 Int_t i = TMath::Min(Int_t(fShaper.size()), time - 1), j(0);
283 while (i >= 0 && j < SHAPER_LUT - 1) {
284 j = (time - i) - 1;
285 out += 1. * fShaper[i] * fgkShaperLUT[j];
286 i--;
287 }
288 if (fTimeFT < 0 && fTimeDY < 0) return out;
289
290 //decay
291 j = time - (fTimeDY - SHAPER_LUT);
292 if (j > SHAPER_LUT - 1) return out;
293 Double_t decay = fFT * fgkDecayLUT[j];
294 if (out > decay) { // reset
295 fTimeDY = -1;
296 fFT = 0;
297 decay = out;
298 }
299 return decay;
300}
301
302//___________________________________________________________________
303void CbmTrdFASP::PhysToRaw(std::vector<std::pair<CbmTrdDigi*, CbmMatch*>>* vdigi)
304{
308 if (fgNeighbour) ScanDigiNE(vdigi);
309 else
310 ScanDigi(vdigi);
311}
312
313//___________________________________________________________________
314void CbmTrdFASP::Print(Option_t* opt) const
315{
318
319 printf("FASP Simulator : Pad[%3d] NeighbourTrigger[%c]\n", fPad, fgNeighbour ? 'y' : 'n');
320 printf(" Main CH : Trigger[V]=%4.2f Flat-Top[ns]=%5.1f\n", fgShaperThr,
322 printf(" Neighbour CH : Trigger[V]=%4.2f Linear-gate[ns]=%5.1f\n", fgNeighbourThr,
324
325 if (strcmp(opt, "all") != 0) return;
326 Int_t time = 0;
327 vector<Float_t>::const_iterator it = fShaper.cbegin(), jt = fShaperNext.cbegin();
328 while (it != fShaper.cend()) {
329 printf("time[ns]=%4d sgn[mV] : sc(%6.1f) sn(%6.1f)\n", time, 1.e3 * (*it), 1.e3 * (*jt));
330
331 it++;
332 jt++;
333 time += 5;
334 }
335}
336
337//___________________________________________________________________
339{
349
350 int id(typ == 'T' ? 0 : 1), chId = fChId[id], asicId = fAsicId[id];
351 if (VERBOSE)
352 printf(" CbmTrdFASP::ProcessShaper(%c) : pad [%3d] asic/ch[%4d/%2d] "
353 "digis[%d] ...\n",
354 typ, fPad, asicId, chId, fNphys[id]);
355 fTimeFT = -1; // length of the FT gate in 5ns bins (see fgkNclkFT)
356 fTimeLG = -1; // length of the LG gate in 5ns bins (see fgkNclkLG)
357 fFT = 0;
358 // digital signals
359 Bool_t ht(0), htf(0), // hit_threshold level/front for current FASP channel
360 ht_next(0), htf_next(0), // hit_threshold level/front for next FASP channel
361 htf_prev(0), // hit_threshold front for previous FASP channel
362 pk_cmd(0), // peak_command level
363 lg_cmd(0), // linear-gate_command level
364 trigger(0); // trigger type [1] = self [0]=neighbour
365 uint n(0), // no of raw digi found in current shaper
366 ht_time(0), // hit time on current channel
367 cs_time(0); // CS time on current channel
368 double out, max(-1), old(-1);
369 for (size_t i = 1; i < fShaper.size() - 1; i++) {
370 // compute hit threshold level/front for current and neighbour channels
371 htf_prev = 0;
372 htf = 0;
373 htf_next = 0;
374 //TODO NE if (/*fShaper[i-1]<fgNeighbourThr && */ fShaper[i] >= fgNeighbourThr) ht_time = i * 5;
375 if (fShaper[i - 1] < fgShaperThr && fShaper[i] >= fgShaperThr) {
376 ht = 1;
377 htf = 1;
378 trigger = 1;
379 ht_time = i * 5;
380 if (VERBOSE) printf(" %4lu : HT 1\n", i * 5);
381 }
382 else if (fShaper[i - 1] >= fgShaperThr && fShaper[i] < fgShaperThr) {
383 ht = 0;
384 if (VERBOSE) printf(" %4lu : HT 0\n", i * 5);
385 }
386 if (fgNeighbour) { // compute neighbour hit threshold if NE selected
387 if (fShaperNext[i - 1] < fgShaperThr && fShaperNext[i] >= fgShaperThr) {
388 ht_next = 1;
389 htf_next = 1;
390 if (VERBOSE) printf(" %4lu : HT_NEXT 1\n", i * 5);
391 }
392 else if (fShaperNext[i - 1] >= fgShaperThr && fShaperNext[i] < fgShaperThr) {
393 ht_next = 0;
394 if (VERBOSE) printf(" %4lu : HT_NEXT 0\n", i * 5);
395 }
396 if (!fHitThPrev[i] && fHitThPrev[i + 1]) {
397 htf_prev = 1;
398 if (VERBOSE) printf(" %4lu : HT_PREV 1\n", i * 5);
399 }
400 if (VERBOSE && fHitThPrev[i] && !fHitThPrev[i + 1]) printf(" %4lu : HT_PREV 0\n", i * 5);
401 }
402
403 // compute linear gate
404 if (fTimeLG < 0) { // check if linear-gate is closed
405 if ((fgNeighbour && (htf_prev || htf || htf_next)) || (!fgNeighbour && htf)) {
406 lg_cmd = 1;
407 fTimeLG =
408 i
409 + 0.2 * fgNclkLG
410 * CbmTrdDigi::Clk(CbmTrdDigi::eCbmTrdAsicType::kFASP); // compute the end time of LG in FASP-sim 5ns clks
411 if (VERBOSE) printf(" %4lu : LG_CMD 1 -> LG_End[ns]=%d\n", i * 5, fTimeLG * 5);
412 }
413 }
414 else if (static_cast<size_t>(fTimeLG) < i && !ht && !pk_cmd && (fgNeighbour && (!fHitThPrev[i] && !ht_next))) {
415 lg_cmd = 0;
416 trigger = 0;
417 fTimeLG = -1;
418 if (VERBOSE) printf(" %4lu : LG_CMD 0\n", i * 5);
419 }
420
421 if (!pk_cmd && lg_cmd && fShaper[i] >= fgNeighbourThr && fShaper[i] >= fShaper[i - 1]
422 && fShaper[i] > fShaper[i + 1]) {
423 pk_cmd = 1;
424 if (VERBOSE) printf(" %4lu : PK_CMD 1\n", i * 5);
425 }
426
427 out = MakeOut(i);
428
429 if (fTimeFT < 0 && lg_cmd && pk_cmd) {
430 if (max > 0 && out > max) {
431 max = -1;
432 if (VERBOSE) printf(" %4lu : RESET MAX\n", i * 5);
433 }
434 if (max < 0 && out < old) {
435 max = out;
436 if (VERBOSE) printf(" %4lu : MAX[V]=%5.2f\n", i * 5, max);
437 }
438 // start digital processing when the analog signal drops 50 mV under the max value.
439 // Adjust this value to the measured data in order to correctly describe the time difference spectra obtained (e.g. mCBM22)
440 if (out < max - 0.02) {
442 fTimeDY = -1;
443 fFT = max + 0.01;
444 cs_time = uint((i + 3) * 5);
445 // save data for digi update
446 fDigiProc.push_back(make_tuple(
447 ht_time, cs_time, UInt_t(4095 * TMath::Min(Float_t(1.), (fgBaseline + fFT) / fgOutGain)), trigger));
448 n++;
449 if (VERBOSE)
450 printf(" %4lu : HT[ns]=%4u FT[V]=%5.2f CS[ns]=%4u CS_End[ns]=%d "
451 "Trig[%c]\n",
452 i * 5, ht_time, fFT, cs_time, fTimeFT * 5, trigger ? 'S' : 'N');
453 }
454 }
455 if (fTimeFT > 0 && static_cast<size_t>(fTimeFT) <= i) {
457 fTimeFT = -1;
458 fTimeLG = -1;
459 pk_cmd = 0;
460 lg_cmd = 0;
461 trigger = 0;
462 max = -1;
463 if (VERBOSE)
464 printf(" %4lu : LG_CMD 0\n"
465 " %4lu : PK_CMD 0\n",
466 i * 5, i * 5);
467 }
468 if (DRAW) fOut[i + 3] = out;
469 old = out;
470 if (fgNeighbour) fHitThPrev[i] = ht; // save hit threshold for next channel
471 }
472
473 // save results for draw
474 if (DRAW) {
475 int gid = fGraphMap[asicId][chId];
476 if (!gid)
477 LOG(warn) << "CbmTrdFASP::ProcessShaper : Draw representation missing for "
478 << Form(" FASP id/ch[%3d/%2d].", asicId, chId);
479 else {
480 double time = fStartTime;
481 //printf("AB :: g%d[%s] size[%lu]\n", gid, fGraph[gid]->GetName(), fShaper.size());
482 for (uint ip(0); ip < fShaper.size(); ip++, time += 5.) {
483 //printf("time[%3d]=%.0f out[%f] shp[%f]\n", ip, time, 0.3+fOut[ip], 0.2+fShaper[ip]);
484 fGraph[gid]->SetPoint(ip, time, 0.3 + fOut[ip]);
485 fGraphShp[gid]->SetPoint(ip, time, 0.2 + fShaper[ip]);
486 }
487 }
488 }
489 // move to the next channel
490 memcpy(fShaper.data(), fShaperNext.data(), fShaper.size() * sizeof(Float_t));
491 memset(fShaperNext.data(), 0, fShaperNext.size() * sizeof(Float_t));
492
493 return n;
494}
495
496//___________________________________________________________________
497void CbmTrdFASP::ScanDigi(std::vector<std::pair<CbmTrdDigi*, CbmMatch*>>* vdigi)
498{
501
502 fPad = -1;
503 fNphys[0] = 0;
504 fNphys[1] = 0;
505
506 // process digi
507 vector<Float_t>::iterator itb;
508 for (auto iv : (*vdigi)) {
509 CbmTrdDigi* digi = iv.first;
510
511 // first time initialize geographic location
512 if (fPad < 0) {
513 fPad = digi->GetAddressChannel();
514
515 if (VERBOSE > 1)
516 printf(" AB :: (%3d) fStartTime = %llu N=%lu ASIC=[%3d | %2d] ASIC=[%3d | %2d]\n", fPad, fStartTime,
517 vdigi->size(), fAsicId[0], fChId[0], fAsicId[1], fChId[1]);
518 }
519 ULong64_t time = (digi->GetTime() - fStartTime) / 5; // get time in 5ns bins
520
521 if (VERBOSE > 1) printf(" AB :: digi->GetTime() = %f sim->time = %llu\n", digi->GetTime(), time * 5);
522 if (time + FASP_WINDOW > fShaper.size()) {
523 LOG(debug) << "CbmTrdFASP::ScanDigi() : Digi @ pad[" << fPad << " time[" << digi->GetTime()
524 << "] dows not fit in the current buffer starting @ " << fStartTime << "ns. Skip this time.";
525 break;
526 }
527
528 int dt;
529 double t, r = digi->GetCharge(t, dt);
530 t /= 10;
531 r /= 10;
532 if (VERBOSE)
533 printf(" time buffer[5ns]=%4llu / phys[ns]=%lu charge[fC]=%5.1f / %5.1f", time, digi->GetTimeDAQ(), t, r);
534 // tilt pad channel
535 if (t > 0 && fChId[0] >= 0) {
536 if (DRAW) {
537 int gid = fGraphMap[fAsicId[0]][fChId[0]];
538 if (!gid) fGraphMap[fAsicId[0]][fChId[0]] = gid = AddGraph('T');
539 if (gid) fGraphPhys[gid]->SetPoint(time, digi->GetTime(), t / 10.);
540 }
542 itb = fShaper.begin();
543 itb += time;
544 for (Int_t it(0); it < FASP_WINDOW && itb != fShaper.end(); it++, itb++)
545 (*itb) += fSignal[it];
546 fNphys[0]++;
547 }
548
549 // rect pad channel
550 if (r > 0 && fChId[1] >= 0) {
551 if (DRAW) {
552 int gid = fGraphMap[fAsicId[1]][fChId[1]];
553 if (!gid) fGraphMap[fAsicId[1]][fChId[1]] = gid = AddGraph('R');
554 fGraphPhys[gid]->SetPoint(time, digi->GetTime(), r / 10.);
555 }
557 itb = fShaperNext.begin();
558 itb += time + dt;
559 for (Int_t it(0); it < FASP_WINDOW && itb != fShaperNext.end(); it++, itb++)
560 (*itb) += fSignal[it];
561 fNphys[1]++;
562 }
563
564 if (VERBOSE) printf("\n");
565 } // end reading the digis
566
567 if (fNphys[0]) fNraw = ProcessShaper('T'); // process tilt
568 if (fNphys[1]) ProcessShaper('R'); // process rect
569
570 fDigi = vdigi;
571 WriteDigi();
572}
573
574//___________________________________________________________________
575void CbmTrdFASP::ScanDigiNE(std::vector<std::pair<CbmTrdDigi*, CbmMatch*>>* /*vdigi*/)
576{
580
581 LOG(info) << "CbmTrdFASP::ScanDigiNE() : **Method under construction**. Use CbmTrdFASP::ScanDigi() instead by "
582 "selecting CbmTrdFASP:: SetNeighbourTrigger(0)";
583 return;
584
585 // CbmTrdDigi* digi(nullptr);
586 // ULong64_t time;
587 // Int_t dt, gid, asic = row * 9 + col / 8; // TODO should come from calibration
588 // Double_t t, r;
589 // std::vector<std::pair<CbmTrdDigi*, CbmMatch*>>::iterator iv;
590 // vector<Float_t>::iterator itb;
591 //
592 // if (fAsicId < 0) fAsicId = asic; // init asic identifier
593 //
594 // // printf("CbmTrdFASP::ScanDigiNE : fCol[%2d] col[%2d] graph[%d]\n", fCol, col, fGraphId);
595 // // No interference from the previous data
596 // if (fCol < 0 || // first data in the module
597 // (fCol >= 0 && col != fCol + 1) || // column jump
598 // (fRow >= 0 && row != fRow)) { // row jump
599 // Clear((fAsicId != asic ? "draw" : ""));
600 //
601 // if (DRAW) {
602 // gid = fGraphId;
603 // if (fGraphPhys[gid]) memset(fGraphPhys[gid]->GetY(), 0, fGraphPhys[gid]->GetN() * sizeof(Double_t));
604 // else {
605 // fGraphPhys[gid] = new TGraph(fOut.size());
606 // fGraphPhys[gid]->SetName(Form("gp%03d", gid));
607 // fGraphPhys[gid]->SetMarkerStyle(20);
608 // fGraphPhys[gid]->SetLineWidth(2);
609 // for (UInt_t ip(0), tm(fStartTime / 5); ip < fOut.size(); ip++, tm += 5)
610 // fGraphPhys[gid]->SetPoint(ip, tm, 0.);
611 // }
612 // }
613 //
614 // // load data from current tilt channel
615 // iv = vdigi->begin();
616 // while (iv != vdigi->end()) {
617 // digi = iv->first;
618 // // TODO:AB GetTimeDAQ() [clk] and fStartTime [ns]
619 // time = (digi->GetTimeDAQ() - fStartTime) / 5; // get time from buffer start in 5ns bins
620 // if (time + FASP_WINDOW > fShaper.size()) {
621 // LOG(debug) << "CbmTrdFASP::ScanDigiNE() : T-Digi @ row[" << row << "] col[" << col << "] time["
622 // << digi->GetTimeDAQ() << "] does not fit in the current buffer starting @ " << fStartTime
623 // << "ns. Skip this time.";
624 // break;
625 // }
626 // digi->GetCharge(t, dt);
627 // t /= 10.;
628 // if (VERBOSE) printf(" [T] time buffer[5ns]=%4llu / phys[ns]=%lu charge[fC]=%5.1f ", time, digi->GetTimeDAQ(), t);
629 // // tilt pad channel
630 // if (t > 0) {
631 // if (DRAW) fGraphPhys[gid]->SetPoint(time, digi->GetTimeDAQ(), t / 100.);
632 // GetShaperSignal(t);
633 // itb = fShaper.begin();
634 // itb += time;
635 // for (Int_t it(0); it < FASP_WINDOW && itb != fShaper.end(); it++, itb++)
636 // (*itb) += fSignal[it];
637 // fNphys[0]++;
638 // }
639 // else if (VERBOSE)
640 // printf("\n");
641 // iv++;
642 // }
643 // }
644 // else {
645 // if (DRAW) {
646 // gid = fGraphId + 1;
647 // if (fGraphPhys[gid]) memset(fGraphPhys[gid]->GetY(), 0, fGraphPhys[gid]->GetN() * sizeof(Double_t));
648 // else {
649 // fGraphPhys[gid] = new TGraph(fOut.size());
650 // fGraphPhys[gid]->SetName(Form("gp%03d", gid));
651 // fGraphPhys[gid]->SetMarkerStyle(20);
652 // fGraphPhys[gid]->SetLineWidth(2);
653 // for (UInt_t ip(0), tm(fStartTime / 5); ip < fOut.size(); ip++, tm += 5)
654 // fGraphPhys[gid]->SetPoint(ip, tm, 0.);
655 // }
656 // }
657 //
658 // // load tilt digi to account for neighbour trigger
659 // iv = vdigi->begin();
660 // while (iv != vdigi->end()) {
661 // digi = iv->first;
662 // time = (digi->GetTimeDAQ() - fStartTime) / 5; // get time from buffer start in 5ns bins
663 // if (time + FASP_WINDOW > fShaper.size()) {
664 // LOG(debug) << "CbmTrdFASP::ScanDigiNE() : T-Digi @ row[" << row << "] col[" << col << "] time["
665 // << digi->GetTimeDAQ() << "] does not fit in the current buffer starting @ " << fStartTime
666 // << "ns. Skip this time.";
667 // break;
668 // }
669 // digi->GetCharge(t, dt);
670 // t /= 10.;
671 //
672 // if (VERBOSE) printf(" [T] time buffer[5ns]=%4llu / phys[ns]=%lu charge[fC]=%5.1f ", time, digi->GetTimeDAQ(), t);
673 // // tilt pad channel
674 // if (t > 0) {
675 // if (DRAW) fGraphPhys[gid]->SetPoint(time, digi->GetTimeDAQ(), t / 100.);
676 // GetShaperSignal(t);
677 // itb = fShaperNext.begin();
678 // itb += time;
679 // for (Int_t it(0); it < FASP_WINDOW && itb != fShaperNext.end(); it++, itb++)
680 // (*itb) += fSignal[it];
681 // fNphys[1]++;
682 // }
683 // else if (VERBOSE)
684 // printf("\n");
685 // iv++;
686 // }
687 // ProcessShaper('R'); // process previous rect channel
688 // WriteDigi(); // finalize fDigi list
689 // }
690 //
691 // fCol = col;
692 // fRow = row;
693 // fAsicId = asic;
694 //
695 // if (DRAW) {
696 // gid = fGraphId + 1;
697 // if (fGraphPhys[gid]) memset(fGraphPhys[gid]->GetY(), 0, fGraphPhys[gid]->GetN() * sizeof(Double_t));
698 // else {
699 // fGraphPhys[gid] = new TGraph(fOut.size());
700 // fGraphPhys[gid]->SetName(Form("gp%03d", gid));
701 // fGraphPhys[gid]->SetMarkerStyle(20);
702 // fGraphPhys[gid]->SetLineWidth(2);
703 // for (UInt_t ip(0), tm(fStartTime / 5); ip < fOut.size(); ip++, tm += 5)
704 // fGraphPhys[gid]->SetPoint(ip, tm, 0.);
705 // }
706 // }
707 // iv = vdigi->begin();
708 // while (iv != vdigi->end()) {
709 // digi = iv->first;
710 // time = (digi->GetTimeDAQ() - fStartTime) / 5; // get time from buffer start in 5ns bins
711 // if (time + FASP_WINDOW > fShaper.size()) {
712 // LOG(debug) << "CbmTrdFASP::ScanDigiNE() : R-Digi @ row[" << row << "] col[" << col << "] time["
713 // << digi->GetTimeDAQ() << "] dows not fit in the current buffer starting @ " << fStartTime
714 // << "ns. Skip this time.";
715 // break;
716 // }
717 // r = digi->GetCharge(t, dt);
718 // r /= 10.;
719 // if (VERBOSE) printf(" [R] time buffer[5ns]=%4llu / phys[ns]=%lu charge[fC]=%5.1f ", time, digi->GetTimeDAQ(), r);
720 //
721 // // rect pad channel
722 // if (r > 0) {
723 // if (DRAW) fGraphPhys[gid]->SetPoint(time, digi->GetTimeDAQ(), r / 100.);
724 // GetShaperSignal(r);
725 // itb = fShaperNext.begin();
726 // itb += time + dt;
727 // for (Int_t it(0); it < FASP_WINDOW && itb != fShaperNext.end(); it++, itb++)
728 // (*itb) += fSignal[it];
729 // fNphys[1]++;
730 // }
731 // else if (VERBOSE)
732 // printf("\n");
733 // iv++;
734 // }
735 // fNraw = ProcessShaper('T'); // process tilt on current channel
736 // fDigi = vdigi; // save link for further processing
737}
738
739//___________________________________________________________________
740void CbmTrdFASP::SetProcTime(ULong64_t t)
741{
742 if (t == 0) fProcTime = 5 * (fShaper.size() - fgkBufferKeep);
743 else
744 fProcTime = t - fStartTime;
745 //printf(" ProcTime[ns] = %d\n", fProcTime);
746}
747
748//___________________________________________________________________
750{
751
752 if (!fDigi) return;
753 if (VERBOSE) printf(" CbmTrdFASP::WriteDigi(T[%d], R[%d]) ...\n", fNraw, Int_t(fDigiProc.size()) - fNraw);
754
755 auto it = fDigiProc.begin(), // [tilt] iterator
756 ir = it + fNraw, // [rect] iterator
757 im = ir; // middle iterator
758 CbmTrdDigi *digi(nullptr), *digi1(nullptr);
759 CbmMatch* dmatch(nullptr);
760 ULong64_t time, // relative time[ns] of PHYSICAL digi wrt the simulator window start
761 dtime; // relative time[ns] of DIGITAL digi wrt the simulator window start
762 int dt, ddt, trigger(0);
763 double t, r, t0, r0;
764 uint ht_time, // hit time [ns]
765 cs_time, // CS time [ns]
766 tADC, // tilt channel signal [ADC]
767 rADC; // rect channel signal [ADC]
768 //Bool_t mask(0), pileup(0);
769 for (auto dd : (*fDigi)) {
770 digi = dd.first;
771 dmatch = dd.second;
772 time = (digi->GetTime() - fStartTime); // digi time in [ns] from buffer start
773 // stop digi finalize
774 if (time > fProcTime) break;
775 r = digi->GetCharge(t, dt);
776 r /= 10.;
777 t /= 10.;
778 if (VERBOSE) cout << " IN : " << digi->ToString() << " " << dmatch->ToString();
779
780 rADC = 0;
781 tADC = 0;
782 dtime = 0; // relative time of the prompt(T and R) CS signals expressed in ns
783 ddt = 0; // time difference between the T and R signals expressed in clks
784 trigger = 0;
785 ht_time = 0;
786 while (it != im) {
787 ht_time = get<0>(*it);
788 if (VERBOSE) printf(" [T] htime[ns]=%d FT[ADC]=%4d ... \n", ht_time, get<2>(*it));
789 if (ht_time > time) break;
790 it++;
791 }
792 if (t > 0. && it != im) {
793 cs_time = get<1>(*it);
794 //printf("match T htime[%d] time[%llu]\n", ht_time, time);
795 if (ht_time - time < 400) { // found converted hit
796 if (VERBOSE)
797 printf(" [T] ht[ns]=%4d CS[ns]=%4d FT[ADC]=%4u Trg[%s]\n", ht_time, cs_time, get<2>(*it),
798 (get<3>(*it) ? "S" : "N"));
799 dtime = cs_time;
800 tADC = get<2>(*it);
801 if (get<3>(*it)) trigger |= 1;
802 //if(csTime-hTime > 350) pileup = kTRUE; // 350 ns max peak time
803 it++;
804 } //else if(t>40.) mask=kTRUE; // hit not converted : possible under threshold, masked
805 // 40fC threshold
806 }
807 //printf("AB :: time[%llu] dtime[%llu] dt[%d]\n", time, dtime, dt);
808
809 time += dt;
810 ht_time = 0;
811 while (ir != fDigiProc.end()) {
812 ht_time = get<0>(*ir);
813 if (VERBOSE) printf(" [R] htime[ns]=%d FT[ADC]=%4d ...\n", ht_time, get<2>(*ir));
814 if (ht_time > time) break;
815 ir++;
816 }
817 if (r > 0. && ir != fDigiProc.end()) {
818 cs_time = get<1>(*ir);
819 //printf("match R htime[%d] time[%llu]\n", ht_time, time);
820 if (ht_time - time < 400) { // found converted hit
821 if (VERBOSE)
822 printf(" [R] ht[ns]=%4d CS[ns]=%4d FT[ADC]=%4u Trg[%s]\n", ht_time, cs_time, get<2>(*ir),
823 (get<3>(*ir) ? "S" : "N"));
824 if (dtime) {
825 if (cs_time > dtime)
826 ddt = TMath::Ceil(Int_t(cs_time - dtime) / CbmTrdDigi::Clk(CbmTrdDigi::eCbmTrdAsicType::kFASP));
827 else
828 ddt = TMath::Floor(Int_t(cs_time - dtime) / CbmTrdDigi::Clk(CbmTrdDigi::eCbmTrdAsicType::kFASP));
829 }
830 else
831 dtime = cs_time;
832 rADC = get<2>(*ir);
833 if (get<3>(*ir)) trigger |= 2;
834 //if(csTime-hTime > 350) pileup = kTRUE; // 350 ns max peak time
835 ir++;
836 } //else if(r>40.) mask=kTRUE; // hit not converted : possible under threshold, masked
837 // 40fC threshold
838 }
839
840 //printf("AB :: time[%llu] dtime[%llu] ddt[%d]\n", time, dtime, ddt);
841 //update digi
842 digi->SetMasked(1); // mark as processed
843 if (dtime) {
844 digi->SetTime(fStartTime + dtime);
845 digi->SetCharge(tADC, rADC, ddt);
846 digi->SetTriggerType(trigger);
847 if (VERBOSE) cout << " OUT: " << digi->ToString();
848 }
849 else
850 digi->SetFlag(0, kTRUE); // mark for deletion
851 if (VERBOSE)
852 cout << " ===================================="
853 << "\n";
854 }
855
856 // try to merge digits
857 ULong64_t time0, time1;
858 Char_t type0, // prompt digi. kTRUE if rectangle
859 type1; // late digi. kTRUE if rectangle
860 for (auto iv = fDigi->begin(); iv != fDigi->end(); iv++) {
861 digi = iv->first;
862 if (digi->IsFlagged(0)) continue; // no output
863 if (!digi->IsMasked()) break; // not finalized
864 r0 = digi->GetCharge(t0, dt);
865 if (r0 > 0. && t0 > 0.) continue; // already complete
866 type0 = (r0 > 0. ? 1 : 0); // mark type for prompt digi
867 time0 = digi->GetTimeDAQ(); // mark time for prompt digi
868 for (auto jv = iv + 1; jv != fDigi->end(); jv++) {
869 digi1 = jv->first;
870 if (digi1->IsFlagged(0)) continue; // no output
871 r = digi1->GetCharge(t, dt);
872 if (r > 0. && t > 0.) continue; // already complete
873 type1 = (r > 0. ? 1 : 0); // mark type for late digi
874 time1 = digi1->GetTimeDAQ(); // mark time for late digi
875
876 // try merge
877 if (type0 == type1) break; // same type
878 dt = time1 - time0;
879 if (dt > 7) break; // digits too far to be merged
880
881 if (VERBOSE) cout << "MERGE: " << digi->ToString() << " " << digi1->ToString();
882 if (type0) { // prompt digi rectangle
883 digi->SetTime(digi1->GetTime());
884 digi->SetCharge(t, r0, -dt);
885 Int_t rtrg(digi->GetTriggerType() & 2), ttrg(digi1->GetTriggerType() & 1);
886 digi->SetTriggerType(rtrg | ttrg); //merge the triggers
887 }
888 else { // prompt digi tilt
889 digi->SetCharge(t0, r, dt);
890 Int_t ttrg(digi->GetTriggerType() & 1), rtrg(digi1->GetTriggerType() & 2);
891 digi->SetTriggerType(rtrg | ttrg); //merge the triggers
892 }
893 if (VERBOSE) cout << "RES : " << digi->ToString();
894
895 digi1->SetFlag(0, kTRUE); // mark for deletion
896 break;
897 }
898 }
899
900 // reset processed digi vector
901 fDigiProc.clear();
902 fNraw = 0;
903 fDigi = nullptr; //remove old link
904}
905
906const Float_t CbmTrdFASP::fgkShaperPar[4] = {0.068, 26.1, 2.61, 50};
908 7.52e-04, 3.79e-03, 9.02e-03, 1.58e-02, 2.33e-02, 3.10e-02, 3.83e-02, 4.48e-02, 5.03e-02, 5.46e-02,
909 5.78e-02, 5.99e-02, 6.10e-02, 6.11e-02, 6.04e-02, 5.90e-02, 5.71e-02, 5.47e-02, 5.20e-02, 4.91e-02,
910 4.60e-02, 4.29e-02, 3.98e-02, 3.67e-02, 3.37e-02, 3.09e-02, 2.81e-02, 2.55e-02, 2.31e-02, 2.08e-02,
911 1.87e-02, 1.68e-02, 1.50e-02, 1.34e-02, 1.20e-02, 1.06e-02, 9.42e-03, 8.34e-03, 7.37e-03, 6.50e-03,
912 5.72e-03, 5.03e-03, 4.42e-03, 3.87e-03, 3.39e-03, 2.97e-03, 2.59e-03, 2.26e-03, 1.97e-03, 1.71e-03,
913 1.49e-03, 1.29e-03, 1.12e-03, 9.73e-04, 8.43e-04, 7.30e-04, 6.31e-04, 5.45e-04, 4.71e-04, 4.06e-04,
914 3.50e-04, 3.01e-04, 2.60e-04, 2.23e-04, 1.92e-04, 1.65e-04, 1.42e-04, 1.22e-04, 1.04e-04, 8.94e-05,
915 7.66e-05, 6.56e-05, 5.61e-05, 4.80e-05, 4.11e-05, 3.51e-05, 3.00e-05, 2.56e-05, 2.19e-05, 1.86e-05};
917 1.00e+00, 9.05e-01, 8.19e-01, 7.41e-01, 6.70e-01, 6.07e-01, 5.49e-01, 4.97e-01, 4.49e-01, 4.07e-01,
918 3.68e-01, 3.33e-01, 3.01e-01, 2.73e-01, 2.47e-01, 2.23e-01, 2.02e-01, 1.83e-01, 1.65e-01, 1.50e-01,
919 1.35e-01, 1.22e-01, 1.11e-01, 1.00e-01, 9.07e-02, 8.21e-02, 7.43e-02, 6.72e-02, 6.08e-02, 5.50e-02,
920 4.98e-02, 4.50e-02, 4.08e-02, 3.69e-02, 3.34e-02, 3.02e-02, 2.73e-02, 2.47e-02, 2.24e-02, 2.02e-02,
921 1.83e-02, 1.66e-02, 1.50e-02, 1.36e-02, 1.23e-02, 1.11e-02, 1.01e-02, 9.10e-03, 8.23e-03, 7.45e-03,
922 6.74e-03, 6.10e-03, 5.52e-03, 4.99e-03, 4.52e-03, 4.09e-03, 3.70e-03, 3.35e-03, 3.03e-03, 2.74e-03,
923 2.48e-03, 2.24e-03, 2.03e-03, 1.84e-03, 1.66e-03, 1.50e-03, 1.36e-03, 1.23e-03, 1.11e-03, 1.01e-03,
924 9.12e-04, 8.25e-04, 7.47e-04, 6.76e-04, 6.11e-04, 5.53e-04, 5.00e-04, 4.53e-04, 4.10e-04, 3.71e-04};
926 6, 12, 18, 24, 30, 36, 42, 48, 54, 60, 66, 72, 78, 84, 90, 96, 102, 108,
927 114, 120, 126, 132, 138, 144, 150, 156, 162, 168, 180, 195, 210, 225, 240, 255, 270, 285,
928 300, 315, 330, 345, 360, 375, 390, 405, 420, 435, 450, 465, 480, 495, 510, 525, 540};
930 {// charge = 6 fC
931 0.000, 0.000, 0.000, 0.001, 0.002, 0.004, 0.005, 0.007, 0.009, 0.011, 0.012, 0.014, 0.015, 0.017,
932 0.018, 0.020, 0.021, 0.022, 0.023, 0.024, 0.025, 0.025, 0.026, 0.027, 0.027, 0.028, 0.028, 0.029,
933 0.029, 0.029, 0.029, 0.029, 0.029, 0.029, 0.029, 0.029, 0.029, 0.029, 0.029, 0.029, 0.028, 0.028,
934 0.028, 0.027, 0.027, 0.027, 0.026, 0.026, 0.025, 0.025, 0.024, 0.024, 0.024, 0.023, 0.022, 0.022,
935 0.021, 0.021, 0.021, 0.020, 0.020, 0.019, 0.018, 0.018, 0.017, 0.017, 0.016, 0.016, 0.016, 0.015,
936 0.015, 0.014, 0.014, 0.013, 0.013, 0.012, 0.012, 0.011, 0.011, 0.011, 0.010, 0.010, 0.009, 0.009,
937 0.009, 0.008, 0.008, 0.008, 0.007, 0.007, 0.007, 0.006, 0.006, 0.006, 0.006, 0.005, 0.005, 0.005,
938 0.005, 0.004, 0.004, 0.004, 0.004, 0.003, 0.003, 0.003, 0.003, 0.003, 0.003, 0.002, 0.002, 0.002,
939 0.002, 0.002, 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
940 0.001, 0.001, 0.001, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
941 0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000,
942 -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000,
943 -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000,
944 -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000,
945 -0.000, -0.000, -0.000, -0.000},
946 {// charge = 12 fC
947 0.000, 0.000, 0.001, 0.002, 0.004, 0.007, 0.011, 0.014, 0.017, 0.021, 0.024, 0.027, 0.031, 0.034,
948 0.036, 0.039, 0.041, 0.044, 0.046, 0.048, 0.049, 0.051, 0.052, 0.053, 0.055, 0.056, 0.056, 0.057,
949 0.058, 0.058, 0.058, 0.059, 0.059, 0.059, 0.059, 0.058, 0.058, 0.058, 0.058, 0.057, 0.057, 0.056,
950 0.055, 0.055, 0.054, 0.053, 0.052, 0.052, 0.051, 0.050, 0.049, 0.048, 0.047, 0.046, 0.045, 0.044,
951 0.043, 0.042, 0.041, 0.040, 0.039, 0.038, 0.037, 0.036, 0.035, 0.034, 0.033, 0.032, 0.031, 0.030,
952 0.029, 0.028, 0.027, 0.026, 0.025, 0.025, 0.024, 0.023, 0.022, 0.021, 0.020, 0.020, 0.019, 0.018,
953 0.017, 0.017, 0.016, 0.015, 0.015, 0.014, 0.013, 0.013, 0.012, 0.012, 0.011, 0.011, 0.010, 0.010,
954 0.009, 0.009, 0.008, 0.008, 0.007, 0.007, 0.007, 0.006, 0.006, 0.006, 0.005, 0.005, 0.005, 0.004,
955 0.004, 0.004, 0.003, 0.003, 0.003, 0.003, 0.003, 0.002, 0.002, 0.002, 0.002, 0.002, 0.002, 0.001,
956 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000,
957 -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.001, -0.001, -0.001,
958 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
959 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
960 -0.001, -0.001, -0.001, -0.001, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.000,
961 -0.000, -0.000, -0.000, -0.000},
962 {// charge = 18 fC
963 0.000, -0.000, 0.001, 0.003, 0.007, 0.011, 0.016, 0.021, 0.026, 0.031, 0.036, 0.041, 0.046, 0.050,
964 0.054, 0.058, 0.062, 0.065, 0.068, 0.071, 0.074, 0.076, 0.078, 0.080, 0.082, 0.083, 0.085, 0.086,
965 0.086, 0.087, 0.088, 0.088, 0.088, 0.088, 0.088, 0.088, 0.087, 0.087, 0.086, 0.086, 0.085, 0.084,
966 0.083, 0.082, 0.081, 0.080, 0.079, 0.077, 0.076, 0.075, 0.073, 0.072, 0.071, 0.069, 0.068, 0.066,
967 0.065, 0.063, 0.062, 0.060, 0.058, 0.057, 0.055, 0.054, 0.052, 0.051, 0.049, 0.048, 0.046, 0.045,
968 0.044, 0.042, 0.041, 0.039, 0.038, 0.037, 0.035, 0.034, 0.033, 0.032, 0.031, 0.029, 0.028, 0.027,
969 0.026, 0.025, 0.024, 0.023, 0.022, 0.021, 0.020, 0.019, 0.018, 0.018, 0.017, 0.016, 0.015, 0.015,
970 0.014, 0.013, 0.012, 0.012, 0.011, 0.011, 0.010, 0.009, 0.009, 0.008, 0.008, 0.007, 0.007, 0.007,
971 0.006, 0.006, 0.005, 0.005, 0.005, 0.004, 0.004, 0.004, 0.003, 0.003, 0.003, 0.003, 0.002, 0.002,
972 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.000, 0.000, 0.000, 0.000, 0.000,
973 -0.000, -0.000, -0.000, -0.000, -0.000, -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
974 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
975 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
976 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
977 -0.001, -0.001, -0.001, -0.001},
978 {// charge = 24 fC
979 0.000, -0.000, 0.001, 0.004, 0.009, 0.015, 0.021, 0.028, 0.035, 0.042, 0.049, 0.055, 0.061, 0.067,
980 0.073, 0.078, 0.083, 0.087, 0.091, 0.095, 0.099, 0.102, 0.104, 0.107, 0.109, 0.111, 0.113, 0.114,
981 0.115, 0.116, 0.117, 0.117, 0.118, 0.118, 0.117, 0.117, 0.117, 0.116, 0.115, 0.114, 0.113, 0.112,
982 0.111, 0.109, 0.108, 0.106, 0.105, 0.103, 0.101, 0.100, 0.098, 0.096, 0.094, 0.092, 0.090, 0.088,
983 0.086, 0.084, 0.082, 0.080, 0.078, 0.076, 0.074, 0.072, 0.070, 0.068, 0.066, 0.064, 0.062, 0.060,
984 0.058, 0.056, 0.054, 0.053, 0.051, 0.049, 0.047, 0.046, 0.044, 0.042, 0.041, 0.039, 0.038, 0.036,
985 0.035, 0.033, 0.032, 0.031, 0.029, 0.028, 0.027, 0.026, 0.025, 0.023, 0.022, 0.021, 0.020, 0.019,
986 0.018, 0.017, 0.016, 0.016, 0.015, 0.014, 0.013, 0.012, 0.012, 0.011, 0.011, 0.010, 0.009, 0.009,
987 0.008, 0.008, 0.007, 0.007, 0.006, 0.006, 0.005, 0.005, 0.004, 0.004, 0.004, 0.003, 0.003, 0.003,
988 0.002, 0.002, 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.001, 0.001, 0.000, 0.000, 0.000, 0.000,
989 -0.000, -0.000, -0.000, -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
990 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
991 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
992 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
993 -0.001, -0.001, -0.001, -0.001},
994 {// charge = 30 fC
995 0.000, -0.000, 0.001, 0.005, 0.011, 0.018, 0.026, 0.035, 0.044, 0.052, 0.061, 0.069, 0.076, 0.084,
996 0.091, 0.097, 0.103, 0.109, 0.114, 0.119, 0.123, 0.127, 0.131, 0.134, 0.137, 0.139, 0.141, 0.143,
997 0.144, 0.145, 0.146, 0.147, 0.147, 0.147, 0.147, 0.146, 0.146, 0.145, 0.144, 0.143, 0.141, 0.140,
998 0.138, 0.137, 0.135, 0.133, 0.131, 0.129, 0.127, 0.125, 0.122, 0.120, 0.118, 0.115, 0.113, 0.110,
999 0.108, 0.105, 0.103, 0.100, 0.097, 0.095, 0.092, 0.090, 0.087, 0.085, 0.082, 0.080, 0.077, 0.075,
1000 0.073, 0.070, 0.068, 0.066, 0.063, 0.061, 0.059, 0.057, 0.055, 0.053, 0.051, 0.049, 0.047, 0.045,
1001 0.044, 0.042, 0.040, 0.038, 0.037, 0.035, 0.034, 0.032, 0.031, 0.029, 0.028, 0.027, 0.025, 0.024,
1002 0.023, 0.022, 0.021, 0.020, 0.019, 0.018, 0.017, 0.016, 0.015, 0.014, 0.013, 0.012, 0.012, 0.011,
1003 0.010, 0.009, 0.009, 0.008, 0.008, 0.007, 0.007, 0.006, 0.006, 0.005, 0.005, 0.004, 0.004, 0.003,
1004 0.003, 0.003, 0.002, 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.001, 0.001, 0.000, 0.000, 0.000,
1005 -0.000, -0.000, -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
1006 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
1007 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
1008 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
1009 -0.001, -0.001, -0.001, -0.001},
1010 {// charge = 36 fC
1011 0.000, -0.000, 0.002, 0.006, 0.013, 0.022, 0.032, 0.042, 0.052, 0.063, 0.073, 0.082, 0.092, 0.101,
1012 0.109, 0.117, 0.124, 0.131, 0.137, 0.143, 0.148, 0.153, 0.157, 0.160, 0.164, 0.167, 0.169, 0.171,
1013 0.173, 0.174, 0.175, 0.176, 0.176, 0.176, 0.176, 0.176, 0.175, 0.174, 0.173, 0.171, 0.170, 0.168,
1014 0.166, 0.164, 0.162, 0.160, 0.157, 0.155, 0.152, 0.150, 0.147, 0.144, 0.141, 0.138, 0.135, 0.132,
1015 0.129, 0.126, 0.123, 0.120, 0.117, 0.114, 0.111, 0.108, 0.105, 0.102, 0.099, 0.096, 0.093, 0.090,
1016 0.087, 0.084, 0.082, 0.079, 0.076, 0.074, 0.071, 0.068, 0.066, 0.064, 0.061, 0.059, 0.057, 0.054,
1017 0.052, 0.050, 0.048, 0.046, 0.044, 0.042, 0.040, 0.039, 0.037, 0.035, 0.034, 0.032, 0.030, 0.029,
1018 0.027, 0.026, 0.025, 0.024, 0.022, 0.021, 0.020, 0.019, 0.018, 0.017, 0.016, 0.015, 0.014, 0.013,
1019 0.012, 0.011, 0.011, 0.010, 0.009, 0.008, 0.008, 0.007, 0.007, 0.006, 0.006, 0.005, 0.005, 0.004,
1020 0.004, 0.003, 0.003, 0.003, 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.001, 0.000, 0.000, 0.000,
1021 -0.000, -0.000, -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.002,
1022 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1023 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.001,
1024 -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
1025 -0.001, -0.001, -0.001, -0.001},
1026 {// charge = 42 fC
1027 0.000, -0.000, 0.002, 0.007, 0.016, 0.026, 0.037, 0.049, 0.061, 0.073, 0.085, 0.096, 0.107, 0.117,
1028 0.127, 0.136, 0.145, 0.152, 0.160, 0.166, 0.173, 0.178, 0.183, 0.187, 0.191, 0.194, 0.197, 0.200,
1029 0.202, 0.203, 0.205, 0.205, 0.206, 0.206, 0.206, 0.205, 0.204, 0.203, 0.202, 0.200, 0.198, 0.196,
1030 0.194, 0.192, 0.189, 0.186, 0.184, 0.181, 0.178, 0.174, 0.171, 0.168, 0.165, 0.161, 0.158, 0.154,
1031 0.151, 0.147, 0.144, 0.140, 0.137, 0.133, 0.129, 0.126, 0.122, 0.119, 0.115, 0.112, 0.109, 0.105,
1032 0.102, 0.099, 0.095, 0.092, 0.089, 0.086, 0.083, 0.080, 0.077, 0.074, 0.071, 0.069, 0.066, 0.063,
1033 0.061, 0.058, 0.056, 0.054, 0.051, 0.049, 0.047, 0.045, 0.043, 0.041, 0.039, 0.037, 0.035, 0.034,
1034 0.032, 0.030, 0.029, 0.027, 0.026, 0.024, 0.023, 0.022, 0.021, 0.020, 0.018, 0.017, 0.016, 0.015,
1035 0.014, 0.013, 0.012, 0.012, 0.011, 0.010, 0.009, 0.008, 0.008, 0.007, 0.007, 0.006, 0.005, 0.005,
1036 0.004, 0.004, 0.003, 0.003, 0.003, 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.001, 0.000, 0.000,
1037 -0.000, -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002,
1038 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1039 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1040 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001,
1041 -0.001, -0.001, -0.001, -0.001},
1042 {// charge = 48 fC
1043 0.000, -0.000, 0.002, 0.008, 0.018, 0.029, 0.042, 0.056, 0.070, 0.083, 0.097, 0.110, 0.122, 0.134,
1044 0.145, 0.155, 0.165, 0.174, 0.183, 0.190, 0.197, 0.203, 0.209, 0.214, 0.218, 0.222, 0.226, 0.228,
1045 0.231, 0.232, 0.234, 0.235, 0.235, 0.235, 0.235, 0.234, 0.233, 0.232, 0.230, 0.229, 0.227, 0.224,
1046 0.222, 0.219, 0.216, 0.213, 0.210, 0.206, 0.203, 0.200, 0.196, 0.192, 0.188, 0.184, 0.180, 0.176,
1047 0.172, 0.168, 0.164, 0.160, 0.156, 0.152, 0.148, 0.144, 0.140, 0.136, 0.132, 0.128, 0.124, 0.120,
1048 0.116, 0.113, 0.109, 0.105, 0.102, 0.098, 0.095, 0.091, 0.088, 0.085, 0.082, 0.078, 0.076, 0.072,
1049 0.070, 0.067, 0.064, 0.061, 0.059, 0.056, 0.054, 0.051, 0.049, 0.047, 0.045, 0.043, 0.041, 0.039,
1050 0.037, 0.035, 0.033, 0.031, 0.030, 0.028, 0.026, 0.025, 0.024, 0.022, 0.021, 0.020, 0.018, 0.017,
1051 0.016, 0.015, 0.014, 0.013, 0.012, 0.011, 0.011, 0.010, 0.009, 0.008, 0.007, 0.007, 0.006, 0.006,
1052 0.005, 0.004, 0.004, 0.003, 0.003, 0.003, 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.000, 0.000,
1053 -0.000, -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1054 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1055 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1056 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.001,
1057 -0.001, -0.001, -0.001, -0.001},
1058 {// charge = 54 fC
1059 0.000, -0.000, 0.002, 0.009, 0.020, 0.033, 0.047, 0.063, 0.078, 0.094, 0.109, 0.123, 0.137, 0.151,
1060 0.163, 0.175, 0.186, 0.196, 0.205, 0.214, 0.222, 0.229, 0.235, 0.241, 0.246, 0.250, 0.254, 0.257,
1061 0.259, 0.261, 0.263, 0.264, 0.264, 0.264, 0.264, 0.263, 0.262, 0.261, 0.259, 0.257, 0.255, 0.252,
1062 0.249, 0.247, 0.243, 0.240, 0.236, 0.232, 0.229, 0.225, 0.220, 0.216, 0.212, 0.207, 0.203, 0.198,
1063 0.194, 0.189, 0.185, 0.180, 0.176, 0.171, 0.166, 0.162, 0.157, 0.153, 0.148, 0.144, 0.140, 0.135,
1064 0.131, 0.127, 0.123, 0.118, 0.114, 0.110, 0.107, 0.103, 0.099, 0.095, 0.092, 0.088, 0.085, 0.082,
1065 0.078, 0.075, 0.072, 0.069, 0.066, 0.063, 0.060, 0.058, 0.055, 0.053, 0.050, 0.048, 0.046, 0.043,
1066 0.041, 0.039, 0.037, 0.035, 0.034, 0.031, 0.030, 0.028, 0.027, 0.025, 0.024, 0.022, 0.021, 0.020,
1067 0.018, 0.017, 0.016, 0.015, 0.014, 0.013, 0.012, 0.011, 0.010, 0.009, 0.008, 0.008, 0.007, 0.006,
1068 0.006, 0.005, 0.004, 0.004, 0.003, 0.003, 0.002, 0.002, 0.002, 0.001, 0.001, 0.001, 0.000, 0.000,
1069 -0.000, -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1070 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1071 -0.003, -0.003, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1072 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1073 -0.002, -0.002, -0.002, -0.002},
1074 {// charge = 60 fC
1075 0.000, -0.000, 0.002, 0.010, 0.022, 0.036, 0.053, 0.070, 0.087, 0.104, 0.121, 0.137, 0.152, 0.167,
1076 0.181, 0.194, 0.206, 0.218, 0.228, 0.238, 0.246, 0.254, 0.261, 0.267, 0.273, 0.278, 0.282, 0.285,
1077 0.288, 0.290, 0.292, 0.293, 0.294, 0.294, 0.293, 0.293, 0.292, 0.290, 0.288, 0.286, 0.283, 0.280,
1078 0.277, 0.274, 0.270, 0.266, 0.262, 0.258, 0.254, 0.249, 0.245, 0.240, 0.235, 0.230, 0.225, 0.220,
1079 0.215, 0.210, 0.205, 0.200, 0.195, 0.190, 0.185, 0.180, 0.175, 0.170, 0.165, 0.160, 0.155, 0.150,
1080 0.146, 0.141, 0.136, 0.132, 0.127, 0.123, 0.118, 0.114, 0.110, 0.106, 0.102, 0.098, 0.094, 0.091,
1081 0.087, 0.083, 0.080, 0.077, 0.073, 0.070, 0.067, 0.064, 0.061, 0.059, 0.056, 0.053, 0.051, 0.048,
1082 0.046, 0.044, 0.041, 0.039, 0.037, 0.035, 0.033, 0.031, 0.030, 0.028, 0.026, 0.025, 0.023, 0.022,
1083 0.020, 0.019, 0.018, 0.016, 0.015, 0.014, 0.013, 0.012, 0.011, 0.010, 0.009, 0.008, 0.008, 0.007,
1084 0.006, 0.006, 0.005, 0.004, 0.004, 0.003, 0.003, 0.002, 0.002, 0.002, 0.001, 0.001, 0.000, 0.000,
1085 -0.000, -0.001, -0.001, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1086 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1087 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.002, -0.002, -0.002, -0.002,
1088 -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1089 -0.002, -0.002, -0.002, -0.002},
1090 {// charge = 66 fC
1091 0.000, -0.000, 0.003, 0.011, 0.024, 0.040, 0.058, 0.076, 0.095, 0.114, 0.133, 0.151, 0.168, 0.184,
1092 0.199, 0.214, 0.227, 0.239, 0.251, 0.261, 0.271, 0.279, 0.287, 0.294, 0.300, 0.305, 0.310, 0.314,
1093 0.317, 0.319, 0.321, 0.322, 0.323, 0.323, 0.323, 0.322, 0.321, 0.319, 0.317, 0.314, 0.312, 0.308,
1094 0.305, 0.301, 0.297, 0.293, 0.289, 0.284, 0.279, 0.274, 0.269, 0.264, 0.259, 0.253, 0.248, 0.243,
1095 0.237, 0.231, 0.226, 0.220, 0.215, 0.209, 0.203, 0.198, 0.192, 0.187, 0.181, 0.176, 0.171, 0.165,
1096 0.160, 0.155, 0.150, 0.145, 0.140, 0.135, 0.130, 0.126, 0.121, 0.117, 0.112, 0.108, 0.104, 0.100,
1097 0.096, 0.092, 0.088, 0.084, 0.081, 0.077, 0.074, 0.071, 0.067, 0.064, 0.062, 0.059, 0.056, 0.053,
1098 0.050, 0.048, 0.045, 0.043, 0.041, 0.039, 0.036, 0.035, 0.032, 0.031, 0.029, 0.027, 0.025, 0.024,
1099 0.022, 0.021, 0.020, 0.018, 0.017, 0.016, 0.014, 0.013, 0.012, 0.011, 0.010, 0.009, 0.008, 0.008,
1100 0.007, 0.006, 0.006, 0.005, 0.004, 0.004, 0.003, 0.003, 0.002, 0.002, 0.001, 0.001, 0.000, 0.000,
1101 -0.000, -0.001, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.003, -0.003,
1102 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1103 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1104 -0.003, -0.003, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1105 -0.002, -0.002, -0.002, -0.002},
1106 {// charge = 72 fC
1107 0.000, -0.000, 0.003, 0.012, 0.026, 0.044, 0.063, 0.083, 0.104, 0.124, 0.145, 0.164, 0.183, 0.201,
1108 0.217, 0.233, 0.247, 0.261, 0.273, 0.285, 0.295, 0.305, 0.313, 0.321, 0.327, 0.333, 0.338, 0.342,
1109 0.346, 0.348, 0.350, 0.352, 0.353, 0.353, 0.352, 0.351, 0.350, 0.348, 0.346, 0.343, 0.340, 0.336,
1110 0.333, 0.329, 0.324, 0.320, 0.315, 0.310, 0.305, 0.299, 0.294, 0.288, 0.282, 0.277, 0.271, 0.265,
1111 0.259, 0.253, 0.247, 0.240, 0.234, 0.228, 0.222, 0.216, 0.210, 0.204, 0.198, 0.192, 0.186, 0.180,
1112 0.175, 0.169, 0.163, 0.158, 0.153, 0.147, 0.142, 0.137, 0.132, 0.127, 0.123, 0.118, 0.113, 0.109,
1113 0.104, 0.100, 0.096, 0.092, 0.088, 0.084, 0.081, 0.077, 0.074, 0.070, 0.067, 0.064, 0.061, 0.058,
1114 0.055, 0.053, 0.049, 0.047, 0.044, 0.042, 0.040, 0.038, 0.035, 0.033, 0.031, 0.030, 0.028, 0.026,
1115 0.024, 0.023, 0.021, 0.020, 0.018, 0.017, 0.016, 0.015, 0.013, 0.012, 0.011, 0.010, 0.009, 0.008,
1116 0.008, 0.007, 0.006, 0.005, 0.005, 0.004, 0.003, 0.003, 0.002, 0.002, 0.001, 0.001, 0.001, 0.000,
1117 -0.000, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003,
1118 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1119 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1120 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002, -0.002,
1121 -0.002, -0.002, -0.002, -0.002},
1122 {// charge = 78 fC
1123 0.000, -0.000, 0.003, 0.013, 0.029, 0.047, 0.068, 0.090, 0.112, 0.135, 0.156, 0.178, 0.198, 0.217,
1124 0.235, 0.252, 0.268, 0.283, 0.296, 0.309, 0.320, 0.330, 0.339, 0.347, 0.355, 0.361, 0.366, 0.371,
1125 0.375, 0.378, 0.380, 0.381, 0.382, 0.382, 0.382, 0.381, 0.379, 0.377, 0.375, 0.372, 0.368, 0.365,
1126 0.361, 0.356, 0.351, 0.346, 0.341, 0.336, 0.330, 0.324, 0.318, 0.312, 0.306, 0.300, 0.293, 0.287,
1127 0.280, 0.274, 0.267, 0.260, 0.254, 0.247, 0.240, 0.234, 0.227, 0.221, 0.214, 0.208, 0.202, 0.195,
1128 0.189, 0.183, 0.177, 0.171, 0.165, 0.160, 0.154, 0.149, 0.143, 0.138, 0.133, 0.128, 0.123, 0.118,
1129 0.113, 0.109, 0.104, 0.100, 0.096, 0.091, 0.087, 0.084, 0.080, 0.076, 0.073, 0.069, 0.066, 0.063,
1130 0.060, 0.057, 0.053, 0.051, 0.048, 0.046, 0.043, 0.041, 0.038, 0.036, 0.034, 0.032, 0.030, 0.028,
1131 0.026, 0.025, 0.023, 0.021, 0.020, 0.018, 0.017, 0.016, 0.015, 0.013, 0.012, 0.011, 0.010, 0.009,
1132 0.008, 0.007, 0.007, 0.006, 0.005, 0.004, 0.004, 0.003, 0.002, 0.002, 0.001, 0.001, 0.001, 0.000,
1133 -0.000, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003, -0.003,
1134 -0.003, -0.003, -0.003, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004,
1135 -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1136 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.002, -0.002, -0.002, -0.002,
1137 -0.002, -0.002, -0.002, -0.002},
1138 {// charge = 84 fC
1139 0.000, -0.000, 0.003, 0.014, 0.031, 0.051, 0.073, 0.097, 0.121, 0.145, 0.168, 0.191, 0.213, 0.234,
1140 0.253, 0.271, 0.288, 0.304, 0.319, 0.332, 0.344, 0.355, 0.365, 0.374, 0.382, 0.389, 0.395, 0.399,
1141 0.403, 0.407, 0.409, 0.410, 0.411, 0.412, 0.411, 0.410, 0.408, 0.406, 0.403, 0.400, 0.397, 0.393,
1142 0.388, 0.384, 0.379, 0.373, 0.368, 0.362, 0.356, 0.349, 0.343, 0.336, 0.330, 0.323, 0.316, 0.309,
1143 0.302, 0.295, 0.288, 0.280, 0.273, 0.266, 0.259, 0.252, 0.245, 0.238, 0.231, 0.224, 0.217, 0.210,
1144 0.204, 0.197, 0.191, 0.184, 0.178, 0.172, 0.166, 0.160, 0.154, 0.149, 0.143, 0.137, 0.132, 0.127,
1145 0.122, 0.117, 0.112, 0.108, 0.103, 0.099, 0.094, 0.090, 0.086, 0.082, 0.078, 0.075, 0.071, 0.067,
1146 0.064, 0.062, 0.058, 0.055, 0.052, 0.049, 0.046, 0.044, 0.041, 0.039, 0.037, 0.035, 0.032, 0.030,
1147 0.028, 0.026, 0.025, 0.023, 0.021, 0.020, 0.018, 0.017, 0.016, 0.014, 0.013, 0.012, 0.011, 0.010,
1148 0.009, 0.008, 0.007, 0.006, 0.005, 0.005, 0.004, 0.003, 0.003, 0.002, 0.002, 0.001, 0.001, 0.000,
1149 -0.000, -0.001, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1150 -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004,
1151 -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.003, -0.003, -0.003,
1152 -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.002,
1153 -0.002, -0.002, -0.002, -0.002},
1154 {// charge = 90 fC
1155 0.000, -0.000, 0.004, 0.015, 0.033, 0.054, 0.078, 0.104, 0.129, 0.155, 0.180, 0.205, 0.228, 0.250,
1156 0.271, 0.291, 0.309, 0.326, 0.341, 0.356, 0.369, 0.381, 0.391, 0.401, 0.409, 0.416, 0.423, 0.428,
1157 0.432, 0.436, 0.438, 0.440, 0.441, 0.441, 0.440, 0.439, 0.438, 0.435, 0.432, 0.429, 0.425, 0.421,
1158 0.416, 0.411, 0.406, 0.400, 0.394, 0.388, 0.381, 0.374, 0.368, 0.360, 0.353, 0.346, 0.339, 0.331,
1159 0.323, 0.316, 0.308, 0.301, 0.293, 0.285, 0.277, 0.270, 0.262, 0.255, 0.247, 0.240, 0.233, 0.225,
1160 0.218, 0.211, 0.204, 0.198, 0.191, 0.184, 0.178, 0.171, 0.165, 0.159, 0.153, 0.147, 0.142, 0.136,
1161 0.131, 0.125, 0.120, 0.115, 0.110, 0.105, 0.101, 0.096, 0.092, 0.088, 0.084, 0.080, 0.076, 0.072,
1162 0.069, 0.066, 0.062, 0.059, 0.056, 0.053, 0.050, 0.047, 0.044, 0.042, 0.039, 0.037, 0.035, 0.032,
1163 0.030, 0.029, 0.026, 0.025, 0.023, 0.021, 0.020, 0.018, 0.017, 0.015, 0.014, 0.013, 0.012, 0.011,
1164 0.009, 0.008, 0.007, 0.007, 0.006, 0.005, 0.004, 0.003, 0.003, 0.002, 0.002, 0.001, 0.001, 0.000,
1165 -0.000, -0.001, -0.001, -0.002, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003, -0.003, -0.004, -0.004,
1166 -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004,
1167 -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004,
1168 -0.004, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1169 -0.003, -0.003, -0.002, -0.002},
1170 {// charge = 96 fC
1171 0.000, -0.000, 0.004, 0.016, 0.035, 0.058, 0.083, 0.110, 0.138, 0.165, 0.192, 0.218, 0.243, 0.266,
1172 0.289, 0.310, 0.329, 0.347, 0.364, 0.379, 0.393, 0.406, 0.417, 0.427, 0.436, 0.444, 0.451, 0.456,
1173 0.461, 0.465, 0.467, 0.469, 0.470, 0.470, 0.470, 0.469, 0.467, 0.464, 0.461, 0.458, 0.453, 0.449,
1174 0.444, 0.439, 0.433, 0.427, 0.420, 0.414, 0.407, 0.399, 0.392, 0.385, 0.377, 0.369, 0.361, 0.353,
1175 0.345, 0.337, 0.329, 0.321, 0.312, 0.304, 0.296, 0.288, 0.280, 0.272, 0.264, 0.256, 0.248, 0.241,
1176 0.233, 0.225, 0.218, 0.211, 0.204, 0.197, 0.190, 0.183, 0.176, 0.170, 0.163, 0.157, 0.151, 0.145,
1177 0.139, 0.134, 0.128, 0.123, 0.118, 0.113, 0.108, 0.103, 0.098, 0.094, 0.089, 0.085, 0.081, 0.077,
1178 0.073, 0.070, 0.066, 0.063, 0.059, 0.056, 0.053, 0.050, 0.047, 0.045, 0.042, 0.039, 0.037, 0.035,
1179 0.032, 0.030, 0.028, 0.026, 0.025, 0.023, 0.021, 0.019, 0.018, 0.016, 0.015, 0.014, 0.012, 0.011,
1180 0.010, 0.009, 0.008, 0.007, 0.006, 0.005, 0.004, 0.004, 0.003, 0.002, 0.002, 0.001, 0.001, 0.000,
1181 -0.000, -0.001, -0.001, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003, -0.003, -0.004, -0.004, -0.004,
1182 -0.004, -0.004, -0.004, -0.004, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005,
1183 -0.005, -0.005, -0.005, -0.005, -0.005, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004,
1184 -0.004, -0.004, -0.004, -0.004, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1185 -0.003, -0.003, -0.003, -0.003},
1186 {// charge = 102 fC
1187 0.000, -0.000, 0.004, 0.017, 0.037, 0.061, 0.088, 0.117, 0.146, 0.175, 0.204, 0.231, 0.258, 0.283,
1188 0.307, 0.329, 0.350, 0.369, 0.387, 0.403, 0.418, 0.431, 0.443, 0.454, 0.464, 0.472, 0.479, 0.485,
1189 0.490, 0.494, 0.497, 0.498, 0.500, 0.500, 0.499, 0.498, 0.496, 0.493, 0.490, 0.486, 0.482, 0.477,
1190 0.472, 0.466, 0.460, 0.453, 0.447, 0.439, 0.432, 0.425, 0.417, 0.409, 0.401, 0.392, 0.384, 0.375,
1191 0.367, 0.358, 0.349, 0.341, 0.332, 0.323, 0.315, 0.306, 0.297, 0.289, 0.280, 0.272, 0.264, 0.256,
1192 0.247, 0.240, 0.232, 0.224, 0.216, 0.209, 0.202, 0.194, 0.187, 0.180, 0.174, 0.167, 0.161, 0.154,
1193 0.148, 0.142, 0.136, 0.131, 0.125, 0.120, 0.114, 0.109, 0.104, 0.100, 0.095, 0.091, 0.086, 0.082,
1194 0.078, 0.075, 0.070, 0.067, 0.063, 0.060, 0.056, 0.053, 0.050, 0.047, 0.045, 0.042, 0.039, 0.037,
1195 0.035, 0.032, 0.030, 0.028, 0.026, 0.024, 0.022, 0.021, 0.019, 0.017, 0.016, 0.015, 0.013, 0.012,
1196 0.011, 0.010, 0.008, 0.007, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.002, 0.001, 0.001, 0.000,
1197 -0.000, -0.001, -0.001, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003, -0.004, -0.004, -0.004, -0.004,
1198 -0.004, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005,
1199 -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.004, -0.004, -0.004, -0.004, -0.004,
1200 -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003, -0.003,
1201 -0.003, -0.003, -0.003, -0.003},
1202 {// charge = 108 fC
1203 0.000, -0.000, 0.004, 0.018, 0.039, 0.065, 0.093, 0.124, 0.155, 0.185, 0.215, 0.245, 0.273, 0.299,
1204 0.324, 0.348, 0.370, 0.390, 0.409, 0.427, 0.442, 0.457, 0.469, 0.481, 0.491, 0.500, 0.507, 0.513,
1205 0.519, 0.523, 0.526, 0.528, 0.529, 0.529, 0.529, 0.527, 0.525, 0.522, 0.519, 0.515, 0.510, 0.505,
1206 0.500, 0.493, 0.487, 0.480, 0.473, 0.465, 0.458, 0.450, 0.441, 0.433, 0.424, 0.415, 0.406, 0.397,
1207 0.388, 0.379, 0.370, 0.361, 0.352, 0.342, 0.333, 0.324, 0.315, 0.306, 0.297, 0.288, 0.279, 0.271,
1208 0.262, 0.254, 0.245, 0.237, 0.229, 0.221, 0.214, 0.206, 0.198, 0.191, 0.184, 0.177, 0.170, 0.163,
1209 0.157, 0.150, 0.144, 0.138, 0.132, 0.127, 0.121, 0.116, 0.111, 0.105, 0.101, 0.096, 0.091, 0.087,
1210 0.083, 0.079, 0.074, 0.070, 0.067, 0.063, 0.060, 0.056, 0.053, 0.050, 0.047, 0.044, 0.042, 0.039,
1211 0.036, 0.034, 0.032, 0.030, 0.027, 0.025, 0.024, 0.022, 0.020, 0.018, 0.017, 0.015, 0.014, 0.013,
1212 0.011, 0.010, 0.009, 0.008, 0.007, 0.006, 0.005, 0.004, 0.003, 0.003, 0.002, 0.001, 0.001, 0.000,
1213 -0.001, -0.001, -0.001, -0.002, -0.002, -0.003, -0.003, -0.003, -0.003, -0.004, -0.004, -0.004, -0.004, -0.005,
1214 -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005,
1215 -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.004,
1216 -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.003, -0.003, -0.003, -0.003,
1217 -0.003, -0.003, -0.003, -0.003},
1218 {// charge = 114 fC
1219 0.000, -0.000, 0.004, 0.019, 0.041, 0.068, 0.098, 0.130, 0.163, 0.195, 0.227, 0.258, 0.288, 0.316,
1220 0.342, 0.367, 0.390, 0.412, 0.432, 0.450, 0.467, 0.482, 0.495, 0.507, 0.518, 0.527, 0.535, 0.542,
1221 0.547, 0.552, 0.555, 0.557, 0.558, 0.559, 0.558, 0.557, 0.554, 0.551, 0.548, 0.544, 0.539, 0.533,
1222 0.527, 0.521, 0.514, 0.507, 0.499, 0.491, 0.483, 0.475, 0.466, 0.457, 0.448, 0.439, 0.429, 0.420,
1223 0.410, 0.400, 0.391, 0.381, 0.371, 0.362, 0.352, 0.342, 0.333, 0.323, 0.314, 0.304, 0.295, 0.286,
1224 0.277, 0.268, 0.259, 0.250, 0.242, 0.234, 0.225, 0.217, 0.209, 0.202, 0.194, 0.187, 0.179, 0.173,
1225 0.166, 0.159, 0.152, 0.146, 0.140, 0.134, 0.128, 0.122, 0.117, 0.112, 0.106, 0.101, 0.096, 0.092,
1226 0.087, 0.083, 0.078, 0.074, 0.071, 0.067, 0.063, 0.060, 0.056, 0.053, 0.050, 0.047, 0.044, 0.041,
1227 0.039, 0.036, 0.034, 0.031, 0.029, 0.027, 0.025, 0.023, 0.021, 0.019, 0.018, 0.016, 0.015, 0.013,
1228 0.012, 0.011, 0.009, 0.008, 0.007, 0.006, 0.005, 0.004, 0.003, 0.003, 0.002, 0.001, 0.001, 0.000,
1229 -0.001, -0.001, -0.002, -0.002, -0.002, -0.003, -0.003, -0.003, -0.004, -0.004, -0.004, -0.004, -0.005, -0.005,
1230 -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006,
1231 -0.006, -0.006, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005,
1232 -0.005, -0.005, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.003, -0.003,
1233 -0.003, -0.003, -0.003, -0.003},
1234 {// charge = 120 fC
1235 0.000, -0.000, 0.005, 0.020, 0.043, 0.072, 0.103, 0.137, 0.171, 0.205, 0.239, 0.271, 0.302, 0.332,
1236 0.360, 0.386, 0.411, 0.433, 0.454, 0.474, 0.491, 0.507, 0.521, 0.534, 0.545, 0.555, 0.563, 0.570,
1237 0.576, 0.581, 0.584, 0.586, 0.588, 0.588, 0.587, 0.586, 0.584, 0.581, 0.577, 0.572, 0.567, 0.561,
1238 0.555, 0.548, 0.541, 0.534, 0.526, 0.517, 0.509, 0.500, 0.491, 0.481, 0.472, 0.462, 0.452, 0.442,
1239 0.432, 0.422, 0.411, 0.401, 0.391, 0.381, 0.371, 0.360, 0.350, 0.340, 0.330, 0.320, 0.311, 0.301,
1240 0.291, 0.282, 0.273, 0.264, 0.255, 0.246, 0.237, 0.229, 0.220, 0.212, 0.204, 0.197, 0.189, 0.182,
1241 0.174, 0.167, 0.160, 0.154, 0.147, 0.141, 0.135, 0.129, 0.123, 0.117, 0.112, 0.107, 0.101, 0.097,
1242 0.092, 0.088, 0.082, 0.078, 0.074, 0.070, 0.066, 0.063, 0.059, 0.056, 0.052, 0.049, 0.046, 0.043,
1243 0.041, 0.038, 0.035, 0.033, 0.031, 0.028, 0.026, 0.024, 0.022, 0.021, 0.019, 0.017, 0.016, 0.014,
1244 0.013, 0.011, 0.010, 0.009, 0.008, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.001, 0.001, 0.000,
1245 -0.001, -0.001, -0.002, -0.002, -0.002, -0.003, -0.003, -0.004, -0.004, -0.004, -0.005, -0.005, -0.005, -0.005,
1246 -0.005, -0.005, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006,
1247 -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005,
1248 -0.005, -0.005, -0.005, -0.005, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004,
1249 -0.003, -0.003, -0.003, -0.003},
1250 {// charge = 126 fC
1251 0.000, -0.000, 0.005, 0.021, 0.045, 0.075, 0.108, 0.143, 0.179, 0.215, 0.250, 0.284, 0.317, 0.348,
1252 0.378, 0.405, 0.431, 0.455, 0.477, 0.497, 0.515, 0.532, 0.547, 0.560, 0.572, 0.583, 0.591, 0.599,
1253 0.605, 0.610, 0.613, 0.616, 0.617, 0.617, 0.617, 0.615, 0.613, 0.610, 0.606, 0.601, 0.596, 0.590,
1254 0.583, 0.576, 0.568, 0.560, 0.552, 0.543, 0.534, 0.525, 0.515, 0.505, 0.495, 0.485, 0.475, 0.464,
1255 0.453, 0.443, 0.432, 0.421, 0.411, 0.400, 0.389, 0.378, 0.368, 0.357, 0.347, 0.336, 0.326, 0.316,
1256 0.306, 0.296, 0.286, 0.277, 0.267, 0.258, 0.249, 0.240, 0.232, 0.223, 0.215, 0.206, 0.198, 0.191,
1257 0.183, 0.176, 0.169, 0.161, 0.155, 0.148, 0.141, 0.135, 0.129, 0.123, 0.118, 0.112, 0.107, 0.101,
1258 0.096, 0.092, 0.086, 0.082, 0.078, 0.074, 0.070, 0.066, 0.062, 0.058, 0.055, 0.052, 0.049, 0.045,
1259 0.043, 0.040, 0.037, 0.035, 0.032, 0.030, 0.027, 0.025, 0.023, 0.021, 0.020, 0.018, 0.016, 0.015,
1260 0.013, 0.012, 0.010, 0.009, 0.008, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.001, 0.001, 0.000,
1261 -0.001, -0.001, -0.002, -0.002, -0.003, -0.003, -0.003, -0.004, -0.004, -0.004, -0.005, -0.005, -0.005, -0.005,
1262 -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006,
1263 -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005, -0.005,
1264 -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004, -0.004,
1265 -0.004, -0.004, -0.003, -0.003},
1266 {// charge = 132 fC
1267 0.000, -0.000, 0.005, 0.022, 0.047, 0.078, 0.113, 0.150, 0.188, 0.225, 0.262, 0.298, 0.332, 0.365,
1268 0.395, 0.424, 0.451, 0.476, 0.499, 0.520, 0.540, 0.557, 0.573, 0.587, 0.599, 0.610, 0.619, 0.627,
1269 0.634, 0.639, 0.642, 0.645, 0.647, 0.647, 0.646, 0.645, 0.642, 0.639, 0.635, 0.630, 0.624, 0.618,
1270 0.611, 0.604, 0.596, 0.587, 0.578, 0.569, 0.560, 0.550, 0.540, 0.529, 0.519, 0.508, 0.497, 0.486,
1271 0.475, 0.464, 0.453, 0.441, 0.430, 0.419, 0.408, 0.397, 0.385, 0.374, 0.363, 0.353, 0.342, 0.331,
1272 0.321, 0.310, 0.300, 0.290, 0.280, 0.271, 0.261, 0.252, 0.243, 0.234, 0.225, 0.216, 0.208, 0.200,
1273 0.192, 0.184, 0.177, 0.169, 0.162, 0.155, 0.148, 0.142, 0.135, 0.129, 0.123, 0.117, 0.112, 0.106,
1274 0.101, 0.097, 0.091, 0.086, 0.082, 0.077, 0.073, 0.069, 0.065, 0.061, 0.058, 0.054, 0.051, 0.048,
1275 0.045, 0.042, 0.039, 0.036, 0.034, 0.031, 0.029, 0.027, 0.025, 0.022, 0.021, 0.019, 0.017, 0.015,
1276 0.014, 0.012, 0.011, 0.010, 0.008, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.002, 0.001, 0.000,
1277 -0.001, -0.001, -0.002, -0.002, -0.003, -0.003, -0.004, -0.004, -0.004, -0.005, -0.005, -0.005, -0.005, -0.006,
1278 -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.007, -0.007, -0.007, -0.007, -0.007, -0.006,
1279 -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.005,
1280 -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.004, -0.004, -0.004, -0.004, -0.004,
1281 -0.004, -0.004, -0.004, -0.004},
1282 {// charge = 138 fC
1283 0.000, -0.000, 0.005, 0.023, 0.049, 0.081, 0.118, 0.156, 0.196, 0.235, 0.274, 0.311, 0.347, 0.381,
1284 0.413, 0.443, 0.471, 0.497, 0.522, 0.544, 0.564, 0.582, 0.599, 0.614, 0.626, 0.638, 0.647, 0.656,
1285 0.662, 0.668, 0.672, 0.674, 0.676, 0.676, 0.676, 0.674, 0.671, 0.668, 0.663, 0.658, 0.653, 0.646,
1286 0.639, 0.631, 0.623, 0.614, 0.605, 0.595, 0.585, 0.575, 0.564, 0.554, 0.543, 0.531, 0.520, 0.509,
1287 0.497, 0.485, 0.473, 0.462, 0.450, 0.438, 0.426, 0.415, 0.403, 0.391, 0.380, 0.369, 0.357, 0.346,
1288 0.335, 0.325, 0.314, 0.303, 0.293, 0.283, 0.273, 0.263, 0.254, 0.244, 0.235, 0.226, 0.218, 0.209,
1289 0.201, 0.193, 0.185, 0.177, 0.169, 0.162, 0.155, 0.148, 0.141, 0.135, 0.129, 0.123, 0.117, 0.111,
1290 0.106, 0.101, 0.095, 0.090, 0.085, 0.081, 0.076, 0.072, 0.068, 0.064, 0.060, 0.057, 0.053, 0.050,
1291 0.047, 0.044, 0.041, 0.038, 0.035, 0.033, 0.030, 0.028, 0.026, 0.024, 0.021, 0.020, 0.018, 0.016,
1292 0.014, 0.013, 0.011, 0.010, 0.009, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.002, 0.001, 0.000,
1293 -0.001, -0.001, -0.002, -0.002, -0.003, -0.003, -0.004, -0.004, -0.005, -0.005, -0.005, -0.005, -0.006, -0.006,
1294 -0.006, -0.006, -0.006, -0.006, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007,
1295 -0.007, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006,
1296 -0.006, -0.006, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.004, -0.004, -0.004,
1297 -0.004, -0.004, -0.004, -0.004},
1298 {// charge = 144 fC
1299 0.000, -0.000, 0.006, 0.024, 0.051, 0.085, 0.123, 0.163, 0.204, 0.245, 0.285, 0.324, 0.361, 0.397,
1300 0.430, 0.462, 0.491, 0.519, 0.544, 0.567, 0.588, 0.607, 0.625, 0.640, 0.654, 0.665, 0.676, 0.684,
1301 0.691, 0.697, 0.701, 0.704, 0.705, 0.706, 0.705, 0.703, 0.701, 0.697, 0.692, 0.687, 0.681, 0.674,
1302 0.667, 0.659, 0.650, 0.641, 0.631, 0.621, 0.611, 0.600, 0.589, 0.578, 0.566, 0.555, 0.543, 0.531,
1303 0.519, 0.506, 0.494, 0.482, 0.470, 0.457, 0.445, 0.433, 0.421, 0.409, 0.397, 0.385, 0.373, 0.362,
1304 0.350, 0.339, 0.328, 0.317, 0.306, 0.295, 0.285, 0.275, 0.265, 0.255, 0.246, 0.236, 0.227, 0.218,
1305 0.210, 0.201, 0.193, 0.185, 0.177, 0.169, 0.162, 0.155, 0.148, 0.141, 0.134, 0.128, 0.122, 0.116,
1306 0.110, 0.105, 0.099, 0.094, 0.089, 0.084, 0.080, 0.075, 0.071, 0.067, 0.063, 0.059, 0.055, 0.052,
1307 0.049, 0.045, 0.042, 0.039, 0.037, 0.034, 0.031, 0.029, 0.027, 0.025, 0.022, 0.020, 0.019, 0.017,
1308 0.015, 0.013, 0.012, 0.011, 0.009, 0.008, 0.007, 0.006, 0.004, 0.003, 0.002, 0.002, 0.001, 0.000,
1309 -0.001, -0.001, -0.002, -0.002, -0.003, -0.003, -0.004, -0.004, -0.005, -0.005, -0.005, -0.006, -0.006, -0.006,
1310 -0.006, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007,
1311 -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006,
1312 -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.004,
1313 -0.004, -0.004, -0.004, -0.004},
1314 {// charge = 150 fC
1315 0.000, -0.000, 0.006, 0.025, 0.053, 0.088, 0.127, 0.169, 0.212, 0.255, 0.296, 0.337, 0.376, 0.413,
1316 0.448, 0.481, 0.512, 0.540, 0.566, 0.591, 0.613, 0.633, 0.650, 0.666, 0.681, 0.693, 0.704, 0.712,
1317 0.720, 0.726, 0.730, 0.733, 0.735, 0.735, 0.734, 0.733, 0.730, 0.726, 0.721, 0.716, 0.709, 0.702,
1318 0.695, 0.686, 0.677, 0.668, 0.658, 0.647, 0.637, 0.625, 0.614, 0.602, 0.590, 0.578, 0.566, 0.553,
1319 0.540, 0.528, 0.515, 0.502, 0.489, 0.476, 0.464, 0.451, 0.438, 0.426, 0.413, 0.401, 0.389, 0.377,
1320 0.365, 0.353, 0.341, 0.330, 0.319, 0.308, 0.297, 0.286, 0.276, 0.266, 0.256, 0.246, 0.237, 0.227,
1321 0.218, 0.209, 0.201, 0.192, 0.184, 0.176, 0.169, 0.161, 0.154, 0.147, 0.140, 0.133, 0.127, 0.121,
1322 0.115, 0.110, 0.103, 0.098, 0.093, 0.088, 0.083, 0.078, 0.074, 0.070, 0.066, 0.062, 0.058, 0.054,
1323 0.051, 0.047, 0.044, 0.041, 0.038, 0.035, 0.033, 0.030, 0.028, 0.026, 0.023, 0.021, 0.019, 0.017,
1324 0.016, 0.014, 0.012, 0.011, 0.009, 0.008, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.001, 0.000,
1325 -0.001, -0.001, -0.002, -0.003, -0.003, -0.004, -0.004, -0.005, -0.005, -0.005, -0.006, -0.006, -0.006, -0.006,
1326 -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007,
1327 -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006,
1328 -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005,
1329 -0.005, -0.004, -0.004, -0.004},
1330 {// charge = 156 fC
1331 0.000, -0.001, 0.006, 0.025, 0.055, 0.091, 0.132, 0.176, 0.220, 0.264, 0.308, 0.350, 0.390, 0.429,
1332 0.465, 0.500, 0.532, 0.561, 0.589, 0.614, 0.637, 0.658, 0.676, 0.693, 0.708, 0.721, 0.732, 0.741,
1333 0.749, 0.755, 0.759, 0.762, 0.764, 0.765, 0.764, 0.762, 0.759, 0.755, 0.750, 0.744, 0.738, 0.731,
1334 0.722, 0.714, 0.704, 0.695, 0.684, 0.673, 0.662, 0.651, 0.639, 0.626, 0.614, 0.601, 0.588, 0.575,
1335 0.562, 0.549, 0.536, 0.522, 0.509, 0.496, 0.482, 0.469, 0.456, 0.443, 0.430, 0.417, 0.404, 0.392,
1336 0.379, 0.367, 0.355, 0.343, 0.332, 0.320, 0.309, 0.298, 0.287, 0.276, 0.266, 0.256, 0.246, 0.237,
1337 0.227, 0.218, 0.209, 0.200, 0.192, 0.183, 0.175, 0.168, 0.160, 0.153, 0.146, 0.139, 0.132, 0.126,
1338 0.119, 0.114, 0.107, 0.102, 0.097, 0.091, 0.086, 0.082, 0.077, 0.072, 0.068, 0.064, 0.060, 0.056,
1339 0.053, 0.049, 0.046, 0.043, 0.040, 0.037, 0.034, 0.031, 0.029, 0.027, 0.024, 0.022, 0.020, 0.018,
1340 0.016, 0.015, 0.013, 0.011, 0.010, 0.008, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.001, 0.000,
1341 -0.001, -0.002, -0.002, -0.003, -0.003, -0.004, -0.004, -0.005, -0.005, -0.006, -0.006, -0.006, -0.006, -0.007,
1342 -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008,
1343 -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.006,
1344 -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005, -0.005, -0.005, -0.005,
1345 -0.005, -0.005, -0.004, -0.004},
1346 {// charge = 162 fC
1347 0.000, -0.001, 0.006, 0.026, 0.057, 0.095, 0.137, 0.182, 0.228, 0.274, 0.319, 0.363, 0.405, 0.445,
1348 0.483, 0.518, 0.552, 0.582, 0.611, 0.637, 0.661, 0.683, 0.702, 0.719, 0.735, 0.748, 0.760, 0.769,
1349 0.777, 0.784, 0.788, 0.792, 0.793, 0.794, 0.793, 0.791, 0.788, 0.784, 0.779, 0.773, 0.766, 0.759,
1350 0.750, 0.741, 0.732, 0.721, 0.711, 0.699, 0.688, 0.676, 0.663, 0.651, 0.638, 0.624, 0.611, 0.598,
1351 0.584, 0.570, 0.556, 0.543, 0.529, 0.515, 0.501, 0.487, 0.474, 0.460, 0.447, 0.433, 0.420, 0.407,
1352 0.394, 0.382, 0.369, 0.357, 0.345, 0.333, 0.321, 0.309, 0.298, 0.287, 0.276, 0.266, 0.256, 0.246,
1353 0.236, 0.226, 0.217, 0.208, 0.199, 0.191, 0.182, 0.174, 0.166, 0.159, 0.151, 0.144, 0.137, 0.131,
1354 0.124, 0.119, 0.112, 0.106, 0.100, 0.095, 0.090, 0.085, 0.080, 0.075, 0.071, 0.067, 0.062, 0.059,
1355 0.055, 0.051, 0.048, 0.044, 0.041, 0.038, 0.035, 0.033, 0.030, 0.028, 0.025, 0.023, 0.021, 0.019,
1356 0.017, 0.015, 0.013, 0.012, 0.010, 0.009, 0.007, 0.006, 0.005, 0.004, 0.003, 0.002, 0.001, 0.000,
1357 -0.001, -0.002, -0.002, -0.003, -0.003, -0.004, -0.005, -0.005, -0.005, -0.006, -0.006, -0.006, -0.007, -0.007,
1358 -0.007, -0.007, -0.007, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008,
1359 -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007,
1360 -0.007, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005, -0.005, -0.005,
1361 -0.005, -0.005, -0.005, -0.005},
1362 {// charge = 168 fC
1363 0.000, -0.001, 0.006, 0.027, 0.059, 0.098, 0.142, 0.188, 0.236, 0.284, 0.330, 0.376, 0.419, 0.461,
1364 0.500, 0.537, 0.572, 0.604, 0.633, 0.660, 0.685, 0.708, 0.728, 0.746, 0.762, 0.776, 0.788, 0.798,
1365 0.806, 0.812, 0.817, 0.821, 0.823, 0.823, 0.823, 0.821, 0.818, 0.813, 0.808, 0.802, 0.795, 0.787,
1366 0.778, 0.769, 0.759, 0.748, 0.737, 0.725, 0.713, 0.701, 0.688, 0.675, 0.661, 0.648, 0.634, 0.620,
1367 0.606, 0.592, 0.577, 0.563, 0.549, 0.534, 0.520, 0.506, 0.491, 0.477, 0.463, 0.450, 0.436, 0.422,
1368 0.409, 0.396, 0.383, 0.370, 0.357, 0.345, 0.333, 0.321, 0.309, 0.298, 0.287, 0.276, 0.265, 0.255,
1369 0.245, 0.235, 0.225, 0.216, 0.207, 0.198, 0.189, 0.181, 0.173, 0.165, 0.157, 0.150, 0.142, 0.136,
1370 0.129, 0.123, 0.116, 0.110, 0.104, 0.099, 0.093, 0.088, 0.083, 0.078, 0.074, 0.069, 0.065, 0.061,
1371 0.057, 0.053, 0.050, 0.046, 0.043, 0.040, 0.037, 0.034, 0.031, 0.029, 0.026, 0.024, 0.022, 0.020,
1372 0.018, 0.016, 0.014, 0.012, 0.011, 0.009, 0.008, 0.006, 0.005, 0.004, 0.003, 0.002, 0.001, 0.000,
1373 -0.001, -0.002, -0.002, -0.003, -0.004, -0.004, -0.005, -0.005, -0.006, -0.006, -0.006, -0.007, -0.007, -0.007,
1374 -0.007, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008,
1375 -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007,
1376 -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005,
1377 -0.005, -0.005, -0.005, -0.005},
1378 {// charge = 180 fC
1379 0.000, -0.001, 0.007, 0.029, 0.062, 0.104, 0.151, 0.201, 0.251, 0.302, 0.353, 0.401, 0.448, 0.493,
1380 0.535, 0.574, 0.611, 0.646, 0.677, 0.707, 0.733, 0.757, 0.779, 0.798, 0.816, 0.831, 0.843, 0.854,
1381 0.863, 0.870, 0.876, 0.879, 0.882, 0.882, 0.881, 0.879, 0.876, 0.872, 0.866, 0.859, 0.852, 0.843,
1382 0.834, 0.824, 0.813, 0.802, 0.790, 0.778, 0.765, 0.751, 0.738, 0.723, 0.709, 0.694, 0.680, 0.665,
1383 0.649, 0.634, 0.619, 0.603, 0.588, 0.573, 0.557, 0.542, 0.527, 0.512, 0.497, 0.482, 0.467, 0.453,
1384 0.439, 0.424, 0.410, 0.397, 0.383, 0.370, 0.357, 0.344, 0.332, 0.320, 0.308, 0.296, 0.284, 0.273,
1385 0.262, 0.252, 0.241, 0.231, 0.221, 0.212, 0.203, 0.194, 0.185, 0.177, 0.168, 0.160, 0.153, 0.145,
1386 0.138, 0.132, 0.124, 0.118, 0.112, 0.106, 0.100, 0.094, 0.089, 0.084, 0.079, 0.074, 0.070, 0.065,
1387 0.061, 0.057, 0.053, 0.049, 0.046, 0.043, 0.039, 0.036, 0.034, 0.031, 0.028, 0.026, 0.023, 0.021,
1388 0.019, 0.017, 0.015, 0.013, 0.011, 0.010, 0.008, 0.007, 0.006, 0.004, 0.003, 0.002, 0.001, 0.000,
1389 -0.001, -0.002, -0.002, -0.003, -0.004, -0.005, -0.005, -0.006, -0.006, -0.006, -0.007, -0.007, -0.007, -0.008,
1390 -0.008, -0.008, -0.008, -0.008, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009,
1391 -0.009, -0.009, -0.009, -0.009, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.007,
1392 -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006, -0.006,
1393 -0.005, -0.005, -0.005, -0.005},
1394 {// charge = 195 fC
1395 0.000, -0.001, 0.007, 0.031, 0.067, 0.112, 0.162, 0.216, 0.271, 0.326, 0.380, 0.433, 0.484, 0.532,
1396 0.578, 0.621, 0.661, 0.698, 0.733, 0.764, 0.793, 0.819, 0.843, 0.864, 0.883, 0.899, 0.913, 0.925,
1397 0.935, 0.943, 0.948, 0.953, 0.955, 0.956, 0.955, 0.953, 0.949, 0.944, 0.938, 0.931, 0.923, 0.914,
1398 0.904, 0.893, 0.882, 0.869, 0.856, 0.843, 0.829, 0.814, 0.800, 0.784, 0.769, 0.753, 0.737, 0.721,
1399 0.704, 0.688, 0.671, 0.654, 0.638, 0.621, 0.604, 0.588, 0.571, 0.555, 0.539, 0.523, 0.507, 0.491,
1400 0.476, 0.460, 0.445, 0.430, 0.416, 0.401, 0.387, 0.373, 0.360, 0.347, 0.334, 0.321, 0.308, 0.296,
1401 0.285, 0.273, 0.262, 0.251, 0.240, 0.230, 0.220, 0.210, 0.201, 0.192, 0.183, 0.174, 0.166, 0.158,
1402 0.150, 0.143, 0.135, 0.128, 0.121, 0.115, 0.108, 0.102, 0.096, 0.091, 0.086, 0.080, 0.076, 0.071,
1403 0.066, 0.062, 0.058, 0.054, 0.050, 0.046, 0.043, 0.039, 0.036, 0.033, 0.030, 0.028, 0.025, 0.023,
1404 0.020, 0.018, 0.016, 0.014, 0.012, 0.011, 0.009, 0.007, 0.006, 0.005, 0.003, 0.002, 0.001, 0.000,
1405 -0.001, -0.002, -0.003, -0.003, -0.004, -0.005, -0.006, -0.006, -0.007, -0.007, -0.007, -0.008, -0.008, -0.008,
1406 -0.009, -0.009, -0.009, -0.009, -0.009, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010,
1407 -0.010, -0.010, -0.010, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009, -0.009, -0.008, -0.008, -0.008,
1408 -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006,
1409 -0.006, -0.006, -0.006, -0.006},
1410 {// charge = 210 fC
1411 0.000, -0.001, 0.008, 0.033, 0.071, 0.119, 0.173, 0.230, 0.290, 0.349, 0.407, 0.464, 0.519, 0.571,
1412 0.620, 0.666, 0.710, 0.750, 0.787, 0.821, 0.853, 0.881, 0.907, 0.930, 0.950, 0.967, 0.983, 0.996,
1413 1.006, 1.015, 1.021, 1.026, 1.028, 1.029, 1.028, 1.026, 1.022, 1.017, 1.011, 1.003, 0.995, 0.985,
1414 0.974, 0.962, 0.950, 0.937, 0.923, 0.908, 0.893, 0.878, 0.862, 0.845, 0.828, 0.811, 0.794, 0.777,
1415 0.759, 0.741, 0.723, 0.705, 0.687, 0.669, 0.652, 0.634, 0.616, 0.598, 0.581, 0.564, 0.546, 0.529,
1416 0.513, 0.496, 0.480, 0.464, 0.448, 0.433, 0.418, 0.403, 0.388, 0.374, 0.360, 0.346, 0.333, 0.320,
1417 0.307, 0.294, 0.282, 0.270, 0.259, 0.248, 0.237, 0.227, 0.216, 0.206, 0.197, 0.188, 0.179, 0.170,
1418 0.161, 0.154, 0.145, 0.138, 0.131, 0.124, 0.117, 0.110, 0.104, 0.098, 0.092, 0.087, 0.081, 0.076,
1419 0.071, 0.067, 0.062, 0.058, 0.054, 0.050, 0.046, 0.043, 0.039, 0.036, 0.033, 0.030, 0.027, 0.025,
1420 0.022, 0.020, 0.017, 0.015, 0.013, 0.011, 0.010, 0.008, 0.006, 0.005, 0.004, 0.002, 0.001, 0.000,
1421 -0.001, -0.002, -0.003, -0.004, -0.005, -0.005, -0.006, -0.007, -0.007, -0.007, -0.008, -0.008, -0.009, -0.009,
1422 -0.009, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.011, -0.011, -0.011, -0.011, -0.011, -0.010,
1423 -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.009, -0.009, -0.009, -0.009, -0.009,
1424 -0.009, -0.008, -0.008, -0.008, -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007, -0.007,
1425 -0.006, -0.006, -0.006, -0.006},
1426 {// charge = 225 fC
1427 0.000, -0.001, 0.008, 0.035, 0.075, 0.126, 0.183, 0.245, 0.308, 0.371, 0.434, 0.495, 0.553, 0.609,
1428 0.662, 0.712, 0.758, 0.801, 0.842, 0.878, 0.912, 0.943, 0.970, 0.995, 1.017, 1.036, 1.052, 1.066,
1429 1.077, 1.087, 1.094, 1.098, 1.101, 1.102, 1.102, 1.099, 1.095, 1.090, 1.083, 1.075, 1.066, 1.055,
1430 1.044, 1.032, 1.018, 1.004, 0.989, 0.974, 0.958, 0.941, 0.924, 0.906, 0.888, 0.870, 0.852, 0.833,
1431 0.814, 0.795, 0.776, 0.756, 0.737, 0.718, 0.699, 0.680, 0.661, 0.642, 0.623, 0.605, 0.586, 0.568,
1432 0.550, 0.532, 0.515, 0.498, 0.481, 0.464, 0.448, 0.432, 0.416, 0.401, 0.386, 0.371, 0.357, 0.343,
1433 0.329, 0.316, 0.303, 0.290, 0.278, 0.266, 0.254, 0.243, 0.232, 0.222, 0.211, 0.201, 0.192, 0.183,
1434 0.173, 0.165, 0.156, 0.148, 0.140, 0.133, 0.126, 0.118, 0.112, 0.105, 0.099, 0.093, 0.087, 0.082,
1435 0.077, 0.072, 0.067, 0.062, 0.058, 0.053, 0.049, 0.046, 0.042, 0.039, 0.035, 0.032, 0.029, 0.026,
1436 0.024, 0.021, 0.019, 0.016, 0.014, 0.012, 0.010, 0.009, 0.007, 0.005, 0.004, 0.002, 0.001, 0.000,
1437 -0.001, -0.002, -0.003, -0.004, -0.005, -0.006, -0.006, -0.007, -0.008, -0.008, -0.009, -0.009, -0.009, -0.010,
1438 -0.010, -0.010, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011,
1439 -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010, -0.010,
1440 -0.009, -0.009, -0.009, -0.009, -0.009, -0.008, -0.008, -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007,
1441 -0.007, -0.007, -0.007, -0.006},
1442 {// charge = 240 fC
1443 0.000, -0.001, 0.009, 0.036, 0.079, 0.132, 0.193, 0.258, 0.325, 0.393, 0.460, 0.524, 0.587, 0.647,
1444 0.703, 0.756, 0.806, 0.852, 0.895, 0.935, 0.971, 1.004, 1.033, 1.060, 1.083, 1.103, 1.121, 1.136,
1445 1.148, 1.158, 1.166, 1.171, 1.174, 1.175, 1.175, 1.172, 1.168, 1.162, 1.155, 1.147, 1.137, 1.126,
1446 1.114, 1.101, 1.087, 1.072, 1.056, 1.039, 1.022, 1.004, 0.986, 0.968, 0.948, 0.929, 0.909, 0.889,
1447 0.869, 0.849, 0.828, 0.808, 0.788, 0.767, 0.746, 0.726, 0.706, 0.686, 0.666, 0.646, 0.626, 0.607,
1448 0.588, 0.569, 0.550, 0.532, 0.514, 0.496, 0.479, 0.462, 0.445, 0.428, 0.412, 0.397, 0.381, 0.366,
1449 0.352, 0.338, 0.324, 0.310, 0.297, 0.284, 0.272, 0.260, 0.248, 0.237, 0.226, 0.215, 0.205, 0.196,
1450 0.185, 0.176, 0.167, 0.158, 0.150, 0.142, 0.134, 0.127, 0.119, 0.113, 0.106, 0.100, 0.093, 0.088,
1451 0.082, 0.077, 0.071, 0.067, 0.062, 0.057, 0.053, 0.049, 0.045, 0.041, 0.038, 0.034, 0.031, 0.028,
1452 0.025, 0.023, 0.020, 0.018, 0.015, 0.013, 0.011, 0.009, 0.007, 0.006, 0.004, 0.003, 0.001, 0.000,
1453 -0.001, -0.002, -0.003, -0.004, -0.005, -0.006, -0.007, -0.007, -0.008, -0.009, -0.009, -0.010, -0.010, -0.010,
1454 -0.011, -0.011, -0.011, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012,
1455 -0.012, -0.012, -0.012, -0.012, -0.012, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.011, -0.010, -0.010,
1456 -0.010, -0.010, -0.010, -0.009, -0.009, -0.009, -0.009, -0.009, -0.008, -0.008, -0.008, -0.008, -0.008, -0.007,
1457 -0.007, -0.007, -0.007, -0.007},
1458 {// charge = 255 fC
1459 0.000, -0.001, 0.009, 0.038, 0.082, 0.139, 0.203, 0.271, 0.343, 0.414, 0.485, 0.554, 0.620, 0.683,
1460 0.744, 0.800, 0.853, 0.903, 0.949, 0.991, 1.029, 1.064, 1.096, 1.124, 1.149, 1.170, 1.189, 1.204,
1461 1.217, 1.227, 1.234, 1.240, 1.243, 1.244, 1.243, 1.241, 1.237, 1.231, 1.224, 1.215, 1.205, 1.194,
1462 1.182, 1.168, 1.153, 1.138, 1.121, 1.104, 1.086, 1.067, 1.048, 1.028, 1.008, 0.988, 0.967, 0.946,
1463 0.924, 0.903, 0.881, 0.860, 0.838, 0.816, 0.795, 0.773, 0.751, 0.730, 0.709, 0.688, 0.667, 0.646,
1464 0.626, 0.606, 0.586, 0.567, 0.548, 0.529, 0.510, 0.492, 0.474, 0.457, 0.440, 0.423, 0.407, 0.391,
1465 0.375, 0.360, 0.346, 0.331, 0.317, 0.304, 0.290, 0.277, 0.265, 0.253, 0.241, 0.230, 0.219, 0.209,
1466 0.198, 0.188, 0.178, 0.169, 0.160, 0.152, 0.143, 0.136, 0.128, 0.120, 0.113, 0.107, 0.100, 0.094,
1467 0.088, 0.082, 0.076, 0.071, 0.066, 0.061, 0.057, 0.053, 0.048, 0.044, 0.041, 0.037, 0.034, 0.030,
1468 0.027, 0.024, 0.022, 0.019, 0.017, 0.014, 0.012, 0.010, 0.008, 0.006, 0.005, 0.003, 0.002, 0.000,
1469 -0.001, -0.002, -0.003, -0.005, -0.006, -0.006, -0.007, -0.008, -0.009, -0.009, -0.010, -0.010, -0.011, -0.011,
1470 -0.011, -0.012, -0.012, -0.012, -0.012, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013,
1471 -0.013, -0.013, -0.013, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012, -0.011, -0.011, -0.011, -0.011,
1472 -0.011, -0.011, -0.010, -0.010, -0.010, -0.010, -0.010, -0.009, -0.009, -0.009, -0.009, -0.008, -0.008, -0.008,
1473 -0.008, -0.008, -0.007, -0.007},
1474 {// charge = 270 fC
1475 0.000, -0.001, 0.010, 0.039, 0.086, 0.145, 0.212, 0.284, 0.359, 0.434, 0.509, 0.582, 0.652, 0.720,
1476 0.783, 0.844, 0.900, 0.953, 1.001, 1.046, 1.087, 1.124, 1.158, 1.187, 1.212, 1.234, 1.252, 1.267,
1477 1.279, 1.289, 1.295, 1.300, 1.302, 1.303, 1.302, 1.299, 1.295, 1.289, 1.282, 1.274, 1.264, 1.253,
1478 1.241, 1.228, 1.213, 1.198, 1.181, 1.164, 1.146, 1.127, 1.107, 1.087, 1.066, 1.045, 1.023, 1.001,
1479 0.979, 0.957, 0.934, 0.911, 0.889, 0.866, 0.843, 0.821, 0.798, 0.776, 0.753, 0.731, 0.709, 0.688,
1480 0.666, 0.645, 0.624, 0.604, 0.584, 0.564, 0.544, 0.525, 0.506, 0.488, 0.470, 0.452, 0.435, 0.418,
1481 0.401, 0.385, 0.370, 0.354, 0.340, 0.325, 0.311, 0.297, 0.284, 0.271, 0.259, 0.247, 0.235, 0.225,
1482 0.213, 0.202, 0.192, 0.182, 0.173, 0.164, 0.155, 0.146, 0.138, 0.130, 0.123, 0.115, 0.108, 0.102,
1483 0.095, 0.089, 0.083, 0.077, 0.072, 0.067, 0.062, 0.057, 0.053, 0.049, 0.044, 0.041, 0.037, 0.033,
1484 0.030, 0.027, 0.024, 0.021, 0.019, 0.016, 0.014, 0.012, 0.009, 0.007, 0.006, 0.004, 0.002, 0.001,
1485 -0.001, -0.002, -0.003, -0.004, -0.005, -0.006, -0.007, -0.008, -0.009, -0.009, -0.010, -0.011, -0.011, -0.012,
1486 -0.012, -0.012, -0.012, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.014, -0.014, -0.014, -0.014, -0.014,
1487 -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.012, -0.012, -0.012, -0.012, -0.012, -0.012,
1488 -0.011, -0.011, -0.011, -0.011, -0.011, -0.010, -0.010, -0.010, -0.010, -0.010, -0.009, -0.009, -0.009, -0.009,
1489 -0.008, -0.008, -0.008, -0.008},
1490 {// charge = 285 fC
1491 0.000, -0.001, 0.010, 0.041, 0.089, 0.150, 0.220, 0.296, 0.374, 0.454, 0.533, 0.609, 0.684, 0.755,
1492 0.822, 0.886, 0.946, 1.002, 1.053, 1.101, 1.144, 1.183, 1.217, 1.246, 1.271, 1.291, 1.308, 1.321,
1493 1.331, 1.338, 1.343, 1.346, 1.348, 1.347, 1.346, 1.342, 1.338, 1.332, 1.326, 1.318, 1.309, 1.299,
1494 1.287, 1.275, 1.261, 1.247, 1.231, 1.214, 1.197, 1.179, 1.159, 1.139, 1.119, 1.097, 1.076, 1.053,
1495 1.031, 1.008, 0.985, 0.962, 0.939, 0.916, 0.893, 0.869, 0.846, 0.823, 0.800, 0.777, 0.754, 0.732,
1496 0.709, 0.687, 0.666, 0.644, 0.623, 0.602, 0.582, 0.562, 0.542, 0.523, 0.504, 0.485, 0.467, 0.449,
1497 0.432, 0.415, 0.398, 0.382, 0.366, 0.351, 0.336, 0.321, 0.307, 0.294, 0.280, 0.268, 0.255, 0.244,
1498 0.231, 0.220, 0.209, 0.199, 0.188, 0.179, 0.169, 0.160, 0.151, 0.143, 0.135, 0.127, 0.119, 0.112,
1499 0.105, 0.098, 0.092, 0.086, 0.080, 0.075, 0.069, 0.064, 0.059, 0.055, 0.050, 0.046, 0.042, 0.038,
1500 0.035, 0.031, 0.028, 0.025, 0.022, 0.019, 0.017, 0.014, 0.012, 0.010, 0.008, 0.006, 0.004, 0.002,
1501 0.001, -0.001, -0.002, -0.003, -0.004, -0.006, -0.007, -0.007, -0.008, -0.009, -0.010, -0.010, -0.011, -0.012,
1502 -0.012, -0.012, -0.013, -0.013, -0.013, -0.013, -0.014, -0.014, -0.014, -0.014, -0.014, -0.014, -0.014, -0.014,
1503 -0.014, -0.014, -0.014, -0.014, -0.014, -0.014, -0.013, -0.013, -0.013, -0.013, -0.013, -0.013, -0.012, -0.012,
1504 -0.012, -0.012, -0.012, -0.011, -0.011, -0.011, -0.011, -0.011, -0.010, -0.010, -0.010, -0.010, -0.010, -0.009,
1505 -0.009, -0.009, -0.009, -0.008},
1506 {// charge = 300 fC
1507 0.000, -0.001, 0.010, 0.042, 0.092, 0.155, 0.228, 0.307, 0.389, 0.472, 0.555, 0.636, 0.714, 0.789,
1508 0.860, 0.928, 0.991, 1.050, 1.104, 1.154, 1.199, 1.238, 1.271, 1.298, 1.320, 1.338, 1.352, 1.362,
1509 1.369, 1.374, 1.377, 1.379, 1.379, 1.378, 1.375, 1.372, 1.368, 1.363, 1.356, 1.349, 1.341, 1.332,
1510 1.322, 1.311, 1.299, 1.285, 1.271, 1.256, 1.240, 1.223, 1.205, 1.186, 1.166, 1.146, 1.125, 1.103,
1511 1.081, 1.058, 1.035, 1.012, 0.989, 0.965, 0.942, 0.918, 0.895, 0.871, 0.848, 0.824, 0.801, 0.778,
1512 0.755, 0.732, 0.710, 0.688, 0.666, 0.644, 0.623, 0.602, 0.581, 0.561, 0.541, 0.522, 0.503, 0.484,
1513 0.466, 0.448, 0.430, 0.413, 0.397, 0.381, 0.365, 0.349, 0.334, 0.320, 0.306, 0.292, 0.279, 0.266,
1514 0.254, 0.242, 0.230, 0.218, 0.208, 0.197, 0.187, 0.177, 0.168, 0.158, 0.150, 0.141, 0.133, 0.125,
1515 0.118, 0.110, 0.104, 0.097, 0.091, 0.085, 0.079, 0.073, 0.068, 0.063, 0.058, 0.053, 0.049, 0.045,
1516 0.041, 0.037, 0.034, 0.030, 0.027, 0.024, 0.021, 0.018, 0.016, 0.013, 0.011, 0.009, 0.007, 0.005,
1517 0.003, 0.002, 0.000, -0.001, -0.003, -0.004, -0.005, -0.006, -0.007, -0.008, -0.009, -0.010, -0.010, -0.011,
1518 -0.012, -0.012, -0.012, -0.013, -0.013, -0.013, -0.014, -0.014, -0.014, -0.014, -0.015, -0.015, -0.015, -0.015,
1519 -0.015, -0.015, -0.015, -0.015, -0.014, -0.014, -0.014, -0.014, -0.014, -0.014, -0.014, -0.013, -0.013, -0.013,
1520 -0.013, -0.013, -0.012, -0.012, -0.012, -0.012, -0.012, -0.011, -0.011, -0.011, -0.011, -0.010, -0.010, -0.010,
1521 -0.010, -0.010, -0.009, -0.009},
1522 {// charge = 315 fC
1523 0.000, -0.001, 0.011, 0.043, 0.094, 0.159, 0.235, 0.317, 0.403, 0.490, 0.576, 0.661, 0.743, 0.822,
1524 0.897, 0.968, 1.035, 1.097, 1.154, 1.205, 1.249, 1.287, 1.317, 1.341, 1.359, 1.373, 1.383, 1.390,
1525 1.395, 1.398, 1.400, 1.400, 1.400, 1.398, 1.396, 1.393, 1.389, 1.384, 1.378, 1.372, 1.365, 1.357,
1526 1.348, 1.339, 1.328, 1.317, 1.304, 1.291, 1.276, 1.261, 1.244, 1.227, 1.209, 1.190, 1.170, 1.149,
1527 1.128, 1.106, 1.084, 1.061, 1.038, 1.015, 0.992, 0.968, 0.944, 0.921, 0.897, 0.873, 0.850, 0.826,
1528 0.803, 0.779, 0.756, 0.734, 0.711, 0.689, 0.667, 0.645, 0.624, 0.603, 0.582, 0.562, 0.542, 0.522,
1529 0.503, 0.485, 0.466, 0.448, 0.431, 0.414, 0.397, 0.381, 0.365, 0.350, 0.335, 0.320, 0.306, 0.292,
1530 0.279, 0.267, 0.253, 0.241, 0.230, 0.218, 0.207, 0.197, 0.187, 0.177, 0.167, 0.158, 0.149, 0.141,
1531 0.133, 0.125, 0.117, 0.110, 0.103, 0.097, 0.090, 0.084, 0.078, 0.073, 0.067, 0.062, 0.058, 0.053,
1532 0.049, 0.044, 0.040, 0.037, 0.033, 0.030, 0.026, 0.023, 0.021, 0.018, 0.015, 0.013, 0.011, 0.008,
1533 0.006, 0.004, 0.003, 0.001, -0.000, -0.002, -0.003, -0.004, -0.006, -0.007, -0.008, -0.008, -0.009, -0.010,
1534 -0.011, -0.011, -0.012, -0.012, -0.013, -0.013, -0.014, -0.014, -0.014, -0.014, -0.015, -0.015, -0.015, -0.015,
1535 -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.014, -0.014, -0.014, -0.014, -0.014,
1536 -0.013, -0.013, -0.013, -0.013, -0.013, -0.012, -0.012, -0.012, -0.012, -0.012, -0.011, -0.011, -0.011, -0.011,
1537 -0.011, -0.010, -0.010, -0.010},
1538 {// charge = 330 fC
1539 0.000, -0.001, 0.011, 0.044, 0.096, 0.163, 0.241, 0.326, 0.415, 0.506, 0.596, 0.685, 0.771, 0.854,
1540 0.933, 1.007, 1.077, 1.142, 1.201, 1.252, 1.294, 1.328, 1.353, 1.373, 1.387, 1.397, 1.404, 1.409,
1541 1.412, 1.414, 1.414, 1.414, 1.414, 1.412, 1.410, 1.407, 1.404, 1.400, 1.395, 1.389, 1.383, 1.377,
1542 1.369, 1.361, 1.352, 1.342, 1.331, 1.320, 1.307, 1.294, 1.279, 1.263, 1.247, 1.230, 1.211, 1.192,
1543 1.173, 1.152, 1.131, 1.109, 1.086, 1.064, 1.041, 1.017, 0.994, 0.970, 0.947, 0.923, 0.899, 0.875,
1544 0.852, 0.828, 0.805, 0.782, 0.758, 0.736, 0.713, 0.691, 0.669, 0.647, 0.626, 0.605, 0.584, 0.563,
1545 0.544, 0.524, 0.505, 0.486, 0.468, 0.450, 0.432, 0.415, 0.398, 0.382, 0.366, 0.351, 0.336, 0.321,
1546 0.307, 0.294, 0.279, 0.267, 0.254, 0.242, 0.230, 0.219, 0.208, 0.197, 0.187, 0.177, 0.168, 0.159,
1547 0.150, 0.141, 0.133, 0.125, 0.118, 0.110, 0.104, 0.097, 0.090, 0.084, 0.078, 0.073, 0.068, 0.062,
1548 0.058, 0.053, 0.049, 0.044, 0.040, 0.037, 0.033, 0.030, 0.026, 0.023, 0.020, 0.018, 0.015, 0.013,
1549 0.010, 0.008, 0.006, 0.004, 0.002, 0.001, -0.001, -0.002, -0.003, -0.005, -0.006, -0.007, -0.008, -0.009,
1550 -0.010, -0.010, -0.011, -0.012, -0.012, -0.013, -0.013, -0.014, -0.014, -0.014, -0.015, -0.015, -0.015, -0.015,
1551 -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.014,
1552 -0.014, -0.014, -0.014, -0.014, -0.013, -0.013, -0.013, -0.013, -0.013, -0.012, -0.012, -0.012, -0.012, -0.012,
1553 -0.011, -0.011, -0.011, -0.011},
1554 {// charge = 345 fC
1555 0.000, -0.001, 0.012, 0.045, 0.098, 0.167, 0.247, 0.335, 0.427, 0.521, 0.615, 0.708, 0.797, 0.884,
1556 0.967, 1.045, 1.118, 1.186, 1.244, 1.293, 1.331, 1.360, 1.380, 1.395, 1.405, 1.412, 1.417, 1.420,
1557 1.422, 1.423, 1.424, 1.424, 1.423, 1.422, 1.420, 1.417, 1.415, 1.411, 1.407, 1.403, 1.398, 1.392,
1558 1.386, 1.379, 1.371, 1.363, 1.354, 1.344, 1.333, 1.322, 1.309, 1.296, 1.281, 1.266, 1.249, 1.232,
1559 1.214, 1.195, 1.175, 1.155, 1.133, 1.111, 1.089, 1.066, 1.043, 1.020, 0.997, 0.973, 0.949, 0.926,
1560 0.902, 0.878, 0.854, 0.831, 0.807, 0.784, 0.761, 0.738, 0.716, 0.693, 0.671, 0.649, 0.628, 0.607,
1561 0.586, 0.566, 0.546, 0.526, 0.507, 0.488, 0.470, 0.452, 0.434, 0.417, 0.400, 0.384, 0.368, 0.352,
1562 0.337, 0.323, 0.308, 0.294, 0.281, 0.268, 0.255, 0.243, 0.231, 0.220, 0.209, 0.198, 0.188, 0.178,
1563 0.169, 0.159, 0.150, 0.142, 0.134, 0.126, 0.118, 0.111, 0.104, 0.097, 0.091, 0.085, 0.079, 0.073,
1564 0.068, 0.063, 0.058, 0.053, 0.049, 0.044, 0.040, 0.037, 0.033, 0.030, 0.026, 0.023, 0.020, 0.018,
1565 0.015, 0.013, 0.010, 0.008, 0.006, 0.004, 0.002, 0.001, -0.001, -0.002, -0.004, -0.005, -0.006, -0.007,
1566 -0.008, -0.009, -0.010, -0.011, -0.011, -0.012, -0.012, -0.013, -0.013, -0.014, -0.014, -0.015, -0.015, -0.015,
1567 -0.015, -0.015, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.015, -0.015, -0.015, -0.015,
1568 -0.015, -0.015, -0.015, -0.014, -0.014, -0.014, -0.014, -0.014, -0.013, -0.013, -0.013, -0.013, -0.012, -0.012,
1569 -0.012, -0.012, -0.012, -0.011},
1570 {// charge = 360 fC
1571 0.000, -0.001, 0.012, 0.046, 0.100, 0.170, 0.252, 0.342, 0.437, 0.534, 0.632, 0.728, 0.822, 0.912,
1572 0.999, 1.081, 1.157, 1.226, 1.282, 1.327, 1.360, 1.383, 1.399, 1.410, 1.417, 1.421, 1.425, 1.427,
1573 1.428, 1.429, 1.430, 1.430, 1.429, 1.428, 1.427, 1.425, 1.423, 1.420, 1.417, 1.413, 1.409, 1.405,
1574 1.400, 1.394, 1.387, 1.380, 1.373, 1.364, 1.355, 1.346, 1.335, 1.323, 1.311, 1.297, 1.283, 1.268,
1575 1.251, 1.234, 1.216, 1.197, 1.178, 1.157, 1.136, 1.114, 1.092, 1.069, 1.046, 1.023, 1.000, 0.976,
1576 0.952, 0.929, 0.905, 0.881, 0.857, 0.834, 0.810, 0.787, 0.764, 0.741, 0.718, 0.696, 0.674, 0.652,
1577 0.631, 0.609, 0.589, 0.568, 0.548, 0.528, 0.509, 0.490, 0.472, 0.454, 0.436, 0.419, 0.402, 0.386,
1578 0.370, 0.355, 0.339, 0.324, 0.310, 0.296, 0.283, 0.269, 0.257, 0.245, 0.233, 0.221, 0.210, 0.200,
1579 0.189, 0.179, 0.170, 0.160, 0.151, 0.143, 0.134, 0.127, 0.119, 0.112, 0.105, 0.098, 0.091, 0.085,
1580 0.079, 0.074, 0.068, 0.063, 0.058, 0.053, 0.049, 0.045, 0.041, 0.037, 0.033, 0.030, 0.026, 0.023,
1581 0.020, 0.018, 0.015, 0.012, 0.010, 0.008, 0.006, 0.004, 0.002, 0.001, -0.001, -0.002, -0.004, -0.005,
1582 -0.006, -0.007, -0.008, -0.009, -0.010, -0.011, -0.012, -0.012, -0.013, -0.013, -0.014, -0.014, -0.015, -0.015,
1583 -0.015, -0.015, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016,
1584 -0.016, -0.015, -0.015, -0.015, -0.015, -0.015, -0.015, -0.014, -0.014, -0.014, -0.014, -0.013, -0.013, -0.013,
1585 -0.013, -0.013, -0.012, -0.012},
1586 {// charge = 375 fC
1587 0.000, -0.001, 0.012, 0.047, 0.102, 0.173, 0.256, 0.348, 0.446, 0.546, 0.647, 0.747, 0.845, 0.939,
1588 1.029, 1.114, 1.193, 1.261, 1.314, 1.353, 1.381, 1.399, 1.411, 1.418, 1.423, 1.427, 1.429, 1.431,
1589 1.432, 1.433, 1.434, 1.434, 1.433, 1.433, 1.432, 1.431, 1.429, 1.427, 1.425, 1.422, 1.419, 1.415,
1590 1.411, 1.406, 1.401, 1.395, 1.389, 1.382, 1.374, 1.366, 1.357, 1.347, 1.336, 1.325, 1.313, 1.299,
1591 1.285, 1.270, 1.253, 1.236, 1.218, 1.200, 1.180, 1.160, 1.139, 1.117, 1.095, 1.072, 1.049, 1.026,
1592 1.003, 0.979, 0.955, 0.932, 0.908, 0.884, 0.860, 0.837, 0.813, 0.790, 0.767, 0.744, 0.721, 0.699,
1593 0.677, 0.655, 0.633, 0.612, 0.591, 0.571, 0.551, 0.531, 0.511, 0.493, 0.474, 0.456, 0.438, 0.421,
1594 0.404, 0.387, 0.371, 0.356, 0.341, 0.326, 0.311, 0.298, 0.284, 0.271, 0.258, 0.246, 0.234, 0.222,
1595 0.211, 0.201, 0.190, 0.180, 0.170, 0.161, 0.152, 0.144, 0.135, 0.127, 0.120, 0.112, 0.105, 0.098,
1596 0.092, 0.086, 0.080, 0.074, 0.069, 0.063, 0.058, 0.054, 0.049, 0.045, 0.041, 0.037, 0.033, 0.030,
1597 0.026, 0.023, 0.020, 0.018, 0.015, 0.012, 0.010, 0.008, 0.006, 0.004, 0.002, 0.000, -0.001, -0.003,
1598 -0.004, -0.005, -0.007, -0.008, -0.009, -0.010, -0.010, -0.011, -0.012, -0.013, -0.013, -0.014, -0.014, -0.015,
1599 -0.015, -0.015, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016,
1600 -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.015, -0.015, -0.015, -0.015, -0.015, -0.014, -0.014, -0.014,
1601 -0.014, -0.013, -0.013, -0.013},
1602 {// charge = 390 fC
1603 0.000, -0.001, 0.012, 0.048, 0.103, 0.175, 0.259, 0.354, 0.454, 0.557, 0.661, 0.764, 0.865, 0.962,
1604 1.056, 1.145, 1.225, 1.290, 1.339, 1.373, 1.395, 1.409, 1.418, 1.424, 1.427, 1.430, 1.432, 1.434,
1605 1.435, 1.436, 1.436, 1.436, 1.436, 1.436, 1.436, 1.435, 1.434, 1.432, 1.430, 1.428, 1.426, 1.423,
1606 1.420, 1.416, 1.412, 1.407, 1.402, 1.396, 1.390, 1.383, 1.375, 1.367, 1.358, 1.348, 1.338, 1.327,
1607 1.314, 1.301, 1.287, 1.272, 1.256, 1.239, 1.221, 1.202, 1.183, 1.163, 1.142, 1.120, 1.098, 1.075,
1608 1.052, 1.029, 1.006, 0.982, 0.958, 0.935, 0.911, 0.887, 0.863, 0.840, 0.816, 0.793, 0.770, 0.747,
1609 0.724, 0.702, 0.679, 0.658, 0.636, 0.615, 0.594, 0.573, 0.553, 0.533, 0.514, 0.495, 0.476, 0.458,
1610 0.440, 0.423, 0.406, 0.389, 0.374, 0.357, 0.342, 0.327, 0.313, 0.299, 0.285, 0.272, 0.260, 0.247,
1611 0.235, 0.224, 0.212, 0.202, 0.191, 0.181, 0.171, 0.162, 0.153, 0.144, 0.136, 0.128, 0.120, 0.113,
1612 0.106, 0.099, 0.092, 0.086, 0.080, 0.074, 0.069, 0.064, 0.059, 0.054, 0.049, 0.045, 0.041, 0.037,
1613 0.033, 0.030, 0.026, 0.023, 0.020, 0.017, 0.015, 0.012, 0.010, 0.008, 0.006, 0.004, 0.002, 0.000,
1614 -0.001, -0.003, -0.004, -0.006, -0.007, -0.008, -0.009, -0.010, -0.011, -0.012, -0.012, -0.013, -0.013, -0.014,
1615 -0.015, -0.015, -0.015, -0.016, -0.016, -0.016, -0.016, -0.016, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017,
1616 -0.017, -0.017, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016, -0.015, -0.015, -0.015, -0.015,
1617 -0.015, -0.014, -0.014, -0.014},
1618 {// charge = 405 fC
1619 0.000, -0.001, 0.013, 0.049, 0.104, 0.177, 0.262, 0.358, 0.460, 0.566, 0.672, 0.778, 0.882, 0.983,
1620 1.080, 1.172, 1.251, 1.313, 1.358, 1.387, 1.405, 1.416, 1.422, 1.427, 1.430, 1.432, 1.434, 1.435,
1621 1.436, 1.437, 1.438, 1.438, 1.438, 1.438, 1.438, 1.438, 1.437, 1.436, 1.435, 1.433, 1.431, 1.429,
1622 1.427, 1.424, 1.421, 1.417, 1.413, 1.408, 1.403, 1.397, 1.391, 1.384, 1.376, 1.368, 1.359, 1.350,
1623 1.339, 1.328, 1.316, 1.303, 1.289, 1.274, 1.258, 1.241, 1.223, 1.205, 1.186, 1.165, 1.145, 1.123,
1624 1.101, 1.078, 1.055, 1.032, 1.009, 0.985, 0.962, 0.938, 0.914, 0.890, 0.867, 0.843, 0.819, 0.796,
1625 0.773, 0.750, 0.727, 0.704, 0.682, 0.660, 0.639, 0.617, 0.596, 0.576, 0.555, 0.535, 0.516, 0.497,
1626 0.478, 0.460, 0.442, 0.425, 0.408, 0.391, 0.375, 0.359, 0.344, 0.329, 0.314, 0.300, 0.287, 0.274,
1627 0.261, 0.248, 0.236, 0.225, 0.213, 0.202, 0.192, 0.182, 0.172, 0.163, 0.154, 0.145, 0.137, 0.128,
1628 0.121, 0.113, 0.106, 0.099, 0.093, 0.086, 0.080, 0.075, 0.069, 0.064, 0.059, 0.054, 0.049, 0.045,
1629 0.041, 0.037, 0.033, 0.030, 0.026, 0.023, 0.020, 0.017, 0.015, 0.012, 0.010, 0.007, 0.005, 0.003,
1630 0.002, -0.000, -0.002, -0.003, -0.005, -0.006, -0.007, -0.008, -0.009, -0.010, -0.011, -0.012, -0.013, -0.013,
1631 -0.014, -0.014, -0.015, -0.015, -0.016, -0.016, -0.016, -0.016, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017,
1632 -0.017, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017, -0.016, -0.016, -0.016, -0.016, -0.016, -0.016,
1633 -0.015, -0.015, -0.015, -0.015},
1634 {// charge = 420 fC
1635 0.000, -0.001, 0.013, 0.050, 0.105, 0.178, 0.265, 0.362, 0.466, 0.573, 0.682, 0.790, 0.897, 1.000,
1636 1.099, 1.193, 1.271, 1.330, 1.370, 1.396, 1.410, 1.419, 1.424, 1.428, 1.431, 1.433, 1.435, 1.436,
1637 1.437, 1.438, 1.439, 1.439, 1.440, 1.440, 1.440, 1.440, 1.439, 1.439, 1.438, 1.437, 1.436, 1.434,
1638 1.432, 1.430, 1.428, 1.425, 1.421, 1.418, 1.414, 1.409, 1.404, 1.398, 1.392, 1.385, 1.378, 1.370,
1639 1.361, 1.351, 1.341, 1.330, 1.318, 1.305, 1.291, 1.276, 1.260, 1.244, 1.226, 1.208, 1.188, 1.168,
1640 1.148, 1.126, 1.104, 1.081, 1.059, 1.035, 1.012, 0.989, 0.965, 0.941, 0.917, 0.893, 0.870, 0.846,
1641 0.822, 0.799, 0.776, 0.753, 0.730, 0.707, 0.685, 0.663, 0.641, 0.620, 0.599, 0.578, 0.558, 0.538,
1642 0.518, 0.499, 0.480, 0.462, 0.444, 0.426, 0.410, 0.392, 0.376, 0.361, 0.345, 0.330, 0.316, 0.302,
1643 0.288, 0.275, 0.262, 0.249, 0.237, 0.225, 0.214, 0.203, 0.193, 0.183, 0.173, 0.163, 0.154, 0.145,
1644 0.137, 0.129, 0.121, 0.113, 0.106, 0.099, 0.093, 0.086, 0.080, 0.075, 0.069, 0.064, 0.059, 0.054,
1645 0.049, 0.045, 0.041, 0.037, 0.033, 0.030, 0.026, 0.023, 0.020, 0.017, 0.014, 0.012, 0.009, 0.007,
1646 0.005, 0.003, 0.001, -0.001, -0.002, -0.004, -0.005, -0.006, -0.008, -0.009, -0.010, -0.011, -0.012, -0.012,
1647 -0.013, -0.014, -0.014, -0.015, -0.015, -0.016, -0.016, -0.016, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017,
1648 -0.018, -0.018, -0.018, -0.018, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017, -0.017, -0.016, -0.016,
1649 -0.016, -0.016, -0.016, -0.016},
1650 {// charge = 435 fC
1651 0.000, -0.001, 0.013, 0.050, 0.107, 0.180, 0.267, 0.365, 0.470, 0.579, 0.690, 0.800, 0.908, 1.013,
1652 1.114, 1.208, 1.285, 1.341, 1.378, 1.400, 1.413, 1.421, 1.425, 1.429, 1.431, 1.433, 1.435, 1.436,
1653 1.438, 1.439, 1.439, 1.440, 1.441, 1.441, 1.441, 1.441, 1.441, 1.441, 1.441, 1.440, 1.439, 1.438,
1654 1.437, 1.435, 1.433, 1.431, 1.429, 1.426, 1.422, 1.419, 1.415, 1.410, 1.405, 1.399, 1.393, 1.387,
1655 1.379, 1.371, 1.363, 1.353, 1.343, 1.332, 1.320, 1.308, 1.294, 1.279, 1.263, 1.247, 1.229, 1.211,
1656 1.192, 1.172, 1.151, 1.130, 1.108, 1.085, 1.062, 1.039, 1.016, 0.992, 0.968, 0.945, 0.921, 0.897,
1657 0.873, 0.849, 0.826, 0.802, 0.779, 0.756, 0.733, 0.710, 0.688, 0.665, 0.644, 0.622, 0.601, 0.580,
1658 0.560, 0.540, 0.520, 0.501, 0.482, 0.464, 0.446, 0.428, 0.411, 0.394, 0.378, 0.362, 0.346, 0.331,
1659 0.317, 0.303, 0.289, 0.275, 0.262, 0.250, 0.238, 0.226, 0.215, 0.204, 0.193, 0.183, 0.173, 0.163,
1660 0.154, 0.145, 0.137, 0.129, 0.121, 0.113, 0.106, 0.099, 0.092, 0.086, 0.080, 0.074, 0.069, 0.063,
1661 0.058, 0.053, 0.049, 0.044, 0.040, 0.036, 0.033, 0.029, 0.026, 0.022, 0.019, 0.016, 0.014, 0.011,
1662 0.009, 0.007, 0.004, 0.002, 0.001, -0.001, -0.003, -0.004, -0.006, -0.007, -0.008, -0.009, -0.010, -0.011,
1663 -0.012, -0.013, -0.014, -0.014, -0.015, -0.016, -0.016, -0.016, -0.017, -0.017, -0.017, -0.018, -0.018, -0.018,
1664 -0.018, -0.018, -0.018, -0.018, -0.018, -0.018, -0.018, -0.018, -0.018, -0.018, -0.018, -0.018, -0.017, -0.017,
1665 -0.017, -0.017, -0.017, -0.016},
1666 {// charge = 450 fC
1667 0.000, -0.001, 0.014, 0.051, 0.108, 0.181, 0.269, 0.368, 0.474, 0.584, 0.696, 0.807, 0.916, 1.022,
1668 1.124, 1.217, 1.292, 1.347, 1.382, 1.398, 1.410, 1.418, 1.424, 1.430, 1.436, 1.439, 1.441, 1.442,
1669 1.443, 1.442, 1.441, 1.439, 1.435, 1.431, 1.425, 1.419, 1.412, 1.404, 1.395, 1.386, 1.376, 1.366,
1670 1.355, 1.343, 1.330, 1.316, 1.301, 1.285, 1.268, 1.250, 1.232, 1.212, 1.192, 1.171, 1.149, 1.127,
1671 1.104, 1.080, 1.056, 1.032, 1.008, 0.983, 0.959, 0.934, 0.910, 0.885, 0.860, 0.836, 0.812, 0.788,
1672 0.764, 0.740, 0.717, 0.694, 0.671, 0.649, 0.627, 0.605, 0.583, 0.562, 0.542, 0.522, 0.502, 0.483,
1673 0.464, 0.445, 0.427, 0.410, 0.393, 0.376, 0.360, 0.344, 0.329, 0.314, 0.299, 0.285, 0.271, 0.259,
1674 0.245, 0.233, 0.221, 0.210, 0.199, 0.188, 0.178, 0.168, 0.158, 0.149, 0.140, 0.131, 0.123, 0.115,
1675 0.108, 0.100, 0.093, 0.087, 0.081, 0.074, 0.069, 0.063, 0.058, 0.053, 0.048, 0.044, 0.039, 0.035,
1676 0.031, 0.027, 0.024, 0.021, 0.018, 0.015, 0.012, 0.009, 0.007, 0.004, 0.002, 0.000, -0.002, -0.003,
1677 -0.005, -0.006, -0.008, -0.009, -0.010, -0.011, -0.012, -0.013, -0.014, -0.015, -0.016, -0.016, -0.017, -0.017,
1678 -0.018, -0.018, -0.018, -0.019, -0.019, -0.019, -0.019, -0.019, -0.019, -0.019, -0.019, -0.019, -0.019, -0.019,
1679 -0.019, -0.019, -0.019, -0.019, -0.018, -0.018, -0.018, -0.018, -0.017, -0.017, -0.017, -0.017, -0.016, -0.016,
1680 -0.016, -0.015, -0.015, -0.015, -0.015, -0.014, -0.014, -0.014, -0.013, -0.013, -0.013, -0.012, -0.012, -0.012,
1681 -0.012, -0.011, -0.011, -0.011},
1682 {// charge = 465 fC
1683 0.000, -0.001, 0.014, 0.052, 0.109, 0.183, 0.271, 0.370, 0.477, 0.588, 0.701, 0.813, 0.922, 1.028,
1684 1.130, 1.223, 1.305, 1.361, 1.393, 1.409, 1.421, 1.430, 1.435, 1.438, 1.439, 1.439, 1.437, 1.433,
1685 1.426, 1.417, 1.405, 1.392, 1.378, 1.362, 1.345, 1.327, 1.308, 1.287, 1.266, 1.244, 1.220, 1.196,
1686 1.172, 1.146, 1.120, 1.094, 1.067, 1.040, 1.013, 0.986, 0.959, 0.932, 0.905, 0.878, 0.851, 0.825,
1687 0.799, 0.773, 0.747, 0.722, 0.697, 0.673, 0.649, 0.625, 0.602, 0.579, 0.557, 0.535, 0.514, 0.493,
1688 0.472, 0.453, 0.433, 0.414, 0.396, 0.378, 0.361, 0.344, 0.328, 0.312, 0.297, 0.282, 0.267, 0.253,
1689 0.240, 0.227, 0.215, 0.203, 0.191, 0.180, 0.169, 0.159, 0.149, 0.139, 0.130, 0.121, 0.113, 0.105,
1690 0.097, 0.090, 0.083, 0.076, 0.069, 0.063, 0.058, 0.052, 0.047, 0.042, 0.037, 0.033, 0.028, 0.024,
1691 0.021, 0.017, 0.014, 0.011, 0.007, 0.005, 0.002, -0.000, -0.003, -0.005, -0.007, -0.009, -0.010, -0.012,
1692 -0.013, -0.015, -0.016, -0.017, -0.018, -0.019, -0.020, -0.020, -0.021, -0.021, -0.022, -0.022, -0.023, -0.023,
1693 -0.023, -0.024, -0.024, -0.024, -0.024, -0.024, -0.024, -0.024, -0.024, -0.023, -0.023, -0.023, -0.023, -0.022,
1694 -0.022, -0.022, -0.022, -0.021, -0.021, -0.021, -0.020, -0.020, -0.020, -0.019, -0.019, -0.018, -0.018, -0.018,
1695 -0.017, -0.017, -0.016, -0.016, -0.016, -0.015, -0.015, -0.015, -0.014, -0.014, -0.013, -0.013, -0.013, -0.012,
1696 -0.012, -0.012, -0.011, -0.011, -0.011, -0.010, -0.010, -0.010, -0.010, -0.009, -0.009, -0.009, -0.008, -0.008,
1697 -0.008, -0.008, -0.007, -0.007},
1698 {// charge = 480 fC
1699 0.000, -0.001, 0.014, 0.052, 0.109, 0.184, 0.272, 0.372, 0.480, 0.591, 0.705, 0.817, 0.927, 1.032,
1700 1.140, 1.246, 1.324, 1.373, 1.399, 1.416, 1.426, 1.433, 1.436, 1.438, 1.437, 1.433, 1.426, 1.415,
1701 1.401, 1.385, 1.366, 1.346, 1.324, 1.300, 1.275, 1.249, 1.222, 1.194, 1.166, 1.137, 1.107, 1.077,
1702 1.047, 1.017, 0.987, 0.957, 0.927, 0.897, 0.868, 0.839, 0.811, 0.782, 0.755, 0.727, 0.700, 0.674,
1703 0.648, 0.622, 0.598, 0.573, 0.549, 0.526, 0.503, 0.481, 0.460, 0.439, 0.419, 0.399, 0.380, 0.361,
1704 0.343, 0.325, 0.308, 0.292, 0.276, 0.261, 0.246, 0.232, 0.218, 0.205, 0.193, 0.180, 0.169, 0.157,
1705 0.147, 0.136, 0.127, 0.117, 0.108, 0.100, 0.091, 0.083, 0.076, 0.068, 0.062, 0.055, 0.049, 0.043,
1706 0.038, 0.033, 0.028, 0.023, 0.019, 0.015, 0.011, 0.007, 0.004, 0.001, -0.002, -0.005, -0.008, -0.010,
1707 -0.012, -0.015, -0.016, -0.018, -0.020, -0.021, -0.023, -0.024, -0.025, -0.026, -0.027, -0.028, -0.028, -0.029,
1708 -0.029, -0.030, -0.030, -0.030, -0.030, -0.030, -0.031, -0.031, -0.030, -0.030, -0.030, -0.030, -0.030, -0.030,
1709 -0.029, -0.029, -0.029, -0.028, -0.028, -0.028, -0.027, -0.027, -0.026, -0.026, -0.025, -0.025, -0.024, -0.024,
1710 -0.023, -0.023, -0.022, -0.022, -0.021, -0.021, -0.020, -0.020, -0.019, -0.019, -0.018, -0.018, -0.017, -0.017,
1711 -0.016, -0.016, -0.015, -0.015, -0.014, -0.014, -0.013, -0.013, -0.013, -0.012, -0.012, -0.012, -0.011, -0.011,
1712 -0.010, -0.010, -0.010, -0.009, -0.009, -0.009, -0.008, -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007,
1713 -0.006, -0.006, -0.006, -0.006},
1714 {// charge = 495 fC
1715 0.000, -0.001, 0.014, 0.053, 0.110, 0.185, 0.274, 0.374, 0.482, 0.594, 0.708, 0.820, 0.929, 1.041,
1716 1.155, 1.259, 1.333, 1.378, 1.403, 1.419, 1.428, 1.434, 1.436, 1.436, 1.433, 1.426, 1.414, 1.398,
1717 1.379, 1.357, 1.333, 1.306, 1.278, 1.248, 1.218, 1.187, 1.155, 1.122, 1.089, 1.056, 1.023, 0.990,
1718 0.958, 0.925, 0.893, 0.861, 0.830, 0.799, 0.769, 0.739, 0.710, 0.682, 0.653, 0.626, 0.599, 0.573,
1719 0.547, 0.522, 0.498, 0.475, 0.452, 0.429, 0.408, 0.387, 0.367, 0.347, 0.328, 0.310, 0.292, 0.275,
1720 0.259, 0.243, 0.228, 0.213, 0.199, 0.185, 0.172, 0.160, 0.148, 0.137, 0.126, 0.115, 0.105, 0.096,
1721 0.087, 0.078, 0.070, 0.063, 0.054, 0.048, 0.041, 0.035, 0.029, 0.024, 0.018, 0.013, 0.009, 0.004,
1722 0.000, -0.003, -0.007, -0.010, -0.013, -0.016, -0.019, -0.021, -0.024, -0.026, -0.028, -0.029, -0.031, -0.032,
1723 -0.033, -0.035, -0.035, -0.036, -0.037, -0.038, -0.038, -0.039, -0.039, -0.039, -0.039, -0.039, -0.039, -0.039,
1724 -0.039, -0.039, -0.039, -0.039, -0.038, -0.038, -0.037, -0.037, -0.036, -0.036, -0.035, -0.035, -0.034, -0.034,
1725 -0.033, -0.033, -0.032, -0.031, -0.030, -0.030, -0.029, -0.029, -0.028, -0.027, -0.026, -0.026, -0.025, -0.025,
1726 -0.024, -0.023, -0.022, -0.022, -0.021, -0.021, -0.020, -0.019, -0.019, -0.018, -0.018, -0.017, -0.017, -0.016,
1727 -0.016, -0.015, -0.015, -0.014, -0.013, -0.013, -0.013, -0.012, -0.012, -0.011, -0.011, -0.011, -0.010, -0.010,
1728 -0.009, -0.009, -0.009, -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.006,
1729 -0.006, -0.005, -0.005, -0.005},
1730 {// charge = 510 fC
1731 0.000, -0.001, 0.015, 0.053, 0.111, 0.186, 0.275, 0.376, 0.484, 0.597, 0.710, 0.823, 0.934, 1.049,
1732 1.165, 1.267, 1.338, 1.381, 1.405, 1.420, 1.429, 1.434, 1.436, 1.435, 1.430, 1.420, 1.405, 1.385,
1733 1.362, 1.335, 1.307, 1.276, 1.243, 1.210, 1.176, 1.141, 1.106, 1.070, 1.034, 0.999, 0.964, 0.929,
1734 0.895, 0.861, 0.828, 0.795, 0.763, 0.731, 0.700, 0.670, 0.640, 0.611, 0.583, 0.556, 0.529, 0.503,
1735 0.478, 0.453, 0.430, 0.407, 0.384, 0.363, 0.342, 0.322, 0.303, 0.284, 0.266, 0.249, 0.232, 0.216,
1736 0.201, 0.186, 0.172, 0.159, 0.146, 0.133, 0.122, 0.110, 0.100, 0.089, 0.080, 0.071, 0.062, 0.053,
1737 0.046, 0.038, 0.031, 0.025, 0.018, 0.012, 0.007, 0.002, -0.003, -0.007, -0.012, -0.016, -0.019, -0.022,
1738 -0.025, -0.028, -0.031, -0.033, -0.036, -0.038, -0.039, -0.041, -0.042, -0.044, -0.045, -0.046, -0.047, -0.047,
1739 -0.048, -0.048, -0.049, -0.049, -0.049, -0.049, -0.049, -0.049, -0.049, -0.048, -0.048, -0.048, -0.047, -0.047,
1740 -0.046, -0.046, -0.045, -0.044, -0.044, -0.043, -0.042, -0.042, -0.041, -0.040, -0.039, -0.038, -0.038, -0.037,
1741 -0.036, -0.035, -0.034, -0.033, -0.032, -0.032, -0.031, -0.030, -0.029, -0.028, -0.028, -0.027, -0.026, -0.025,
1742 -0.024, -0.024, -0.023, -0.022, -0.021, -0.021, -0.020, -0.019, -0.019, -0.018, -0.017, -0.017, -0.016, -0.016,
1743 -0.015, -0.015, -0.014, -0.014, -0.013, -0.013, -0.012, -0.012, -0.011, -0.011, -0.010, -0.010, -0.010, -0.009,
1744 -0.009, -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005,
1745 -0.005, -0.005, -0.005, -0.004},
1746 {// charge = 525 fC
1747 0.000, -0.001, 0.015, 0.054, 0.112, 0.187, 0.277, 0.377, 0.486, 0.599, 0.713, 0.825, 0.939, 1.056,
1748 1.172, 1.272, 1.342, 1.383, 1.406, 1.421, 1.430, 1.435, 1.436, 1.435, 1.428, 1.416, 1.398, 1.376,
1749 1.350, 1.321, 1.289, 1.255, 1.220, 1.184, 1.147, 1.110, 1.072, 1.035, 0.997, 0.960, 0.924, 0.888,
1750 0.852, 0.817, 0.783, 0.749, 0.717, 0.684, 0.653, 0.622, 0.593, 0.564, 0.535, 0.508, 0.481, 0.455,
1751 0.430, 0.406, 0.383, 0.360, 0.338, 0.317, 0.297, 0.277, 0.259, 0.241, 0.223, 0.207, 0.191, 0.175,
1752 0.161, 0.147, 0.134, 0.121, 0.109, 0.097, 0.086, 0.076, 0.066, 0.057, 0.048, 0.039, 0.032, 0.024,
1753 0.017, 0.011, 0.004, -0.001, -0.007, -0.012, -0.017, -0.021, -0.025, -0.029, -0.032, -0.036, -0.039, -0.041,
1754 -0.044, -0.046, -0.048, -0.050, -0.051, -0.053, -0.054, -0.055, -0.056, -0.056, -0.057, -0.058, -0.058, -0.058,
1755 -0.058, -0.058, -0.058, -0.058, -0.058, -0.057, -0.057, -0.056, -0.056, -0.055, -0.054, -0.054, -0.053, -0.052,
1756 -0.051, -0.050, -0.050, -0.049, -0.048, -0.047, -0.046, -0.045, -0.044, -0.043, -0.042, -0.041, -0.040, -0.039,
1757 -0.038, -0.037, -0.036, -0.035, -0.034, -0.033, -0.032, -0.031, -0.030, -0.029, -0.028, -0.028, -0.027, -0.026,
1758 -0.025, -0.024, -0.023, -0.022, -0.022, -0.021, -0.020, -0.020, -0.019, -0.018, -0.017, -0.017, -0.016, -0.016,
1759 -0.015, -0.015, -0.014, -0.013, -0.013, -0.012, -0.012, -0.011, -0.011, -0.011, -0.010, -0.010, -0.009, -0.009,
1760 -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005,
1761 -0.005, -0.005, -0.004, -0.004},
1762 {// charge = 540 fC
1763 0.000, -0.001, 0.015, 0.054, 0.113, 0.188, 0.278, 0.379, 0.487, 0.601, 0.715, 0.827, 0.944, 1.061,
1764 1.177, 1.276, 1.344, 1.384, 1.407, 1.422, 1.430, 1.435, 1.437, 1.435, 1.428, 1.415, 1.395, 1.371,
1765 1.343, 1.312, 1.278, 1.242, 1.205, 1.167, 1.129, 1.090, 1.051, 1.012, 0.974, 0.935, 0.898, 0.861,
1766 0.825, 0.789, 0.754, 0.720, 0.686, 0.654, 0.622, 0.591, 0.561, 0.532, 0.503, 0.476, 0.449, 0.423,
1767 0.399, 0.374, 0.351, 0.329, 0.307, 0.286, 0.266, 0.247, 0.229, 0.211, 0.194, 0.178, 0.163, 0.148,
1768 0.134, 0.120, 0.108, 0.095, 0.084, 0.073, 0.062, 0.052, 0.043, 0.034, 0.026, 0.018, 0.011, 0.004,
1769 -0.003, -0.009, -0.014, -0.019, -0.025, -0.029, -0.033, -0.037, -0.041, -0.044, -0.047, -0.050, -0.052, -0.054,
1770 -0.056, -0.058, -0.060, -0.061, -0.062, -0.063, -0.064, -0.065, -0.065, -0.066, -0.066, -0.066, -0.066, -0.066,
1771 -0.066, -0.065, -0.065, -0.064, -0.064, -0.063, -0.062, -0.062, -0.061, -0.060, -0.059, -0.058, -0.057, -0.056,
1772 -0.055, -0.054, -0.053, -0.052, -0.051, -0.050, -0.049, -0.048, -0.046, -0.045, -0.044, -0.043, -0.042, -0.041,
1773 -0.040, -0.039, -0.037, -0.036, -0.035, -0.034, -0.033, -0.032, -0.031, -0.030, -0.029, -0.028, -0.027, -0.026,
1774 -0.025, -0.025, -0.024, -0.023, -0.022, -0.021, -0.021, -0.020, -0.019, -0.018, -0.018, -0.017, -0.016, -0.016,
1775 -0.015, -0.015, -0.014, -0.013, -0.013, -0.012, -0.012, -0.011, -0.011, -0.010, -0.010, -0.010, -0.009, -0.009,
1776 -0.008, -0.008, -0.008, -0.007, -0.007, -0.007, -0.007, -0.006, -0.006, -0.006, -0.005, -0.005, -0.005, -0.005,
1777 -0.005, -0.004, -0.004, -0.004}};
ClassImp(CbmConverterManager)
#define VERBOSE
#define DRAW
#define FASP_WINDOW
Definition CbmTrdFASP.h:25
#define SHAPER_LUT
Definition CbmTrdFASP.h:26
#define NGRAPH
Definition CbmTrdFASP.h:27
friend fscal max(fscal x, fscal y)
float Float_t
int Int_t
bool Bool_t
Generates beam ions for transport simulation.
virtual std::string ToString() const
Return string representation of the object.
Definition CbmMatch.cxx:25
int32_t GetTriggerType() const
Channel trigger type. SPADIC specific see CbmTrdTriggerType.
Definition CbmTrdDigi.h:160
void SetFlag(const int32_t iflag, bool set=true)
Generic flag status setter.
std::string ToString() const
String representation of a TRD digi. Account for digi type and specific information.
uint64_t GetTimeDAQ() const
Getter for global DAQ time [clk]. Differs for each ASIC. [In FASP case DAQ time is already stored in ...
Definition CbmTrdDigi.h:158
int32_t GetAddressChannel() const
Getter read-out id.
bool IsFlagged(const int32_t iflag) const
Query flag status (generic)
static float Clk(eCbmTrdAsicType ty)
DAQ clock accessor for each ASIC.
Definition CbmTrdDigi.h:110
double GetTime() const
Getter for physical time [ns]. Accounts for clock representation of each ASIC. In SPADIC case physica...
Definition CbmTrdDigi.h:155
double GetCharge() const
Common purpose charge getter.
FASP channel simulator.
Definition CbmTrdFASP.h:30
void GetShaperSignal(Double_t charge)
Retrive linear interpolation of CADENCE for signal.
int fAsicId[2]
identifier of FASP(s) in module
Definition CbmTrdFASP.h:127
std::vector< Float_t > fShaperNext
next channel shaper analog
Definition CbmTrdFASP.h:112
TLine * fGthr
graph representation of various thresholds
Definition CbmTrdFASP.h:135
static const Int_t fgkNDB
DB shaper size.
Definition CbmTrdFASP.h:139
Int_t fNphys[2]
number of physical digi in the current [0] and next [1] shaper
Definition CbmTrdFASP.h:105
Int_t ProcessShaper(Char_t typ='T')
Calculate output FASP signal and CS timming for the signal array stored in fShaper.
UInt_t fProcTime
time window [ns] for actual digi processing (excluded fgkBufferKeep)
Definition CbmTrdFASP.h:103
static Bool_t fgNeighbour
Neighbour enable flag.
Definition CbmTrdFASP.h:149
void SetProcTime(ULong64_t t=0)
Set limit in time for processing digis.
Int_t fTimeDY
Time decay from FT [5*ns].
Definition CbmTrdFASP.h:122
Int_t fNraw
number of raw digi for the tilt channel
Definition CbmTrdFASP.h:106
virtual void InitChannel(int id, const CbmTrdParFaspChannel *par, int asicId=-1, int chId=-1)
[Re]Initialize one of the two FASP channels
static const Int_t fgkBufferKeep
length of buffer time in 5ns which is kept between cycles
Definition CbmTrdFASP.h:156
TGraph * fGraphShp[NGRAPH]
graph representations of FASP shaper
Definition CbmTrdFASP.h:133
virtual Bool_t Go(ULong64_t time)
Check if there is enough time elapsed from fStartTime to run simulator.
static const Int_t fgkNclkFT
length of flat top in FASP clocks
Definition CbmTrdFASP.h:147
Float_t fFT
Flat Top value [V].
Definition CbmTrdFASP.h:123
static Float_t fgShaperThr
shaper threshold [V]
Definition CbmTrdFASP.h:151
virtual void Clear(Option_t *opt="")
Finalize currently stored data.
void ScanDigi(std::vector< std::pair< CbmTrdDigi *, CbmMatch * > > *digi)
Read digi array for one pair T/R defined by the pad column.
static const Float_t fgkShaperLUT[SHAPER_LUT]
shaper LUT
Definition CbmTrdFASP.h:907
Double_t MakeOut(Int_t t)
Make convolution of shaper1 superposition and theoretic shaping model (see fgkShaperPar)
TGraph * fGraphPhys[NGRAPH]
graph representations of physics digi
Definition CbmTrdFASP.h:134
virtual void Print(Option_t *opt="") const
Print-out FASP analog/digital response to currently stored data.
virtual void Draw(Option_t *opt="")
Graphical representation of FASP analog/digital response to currently stored data.
TCanvas * fMonitor
monitor canvas when drawing
Definition CbmTrdFASP.h:136
int fChId[2]
FASP channels being processed.
Definition CbmTrdFASP.h:128
void WriteDigi()
Write processed digi to output array.
std::vector< bool > fHitThPrev
previous channel hit threshold
Definition CbmTrdFASP.h:110
static const Float_t fgkShaperPar[4]
shaper parameters
Definition CbmTrdFASP.h:906
void ScanDigiNE(std::vector< std::pair< CbmTrdDigi *, CbmMatch * > > *digi)
Read digi array for neighbour trigger processing.
CbmTrdFASP(UInt_t uslice=1000)
Constructor of FASP simulator.
std::vector< std::pair< CbmTrdDigi *, CbmMatch * > > * fDigi
link to digi vector to be transformed
Definition CbmTrdFASP.h:107
static Float_t fgBaseline
FASP baseline [V].
Definition CbmTrdFASP.h:152
int fPad
current pad as defined by CbmTrdModuleAbstract::GetPadAddress()
Definition CbmTrdFASP.h:104
virtual void PhysToRaw(std::vector< std::pair< CbmTrdDigi *, CbmMatch * > > *digi)
Convert physics information in digi to the raw format.
Int_t fTimeLG
Linear gate time length [5*ns].
Definition CbmTrdFASP.h:120
static Float_t fgNeighbourThr
neighbour threshold [V] for fgNeighbour=kTRUE
Definition CbmTrdFASP.h:150
virtual ~CbmTrdFASP()
Float_t fSignal[FASP_WINDOW]
temporary array to store shaper analog signal for current charge interpolation
Definition CbmTrdFASP.h:115
static Int_t fgNclkLG
length of linear-gate command in FASP clocks
Definition CbmTrdFASP.h:148
int AddGraph(char typ='T')
std::vector< Float_t > fOut
analog output for the current channel
Definition CbmTrdFASP.h:129
ULong64_t fStartTime
time offset [ns] for the current simulation
Definition CbmTrdFASP.h:102
std::vector< Float_t > fShaper
current channel shaper analog
Definition CbmTrdFASP.h:111
TGraph * fGraph[NGRAPH]
graph representations of analog FASP response
Definition CbmTrdFASP.h:132
std::map< int, std::array< int, NFASPCH > > fGraphMap
map of ASIC_id and (ch_id, output_id of FASP signals graphs) pairs
Definition CbmTrdFASP.h:131
static const Float_t fgkShaper[fgkNDB][FASP_WINDOW]
DB shaper signals for each input charge discretization.
Definition CbmTrdFASP.h:929
int fNgraph
No of graphs generated.
Definition CbmTrdFASP.h:126
const CbmTrdParFaspChannel * fPar[2]
current FASP ASIC parametrization
Definition CbmTrdFASP.h:119
static Float_t fgOutGain
FASP -> ADC gain [V/4095 ADC].
Definition CbmTrdFASP.h:153
static const Float_t fgkCharge[fgkNDB]
DB input charge discretization.
Definition CbmTrdFASP.h:925
Int_t fTimeFT
Chip Select time legth [5*ns].
Definition CbmTrdFASP.h:121
static const Float_t fgkDecayLUT[SHAPER_LUT]
forced discharged of FASP LUT
Definition CbmTrdFASP.h:916
std::vector< std::tuple< UInt_t, UInt_t, UInt_t, Bool_t > > fDigiProc
proccessed info wrt fStartTime <hit_time[ns], CS_time[ns], OUT[ADC], trigger>
Definition CbmTrdFASP.h:113
Definition of FASP channel calibration container.
Hash for CbmL1LinkKey.
#define NFASPCH