SND@LHC Software
Loading...
Searching...
No Matches
convertNoisyMap.py
Go to the documentation of this file.
1#!/bin/python
2
3# Python script to convert a B field map from an ascii text file into
4# a ROOT file for FairShip. The field map needs to use a regular
5# x,y,z binning structure where the co-ordinates are assumed to
6# be in ascending z, y and x, in that order. Each data line needs to
7# contain the values "x y z Bx By Bz", with the field components in Tesla.
8# Comment lines in the datafile starting with the hash symbol # are ignored.
9# Distances for FairShip need to be in cm, so we use a scaling factor
10# cmScale (default = 1.0) to convert the text file distances into cm.
11# For example, if the input data uses mm for lengths, cmScale = 0.1.
12
13import ROOT
14import pandas as pd
15import os
16
17# Struct for the ROOT file TTree data: coord range and field binning
18
19ROOT.gROOT.ProcessLine(
20 "struct rangeStruct{\
21 float xMin;\
22 float xMax;\
23 float dx;\
24 float yMin;\
25 float yMax;\
26 float dy;\
27 float zMin;\
28 float zMax;\
29 float dz;\
30 };");
31
32# The field map is assumed to obey the following co-ordinate bin ordering:
33# z is increased first, y is increased 2nd, x is increased last.
34# For the coordinate bin (iX, iY, iZ), the field bin = (iX*Ny + iY)*Nz + iZ,
35# where Ny and Nz are the number of y and z bins
36
37ROOT.gROOT.ProcessLine(
38 "struct dataStruct{\
39 float x;\
40 float y;\
41 float z;\
42 float Bx;\
43 float By;\
44 float Bz;\
45 };");
46
47
48def run(inFileName='FieldTest.txt', rootFileName='BFieldTest.root',
49 cmScale=1.0, storeCoords=False):
50 createRootMap(inFileName, rootFileName, cmScale, storeCoords)
51
52
53def createRootMap(inFileName, rootFileName, cmScale, storeCoords):
54 print ('Create map {0} from {1} using cmScale = {2}'.format(rootFileName,
55 inFileName, cmScale))
56 if storeCoords is True:
57 print ('We will also store the x,y,z field coordinates in {0}'.format(rootFileName))
58
59 rangeInfo = findRanges(inFileName, cmScale)
60
61 # Define ROOT file and its TTree
62 theFile = ROOT.TFile.Open(rootFileName, 'recreate')
63
64 rangeTree = ROOT.TTree('Range', 'Range')
65 rangeTree.SetDirectory(theFile)
66
67 # Co-ordinate ranges
68 rStruct = ROOT.rangeStruct()
69 rangeTree.Branch('xMin', ROOT.AddressOf(rStruct, 'xMin'), 'xMin/F')
70 rangeTree.Branch('xMax', ROOT.AddressOf(rStruct, 'xMax'), 'xMax/F')
71 rangeTree.Branch('dx', ROOT.AddressOf(rStruct, 'dx'), 'dx/F')
72 rangeTree.Branch('yMin', ROOT.AddressOf(rStruct, 'yMin'), 'yMin/F')
73 rangeTree.Branch('yMax', ROOT.AddressOf(rStruct, 'yMax'), 'yMax/F')
74 rangeTree.Branch('dy', ROOT.AddressOf(rStruct, 'dy'), 'dy/F')
75 rangeTree.Branch('zMin', ROOT.AddressOf(rStruct, 'zMin'), 'zMin/F')
76 rangeTree.Branch('zMax', ROOT.AddressOf(rStruct, 'zMax'), 'zMax/F')
77 rangeTree.Branch('dz', ROOT.AddressOf(rStruct, 'dz'), 'dz/F')
78
79 rStruct.xMin = rangeInfo['xMin']
80 rStruct.xMax = rangeInfo['xMax']
81 rStruct.dx = rangeInfo['dx']
82 rStruct.yMin = rangeInfo['yMin']
83 rStruct.yMax = rangeInfo['yMax']
84 rStruct.dy = rangeInfo['dy']
85 rStruct.zMin = rangeInfo['zMin']
86 rStruct.zMax = rangeInfo['zMax']
87 rStruct.dz = rangeInfo['dz']
88
89 # Centre the field map on the local origin (cm)
90 x0 = 0#.5 * (rStruct.xMin + rStruct.xMax)
91 y0 = 0#.5 * (rStruct.yMin + rStruct.yMax)
92 z0 = 0.5 * (rStruct.zMin + rStruct.zMax)
93
94 # Use this if we don't want to centre the field map
95 # x0 = 0.0
96 # y0 = 0.0
97 # z0 = 0.0
98
99 print ('Centering field map using co-ordinate shift {0} {1} {2} cm'.format(x0, y0, z0))
100
101 # Center co-ordinate range limits (cm)
102 rStruct.xMin = rStruct.xMin - x0
103 rStruct.xMax = rStruct.xMax - x0
104
105 rStruct.yMin = rStruct.yMin - y0
106 rStruct.yMax = rStruct.yMax - y0
107
108 rStruct.zMin = rStruct.zMin - z0
109 rStruct.zMax = rStruct.zMax - z0
110
111 print ('x range = {0} to {1}'.format(rStruct.xMin, rStruct.xMax))
112 print ('y range = {0} to {1}'.format(rStruct.yMin, rStruct.yMax))
113 print ('z range = {0} to {1}'.format(rStruct.zMin, rStruct.zMax))
114
115 # Fill info into range tree
116 rangeTree.Fill()
117
118 # Store field data components
119 dataTree = ROOT.TTree('Data', 'Data')
120 dataTree.SetDirectory(theFile)
121
122 # Field components with (x,y,z) coordinate binning ordered such that
123 # z, then y, then x is increased. For the coordinate bin (iX, iY, iZ),
124 # the field bin = (iX*Ny + iY)*Nz + iZ, where Ny and Nz are the number
125 # of y and z bins
126 dStruct = ROOT.dataStruct()
127 if storeCoords is True:
128 dataTree.Branch('x', ROOT.AddressOf(dStruct, 'x'), 'x/F')
129 dataTree.Branch('y', ROOT.AddressOf(dStruct, 'y'), 'y/F')
130 dataTree.Branch('z', ROOT.AddressOf(dStruct, 'z'), 'z/F')
131
132 dataTree.Branch('Bx', ROOT.AddressOf(dStruct, 'Bx'), 'Bx/F')
133 dataTree.Branch('By', ROOT.AddressOf(dStruct, 'By'), 'By/F')
134 dataTree.Branch('Bz', ROOT.AddressOf(dStruct, 'Bz'), 'Bz/F')
135
136 # Reopen the file and store the information in the ROOT file
137
138 inData = pd.read_csv(inFileName, delim_whitespace=True, header=None)
139 inData.columns=["x", "y", "z", "bx", "by", "bz"]
140 inData = inData.sort_values(by=["x","y","z"])
141 inData = inData.astype(float)
142
143 count = 0.
144 data_shape = float(inData.shape[0])
145 for row in inData.itertuples():
146 if row.Index / data_shape >= count:
147 print("Processed: {} %".format(count * 100))
148 count += 0.1
149
150 # Bin centre coordinates with origin shift (all in cm)
151 if storeCoords is True:
152 dStruct.x = row.x * cmScale - x0
153 dStruct.y = row.y * cmScale - y0
154 dStruct.z = row.z * cmScale - z0
155
156 # B field components (Tesla)
157 dStruct.Bx = row.bx
158 dStruct.By = row.by
159 dStruct.Bz = row.bz
160
161 dataTree.Fill()
162
163 theFile.cd()
164 rangeTree.Write()
165 dataTree.Write()
166 theFile.Close()
167
168
169def findRanges(inFileName, cmScale):
170 # First read the data file to find the binning and coordinate ranges.
171 # Store the unique (ordered) x, y and z values so we can then find the
172 # bin widths, min/max limits and central offset
173
174 xArray = []
175 yArray = []
176 zArray = []
177 x_set = set()
178 y_set = set()
179 z_set = set()
180
181 with open(inFileName, 'r') as f:
182
183 # Read each line
184 for line in f:
185
186 # Ignore comment lines which begin with "#"
187 if '#' not in line:
188
189 sLine = line.split()
190
191 x = float(sLine[0]) * cmScale
192 y = float(sLine[1]) * cmScale
193 z = float(sLine[2]) * cmScale
194
195 if x not in x_set:
196 xArray.append(x)
197 x_set.add(x)
198
199 if y not in y_set:
200 yArray.append(y)
201 y_set.add(y)
202
203 if z not in z_set:
204 zArray.append(z)
205 z_set.add(z)
206
207 Nx = len(xArray)
208 Ny = len(yArray)
209 Nz = len(zArray)
210
211 xMin = 0.0
212 xMax = 0.0
213 dx = 0.0
214
215 if Nx > 0:
216 xMin = xArray[0]
217 Nx1 = Nx - 1
218 xMax = xArray[Nx1]
219 dx = (xMax - xMin) / (Nx1 * 1.0)
220
221 if Ny > 0:
222 yMin = yArray[0]
223 Ny1 = Ny - 1
224 yMax = yArray[Ny1]
225 dy = (yMax - yMin) / (Ny1 * 1.0)
226
227 if Nz > 0:
228 zMin = zArray[0]
229 Nz1 = Nz - 1
230 zMax = zArray[Nz1]
231 dz = (zMax - zMin) / (Nz1 * 1.0)
232
233 rangeInfo = {'Nx': Nx, 'xMin': xMin, 'xMax': xMax, 'dx': dx,
234 'Ny': Ny, 'yMin': yMin, 'yMax': yMax, 'dy': dy,
235 'Nz': Nz, 'zMin': zMin, 'zMax': zMax, 'dz': dz}
236
237 print ('rangeInfo = {0}'.format(rangeInfo))
238
239 return rangeInfo
240
241
242if __name__ == "__main__":
243 run(os.path.expandvars("$FAIRSHIP/files/noisy_fieldMap.csv"),
244 os.path.expandvars("$FAIRSHIP/files/MuonShieldField.root"), 1, False)
set(INCLUDE_DIRECTORIES ${SYSTEM_INCLUDE_DIRECTORIES} ${VMC_INCLUDE_DIRS} ${CMAKE_SOURCE_DIR}/shipdata ${CMAKE_SOURCE_DIR}/shipLHC ${CMAKE_SOURCE_DIR}/analysis/cuts ${CMAKE_SOURCE_DIR}/analysis/tools ${FMT_INCLUDE_DIR}) include_directories($
Definition CMakeLists.txt:1
findRanges(inFileName, cmScale)
run(inFileName='FieldTest.txt', rootFileName='BFieldTest.root', cmScale=1.0, storeCoords=False)
createRootMap(inFileName, rootFileName, cmScale, storeCoords)