[831] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
| 26 | // |
---|
| 27 | // $Id: G4ParameterisationTrd.cc,v 1.14 2006/06/29 18:18:48 gunter Exp $ |
---|
[850] | 28 | // GEANT4 tag $Name: HEAD $ |
---|
[831] | 29 | // |
---|
| 30 | // class G4ParameterisationTrd Implementation file |
---|
| 31 | // |
---|
| 32 | // 26.05.03 - P.Arce, Initial version |
---|
| 33 | // 08.04.04 - I.Hrivnacova, Implemented reflection |
---|
| 34 | // -------------------------------------------------------------------- |
---|
| 35 | |
---|
| 36 | #include "G4ParameterisationTrd.hh" |
---|
| 37 | |
---|
| 38 | #include <iomanip> |
---|
| 39 | #include "G4ThreeVector.hh" |
---|
| 40 | #include "G4RotationMatrix.hh" |
---|
| 41 | #include "G4VPhysicalVolume.hh" |
---|
| 42 | #include "G4LogicalVolume.hh" |
---|
| 43 | #include "G4ReflectedSolid.hh" |
---|
| 44 | #include "G4Trd.hh" |
---|
| 45 | #include "G4Trap.hh" |
---|
| 46 | |
---|
| 47 | //-------------------------------------------------------------------------- |
---|
| 48 | G4VParameterisationTrd:: |
---|
| 49 | G4VParameterisationTrd( EAxis axis, G4int nDiv, G4double width, |
---|
| 50 | G4double offset, G4VSolid* msolid, |
---|
| 51 | DivisionType divType ) |
---|
| 52 | : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid ) |
---|
| 53 | { |
---|
| 54 | G4Trd* msol = (G4Trd*)(msolid); |
---|
| 55 | if (msolid->GetEntityType() == "G4ReflectedSolid") |
---|
| 56 | { |
---|
| 57 | // Get constituent solid |
---|
| 58 | G4VSolid* mConstituentSolid |
---|
| 59 | = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid(); |
---|
| 60 | msol = (G4Trd*)(mConstituentSolid); |
---|
| 61 | |
---|
| 62 | // Create a new solid with inversed parameters |
---|
| 63 | G4Trd* newSolid |
---|
| 64 | = new G4Trd(msol->GetName(), |
---|
| 65 | msol->GetXHalfLength2(), msol->GetXHalfLength1(), |
---|
| 66 | msol->GetYHalfLength2(), msol->GetYHalfLength1(), |
---|
| 67 | msol->GetZHalfLength()); |
---|
| 68 | msol = newSolid; |
---|
| 69 | fmotherSolid = newSolid; |
---|
| 70 | fReflectedSolid = true; |
---|
| 71 | fDeleteSolid = true; |
---|
| 72 | } |
---|
| 73 | } |
---|
| 74 | |
---|
| 75 | //------------------------------------------------------------------------ |
---|
| 76 | G4VParameterisationTrd::~G4VParameterisationTrd() |
---|
| 77 | { |
---|
| 78 | } |
---|
| 79 | |
---|
| 80 | //------------------------------------------------------------------------ |
---|
| 81 | G4ParameterisationTrdX:: |
---|
| 82 | G4ParameterisationTrdX( EAxis axis, G4int nDiv, |
---|
| 83 | G4double width, G4double offset, |
---|
| 84 | G4VSolid* msolid, DivisionType divType ) |
---|
| 85 | : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType ) |
---|
| 86 | { |
---|
| 87 | CheckParametersValidity(); |
---|
| 88 | SetType( "DivisionTrdX" ); |
---|
| 89 | |
---|
| 90 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 91 | if( divType == DivWIDTH ) |
---|
| 92 | { |
---|
| 93 | fnDiv = CalculateNDiv( 2*msol->GetXHalfLength1(), |
---|
| 94 | width, offset ); |
---|
| 95 | } |
---|
| 96 | else if( divType == DivNDIV ) |
---|
| 97 | { |
---|
| 98 | fwidth = CalculateWidth( 2*msol->GetXHalfLength1(), nDiv, offset ); |
---|
| 99 | } |
---|
| 100 | |
---|
| 101 | #ifdef G4DIVDEBUG |
---|
| 102 | if( verbose >= -1 ) |
---|
| 103 | { |
---|
| 104 | G4cout << " G4ParameterisationTrdX - ## divisions " << fnDiv << " = " |
---|
| 105 | << nDiv << G4endl |
---|
| 106 | << " Offset " << foffset << " = " << offset << G4endl |
---|
| 107 | << " Width " << fwidth << " = " << width << G4endl; |
---|
| 108 | } |
---|
| 109 | #endif |
---|
| 110 | } |
---|
| 111 | |
---|
| 112 | //------------------------------------------------------------------------ |
---|
| 113 | G4ParameterisationTrdX::~G4ParameterisationTrdX() |
---|
| 114 | { |
---|
| 115 | } |
---|
| 116 | |
---|
| 117 | //------------------------------------------------------------------------ |
---|
| 118 | G4double G4ParameterisationTrdX::GetMaxParameter() const |
---|
| 119 | { |
---|
| 120 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 121 | return 2*msol->GetXHalfLength1(); |
---|
| 122 | } |
---|
| 123 | |
---|
| 124 | //------------------------------------------------------------------------ |
---|
| 125 | void |
---|
| 126 | G4ParameterisationTrdX:: |
---|
| 127 | ComputeTransformation( const G4int copyNo, |
---|
| 128 | G4VPhysicalVolume *physVol ) const |
---|
| 129 | { |
---|
| 130 | G4Trd* msol = (G4Trd*)(fmotherSolid ); |
---|
| 131 | G4double mdx = msol->GetXHalfLength1(); |
---|
| 132 | |
---|
| 133 | //----- translation |
---|
| 134 | G4ThreeVector origin(0.,0.,0.); |
---|
| 135 | G4double posi = -mdx + foffset + (copyNo+0.5)*fwidth; |
---|
| 136 | if( faxis == kXAxis ) |
---|
| 137 | { |
---|
| 138 | origin.setX( posi ); |
---|
| 139 | } |
---|
| 140 | else |
---|
| 141 | { |
---|
| 142 | G4Exception("G4ParameterisationTrdX::ComputeTransformation()", |
---|
| 143 | "IllegalConstruct", FatalException, |
---|
| 144 | "Only axes along X are allowed, and axis is: "+faxis); |
---|
| 145 | } |
---|
| 146 | |
---|
| 147 | #ifdef G4DIVDEBUG |
---|
| 148 | if( verbose >= -2 ) |
---|
| 149 | { |
---|
| 150 | G4cout << std::setprecision(8) << " G4ParameterisationTrdX: " |
---|
| 151 | << copyNo << G4endl |
---|
| 152 | << " Position: " << origin << " - Axis: " << faxis << G4endl; |
---|
| 153 | } |
---|
| 154 | #endif |
---|
| 155 | |
---|
| 156 | //----- set translation |
---|
| 157 | physVol->SetTranslation( origin ); |
---|
| 158 | } |
---|
| 159 | |
---|
| 160 | //-------------------------------------------------------------------------- |
---|
| 161 | void |
---|
| 162 | G4ParameterisationTrdX:: |
---|
| 163 | ComputeDimensions( G4Trd& trd, const G4int, const G4VPhysicalVolume* ) const |
---|
| 164 | { |
---|
| 165 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 166 | |
---|
| 167 | G4double pDy1 = msol->GetYHalfLength1(); |
---|
| 168 | G4double pDy2 = msol->GetYHalfLength2(); |
---|
| 169 | G4double pDz = msol->GetZHalfLength(); |
---|
| 170 | G4double pDx = fwidth/2.; |
---|
| 171 | |
---|
| 172 | trd.SetAllParameters ( pDx, pDx, pDy1, pDy2, pDz ); |
---|
| 173 | |
---|
| 174 | #ifdef G4DIVDEBUG |
---|
| 175 | if( verbose >= -2 ) |
---|
| 176 | { |
---|
| 177 | G4cout << " G4ParameterisationTrdX::ComputeDimensions():" << G4endl; |
---|
| 178 | trd.DumpInfo(); |
---|
| 179 | } |
---|
| 180 | #endif |
---|
| 181 | } |
---|
| 182 | |
---|
| 183 | //-------------------------------------------------------------------------- |
---|
| 184 | void G4ParameterisationTrdX::CheckParametersValidity() |
---|
| 185 | { |
---|
| 186 | G4VDivisionParameterisation::CheckParametersValidity(); |
---|
| 187 | |
---|
| 188 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 189 | |
---|
| 190 | G4double mpDx1 = msol->GetXHalfLength1(); |
---|
| 191 | G4double mpDx2 = msol->GetXHalfLength2(); |
---|
| 192 | |
---|
| 193 | if( std::fabs(mpDx1 - mpDx2) > kCarTolerance ) |
---|
| 194 | { |
---|
| 195 | G4cerr << "ERROR - G4ParameterisationTrdX::CheckParametersValidity()" |
---|
| 196 | << G4endl |
---|
| 197 | << " Making a division of a TRD along axis X," << G4endl |
---|
| 198 | << " while the X half lengths are not equal," << G4endl |
---|
| 199 | << " is not (yet) supported. It will result" << G4endl |
---|
| 200 | << " in non-equal division solids." << G4endl; |
---|
| 201 | G4Exception("G4ParameterisationTrdX::CheckParametersValidity()", |
---|
| 202 | "IllegalConstruct", FatalException, |
---|
| 203 | "Invalid solid specification. NOT supported."); |
---|
| 204 | } |
---|
| 205 | } |
---|
| 206 | |
---|
| 207 | //-------------------------------------------------------------------------- |
---|
| 208 | G4ParameterisationTrdY:: |
---|
| 209 | G4ParameterisationTrdY( EAxis axis, G4int nDiv, |
---|
| 210 | G4double width, G4double offset, |
---|
| 211 | G4VSolid* msolid, DivisionType divType) |
---|
| 212 | : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType ) |
---|
| 213 | { |
---|
| 214 | CheckParametersValidity(); |
---|
| 215 | SetType( "DivisionTrdY" ); |
---|
| 216 | |
---|
| 217 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 218 | if( divType == DivWIDTH ) |
---|
| 219 | { |
---|
| 220 | fnDiv = CalculateNDiv( 2*msol->GetYHalfLength1(), |
---|
| 221 | width, offset ); |
---|
| 222 | } |
---|
| 223 | else if( divType == DivNDIV ) |
---|
| 224 | { |
---|
| 225 | fwidth = CalculateWidth( 2*msol->GetYHalfLength1(), |
---|
| 226 | nDiv, offset ); |
---|
| 227 | } |
---|
| 228 | |
---|
| 229 | #ifdef G4DIVDEBUG |
---|
| 230 | if( verbose >= 1 ) |
---|
| 231 | { |
---|
| 232 | G4cout << " G4ParameterisationTrdY no divisions " << fnDiv |
---|
| 233 | << " = " << nDiv << G4endl |
---|
| 234 | << " Offset " << foffset << " = " << offset << G4endl |
---|
| 235 | << " width " << fwidth << " = " << width << G4endl; |
---|
| 236 | } |
---|
| 237 | #endif |
---|
| 238 | } |
---|
| 239 | |
---|
| 240 | //------------------------------------------------------------------------ |
---|
| 241 | G4ParameterisationTrdY::~G4ParameterisationTrdY() |
---|
| 242 | { |
---|
| 243 | } |
---|
| 244 | |
---|
| 245 | //------------------------------------------------------------------------ |
---|
| 246 | G4double G4ParameterisationTrdY::GetMaxParameter() const |
---|
| 247 | { |
---|
| 248 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 249 | return 2*msol->GetYHalfLength1(); |
---|
| 250 | } |
---|
| 251 | |
---|
| 252 | //-------------------------------------------------------------------------- |
---|
| 253 | void |
---|
| 254 | G4ParameterisationTrdY:: |
---|
| 255 | ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol ) const |
---|
| 256 | { |
---|
| 257 | G4Trd* msol = (G4Trd*)(fmotherSolid ); |
---|
| 258 | G4double mdy = msol->GetYHalfLength1(); |
---|
| 259 | |
---|
| 260 | //----- translation |
---|
| 261 | G4ThreeVector origin(0.,0.,0.); |
---|
| 262 | G4double posi = -mdy + foffset + (copyNo+0.5)*fwidth; |
---|
| 263 | if( faxis == kYAxis ) |
---|
| 264 | { |
---|
| 265 | origin.setY( posi ); |
---|
| 266 | } |
---|
| 267 | else |
---|
| 268 | { |
---|
| 269 | G4Exception("G4ParameterisationTrdY::ComputeTransformation()", |
---|
| 270 | "IllegalConstruct", FatalException, |
---|
| 271 | "Only axes along Y are allowed !"); |
---|
| 272 | } |
---|
| 273 | |
---|
| 274 | #ifdef G4DIVDEBUG |
---|
| 275 | if( verbose >= 2 ) |
---|
| 276 | { |
---|
| 277 | G4cout << std::setprecision(8) << " G4ParameterisationTrdY " << copyNo |
---|
| 278 | << " pos " << origin << " rot mat " << " axis " << faxis << G4endl; |
---|
| 279 | } |
---|
| 280 | #endif |
---|
| 281 | |
---|
| 282 | //----- set translation |
---|
| 283 | physVol->SetTranslation( origin ); |
---|
| 284 | } |
---|
| 285 | |
---|
| 286 | //-------------------------------------------------------------------------- |
---|
| 287 | void |
---|
| 288 | G4ParameterisationTrdY:: |
---|
| 289 | ComputeDimensions(G4Trd& trd, const G4int, const G4VPhysicalVolume*) const |
---|
| 290 | { |
---|
| 291 | //---- The division along Y of a Trd will result a Trd, only |
---|
| 292 | //--- if Y at -Z and +Z are equal, else use the G4Trap version |
---|
| 293 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 294 | |
---|
| 295 | G4double pDx1 = msol->GetXHalfLength1(); |
---|
| 296 | G4double pDx2 = msol->GetXHalfLength2(); |
---|
| 297 | G4double pDz = msol->GetZHalfLength(); |
---|
| 298 | G4double pDy = fwidth/2.; |
---|
| 299 | |
---|
| 300 | trd.SetAllParameters ( pDx1, pDx2, pDy, pDy, pDz ); |
---|
| 301 | |
---|
| 302 | #ifdef G4DIVDEBUG |
---|
| 303 | if( verbose >= 2 ) |
---|
| 304 | { |
---|
| 305 | G4cout << " G4ParameterisationTrdY::ComputeDimensions():" << G4endl; |
---|
| 306 | trd.DumpInfo(); |
---|
| 307 | } |
---|
| 308 | #endif |
---|
| 309 | } |
---|
| 310 | |
---|
| 311 | //-------------------------------------------------------------------------- |
---|
| 312 | void G4ParameterisationTrdY::CheckParametersValidity() |
---|
| 313 | { |
---|
| 314 | G4VDivisionParameterisation::CheckParametersValidity(); |
---|
| 315 | |
---|
| 316 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 317 | |
---|
| 318 | G4double mpDy1 = msol->GetYHalfLength1(); |
---|
| 319 | G4double mpDy2 = msol->GetYHalfLength2(); |
---|
| 320 | |
---|
| 321 | if( std::fabs(mpDy1 - mpDy2) > kCarTolerance ) |
---|
| 322 | { |
---|
| 323 | G4cerr << "ERROR - G4ParameterisationTrdY::CheckParametersValidity()" |
---|
| 324 | << G4endl |
---|
| 325 | << " Making a division of a TRD along axis Y while" << G4endl |
---|
| 326 | << " the Y half lengths are not equal is not (yet)" << G4endl |
---|
| 327 | << " supported. It will result in non-equal" << G4endl |
---|
| 328 | << " division solids." << G4endl; |
---|
| 329 | G4Exception("G4ParameterisationTrdY::CheckParametersValidity()", |
---|
| 330 | "IllegalConstruct", FatalException, |
---|
| 331 | "Invalid solid specification. NOT supported."); |
---|
| 332 | } |
---|
| 333 | } |
---|
| 334 | |
---|
| 335 | //-------------------------------------------------------------------------- |
---|
| 336 | G4ParameterisationTrdZ:: |
---|
| 337 | G4ParameterisationTrdZ( EAxis axis, G4int nDiv, |
---|
| 338 | G4double width, G4double offset, |
---|
| 339 | G4VSolid* msolid, DivisionType divType ) |
---|
| 340 | : G4VParameterisationTrd( axis, nDiv, width, offset, msolid, divType ) |
---|
| 341 | { |
---|
| 342 | CheckParametersValidity(); |
---|
| 343 | SetType( "DivTrdZ" ); |
---|
| 344 | |
---|
| 345 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 346 | if( divType == DivWIDTH ) |
---|
| 347 | { |
---|
| 348 | fnDiv = CalculateNDiv( 2*msol->GetZHalfLength(), |
---|
| 349 | width, offset ); |
---|
| 350 | } |
---|
| 351 | else if( divType == DivNDIV ) |
---|
| 352 | { |
---|
| 353 | fwidth = CalculateWidth( 2*msol->GetZHalfLength(), |
---|
| 354 | nDiv, offset ); |
---|
| 355 | } |
---|
| 356 | |
---|
| 357 | #ifdef G4DIVDEBUG |
---|
| 358 | if( verbose >= 1 ) |
---|
| 359 | { |
---|
| 360 | G4cout << " G4ParameterisationTrdZ no divisions " << fnDiv |
---|
| 361 | << " = " << nDiv << G4endl |
---|
| 362 | << " Offset " << foffset << " = " << offset << G4endl |
---|
| 363 | << " Width " << fwidth << " = " << width << G4endl; |
---|
| 364 | } |
---|
| 365 | #endif |
---|
| 366 | } |
---|
| 367 | |
---|
| 368 | //------------------------------------------------------------------------ |
---|
| 369 | G4ParameterisationTrdZ::~G4ParameterisationTrdZ() |
---|
| 370 | { |
---|
| 371 | } |
---|
| 372 | |
---|
| 373 | //------------------------------------------------------------------------ |
---|
| 374 | G4double G4ParameterisationTrdZ::GetMaxParameter() const |
---|
| 375 | { |
---|
| 376 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 377 | return 2*msol->GetZHalfLength(); |
---|
| 378 | } |
---|
| 379 | |
---|
| 380 | //-------------------------------------------------------------------------- |
---|
| 381 | void |
---|
| 382 | G4ParameterisationTrdZ:: |
---|
| 383 | ComputeTransformation(const G4int copyNo, G4VPhysicalVolume *physVol) const |
---|
| 384 | { |
---|
| 385 | G4Trd* msol = (G4Trd*)(fmotherSolid ); |
---|
| 386 | G4double mdz = msol->GetZHalfLength(); |
---|
| 387 | |
---|
| 388 | //----- translation |
---|
| 389 | G4ThreeVector origin(0.,0.,0.); |
---|
| 390 | G4double posi = -mdz + OffsetZ() + (copyNo+0.5)*fwidth; |
---|
| 391 | if( faxis == kZAxis ) |
---|
| 392 | { |
---|
| 393 | origin.setZ( posi ); |
---|
| 394 | } |
---|
| 395 | else |
---|
| 396 | { |
---|
| 397 | G4Exception("G4ParameterisationTrdZ::ComputeTransformation()", |
---|
| 398 | "IllegalConstruct", FatalException, |
---|
| 399 | "Only axes along Z are allowed !"); |
---|
| 400 | } |
---|
| 401 | |
---|
| 402 | #ifdef G4DIVDEBUG |
---|
| 403 | if( verbose >= 1 ) |
---|
| 404 | { |
---|
| 405 | G4cout << std::setprecision(8) << " G4ParameterisationTrdZ: " |
---|
| 406 | << copyNo << G4endl |
---|
| 407 | << " Position: " << origin << " - Offset: " << foffset |
---|
| 408 | << " - Width: " << fwidth << " Axis " << faxis << G4endl; |
---|
| 409 | } |
---|
| 410 | #endif |
---|
| 411 | |
---|
| 412 | //----- set translation |
---|
| 413 | physVol->SetTranslation( origin ); |
---|
| 414 | } |
---|
| 415 | |
---|
| 416 | //-------------------------------------------------------------------------- |
---|
| 417 | void |
---|
| 418 | G4ParameterisationTrdZ:: |
---|
| 419 | ComputeDimensions(G4Trd& trd, const G4int copyNo, |
---|
| 420 | const G4VPhysicalVolume*) const |
---|
| 421 | { |
---|
| 422 | //---- The division along Z of a Trd will result a Trd |
---|
| 423 | G4Trd* msol = (G4Trd*)(fmotherSolid); |
---|
| 424 | |
---|
| 425 | G4double pDx1 = msol->GetXHalfLength1(); |
---|
| 426 | G4double DDx = (msol->GetXHalfLength2() - msol->GetXHalfLength1() ); |
---|
| 427 | G4double pDy1 = msol->GetYHalfLength1(); |
---|
| 428 | G4double DDy = (msol->GetYHalfLength2() - msol->GetYHalfLength1() ); |
---|
| 429 | G4double pDz = fwidth/2.; |
---|
| 430 | G4double zLength = 2*msol->GetZHalfLength(); |
---|
| 431 | |
---|
| 432 | trd.SetAllParameters( pDx1+DDx*(OffsetZ()+copyNo*fwidth)/zLength, |
---|
| 433 | pDx1+DDx*(OffsetZ()+(copyNo+1)*fwidth)/zLength, |
---|
| 434 | pDy1+DDy*(OffsetZ()+copyNo*fwidth)/zLength, |
---|
| 435 | pDy1+DDy*(OffsetZ()+(copyNo+1)*fwidth)/zLength, pDz ); |
---|
| 436 | |
---|
| 437 | #ifdef G4DIVDEBUG |
---|
| 438 | if( verbose >= 1 ) |
---|
| 439 | { |
---|
| 440 | G4cout << " G4ParameterisationTrdZ::ComputeDimensions()" |
---|
| 441 | << " - Mother TRD " << G4endl; |
---|
| 442 | msol->DumpInfo(); |
---|
| 443 | G4cout << " - Parameterised TRD: " |
---|
| 444 | << copyNo << G4endl; |
---|
| 445 | trd.DumpInfo(); |
---|
| 446 | } |
---|
| 447 | #endif |
---|
| 448 | } |
---|