BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSFieldQueryForVis.hh
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#ifndef BDSFIELDQUERYFORVIS_H
20#define BDSFIELDQUERYFORVIS_H
21
22#include "BDSFieldQuery.hh"
23
24#include "G4ThreeVector.hh"
25#include "G4Types.hh"
26
27#include <array>
28#include <vector>
29
31
38{
39public:
41 virtual ~BDSFieldQueryForVis();
42
43 inline const std::vector<std::array<G4double, 9>> Values() const {return values;}
44 inline G4double MaxFieldE() const {return maxFieldE;}
45 inline G4double MaxFieldB() const {return maxFieldB;}
46
48 virtual void QueryField(const BDSFieldQueryInfo* query);
49
51 virtual void CleanUp();
52
54 static G4double SimpleMag(G4double x, G4double y, G4double z);
55
57 virtual void OpenFiles(const BDSFieldQueryInfo*){;}
58 virtual void WriteHeader(std::ofstream&, const BDSFieldQueryInfo*) const {;}
59 virtual void CloseFiles() {;}
60 virtual void WriteFieldValue(const G4ThreeVector&,
61 G4double,
62 const G4double[6]) {;}
63 virtual void PrintBAndEInfo(const BDSFieldQueryInfo*) const {;}
65
66protected:
68 virtual void GetFieldValue(const G4ThreeVector& globalXYZ,
69 const G4ThreeVector& globalDirection,
70 G4double tGlobal,
71 G4double fieldValue[6]);
72
73private:
74 std::vector<std::array<G4double, 9>> values;
75 G4double maxFieldB;
76 G4double maxFieldE;
77 G4bool drawZeroValuePoints;
78};
79
80#endif
Wrapper class of BDSFieldQuery that accumulates values for visualisation.
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 CloseFiles()
Replace with empty implementation so no file handling is done for visualisation.
virtual void WriteFieldValue(const G4ThreeVector &, G4double, const G4double[6])
Replace with empty implementation so no file handling is done for visualisation.
virtual void CleanUp()
Call the base class method and also clear the vector of values in this class.
virtual void WriteHeader(std::ofstream &, const BDSFieldQueryInfo *) const
Replace with empty implementation so no file handling is done for visualisation.
virtual void PrintBAndEInfo(const BDSFieldQueryInfo *) const
Replace with empty implementation so no file handling is done for visualisation.
virtual void OpenFiles(const BDSFieldQueryInfo *)
Replace with empty implementation so no file handling is done for visualisation.
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.
Class for querying the Geant4 model for field at any point.