]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveFMDLoader.cxx
0a53e57c0feb5442b3d0d246bdac058774d6540a
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveFMDLoader.cxx
1 /**************************************************************************
2  * Copyright(c) 2008, Christian Holm Christensen                          *
3  *                                                                        *
4  * Author: Christian Holm Christensen.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation for any purposes is hereby granted without fee,          * 
9  * provided that the above copyright notice appears in all copies and     *
10  * that both the copyright notice and this permission notice remains      *
11  * intact.  The authors make no claims about the suitability of this      * 
12  * software for any purpose. It is provided "as is" without express or    *
13  * implied warranty.                                                      * 
14  **************************************************************************/
15 /* $Id$ */
16 /** @file    AliEveFMDLoader.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Sun Mar 26 17:59:18 2006
19     @brief   Implementation of AliEveFMDLoader singleton class 
20 */
21 //____________________________________________________________________
22 //                                                                          
23 // Forward Multiplicity Detector based on Silicon wafers. This class
24 // is the loader for the event display. 
25 //
26 // This class is a singleton, meaning that there's only one instance
27 // of this.  This is done to speed up the processing by putting all
28 // things that are needed every time into the constructor. 
29 //
30 #include "AliRunLoader.h"
31 #include "EveBase/AliEveEventManager.h"
32 #include "AliEveFMDLoader.h"
33 #include "../FMD/AliFMDUShortMap.h"
34 #include "../FMD/AliFMDBoolMap.h"
35 #include "../FMD/AliFMDGeometry.h"
36 #include "../FMD/AliFMDParameters.h"
37 #include "../FMD/AliFMDDetector.h"
38 #include "../FMD/AliFMDRing.h"
39 #include "../FMD/AliFMDBaseDigit.h"
40 #include "../FMD/AliFMDDigit.h"
41 #include "../FMD/AliFMDRawReader.h"
42 #include "../FMD/AliFMDHit.h"
43 #include "AliESDEvent.h"
44 #include "AliESDFMD.h"
45 #include "AliLog.h"
46 #include "AliRawReader.h"
47 #include <TClonesArray.h>
48 #include <TTree.h>
49 #include <TGeoShape.h>
50 #include <TGeoManager.h>
51 #include <TEveGeoNode.h>
52 #include <TEveBoxSet.h>
53 #include <TEveQuadSet.h>
54 #include <TEveManager.h>
55 #include <TEveUtil.h>
56 #include <TStyle.h>
57 #include <TMath.h>
58 #include <iostream>
59
60 // Some very private variables. 
61 namespace 
62 {
63   const Char_t* kDetector = "FMD%d";
64   const Char_t* kRing     = "FMD%d%c";
65   const Char_t* kModule   = "FMD%d%c[%02d-%02d]";
66   const Char_t* kSector   = "FMD%d%c[%02d] %s";
67
68   const Char_t* kHits     = "Hits";
69   const Char_t* kDigits   = "Digits";
70   const Char_t* kRaw      = "Raw";
71   const Char_t* kESD      = "ESD";
72   
73 }
74
75 //____________________________________________________________________
76 AliEveFMDLoader* AliEveFMDLoader::fgInstance = 0;
77
78 //____________________________________________________________________
79 AliEveFMDLoader* AliEveFMDLoader::Instance()
80 {
81   // Get the singleton instance.  If the instance has not been
82   // instantised yet, it will be after this call.
83   if (!fgInstance) 
84     fgInstance = new AliEveFMDLoader();
85   return fgInstance;
86 }
87
88 //____________________________________________________________________
89 AliEveFMDLoader::AliEveFMDLoader(const char* name, Bool_t useBoxes, 
90                                  Bool_t old)
91   : TEveElementList(name, 0), 
92     fHitPalette(0, 1000),
93     fDigitPalette(0, 1023), 
94     fMultPalette(0, 20),
95     fUseBoxDigits(useBoxes), 
96     fHitCache("AliFMDHit",0),
97     fDigitCache("AliFMDDigit", 0),
98     fRawCache("AliFMDDigit", 0)
99 {
100   // Constructor 
101   // Parameters:
102   // @param name     Name of the folder. 
103   // @param useBoxes Whether to use boxes or Quads for the signals 
104   // @param old      Whether to enable reading old RCU data format
105
106   // increase reference count 
107   IncDenyDestroy();
108   
109   // Increase reference counts on palettes 
110   fHitPalette.IncRefCount();
111   fDigitPalette.IncRefCount();
112   fMultPalette.IncRefCount();
113   
114   
115   // Initialize the FMD geometry manager 
116   AliEveEventManager::AssertGeometry();
117   AliFMDGeometry* geom = AliFMDGeometry::Instance();
118   geom->Init();
119   geom->InitTransformations();
120
121   AliFMDParameters* pars = AliFMDParameters::Instance();
122   // pars->UseRcuTrailer(!old);
123   pars->UseCompleteHeader(old);
124   pars->SetSampleRate(4);
125   pars->Init(kFALSE, 0);
126
127   // Get shapes
128   TGeoShape* inner = static_cast<TGeoShape*>(gGeoManager->GetListOfShapes()
129                                              ->FindObject("FMDI_physical_sensor"));
130   if (!inner) throw TEveException("Shape of inner type sensors not found");
131   TGeoShape* outer = static_cast<TGeoShape*>(gGeoManager->GetListOfShapes()
132                                              ->FindObject("FMDO_physical_sensor"));
133   if (!outer) throw TEveException("Shape of outer type sensors not found");
134
135   // Emulate reference counting 
136   inner->SetUniqueID(1000);
137   outer->SetUniqueID(1000);
138   
139   // Loop over detectors 
140   for (UShort_t d = 1; d <= 3; d++) { 
141     AliFMDDetector*  detector = geom->GetDetector(d);
142     if (!detector) continue;
143     TEveElementList* ed       = new TEveElementList(Form(kDetector, 
144                                                          detector->GetId()));
145     AddElement(ed);
146     ed->IncDenyDestroy();
147     ed->SetUserData(detector);
148
149     // Loop over rings 
150     Char_t           rings[]  = { 'I', 'O', 0 };
151     Char_t*          pr       = &(rings[0]);
152     while (*pr) { 
153       AliFMDRing* ring = detector->GetRing(*pr);
154       pr++;
155       if (!ring) continue;
156       TEveElementList* er      = new TEveElementList(Form(kRing, 
157                                                           detector->GetId(),
158                                                           ring->GetId()));
159       ed->AddElement(er);
160       er->IncDenyDestroy();
161       er->SetUserData(ring);
162       
163       // UShort_t      nsec    = ring->GetNSectors();
164       // UShort_t      nstr    = ring->GetNStrips();
165       UShort_t         nmod    = ring->GetNModules();
166       // Loop over modules 
167       for (UShort_t m = 0; m < nmod; m++) { 
168         TEveGeoShape* em = new TEveGeoShape(Form(kModule, 
169                                                  detector->GetId(),
170                                                  ring->GetId(), 
171                                                  2*m, 2*m+1));
172         er->AddElement(em);
173         em->SetTransMatrix(*(detector->FindTransform(ring->GetId(), 2*m)));
174         em->SetShape(ring->GetId() == 'I' ? inner : outer);
175         em->SetMainColor(Color_t(kRed));
176         em->SetMainTransparency(32);
177         em->IncDenyDestroy();
178
179 #if 0
180         for (UShort_t s = 2*m; s < 2*m+2 && s < nsec; s++) { 
181           TEveDigitSet* eb = MakeDigitSet(Form(kSector,
182                                                detector->GetId(), 
183                                                ring->GetId(), s), nstr);
184           em->AddElement(eb);
185           eb->SetEmitSignals(kFALSE);
186           eb->SetPickable(kTRUE);
187           // eb->SetOwnIds(kTRUE);
188         } // for (UShort_t s ...)
189 #endif
190       }  // for (UShort_t m ...)
191     }  // while (pr)
192   }  // for (UShort_t d ...)
193 }
194
195 //____________________________________________________________________
196 AliEveFMDLoader::~AliEveFMDLoader()
197 {
198   // Destructor
199   AliWarning("AliEveFMDLoader being destroyed!");
200 }
201
202 //____________________________________________________________________
203 Int_t
204 AliEveFMDLoader::RemoveFromListTrees(TEveElement* el)
205 {
206   // Called when the element should be removed from the list. We
207   // overload this to allow clearing of signals.
208   // Parameters:
209   // @param el Tree to remove from.
210
211   // Since we're most likely setting up for a new event, we clear all
212   // signals here - a little tricky, but it works(tm). 
213   ClearDigitSets("All");
214   
215   // Do normal TEveElement::RemoveElement
216   return TEveElementList::RemoveFromListTrees(el);
217 }
218 //____________________________________________________________________
219 void
220 AliEveFMDLoader::RemoveParent(TEveElement* el)
221 {
222   // Called when the element should be removed from the list. We
223   // overload this to allow clearing of signals.
224   // Parameters:
225   // @param el Parent to remove from.
226   TEveElementList::RemoveParent(el);
227 }
228
229 //____________________________________________________________________
230 TEveDigitSet* 
231 AliEveFMDLoader::MakeDigitSet(const char* name, UShort_t nstr)
232 {
233   // Make a digit set.  The type of digit set depends on the setting
234   // of fUseBoxDigits.  If this is true, we return a TEveBoxSet,
235   // otherwise a TEveQuadSet   
236   // Parameters: 
237   //    name    The name 
238   //    nstr    The number of strips 
239   // Return 
240   //    newly allocated digit set
241   TEveDigitSet* ret = 0;
242   if (fUseBoxDigits) { 
243     TEveBoxSet* boxes = new TEveBoxSet(name);
244     // boxes->Reset(TEveBoxSet::kBT_AABox, kFALSE, nstr);
245     boxes->Reset(TEveBoxSet::kBT_FreeBox, kFALSE, nstr);
246     ret = boxes;
247   }
248   else { 
249     TEveQuadSet* quads = new TEveQuadSet(name);
250     quads->Reset(TEveQuadSet::kQT_RectangleXY, kFALSE, nstr);
251     ret = quads;
252   }
253   return ret;
254 }
255
256 //____________________________________________________________________
257 void
258 AliEveFMDLoader::ClearDigitSets(const char* type)
259 {
260   // Clear signals of some type.   
261   // Parameters:
262   // @param type Type of signals to clear 
263   //     Type can be one of 
264   //     - All    All signals 
265   //     - Hits   Hits 
266   //     - Digits Digits 
267   //     - Raw    Raw 
268   //     - ESD    ESD 
269   TString stype(type);
270
271   for (TEveElement::List_i di = BeginChildren(); 
272        di != EndChildren(); ++di) { 
273     for (TEveElement::List_i ri = (*di)->BeginChildren(); 
274          ri != (*di)->EndChildren(); ++ri) { 
275       for (TEveElement::List_i mi = (*ri)->BeginChildren();
276            mi != (*ri)->EndChildren(); ++mi) { 
277         if (stype == "All") {
278           (*mi)->RemoveElements();
279           continue;
280         }
281         for (TEveElement::List_i si = (*mi)->BeginChildren(); 
282              si != (*mi)->EndChildren(); ++si) { 
283           TEveDigitSet* signals = static_cast<TEveDigitSet*>((*si));
284           if (!signals) continue;
285           TString s(signals->GetName());
286           if (!s.Contains(type)) continue;
287           (*mi)->RemoveElement(signals);
288         }
289       }
290     }
291   }
292 }
293
294 //____________________________________________________________________
295 TEveDigitSet*
296 AliEveFMDLoader::FindDigitSet(const char* t, UShort_t d, Char_t r, UShort_t s)
297 {
298   // Find a digit set corresponding to the passed parameters.  If it
299   // is not found, one is created
300   // Parameters:
301   // @param type   Type of data 
302   // @param d      Detector 
303   // @param r      Ring 
304   // @param s      Sector 
305   // @return a digit set
306   TEveElement* detector = FindChild(Form(kDetector, d));
307   if (!detector) { 
308     AliError(Form("Detector %s not found", Form(kDetector, d)));
309     return 0;
310   }
311   
312   TEveElement* ring = detector->FindChild(Form(kRing, d, r));
313   if (!ring) { 
314     AliError(Form("Ring %s not found", Form(kRing, d, r)));
315     return 0;
316   }
317   
318   Int_t mod = 2*(s/2);
319   TEveElement* module = ring->FindChild(Form(kModule, d, r, mod, mod+1));
320   if (!module) { 
321     AliError(Form("Module %s not found", Form(kModule, d, r, s, s+1)));
322     return 0;
323   }
324   
325   TEveElement*  sector = module->FindChild(Form(kSector, d, r, s, t));
326   TEveDigitSet* signal = static_cast<TEveDigitSet*>(sector);
327   if (!sector) { 
328     AliFMDRing* rng = AliFMDGeometry::Instance()->GetRing(r);
329     signal = MakeDigitSet(Form(kSector, d, r, s, t), rng->GetNStrips());
330     module->AddElement(signal);
331     signal->SetEmitSignals(kFALSE);
332     signal->SetPickable(kTRUE);
333     TString st(t);
334     if      (t == kHits)   signal->SetPalette(&fHitPalette);
335     else if (t == kDigits) signal->SetPalette(&fDigitPalette);
336     else if (t == kRaw)    signal->SetPalette(&fDigitPalette);
337     else if (t == kESD) {
338       signal->SetPalette(&fMultPalette);
339       signal->SetOwnIds(kTRUE);
340     }
341   }
342   return signal;
343 }
344
345 //____________________________________________________________________
346 void
347 AliEveFMDLoader::AddSignal(const char* t, 
348                            UShort_t det, Char_t rng, UShort_t sec, 
349                            UShort_t str, Float_t signal, Float_t min, 
350                            Float_t  max, TObject* ref)
351 {
352   // Add a signal to a digit set
353   // Parameters:
354   // @param type   Type of data 
355   // @param det    Detector 
356   // @param rng    Ring 
357   // @param sec    Sector 
358   // @param str    Strip
359   // @param signal Signal value 
360   // @param min    Minimum of this kind of signal 
361   // @param max    Maximum of this kind of signal 
362   // @param ref    Reference object 
363   AliFMDGeometry* geom = AliFMDGeometry::Instance();
364   Double_t x, y, z;
365   geom->Detector2XYZ(det, rng, sec, str, x, y, z);
366   AddSignal(t, det, rng, sec, str, x, y, z, signal, min, max, ref);
367 }
368
369 //____________________________________________________________________
370 void
371 AliEveFMDLoader::AddSignal(const char* t, 
372                            UShort_t det, Char_t rng, UShort_t sec, 
373                            UShort_t str, Double_t x, Double_t y, Double_t z, 
374                            Float_t signal, Float_t min, Float_t  max, 
375                            TObject* ref)
376 {
377   // Add a signal to a digit set, with known (x,y,z) coordinates
378   // (this is for hits)
379   // Parameters:
380   // @param type   Type of data 
381   // @param det    Detector 
382   // @param rng    Ring 
383   // @param sec    Sector 
384   // @param str    Strip
385   // @param x      X coordinate  
386   // @param y      Y coordinate  
387   // @param z      Z coordinate  
388   // @param signal Signal value 
389   // @param min    Minimum of this kind of signal 
390   // @param max    Maximum of this kind of signal 
391   // @param ref    Reference object 
392   AliFMDGeometry* geom = AliFMDGeometry::Instance();
393   AliFMDRing*     ring = geom->GetRing(rng);
394   if (!ring) return;
395
396   TEveDigitSet* signals = FindDigitSet(t, det, rng, sec);
397   if (!signals) { 
398     AliWarning(Form("No signal (%s) found for FMD%d%c[%02d,%03d]", 
399                     t, det, rng, sec, str));
400     return;
401   }
402   
403   Float_t  scaled = TMath::Min((signal - min) / (max - min) * 10., 10.);
404   Double_t w      = 2*ring->GetPitch();
405   Int_t    value  = int(TMath::Nint(signal));
406   AliDebug(1, Form("New signal at FMD%d%c[%02d,%03d]=%f (v=%d, s=%f)", 
407                    det, rng, sec, str, signal, value, scaled));
408   AddDigit(signals, x, y, z, w, scaled, value, ref);
409 }
410
411 //____________________________________________________________________
412 void
413 AliEveFMDLoader::AddDigit(TEveDigitSet* signals, 
414                           Double_t x, Double_t y, Double_t z, 
415                           Double_t w, Float_t scaled, Int_t value, 
416                           TObject* ref)
417 {
418   // Add a digit to a digit set. 
419   // Parameters:
420   // @param signals Digit set. 
421   // @param x      X coordinate  
422   // @param y      Y coordinate  
423   // @param z      Z coordinate  
424   // @param w      strip pitch 
425   // @param scaled Scaled value 
426   // @param value  Signal value 
427   // @param ref    Reference object
428   if (fUseBoxDigits) { 
429     TEveBoxSet* boxes = static_cast<TEveBoxSet*>(signals);
430     Float_t zc   = (z > 0 ? -1 : 1) * scaled + z;
431     Float_t vs[] = { -w, -5*w, zc-scaled,   // Lower back  left
432                      +w, -5*w, zc-scaled,   // Lower back  right
433                      +w, +5*w, zc-scaled,   // Lower front right
434                      -w, +5*w, zc-scaled,   // Lower front left 
435                      -w, -5*w, zc+scaled,   // Upper back  left
436                      +w, -5*w, zc+scaled,   // Upper back  right
437                      +w, +5*w, zc+scaled,   // Upper front right
438                      -w, +5*w, zc+scaled }; // Upper front left
439     Float_t ang  = TMath::ATan2(y,x);
440     for (size_t i = 0; i < 8; i++) { 
441       Float_t bx = vs[3*i+0];
442       Float_t by = vs[3*i+1];
443       Float_t ca = TMath::Cos(ang);
444       Float_t sa = TMath::Sin(ang);
445       vs[3*i+0]  = bx * ca - by * sa + x;
446       vs[3*i+1]  = bx * sa + by * ca + y;
447     }
448     // boxes->AddBox(x, y, (z > 0 ? -scaled : 0) + z , 5*w, w, scaled);
449     boxes->AddBox(vs);
450     boxes->DigitValue(value);
451     if (ref) boxes->DigitId(ref);
452   }
453   else { 
454     TEveQuadSet* quads = static_cast<TEveQuadSet*>(signals);
455     quads->AddQuad(x,y,z,w,w);
456     quads->QuadValue(value);
457     if (ref) quads->QuadId(ref);
458   }
459 }
460
461
462 //____________________________________________________________________
463 void
464 AliEveFMDLoader::CheckAdd()
465 {
466   // check if we shoul re-add ourselves to the current event node 
467   TEveElement* event = gEve->GetCurrentEvent();
468   if (event && event->FindChild(GetName())) return;
469   gEve->AddElement(this);
470 }
471
472
473 //____________________________________________________________________
474 void
475 AliEveFMDLoader::LoadHits()
476 {
477   // Load and display hits
478   ClearDigitSets(kHits);
479
480   AliRunLoader* rl =  AliEveEventManager::AssertRunLoader();
481   if (!rl) { 
482     AliError("No run loader");
483     return;
484   }
485   
486   rl->LoadHits("FMD");
487   TTree* ht = rl->GetTreeH("FMD", false);
488   if (!ht) { 
489     AliError("No FMD tree");
490     return;
491   }
492   
493   TClonesArray* hits = &fHitCache;
494   fHitCache.Clear();
495   ht->SetBranchAddress("FMD", &hits);
496
497   Float_t min = fHitPalette.GetMinVal();
498   Float_t max = fHitPalette.GetMaxVal();
499
500   Int_t nTracks = ht->GetEntriesFast();
501   for (Int_t i = 0; i < nTracks; i++) {
502     Int_t hitRead  = ht->GetEntry(i);
503     if (hitRead <= 0) continue;
504
505     Int_t nHit = hits->GetEntriesFast();
506     if (nHit <= 0) continue;
507   
508     for (Int_t j = 0; j < nHit; j++) {
509       AliFMDHit* hit = static_cast<AliFMDHit*>(hits->At(j));
510       if (!hit) continue;
511   
512       AddSignal(kHits, 
513                 hit->Detector(), hit->Ring(), hit->Sector(), hit->Strip(),
514                 hit->X(), hit->Y(), hit->Z(), int(hit->Edep()*1000), 
515                 min, max, hit);
516     }
517   }
518   CheckAdd();
519 }
520
521 //____________________________________________________________________
522 void
523 AliEveFMDLoader::DoLoadDigits(const char* t, TClonesArray* digits)
524 {
525   // Do the actual display of digits 
526   // Parameters:
527   // @param type What to show 
528   // @param digits The digits
529   Float_t min = fDigitPalette.GetMinVal();
530   Float_t max = fDigitPalette.GetMaxVal();
531
532   Int_t n = digits->GetEntriesFast();
533   for (Int_t i = 0; i < n; i++) { 
534     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
535     if (!digit) return;
536     AddSignal(t, digit->Detector(), digit->Ring(), digit->Sector(),
537               digit->Strip(), digit->Counts(), min, max, digit);
538   }
539   CheckAdd();
540 }
541
542 //____________________________________________________________________
543 void
544 AliEveFMDLoader::LoadDigits()
545 {
546   // Load and display simulated digits 
547   ClearDigitSets(kDigits);
548
549   AliRunLoader* rl =  AliEveEventManager::AssertRunLoader();
550   if (!rl) { 
551     AliError("No run-loader");
552     return;
553   }
554   
555   rl->LoadDigits("FMD");
556   TTree* dt = rl->GetTreeD("FMD", false);
557   if (!dt) { 
558     AliError("No FMD tree");
559     return;
560   }
561   
562   TClonesArray* digits = &fDigitCache;
563   fDigitCache.Clear();
564   dt->SetBranchAddress("FMD", &digits);
565
566   Int_t read = dt->GetEntry(0);
567   if (read <= 0) { 
568     AliWarning("Nothing read");
569     return;
570   }
571   DoLoadDigits(kDigits, digits);
572 }
573
574
575 //____________________________________________________________________
576 void
577 AliEveFMDLoader::LoadRaw()
578 {
579   // Load and display raw digits 
580   ClearDigitSets(kRaw);
581
582   AliRawReader* rr =  AliEveEventManager::AssertRawReader();
583   if (!rr) { 
584     AliError("No raw-reader");
585     return;
586   }
587   rr->Reset();
588   std::cout<<"Now in event # " << *(rr->GetEventId()) << std::endl;
589   AliFMDRawReader* fr = new AliFMDRawReader(rr, 0);
590   TClonesArray* digits = &fRawCache;
591   fRawCache.Clear();
592   
593   fr->ReadAdcs(digits);
594   
595   DoLoadDigits(kRaw, digits);
596 }
597
598 //____________________________________________________________________
599 void
600 AliEveFMDLoader::LoadESD()
601 {
602   // Load and display ESD information 
603   ClearDigitSets(kESD);
604
605   AliESDEvent* esd =  AliEveEventManager::AssertESD();
606   if (!esd) { 
607     AliError("No ESD");
608     return;
609   }
610
611   AliESDFMD* fmd = esd->GetFMDData();
612   if (!fmd) { 
613     AliError("No FMD ESD data");
614     return;
615   }
616
617   Float_t min = fMultPalette.GetMinVal();
618   Float_t max = fMultPalette.GetMaxVal();
619   
620   for (UShort_t det = 1; det <= 3; det++) {
621     Char_t rings[] = { 'I', (det == 1 ? '\0' : 'O'), '\0' };
622     for (Char_t* rng = rings; *rng != '\0'; rng++) {
623       UShort_t nsec = (*rng == 'I' ?  20 :  40);
624       UShort_t nstr = (*rng == 'I' ? 512 : 256);
625       for (UShort_t sec = 0; sec < nsec; sec++) {
626         for (UShort_t str = 0; str < nstr; str++) {
627           Float_t mult = fmd->Multiplicity(det,*rng,sec,str);
628           if (mult == AliESDFMD::kInvalidMult) continue;
629           Float_t eta  = fmd->Eta(det,*rng,sec,str);
630           AddSignal(kESD, det, *rng, sec, str, mult, min, max, 
631                     new TNamed(Form("FMD%d%c[%02d,%03d]", det, *rng, sec, str), 
632                                Form("Mch=%f, eta=%f", mult, eta)));
633         }
634       }
635     }
636   }
637   CheckAdd();
638 }
639
640
641         
642       
643
644   
645 //____________________________________________________________________
646 //
647 // EOF
648 //