]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSPreprocessor.cxx
move minuit initialization to the unfolding method since it leaks and it is not used...
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPreprocessor.cxx
CommitLineData
1ab07e55 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/* $Id$ */
17
18///////////////////////////////////////////////////////////////////////////////
19// PHOS Preprocessor class. It runs by Shuttle at the end of the run,
20// calculates calibration coefficients and dead/bad channels
21// to be posted in OCDB
22//
23// Author: Boris Polichtchouk, 4 October 2006
24///////////////////////////////////////////////////////////////////////////////
25
26#include "AliPHOSPreprocessor.h"
27#include "AliLog.h"
28#include "AliCDBMetaData.h"
30a0d5a2 29#include "AliCDBEntry.h"
1ab07e55 30#include "AliPHOSEmcCalibData.h"
31#include "TFile.h"
32#include "TH1.h"
450717bc 33#include "TF1.h"
34#include "TH2.h"
1ab07e55 35#include "TMap.h"
36#include "TRandom.h"
f97eb136 37#include "TKey.h"
516e7ab3 38#include "TList.h"
39#include "TObjString.h"
30a0d5a2 40#include "TObjArray.h"
450717bc 41#include "TObject.h"
42#include "TString.h"
43#include "TMath.h"
306e5917 44#include "TAxis.h"
3c0265c1 45#include "AliPHOSEmcBadChannelsMap.h"
1ab07e55 46
47ClassImp(AliPHOSPreprocessor)
48
450717bc 49 Double_t rgaus(Double_t *x, Double_t *par)
50{
51 Double_t gaus = par[0] * TMath::Exp( -(x[0]-par[1])*(x[0]-par[1]) /
52 (2*par[2]*par[2]) );
53 return gaus;
54}
55
1ab07e55 56//_______________________________________________________________________________________
57AliPHOSPreprocessor::AliPHOSPreprocessor() :
81587b1e 58AliPreprocessor("PHS",0)
1ab07e55 59{
60 //default constructor
61}
62
63//_______________________________________________________________________________________
81587b1e 64AliPHOSPreprocessor::AliPHOSPreprocessor(AliShuttleInterface* shuttle):
65AliPreprocessor("PHS",shuttle)
1ab07e55 66{
67 // Constructor
ce2c6758 68
69 AddRunType("PHYSICS");
71eef542 70 AddRunType("LED");
1ab07e55 71}
72
73//_______________________________________________________________________________________
74UInt_t AliPHOSPreprocessor::Process(TMap* /*valueSet*/)
75{
76 // process data retrieved by the Shuttle
516e7ab3 77
3c0265c1 78 TString runType = GetRunType();
79 Log(Form("Run type: %s",runType.Data()));
80
71eef542 81 if(runType=="LED") {
eba66a50 82 Bool_t ledOK = ProcessLEDRun();
4a7b5a26 83 Bool_t badmap_OK = FindBadChannelsEmc();
84 if(!badmap_OK) Log(Form("WARNING! FindBadChannels() completed with BAD status!"));
eba66a50 85 if(ledOK) return 0;
86 else
87 return 1;
88 }
9f43d7da 89
90 if(runType=="PHYSICS") {
91
9f43d7da 92 Bool_t calibEmc_OK = CalibrateEmc();
5db29457 93
4a7b5a26 94 if(calibEmc_OK) return 0;
9f43d7da 95 else
96 return 1;
97 }
98
99 Log(Form("Unknown run type %s. Do nothing and return OK.",runType.Data()));
100 return 0;
101
102}
103
104
105Bool_t AliPHOSPreprocessor::ProcessLEDRun()
106{
12661472 107 //Process LED run, update High Gain/Low Gain ratios.
9f43d7da 108
12661472 109 AliPHOSEmcCalibData calibData;
9f43d7da 110
111 TList* list = GetFileSources(kDAQ, "LED");
112 if(!list) {
113 Log("Sources list for LED run not found, exit.");
114 return kFALSE;
115 }
116
12661472 117 if(!list->GetEntries()) {
118 Log("Sources list for LED run is empty, exit.");
119 return kFALSE;
120 }
121
122 //Retrieve the last EMC calibration object
123 const AliPHOSEmcCalibData* clb=0;
124 AliCDBEntry* entryCalib = GetFromOCDB("Calib", "EmcGainPedestals");
125
126 if(!entryCalib)
127 Log(Form("Cannot find any AliCDBEntry for [Calib, EmcGainPedestals]!"));
128 else
129 clb = (AliPHOSEmcCalibData*)entryCalib->GetObject();
130
9f43d7da 131 TIter iter(list);
132 TObjString *source;
9f43d7da 133
134 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
12661472 135
9f43d7da 136 AliInfo(Form("found source %s", source->String().Data()));
137
138 TString fileName = GetFile(kDAQ, "LED", source->GetName());
139 AliInfo(Form("Got filename: %s",fileName.Data()));
140
141 TFile f(fileName);
142
143 if(!f.IsOpen()) {
144 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
145 return kFALSE;
146 }
5eccbe19 147
7708b003 148 TH1I* fFiredCells = (TH1I*)f.Get("fFiredCells");
5eccbe19 149
150 if(fFiredCells) {
7708b003 151 const Double_t nFiredCells = fFiredCells->GetMean();
5eccbe19 152 Log(Form("Number of fired cells per event is %.1f",nFiredCells));
7708b003 153 }
5eccbe19 154
9f43d7da 155 const Int_t nMod=5; // 1:5 modules
156 const Int_t nCol=56; //1:56 columns in each module
157 const Int_t nRow=64; //1:64 rows in each module
12661472 158
9f43d7da 159 for(Int_t mod=0; mod<nMod; mod++) {
160 for(Int_t col=0; col<nCol; col++) {
161 for(Int_t row=0; row<nRow; row++) {
a60e0e07 162
163 if(clb) {
0aa96a7a 164 Float_t hg2lg = clb->GetHighLowRatioEmc(5-mod,col+1,row+1);
165 Double_t coeff = clb->GetADCchannelEmc(5-mod,col+1,row+1);
166 calibData.SetADCchannelEmc(5-mod,col+1,row+1,coeff);
167 calibData.SetHighLowRatioEmc(5-mod,col+1,row+1,hg2lg);
a60e0e07 168 }
169
7708b003 170 //High Gain to Low Gain ratio
5eccbe19 171 Float_t ratio = HG2LG(mod,row,col,&f);
6d58d316 172 if(ratio != 16.) {
173 calibData.SetHighLowRatioEmc(5-mod,col+1,row+1,ratio);
5eccbe19 174 AliInfo(Form("mod %d iX %d iZ %d ratio %.3f\n",mod,row,col,ratio));
6d58d316 175 }
9f43d7da 176 }
177 }
178 }
179
12661472 180 } // end of loop over files
9f43d7da 181
4a7b5a26 182 //Store the updated High Gain/Low Gain ratios
12661472 183 AliCDBMetaData emcMetaData;
9f43d7da 184
27bf4542 185 //Data valid from current run until updated (validityInfinite=kTRUE)
186 Bool_t result = Store("Calib","EmcGainPedestals",&calibData,&emcMetaData,0,kTRUE);
9f43d7da 187 return result;
188
189}
190
191Float_t AliPHOSPreprocessor::HG2LG(Int_t mod, Int_t X, Int_t Z, TFile* f)
192{
193 //Calculates High gain to Low gain ratio
194 //for crystal at the position (X,Z) in the PHOS module mod.
195
196 char hname[128];
197 sprintf(hname,"%d_%d_%d",mod,X,Z);
198
199 TH1F* h1 = (TH1F*)f->Get(hname);
200 if(!h1) return 16.;
1eff11b2 201
6d58d316 202 if(h1->GetEntries()<2000.) return 16.;
1eff11b2 203
ff8d7889 204 if(h1->GetMaximum()<10.) h1->Rebin(4);
9f43d7da 205 if(h1->GetMaximum()<10.) return 16.;
ff8d7889 206
9f43d7da 207 Double_t max = h1->GetBinCenter(h1->GetMaximumBin()); // peak
208 Double_t xmin = max - (h1->GetRMS()/3);
209 Double_t xmax = max + (h1->GetRMS()/2);
210 // Double_t xmin = max - (h1->GetRMS());
211 // Double_t xmax = max + (h1->GetRMS());
212
213 TF1* gaus1 = new TF1("gaus1",rgaus,xmin,xmax,3);
214 gaus1->SetParNames("Constant","Mean","Sigma");
215 gaus1->SetParameter("Constant",h1->GetMaximum());
216 gaus1->SetParameter("Mean",max);
217 gaus1->SetParameter("Sigma",1.);
218 gaus1->SetLineColor(kBlue);
efa56161 219
220 Double_t mean_min = h1->GetXaxis()->GetXmin();
221 Double_t mean_max = h1->GetXaxis()->GetXmax();
222 gaus1->SetParLimits(1,mean_min,mean_max);
223
ff8d7889 224 h1->Fit(gaus1,"RQ+");
225 Double_t hg2lg = gaus1->GetParameter("Mean");
226 if( (hg2lg-mean_min<0.001) || (mean_max-hg2lg<0.001)) hg2lg=max;
227
4d38201a 228 AliInfo(Form("%s: %.1f entries, mean=%.3f, peak=%.3f, rms= %.3f. HG/LG = %.3f\n",
ff8d7889 229 h1->GetTitle(),h1->GetEntries(),h1->GetMean(),max,h1->GetRMS(),hg2lg));
9f43d7da 230
ff8d7889 231 return hg2lg;
232
9f43d7da 233}
234
235Bool_t AliPHOSPreprocessor::FindBadChannelsEmc()
236{
5db29457 237 //Loop over two systems: DAQ and HLT.
238 //For each system the same algorithm implemented in DoFindBadChannelsEmc() invokes.
9f43d7da 239
5db29457 240 TList* list=0;
241 TString path;
242
243 Int_t system[2] = { kDAQ, kHLT };
244 const char* sysn[] = { "DAQ","HLT" };
245 Bool_t result[2] = { kTRUE, kTRUE };
9f43d7da 246
5db29457 247 for (Int_t i=0; i<2; i++) {
e9122c68 248
249 if(system[i] == kHLT) continue;
250
5db29457 251 AliPHOSEmcBadChannelsMap badMap;
252 list = GetFileSources(system[i], "BAD_CHANNELS");
9f43d7da 253
5db29457 254 if(!list) {
255 Log(Form("%s sources list for BAD_CHANNELS not found!",sysn[i]));
256 result[i] = kFALSE;
257 continue;
258 }
9f43d7da 259
5db29457 260 if(!list->GetEntries()) {
261 Log(Form("Got empty sources list. It seems %s DA2 did not produce any files!",sysn[i]));
262 result[i] = kFALSE;
263 continue;
264 }
265
266 result[i] *= DoFindBadChannelsEmc(system[i],list,badMap);
267
268 // Store the bad channels map.
269
270 AliCDBMetaData md;
271 md.SetResponsible("Boris Polishchuk");
272
273 if(system[i] == kDAQ)
274 path = "Calib";
275 else
276 path = "HLT";
277
27bf4542 278 // Data valid from current run until being updated (validityInfinite=kTRUE)
279 result[i] *= Store(path.Data(), "EmcBadChannels", &badMap, &md, 0, kTRUE);
5db29457 280
9f43d7da 281 }
5db29457 282
283 if(result[0] || result[1]) return kTRUE;
284 else return kFALSE;
285}
286
287Bool_t AliPHOSPreprocessor::DoFindBadChannelsEmc(Int_t system, TList* list, AliPHOSEmcBadChannelsMap& badMap)
288{
289 //Creates the bad channels map for PHOS EMC.
290
291 // The file fileName contains histograms which have been produced by DA2 detector algorithm.
292 // It is a responsibility of the SHUTTLE framework to form the fileName.
9f43d7da 293
294 TIter iter(list);
295 TObjString *source;
296 char hnam[80];
297 TH1F* h1=0;
298
299 const Float_t fQualityCut = 1.;
300 Int_t nGoods[5] = {0,0,0,0,0};
301
302 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
303
304 AliInfo(Form("found source %s", source->String().Data()));
305
5db29457 306 TString fileName = GetFile(system, "BAD_CHANNELS", source->GetName());
9f43d7da 307 AliInfo(Form("Got filename: %s",fileName.Data()));
308
309 TFile f(fileName);
310
311 if(!f.IsOpen()) {
312 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
313 return kFALSE;
314 }
7708b003 315
5eccbe19 316 Log(Form("Begin check for bad channels."));
7708b003 317
9f43d7da 318 for(Int_t mod=0; mod<5; mod++) {
319 for(Int_t iX=0; iX<64; iX++) {
320 for(Int_t iZ=0; iZ<56; iZ++) {
5eccbe19 321
9f43d7da 322 sprintf(hnam,"%d_%d_%d_%d",mod,iX,iZ,1); // high gain
323 h1 = (TH1F*)f.Get(hnam);
324
325 if(h1) {
326 Double_t mean = h1->GetMean();
327
328 if(mean)
329 Log(Form("iX=%d iZ=%d gain=%d mean=%.3f\n",iX,iZ,1,mean));
330
331 if( mean>0 && mean<fQualityCut ) {
332 nGoods[mod]++;
333 }
334 else
335 badMap.SetBadChannel(mod+1,iZ+1,iX+1); //module, col,row
336 }
337
338 }
339 }
340
341 if(nGoods[mod])
342 Log(Form("Module %d: %d good channels.",mod,nGoods[mod]));
343 }
344
345
346 } // end of loop over sources
347
5db29457 348 return kTRUE;
9f43d7da 349}
350
351Bool_t AliPHOSPreprocessor::CalibrateEmc()
352{
5db29457 353 //Loop over two systems: DAQ and HLT.
354 //For each system the same algorithm implemented in DoCalibrateEmc() invokes.
9f43d7da 355
dd7c100d 356 AliPHOSEmcCalibData* lastCalib=0;
5db29457 357 const AliPHOSEmcBadChannelsMap* badMap=0;
358 AliCDBEntry* entryBCM=0;
dd7c100d 359 AliCDBEntry* entryEmc=0;
5db29457 360 TList* list=0;
361 TString path;
362
363 Int_t system[2] = { kDAQ, kHLT };
364 const char* sysn[] = { "DAQ","HLT" };
365 Bool_t result[2] = { kTRUE, kTRUE };
eba66a50 366
5db29457 367 for (Int_t i=0; i<2; i++) {
368
9bf526d6 369 if(system[i] == kHLT) continue;
370
5db29457 371 AliPHOSEmcCalibData calibData;
372 list = GetFileSources(system[i], "AMPLITUDES");
9f43d7da 373
5db29457 374 if(!list) {
375 Log(Form("%s sources list not found!",sysn[i]));
376 result[i] = kFALSE;
377 continue;
378 }
379
380 if(!list->GetEntries()) {
381 Log(Form("Got empty sources list. It seems %s DA1 did not produce any files!",sysn[i]));
382 result[i] = kFALSE;
383 continue;
384 }
385
386 // Retrieve the Bad Channels Map (BCM)
387
388 if(system[i] == kDAQ)
389 path = "Calib";
390 else
391 path = "HLT";
9f43d7da 392
5db29457 393 entryBCM = GetFromOCDB(path.Data(), "EmcBadChannels");
1ab07e55 394
5db29457 395 if(!entryBCM)
396 Log(Form("WARNING!! Cannot find any AliCDBEntry for [%s, EmcBadChannels]!",path.Data()));
397 else
398 badMap = (AliPHOSEmcBadChannelsMap*)entryBCM->GetObject();
399
400 if(!badMap)
401 Log(Form("WARNING!! Nothing for %s in AliCDBEntry. All cells considered GOOD!",sysn[i]));
402
dd7c100d 403 // Retrieve the last EMC calibration object
404
405 entryEmc = GetFromOCDB(path.Data(), "EmcGainPedestals");
406
407 if(!entryEmc)
408 Log(Form("Cannot find any EmcGainPedestals entry for this run and path %s",path.Data()));
409 else
410 lastCalib = (AliPHOSEmcCalibData*)entryEmc->GetObject();
411
af1c70b8 412 if(lastCalib)
413 result[i] *= DoCalibrateEmc(system[i],list,badMap,*lastCalib);
414 else
dd7c100d 415 result[i] *= DoCalibrateEmc(system[i],list,badMap,calibData);
af1c70b8 416
5db29457 417 //Store EMC calibration data
5db29457 418 AliCDBMetaData emcMetaData;
dd7c100d 419
420 if(lastCalib)
af1c70b8 421 result[i] *= Store(path.Data(), "EmcGainPedestals", lastCalib, &emcMetaData, 0, kFALSE);
dd7c100d 422 else
af1c70b8 423 result[i] *= Store(path.Data(), "EmcGainPedestals", &calibData, &emcMetaData, 0, kFALSE);
424
425 //Store reference data
d70d764a 426 Bool_t refOK = StoreReferenceEmc(system[i],list);
427 if(refOK) Log(Form("Reference data for %s amplitudes successfully stored.",sysn[i]));
dd7c100d 428
9f43d7da 429 }
5db29457 430
431 if(result[0] || result[1]) return kTRUE;
432 else return kFALSE;
433}
9f43d7da 434
d70d764a 435Bool_t AliPHOSPreprocessor::StoreReferenceEmc(Int_t system, TList* list)
af1c70b8 436{
437 //Put 2D calibration histograms (E vs Time) prepared by DAQ/HLT to the reference storage.
438 //system is DAQ or HLT, TList is the list of FES sources.
439
d70d764a 440 if(system!=kDAQ) return kFALSE;
af1c70b8 441
d70d764a 442 TObjString *source = dynamic_cast<TObjString *> (list->First());
443 if(!source) return kFALSE;
444
445 TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
446
447 Bool_t resultRef = StoreReferenceFile(fileName.Data(),"CalibRefPHOS.root");
af1c70b8 448 return resultRef;
449
450}
451
30a0d5a2 452
5db29457 453Bool_t AliPHOSPreprocessor::DoCalibrateEmc(Int_t system, TList* list, const AliPHOSEmcBadChannelsMap* badMap, AliPHOSEmcCalibData& calibData)
454{
455 // Return kTRUE if OK.
456 // I. Calculates the set of calibration coefficients to equalyze the mean energies deposited at high gain.
457 // II. Extracts High_Gain/Low_Gain ratio for each channel.
30a0d5a2 458
5db29457 459 // The file fileName contains histograms which have been produced by DA1 detector algorithm.
460 // It is a responsibility of the SHUTTLE framework to form the fileName.
30a0d5a2 461
5db29457 462 gRandom->SetSeed(0); //the seed is set to the current machine clock!
af1c70b8 463 Int_t minEntries=1000; // recalculate calibration coeff. if Nentries > minEntries.
30a0d5a2 464
516e7ab3 465 TIter iter(list);
466 TObjString *source;
3c0265c1 467
516e7ab3 468 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
469 AliInfo(Form("found source %s", source->String().Data()));
f97eb136 470
5db29457 471 TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
516e7ab3 472 AliInfo(Form("Got filename: %s",fileName.Data()));
f97eb136 473
516e7ab3 474 TFile f(fileName);
f97eb136 475
516e7ab3 476 if(!f.IsOpen()) {
477 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
5db29457 478 return kFALSE;
516e7ab3 479 }
9f43d7da 480
516e7ab3 481 const Int_t nMod=5; // 1:5 modules
482 const Int_t nCol=56; //1:56 columns in each module
483 const Int_t nRow=64; //1:64 rows in each module
484
485 Double_t coeff;
486 char hnam[80];
450717bc 487 TH2F* h2=0;
488 TH1D* h1=0;
3c0265c1 489
516e7ab3 490 //Get the reference histogram
491 //(author: Gustavo Conesa Balbastre)
492
493 TList * keylist = f.GetListOfKeys();
494 Int_t nkeys = f.GetNkeys();
495 Bool_t ok = kFALSE;
496 TKey *key;
516e7ab3 497 Int_t ikey = 0;
498 Int_t counter = 0;
450717bc 499 TH1D* hRef = 0;
500
516e7ab3 501 //Check if the file contains any histogram
502
503 if(nkeys< 2){
504 Log(Form("Not enough histograms (%d) for calibration.",nkeys));
4a7b5a26 505 return 1; // it's not fatal! May be short run..
f97eb136 506 }
450717bc 507
516e7ab3 508 while(!ok){
509 ikey = gRandom->Integer(nkeys);
510 key = (TKey*)keylist->At(ikey);
450717bc 511 TObject* obj = f.Get(key->GetName());
512 TString cname(obj->ClassName());
513 if(cname == "TH2F") {
514 h2 = (TH2F*)obj;
515 TString htitl = h2->GetTitle();
516 if(htitl.Contains("and gain 1")) {
517 hRef = h2->ProjectionX();
dd7c100d 518 hRef->GetXaxis()->SetRangeUser(10.,1000.); // to cut off saturation peak and noise
450717bc 519 // Check if the reference histogram has too little statistics
dd7c100d 520 if(hRef->GetMean() && hRef->GetEntries()>minEntries) ok=kTRUE;
30a0d5a2 521
522 const TString delim = "_";
523 TString str = hRef->GetName();
524 TObjArray* tks = str.Tokenize(delim);
525 const Int_t md = ((TObjString*)tks->At(0))->GetString().Atoi();
526 const Int_t X = ((TObjString*)tks->At(1))->GetString().Atoi();
527 const Int_t Z = ((TObjString*)tks->At(2))->GetString().Atoi();
528
529 if(badMap) {
6d58d316 530 if(badMap->IsBadChannel(5-md,Z+1,X+1)) {
30a0d5a2 531 AliInfo(Form("Cell mod=%d col=%d row=%d is bad. Histogram %s rejected.",
6d58d316 532 5-md,Z+1,X+1,hRef->GetName()));
30a0d5a2 533 ok=kFALSE;
534 }
535 }
536
450717bc 537 }
516e7ab3 538 }
b8093789 539
450717bc 540 counter++;
b8093789 541
542 if(!ok && counter > nkeys){
543 Log("No histogram with enough statistics for reference. Exit.");
4a7b5a26 544 return 1; // Not fatal, just wait..
b8093789 545 }
516e7ab3 546 }
3c0265c1 547
516e7ab3 548 Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
549 hRef->GetName(),hRef->GetEntries(),
550 hRef->GetMean(),hRef->GetRMS()));
551
552 Double_t refMean=hRef->GetMean();
3c0265c1 553
516e7ab3 554 // Calculates relative calibration coefficients for all non-zero channels
3c0265c1 555
516e7ab3 556 for(Int_t mod=0; mod<nMod; mod++) {
557 for(Int_t col=0; col<nCol; col++) {
558 for(Int_t row=0; row<nRow; row++) {
450717bc 559
560 sprintf(hnam,"%d_%d_%d_1",mod,row,col); // high gain!
561 h2 = (TH2F*)f.Get(hnam);
4a7b5a26 562
516e7ab3 563 //TODO: dead channels exclusion!
450717bc 564 if(h2) {
565 h1 = h2->ProjectionX();
dd7c100d 566 h1->GetXaxis()->SetRangeUser(10.,1000.); //to cut off saturation peak and noise
450717bc 567 coeff = h1->GetMean()/refMean;
dd7c100d 568 if(coeff>0 && h1->GetEntries()>minEntries) {
6d58d316 569 calibData.SetADCchannelEmc(5-mod,col+1,row+1,0.005/coeff);
dd7c100d 570 AliInfo(Form("mod %d col %d row %d coeff %f\n",mod,col,row,coeff));
571 }
516e7ab3 572 }
1ab07e55 573 }
1ab07e55 574 }
575 }
516e7ab3 576
577 f.Close();
1ab07e55 578 }
516e7ab3 579
5db29457 580 return 1;
1ab07e55 581}
450717bc 582
9f43d7da 583
450717bc 584