SND@LHC Software
Loading...
Searching...
No Matches
convertMap.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
13from __future__ import print_function
14import ROOT
15
16# Struct for the ROOT file TTree data: coord range and field binning
17
18ROOT.gROOT.ProcessLine(
19"struct rangeStruct{\
20 float xMin;\
21 float xMax;\
22 float dx;\
23 float yMin;\
24 float yMax;\
25 float dy;\
26 float zMin;\
27 float zMax;\
28 float dz;\
29};");
30
31# The field map is assumed to obey the following co-ordinate bin ordering:
32# z is increased first, y is increased 2nd, x is increased last.
33# For the coordinate bin (iX, iY, iZ), the field bin = (iX*Ny + iY)*Nz + iZ,
34# where Ny and Nz are the number of y and z bins
35
36ROOT.gROOT.ProcessLine(
37"struct dataStruct{\
38 float x;\
39 float y;\
40 float z;\
41 float Bx;\
42 float By;\
43 float Bz;\
44};");
45
46
47def run(inFileName = 'FieldTest.txt', rootFileName = 'BFieldTest.root',
48 cmScale = 1.0, storeCoords = False):
49
50 createRootMap(inFileName, rootFileName, cmScale, storeCoords)
51
52
53def createRootMap(inFileName, rootFileName, cmScale, storeCoords):
54
55 print('Create map {0} from {1} using cmScale = {2}'.format(rootFileName,
56 inFileName, cmScale))
57 if storeCoords is True:
58 print('We will also store the x,y,z field coordinates in {0}'.format(rootFileName))
59
60 rangeInfo = findRanges(inFileName, cmScale)
61
62 # Define ROOT file and its TTree
63 theFile = ROOT.TFile.Open(rootFileName, 'recreate')
64
65 rangeTree = ROOT.TTree('Range', 'Range')
66 rangeTree.SetDirectory(theFile)
67
68 # Co-ordinate ranges
69 rStruct = ROOT.rangeStruct()
70 rangeTree.Branch('xMin', ROOT.AddressOf(rStruct, 'xMin'), 'xMin/F')
71 rangeTree.Branch('xMax', ROOT.AddressOf(rStruct, 'xMax'), 'xMax/F')
72 rangeTree.Branch('dx', ROOT.AddressOf(rStruct, 'dx'), 'dx/F')
73 rangeTree.Branch('yMin', ROOT.AddressOf(rStruct, 'yMin'), 'yMin/F')
74 rangeTree.Branch('yMax', ROOT.AddressOf(rStruct, 'yMax'), 'yMax/F')
75 rangeTree.Branch('dy', ROOT.AddressOf(rStruct, 'dy'), 'dy/F')
76 rangeTree.Branch('zMin', ROOT.AddressOf(rStruct, 'zMin'), 'zMin/F')
77 rangeTree.Branch('zMax', ROOT.AddressOf(rStruct, 'zMax'), 'zMax/F')
78 rangeTree.Branch('dz', ROOT.AddressOf(rStruct, 'dz'), 'dz/F')
79
80 rStruct.xMin = rangeInfo['xMin']
81 rStruct.xMax = rangeInfo['xMax']
82 rStruct.dx = rangeInfo['dx']
83 rStruct.yMin = rangeInfo['yMin']
84 rStruct.yMax = rangeInfo['yMax']
85 rStruct.dy = rangeInfo['dy']
86 rStruct.zMin = rangeInfo['zMin']
87 rStruct.zMax = rangeInfo['zMax']
88 rStruct.dz = rangeInfo['dz']
89
90 # Centre the field map on the local origin (cm)
91 x0 = 0.5*(rStruct.xMin + rStruct.xMax)
92 y0 = 0.5*(rStruct.yMin + rStruct.yMax)
93 z0 = 0.5*(rStruct.zMin + rStruct.zMax)
94
95 # Use this if we don't want to centre the field map
96 #x0 = 0.0
97 #y0 = 0.0
98 #z0 = 0.0
99
100 print('Centering field map using co-ordinate shift {0} {1} {2} cm'.format(x0, y0, z0))
101
102 # Center co-ordinate range limits (cm)
103 rStruct.xMin = rStruct.xMin - x0
104 rStruct.xMax = rStruct.xMax - x0
105
106 rStruct.yMin = rStruct.yMin - y0
107 rStruct.yMax = rStruct.yMax - y0
108
109 rStruct.zMin = rStruct.zMin - z0
110 rStruct.zMax = rStruct.zMax - z0
111
112 print('x range = {0} to {1}'.format(rStruct.xMin, rStruct.xMax))
113 print('y range = {0} to {1}'.format(rStruct.yMin, rStruct.yMax))
114 print('z range = {0} to {1}'.format(rStruct.zMin, rStruct.zMax))
115
116 # Fill info into range tree
117 rangeTree.Fill()
118
119 # Store field data components
120 dataTree = ROOT.TTree('Data', 'Data')
121 dataTree.SetDirectory(theFile)
122
123 # Field components with (x,y,z) coordinate binning ordered such that
124 # z, then y, then x is increased. For the coordinate bin (iX, iY, iZ),
125 # the field bin = (iX*Ny + iY)*Nz + iZ, where Ny and Nz are the number
126 # of y and z bins
127 dStruct = ROOT.dataStruct()
128 if storeCoords is True:
129 dataTree.Branch('x', ROOT.AddressOf(dStruct, 'x'), 'x/F')
130 dataTree.Branch('y', ROOT.AddressOf(dStruct, 'y'), 'y/F')
131 dataTree.Branch('z', ROOT.AddressOf(dStruct, 'z'), 'z/F')
132
133 dataTree.Branch('Bx', ROOT.AddressOf(dStruct, 'Bx'), 'Bx/F')
134 dataTree.Branch('By', ROOT.AddressOf(dStruct, 'By'), 'By/F')
135 dataTree.Branch('Bz', ROOT.AddressOf(dStruct, 'Bz'), 'Bz/F')
136
137 # Reopen the file and store the information in the ROOT file
138 with open(inFileName, 'r') as f:
139
140 # Read each line
141 for line in f:
142
143 # Ignore comment lines which begin with "#"
144 if '#' not in line:
145
146 sLine = line.split()
147
148 # Bin centre coordinates with origin shift (all in cm)
149 if storeCoords is True:
150 dStruct.x = float(sLine[0])*cmScale - x0
151 dStruct.y = float(sLine[1])*cmScale - y0
152 dStruct.z = float(sLine[2])*cmScale - z0
153
154 # B field components (Tesla)
155 dStruct.Bx = float(sLine[3])
156 dStruct.By = float(sLine[4])
157 dStruct.Bz = float(sLine[5])
158
159 dataTree.Fill()
160
161
162 theFile.cd()
163 rangeTree.Write()
164 dataTree.Write()
165 theFile.Close()
166
167
168def findRanges(inFileName, cmScale):
169
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
178 with open(inFileName, 'r') as f:
179
180 # Read each line
181 for line in f:
182
183 # Ignore comment lines which begin with "#"
184 if '#' not in line:
185
186 sLine = line.split()
187
188 x = float(sLine[0])*cmScale
189 y = float(sLine[1])*cmScale
190 z = float(sLine[2])*cmScale
191
192 if x not in xArray:
193 xArray.append(x)
194
195 if y not in yArray:
196 yArray.append(y)
197
198 if z not in zArray:
199 zArray.append(z)
200
201 Nx = len(xArray)
202 Ny = len(yArray)
203 Nz = len(zArray)
204
205 xMin = 0.0
206 xMax = 0.0
207 dx = 0.0
208
209 if Nx > 0:
210 xMin = xArray[0]
211 Nx1 = Nx - 1
212 xMax = xArray[Nx1]
213 dx = (xMax - xMin)/(Nx1*1.0)
214
215 if Ny > 0:
216 yMin = yArray[0]
217 Ny1 = Ny - 1
218 yMax = yArray[Ny1]
219 dy = (yMax - yMin)/(Ny1*1.0)
220
221 if Nz > 0:
222 zMin = zArray[0]
223 Nz1 = Nz - 1
224 zMax = zArray[Nz1]
225 dz = (zMax - zMin)/(Nz1*1.0)
226
227 rangeInfo = {'Nx': Nx, 'xMin': xMin, 'xMax': xMax, 'dx': dx,
228 'Ny': Ny, 'yMin': yMin, 'yMax': yMax, 'dy': dy,
229 'Nz': Nz, 'zMin': zMin, 'zMax': zMax, 'dz': dz}
230
231 #print 'rangeInfo = {0}'.format(rangeInfo)
232
233 return rangeInfo
234
235
236if __name__ == "__main__":
237
238 run('GoliathFieldMap.txt', 'GoliathFieldMap.root', 0.1, True)
239 #run('BFieldTest.txt', 'BFieldTest.root', 1.0)
createRootMap(inFileName, rootFileName, cmScale, storeCoords)
Definition convertMap.py:53
run(inFileName='FieldTest.txt', rootFileName='BFieldTest.root', cmScale=1.0, storeCoords=False)
Definition convertMap.py:48
findRanges(inFileName, cmScale)