source: trunk/source/visualization/HepRep/src/G4HepRepFileXMLWriter.cc@ 1219

Last change on this file since 1219 was 834, checked in by garnier, 17 years ago

import all except CVS

File size: 10.3 KB
Line 
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: G4HepRepFileXMLWriter.cc,v 1.3 2006/06/29 21:17:28 gunter Exp $
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
43#include "G4HepRepFileXMLWriter.hh"
44
45#include "G4ios.hh"
46
47G4HepRepFileXMLWriter::G4HepRepFileXMLWriter()
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 prevTypeName[i] = new char[1];
60 strcpy(prevTypeName[i],"");
61
62 inType[i] = false;
63 inInstance[i] = false;
64 }
65
66 inPrimitive = false;
67 inPoint = false;
68}
69
70void G4HepRepFileXMLWriter::addType(const char* name,int newTypeDepth)
71{
72 if (fout.good())
73 {
74 // Flatten structure if it exceeds maximum allowed typeDepth of 49.
75 if (newTypeDepth > 49)
76 newTypeDepth = 49;
77
78 if (newTypeDepth < 0)
79 newTypeDepth = 0;
80
81 // Insert any layers that are missing from the hierarchy (protects against
82 // callers that skip from, say, layer 1 to layer 3 with no layer 2).
83 while (typeDepth < (newTypeDepth-1)) {
84 addType("Layer Inserted by G4HepRepFileXMLWriter", typeDepth + 1);
85 addInstance();
86 }
87
88 // If moving closer to the root, close previously open types.
89 while (newTypeDepth<typeDepth)
90 endType();
91
92 // Close any remaining primitives of the current instance.
93 endPrimitive();
94
95 // If this is a new type name for the current depth, declare the
96 // new Type. Otherwise, it is just another Instance of the current Type.
97 if (strcmp(name,prevTypeName[newTypeDepth])!=0)
98 {
99 if (inType[newTypeDepth])
100 endType();
101
102 prevTypeName[newTypeDepth] = new char[strlen(name)+1];
103 strcpy(prevTypeName[newTypeDepth],name);
104
105 inType[newTypeDepth] = true;
106 indent();
107 fout << "<heprep:type version=\"null\" name=\"" << name << "\">"
108 << G4endl;
109
110 typeDepth = newTypeDepth;
111 }
112 } else {
113#ifdef G4HEPREPFILEDEBUG
114 G4cout << "G4HepRepFileXMLWriter:addType No file is currently open." << G4endl;
115#endif
116 }
117}
118
119void G4HepRepFileXMLWriter::addInstance()
120{
121 if (fout.good())
122 {
123 if (inType[typeDepth])
124 {
125 endInstance();
126 inInstance[typeDepth] = true;
127 indent();
128 fout << "<heprep:instance>" << G4endl;
129 } else {
130#ifdef G4HEPREPFILEDEBUG
131 G4cout << "G4HepRepFileXMLWriter:addInstance No HepRep Type is currently open" << G4endl;
132#endif
133 }
134 } else {
135#ifdef G4HEPREPFILEDEBUG
136 G4cout << "G4HepRepFileXMLWriter:addInstance No file is currently open" << G4endl;
137#endif
138 }
139}
140
141void G4HepRepFileXMLWriter::addPrimitive()
142{
143 if (fout.good())
144 {
145 if (inInstance[typeDepth])
146 {
147 endPrimitive();
148 inPrimitive = true;
149 indent();
150 fout << "<heprep:primitive>" << G4endl;
151 } else {
152#ifdef G4HEPREPFILEDEBUG
153 G4cout << "G4HepRepFileXMLWriter:addPrimitive No HepRep Instance is currently open" << G4endl;
154#endif
155 }
156 } else {
157#ifdef G4HEPREPFILEDEBUG
158 G4cout << "G4HepRepFileXMLWriter:addPrimitive No file is currently open" << G4endl;
159#endif
160 }
161}
162
163void G4HepRepFileXMLWriter::addPoint(double x, double y, double z)
164{
165 if (fout.good())
166 {
167 if (inPrimitive)
168 {
169 endPoint();
170 inPoint = true;
171 indent();
172 fout << "<heprep:point x=\"" << x << "\" y=\"" << y << "\" z=\"" << z << "\">" << G4endl;
173 } else {
174#ifdef G4HEPREPFILEDEBUG
175 G4cout << "G4HepRepFileXMLWriter:addPoint No HepRep Primitive is currently open" << G4endl;
176#endif
177 }
178 } else {
179#ifdef G4HEPREPFILEDEBUG
180 G4cout << "G4HepRepFileXMLWriter:addPoint No file is currently open" << G4endl;
181#endif
182 }
183}
184
185void G4HepRepFileXMLWriter::addAttDef(const char* name,
186 const char* desc,
187 const char* type,
188 const char* extra)
189{
190 if (fout.good())
191 {
192 indent();
193 fout << " <heprep:attdef extra=\"" << extra << "\" name=\"" << name << "\" type=\"" << type << "\"" << G4endl;
194 indent();
195 fout << " desc=\"" << desc << "\"/>" << G4endl;
196 } else {
197#ifdef G4HEPREPFILEDEBUG
198 G4cout << "G4HepRepFileXMLWriter:addAttDef No file is currently open" << G4endl;
199#endif
200 }
201}
202
203// Four methods to fill attValues
204void G4HepRepFileXMLWriter::addAttValue (const char* name,
205 const char* value)
206{
207 if (fout.good())
208 {
209 indent();
210 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
211 indent();
212 fout << " value=\"" << value << "\"/>" << G4endl;
213 } else {
214#ifdef G4HEPREPFILEDEBUG
215 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
216#endif
217 }
218}
219
220void G4HepRepFileXMLWriter::addAttValue (const char* name,
221 double value)
222{
223 if (fout.good())
224 {
225 indent();
226 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
227 indent();
228 fout << " value=\"" << value << "\"/>" << G4endl;
229 } else {
230#ifdef G4HEPREPFILEDEBUG
231 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
232#endif
233 }
234}
235
236void G4HepRepFileXMLWriter::addAttValue (const char* name,
237 int value)
238{
239 if (fout.good())
240 {
241 indent();
242 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
243 indent();
244 fout << " value=\"" << value << "\"/>" << G4endl;
245 } else {
246#ifdef G4HEPREPFILEDEBUG
247 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
248#endif
249 }
250}
251
252void G4HepRepFileXMLWriter::addAttValue (const char* name,
253 bool value)
254{
255 if (fout.good())
256 {
257 indent();
258 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
259 indent();
260 if (value)
261 fout << " value=\"True\"/>" << G4endl;
262 else
263 fout << " value=\"False\"/>" << G4endl;
264 } else {
265#ifdef G4HEPREPFILEDEBUG
266 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
267#endif
268 }
269}
270
271void G4HepRepFileXMLWriter::addAttValue (const char* name,
272 double value1,
273 double value2,
274 double value3)
275{
276 if (fout.good())
277 {
278 int redness = int(value1*255.);
279 int greenness = int(value2*255.);
280 int blueness = int(value3*255.);
281 indent();
282 fout << " <heprep:attvalue showLabel=\"NONE\" name=\"" << name << "\"" << G4endl;
283 indent();
284 fout << " value=\"" << redness << "," << greenness << "," << blueness << "\"/>" << G4endl;
285 } else {
286#ifdef G4HEPREPFILEDEBUG
287 G4cout << "G4HepRepFileXMLWriter:addAttValue No file is currently open" << G4endl;
288#endif
289 }
290}
291
292void G4HepRepFileXMLWriter::open(const char* fileSpec)
293{
294 if (isOpen)
295 close();
296
297 fout.open(fileSpec);
298
299 if (fout.good()) {
300 fout << "<?xml version=\"1.0\" ?>" << G4endl;
301 fout << "<heprep:heprep xmlns:heprep=\"http://www.slac.stanford.edu/~perl/heprep/\"" << G4endl;
302 fout << " xmlns:xsi=\"http://www.w3.org/1999/XMLSchema-instance\" xsi:schemaLocation=\"HepRep.xsd\">" << G4endl;
303
304 isOpen = true;
305 init();
306 } else {
307 G4cout << "G4HepRepFileXMLWriter:open Unable to write to file " << fileSpec << G4endl;
308 }
309}
310
311void G4HepRepFileXMLWriter::close()
312{
313 // Close any remaining open Types
314 endTypes();
315
316 if (fout.good()) {
317 fout << "</heprep:heprep>" << G4endl;
318 fout.close( );
319 isOpen = false;
320 } else {
321 G4cout << "G4HepRepFileXMLWriter:close No file is currently open" << G4endl;
322 }
323}
324
325void G4HepRepFileXMLWriter::endTypes()
326{
327 // Close any remaining open Types
328 while(typeDepth>-1)
329 endType();
330}
331
332void G4HepRepFileXMLWriter::endType()
333{
334 endInstance();
335 indent();
336 fout << "</heprep:type>" << G4endl;
337 inType[typeDepth] = false;
338 delete [] prevTypeName[typeDepth];
339 prevTypeName[typeDepth] = new char[1];
340 strcpy(prevTypeName[typeDepth],"");
341 typeDepth--;
342}
343
344void G4HepRepFileXMLWriter::endInstance()
345{
346 if (inInstance[typeDepth])
347 {
348 endPrimitive();
349 indent();
350 fout << "</heprep:instance>" << G4endl;
351 inInstance[typeDepth] = false;
352 }
353}
354
355void G4HepRepFileXMLWriter::endPrimitive()
356{
357 if (inPrimitive)
358 {
359 endPoint();
360 indent();
361 fout << "</heprep:primitive>" << G4endl;
362 inPrimitive = false;
363 }
364}
365
366void G4HepRepFileXMLWriter::endPoint()
367{
368 if (inPoint)
369 {
370 indent();
371 fout << "</heprep:point>" << G4endl;
372 inPoint = false;
373 }
374}
375
376void G4HepRepFileXMLWriter::indent()
377{
378 if (fout.good())
379 {
380 int i = 0;
381 while (inType[i] && i<12) {
382 fout << " ";
383 if (inInstance[i])
384 fout << " ";
385 i++;
386 }
387
388 if (inPrimitive)
389 fout << " ";
390 if (inPoint)
391 fout << " ";
392 }
393}
Note: See TracBrowser for help on using the repository browser.