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