BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldQueryForVis.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSFieldQueryForVis.hh"
20#include "BDSFieldQueryInfo.hh"
21
22#include "G4ThreeVector.hh"
23#include "G4Types.hh"
24
25#include <algorithm>
26#include <array>
27#include <cmath>
28#include <limits>
29#include <vector>
30
31
32BDSFieldQueryForVis::BDSFieldQueryForVis():
33 maxFieldB(0),
34 maxFieldE(0),
35 drawZeroValuePoints(false)
36{;}
37
38BDSFieldQueryForVis::~BDSFieldQueryForVis()
39{;}
40
42{
44 values.clear();
45 maxFieldB = -std::numeric_limits<G4double>::max();
46 maxFieldE = -std::numeric_limits<G4double>::max();
47}
48
50{
51 G4int totalNPoints = query->xInfo.n * query->yInfo.n * query->zInfo.n * query->tInfo.n;
52 values.reserve(totalNPoints);
53 drawZeroValuePoints = query->drawZeroValuePoints;
55}
56
57void BDSFieldQueryForVis::GetFieldValue(const G4ThreeVector& globalXYZ,
58 const G4ThreeVector& globalDirection,
59 G4double tGlobal,
60 G4double fieldValue[6])
61{
62 BDSFieldQuery::GetFieldValue(globalXYZ, globalDirection, tGlobal, fieldValue);
63 G4double bFieldMag = SimpleMag(fieldValue[0], fieldValue[1], fieldValue[2]);
64 G4double eFieldMag = SimpleMag(fieldValue[3], fieldValue[4], fieldValue[5]);
65 maxFieldB = std::max(maxFieldB, bFieldMag);
66 maxFieldE = std::max(maxFieldE, eFieldMag);
67 if ( (bFieldMag == 0) && (eFieldMag == 0) && !drawZeroValuePoints)
68 {return;} // don't store 0 field points - we don't need to retain the structure or complete array shape
69 values.emplace_back( std::array<G4double, 9>({globalXYZ.x(), globalXYZ.y(), globalXYZ.z(),
70 fieldValue[0], fieldValue[1], fieldValue[2],
71 fieldValue[3], fieldValue[4], fieldValue[5]}) );
72}
73
74G4double BDSFieldQueryForVis::SimpleMag(G4double x, G4double y, G4double z)
75{
76 return std::sqrt(x*x + y*y + z*z);
77}
virtual void GetFieldValue(const G4ThreeVector &globalXYZ, const G4ThreeVector &globalDirection, G4double tGlobal, G4double fieldValue[6])
Call base class method but make a copy of the global position and fields in this class.
virtual void QueryField(const BDSFieldQueryInfo *query)
Reserve the vector size we need first, then proceed as normal. Avoids big copies.
virtual void CleanUp()
Call the base class method and also clear the vector of values in this class.
static G4double SimpleMag(G4double x, G4double y, G4double z)
Calculate magnitude of 3-vector without having to construct one.
Holder class for all information required for a field query.
virtual void CleanUp()
Reset any member variables used during a query. Closes any files if open.
virtual void QueryField(const BDSFieldQueryInfo *query)
Query the field in the Geant4 model according to information in query.
virtual void GetFieldValue(const G4ThreeVector &globalXYZ, const G4ThreeVector &globalDirection, G4double tGlobal, G4double fieldValue[6])