BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamlineEndPieceBuilder.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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 "BDSAcceleratorComponent.hh"
20#include "BDSBeamline.hh"
21#include "BDSBeamlineEndPieceBuilder.hh"
22#include "BDSBeamlineElement.hh"
23#include "BDSExtent.hh"
24#include "BDSGlobalConstants.hh"
25#include "BDSSimpleComponent.hh"
26#include "BDSTiltOffset.hh"
27#include "BDSUtilities.hh"
28
29#include "globals.hh" // geant4 types / globals
30#include "G4ThreeVector.hh"
31
33 const G4bool circularMachine)
34{
35 const G4double lengthSafetyLarge = BDSGlobalConstants::Instance()->LengthSafetyLarge();
36
37 // the beamline of end pieces to be placed.
38 BDSBeamline* endPieces = new BDSBeamline();
39 if (beamline->empty())
40 {return endPieces;}
41
42 // references to first and last item for checking
43 const BDSBeamlineElement* firstItem = beamline->GetFirstItem();
44 const BDSBeamlineElement* lastItem = beamline->GetLastItem();
45
46 // loop over main beamline
47 G4int currentIndex = 0;
48 for (auto element : *beamline)
49 {
50 const auto accComponent = element->GetAcceleratorComponent();
51 const auto accComponentTO = element->GetTiltOffset();
52 const auto endPieceBefore = accComponent->EndPieceBefore();
53 const auto endPieceAfter = accComponent->EndPieceAfter();
54
55 G4bool placeBefore = true;
56 G4bool placeAfter = true;
57 if (!endPieceBefore && !endPieceAfter)
58 {// no end piece on either side - skip this element
59 currentIndex++;
60 continue;
61 }
62 else if (!endPieceBefore)
63 {placeBefore = false;}
64 else if (!endPieceAfter)
65 {placeAfter = false;}
66
67 G4double requiredBeforeLength = 0; // required length for previous end piece
68 G4double requiredAfterLength = 0; // required length for next end piece
69 if (endPieceBefore)
70 {requiredBeforeLength = endPieceBefore->GetChordLength() + 2*lengthSafetyLarge;}
71 if (endPieceAfter)
72 {requiredAfterLength = endPieceAfter->GetChordLength() + 2*lengthSafetyLarge;}
73
74 if (placeBefore)
75 {
76 // search behind current element until we hit a previous component that's not a drift.
77 // ie work out the amount of drift space available. If we hit a previous element, check
78 // if it has an end piece and take that length out of the available space. if there's then
79 // available space for this end piece we place it.
80 G4bool keepGoing = true;
81 const BDSBeamlineElement* inspectedElement = element; // start with current element
82 G4double availableLength = 0;
83 G4double previousNonDriftEndPieceL = 0;
84 G4bool driftIsFirstItem = false;
85 BDSExtent endPieceInnerExtent = endPieceBefore->GetInnerExtent();
86 G4bool driftsAreTooBigBefore = false;
87 while (keepGoing)
88 {
89 inspectedElement = beamline->GetPrevious(inspectedElement);
90 if (inspectedElement)
91 { // there is a previous element - inspect it
92 G4String inspectedElementType = inspectedElement->GetType();
93 if (inspectedElementType == "drift" || inspectedElementType == "dipolefringe" || inspectedElementType == "element")
94 {// leave keepGoing true here to keep going
95 // check extents first
96 BDSExtent extPipe = inspectedElement->GetAcceleratorComponent()->GetExtent();
97 G4double tiltPipe = inspectedElement->GetTilt();
98 if (accComponentTO) // could be nullptr
99 {extPipe = extPipe.Tilted(accComponentTO->GetTilt() - tiltPipe);}
100 if (extPipe.TransverselyGreaterThan(endPieceInnerExtent))
101 {
102 keepGoing = false;
103 if (inspectedElementType != "element")
104 {driftsAreTooBigBefore = true;}
105 }
106 if (keepGoing)
107 {
108 availableLength += inspectedElement->GetChordLength();
109 driftIsFirstItem = inspectedElement == firstItem;
110 }
111 }
112 else
113 {
114 keepGoing = false;
115 auto previousNonDriftAccCompEndPiece = inspectedElement->GetAcceleratorComponent()->EndPieceAfter();
116 if (previousNonDriftAccCompEndPiece)
117 {previousNonDriftEndPieceL = previousNonDriftAccCompEndPiece->GetChordLength();}
118 } // not a drift - can't be used
119 }
120 else
121 {keepGoing = false;} // quit the loop - no previous element
122 }
123
124 // don't place if there isn't enough space for BOTH endpieces - the previous one and this one
125 if (requiredBeforeLength > (availableLength - previousNonDriftEndPieceL))
126 {placeBefore = false;}
127
128 // don't place if the drifts on either side are too big for the end pieces to fit around
129 if (driftsAreTooBigBefore)
130 {placeBefore = false;}
131
132 // place the last item anyway even if there technically isn't space - either too short
133 // a drift or the magnet is first in the line
134 if ( (element == firstItem) || (driftIsFirstItem && !driftsAreTooBigBefore) )
135 {placeBefore = true;}
136
137 // If it's a circular machine the end piece will overlap with the teleporter or
138 // terminator or previous element which we don't check against, so play it safe
139 // and don't build it.
140 if ( (element == firstItem) && circularMachine )
141 {placeBefore = false;}
142
143 // Now check if the coil volumes will overlap
144 // This check is only done on before, as something can always be naturally added
145 // afterwards - this just sets a precedence of after over before.
146 if (placeBefore && !endPieces->empty() ) // only look backwards if we have already placed an end piece
147 {
148 // endPieces is a beam line of only end pieces
149 BDSAcceleratorComponent* prevEndPiece = (endPieces->back())->GetAcceleratorComponent();
150 // OF = outgoing face, IF = incoming face
151 // outgoing is end of next element being added, incoming is front face previous already added
152 BDSExtent extOF = endPieceBefore->GetExtent();
153 BDSExtent extIF = prevEndPiece->GetExtent();
154 G4ThreeVector oFNormal = endPieceBefore->InputFaceNormal();
155 G4ThreeVector iFNormal = prevEndPiece->OutputFaceNormal();
156 G4double lastEndPieceL = prevEndPiece->GetChordLength();
157 G4double zSeparation = availableLength - requiredBeforeLength - lastEndPieceL;
158 G4bool willIntersect = BDS::WillIntersect(iFNormal, oFNormal, zSeparation, extIF, extOF);
159 if (willIntersect)
160 {placeBefore = false;}
161
162 /*
163 // 3d check
164 auto thisEl = beamline->ProvideEndPieceElementBefore(endPieceBefore,
165 currentIndex);
166
167 G4bool willIntersect3D = endPieces->back()->Overlaps(thisEl);
168 if (willIntersect3D)
169 {
170 placeBefore = false;
171 delete thisEl;
172 }
173 */
174 }
175 if (placeBefore)
176 { // provide a BDSBeamlineElement for the end piece w.r.t. the original beam line
177 auto beforeEl = beamline->ProvideEndPieceElementBefore(endPieceBefore,
178 currentIndex);
179 endPieces->AddBeamlineElement(beforeEl); // append to end piece beam line
180 }
181 }
182
183 if (placeAfter)
184 {
185 // search after current element until we hit another component that's not a drift.
186 // ie work out the amount of drift space available. If we hit another element, check
187 // if it has an end piece and take that length out of the available space. if
188 // there's then available space for this end piece we place it.
189 G4bool keepGoing = true;
190 const BDSBeamlineElement* inspectedElement = element; // start with current element
191 G4double availableLength = 0;
192 G4double nextNonDriftEndPieceL = 0;
193 G4bool driftIsLastItem = false;
194 BDSExtent endPieceInnerExtent = endPieceAfter->GetInnerExtent();
195 G4bool driftsAreTooBigAfter = false;
196 while (keepGoing)
197 {
198 inspectedElement = beamline->GetNext(inspectedElement);
199 if (inspectedElement)
200 { // there is a next element - inspect it
201 G4String inspectedElementType = inspectedElement->GetType();
202 if (inspectedElementType == "drift" || inspectedElementType == "dipolefringe" || inspectedElementType == "element")
203 {// leave keepGoing true here to keep going
204 // check extents first
205 BDSExtent extPipe = inspectedElement->GetAcceleratorComponent()->GetExtent();
206 G4double tiltPipe = inspectedElement->GetTilt();
207 if (accComponentTO) // could be nullptr
208 {extPipe = extPipe.Tilted(accComponentTO->GetTilt() - tiltPipe);}
209 if (extPipe.TransverselyGreaterThan(endPieceInnerExtent))
210 {
211 keepGoing = false;
212 if (inspectedElementType != "element")
213 {driftsAreTooBigAfter = true;}
214 }
215 if (keepGoing)
216 {
217 availableLength += inspectedElement->GetChordLength();
218 driftIsLastItem = inspectedElement == lastItem;
219 }
220 }
221 else
222 {
223 keepGoing = false;
224 auto nextNonDriftAccCompEndPiece = inspectedElement->GetAcceleratorComponent()->EndPieceBefore();
225 if (nextNonDriftAccCompEndPiece)
226 {nextNonDriftEndPieceL = nextNonDriftAccCompEndPiece->GetChordLength();}
227 } // not a drift - can't be used
228 }
229 else
230 {keepGoing = false;} // quit the loop - no previous element
231 }
232
233 // don't place if there isn't enough space for BOTH endpieces - this one and the upcoming one
234 if (requiredAfterLength > (availableLength - nextNonDriftEndPieceL))
235 {placeAfter = false;}
236
237 // don't place if the drifts on either side are too big for the end pieces to fit around
238 if (driftsAreTooBigAfter)
239 {placeAfter = false;}
240
241 // place the last item anyway even if there technically isn't space - either too short
242 // a drift or the magnet is last in the line
243 if ( (element == lastItem) || (driftIsLastItem && !driftsAreTooBigAfter) )
244 {placeAfter = true;}
245
246 if (placeAfter)
247 { // provide a BDSBeamlineElement for the end piece w.r.t. the original beam line
248 G4bool flip = endPieceAfter == endPieceBefore; // sometimes the after piece is a copy of the before one
249 auto afterEl = beamline->ProvideEndPieceElementAfter(endPieceAfter,
250 currentIndex,
251 flip); // flip it round for placement
252 endPieces->AddBeamlineElement(afterEl); // append to end piece beam line
253 }
254 }
255 currentIndex++; // increment iterator index on beamline
256 }
257
258 return endPieces;
259}
Abstract class that represents a component of an accelerator.
BDSSimpleComponent * EndPieceBefore() const
virtual G4double GetChordLength() const
G4ThreeVector OutputFaceNormal() const
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4double GetChordLength() const
Accessor.
BDSAcceleratorComponent * GetAcceleratorComponent() const
Accessor.
G4double GetTilt() const
Convenience accessor.
G4String GetType() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
const BDSBeamlineElement * GetPrevious(const BDSBeamlineElement *element) const
Definition: BDSBeamline.cc:686
BDSBeamlineElement * back() const
Return a reference to the last element.
Definition: BDSBeamline.hh:211
const BDSBeamlineElement * GetLastItem() const
Return a reference to the last element.
Definition: BDSBeamline.hh:126
BDSBeamlineElement * ProvideEndPieceElementAfter(BDSSimpleComponent *endPiece, G4int index, G4bool flip=false) const
Definition: BDSBeamline.cc:853
const BDSBeamlineElement * GetFirstItem() const
Return a reference to the first element.
Definition: BDSBeamline.hh:123
void AddBeamlineElement(BDSBeamlineElement *element)
Definition: BDSBeamline.cc:571
G4bool empty() const
Iterator mechanics.
Definition: BDSBeamline.hh:175
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
BDSExtent Tilted(G4double angle) const
Provide a new copy of this extent but rotated along Z by a given tilt angle.
Definition: BDSExtent.cc:121
G4bool TransverselyGreaterThan(const BDSExtent &r) const
Comparison operator for x,y only. Ignores z (length).
Definition: BDSExtent.hh:102
BDSExtent GetExtent() const
Accessor - see member for more info.
static BDSGlobalConstants * Instance()
Access method.
BDSBeamline * BuildEndPieceBeamline(const BDSBeamline *beamline, const G4bool circularMachine)
G4bool WillIntersect(const G4ThreeVector &incomingNormal, const G4ThreeVector &outgoingNormal, const G4double &zSeparation, const BDSExtent &incomingExtent, const BDSExtent &outgoingExtent)