Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009 #include "MagneticField/VolumeBasedEngine/interface/MagGeometry.h"
00010 #include "MagneticField/VolumeGeometry/interface/MagVolume.h"
00011 #include "MagneticField/VolumeGeometry/interface/MagVolume6Faces.h"
00012 #include "MagneticField/Layers/interface/MagBLayer.h"
00013 #include "MagneticField/Layers/interface/MagESector.h"
00014
00015 #include "Utilities/BinningTools/interface/PeriodicBinFinderInPhi.h"
00016
00017 #include <FWCore/ParameterSet/interface/ParameterSet.h>
00018 #include "FWCore/Utilities/interface/isFinite.h"
00019
00020 #include "MagneticField/Layers/interface/MagVerbosity.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022
00023 using namespace std;
00024 using namespace edm;
00025
00026 MagGeometry::MagGeometry(const edm::ParameterSet& config, std::vector<MagBLayer *> tbl,
00027 std::vector<MagESector *> tes,
00028 std::vector<MagVolume6Faces*> tbv,
00029 std::vector<MagVolume6Faces*> tev) :
00030 lastVolume(0), theBLayers(tbl), theESectors(tes), theBVolumes(tbv), theEVolumes(tev), geometryVersion(0)
00031 {
00032
00033 cacheLastVolume = config.getUntrackedParameter<bool>("cacheLastVolume", true);
00034 geometryVersion = config.getParameter<int>("geometryVersion");
00035
00036 vector<double> rBorders;
00037
00038 for (vector<MagBLayer *>::const_iterator ilay = theBLayers.begin();
00039 ilay != theBLayers.end(); ++ilay) {
00040 if (verbose::debugOut) cout << " Barrel layer at " << (*ilay)->minR() <<endl;
00041
00042 rBorders.push_back((*ilay)->minR());
00043 }
00044
00045 theBarrelBinFinder = new MagBinFinders::GeneralBinFinderInR<double>(rBorders);
00046
00047 if (verbose::debugOut) {
00048 for (vector<MagESector *>::const_iterator isec = theESectors.begin();
00049 isec != theESectors.end(); ++isec) {
00050 cout << " Endcap sector at " << (*isec)->minPhi() << endl;
00051 }
00052 }
00053
00054
00055
00056 int nEBins = theESectors.size();
00057 theEndcapBinFinder = new PeriodicBinFinderInPhi<float>(theESectors.front()->minPhi()+Geom::pi()/nEBins, nEBins);
00058
00059 }
00060
00061 MagGeometry::~MagGeometry(){
00062 delete theBarrelBinFinder;
00063 delete theEndcapBinFinder;
00064
00065 for (vector<MagBLayer *>::const_iterator ilay = theBLayers.begin();
00066 ilay != theBLayers.end(); ++ilay) {
00067 delete (*ilay);
00068 }
00069
00070 for (vector<MagESector *>::const_iterator ilay = theESectors.begin();
00071 ilay != theESectors.end(); ++ilay) {
00072 delete (*ilay);
00073 }
00074 }
00075
00076
00077
00078 GlobalVector MagGeometry::fieldInTesla(const GlobalPoint & gp) const {
00079 MagVolume * v = 0;
00080
00081
00082 v = findVolume(gp);
00083 if (v!=0) {
00084 return v->fieldInTesla(gp);
00085 }
00086
00087
00088
00089 if (edm::isNotFinite(gp.mag())) {
00090 LogWarning("InvalidInput") << "Input value invalid (not a number): " << gp << endl;
00091
00092 } else {
00093 LogWarning("MagneticField") << "MagGeometry::fieldInTesla: failed to find volume for " << gp << endl;
00094 }
00095 return GlobalVector();
00096 }
00097
00098
00099
00100 MagVolume*
00101 MagGeometry::findVolume1(const GlobalPoint & gp, double tolerance) const {
00102
00103 MagVolume6Faces * found = 0;
00104
00105 if (inBarrel(gp)) {
00106 for (vector<MagVolume6Faces*>::const_iterator v = theBVolumes.begin();
00107 v!=theBVolumes.end(); ++v){
00108 if ((*v)==0) {
00109 cout << endl << "***ERROR: MagGeometry::findVolume: MagVolume not set" << endl;
00110 continue;
00111 }
00112 if ((*v)->inside(gp,tolerance)) {
00113 found = (*v);
00114 break;
00115 }
00116 }
00117
00118 } else {
00119 for (vector<MagVolume6Faces*>::const_iterator v = theEVolumes.begin();
00120 v!=theEVolumes.end(); ++v){
00121 if ((*v)==0) {
00122 cout << endl << "***ERROR: MagGeometry::findVolume: MagVolume not set" << endl;
00123 continue;
00124 }
00125 if ((*v)->inside(gp,tolerance)) {
00126 found = (*v);
00127 break;
00128 }
00129 }
00130 }
00131
00132 return found;
00133 }
00134
00135
00136 MagVolume*
00137 MagGeometry::findVolume(const GlobalPoint & gp, double tolerance) const{
00138
00139 if (lastVolume!=0 && lastVolume->inside(gp)){
00140 return lastVolume;
00141 }
00142
00143 MagVolume * result=0;
00144 if (inBarrel(gp)) {
00145 double R = gp.perp();
00146 int bin = theBarrelBinFinder->binIndex(R);
00147
00148 for (int bin1 = bin; bin1 >= max(0,bin-2); --bin1) {
00149 if (verbose::debugOut) cout << "Trying layer at R " << theBLayers[bin1]->minR()
00150 << " " << R << endl ;
00151 result = theBLayers[bin1]->findVolume(gp, tolerance);
00152 if (verbose::debugOut) cout << "***In blayer " << bin1-bin << " "
00153 << (result==0? " failed " : " OK ") <<endl;
00154 if (result != 0) break;
00155 }
00156
00157 } else {
00158 Geom::Phi<float> phi = gp.phi();
00159 int bin = theEndcapBinFinder->binIndex(phi);
00160 if (verbose::debugOut) cout << "Trying endcap sector at phi "
00161 << theESectors[bin]->minPhi() << " " << phi << endl ;
00162 result = theESectors[bin]->findVolume(gp, tolerance);
00163 if (verbose::debugOut) cout << "***In guessed esector "
00164 << (result==0? " failed " : " OK ") <<endl;
00165 }
00166
00167
00168 if (result==0 && tolerance < 0.0001) {
00169
00170
00171
00172 if (verbose::debugOut) cout << "Increasing the tolerance to 0.03" <<endl;
00173 result = findVolume(gp, 0.03);
00174 }
00175
00176 if (cacheLastVolume) lastVolume = result;
00177
00178 return result;
00179 }
00180
00181
00182
00183
00184 bool MagGeometry::inBarrel(const GlobalPoint& gp) const {
00185 float Z = fabs(gp.z());
00186 float R = gp.perp();
00187
00188
00189 if (geometryVersion>=120812) {
00190 return (Z<350. ||
00191 (R>172.4 && Z<633.29) ||
00192 (R>308.7345 && Z<662.01));
00193 } else if (geometryVersion>=90812) {
00194 return (Z<350. ||
00195 (R>172.4 && Z<633.89) ||
00196 (R>308.755 && Z<662.01));
00197 } else {
00198 return (Z<350. ||
00199 (R>172.4 && Z<633.29) ||
00200 (R>308.755 && Z<661.01));
00201 }
00202 }
00203