Line data Source code
1 : // $Id$
2 :
3 : //**************************************************************************
4 : //* This file is property of and copyright by the ALICE HLT Project *
5 : //* ALICE Experiment at CERN, All rights reserved. *
6 : //* *
7 : //* Primary Authors: Timur Pocheptsov <Timur.Pocheptsov@cern.ch> *
8 : //* Matthias Richter <Matthias.Richter@cern.ch>
9 : //* for The ALICE HLT Project. *
10 : //* *
11 : //* Permission to use, copy, modify and distribute this software and its *
12 : //* documentation strictly for non-commercial purposes is hereby granted *
13 : //* without fee, provided that the above copyright notice appears in all *
14 : //* copies and that both the copyright notice and this permission notice *
15 : //* appear in the supporting documentation. The authors make no claims *
16 : //* about the suitability of this software for any purpose. It is *
17 : //* provided "as is" without express or implied warranty. *
18 : //**************************************************************************
19 :
20 : /// @file AliHLTTTreeProcessor.cxx
21 : /// @author Timur Pocheptsov, Matthias Richter
22 : /// @date 05.07.2010
23 : /// @brief Generic component for data collection in a TTree
24 :
25 : #include <cerrno>
26 : #include <memory>
27 :
28 : #include "AliHLTTTreeProcessor.h"
29 : #include "AliHLTErrorGuard.h"
30 : #include "TDirectory.h"
31 : #include "TDatime.h"
32 : #include "TString.h"
33 : #include "TTree.h"
34 : #include "TH1.h"
35 : #include "TStopwatch.h"
36 : #include "TUUID.h"
37 : #include "TSystem.h"
38 : #include "TRandom3.h"
39 :
40 : /** ROOT macro for the implementation of ROOT specific class methods */
41 126 : ClassImp(AliHLTTTreeProcessor)
42 :
43 : AliHLTTTreeProcessor::AliHLTTTreeProcessor()
44 0 : : AliHLTProcessor(),
45 0 : fDefinitions(),
46 0 : fTree(0),
47 0 : fMaxEntries(kMaxEntries),
48 0 : fPublishInterval(kInterval),
49 0 : fLastTime(0),
50 0 : fpEventTimer(NULL),
51 0 : fpCycleTimer(NULL),
52 0 : fMaxMemory(700000),
53 0 : fMaxEventTime(0),
54 0 : fNofEventsForce(0),
55 0 : fForcedEventsCount(0),
56 0 : fSkippedEventsCount(0),
57 0 : fNewEventsCount(0),
58 0 : fUniqueId(0),
59 0 : fIgnoreCycleTime(10),
60 0 : fCycleTimeFactor(1.0)
61 0 : {
62 : // see header file for class documentation
63 : // or
64 : // refer to README to build package
65 : // or
66 : // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
67 0 : }
68 :
69 : const AliHLTUInt32_t AliHLTTTreeProcessor::fgkTimeScale=1000000; // ticks per second
70 :
71 : AliHLTTTreeProcessor::~AliHLTTTreeProcessor()
72 0 : {
73 : // see header file for class documentation
74 0 : }
75 :
76 : AliHLTComponentDataType AliHLTTTreeProcessor::GetOutputDataType()
77 : {
78 : // get the component output data type
79 0 : return kAliHLTDataTypeHistogram;
80 : }
81 :
82 : void AliHLTTTreeProcessor::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
83 : {
84 : // get the output size estimator
85 : //
86 0 : if (!fDefinitions.size()) {
87 0 : HLTError("Can not calculate output data size, no histogram definitions were provided");
88 : return;
89 : }
90 :
91 0 : constBase = 0;
92 0 : for (list_const_iterator i = fDefinitions.begin(); i != fDefinitions.end(); ++i)
93 0 : constBase += i->GetSize();
94 :
95 0 : inputMultiplier = 1.;
96 0 : }
97 :
98 : int AliHLTTTreeProcessor::DoInit(int argc, const char** argv)
99 : {
100 : // init component
101 : // ask child to create the tree.
102 : int iResult = 0;
103 :
104 : // component configuration
105 : //Stage 1: default initialization.
106 : //"Default" (for derived component) histograms.
107 0 : FillHistogramDefinitions();
108 : //Default values.
109 0 : fMaxEntries = kMaxEntries;
110 0 : fPublishInterval = kInterval;
111 0 : fLastTime = 0;
112 : //Stage 2: OCDB.
113 0 : TString cdbPath("HLT/ConfigHLT/");
114 0 : cdbPath += GetComponentID();
115 : //
116 0 : iResult = ConfigureFromCDBTObjString(cdbPath);
117 : //
118 0 : if (iResult < 0)
119 0 : return iResult;
120 : //Stage 3: command line arguments.
121 0 : if (argc && (iResult = ConfigureFromArgumentString(argc, argv)) < 0)
122 0 : return iResult;
123 :
124 : // calculating a unique id from the hostname and process id
125 : // used for identifying output of multiple components
126 0 : TUUID guid = GenerateGUID();
127 0 : union
128 : {
129 : UChar_t buf[16];
130 : UInt_t bufAsInt[4];
131 : };
132 0 : guid.GetUUID(buf);
133 0 : fUniqueId = bufAsInt[0];
134 :
135 0 : if (!fTree) {
136 : // originally foreseen to pass the arguments to the function, however
137 : // this is not appropriate. Argument scan via overloaded function
138 : // ScanConfigurationArgument
139 0 : std::auto_ptr<TTree> ptr(CreateTree(0, NULL));
140 0 : if (ptr.get()) {
141 0 : ptr->SetDirectory(0);
142 0 : ptr->SetCircular(fMaxEntries);
143 0 : fTree = ptr.release();
144 : } else //No way to process error correctly - error is unknown here.
145 0 : return -EINVAL;
146 0 : } else {
147 0 : HLTError("fTree pointer must be null before DoInit call");
148 0 : return -EINVAL;
149 : }
150 :
151 0 : if (iResult>=0 && fMaxEventTime>0) {
152 0 : fpEventTimer=new TStopwatch;
153 0 : if (fpEventTimer) {
154 0 : fpEventTimer->Reset();
155 : }
156 0 : fpCycleTimer=new TStopwatch;
157 0 : if (fpCycleTimer) {
158 0 : fpCycleTimer->Reset();
159 : }
160 : }
161 0 : fSkippedEventsCount=0;
162 :
163 0 : return iResult;
164 0 : }
165 :
166 : int AliHLTTTreeProcessor::DoDeinit()
167 : {
168 : // cleanup component
169 0 : delete fTree;
170 0 : fTree = 0;
171 0 : fDefinitions.clear();
172 :
173 0 : if (fpEventTimer) delete fpEventTimer;
174 0 : fpEventTimer=NULL;
175 0 : if (fpCycleTimer) delete fpCycleTimer;
176 0 : fpCycleTimer=NULL;
177 :
178 0 : return 0;
179 : }
180 :
181 : int AliHLTTTreeProcessor::DoEvent(const AliHLTComponentEventData& evtData, AliHLTComponentTriggerData& trigData)
182 : {
183 : //Process event and publish histograms.
184 0 : AliHLTUInt32_t eventType=0;
185 0 : if (!IsDataEvent(&eventType) && eventType!=gkAliEventTypeEndOfRun) return 0;
186 :
187 : //I'm pretty sure, that if fTree == 0 (DoInit failed) DoEvent is not called.
188 : //But interface itself does not force you to call DoInit before DoEvent, so,
189 : //I make this check explicit.
190 0 : if (!fTree) {
191 0 : HLTError("fTree is a null pointer, try to call AliHLTTTreeProcessor::DoInit first.");
192 0 : return -EINVAL;//-ENULLTREE? :)
193 : }
194 :
195 : AliHLTUInt32_t averageEventTime=0;
196 : AliHLTUInt32_t averageCycleTime=0;
197 :
198 : int fillingtime=0;
199 : int publishtime=0;
200 : bool bDoFilling=true;
201 : bool bDoPublishing=false;
202 : const int cycleResetInterval=1000;
203 0 : if (fpEventTimer && fpCycleTimer) {
204 0 : averageEventTime=AliHLTUInt32_t(fpEventTimer->RealTime()*fgkTimeScale)/(GetEventCount()+1);
205 0 : fillingtime=int(fpEventTimer->RealTime()*fgkTimeScale);
206 : publishtime=fillingtime;
207 0 : fpEventTimer->Start(kFALSE);
208 0 : fpCycleTimer->Stop();
209 0 : averageCycleTime=AliHLTUInt32_t(fpCycleTimer->RealTime()*fgkTimeScale)/((GetEventCount()%cycleResetInterval)+1);
210 : // adapt processing to 3/4 of the max time
211 0 : bDoFilling=4*averageEventTime<3*fMaxEventTime ||
212 0 : (averageEventTime<fCycleTimeFactor*averageCycleTime && fpCycleTimer->RealTime()>fIgnoreCycleTime);
213 0 : if (fNofEventsForce>0 && fForcedEventsCount<fNofEventsForce) {
214 0 : fForcedEventsCount++;
215 : bDoFilling=true;
216 0 : }
217 : }
218 :
219 : // FIXME: there is still an unclear increase in memory consumption, even if the number of entries
220 : // in the tree is restricted. Valgrind studies did not show an obvious memory leak. This is likely
221 : // to be caused by something deep in the Root TTree functionality and needs to be studied in detail.
222 0 : ProcInfo_t ProcInfo;
223 0 : gSystem->GetProcInfo(&ProcInfo);
224 0 : if (ProcInfo.fMemResident>fMaxMemory) bDoFilling=false;
225 :
226 : // process input data blocks and fill the tree
227 : int iResult = 0;
228 0 : if (eventType!=gkAliEventTypeEndOfRun) {
229 0 : if (bDoFilling) {iResult=FillTree(fTree, evtData, trigData); fNewEventsCount++;}
230 0 : else fSkippedEventsCount++;
231 : }
232 0 : if (fpEventTimer) {
233 0 : fpEventTimer->Stop();
234 0 : fillingtime=int(fpEventTimer->RealTime()*fgkTimeScale)-fillingtime;
235 0 : if (fillingtime<0) fillingtime=0;
236 0 : fpEventTimer->Start(kFALSE);
237 : }
238 :
239 0 : if (iResult < 0) {
240 0 : ALIHLTERRORGUARD(5, "FillTree failed with %d, first event %d", iResult, GetEventCount());
241 0 : return iResult;
242 : }
243 :
244 0 : const TDatime time;
245 :
246 0 : if (( time.Get() - fLastTime > fPublishInterval && fNewEventsCount>0) ||
247 0 : eventType==gkAliEventTypeEndOfRun) {
248 0 : if ((bDoPublishing=fLastTime>0)) { // publish earliest after the first interval but set the timer
249 :
250 0 : for (list_const_iterator i = fDefinitions.begin(); i != fDefinitions.end(); ++i) {
251 0 : if (TH1* h = CreateHistogram(*i)) {
252 : //I do not care about errors here - since I'm not able
253 : //to rollback changes.
254 : // TODO: in case of -ENOSPC et the size of the last object by calling
255 : // GetLastObjectSize() and accumulate the necessary output buffer size
256 0 : PushBack(h, GetOriginDataType(), GetDataSpec());
257 0 : delete h;
258 : }
259 : }
260 0 : unsigned eventcount=GetEventCount()+1;
261 0 : HLTBenchmark("publishing %d histograms, %d entries in tree, %d new events since last publishing, accumulated %d of %d events (%.1f%%)", fDefinitions.size(), fTree->GetEntriesFast(), fNewEventsCount, eventcount-fSkippedEventsCount, eventcount, eventcount>0?(100*float(eventcount-fSkippedEventsCount)/eventcount):0);
262 0 : fNewEventsCount=0;
263 0 : HLTBenchmark("current memory usage %d %d", ProcInfo.fMemResident, ProcInfo.fMemVirtual);
264 0 : }
265 :
266 0 : fLastTime=time.Get();
267 0 : if (fLastTime==0) {
268 : // choose a random offset at beginning to equalize traffic for multiple instances
269 : // of the component
270 0 : gRandom->SetSeed(fUniqueId);
271 0 : fLastTime-=gRandom->Integer(fPublishInterval);
272 0 : }
273 : }
274 :
275 0 : if (fpEventTimer) {
276 0 : fpEventTimer->Stop();
277 0 : publishtime=int(fpEventTimer->RealTime()*fgkTimeScale)-publishtime;
278 0 : if (publishtime>fillingtime) publishtime-=fillingtime;
279 : else publishtime=0;
280 :
281 0 : averageEventTime=AliHLTUInt32_t(fpEventTimer->RealTime()*fgkTimeScale)/(GetEventCount()+1);
282 :
283 : // info output once every 5 seconds
284 : static UInt_t lastTime=0;
285 0 : if (time.Get()-lastTime>5 ||
286 0 : eventType==gkAliEventTypeEndOfRun ||
287 0 : bDoPublishing) {
288 0 : lastTime=time.Get();
289 0 : unsigned eventcount=GetEventCount()+1;
290 0 : HLTBenchmark("filling time %d us, publishing time %d, average total processing time %d us, cycle time %d us, accumulated %d of %d events (%.1f%%)", fillingtime, publishtime, averageEventTime, averageCycleTime, eventcount-fSkippedEventsCount, eventcount, eventcount>0?(100*float(eventcount-fSkippedEventsCount)/eventcount):0);
291 0 : }
292 : }
293 0 : if (fpCycleTimer) {
294 0 : bool bReset=(GetEventCount()%cycleResetInterval)==0;
295 0 : fpCycleTimer->Start(bReset);
296 0 : }
297 :
298 : return iResult;
299 0 : }
300 :
301 : int AliHLTTTreeProcessor::ScanConfigurationArgument(int argc, const char** argv)
302 : {
303 : // scan one argument and its parameters from the list
304 : // return number of processed entries.
305 : // possible arguments:
306 : // -maxentries number
307 : // -interval number
308 : // -histogram name -size number -expression expression [-title expression ] -cut expression ][-opt option]
309 : // As soon as "-histogram" found, -size and -expression and -outtype are required,
310 : // cut and option can be omitted.
311 0 : if (argc <= 0)
312 0 : return 0;
313 :
314 0 : std::list<AliHLTHistogramDefinition> newDefs;
315 0 : AliHLTHistogramDefinition def;
316 :
317 : int i = 0;
318 0 : int maxEntries = fMaxEntries;
319 :
320 0 : while (i < argc) {
321 0 : const TString argument(argv[i]);
322 :
323 0 : if (argument.CompareTo("-maxentries") == 0) { //1. Max entries argument for TTree.
324 0 : if (i + 1 == argc) {
325 0 : HLTError("Numeric value for '-maxentries' is expected");
326 0 : return -EPROTO;
327 : }
328 : //Next must be a number.
329 : //TString returns 0 (number) even if string contains non-numeric symbols.
330 0 : maxEntries = TString(argv[i + 1]).Atoi();
331 0 : if (maxEntries <= 0) {
332 0 : HLTError("Bad value for '-maxentries': %d", maxEntries);
333 0 : return -EPROTO;
334 : }
335 :
336 0 : i += 2;
337 0 : } else if (argument.CompareTo("-interval") == 0) { //2. Interval argument for publishing.
338 0 : if (i + 1 == argc) {
339 0 : HLTError("Numeric value for '-interval' is expected");
340 0 : return -EPROTO;
341 : }
342 :
343 0 : const Int_t interval = TString(argv[i + 1]).Atoi();
344 0 : if (interval < 0) {
345 0 : HLTError("Bad value for '-interval' argument: %d", interval);
346 0 : return -EPROTO;
347 : }
348 :
349 0 : fPublishInterval = interval;
350 :
351 0 : i += 2;
352 0 : } else if (argument.CompareTo("-maxeventtime") == 0) { // max average processing time in us
353 0 : if (i + 1 == argc) {
354 0 : HLTError("Numeric value for '-maxeventtime' is expected");
355 0 : return -EPROTO;
356 : }
357 :
358 0 : const Int_t time = TString(argv[i + 1]).Atoi();
359 0 : if (time < 0) {
360 0 : HLTError("Bad value for '-maxeventtime' argument: %d", time);
361 0 : return -EPROTO;
362 : }
363 :
364 0 : fMaxEventTime = time;
365 :
366 0 : i += 2;
367 0 : } else if (argument.CompareTo("-forced-events") == 0) { // number of forced events
368 0 : if (i + 1 == argc) {
369 0 : HLTError("Numeric value for '-forced-events' is expected");
370 0 : return -EPROTO;
371 : }
372 :
373 0 : const Int_t count = TString(argv[i + 1]).Atoi();
374 0 : if (count < 0) {
375 0 : HLTError("Bad value for '-forced-events' argument: %d", count);
376 0 : return -EPROTO;
377 : }
378 :
379 0 : fNofEventsForce = count;
380 0 : fForcedEventsCount=0;
381 :
382 0 : i += 2;
383 0 : } else if (argument.CompareTo("-ignore-cycletime") == 0) { // ignore cycle time for n sec
384 0 : if (i + 1 == argc) {
385 0 : HLTError("Numeric value for '-ignore-cycletime' is expected");
386 0 : return -EPROTO;
387 : }
388 :
389 0 : const Int_t time = TString(argv[i + 1]).Atoi();
390 0 : if (time < 0) {
391 0 : HLTError("Bad value for '-ignore-cycletime' argument: %d", time);
392 0 : return -EPROTO;
393 : }
394 :
395 0 : fIgnoreCycleTime = time;
396 0 : i += 2;
397 0 : } else if (argument.CompareTo("-maxmemory") == 0) { // maximum of memory in kByte to be used by the component
398 0 : if (i + 1 == argc) {
399 0 : HLTError("Numeric value for '-maxmemory' is expected");
400 0 : return -EPROTO;
401 : }
402 :
403 0 : const Int_t mem = TString(argv[i + 1]).Atoi();
404 0 : if (mem < 0) {
405 0 : HLTError("Bad value for '-maxmemory' argument: %d", time);
406 0 : return -EPROTO;
407 : }
408 :
409 0 : fMaxMemory = mem;
410 0 : i += 2;
411 0 : } else if (argument.CompareTo("-cycletime-factor") == 0) { // weight factor for cycle time
412 0 : if (i + 1 == argc) {
413 0 : HLTError("Numeric value for '-cycletime-factor' is expected");
414 0 : return -EPROTO;
415 : }
416 :
417 0 : const Float_t factor = TString(argv[i + 1]).Atof();
418 0 : if (factor < 0) {
419 0 : HLTError("Bad value for '-cycletime-factor' argument: %f", factor);
420 0 : return -EPROTO;
421 : }
422 :
423 0 : fCycleTimeFactor = factor;
424 0 : i += 2;
425 0 : } else if (argument.CompareTo("-histogram") == 0) { //3. Histogramm definition.
426 0 : const int nParsed = ParseHistogramDefinition(argc, argv, i, def);
427 0 : if (!nParsed)
428 0 : return -EPROTO;
429 :
430 0 : newDefs.push_back(def);
431 :
432 0 : i += nParsed;
433 0 : } else {
434 0 : HLTError("Unknown argument %s", argument.Data());
435 0 : return -EPROTO;
436 : }
437 0 : }
438 :
439 0 : if (maxEntries != fMaxEntries) {
440 0 : fMaxEntries = maxEntries;
441 0 : if (fTree) {
442 0 : fTree->Reset();
443 0 : fTree->SetCircular(fMaxEntries);
444 : }
445 : }
446 :
447 0 : if (newDefs.size())
448 0 : fDefinitions.swap(newDefs);
449 :
450 0 : return i;
451 0 : }
452 :
453 : TH1* AliHLTTTreeProcessor::CreateHistogram(const AliHLTHistogramDefinition& d)
454 : {
455 :
456 : // create a histogram from the tree
457 0 : if (!fTree) {
458 0 : HLTError("fTree is a null pointer, try to call AliHLTTTreeProcessor::DoInit first.");
459 0 : return 0;
460 : }
461 :
462 0 : TString histName(d.GetName());
463 0 : if (!histName.Contains("(")) {
464 : //Without number of bins, the histogram will be "fixed"
465 : //and most of values can go to underflow/overflow bins,
466 : //since kCanRebin will be false.
467 0 : histName += TString::Format("(%d)", Int_t(kDefaultNBins));
468 0 : }
469 :
470 0 : const Long64_t rez = fTree->Project(histName.Data(), d.GetExpression().Data(), d.GetCut().Data(), d.GetDrawOption().Data());
471 :
472 0 : if (rez == -1) {
473 0 : HLTError("TTree::Project failed");
474 0 : return 0;
475 : }
476 :
477 : //Now, cut off the binning part of a name
478 0 : histName = histName(0, histName.Index("("));
479 0 : TH1 * hist = dynamic_cast<TH1*>(gDirectory->Get(histName.Data()));
480 0 : if (!hist) {
481 0 : const TString msg(Form("Hist %s is a null pointer, selection was %s, strange name or hist's type\n", histName.Data(), d.GetExpression().Data()));
482 0 : HLTError(msg.Data());
483 0 : }else if (d.GetDrawOption().Length()) {
484 0 : hist->SetOption(d.GetDrawOption().Data());
485 : }
486 :
487 : //Reformatting the histogram name
488 0 : TString str2=d.GetCut().Data();
489 0 : str2.ReplaceAll("Track_", "");
490 0 : str2.ReplaceAll("&&", " ");
491 0 : str2 = histName+" "+str2;
492 0 : hist->SetTitle(str2);
493 :
494 0 : if(d.GetTitle().Length()){
495 :
496 : //removing underscore
497 0 : TString axis=d.GetTitle().Data();
498 0 : axis.ReplaceAll("_{T}", "underscore{T}");
499 0 : axis.ReplaceAll("_", " ");
500 0 : axis.ReplaceAll("underscore{T}", "_{T}");
501 :
502 0 : hist->SetXTitle(axis);
503 0 : hist->GetXaxis()->CenterTitle();
504 0 : }
505 : return hist;
506 0 : }
507 :
508 : int AliHLTTTreeProcessor::ParseHistogramDefinition(int argc, const char** argv, int pos, AliHLTHistogramDefinition& dst)const
509 : {
510 : //Histogram-definition:
511 : // -histogram name -size number -expression expression [-title expression][-cut expression][-opt option]
512 :
513 : //at pos we have '-histogram', at pos + 1 must be the name.
514 0 : if (pos + 1 == argc) {
515 0 : HLTError("Bad histogram definition, histogram name is expected");
516 0 : return 0;
517 : }
518 :
519 0 : dst.SetName(argv[pos + 1]);
520 0 : pos += 2;
521 :
522 : //At pos must be '-size', and number at pos + 1.
523 0 : if (pos == argc || TString(argv[pos]).CompareTo("-size")) {
524 0 : HLTError("Bad histogram definition, '-size' is expected");
525 0 : return 0;
526 : }
527 :
528 0 : if (pos + 1 == argc) {
529 0 : HLTError("Bad histogram definition, size is expected");
530 0 : return 0;
531 : }
532 :
533 0 : dst.SetSize(TString(argv[pos + 1]).Atoi());
534 0 : if (dst.GetSize() <= 0) {
535 0 : HLTError("Bad histogram definition, positive size is required");
536 0 : return 0;
537 : }
538 :
539 0 : pos += 2;
540 : //At pos must be '-expression', and expression at pos + 1.
541 0 : if (pos == argc || TString(argv[pos]).CompareTo("-expression")) {
542 0 : HLTError("Bad histogram definition, '-expression' is expected");
543 0 : return 0;
544 : }
545 :
546 0 : if (pos + 1 == argc) {
547 0 : HLTError("Bad histogram definition, expression is expected");
548 0 : return 0;
549 : }
550 :
551 0 : dst.SetExpression(argv[pos + 1]);
552 0 : pos += 2;
553 :
554 : int processed = 6;
555 0 : dst.SetTitle("");
556 0 : dst.SetCut("");
557 0 : dst.SetDrawOption("");
558 :
559 : //remaining options can be the title, cut and Draw option.
560 : //title must be first
561 0 : if (pos + 1 >= argc){
562 0 : return processed;
563 : }
564 0 : if (TString(argv[pos]).CompareTo("-title") == 0) {
565 0 : dst.SetTitle(argv[pos + 1]);
566 0 : pos += 2;
567 : processed += 2;
568 0 : }
569 :
570 : //cut must be second.
571 0 : if (pos + 1 >= argc)
572 0 : return processed;
573 :
574 0 : if (TString(argv[pos]).CompareTo("-cut") == 0) {
575 0 : dst.SetCut(argv[pos + 1]);
576 0 : pos += 2;
577 0 : processed += 2;
578 0 : }
579 :
580 0 : if (pos + 1 >= argc)
581 0 : return processed;
582 :
583 0 : if (TString(argv[pos]).CompareTo("-opt") == 0) {
584 0 : dst.SetDrawOption(argv[pos + 1]);
585 0 : processed += 2;
586 0 : }
587 :
588 0 : return processed;
589 0 : }
|