source: trunk/source/processes/electromagnetic/utils/include/G4VEmProcess.hh@ 991

Last change on this file since 991 was 991, checked in by garnier, 17 years ago

update

File size: 21.1 KB
Line 
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: G4VEmProcess.hh,v 1.47 2008/07/31 13:01:26 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name: G4VEmProcess
35//
36// Author: Vladimir Ivanchenko
37//
38// Creation date: 01.10.2003
39//
40// Modifications:
41// 30-06-04 make destructor virtual (V.Ivanchenko)
42// 09-08-04 optimise integral option (V.Ivanchenko)
43// 11-08-04 add protected methods to access cuts (V.Ivanchenko)
44// 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko)
45// 16-09-04 Add flag for LambdaTable and method RecalculateLambda (VI)
46// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
47// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
48// 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
49// 09-05-05 Fix problem in logic when path boundary between materials (VI)
50// 11-01-06 add A to parameters of ComputeCrossSectionPerAtom (VI)
51// 01-02-06 put default value A=0. to keep compatibility with v5.2 (mma)
52// 13-05-06 Add method to access model by index (V.Ivanchenko)
53// 12-09-06 add SetModel() (mma)
54// 25-09-07 More accurate handling zero xsect in
55// PostStepGetPhysicalInteractionLength (V.Ivanchenko)
56// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
57// 15-07-08 Reorder class members for further multi-thread development (VI)
58//
59// Class Description:
60//
61// It is the unified Discrete process
62
63// -------------------------------------------------------------------
64//
65
66#ifndef G4VEmProcess_h
67#define G4VEmProcess_h 1
68
69#include "G4VDiscreteProcess.hh"
70#include "globals.hh"
71#include "G4Material.hh"
72#include "G4MaterialCutsCouple.hh"
73#include "G4Track.hh"
74#include "G4EmModelManager.hh"
75#include "G4UnitsTable.hh"
76#include "G4ParticleDefinition.hh"
77#include "G4ParticleChangeForGamma.hh"
78
79class G4Step;
80class G4VEmModel;
81class G4DataVector;
82class G4VParticleChange;
83class G4PhysicsTable;
84class G4PhysicsVector;
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87
88class G4VEmProcess : public G4VDiscreteProcess
89{
90public:
91
92 G4VEmProcess(const G4String& name,
93 G4ProcessType type = fElectromagnetic);
94
95 virtual ~G4VEmProcess();
96
97 //------------------------------------------------------------------------
98 // Virtual methods to be implemented in concrete processes
99 //------------------------------------------------------------------------
100
101 virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
102
103 virtual void PrintInfo() = 0;
104
105protected:
106
107 virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
108
109 //------------------------------------------------------------------------
110 // Methods with standard implementation; may be overwritten if needed
111 //------------------------------------------------------------------------
112
113 inline G4double RecalculateLambda(G4double kinEnergy,
114 const G4MaterialCutsCouple* couple);
115
116 //------------------------------------------------------------------------
117 // Generic methods common to all Discrete processes
118 //------------------------------------------------------------------------
119
120public:
121
122 void PrintInfoDefinition();
123
124 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
125
126 // Initialise for build of tables
127 void PreparePhysicsTable(const G4ParticleDefinition&);
128
129 // Build physics table during initialisation
130 void BuildPhysicsTable(const G4ParticleDefinition&);
131
132 // Store PhysicsTable in a file.
133 // Return false in case of failure at I/O
134 G4bool StorePhysicsTable(const G4ParticleDefinition*,
135 const G4String& directory,
136 G4bool ascii = false);
137
138 // Retrieve Physics from a file.
139 // (return true if the Physics Table can be build by using file)
140 // (return false if the process has no functionality or in case of failure)
141 // File name should is constructed as processName+particleName and the
142 // should be placed under the directory specifed by the argument.
143 G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
144 const G4String& directory,
145 G4bool ascii);
146
147 //------------------------------------------------------------------------
148 // Specific methods for Discrete EM post step simulation
149 //------------------------------------------------------------------------
150
151 // It returns the cross section per volume for energy/ material
152 G4double CrossSectionPerVolume(G4double kineticEnergy,
153 const G4MaterialCutsCouple* couple);
154
155 // implementation of virtual method
156 virtual G4double PostStepGetPhysicalInteractionLength(
157 const G4Track& track,
158 G4double previousStepSize,
159 G4ForceCondition* condition
160 );
161
162 // It returns the cross section of the process per atom
163 inline G4double ComputeCrossSectionPerAtom(G4double kineticEnergy,
164 G4double Z, G4double A=0.,
165 G4double cut=0.0);
166
167 inline G4double MeanFreePath(const G4Track& track);
168
169 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
170 size_t& idxRegion) const;
171
172 // It returns cross section per volume
173 inline G4double GetLambda(G4double& kinEnergy,
174 const G4MaterialCutsCouple* couple);
175
176 //------------------------------------------------------------------------
177 // Specific methods to build and access Physics Tables
178 //------------------------------------------------------------------------
179
180 // Binning for lambda table
181 inline void SetLambdaBinning(G4int nbins);
182 inline G4int LambdaBinning() const;
183
184 // Min kinetic energy for tables
185 inline void SetMinKinEnergy(G4double e);
186 inline G4double MinKinEnergy() const;
187
188 // Max kinetic energy for tables
189 inline void SetMaxKinEnergy(G4double e);
190 inline G4double MaxKinEnergy() const;
191
192 inline void SetPolarAngleLimit(G4double a);
193 inline G4double PolarAngleLimit() const;
194
195 inline const G4PhysicsTable* LambdaTable() const;
196
197 //------------------------------------------------------------------------
198 // Define and access particle type
199 //------------------------------------------------------------------------
200
201 inline const G4ParticleDefinition* Particle() const;
202 inline const G4ParticleDefinition* SecondaryParticle() const;
203
204 //------------------------------------------------------------------------
205 // Specific methods to set, access, modify models
206 //------------------------------------------------------------------------
207
208 // Add EM model coupled for the region
209 inline void AddEmModel(G4int, G4VEmModel*, const G4Region* region = 0);
210
211 // Assign a model to a process
212 inline void SetModel(G4VEmModel*);
213
214 // return the assigned model
215 inline G4VEmModel* Model();
216
217 // Define new energy range for the model identified by the name
218 inline void UpdateEmModel(const G4String&, G4double, G4double);
219
220 // Access to models
221 inline G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
222
223 //------------------------------------------------------------------------
224 // Get/set parameters used for simulation of energy loss
225 //------------------------------------------------------------------------
226
227 inline void ActivateDeexcitation(G4bool, const G4Region* r = 0);
228
229 inline void SetLambdaFactor(G4double val);
230
231 inline void SetIntegral(G4bool val);
232 inline G4bool IsIntegral() const;
233
234 inline void SetApplyCuts(G4bool val);
235
236protected:
237
238 G4double GetMeanFreePath(const G4Track& track,
239 G4double previousStepSize,
240 G4ForceCondition* condition);
241
242 G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*);
243
244 inline G4ParticleChangeForGamma* GetParticleChange();
245
246 inline void SetParticle(const G4ParticleDefinition* p);
247
248 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
249
250 inline G4VEmModel* SelectModel(G4double& kinEnergy);
251
252 inline size_t CurrentMaterialCutsCoupleIndex() const;
253
254 inline G4double GetGammaEnergyCut();
255
256 inline G4double GetElectronEnergyCut();
257
258 inline void SetBuildTableFlag(G4bool val);
259
260 inline void SetStartFromNullFlag(G4bool val);
261
262private:
263
264 void Clear();
265
266 void BuildLambdaTable();
267
268 void FindLambdaMax();
269
270 inline void InitialiseStep(const G4Track&);
271
272 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
273
274 inline void ComputeIntegralLambda(G4double kinEnergy);
275
276 inline G4double GetLambdaFromTable(G4double kinEnergy);
277
278 inline G4double GetCurrentLambda(G4double kinEnergy);
279
280 inline G4double ComputeCurrentLambda(G4double kinEnergy);
281
282 // hide assignment operator
283
284 G4VEmProcess(G4VEmProcess &);
285 G4VEmProcess & operator=(const G4VEmProcess &right);
286
287 // ======== Parameters of the class fixed at construction =========
288
289 G4EmModelManager* modelManager;
290 const G4ParticleDefinition* theGamma;
291 const G4ParticleDefinition* theElectron;
292 const G4ParticleDefinition* thePositron;
293 const G4ParticleDefinition* secondaryParticle;
294
295 G4bool buildLambdaTable;
296
297 // ======== Parameters of the class fixed at initialisation =======
298
299 // tables and vectors
300 G4PhysicsTable* theLambdaTable;
301 G4double* theEnergyOfCrossSectionMax;
302 G4double* theCrossSectionMax;
303
304 const std::vector<G4double>* theCuts;
305 const std::vector<G4double>* theCutsGamma;
306 const std::vector<G4double>* theCutsElectron;
307 const std::vector<G4double>* theCutsPositron;
308
309 G4int nLambdaBins;
310
311 G4double minKinEnergy;
312 G4double maxKinEnergy;
313 G4double lambdaFactor;
314 G4double polarAngleLimit;
315
316 G4bool integral;
317 G4bool applyCuts;
318 G4bool startFromNull;
319
320 G4int nRegions;
321 std::vector<G4Region*> regions;
322 std::vector<G4bool> flagsDeexcitation;
323
324 // ======== Cashed values - may be state dependent ================
325
326protected:
327
328 G4ParticleChangeForGamma fParticleChange;
329
330private:
331
332 std::vector<G4DynamicParticle*> secParticles;
333
334 G4VEmModel* selectedModel;
335
336 const G4ParticleDefinition* particle;
337
338 // cash
339 const G4Material* currentMaterial;
340 const G4MaterialCutsCouple* currentCouple;
341 size_t currentMaterialIndex;
342
343 G4double mfpKinEnergy;
344 G4double preStepKinEnergy;
345 G4double preStepLambda;
346
347};
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
350//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
351
352inline void G4VEmProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
353{
354 if(couple != currentCouple) {
355 currentCouple = couple;
356 currentMaterial = couple->GetMaterial();
357 currentMaterialIndex = couple->GetIndex();
358 mfpKinEnergy = DBL_MAX;
359 }
360}
361
362//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
363
364inline void G4VEmProcess::InitialiseStep(const G4Track& track)
365{
366 preStepKinEnergy = track.GetKineticEnergy();
367 DefineMaterial(track.GetMaterialCutsCouple());
368 if (theNumberOfInteractionLengthLeft < 0.0) mfpKinEnergy = DBL_MAX;
369}
370
371//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
372
373inline G4double G4VEmProcess::GetLambda(G4double& kineticEnergy,
374 const G4MaterialCutsCouple* couple)
375{
376 DefineMaterial(couple);
377 return GetCurrentLambda(kineticEnergy);
378}
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
381
382inline G4double G4VEmProcess::GetCurrentLambda(G4double e)
383{
384 G4double x = 0.0;
385 if(theLambdaTable) x = GetLambdaFromTable(e);
386 else x = ComputeCurrentLambda(e);
387 return x;
388}
389
390//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
391
392inline G4double G4VEmProcess::RecalculateLambda(G4double e,
393 const G4MaterialCutsCouple* couple)
394{
395 DefineMaterial(couple);
396 return ComputeCurrentLambda(e);
397}
398
399//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
400
401inline G4double G4VEmProcess::ComputeCurrentLambda(G4double e)
402{
403 G4VEmModel* currentModel = SelectModel(e);
404 G4double x = 0.0;
405 if(currentModel)
406 x = currentModel->CrossSectionPerVolume(currentMaterial,particle,
407 e,(*theCuts)[currentMaterialIndex]);
408 return x;
409}
410
411//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
412
413inline G4double G4VEmProcess::GetLambdaFromTable(G4double e)
414{
415 G4bool b;
416 return (((*theLambdaTable)[currentMaterialIndex])->GetValue(e, b));
417}
418
419//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
420
421inline void G4VEmProcess::ComputeIntegralLambda(G4double e)
422{
423 mfpKinEnergy = theEnergyOfCrossSectionMax[currentMaterialIndex];
424 if (e <= mfpKinEnergy) {
425 preStepLambda = GetLambdaFromTable(e);
426
427 } else {
428 G4double e1 = e*lambdaFactor;
429 if(e1 > mfpKinEnergy) {
430 preStepLambda = GetLambdaFromTable(e);
431 G4double preStepLambda1 = GetLambdaFromTable(e1);
432 if(preStepLambda1 > preStepLambda) {
433 mfpKinEnergy = e1;
434 preStepLambda = preStepLambda1;
435 }
436 } else {
437 preStepLambda = theCrossSectionMax[currentMaterialIndex];
438 }
439 }
440}
441
442//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
443
444inline G4double G4VEmProcess::MeanFreePath(const G4Track& track)
445{
446 DefineMaterial(track.GetMaterialCutsCouple());
447 preStepLambda = GetCurrentLambda(track.GetKineticEnergy());
448 G4double x = DBL_MAX;
449 if(DBL_MIN < preStepLambda) x = 1.0/preStepLambda;
450 return x;
451}
452
453//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
454
455inline G4VEmModel* G4VEmProcess::SelectModel(G4double& kinEnergy)
456{
457 return modelManager->SelectModel(kinEnergy, currentMaterialIndex);
458}
459
460//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
461
462inline G4VEmModel* G4VEmProcess::SelectModelForMaterial(
463 G4double kinEnergy, size_t& idxRegion) const
464{
465 return modelManager->SelectModel(kinEnergy, idxRegion);
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
469
470inline const G4ParticleDefinition* G4VEmProcess::Particle() const
471{
472 return particle;
473}
474
475//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
476
477inline const G4ParticleDefinition* G4VEmProcess::SecondaryParticle() const
478{
479 return secondaryParticle;
480}
481
482//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
483
484inline G4double G4VEmProcess::GetGammaEnergyCut()
485{
486 return (*theCutsGamma)[currentMaterialIndex];
487}
488
489//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
490
491inline G4double G4VEmProcess::GetElectronEnergyCut()
492{
493 return (*theCutsElectron)[currentMaterialIndex];
494}
495
496//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
497
498inline void G4VEmProcess::SetLambdaFactor(G4double val)
499{
500 if(val > 0.0 && val <= 1.0) lambdaFactor = val;
501}
502
503//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
504
505inline G4VEmModel* G4VEmProcess::GetModelByIndex(G4int idx, G4bool ver)
506{
507 return modelManager->GetModel(idx, ver);
508}
509
510//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
511
512inline G4ParticleChangeForGamma* G4VEmProcess::GetParticleChange()
513{
514 return &fParticleChange;
515}
516
517//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
518
519inline void G4VEmProcess::SetParticle(const G4ParticleDefinition* p)
520{
521 particle = p;
522}
523
524//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
525
526inline void G4VEmProcess::SetSecondaryParticle(const G4ParticleDefinition* p)
527{
528 secondaryParticle = p;
529}
530
531//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
532
533inline void G4VEmProcess::AddEmModel(G4int order, G4VEmModel* p,
534 const G4Region* region)
535{
536 G4VEmFluctuationModel* fm = 0;
537 modelManager->AddEmModel(order, p, fm, region);
538 if(p) p->SetParticleChange(pParticleChange);
539}
540
541//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
542
543inline void G4VEmProcess::SetModel(G4VEmModel* model)
544{
545 selectedModel = model;
546}
547
548//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
549
550inline G4VEmModel* G4VEmProcess::Model()
551{
552 return selectedModel;
553}
554
555//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
556
557inline void G4VEmProcess::UpdateEmModel(const G4String& nam,
558 G4double emin, G4double emax)
559{
560 modelManager->UpdateEmModel(nam, emin, emax);
561}
562
563//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
564
565inline G4double G4VEmProcess::ComputeCrossSectionPerAtom(
566 G4double kineticEnergy, G4double Z, G4double A, G4double cut)
567{
568 G4VEmModel* model = SelectModel(kineticEnergy);
569 return model->ComputeCrossSectionPerAtom(particle,kineticEnergy,Z,A,cut);
570}
571
572//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
573
574inline void G4VEmProcess::SetLambdaBinning(G4int nbins)
575{
576 nLambdaBins = nbins;
577}
578
579//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
580
581inline G4int G4VEmProcess::LambdaBinning() const
582{
583 return nLambdaBins;
584}
585
586//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
587
588inline void G4VEmProcess::SetMinKinEnergy(G4double e)
589{
590 minKinEnergy = e;
591}
592
593//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
594
595inline G4double G4VEmProcess::MinKinEnergy() const
596{
597 return minKinEnergy;
598}
599
600//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
601
602inline void G4VEmProcess::SetMaxKinEnergy(G4double e)
603{
604 maxKinEnergy = e;
605}
606
607//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
608
609inline G4double G4VEmProcess::MaxKinEnergy() const
610{
611 return maxKinEnergy;
612}
613
614//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
615
616inline void G4VEmProcess::SetPolarAngleLimit(G4double val)
617{
618 if(val < 0.0) polarAngleLimit = 0.0;
619 else if(val > pi) polarAngleLimit = pi;
620 else polarAngleLimit = val;
621}
622
623//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
624
625inline G4double G4VEmProcess::PolarAngleLimit() const
626{
627 return polarAngleLimit;
628}
629
630//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
631
632inline void G4VEmProcess::ActivateDeexcitation(G4bool, const G4Region*)
633{}
634
635//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
636
637inline const G4PhysicsTable* G4VEmProcess::LambdaTable() const
638{
639 return theLambdaTable;
640}
641
642//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
643
644inline void G4VEmProcess::SetIntegral(G4bool val)
645{
646 if(particle && particle != theGamma) integral = val;
647 if(integral) buildLambdaTable = true;
648}
649
650//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
651
652inline G4bool G4VEmProcess::IsIntegral() const
653{
654 return integral;
655}
656
657//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
658
659inline void G4VEmProcess::SetBuildTableFlag(G4bool val)
660{
661 buildLambdaTable = val;
662 if(!val) integral = false;
663}
664
665//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
666
667inline void G4VEmProcess::SetStartFromNullFlag(G4bool val)
668{
669 startFromNull = val;
670}
671
672//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
673
674inline void G4VEmProcess::SetApplyCuts(G4bool val)
675{
676 applyCuts = val;
677}
678
679//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
680
681inline size_t G4VEmProcess::CurrentMaterialCutsCoupleIndex() const
682{
683 return currentMaterialIndex;
684}
685
686//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
687
688#endif
Note: See TracBrowser for help on using the repository browser.