ATLAS Offline Software
Loading...
Searching...
No Matches
AtlasFieldMap.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
11
12#ifndef MAGFIELDCONDITIONS_ATLASFIELDMAP_H
13#define MAGFIELDCONDITIONS_ATLASFIELDMAP_H 1
14
15// MagField includes
16#include "CxxUtils/restrict.h"
17#include "GaudiKernel/ServiceHandle.h"
23
24#include <iostream>
25#include <memory>
26
27// forward declarations
28class TFile;
29
30namespace MagField {
31
39{
40public:
41 AtlasFieldMap() = default;
43
44 // initialize map from root file
45 bool initializeMap(TFile* rootfile,
46 float solenoidCurrent,
47 float toroidCurrent);
48
49 // Functions used by getField[ZR] in AtlasFieldCache
50 // search for a "zone" to which the point (z,r,phi) belongs
51 const BFieldZone* findBFieldZone(double z, double r, double phi) const;
52
53 // fast 2d map (made of one zone)
55
57 bool solenoidOn() const { return solenoidCurrent() > 0.0; }
58 bool toroidOn() const { return toroidCurrent() > 0.0; }
59 // magnet currents read with map - needed for scaling
60 float solenoidCurrent() const { return m_solenoidCurrent; }
61 float toroidCurrent() const { return m_toroidCurrent; }
62 int solenoidZoneId() const { return m_solenoidZoneId; }
63
64private:
66 AtlasFieldMap(const AtlasFieldMap& other) = delete;
67 AtlasFieldMap& operator=(const AtlasFieldMap& other) = delete;
68 AtlasFieldMap(AtlasFieldMap&& other) = delete;
69
70 // slow zone search is used during initialization to build the LUT
71 BFieldZone* findZoneSlow(double z, double r, double phi);
72
73 // utility functions used by readMap
74 int read_packed_data(std::istream& input, std::vector<int>& data) const;
75 int read_packed_int(std::istream& input, int& n) const;
76 void buildLUT();
77 void buildZR();
78
80 int memSize() const;
81
83
84 // field map name
85 std::string m_filename;
86
87 // currents read in with map
88 float m_solenoidCurrent{ 0 }; // solenoid current in ampere
89 float m_toroidCurrent{ 0 }; // toroid current in ampere
90 int m_solenoidZoneId{ -1 }; // solenoid zone id
91
92 // full 3d map (made of multiple zones)
93 std::vector<BFieldZone> m_zone;
94
95 // fast 2d map (made of one zone)
97
98 // data members used in zone-finding
99 std::vector<double> m_edge[3]; // zone boundaries in z, r, phi
100 std::vector<int> m_edgeLUT[3]; // look-up table for zone edges
101 double m_invq[3]{}; // 1/stepsize in m_edgeLUT
102 std::vector<const BFieldZone*> m_zoneLUT; // look-up table for zones
103 // more data members to speed up zone-finding
104 double m_zmin{ 0 }; // minimum z
105 double m_zmax{ 0 }; // maximum z
106 int m_nz{ 0 }; // number of z bins in zoneLUT
107 double m_rmax{ 0 }; // maximum r
108 int m_nr{ 0 }; // number of r bins in zoneLUT
109 int m_nphi{ 0 }; // number of phi bins in zoneLUT
110 bool m_mapIsInitialized{ false };
111};
112
113} // namespace MagField
114
115#include "AtlasFieldMap.icc"
116#endif // MAGFIELDCONDITIONS_ATLASFIELDMAP_H
Scalar phi() const
phi method
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
#define z
AtlasFieldMap & operator=(AtlasFieldMap &&other)=delete
std::string m_filename
Data Members.
float toroidCurrent() const
int read_packed_data(std::istream &input, std::vector< int > &data) const
int read_packed_int(std::istream &input, int &n) const
std::vector< BFieldZone > m_zone
std::vector< double > m_edge[3]
BFieldMeshZR * m_meshZR
const BFieldMeshZR * getBFieldMesh() const
float solenoidCurrent() const
bool solenoidOn() const
status of the magnets
const BFieldZone * findBFieldZone(double z, double r, double phi) const
BFieldZone * findZoneSlow(double z, double r, double phi)
AtlasFieldMap(const AtlasFieldMap &other)=delete
std::vector< const BFieldZone * > m_zoneLUT
AtlasFieldMap(AtlasFieldMap &&other)=delete
AtlasFieldMap & operator=(const AtlasFieldMap &other)=delete
std::vector< int > m_edgeLUT[3]
int memSize() const
approximate memory footprint in bytes
bool initializeMap(TFile *rootfile, float solenoidCurrent, float toroidCurrent)
int r
Definition globals.cxx:22
static std::vector< std::string > rootfile
Definition iLumiCalc.h:30
Local cache for magnetic field (based on MagFieldServices/AtlasFieldSvcTLS.h)
Macro wrapping the nonstandard restrict keyword.