00001 
00002 
00003 
00004 
00005 
00006 
00007 
00008 
00009 
00010 
00011 
00012 
00013 
00014 
00015 
00016 
00017 
00018 
00019 
00020 
00021 
00022 
00023 
00024 
00025 
00026 
00027 
00028 
00029 
00030 
00031 
00032 
00033 
00034 #include "G4PVReplica.hh"
00035 #include "G4LogicalVolume.hh"
00036 
00037 G4PVReplica::G4PVReplica( const G4String& pName,
00038                                 G4LogicalVolume* pLogical,
00039                                 G4VPhysicalVolume* pMother,
00040                           const EAxis pAxis,
00041                           const G4int nReplicas,
00042                           const G4double width,
00043                           const G4double offset )
00044   : G4VPhysicalVolume(0, G4ThreeVector(), pName, pLogical, pMother),
00045     fcopyNo(-1), fRegularVolsId(0)
00046 {
00047   if ((!pMother) || (!pMother->GetLogicalVolume()))
00048   {
00049     std::ostringstream message;
00050     message << "NULL pointer specified as mother volume." << G4endl
00051             << "The world volume cannot be sliced or parameterised !";
00052     G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002",
00053                 FatalException, message);
00054     return;
00055   }
00056   G4LogicalVolume* motherLogical = pMother->GetLogicalVolume();
00057   if (pLogical == motherLogical)
00058   {
00059     G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002",
00060                 FatalException, "Cannot place a volume inside itself!");
00061     return;
00062   }
00063   SetMotherLogical(motherLogical);
00064   motherLogical->AddDaughter(this);
00065   if (motherLogical->GetNoDaughters() != 1)
00066   {
00067     std::ostringstream message;
00068     message << "Replica or parameterised volume must be the only daughter !"
00069             << G4endl
00070             << "     Mother physical volume: " << pMother->GetName() << G4endl
00071             << "     Replicated volume: " << pName;
00072     G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002",
00073                 FatalException, message);
00074     return;
00075   }
00076   CheckAndSetParameters (pAxis, nReplicas, width, offset);
00077 }
00078 
00079 G4PVReplica::G4PVReplica( const G4String& pName,
00080                                 G4LogicalVolume* pLogical,
00081                                 G4LogicalVolume* pMotherLogical,
00082                           const EAxis pAxis,
00083                           const G4int nReplicas,
00084                           const G4double width,
00085                           const G4double offset )
00086   : G4VPhysicalVolume(0,G4ThreeVector(),pName,pLogical,0),
00087     fcopyNo(-1), fRegularVolsId(0)
00088 {
00089   if (!pMotherLogical)
00090   {
00091     std::ostringstream message;
00092     message << "NULL pointer specified as mother volume for "
00093             << pName << ".";
00094     G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002",
00095                 FatalException, message);
00096     return;
00097   }
00098   if (pLogical == pMotherLogical)
00099   {
00100     G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002",
00101                 FatalException, "Cannot place a volume inside itself!");
00102     return;
00103   }
00104   pMotherLogical->AddDaughter(this);
00105   SetMotherLogical(pMotherLogical);
00106   if (pMotherLogical->GetNoDaughters() != 1)
00107   {
00108     std::ostringstream message;
00109     message << "Replica or parameterised volume must be the only daughter !"
00110             << G4endl
00111             << "     Mother logical volume: " << pMotherLogical->GetName()
00112             << G4endl
00113             << "     Replicated volume: " << pName;
00114     G4Exception("G4PVReplica::G4PVReplica()", "GeomVol0002",
00115                 FatalException, message);
00116     return;
00117   }
00118   CheckAndSetParameters (pAxis, nReplicas, width, offset);
00119 }
00120 
00121 void G4PVReplica::CheckAndSetParameters( const EAxis pAxis,
00122                                          const G4int nReplicas,
00123                                          const G4double width,
00124                                          const G4double offset)
00125 {    
00126   if (nReplicas<1)
00127   {
00128     G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0002",
00129                 FatalException, "Illegal number of replicas.");
00130   }
00131   fnReplicas=nReplicas;
00132   if (width<0)
00133   {
00134     G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0002",
00135                 FatalException, "Width must be positive.");
00136   }
00137   fwidth  = width;
00138   foffset = offset;
00139   faxis   = pAxis;
00140 
00141   
00142   
00143   G4RotationMatrix* pRMat=0;
00144   switch (faxis)
00145   {
00146     case kPhi:
00147       pRMat=new G4RotationMatrix();
00148       if (!pRMat)
00149       {
00150         G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0003",
00151                     FatalException, "Rotation matrix allocation failed.");
00152       }
00153       SetRotation(pRMat);
00154       break;
00155     case kRho:
00156     case kXAxis:
00157     case kYAxis:
00158     case kZAxis:
00159     case kUndefined:
00160       break;
00161     default:
00162       G4Exception("G4PVReplica::CheckAndSetParameters()", "GeomVol0002",
00163                   FatalException, "Unknown axis of replication.");
00164       break;
00165   }
00166 }
00167 
00168 G4PVReplica::G4PVReplica( __void__& a )
00169   : G4VPhysicalVolume(a), faxis(kZAxis), fnReplicas(0), fwidth(0.),
00170     foffset(0.), fcopyNo(-1), fRegularStructureCode(0), fRegularVolsId(0)
00171 {
00172 }
00173 
00174 G4PVReplica::~G4PVReplica()
00175 {
00176   if ( faxis==kPhi )
00177   {
00178     delete GetRotation();
00179   }
00180 }
00181 
00182 G4bool G4PVReplica::IsMany() const
00183 {
00184   return false; 
00185 }
00186 
00187 G4int G4PVReplica::GetCopyNo() const
00188 {
00189   return fcopyNo;
00190 }
00191 
00192 void  G4PVReplica::SetCopyNo(G4int newCopyNo)
00193 {
00194   fcopyNo = newCopyNo;
00195 }
00196 
00197 G4bool G4PVReplica::IsReplicated() const
00198 {
00199   return true;
00200 }
00201 
00202 G4bool G4PVReplica::IsParameterised() const
00203 {
00204   return false;
00205 }
00206 
00207 G4VPVParameterisation* G4PVReplica::GetParameterisation() const
00208 {
00209   return 0;
00210 }
00211 
00212 G4int G4PVReplica::GetMultiplicity() const
00213 {
00214   return fnReplicas;
00215 }
00216 
00217 
00218 
00219 void G4PVReplica::GetReplicationData( EAxis& axis,
00220                                       G4int& nReplicas,
00221                                       G4double& width,
00222                                       G4double& offset,
00223                                       G4bool& consuming ) const
00224 {
00225   axis = faxis;
00226   nReplicas = fnReplicas;
00227   width = fwidth;
00228   offset = foffset;
00229   consuming = true;
00230 }
00231 
00232 G4bool G4PVReplica::IsRegularStructure() const
00233 {
00234   return (fRegularVolsId!=0); 
00235 }
00236 
00237 G4int  G4PVReplica::GetRegularStructureId() const
00238 {
00239   return fRegularVolsId; 
00240 }
00241 
00242 void   G4PVReplica::SetRegularStructureId( G4int Code )
00243 {
00244   fRegularVolsId= Code; 
00245 }