00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "MagneticField/GeomBuilder/src/MagGeoBuilderFromDDD.h"
00010 #include "MagneticField/GeomBuilder/src/volumeHandle.h"
00011 #include "MagneticField/GeomBuilder/src/bSlab.h"
00012 #include "MagneticField/GeomBuilder/src/bRod.h"
00013 #include "MagneticField/GeomBuilder/src/bSector.h"
00014 #include "MagneticField/GeomBuilder/src/bLayer.h"
00015 #include "MagneticField/GeomBuilder/src/eSector.h"
00016 #include "MagneticField/GeomBuilder/src/eLayer.h"
00017 #include "MagneticField/GeomBuilder/src/FakeInterpolator.h"
00018
00019 #include "MagneticField/Layers/interface/MagBLayer.h"
00020 #include "MagneticField/Layers/interface/MagESector.h"
00021
00022 #include "FWCore/ParameterSet/interface/FileInPath.h"
00023
00024 #include "DetectorDescription/Core/interface/DDCompactView.h"
00025 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00026 #include "DetectorDescription/Core/interface/DDFilter.h"
00027
00028 #include "Utilities/BinningTools/interface/ClusterizingHistogram.h"
00029 #include "CLHEP/Units/SystemOfUnits.h"
00030
00031 #include "MagneticField/Interpolation/interface/MagProviderInterpol.h"
00032 #include "MagneticField/Interpolation/interface/MFGridFactory.h"
00033 #include "MagneticField/Interpolation/interface/MFGrid3D.h"
00034
00035 #include "MagneticField/VolumeGeometry/interface/MagVolume6Faces.h"
00036 #include "MagneticField/VolumeGeometry/interface/MagExceptions.h"
00037 #include "MagneticField/Layers/interface/MagVerbosity.h"
00038
00039 #include "DataFormats/GeometryVector/interface/Pi.h"
00040
00041 #include "FWCore/Utilities/interface/Exception.h"
00042 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00043
00044 #include <string>
00045 #include <vector>
00046 #include <iostream>
00047 #include <sstream>
00048 #include <algorithm>
00049 #include <iterator>
00050 #include <map>
00051 #include <set>
00052 #include <iomanip>
00053 #include "Utilities/General/interface/precomputed_value_sort.h"
00054
00055
00056 bool MagGeoBuilderFromDDD::debug;
00057
00058 using namespace std;
00059
00060 MagGeoBuilderFromDDD::MagGeoBuilderFromDDD(string version_, bool debug_, bool overrideMasterSector_) :
00061 version (version_),
00062 overrideMasterSector(overrideMasterSector_)
00063 {
00064 debug = debug_;
00065 if (debug) cout << "Constructing a MagGeoBuilderFromDDD" <<endl;
00066 }
00067
00068 MagGeoBuilderFromDDD::~MagGeoBuilderFromDDD(){
00069 for (handles::const_iterator i=bVolumes.begin();
00070 i!=bVolumes.end(); ++i){
00071 delete (*i);
00072 }
00073
00074 for (handles::const_iterator i=eVolumes.begin();
00075 i!=eVolumes.end(); ++i){
00076 delete (*i);
00077 }
00078 }
00079
00080
00081 void MagGeoBuilderFromDDD::summary(handles & volumes){
00082
00083 int ivolumes = volumes.size();
00084 int isurfaces = ivolumes*6;
00085 int iassigned = 0;
00086 int iunique = 0;
00087 int iref_ass = 0;
00088 int iref_nass = 0;
00089
00090 set<const void *> ptrs;
00091
00092 handles::const_iterator first = volumes.begin();
00093 handles::const_iterator last = volumes.end();
00094
00095 for (handles::const_iterator i=first; i!=last; ++i){
00096 if (int((*i)->shape())>4) continue;
00097 for (int side = 0; side < 6; ++side) {
00098 int references = (*i)->references(side);
00099 if ((*i)->isPlaneMatched(side)) {
00100 ++iassigned;
00101 bool firstOcc = (ptrs.insert(&((*i)->surface(side)))).second;
00102 if (firstOcc) iref_ass+=references;
00103 if (references<2){
00104 cout << "*** Only 1 ref, vol: " << (*i)->name << " # "
00105 << (*i)->copyno << " side: " << side << endl;
00106 }
00107 } else {
00108 iref_nass+=references;
00109 if (references>1){
00110 cout << "*** Ref_nass >1 " <<endl;
00111 }
00112 }
00113 }
00114 }
00115 iunique = ptrs.size();
00116
00117 cout << " volumes " << ivolumes << endl
00118 << " surfaces " << isurfaces << endl
00119 << " assigned " << iassigned << endl
00120 << " unique " << iunique << endl
00121 << " iref_ass " << iref_ass << endl
00122 << " iref_nass " << iref_nass << endl;
00123 }
00124
00125
00126 void MagGeoBuilderFromDDD::build(const DDCompactView & cpva)
00127 {
00128
00129 DDExpandedView fv(cpva);
00130
00131 if (debug) cout << "**********************************************************" <<endl;
00132
00133
00134 map<string, MagProviderInterpol*> bInterpolators;
00135 map<string, MagProviderInterpol*> eInterpolators;
00136
00137
00138 int bVolCount = 0;
00139 int eVolCount = 0;
00140
00141 if (fv.logicalPart().name().name()!="MAGF") {
00142 std::string topNodeName(fv.logicalPart().name().name());
00143
00144
00145 bool doSubDets = fv.firstChild();
00146
00147 bool go=true;
00148 while(go&& doSubDets) {
00149 if (fv.logicalPart().name().name()=="MAGF")
00150 break;
00151 else
00152 go = fv.nextSibling();
00153 }
00154 if (!go) {
00155 throw cms::Exception("NoMAGFinDDD")<<" Neither he top node, nor any child node of the DDCompactView is \"MAGF\" but the top node is instead \""<<topNodeName<<"\"";
00156 }
00157 }
00158
00159 if (debug) { cout << endl << "*** MAGF: " << fv.geoHistory() << endl
00160 << "translation: " << fv.translation() << endl
00161 << " rotation: " << fv.rotation() << endl;
00162 }
00163
00164 bool doSubDets = fv.firstChild();
00165 while (doSubDets){
00166
00167 string name = fv.logicalPart().name().name();
00168 if (debug) cout << endl << "Name: " << name << endl
00169 << " " << fv.geoHistory() <<endl;
00170
00171
00172
00173
00174
00175 if (name.substr(2,2)=="ZP") {
00176 doSubDets = fv.nextSibling();
00177 continue;
00178 }
00179
00180
00181
00182
00183
00184
00185
00186 bool expand = false;
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204 volumeHandle* v = new volumeHandle(fv, expand);
00205
00206
00207
00208
00209 if (overrideMasterSector && v->masterSector!=1) {
00210 v->masterSector=1;
00211 std::string::size_type ibeg, iend;
00212 ibeg = (v->magFile).rfind('/');
00213 iend = (v->magFile).size();
00214 v->magFile = (v->magFile).substr(ibeg+1,iend-ibeg-1);
00215 }
00216
00217
00218 float Z = v->center().z();
00219 float R = v->center().perp();
00220
00221
00222
00223
00224
00225
00226
00227 if ((fabs(Z)<647. || (R>350. && fabs(Z)<662.)) &&
00228 !(fabs(Z)>480 && R<172)
00229
00230
00231
00232 ) {
00233 if (debug) cout << " (Barrel)" <<endl;
00234 bVolumes.push_back(v);
00235
00236
00237
00238
00239
00240
00241
00242
00243 if (v->copyno==v->masterSector) {
00244 buildInterpolator(v, bInterpolators);
00245 ++bVolCount;
00246 }
00247 } else {
00248 if (debug) cout << " (Endcaps)" <<endl;
00249 eVolumes.push_back(v);
00250 if (v->copyno==v->masterSector) {
00251 buildInterpolator(v, eInterpolators);
00252 ++eVolCount;
00253 }
00254 }
00255
00256 doSubDets = fv.nextSibling();
00257 }
00258
00259 if (debug) {
00260 cout << "Number of volumes (barrel): " << bVolumes.size() <<endl
00261 << "Number of volumes (endcap): " << eVolumes.size() <<endl;
00262 cout << "**********************************************************" <<endl;
00263 }
00264
00265
00266
00267
00268
00269
00270
00271 if (debug) {
00272 cout << "-----------------------" << endl;
00273 cout << "SUMMARY: Barrel " << endl;
00274 summary(bVolumes);
00275
00276 cout << endl << "SUMMARY: Endcaps " << endl;
00277 summary(eVolumes);
00278 cout << "-----------------------" << endl;
00279 }
00280
00281
00282
00283
00284
00285 vector<bLayer> layers;
00286 precomputed_value_sort(bVolumes.begin(), bVolumes.end(), ExtractRN());
00287
00288
00289 const float resolution = 1.;
00290 float rmin = bVolumes.front()->RN()-resolution;
00291 float rmax = bVolumes.back()->RN()+resolution;
00292 ClusterizingHistogram hisR( int((rmax-rmin)/resolution) + 1, rmin, rmax);
00293
00294 if (debug) cout << " R layers: " << rmin << " " << rmax << endl;
00295
00296 handles::const_iterator first = bVolumes.begin();
00297 handles::const_iterator last = bVolumes.end();
00298
00299 for (handles::const_iterator i=first; i!=last; ++i){
00300 hisR.fill((*i)->RN());
00301 }
00302 vector<float> rClust = hisR.clusterize(resolution);
00303
00304 handles::const_iterator ringStart = first;
00305 handles::const_iterator separ = first;
00306
00307 for (unsigned int i=0; i<rClust.size() - 1; ++i) {
00308 if (debug) cout << " Layer at RN = " << rClust[i];
00309 float rSepar = (rClust[i] + rClust[i+1])/2.f;
00310 while ((*separ)->RN() < rSepar) ++separ;
00311
00312 bLayer thislayer(ringStart, separ);
00313 layers.push_back(thislayer);
00314 ringStart = separ;
00315 }
00316 {
00317 if (debug) cout << " Layer at RN = " << rClust.back();
00318 bLayer thislayer(separ, last);
00319 layers.push_back(thislayer);
00320 }
00321
00322 if (debug) cout << "Barrel: Found " << rClust.size() << " clusters in R, "
00323 << layers.size() << " layers " << endl << endl;
00324
00325
00326
00327
00328
00329 vector<eSector> sectors;
00330 precomputed_value_sort(eVolumes.begin(), eVolumes.end(), ExtractPhi());
00331
00332
00333
00334 for (int i = 0; i<12; ++i) {
00335 int offset = eVolumes.size()/12;
00336
00337 if (debug) cout << " Sector at phi = "
00338 << (*(eVolumes.begin()+((i)*offset)))->center().phi()
00339 << endl;
00340 sectors.push_back(eSector(eVolumes.begin()+((i)*offset),
00341 eVolumes.begin()+((i+1)*offset)));
00342 }
00343
00344 if (debug) cout << "Endcap: Found "
00345 << sectors.size() << " sectors " << endl;
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387 buildMagVolumes(bVolumes, bInterpolators);
00388
00389
00390 for (vector<bLayer>::const_iterator ilay = layers.begin();
00391 ilay!= layers.end(); ++ilay) {
00392 mBLayers.push_back((*ilay).buildMagBLayer());
00393 }
00394
00395 if (debug) {
00396 cout << "*** BARREL ********************************************" << endl
00397 << "Number of different volumes = " << bVolCount << endl
00398 << "Number of interpolators built = " << bInterpolators.size() << endl
00399 << "Number of MagBLayers built = " << mBLayers.size() << endl;
00400
00401 testInside(bVolumes);
00402 }
00403
00404
00405
00406 buildMagVolumes(eVolumes, eInterpolators);
00407
00408
00409 for (vector<eSector>::const_iterator isec = sectors.begin();
00410 isec!= sectors.end(); ++isec) {
00411 mESectors.push_back((*isec).buildMagESector());
00412 }
00413
00414 if (debug) {
00415 cout << "*** ENDCAP ********************************************" << endl
00416 << "Number of different volumes = " << eVolCount << endl
00417 << "Number of interpolators built = " << eInterpolators.size() << endl
00418 << "Number of MagESector built = " << mESectors.size() << endl;
00419
00420 testInside(eVolumes);
00421 }
00422 }
00423
00424
00425 void MagGeoBuilderFromDDD::buildMagVolumes(const handles & volumes, map<string, MagProviderInterpol*> & interpolators) {
00426
00427 for (handles::const_iterator vol=volumes.begin(); vol!=volumes.end();
00428 ++vol){
00429 const MagProviderInterpol* mp = 0;
00430 if (interpolators.find((*vol)->magFile)!=interpolators.end()) {
00431 mp = interpolators[(*vol)->magFile];
00432 } else {
00433 cout << "No interpolator found for file " << (*vol)->magFile
00434 << " vol: " << (*vol)->name << endl;
00435 cout << interpolators.size() <<endl;
00436 }
00437
00438
00439
00440
00441 int volNum;
00442 string name = (*vol)->name;
00443 name.erase(0,name.rfind('_')+1);
00444 stringstream str;
00445 str << name;
00446 str >> volNum;
00447
00448
00449
00450
00451
00452
00453 int key = volNum*100+(*vol)->copyno;
00454 map<int, double>::const_iterator isf = theScalingFactors.find(key);
00455 if (isf == theScalingFactors.end()) {
00456 key = volNum*100;
00457 isf = theScalingFactors.find(key);
00458 }
00459
00460 double sf = 1.;
00461 if (isf != theScalingFactors.end()) {
00462 sf = (*isf).second;
00463
00464 edm::LogInfo("MagneticField|VolumeBasedMagneticFieldESProducer") << "Applying scaling factor " << sf << " to "<< (*vol)->name << "["<< (*vol)->copyno << "] (key:" << key << ")" << endl;
00465 }
00466
00467 const GloballyPositioned<float> * gpos = (*vol)->placement();
00468 (*vol)->magVolume = new MagVolume6Faces(gpos->position(),
00469 gpos->rotation(),
00470 (*vol)->shape(),
00471 (*vol)->sides(),
00472 mp, sf);
00473
00474 (*vol)->magVolume->setIsIron((*vol)->isIron());
00475
00476
00477 (*vol)->magVolume->name = (*vol)->name;
00478 (*vol)->magVolume->copyno = (*vol)->copyno;
00479 }
00480 }
00481
00482
00483 void MagGeoBuilderFromDDD::buildInterpolator(const volumeHandle * vol, map<string, MagProviderInterpol*> & interpolators){
00484
00485
00486
00487
00488
00489
00490
00491
00492 if (version=="grid_85l_030919" && vol->center().z()>0) return;
00493
00494
00495 double masterSectorPhi = (vol->masterSector-1)*Geom::pi()/6.;
00496
00497
00498
00499 if (version=="grid_85l_030919") {
00500 masterSectorPhi=Geom::pi()/2.;
00501 }
00502
00503 if (debug) {
00504 cout << "Building interpolator from "
00505 << vol->name << " copyno " << vol->copyno
00506 << " at " << vol->center()
00507 << " phi: " << vol->center().phi()
00508 << " file: " << vol->magFile
00509 << endl;
00510
00511 if ( fabs(vol->center().phi() - masterSectorPhi) > Geom::pi()/9.) {
00512 cout << "***WARNING wrong sector? " << endl;
00513 }
00514 }
00515
00516 if (version == "fake") {
00517 interpolators[vol->magFile] = new magneticfield::FakeInterpolator();
00518 return;
00519 }
00520
00521 string fullPath;
00522
00523 try {
00524 edm::FileInPath mydata("MagneticField/Interpolation/data/"+version+"/"+vol->magFile);
00525 fullPath = mydata.fullPath();
00526 } catch (edm::Exception& exc) {
00527 cerr << "MagGeoBuilderFromDDD: exception in reading table; " << exc.what() << endl;
00528 if (!debug) throw;
00529 return;
00530 }
00531
00532
00533 try{
00534 if (vol->toExpand()){
00535
00536
00537
00538 } else {
00539
00540
00541
00542 GloballyPositioned<float> rf = *(vol->placement());
00543
00544 if (vol->masterSector != 1) {
00545 typedef Basic3DVector<float> Vector;
00546
00547 GloballyPositioned<float>::RotationType rot(Vector(0,0,1), -masterSectorPhi);
00548 Vector vpos(vol->placement()->position());
00549
00550
00551 rf = GloballyPositioned<float>(GloballyPositioned<float>::PositionType(rot.multiplyInverse(vpos)), vol->placement()->rotation()*rot);
00552 }
00553
00554 interpolators[vol->magFile] =
00555 MFGridFactory::build( fullPath, rf);
00556 }
00557 } catch (MagException& exc) {
00558 cout << exc.what() << endl;
00559 interpolators.erase(vol->magFile);
00560 if (!debug) throw;
00561 return;
00562 }
00563
00564
00565 if (debug) {
00566
00567 const MagVolume6Faces tempVolume(vol->placement()->position(),
00568 vol->placement()->rotation(),
00569 vol->shape(),
00570 vol->sides(),
00571 interpolators[vol->magFile]);
00572
00573 const MFGrid3D* grid = dynamic_cast<const MFGrid3D*>(interpolators[vol->magFile]);
00574 if (grid!=0) {
00575
00576 vector<int> sizes = grid->dimensions();
00577 cout << "Grid has " << sizes.size() << " dimensions "
00578 << " number of nodes is " << sizes[0] << " " << sizes[1]
00579 << " " << sizes[2] << endl;
00580
00581 const double tolerance = 0.03;
00582
00583
00584 int dumpCount = 0;
00585 for (int j=0; j < sizes[1]; j++) {
00586 for (int k=0; k < sizes[2]; k++) {
00587 for (int i=0; i < sizes[0]; i++) {
00588 MFGrid::LocalPoint lp = grid->nodePosition( i, j, k);
00589 if (! tempVolume.inside(lp, tolerance)) {
00590 if (++dumpCount < 2) {
00591 MFGrid::GlobalPoint gp = tempVolume.toGlobal(lp);
00592 cout << "GRID ERROR: " << i << " " << j << " " << k
00593 << " local: " << lp
00594 << " global: " << gp
00595 << " R= " << gp.perp() << " phi=" << gp.phi() << endl;
00596 }
00597 }
00598 }
00599 }
00600 }
00601
00602 cout << vol->name << " : Number of grid points outside the MagVolume: " << dumpCount << "/" << sizes[0]*sizes[1]*sizes[2] << endl;
00603 }
00604 }
00605 }
00606
00607
00608
00609 void MagGeoBuilderFromDDD::testInside(handles & volumes) {
00610
00611 cout << "--------------------------------------------------" << endl;
00612 cout << " inside(center) test" << endl;
00613 for (handles::const_iterator vol=volumes.begin(); vol!=volumes.end();
00614 ++vol){
00615 for (handles::const_iterator i=volumes.begin(); i!=volumes.end();
00616 ++i){
00617 if ((*i)==(*vol)) continue;
00618
00619 if ((*i)->magVolume->inside((*vol)->center())) {
00620 cout << "*** ERROR: center of " << (*vol)->name << " is inside "
00621 << (*i)->name <<endl;
00622 }
00623 }
00624
00625 if ((*vol)->magVolume->inside((*vol)->center())) {
00626 cout << (*vol)->name << " OK " << endl;
00627 } else {
00628 cout << "*** ERROR: center of volume is not inside it, "
00629 << (*vol)->name << endl;
00630 }
00631 }
00632 cout << "--------------------------------------------------" << endl;
00633 }
00634
00635
00636 vector<MagBLayer*> MagGeoBuilderFromDDD::barrelLayers() const{
00637 return mBLayers;
00638 }
00639
00640 vector<MagESector*> MagGeoBuilderFromDDD::endcapSectors() const{
00641 return mESectors;
00642 }
00643
00644 vector<MagVolume6Faces*> MagGeoBuilderFromDDD::barrelVolumes() const{
00645 vector<MagVolume6Faces*> v;
00646 v.reserve(bVolumes.size());
00647 for (handles::const_iterator i=bVolumes.begin();
00648 i!=bVolumes.end(); ++i){
00649 v.push_back((*i)->magVolume);
00650 }
00651 return v;
00652 }
00653
00654 vector<MagVolume6Faces*> MagGeoBuilderFromDDD::endcapVolumes() const{
00655 vector<MagVolume6Faces*> v;
00656 v.reserve(eVolumes.size());
00657 for (handles::const_iterator i=eVolumes.begin();
00658 i!=eVolumes.end(); ++i){
00659 v.push_back((*i)->magVolume);
00660 }
00661 return v;
00662 }
00663
00664
00665 float MagGeoBuilderFromDDD::maxR() const{
00666
00667
00668
00669 if (version=="grid_85l_030919") return 1000.;
00670 else return 900.;
00671 }
00672
00673 float MagGeoBuilderFromDDD::maxZ() const{
00674
00675 return 1600.;
00676 }
00677
00678
00679 void MagGeoBuilderFromDDD::setScaling(std::vector<int> keys,
00680 std::vector<double> values)
00681 {
00682 if (keys.size() != values.size()) {
00683 throw cms::Exception("InvalidParameter") << "Invalid field scaling parameters 'scalingVolumes' and 'scalingFactors' ";
00684 }
00685 for (unsigned int i=0; i<keys.size(); ++i) {
00686 theScalingFactors[keys[i]] = values[i];
00687 }
00688 }
00689
00690
00691