ATLAS Offline Software
GeoPrimitivesToStringConverter.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // GeoPrimitvesToStringConverter.h, (c) ATLAS Detector software
8 
9 #ifndef GEOPRIMITIVESTOSTRINGCONVERTER_H_
10 #define GEOPRIMITIVESTOSTRINGCONVERTER_H_
11 
16 
17 #ifndef XAOD_STANDALONE
18 # include "CLHEP/Geometry/Transform3D.h"
19 # include "CLHEP/Geometry/Point3D.h"
20 # include "CLHEP/Vector/TwoVector.h"
21 #endif // not XAOD_STANDALONE
22 
23 namespace Amg {
24 
40  inline std::string toString( const Translation3D& translation, int precision = 4 ){
41  Vector3D trans{translation.x(), translation.y(), translation.z()};
42  return toString( trans, precision );
43  }
44 
45 
46  inline std::string toString( const Transform3D& transform, int precision = 4, const std::string& rotOffSet = "") {
47  std::stringstream sstr{};
48  bool printed{false};
49  if (transform.translation().mag() > std::numeric_limits<float>::epsilon()) {
50  sstr<<"translation: "<<toString(transform.translation(), precision);
51  printed = true;
52  }
54  if (printed) sstr<<rotOffSet<<", ";
55  sstr<<"rotation: {"<<toString(transform.linear()*Vector3D::UnitX(), precision)<<",";
56  sstr<<toString(transform.linear()*Vector3D::UnitY(), precision)<<",";
57  sstr<<toString(transform.linear()*Vector3D::UnitZ(), precision)<<"}";
58  printed = true;
59  }
60  if (!printed) sstr<<"Identity matrix ";
61  return sstr.str();
62 }
63 
64 #ifndef XAOD_STANDALONE
65 
66  inline std::string toString( const CLHEP::HepRotation& rot, int precision = 4, const std::string& offset="" ){
67  std::ostringstream sout;
68 
69  sout << std::setiosflags(std::ios::fixed) << std::setprecision(precision);
70  for( int i=0;i<3;++i ){
71  for( int j=0;j<3;++j ){
72  if( j == 0 ) sout << "(";
73  double val = roundWithPrecision(rot(i,j),precision);
74  sout << val;
75  if( j == 2 ) sout << ")";
76  else sout << ", ";
77  }
78  if( i != 2 ) {
79  sout << std::endl;
80  sout << offset;
81  }
82  }
83  return sout.str();
84  }
85 
86 
87  inline std::string toString( const CLHEP::Hep3Vector& translation, int precision = 4){
88  std::ostringstream sout;
89  sout << std::setiosflags(std::ios::fixed) << std::setprecision(precision);
90  for( int j=0;j<3;++j ){
91  if( j == 0 ) sout << "(";
92  double val = roundWithPrecision( translation[j],precision);
93  sout << val;
94  if( j == 2 ) sout << ")";
95  else sout << ", ";
96  }
97  return sout.str();
98  }
99 
100  inline std::string toString( const CLHEP::Hep2Vector& translation, int precision = 4){
101  std::ostringstream sout;
102  sout << std::setiosflags(std::ios::fixed) << std::setprecision(precision);
103  for( int j=0;j<2;++j ){
104  if( j == 0 ) sout << "(";
105  double val = roundWithPrecision( translation[j],precision);
106  sout << val;
107  if( j == 1 ) sout << ")";
108  else sout << ", ";
109  }
110  return sout.str();
111  }
112 
113  inline std::string toString( const HepGeom::Transform3D& transf, int precision = 4, const std::string& offset=""){
114  std::ostringstream sout;
115  sout << "Translation : " << toString( transf.getTranslation(), precision ) << std::endl;
116  std::string rotationOffset = offset + " ";
117  sout << offset << "Rotation : " << toString( transf.getRotation(), precision+2, rotationOffset );
118  return sout.str();
119  }
120 
121 #endif // not XAOD_STANDALONE
122 
123 }
124 
125 #endif /* GEOPRIMITIVESTOSTRINGCONVERTER_H_ */
GeoPrimitives.h
Amg::toString
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Definition: GeoPrimitivesToStringConverter.h:40
EventPrimitivesToStringConverter.h
lumiFormat.i
int i
Definition: lumiFormat.py:85
Amg::Transform3D
Eigen::Affine3d Transform3D
Definition: GeoPrimitives.h:46
Amg::doesNotDeform
bool doesNotDeform(const Amg::Transform3D &trans)
Checks whether the linear part of the transformation rotates or stetches any of the basis vectors.
Definition: GeoPrimitivesHelpers.h:361
Amg::transform
Amg::Vector3D transform(Amg::Vector3D &v, Amg::Transform3D &tr)
Transform a point from a Trasformation3D.
Definition: GeoPrimitivesHelpers.h:156
CLHEPtoEigenConverter.h
Amg
Definition of ATLAS Math & Geometry primitives (Amg)
Definition: AmgStringHelpers.h:19
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
GeoPrimitivesHelpers.h
Amg::Translation3D
Eigen::Translation< double, 3 > Translation3D
Definition: GeoPrimitives.h:44
convertTimingResiduals.offset
offset
Definition: convertTimingResiduals.py:71
Amg::roundWithPrecision
double roundWithPrecision(double val, int precision)
EventPrimitvesToStringConverter.
Definition: EventPrimitivesToStringConverter.h:35