1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215
|
/*
* PlacedBooleanVolume.cpp
*
* Created on: Nov 7, 2014
* Author: swenzel
*/
#include "VecGeom/volumes/PlacedBooleanVolume.h"
#include "VecGeom/volumes/SpecializedBooleanVolume.h"
#include "VecGeom/volumes/UnplacedBooleanVolume.h"
#include "VecGeom/volumes/LogicalVolume.h"
#include "VecGeom/base/Vector3D.h"
#include "VecGeom/base/RNG.h"
#include <map>
#ifdef VECGEOM_ROOT
#include "TGeoShape.h"
#include "TGeoVolume.h"
#include "TGeoCompositeShape.h"
#include "TGeoBoolNode.h"
#include "TGeoMatrix.h"
#include "TGeoManager.h"
#endif
#ifdef VECGEOM_GEANT4
#include "G4SubtractionSolid.hh"
#include "G4UnionSolid.hh"
#include "G4IntersectionSolid.hh"
#include "G4ThreeVector.hh"
#include "G4RotationMatrix.hh"
#endif
#include <iostream>
namespace vecgeom {
template <>
VECCORE_ATT_HOST_DEVICE void PlacedBooleanVolume<kUnion>::PrintType() const
{
VPlacedVolume const *pleft = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *pright = GetUnplacedVolume()->GetRight();
printf("PlacedBooleanVolume<kUnion>(");
pleft->PrintType();
printf(",");
pright->PrintType();
printf(")");
}
template <>
VECCORE_ATT_HOST_DEVICE void PlacedBooleanVolume<kIntersection>::PrintType() const
{
VPlacedVolume const *pleft = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *pright = GetUnplacedVolume()->GetRight();
printf("PlacedBooleanVolume<kIntersection>(");
pleft->PrintType();
printf(",");
pright->PrintType();
printf(")");
}
template <>
VECCORE_ATT_HOST_DEVICE void PlacedBooleanVolume<kSubtraction>::PrintType() const
{
VPlacedVolume const *pleft = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *pright = GetUnplacedVolume()->GetRight();
printf("PlacedBooleanVolume<kSubtraction>(");
pleft->PrintType();
printf(",");
pright->PrintType();
printf(")");
}
template <>
void PlacedBooleanVolume<kUnion>::PrintType(std::ostream &s) const
{
VPlacedVolume const *pleft = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *pright = GetUnplacedVolume()->GetRight();
s << "PlacedBooleanVolume<kUnion>(";
pleft->PrintType(s);
s << ",";
pright->PrintType(s);
s << ")";
}
template <>
void PlacedBooleanVolume<kIntersection>::PrintType(std::ostream &s) const
{
VPlacedVolume const *pleft = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *pright = GetUnplacedVolume()->GetRight();
s << "PlacedBooleanVolume<kIntersection>(";
pleft->PrintType(s);
s << ",";
pright->PrintType(s);
s << ")";
}
template <>
void PlacedBooleanVolume<kSubtraction>::PrintType(std::ostream &s) const
{
VPlacedVolume const *pleft = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *pright = GetUnplacedVolume()->GetRight();
s << "PlacedBooleanVolume<kSubtraction>(";
pleft->PrintType();
s << ",";
pright->PrintType();
s << ")";
}
#ifdef VECGEOM_ROOT
template <>
TGeoShape const *PlacedBooleanVolume<kUnion>::ConvertToRoot() const
{
VPlacedVolume const *left = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *right = GetUnplacedVolume()->GetRight();
Transformation3D const *leftm = left->GetTransformation();
Transformation3D const *rightm = right->GetTransformation();
TGeoUnion *node =
new TGeoUnion(const_cast<TGeoShape *>(left->ConvertToRoot()), const_cast<TGeoShape *>(right->ConvertToRoot()),
Transformation3D::ConvertToTGeoMatrix(*leftm), Transformation3D::ConvertToTGeoMatrix(*rightm));
return new TGeoCompositeShape("RootComposite", node);
}
template <>
TGeoShape const *PlacedBooleanVolume<kIntersection>::ConvertToRoot() const
{
VPlacedVolume const *left = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *right = GetUnplacedVolume()->GetRight();
Transformation3D const *leftm = left->GetTransformation();
Transformation3D const *rightm = right->GetTransformation();
TGeoIntersection *node = new TGeoIntersection(const_cast<TGeoShape *>(left->ConvertToRoot()),
const_cast<TGeoShape *>(right->ConvertToRoot()),
Transformation3D::ConvertToTGeoMatrix(*leftm), Transformation3D::ConvertToTGeoMatrix(*rightm));
return new TGeoCompositeShape("RootComposite", node);
}
template <>
TGeoShape const *PlacedBooleanVolume<kSubtraction>::ConvertToRoot() const
{
VPlacedVolume const *left = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *right = GetUnplacedVolume()->GetRight();
Transformation3D const *leftm = left->GetTransformation();
Transformation3D const *rightm = right->GetTransformation();
TGeoSubtraction *node = new TGeoSubtraction(
const_cast<TGeoShape *>(left->ConvertToRoot()), const_cast<TGeoShape *>(right->ConvertToRoot()),
Transformation3D::ConvertToTGeoMatrix(*leftm), Transformation3D::ConvertToTGeoMatrix(*rightm));
return new TGeoCompositeShape("RootComposite", node);
}
#endif
#ifdef VECGEOM_GEANT4
template <>
G4VSolid const *PlacedBooleanVolume<kUnion>::ConvertToGeant4() const
{
VPlacedVolume const *left = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *right = GetUnplacedVolume()->GetRight();
if (!left->GetTransformation()->IsIdentity()) {
std::cerr << "WARNING : For the moment left transformations are not implemented\n";
}
Transformation3D const *rightm = right->GetTransformation();
G4RotationMatrix *g4rot = new G4RotationMatrix();
auto rot = rightm->Rotation();
// HepRep3x3 seems? to be column major order:
g4rot->set(CLHEP::HepRep3x3(rot[0], rot[3], rot[6], rot[1], rot[4], rot[7], rot[2], rot[5], rot[8]));
return new G4UnionSolid(GetLabel(), const_cast<G4VSolid *>(left->ConvertToGeant4()),
const_cast<G4VSolid *>(right->ConvertToGeant4()), g4rot,
G4ThreeVector(rightm->Translation(0), rightm->Translation(1), rightm->Translation(2)));
}
template <>
G4VSolid const *PlacedBooleanVolume<kIntersection>::ConvertToGeant4() const
{
VPlacedVolume const *left = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *right = GetUnplacedVolume()->GetRight();
if (!left->GetTransformation()->IsIdentity()) {
std::cerr << "WARNING : For the moment left transformations are not implemented\n";
}
Transformation3D const *rightm = right->GetTransformation();
G4RotationMatrix *g4rot = new G4RotationMatrix();
auto rot = rightm->Rotation();
// HepRep3x3 seems? to be column major order:
g4rot->set(CLHEP::HepRep3x3(rot[0], rot[3], rot[6], rot[1], rot[4], rot[7], rot[2], rot[5], rot[8]));
return new G4IntersectionSolid(GetLabel(), const_cast<G4VSolid *>(left->ConvertToGeant4()),
const_cast<G4VSolid *>(right->ConvertToGeant4()), g4rot,
G4ThreeVector(rightm->Translation(0), rightm->Translation(1), rightm->Translation(2)));
}
template <>
G4VSolid const *PlacedBooleanVolume<kSubtraction>::ConvertToGeant4() const
{
VPlacedVolume const *left = GetUnplacedVolume()->GetLeft();
VPlacedVolume const *right = GetUnplacedVolume()->GetRight();
if (!left->GetTransformation()->IsIdentity()) {
std::cerr << "WARNING : For the moment left transformations are not implemented\n";
}
Transformation3D const *rightm = right->GetTransformation();
G4RotationMatrix *g4rot = new G4RotationMatrix();
return new G4SubtractionSolid(GetLabel(), const_cast<G4VSolid *>(left->ConvertToGeant4()),
const_cast<G4VSolid *>(right->ConvertToGeant4()), g4rot,
G4ThreeVector(rightm->Translation(0), rightm->Translation(1), rightm->Translation(2)));
}
#endif
#ifdef VECCORE_CUDA
VECGEOM_DEVICE_INST_PLACED_VOLUME_ALLSPEC_BOOLEAN(SpecializedBooleanVolume, kUnion)
VECGEOM_DEVICE_INST_PLACED_VOLUME_ALLSPEC_BOOLEAN(SpecializedBooleanVolume, kIntersection)
VECGEOM_DEVICE_INST_PLACED_VOLUME_ALLSPEC_BOOLEAN(SpecializedBooleanVolume, kSubtraction)
#endif // VECCORE_CUDA
} // End namespace vecgeom
|