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