00001 // 00002 // ******************************************************************** 00003 // * License and Disclaimer * 00004 // * * 00005 // * The Geant4 software is copyright of the Copyright Holders of * 00006 // * the Geant4 Collaboration. It is provided under the terms and * 00007 // * conditions of the Geant4 Software License, included in the file * 00008 // * LICENSE and available at http://cern.ch/geant4/license . These * 00009 // * include a list of copyright holders. * 00010 // * * 00011 // * Neither the authors of this software system, nor their employing * 00012 // * institutes,nor the agencies providing financial support for this * 00013 // * work make any representation or warranty, express or implied, * 00014 // * regarding this software system or assume any liability for its * 00015 // * use. Please see the license in the file LICENSE and URL above * 00016 // * for the full disclaimer and the limitation of liability. * 00017 // * * 00018 // * This code implementation is the result of the scientific and * 00019 // * technical work of the GEANT4 collaboration. * 00020 // * By using, copying, modifying or distributing the software (or * 00021 // * any work based on the software) you agree to acknowledge its * 00022 // * use in resulting scientific publications, and indicate your * 00023 // * acceptance of all terms of the Geant4 Software license. * 00024 // ******************************************************************** 00025 // 00026 // 00027 // $Id: G4LineSection.cc 69786 2013-05-15 09:38:51Z gcosmo $ 00028 // 00029 // -------------------------------------------------------------------- 00030 00031 #include "G4LineSection.hh" 00032 00033 G4LineSection::G4LineSection( const G4ThreeVector& PntA, 00034 const G4ThreeVector& PntB ) 00035 : EndpointA(PntA), 00036 VecAtoB(PntB-PntA) 00037 { 00038 fABdistanceSq = VecAtoB.mag2() ; 00039 /* 00040 G4double distABsquared = VecAtoB.mag2() ; 00041 if ( distABsquared == 0.0) 00042 { 00043 G4Exception("G4LineSection::G4LineSection()", "GeomField0002", 00044 FatalException, "Equal points in input (line->point) ?") ; 00045 } 00046 else 00047 { 00048 inverse_square_distAB=1.0 / distABsquared ; 00049 } 00050 */ 00051 } 00052 00053 G4double G4LineSection::Dist( G4ThreeVector OtherPnt ) const 00054 { 00055 G4double dist_sq; 00056 G4ThreeVector VecAZ; 00057 G4double sq_VecAZ, inner_prod, unit_projection ; 00058 00059 VecAZ= OtherPnt - EndpointA; 00060 sq_VecAZ = VecAZ.mag2(); 00061 00062 inner_prod= VecAtoB.dot( VecAZ ); 00063 00064 // Determine Projection(AZ on AB) / Length(AB) 00065 // 00066 if( fABdistanceSq != 0.0 ) 00067 { 00068 // unit_projection= inner_prod * InvsqDistAB(); 00069 unit_projection = inner_prod/fABdistanceSq; 00070 00071 if( (0. <= unit_projection ) && (unit_projection <= 1.0 ) ) 00072 { 00073 dist_sq= sq_VecAZ - unit_projection * inner_prod ; 00074 } 00075 else 00076 { 00077 // The perpendicular from the point to the line AB meets the line 00078 // in a point outside the line segment! 00079 00080 if( unit_projection < 0. ) // A is the closest point 00081 { 00082 dist_sq= sq_VecAZ; 00083 } 00084 else // B is the closest point 00085 { 00086 G4ThreeVector EndpointB = EndpointA + VecAtoB; 00087 G4ThreeVector VecBZ = OtherPnt - EndpointB; 00088 dist_sq = VecBZ.mag2(); 00089 } 00090 } 00091 } 00092 else 00093 { 00094 dist_sq = (OtherPnt - EndpointA).mag2() ; 00095 } 00096 if( dist_sq < 0.0 ) dist_sq = 0.0 ; 00097 00098 return std::sqrt(dist_sq) ; 00099 } 00100 00101 G4double G4LineSection::Distline( const G4ThreeVector& OtherPnt, 00102 const G4ThreeVector& LinePntA, 00103 const G4ThreeVector& LinePntB ) 00104 { 00105 G4LineSection LineAB( LinePntA, LinePntB ); // Line from A to B 00106 00107 return LineAB.Dist( OtherPnt ); 00108 }