]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TGeant4/TG4XMLConvertor.h
updated to introduction of new TG4SDServices class; added comment lines separating...
[u/mrichter/AliRoot.git] / TGeant4 / TG4XMLConvertor.h
1 // $Id$
2 // Category: geometry
3 // by I. Hrivnacova, 27.07.2000 
4 //
5 // XML convertor that converts G4 basic geometry objects 
6 // to XML defined by AGDD.dtd
7 // (ATLAS Generic Detector Description)
8
9 #ifndef TG4_XML_CONVERTOR_H
10 #define TG4_XML_CONVERTOR_H
11
12 #include "TG4VXMLConvertor.h"
13 #include "TG4Globals.h"
14
15 #include <globals.hh>
16 #include <g4std/fstream>
17
18 class G4Material;
19 class G4VSolid;
20 class G4LogicalVolume;
21 class G4PVReplica;
22 class G4Box;
23 class G4Tubs;
24 class G4Cons;
25 class G4Trd;
26 class G4Trap;
27 class G4Para;
28 class G4Polycone;
29 class G4Polyhedra;
30
31 class TG4XMLConvertor : public TG4VXMLConvertor
32 {
33   public:
34     TG4XMLConvertor(G4std::ofstream& outFile);
35     virtual ~TG4XMLConvertor();
36
37     // methods
38     virtual void OpenMaterials(const G4String& version, const G4String& date, 
39             const G4String& author, const G4String dtdVersion);
40     virtual void OpenSection(const G4String& name, const G4String& version,
41             const G4String& date, const G4String& author,
42             const G4String& topVolume);
43     virtual void OpenComposition(const G4String& name);
44     virtual void CloseMaterials();
45     virtual void CloseSection();
46     virtual void CloseComposition();
47
48     virtual void WriteMaterial(const G4Material* material); 
49     virtual void WriteSolid(G4String lvName, const G4VSolid* solid, 
50                             G4String materialName); 
51     virtual void WriteRotation(const G4RotationMatrix* rotation); 
52     virtual void WritePosition(G4String lvName, G4ThreeVector position); 
53     virtual void WritePositionWithRotation(
54                                G4String lvName, G4ThreeVector position,
55                                const G4RotationMatrix* rotation); 
56     virtual void WriteReplica(G4String lvName, G4PVReplica* pvr);                              
57     virtual void WriteEmptyLine();
58     virtual void IncreaseIndention();
59     virtual void DecreaseIndention();
60     
61     // set methods
62     void SetNumWidth(G4int width);
63     void SetNumPrecision(G4int precision);
64
65   private:
66     //methods
67     void CutName(G4String& name) const;
68     void CutName(G4String& name, G4int size) const;
69     void PutName(G4String& element, G4String name, G4String templ) const;
70     
71          // writing solids
72     void WriteBox (G4String lvName, const G4Box*  box,  G4String materialName); 
73     void WriteTubs(G4String lvName, const G4Tubs* tubs, G4String materialName); 
74     void WriteCons(G4String lvName, const G4Cons* cons, G4String materialName); 
75     void WriteTrd (G4String lvName, const G4Trd*  trd,  G4String materialName); 
76     void WriteTrap(G4String lvName, const G4Trap* trap, G4String materialName); 
77     void WritePara(G4String lvName, const G4Para* para, G4String materialName); 
78     void WritePolycone(G4String lvName, const G4Polycone* polycone, 
79                    G4String materialName); 
80     void WritePolyhedra(G4String lvName, const G4Polyhedra* polyhedra, 
81                    G4String materialName); 
82   
83     // static data members
84     static const G4int fgkMaxVolumeNameLength;  //maximal volume name length
85     static const G4int fgkMaxMaterialNameLength;//maximal material name length
86     static const G4int fgkDefaultNumWidth;      //default output numbers width
87     static const G4int fgkDefaultNumPrecision;  //default output numbers precision 
88
89     // data members
90     G4std::ofstream&  fOutFile;          //output file
91     const G4String    fkBasicIndention;  //basic indention 
92     G4String          fIndention;        //indention string
93     G4int             fNW;               //output numbers width
94     G4int             fNP;               //output numbers precision 
95     G4int             fRotationCounter;  //counter of rotations
96     TG4RotationMatrixVector  fRotations; //vector of rot matrices
97 };
98
99 inline void TG4XMLConvertor::SetNumWidth(G4int width)
100 { fNW = width; }
101
102 inline void TG4XMLConvertor::SetNumPrecision(G4int precision)
103 { fNP = precision; }
104
105 #endif //TG4_XML_CONVERTOR_H
106