CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/CondTools/DT/plugins/DTHVCheckByAbsoluteValues.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2010/09/14 13:54:22 $
00005  *  $Revision: 1.5 $
00006  *  \author Paolo Ronchese INFN Padova
00007  *
00008  */
00009 
00010 //-----------------------
00011 // This Class' Header --
00012 //-----------------------
00013 #include "CondTools/DT/plugins/DTHVCheckByAbsoluteValues.h"
00014 
00015 //-------------------------------
00016 // Collaborating Class Headers --
00017 //-------------------------------
00018 #include "DataFormats/MuonDetId/interface/DTWireId.h"
00019 #include "FWCore/PluginManager/interface/ModuleDef.h"
00020 #include "FWCore/ServiceRegistry/interface/ServiceMaker.h"
00021 
00022 //---------------
00023 // C++ Headers --
00024 //---------------
00025 #include <iostream>
00026 
00027 namespace cond { namespace service {
00028 
00029 //-------------------
00030 // Initializations --
00031 //-------------------
00032 
00033 
00034 //----------------
00035 // Constructors --
00036 //----------------
00037 DTHVCheckByAbsoluteValues::DTHVCheckByAbsoluteValues(
00038                            const edm::ParameterSet & iConfig, 
00039                            edm::ActivityRegistry & iAR ) {
00040   if ( instance == 0 ) {
00041     std::cout << "create DTHVCheckByAbsoluteValues" << std::endl;
00042     minHV = new float[4];
00043     maxHV = new float[4];
00044 //    minHV[0] = 3500.0;
00045 //    minHV[1] = 3500.0;
00046 //    minHV[2] = 1400.0;
00047 //    minHV[3] =  800.0;
00048     minHV[0] = 3500.0;
00049     minHV[1] = 3500.0;
00050     minHV[2] = 1700.0;
00051     minHV[3] = 1100.0;
00052     maxHV[0] = 4000.0;
00053     maxHV[1] = 4000.0;
00054     maxHV[2] = 2200.0;
00055     maxHV[3] = 1600.0;
00056     maxCurrent = 30.0;
00057     instance = this;
00058   }
00059 }
00060 
00061 //--------------
00062 // Destructor --
00063 //--------------
00064 DTHVCheckByAbsoluteValues::~DTHVCheckByAbsoluteValues() {
00065 }
00066 
00067 //--------------
00068 // Operations --
00069 //--------------
00070 DTHVAbstractCheck::flag DTHVCheckByAbsoluteValues::checkCurrentStatus(
00071                         int rawId, int type,
00072                         float valueA, float valueC, float valueS,
00073                         const std::map<int,timedMeasurement>& snapshotValues,
00074                         const std::map<int,int>& aliasMap,
00075                         const std::map<int,int>& layerMap ) {
00076 
00077 // find all values for this channel
00078 //  ind dpid = 0;
00079 //  std::map<int,int>::const_iterator lpartIter;
00080 //  std::map<int,int>::const_iterator lpartIend = layerMap.end();
00081 //  if ( ( layerIter = layerMap.find( chp0 ) ) != layerIend ) 
00082 //      dpid = layerIter.second;
00083 //  std::map<int,timedMeasurement>::const_iterator snapIter;
00084 //  std::map<int,timedMeasurement>::const_iterator snapIend =
00085 //                                                 snapshotValues.end();
00086 //  float val1 = -999999.0;
00087 //  float val2 = -999999.0;
00088 //  int chan = dpId * 10;
00089 //  if ( ( snapIter = snapshotValues.find( chan + 1 ) ) != snapIend )
00090 //      val1 = snapIter->second.second;
00091 //  if ( ( snapIter = snapshotValues.find( chan + 2 ) ) != snapIend )
00092 //      val2 = snapIter->second.second;
00093 
00094 // find dp identifier for all channels in this layer
00095 //  DTLayerId lay = chlId.layerId();
00096 //  int chp0 = DTWireId( lay, 10 ).rawId();
00097 //  int chp1 = DTWireId( lay, 11 ).rawId();
00098 //  int chp2 = DTWireId( lay, 12 ).rawId();
00099 //  int chp3 = DTWireId( lay, 13 ).rawId();
00100 //  ind dpi0 = 0;
00101 //  ind dpi1 = 0;
00102 //  ind dpi2 = 0;
00103 //  ind dpi3 = 0;
00104 //  std::map<int,int>::const_iterator layerIter;
00105 //  std::map<int,int>::const_iterator layerIend = layerMap.end();
00106 //  if ( ( layerIter = layerMap.find( chp0 ) ) != layerIend ) 
00107 //      dpi0 = layerIter.second;
00108 //  if ( ( layerIter = layerMap.find( chp1 ) ) != layerIend ) 
00109 //      dpi1 = layerIter.second;
00110 //  if ( ( layerIter = layerMap.find( chp2 ) ) != layerIend ) 
00111 //      dpi2 = layerIter.second;
00112 //  if ( ( layerIter = layerMap.find( chp3 ) ) != layerIend ) 
00113 //      dpi3 = layerIter.second;
00114 
00115   DTWireId chlId( rawId );
00116   int part = chlId.wire() - 10;
00117   DTHVAbstractCheck::flag flag;
00118   flag.a = flag.c = flag.s = 0;
00119   if ( type == 1 ) {
00120     if ( valueA < minHV[part] ) flag.a = 1;
00121     if ( valueA > maxHV[part] ) flag.a = 2;
00122     if ( valueS < minHV[   2] ) flag.s = 1;
00123     if ( valueS > maxHV[   2] ) flag.s = 2;
00124     if ( valueC < minHV[   3] ) flag.c = 1;
00125     if ( valueC > maxHV[   3] ) flag.c = 2;
00126   }
00127   if ( type == 2 ) {
00128     float voltA = 0.0;
00129     float voltS = 0.0;
00130     float voltC = 0.0;
00131     DTLayerId lay = chlId.layerId();
00132     int l_p = chlId.wire();
00133     DTWireId chA( lay, l_p );
00134     DTWireId chS( lay, 12 );
00135     DTWireId chC( lay, 13 );
00136     std::map<int,int>::const_iterator layerIter;
00137     std::map<int,int>::const_iterator layerIend = layerMap.end();
00138     std::map<int,timedMeasurement>::const_iterator snapIter;
00139     std::map<int,timedMeasurement>::const_iterator snapIend =
00140                                                    snapshotValues.end();
00141     int chan;
00142     if ( ( layerIter = layerMap.find( chA.rawId() ) ) != layerIend ) {
00143       chan = ( layerIter->second * 10 ) + l_p;
00144       if ( ( snapIter = snapshotValues.find( chan ) ) != snapIend ) {
00145         voltA = snapIter->second.second;
00146       }
00147     }
00148     if ( ( layerIter = layerMap.find( chS.rawId() ) ) != layerIend ) {
00149       chan = ( layerIter->second * 10 ) + 2;
00150       if ( ( snapIter = snapshotValues.find( chan ) ) != snapIend ) {
00151         voltS = snapIter->second.second;
00152       }
00153     }
00154     if ( ( layerIter = layerMap.find( chC.rawId() ) ) != layerIend ) {
00155       chan = ( layerIter->second * 10 ) + 3;
00156       if ( ( snapIter = snapshotValues.find( chan ) ) != snapIend ) {
00157         voltC = snapIter->second.second;
00158       }
00159     }
00160     if ( ( valueA > maxCurrent  ) &&
00161          ( voltA >= minHV[part] ) ) flag.a = 4;
00162     if ( ( valueS > maxCurrent  ) &&
00163          ( voltS >= minHV[   2] ) ) flag.s = 4;
00164     if ( ( valueC > maxCurrent  ) &&
00165          ( voltC >= minHV[   3] ) ) flag.c = 4;
00166   }
00167   return flag;
00168 
00169 }
00170 
00171 
00172 DEFINE_FWK_SERVICE( DTHVCheckByAbsoluteValues );
00173 } }
00174 
00175