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 | // $Id: G4BetheBlochModel.cc,v 1.24 2008/10/22 16:00:57 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class header file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4BetheBlochModel |
---|
35 | // |
---|
36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code |
---|
37 | // |
---|
38 | // Creation date: 03.01.2002 |
---|
39 | // |
---|
40 | // Modifications: |
---|
41 | // |
---|
42 | // 04-12-02 Fix problem of G4DynamicParticle constructor (V.Ivanchenko) |
---|
43 | // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko) |
---|
44 | // 27-01-03 Make models region aware (V.Ivanchenko) |
---|
45 | // 13-02-03 Add name (V.Ivanchenko) |
---|
46 | // 24-03-05 Add G4EmCorrections (V.Ivanchenko) |
---|
47 | // 11-04-05 Major optimisation of internal interfaces (V.Ivanchenko) |
---|
48 | // 11-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma) |
---|
49 | // 12-02-06 move G4LossTableManager::Instance()->EmCorrections() |
---|
50 | // in constructor (mma) |
---|
51 | // 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio, |
---|
52 | // CorrectionsAlongStep needed for ions(V.Ivanchenko) |
---|
53 | // |
---|
54 | // ------------------------------------------------------------------- |
---|
55 | // |
---|
56 | |
---|
57 | |
---|
58 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
59 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
60 | |
---|
61 | #include "G4BetheBlochModel.hh" |
---|
62 | #include "Randomize.hh" |
---|
63 | #include "G4Electron.hh" |
---|
64 | #include "G4LossTableManager.hh" |
---|
65 | #include "G4EmCorrections.hh" |
---|
66 | #include "G4ParticleChangeForLoss.hh" |
---|
67 | #include "G4NistManager.hh" |
---|
68 | |
---|
69 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
70 | |
---|
71 | using namespace std; |
---|
72 | |
---|
73 | G4BetheBlochModel::G4BetheBlochModel(const G4ParticleDefinition* p, |
---|
74 | const G4String& nam) |
---|
75 | : G4VEmModel(nam), |
---|
76 | particle(0), |
---|
77 | tlimit(DBL_MAX), |
---|
78 | twoln10(2.0*log(10.0)), |
---|
79 | bg2lim(0.0169), |
---|
80 | taulim(8.4146e-3), |
---|
81 | isIon(false), |
---|
82 | isInitialised(false) |
---|
83 | { |
---|
84 | fParticleChange = 0; |
---|
85 | if(p) SetParticle(p); |
---|
86 | theElectron = G4Electron::Electron(); |
---|
87 | corr = G4LossTableManager::Instance()->EmCorrections(); |
---|
88 | nist = G4NistManager::Instance(); |
---|
89 | SetLowEnergyLimit(2.0*MeV); |
---|
90 | } |
---|
91 | |
---|
92 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
93 | |
---|
94 | G4BetheBlochModel::~G4BetheBlochModel() |
---|
95 | {} |
---|
96 | |
---|
97 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
98 | |
---|
99 | G4double G4BetheBlochModel::MinEnergyCut(const G4ParticleDefinition*, |
---|
100 | const G4MaterialCutsCouple* couple) |
---|
101 | { |
---|
102 | return couple->GetMaterial()->GetIonisation()->GetMeanExcitationEnergy(); |
---|
103 | } |
---|
104 | |
---|
105 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
106 | |
---|
107 | void G4BetheBlochModel::Initialise(const G4ParticleDefinition* p, |
---|
108 | const G4DataVector&) |
---|
109 | { |
---|
110 | if (!particle) SetParticle(p); |
---|
111 | |
---|
112 | corrFactor = chargeSquare; |
---|
113 | |
---|
114 | if(!isInitialised) { |
---|
115 | isInitialised = true; |
---|
116 | |
---|
117 | if(!fParticleChange) { |
---|
118 | if (pParticleChange) { |
---|
119 | fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*> |
---|
120 | (pParticleChange); |
---|
121 | } else { |
---|
122 | fParticleChange = new G4ParticleChangeForLoss(); |
---|
123 | } |
---|
124 | } |
---|
125 | } |
---|
126 | } |
---|
127 | |
---|
128 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
129 | |
---|
130 | void G4BetheBlochModel::SetParticle(const G4ParticleDefinition* p) |
---|
131 | { |
---|
132 | if(particle != p) { |
---|
133 | particle = p; |
---|
134 | G4String pname = particle->GetParticleName(); |
---|
135 | if (particle->GetParticleType() == "nucleus" && |
---|
136 | pname != "deuteron" && pname != "triton") { |
---|
137 | isIon = true; |
---|
138 | } |
---|
139 | |
---|
140 | mass = particle->GetPDGMass(); |
---|
141 | spin = particle->GetPDGSpin(); |
---|
142 | G4double q = particle->GetPDGCharge()/eplus; |
---|
143 | chargeSquare = q*q; |
---|
144 | ratio = electron_mass_c2/mass; |
---|
145 | G4double magmom = particle->GetPDGMagneticMoment() |
---|
146 | *mass/(0.5*eplus*hbar_Planck*c_squared); |
---|
147 | magMoment2 = magmom*magmom - 1.0; |
---|
148 | formfact = 0.0; |
---|
149 | if(particle->GetLeptonNumber() == 0) { |
---|
150 | G4double x = 0.8426*GeV; |
---|
151 | if(spin == 0.0 && mass < GeV) {x = 0.736*GeV;} |
---|
152 | else if(mass > GeV) { |
---|
153 | x /= nist->GetZ13(mass/proton_mass_c2); |
---|
154 | // tlimit = 51.2*GeV*A13[iz]*A13[iz]; |
---|
155 | } |
---|
156 | formfact = 2.0*electron_mass_c2/(x*x); |
---|
157 | tlimit = 2.0/formfact; |
---|
158 | } |
---|
159 | } |
---|
160 | } |
---|
161 | |
---|
162 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
163 | |
---|
164 | G4double G4BetheBlochModel::GetChargeSquareRatio(const G4ParticleDefinition* p, |
---|
165 | const G4Material* mat, |
---|
166 | G4double kineticEnergy) |
---|
167 | { |
---|
168 | // this method is called only for ions |
---|
169 | G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,kineticEnergy); |
---|
170 | corrFactor = q2*corr->EffectiveChargeCorrection(p,mat,kineticEnergy); |
---|
171 | return corrFactor; |
---|
172 | } |
---|
173 | |
---|
174 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
175 | |
---|
176 | G4double G4BetheBlochModel::GetParticleCharge(const G4ParticleDefinition* p, |
---|
177 | const G4Material* mat, |
---|
178 | G4double kineticEnergy) |
---|
179 | { |
---|
180 | // this method is called only for ions |
---|
181 | return corr->GetParticleCharge(p,mat,kineticEnergy); |
---|
182 | } |
---|
183 | |
---|
184 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
185 | |
---|
186 | G4double |
---|
187 | G4BetheBlochModel::ComputeCrossSectionPerElectron(const G4ParticleDefinition* p, |
---|
188 | G4double kineticEnergy, |
---|
189 | G4double cutEnergy, |
---|
190 | G4double maxKinEnergy) |
---|
191 | { |
---|
192 | G4double cross = 0.0; |
---|
193 | G4double tmax = MaxSecondaryEnergy(p, kineticEnergy); |
---|
194 | G4double maxEnergy = min(tmax,maxKinEnergy); |
---|
195 | if(cutEnergy < maxEnergy) { |
---|
196 | |
---|
197 | G4double totEnergy = kineticEnergy + mass; |
---|
198 | G4double energy2 = totEnergy*totEnergy; |
---|
199 | G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2; |
---|
200 | |
---|
201 | cross = 1.0/cutEnergy - 1.0/maxEnergy |
---|
202 | - beta2*log(maxEnergy/cutEnergy)/tmax; |
---|
203 | |
---|
204 | // +term for spin=1/2 particle |
---|
205 | if( 0.5 == spin ) cross += 0.5*(maxEnergy - cutEnergy)/energy2; |
---|
206 | |
---|
207 | // High order correction different for hadrons and ions |
---|
208 | // nevetheless they are applied to reduce high energy transfers |
---|
209 | // if(!isIon) |
---|
210 | //cross += corr->FiniteSizeCorrectionXS(p,currentMaterial, |
---|
211 | // kineticEnergy,cutEnergy); |
---|
212 | |
---|
213 | cross *= twopi_mc2_rcl2*chargeSquare/beta2; |
---|
214 | } |
---|
215 | |
---|
216 | // G4cout << "BB: e= " << kineticEnergy << " tmin= " << cutEnergy |
---|
217 | // << " tmax= " << tmax << " cross= " << cross << G4endl; |
---|
218 | |
---|
219 | return cross; |
---|
220 | } |
---|
221 | |
---|
222 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
223 | |
---|
224 | G4double G4BetheBlochModel::ComputeCrossSectionPerAtom( |
---|
225 | const G4ParticleDefinition* p, |
---|
226 | G4double kineticEnergy, |
---|
227 | G4double Z, G4double, |
---|
228 | G4double cutEnergy, |
---|
229 | G4double maxEnergy) |
---|
230 | { |
---|
231 | G4double cross = Z*ComputeCrossSectionPerElectron |
---|
232 | (p,kineticEnergy,cutEnergy,maxEnergy); |
---|
233 | return cross; |
---|
234 | } |
---|
235 | |
---|
236 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
237 | |
---|
238 | G4double G4BetheBlochModel::CrossSectionPerVolume( |
---|
239 | const G4Material* material, |
---|
240 | const G4ParticleDefinition* p, |
---|
241 | G4double kineticEnergy, |
---|
242 | G4double cutEnergy, |
---|
243 | G4double maxEnergy) |
---|
244 | { |
---|
245 | currentMaterial = material; |
---|
246 | G4double eDensity = material->GetElectronDensity(); |
---|
247 | G4double cross = eDensity*ComputeCrossSectionPerElectron |
---|
248 | (p,kineticEnergy,cutEnergy,maxEnergy); |
---|
249 | return cross; |
---|
250 | } |
---|
251 | |
---|
252 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
253 | |
---|
254 | G4double G4BetheBlochModel::ComputeDEDXPerVolume(const G4Material* material, |
---|
255 | const G4ParticleDefinition* p, |
---|
256 | G4double kineticEnergy, |
---|
257 | G4double cut) |
---|
258 | { |
---|
259 | G4double tmax = MaxSecondaryEnergy(p, kineticEnergy); |
---|
260 | G4double cutEnergy = min(cut,tmax); |
---|
261 | |
---|
262 | G4double tau = kineticEnergy/mass; |
---|
263 | G4double gam = tau + 1.0; |
---|
264 | G4double bg2 = tau * (tau+2.0); |
---|
265 | G4double beta2 = bg2/(gam*gam); |
---|
266 | |
---|
267 | G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy(); |
---|
268 | G4double eexc2 = eexc*eexc; |
---|
269 | G4double cden = material->GetIonisation()->GetCdensity(); |
---|
270 | G4double mden = material->GetIonisation()->GetMdensity(); |
---|
271 | G4double aden = material->GetIonisation()->GetAdensity(); |
---|
272 | G4double x0den = material->GetIonisation()->GetX0density(); |
---|
273 | G4double x1den = material->GetIonisation()->GetX1density(); |
---|
274 | |
---|
275 | G4double eDensity = material->GetElectronDensity(); |
---|
276 | |
---|
277 | G4double dedx = log(2.0*electron_mass_c2*bg2*cutEnergy/eexc2) |
---|
278 | - (1.0 + cutEnergy/tmax)*beta2; |
---|
279 | |
---|
280 | if(0.5 == spin) { |
---|
281 | G4double del = 0.5*cutEnergy/(kineticEnergy + mass); |
---|
282 | dedx += del*del; |
---|
283 | } |
---|
284 | |
---|
285 | // density correction |
---|
286 | G4double x = log(bg2)/twoln10; |
---|
287 | if ( x >= x0den ) { |
---|
288 | dedx -= twoln10*x - cden ; |
---|
289 | if ( x < x1den ) dedx -= aden*pow((x1den-x),mden) ; |
---|
290 | } |
---|
291 | |
---|
292 | // shell correction |
---|
293 | dedx -= 2.0*corr->ShellCorrection(p,material,kineticEnergy); |
---|
294 | |
---|
295 | // now compute the total ionization loss |
---|
296 | |
---|
297 | if (dedx < 0.0) dedx = 0.0 ; |
---|
298 | |
---|
299 | dedx *= twopi_mc2_rcl2*chargeSquare*eDensity/beta2; |
---|
300 | |
---|
301 | //High order correction different for hadrons and ions |
---|
302 | if(isIon) { |
---|
303 | dedx += corr->IonBarkasCorrection(p,material,kineticEnergy); |
---|
304 | } else { |
---|
305 | dedx += corr->HighOrderCorrections(p,material,kineticEnergy,cutEnergy); |
---|
306 | } |
---|
307 | return dedx; |
---|
308 | } |
---|
309 | |
---|
310 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
311 | |
---|
312 | void G4BetheBlochModel::CorrectionsAlongStep(const G4MaterialCutsCouple* couple, |
---|
313 | const G4DynamicParticle* dp, |
---|
314 | G4double& eloss, |
---|
315 | G4double&, |
---|
316 | G4double length) |
---|
317 | { |
---|
318 | const G4ParticleDefinition* p = dp->GetDefinition(); |
---|
319 | const G4Material* mat = couple->GetMaterial(); |
---|
320 | G4double preKinEnergy = dp->GetKineticEnergy(); |
---|
321 | G4double e = preKinEnergy - eloss*0.5; |
---|
322 | if(e < 0.0) e = preKinEnergy*0.5; |
---|
323 | |
---|
324 | if(isIon) { |
---|
325 | G4double q2 = corr->EffectiveChargeSquareRatio(p,mat,e); |
---|
326 | GetModelOfFluctuations()->SetParticleAndCharge(p, q2); |
---|
327 | eloss *= q2*corr->EffectiveChargeCorrection(p,mat,e)/corrFactor; |
---|
328 | eloss += length*corr->IonHighOrderCorrections(p,couple,e); |
---|
329 | } |
---|
330 | |
---|
331 | if(nuclearStopping && preKinEnergy*proton_mass_c2/mass < chargeSquare*100.*MeV) { |
---|
332 | |
---|
333 | G4double nloss = length*corr->NuclearDEDX(p,mat,e,false); |
---|
334 | |
---|
335 | // too big energy loss |
---|
336 | if(eloss + nloss > preKinEnergy) { |
---|
337 | nloss *= (preKinEnergy/(eloss + nloss)); |
---|
338 | eloss = preKinEnergy; |
---|
339 | } else { |
---|
340 | eloss += nloss; |
---|
341 | } |
---|
342 | /* |
---|
343 | G4cout << "G4ionIonisation::CorrectionsAlongStep: e= " << preKinEnergy |
---|
344 | << " de= " << eloss << " NIEL= " << nloss |
---|
345 | << " dynQ= " << dp->GetCharge()/eplus << G4endl; |
---|
346 | */ |
---|
347 | fParticleChange->ProposeNonIonizingEnergyDeposit(nloss); |
---|
348 | } |
---|
349 | |
---|
350 | } |
---|
351 | |
---|
352 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
353 | |
---|
354 | void G4BetheBlochModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp, |
---|
355 | const G4MaterialCutsCouple*, |
---|
356 | const G4DynamicParticle* dp, |
---|
357 | G4double minKinEnergy, |
---|
358 | G4double maxEnergy) |
---|
359 | { |
---|
360 | G4double kineticEnergy = dp->GetKineticEnergy(); |
---|
361 | G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy); |
---|
362 | |
---|
363 | G4double maxKinEnergy = std::min(maxEnergy,tmax); |
---|
364 | if(minKinEnergy >= maxKinEnergy) return; |
---|
365 | |
---|
366 | G4double totEnergy = kineticEnergy + mass; |
---|
367 | G4double etot2 = totEnergy*totEnergy; |
---|
368 | G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2; |
---|
369 | |
---|
370 | G4double deltaKinEnergy, f; |
---|
371 | G4double f1 = 0.0; |
---|
372 | G4double fmax = 1.0; |
---|
373 | if( 0.5 == spin ) fmax += 0.5*maxKinEnergy*maxKinEnergy/etot2; |
---|
374 | |
---|
375 | // sampling without nuclear size effect |
---|
376 | do { |
---|
377 | G4double q = G4UniformRand(); |
---|
378 | deltaKinEnergy = minKinEnergy*maxKinEnergy |
---|
379 | /(minKinEnergy*(1.0 - q) + maxKinEnergy*q); |
---|
380 | |
---|
381 | f = 1.0 - beta2*deltaKinEnergy/tmax; |
---|
382 | if( 0.5 == spin ) { |
---|
383 | f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2; |
---|
384 | f += f1; |
---|
385 | } |
---|
386 | |
---|
387 | } while( fmax*G4UniformRand() > f); |
---|
388 | |
---|
389 | // projectile formfactor - suppresion of high energy |
---|
390 | // delta-electron production at high energy |
---|
391 | |
---|
392 | G4double x = formfact*deltaKinEnergy; |
---|
393 | if(x > 1.e-6) { |
---|
394 | |
---|
395 | G4double x1 = 1.0 + x; |
---|
396 | G4double g = 1.0/(x1*x1); |
---|
397 | if( 0.5 == spin ) { |
---|
398 | G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass); |
---|
399 | g *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2)); |
---|
400 | } |
---|
401 | if(g > 1.0) { |
---|
402 | G4cout << "### G4BetheBlochModel WARNING: g= " << g |
---|
403 | << dp->GetDefinition()->GetParticleName() |
---|
404 | << " Ekin(MeV)= " << kineticEnergy |
---|
405 | << " delEkin(MeV)= " << deltaKinEnergy |
---|
406 | << G4endl; |
---|
407 | } |
---|
408 | if(G4UniformRand() > g) return; |
---|
409 | } |
---|
410 | |
---|
411 | // delta-electron is produced |
---|
412 | G4double totMomentum = totEnergy*sqrt(beta2); |
---|
413 | G4double deltaMomentum = |
---|
414 | sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2)); |
---|
415 | G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) / |
---|
416 | (deltaMomentum * totMomentum); |
---|
417 | if(cost > 1.0) { |
---|
418 | G4cout << "### G4BetheBlochModel WARNING: cost= " |
---|
419 | << cost << " > 1 for " |
---|
420 | << dp->GetDefinition()->GetParticleName() |
---|
421 | << " Ekin(MeV)= " << kineticEnergy |
---|
422 | << " p(MeV/c)= " << totMomentum |
---|
423 | << " delEkin(MeV)= " << deltaKinEnergy |
---|
424 | << " delMom(MeV/c)= " << deltaMomentum |
---|
425 | << " tmin(MeV)= " << minKinEnergy |
---|
426 | << " tmax(MeV)= " << maxKinEnergy |
---|
427 | << " dir= " << dp->GetMomentumDirection() |
---|
428 | << G4endl; |
---|
429 | cost = 1.0; |
---|
430 | } |
---|
431 | G4double sint = sqrt((1.0 - cost)*(1.0 + cost)); |
---|
432 | |
---|
433 | G4double phi = twopi * G4UniformRand() ; |
---|
434 | |
---|
435 | |
---|
436 | G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost); |
---|
437 | G4ThreeVector direction = dp->GetMomentumDirection(); |
---|
438 | deltaDirection.rotateUz(direction); |
---|
439 | |
---|
440 | // create G4DynamicParticle object for delta ray |
---|
441 | G4DynamicParticle* delta = new G4DynamicParticle(theElectron, |
---|
442 | deltaDirection,deltaKinEnergy); |
---|
443 | |
---|
444 | vdp->push_back(delta); |
---|
445 | |
---|
446 | // Change kinematics of primary particle |
---|
447 | kineticEnergy -= deltaKinEnergy; |
---|
448 | G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum; |
---|
449 | finalP = finalP.unit(); |
---|
450 | |
---|
451 | fParticleChange->SetProposedKineticEnergy(kineticEnergy); |
---|
452 | fParticleChange->SetProposedMomentumDirection(finalP); |
---|
453 | } |
---|
454 | |
---|
455 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|