Skip to content
Snippets Groups Projects
Commit e5b8baf1 authored by MANZANILLAS Luis's avatar MANZANILLAS Luis
Browse files

add second angle in stokes vector for circular polarization

parent e2b4a488
No related branches found
No related tags found
No related merge requests found
......@@ -36,6 +36,7 @@ public:
void SetSourceDiameter(G4double newValue);
void SetSourcePolarizationDegree(G4double newValue);
void SetSourcePolarizationAngle(G4double newValue);
void SetSourceCircularPolarizationAngle(G4double newValue);
G4int GetSourceType(void){return fSourceType;};
G4int GetSourceGeometry(void){return fSourceGeometry;};
G4double GetSourceEnergy(void){return fSourceEnergy;};
......@@ -45,6 +46,7 @@ public:
G4double GetSourcePositionZ(void){return fPositionZ;};
G4double GetSourceDiameter(void){return fDiameter;};
G4double GetSourcePolarizationAngle(void){return fTheta_polar;};
G4double GetSourceCircularPolarizationAngle(void){return fOmega_polar;};
G4double GetSourcePolarizationDegree(void){return fPolarization_degree;};
private:
G4double fPositionX;
......@@ -59,6 +61,7 @@ private:
G4int fSourceGeometry;
G4double fSourceEnergy;
G4double fTheta_polar;
G4double fOmega_polar;
G4double fPolarization_degree;
DetectorConstruction* fDetector;
G4double fPhotonWavelength;
......
......@@ -28,6 +28,7 @@ class PrimaryGeneratorMessenger: public G4UImessenger
G4UIcmdWithAnInteger* fSourceGeometry;
G4UIcmdWithADoubleAndUnit* fSourceEnergy;
G4UIcmdWithADouble* fSourcePolarizationAngle;
G4UIcmdWithADouble* fSourceCircularPolarizationAngle;
G4UIcmdWithADouble* fSourcePolarizationDegree;
G4UIcmdWithADoubleAndUnit* fSourcePositionX;
G4UIcmdWithADoubleAndUnit* fSourcePositionY;
......
#!/bin/sh
x=0.
y=0.
for angle in 0 45 90
do
echo $angle
for job in {1..4..1}
do
for polarizationDegree in 0.5 0.6 0.7 0.8 0.9 1.0
do
FILE="run_polarized_at_${polarizationDegree}_angle_${angle}mm_job_${job}.mac"
/bin/cat <<EOM >$FILE
/process/em/fluo true
/process/em/auger true
/process/em/augerCascade true
/process/em/pixe true
/PoSOS/det/setOutputDirectory /nfs/tegile/work/experiences/detecteurs/manzanillas/Polarimetry/Geant4Output/
/PoSOS/det/setSetupName angle_${angle}_at_${polarizationDegree}0_pc_run_${job}
#select output format, options are: csv root hdf5
/PoSOS/det/setDataType csv
#Set the position of the collimator, a values between -37 and 37 mm (size of samples)
#Choise of source type: 0 gamma,1 Fe-55,2 Cs-137; 1 Bi-207; 2 Sr-90; 3 e-
/run/initialize
/PoSOS/gun/sourceType 0
/PoSOS/gun/sourcePositionX ${x} mm
/PoSOS/gun/sourcePositionY ${y} mm
/PoSOS/gun/sourceDiameter 0.1 mm
/PoSOS/gun/sourceEnergy 7.0 keV
#the default is linear polarization
#/PoSOS/gun/sourceGammaCircularPolarizationAngle ${angle}.
/PoSOS/gun/sourceGammaPolarizationAngle ${angle}.
/PoSOS/gun/sourceGammaPolarizationDegree ${polarizationDegree}
/run/beamOn 80000000
EOM
FILEJOB="job_polarization_${polarizationDegree}_${angle}degrees_job_${job}.sh"
/bin/cat <<EOM >$FILEJOB
#!/bin/bash
#SBATCH -n 1
#SBATCH --qos=parallel
#SBATCH --partition=sumo
#SBATCH --time=2:00:00
#SBATCH --cpus-per-task=2
#
source ~/environment_GeSOS.sh
./PoSOS -m $FILE
exit 0
EOM
sbatch $FILEJOB
done
done
done
......@@ -36,6 +36,7 @@ PrimaryGeneratorAction::PrimaryGeneratorAction(DetectorConstruction* det)
fSourceGeometry = 0;
fSourceEnergy = 10*keV;
fTheta_polar = 0; //in degrees, so the default is horizontal polarization
fOmega_polar = 90; //in degrees, so the default is linear polarization
fPolarization_degree = 0.95; // values from 0 to 1, so the default is 95% polarization
fPhotonWavelength = 0;
fParticleName = "void";
......@@ -78,9 +79,11 @@ void PrimaryGeneratorAction::GeneratePrimaries(G4Event* anEvent)
//Polarazation of gammas
G4double theta_polar_radians = fTheta_polar * pi * rad / 180;
G4double stoke_vector_e1 = cos(theta_polar_radians);
G4double stoke_vector_e2 = sin(theta_polar_radians);
G4ThreeVector stokes_vector = G4ThreeVector(stoke_vector_e1,stoke_vector_e2,0);
G4double omega_polar_radians = fOmega_polar * pi * rad / 180;
G4double stoke_vector_e1 = cos(theta_polar_radians)*cos(omega_polar_radians);
G4double stoke_vector_e2 = sin(theta_polar_radians)*cos(omega_polar_radians);
G4double stoke_vector_e3 = sin(omega_polar_radians);
G4ThreeVector stokes_vector = G4ThreeVector(stoke_vector_e1,stoke_vector_e2,stoke_vector_e3);
if( fPolarization_degree < G4UniformRand()){
......@@ -201,6 +204,10 @@ void PrimaryGeneratorAction::SetSourcePolarizationAngle(G4double newAngle){
fTheta_polar = newAngle;
}
void PrimaryGeneratorAction::SetSourceCircularPolarizationAngle(G4double newAngle){
fOmega_polar = newAngle;
}
void PrimaryGeneratorAction::SetSourcePolarizationDegree(G4double newPolarDegree){
fPolarization_degree = newPolarDegree;
}
......
......@@ -46,6 +46,12 @@ PrimaryGeneratorMessenger::PrimaryGeneratorMessenger(PrimaryGeneratorAction* Gun
fSourcePolarizationAngle->SetDefaultValue(0.);
fSourcePolarizationAngle->AvailableForStates(G4State_Idle);
fSourceCircularPolarizationAngle = new G4UIcmdWithADouble("/PoSOS/gun/sourceGammaCircularPolarizationAngle", this);
fSourceCircularPolarizationAngle->SetGuidance("Choose circular angle of polarization in degrees");
fSourceCircularPolarizationAngle->SetParameterName("sourceCircularPolarAngle",true);
fSourceCircularPolarizationAngle->SetDefaultValue(0.);
fSourceCircularPolarizationAngle->AvailableForStates(G4State_Idle);
fSourcePolarizationDegree = new G4UIcmdWithADouble("/PoSOS/gun/sourceGammaPolarizationDegree", this);
fSourcePolarizationDegree->SetGuidance("Choose degree of polarization o to 1");
fSourcePolarizationDegree->SetParameterName("sourcePolarDegree",true);
......@@ -87,6 +93,7 @@ PrimaryGeneratorMessenger::~PrimaryGeneratorMessenger()
delete fSourceGeometry;
delete fSourceEnergy;
delete fSourcePolarizationAngle;
delete fSourceCircularPolarizationAngle;
delete fSourcePolarizationDegree;
delete fGunDir;
delete fSourcePositionX;
......@@ -115,6 +122,9 @@ void PrimaryGeneratorMessenger::SetNewValue(G4UIcommand* command, G4String newVa
if(command == fSourcePolarizationAngle){
fAction->SetSourcePolarizationAngle(fSourcePolarizationAngle->GetNewDoubleValue(newValue));
}
if(command == fSourceCircularPolarizationAngle){
fAction->SetSourceCircularPolarizationAngle(fSourceCircularPolarizationAngle->GetNewDoubleValue(newValue));
}
if(command == fSourcePolarizationDegree){
fAction->SetSourcePolarizationDegree(fSourcePolarizationDegree->GetNewDoubleValue(newValue));
}
......
......@@ -5,7 +5,7 @@
/process/em/pixe true
#/run/numberOfThreads 4
#/PoSOS/det/setOutputDirectory /nfs/tegile/work/experiences/detecteurs/manzanillas/LEAPS_INNOV/Geant4_output/
/PoSOS/det/setOutputDirectory /nfs/tegile/work/experiences/detecteurs/manzanillas/Polarimetry/Geant4Output/
/PoSOS/det/setSetupName run_1
#select output format, options are: csv root hdf5
/PoSOS/det/setDataType csv
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment