Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
29//
30// Description:
31// Create a HepRep XML File (HepRep version 1).
32//
33// Environment:
34// Software developed for the general High Energy Physics community.
35//
36// Author :
37// J. Perl Original Author
38//
39// Copyright Information:
40// Copyright (C) 2001 Stanford Linear Accelerator Center
41//------------------------------------------------------------------------
42
44
45#include "G4HepRepMessenger.hh"
46#include "G4ios.hh"
47
49{
50 isOpen = false;
51 init();
52}
53
54void G4HepRepFileXMLWriter::init()
55{
56 typeDepth = -1;
57
58 int i = -1;
59 while (i++<49) {
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 addType("Layer Inserted by G4HepRepFileXMLWriter", typeDepth + 1);
87 }
88
89 // If moving closer to the root, close previously open types.
90 while (newTypeDepth<typeDepth)
91 endType();
92
93 // Close any remaining primitives of the current instance.
94 endPrimitive();
95
96 // If this is a new type name for the current depth, declare the
97 // new Type. Otherwise, it is just another Instance of the current Type.
98 if (strcmp(name,prevTypeName[newTypeDepth])!=0)
99 {
100 if (inType[newTypeDepth])
101 endType();
102
103 prevTypeName[newTypeDepth] = new char[strlen(name)+1];
104 strcpy(prevTypeName[newTypeDepth],name);
105
106 inType[newTypeDepth] = true;
107 indent();
108 fout << "<heprep:type version=\"null\" name=\"" << name << "\">"
109 << G4endl;
110
111 typeDepth = newTypeDepth;
112 }
113 } else {
114#ifdef G4HEPREPFILEDEBUG
115 G4cout << "G4HepRepFileXMLWriter:addType No file is currently open." << G4endl;
116#endif
117 }
118}
119
121{
122 if (fout.good())
123 {
124 if (inType[typeDepth])
125 {
126 endInstance();
127 inInstance[typeDepth] = true;
128 indent();
129 fout << "<heprep:instance>" << G4endl;
130 } else {
131#ifdef G4HEPREPFILEDEBUG
132 G4cout << "G4HepRepFileXMLWriter:addInstance No HepRep Type is currently open" << G4endl;
133#endif
134 }
135 } else {
136#ifdef G4HEPREPFILEDEBUG
137 G4cout << "G4HepRepFileXMLWriter:addInstance No file is currently open" << G4endl;
138#endif
139 }
140}
141
143{
144 if (fout.good())
145 {
147 {
148 endPrimitive();
149 inPrimitive = true;
150 indent();
151 fout << "<heprep:primitive>" << G4endl;
152 } else {
153#ifdef G4HEPREPFILEDEBUG
154 G4cout << "G4HepRepFileXMLWriter:addPrimitive No HepRep Instance is currently open" << G4endl;
155#endif
156 }
157 } else {
158#ifdef G4HEPREPFILEDEBUG
159 G4cout << "G4HepRepFileXMLWriter:addPrimitive No file is currently open" << G4endl;
160#endif
161 }
162}
163
164void G4HepRepFileXMLWriter::addPoint(double x, double y, double z)
165{
166 if (fout.good())
167 {
168 if (inPrimitive)
169 {
170 endPoint();
171 inPoint = true;
172 indent();
173
174 // Include scale and center values
176 G4double scale = messenger->getScale();
177 G4ThreeVector center = messenger->getCenter();
178 G4double xNew = scale * ( x - center.x());
179 G4double yNew = scale * ( y - center.y());
180 G4double zNew = scale * ( z - center.z());
181
182 fout << "<heprep:point x=\"" << xNew << "\" y=\"" << yNew << "\" z=\"" << zNew << "\">" << G4endl;
183 } else {
184#ifdef G4HEPREPFILEDEBUG
185 G4cout << "G4HepRepFileXMLWriter:addPoint No HepRep Primitive is currently open" << G4endl;
186#endif
187 }
188 } else {
189#ifdef G4HEPREPFILEDEBUG
190 G4cout << "G4HepRepFileXMLWriter:addPoint No file is currently open" << G4endl;
191#endif
192 }
193}
194
196 const char* desc,
197 const char* type,
198 const char* extra)
199{
200 if (fout.good())
201 {
202 indent();
203 fout << " <heprep:attdef extra=\"" << extra << "\" name=\"" << name << "\" type=\"" << type << "\"" << G4endl;
204 indent();
205 fout << " desc=\"" << desc << "\"/>" << G4endl;
206 } else {
207#ifdef G4HEPREPFILEDEBUG
208 G4cout << "G4HepRepFileXMLWriter:addAttDef No file is currently open" << G4endl;
209#endif
210 }
211}
212
213// Four methods to fill attValues
215 const char* value)
216{
217 if (fout.good())
218 {
219 indent();
220 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
221 indent();
222 fout << " value=\"" << value << "\"/>" << G4endl;
223 } else {
224#ifdef G4HEPREPFILEDEBUG
225 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
226#endif
227 }
228}
229
231 double value)
232{
233 if (fout.good())
234 {
235 indent();
236 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
237 indent();
238 fout << " value=\"" << value << "\"/>" << G4endl;
239 } else {
240#ifdef G4HEPREPFILEDEBUG
241 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
242#endif
243 }
244}
245
247 int value)
248{
249 if (fout.good())
250 {
251 indent();
252 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
253 indent();
254 fout << " value=\"" << value << "\"/>" << G4endl;
255 } else {
256#ifdef G4HEPREPFILEDEBUG
257 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
258#endif
259 }
260}
261
263 bool value)
264{
265 if (fout.good())
266 {
267 indent();
268 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
269 indent();
270 if (value)
271 fout << " value=\"True\"/>" << G4endl;
272 else
273 fout << " value=\"False\"/>" << G4endl;
274 } else {
275#ifdef G4HEPREPFILEDEBUG
276 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
277#endif
278 }
279}
280
282 double value1,
283 double value2,
284 double value3)
285{
286 if (fout.good())
287 {
288 int redness = int(value1*255.);
289 int greenness = int(value2*255.);
290 int blueness = int(value3*255.);
291 indent();
292 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
293 indent();
294 fout << " value=\"" << redness << "," << greenness << "," << blueness << "\"/>" << G4endl;
295 } else {
296#ifdef G4HEPREPFILEDEBUG
297 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
298#endif
299 }
300}
301
302void G4HepRepFileXMLWriter::open(const char* fileSpec)
303{
304 if (isOpen)
305 close();
306
307 fout.open(fileSpec);
308
309 if (fout.good()) {
310 fout << "<?xml version=\"1.0\" ?>" << G4endl;
311 fout << "<heprep:heprep xmlns:heprep=\"http://www.slac.stanford.edu/~perl/heprep/\"" << G4endl;
312 fout << " xmlns:xsi=\"http://www.w3.org/1999/XMLSchema-instance\" xsi:schemaLocation=\"HepRep.xsd\">" << G4endl;
313
314 isOpen = true;
315 init();
316 } else {
317 G4cout << "G4HepRepFileXMLWriter:open Unable to write to file " << fileSpec << G4endl;
318 }
319}
320
322{
323 // Close any remaining open Types
324 endTypes();
325
326 if (fout.good()) {
327 fout << "</heprep:heprep>" << G4endl;
328 fout.close( );
329 isOpen = false;
330 } else {
331 G4cout << "G4HepRepFileXMLWriter:close No file is currently open" << G4endl;
332 }
333}
334
336{
337 // Close any remaining open Types
338 while(typeDepth>-1)
339 endType();
340}
341
342void G4HepRepFileXMLWriter::endType()
343{
344 endInstance();
345 indent();
346 fout << "</heprep:type>" << G4endl;
347 inType[typeDepth] = false;
348 delete [] prevTypeName[typeDepth];
349 prevTypeName[typeDepth] = new char[1];
350 strcpy(prevTypeName[typeDepth],"");
351 typeDepth--;
352}
353
354void G4HepRepFileXMLWriter::endInstance()
355{
357 {
358 endPrimitive();
359 indent();
360 fout << "</heprep:instance>" << G4endl;
361 inInstance[typeDepth] = false;
362 }
363}
364
365void G4HepRepFileXMLWriter::endPrimitive()
366{
367 if (inPrimitive)
368 {
369 endPoint();
370 indent();
371 fout << "</heprep:primitive>" << G4endl;
372 inPrimitive = false;
373 }
374}
375
376void G4HepRepFileXMLWriter::endPoint()
377{
378 if (inPoint)
379 {
380 indent();
381 fout << "</heprep:point>" << G4endl;
382 inPoint = false;
383 }
384}
385
386void G4HepRepFileXMLWriter::indent()
387{
388 if (fout.good())
389 {
390 int i = 0;
391 while (inType[i] && i<12) {
392 fout << " ";
393 if (inInstance[i])
394 fout << " ";
395 i++;
396 }
397
398 if (inPrimitive)
399 fout << " ";
400 if (inPoint)
401 fout << " ";
402 }
403}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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()