Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4HepRepFileXMLWriter.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//--------------------------------------------------------------------------
27// File and Version Information:
28//
29// Description:
30// Create a HepRep XML File (HepRep version 1).
31//
32// Environment:
33// Software developed for the general High Energy Physics community.
34//
35// Author :
36// J. Perl Original Author
37//
38// Copyright Information:
39// Copyright (C) 2001 Stanford Linear Accelerator Center
40//------------------------------------------------------------------------
41
43
44#include "G4HepRepMessenger.hh"
45#include "G4ios.hh"
46
48{
49 isOpen = false;
50 init();
51}
52
53void G4HepRepFileXMLWriter::init()
54{
55 typeDepth = -1;
56
57 int i = -1;
58 while(i++ < 49)
59 {
60 prevTypeName[i] = new char[1];
61 strcpy(prevTypeName[i], "");
62
63 inType[i] = false;
64 inInstance[i] = false;
65 }
66
67 inPrimitive = false;
68 inPoint = false;
69}
70
71void G4HepRepFileXMLWriter::addType(const char* name, int newTypeDepth)
72{
73 if(fout.good())
74 {
75 // Flatten structure if it exceeds maximum allowed typeDepth of 49.
76 if(newTypeDepth > 49)
77 newTypeDepth = 49;
78
79 if(newTypeDepth < 0)
80 newTypeDepth = 0;
81
82 // Insert any layers that are missing from the hierarchy (protects against
83 // callers that skip from, say, layer 1 to layer 3 with no layer 2).
84 while(typeDepth < (newTypeDepth - 1))
85 {
86 addType("Layer Inserted by G4HepRepFileXMLWriter", typeDepth + 1);
88 }
89
90 // If moving closer to the root, close previously open types.
91 while(newTypeDepth < typeDepth)
92 endType();
93
94 // Close any remaining primitives of the current instance.
95 endPrimitive();
96
97 // If this is a new type name for the current depth, declare the
98 // new Type. Otherwise, it is just another Instance of the current Type.
99 if(strcmp(name, prevTypeName[newTypeDepth]) != 0)
100 {
101 if(inType[newTypeDepth])
102 endType();
103
104 prevTypeName[newTypeDepth] = new char[strlen(name) + 1];
105 strcpy(prevTypeName[newTypeDepth], name);
106
107 inType[newTypeDepth] = true;
108 indent();
109 fout << "<heprep:type version=\"null\" name=\"" << name << "\">"
110 << G4endl;
111
112 typeDepth = newTypeDepth;
113 }
114 }
115 else
116 {
117#ifdef G4HEPREPFILEDEBUG
118 G4cout << "G4HepRepFileXMLWriter:addType No file is currently open."
119 << G4endl;
120#endif
121 }
122}
123
125{
126 if(fout.good())
127 {
128 if(inType[typeDepth])
129 {
130 endInstance();
131 inInstance[typeDepth] = true;
132 indent();
133 fout << "<heprep:instance>" << G4endl;
134 }
135 else
136 {
137#ifdef G4HEPREPFILEDEBUG
138 G4cout
139 << "G4HepRepFileXMLWriter:addInstance No HepRep Type is currently open"
140 << G4endl;
141#endif
142 }
143 }
144 else
145 {
146#ifdef G4HEPREPFILEDEBUG
147 G4cout << "G4HepRepFileXMLWriter:addInstance No file is currently open"
148 << G4endl;
149#endif
150 }
151}
152
154{
155 if(fout.good())
156 {
158 {
159 endPrimitive();
160 inPrimitive = true;
161 indent();
162 fout << "<heprep:primitive>" << G4endl;
163 }
164 else
165 {
166#ifdef G4HEPREPFILEDEBUG
167 G4cout << "G4HepRepFileXMLWriter:addPrimitive No HepRep Instance is "
168 "currently open"
169 << G4endl;
170#endif
171 }
172 }
173 else
174 {
175#ifdef G4HEPREPFILEDEBUG
176 G4cout << "G4HepRepFileXMLWriter:addPrimitive No file is currently open"
177 << G4endl;
178#endif
179 }
180}
181
182void G4HepRepFileXMLWriter::addPoint(double x, double y, double z)
183{
184 if(fout.good())
185 {
186 if(inPrimitive)
187 {
188 endPoint();
189 inPoint = true;
190 indent();
191
192 // Include scale and center values
194 G4double scale = messenger->getScale();
195 G4ThreeVector center = messenger->getCenter();
196 G4double xNew = scale * (x - center.x());
197 G4double yNew = scale * (y - center.y());
198 G4double zNew = scale * (z - center.z());
199
200 fout << "<heprep:point x=\"" << xNew << "\" y=\"" << yNew << "\" z=\""
201 << zNew << "\">" << G4endl;
202 }
203 else
204 {
205#ifdef G4HEPREPFILEDEBUG
206 G4cout << "G4HepRepFileXMLWriter:addPoint No HepRep Primitive is "
207 "currently open"
208 << G4endl;
209#endif
210 }
211 }
212 else
213 {
214#ifdef G4HEPREPFILEDEBUG
215 G4cout << "G4HepRepFileXMLWriter:addPoint No file is currently open"
216 << G4endl;
217#endif
218 }
219}
220
221void G4HepRepFileXMLWriter::addAttDef(const char* name, const char* desc,
222 const char* type, const char* extra)
223{
224 if(fout.good())
225 {
226 indent();
227 fout << " <heprep:attdef extra=\"" << extra << "\" name=\"" << name
228 << "\" type=\"" << type << "\"" << G4endl;
229 indent();
230 fout << " desc=\"" << desc << "\"/>" << G4endl;
231 }
232 else
233 {
234#ifdef G4HEPREPFILEDEBUG
235 G4cout << "G4HepRepFileXMLWriter:addAttDef No file is currently open"
236 << G4endl;
237#endif
238 }
239}
240
241// Four methods to fill attValues
242void G4HepRepFileXMLWriter::addAttValue(const char* name, const char* value)
243{
244 if(fout.good())
245 {
246 indent();
247 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\""
248 << G4endl;
249 indent();
250 fout << " value=\"" << value << "\"/>" << G4endl;
251 }
252 else
253 {
254#ifdef G4HEPREPFILEDEBUG
255 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open"
256 << G4endl;
257#endif
258 }
259}
260
261void G4HepRepFileXMLWriter::addAttValue(const char* name, double value)
262{
263 if(fout.good())
264 {
265 indent();
266 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\""
267 << G4endl;
268 indent();
269 fout << " value=\"" << value << "\"/>" << G4endl;
270 }
271 else
272 {
273#ifdef G4HEPREPFILEDEBUG
274 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open"
275 << G4endl;
276#endif
277 }
278}
279
280void G4HepRepFileXMLWriter::addAttValue(const char* name, int value)
281{
282 if(fout.good())
283 {
284 indent();
285 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\""
286 << G4endl;
287 indent();
288 fout << " value=\"" << value << "\"/>" << G4endl;
289 }
290 else
291 {
292#ifdef G4HEPREPFILEDEBUG
293 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open"
294 << G4endl;
295#endif
296 }
297}
298
299void G4HepRepFileXMLWriter::addAttValue(const char* name, bool value)
300{
301 if(fout.good())
302 {
303 indent();
304 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\""
305 << G4endl;
306 indent();
307 if(value)
308 fout << " value=\"True\"/>" << G4endl;
309 else
310 fout << " value=\"False\"/>" << G4endl;
311 }
312 else
313 {
314#ifdef G4HEPREPFILEDEBUG
315 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open"
316 << G4endl;
317#endif
318 }
319}
320
321void G4HepRepFileXMLWriter::addAttValue(const char* name, double value1,
322 double value2, double value3)
323{
324 if(fout.good())
325 {
326 int redness = int(value1 * 255.);
327 int greenness = int(value2 * 255.);
328 int blueness = int(value3 * 255.);
329 indent();
330 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\""
331 << G4endl;
332 indent();
333 fout << " value=\"" << redness << "," << greenness << "," << blueness
334 << "\"/>" << G4endl;
335 }
336 else
337 {
338#ifdef G4HEPREPFILEDEBUG
339 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open"
340 << G4endl;
341#endif
342 }
343}
344
345void G4HepRepFileXMLWriter::open(const char* fileSpec)
346{
347 if(isOpen)
348 close();
349
350 fout.open(fileSpec);
351
352 if(fout.good())
353 {
354 fout << "<?xml version=\"1.0\" ?>" << G4endl;
355 fout << "<heprep:heprep "
356 "xmlns:heprep=\"http://www.slac.stanford.edu/~perl/heprep/\""
357 << G4endl;
358 fout << " xmlns:xsi=\"http://www.w3.org/1999/XMLSchema-instance\" "
359 "xsi:schemaLocation=\"HepRep.xsd\">"
360 << G4endl;
361
362 isOpen = true;
363 init();
364 }
365 else
366 {
367 G4cout << "G4HepRepFileXMLWriter:open Unable to write to file " << fileSpec
368 << G4endl;
369 }
370}
371
373{
374 // Close any remaining open Types
375 endTypes();
376
377 if(fout.good())
378 {
379 fout << "</heprep:heprep>" << G4endl;
380 fout.close();
381 isOpen = false;
382 }
383 else
384 {
385 G4cout << "G4HepRepFileXMLWriter:close No file is currently open" << G4endl;
386 }
387}
388
390{
391 // Close any remaining open Types
392 while(typeDepth > -1)
393 endType();
394}
395
396void G4HepRepFileXMLWriter::endType()
397{
398 endInstance();
399 indent();
400 fout << "</heprep:type>" << G4endl;
401 inType[typeDepth] = false;
402 delete[] prevTypeName[typeDepth];
403 prevTypeName[typeDepth] = new char[1];
404 strcpy(prevTypeName[typeDepth], "");
405 typeDepth--;
406}
407
408void G4HepRepFileXMLWriter::endInstance()
409{
411 {
412 endPrimitive();
413 indent();
414 fout << "</heprep:instance>" << G4endl;
415 inInstance[typeDepth] = false;
416 }
417}
418
419void G4HepRepFileXMLWriter::endPrimitive()
420{
421 if(inPrimitive)
422 {
423 endPoint();
424 indent();
425 fout << "</heprep:primitive>" << G4endl;
426 inPrimitive = false;
427 }
428}
429
430void G4HepRepFileXMLWriter::endPoint()
431{
432 if(inPoint)
433 {
434 indent();
435 fout << "</heprep:point>" << G4endl;
436 inPoint = false;
437 }
438}
439
440void G4HepRepFileXMLWriter::indent()
441{
442 if(fout.good())
443 {
444 int i = 0;
445 while(inType[i] && i < 12)
446 {
447 fout << " ";
448 if(inInstance[i])
449 fout << " ";
450 i++;
451 }
452
453 if(inPrimitive)
454 fout << " ";
455 if(inPoint)
456 fout << " ";
457 }
458}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
void addAttValue(const char *name, const char *value)
void addPoint(double x, double y, double z)
void addType(const char *name, int newTypeDepth)
void open(const char *filespec)
void addAttDef(const char *name, const char *desc, const char *type, const char *extra)
virtual G4double getScale()
static G4HepRepMessenger * GetInstance()
virtual G4ThreeVector getCenter()