Laser Profiles

Gaussian Beam

template<typename T_Params>
struct GaussianBeam : public picongpu::fields::laserProfiles::gaussianBeam::Unitless<T_Params>

Gaussian Beam laser profile with finite pulse length.

Template Parameters
  • T_Params: class parameter to configure the Gaussian Beam profile, see members of gaussianBeam::default::GaussianBeamParam for required members

    //! Use only the 0th Laguerremode for a standard Gaussian
    static constexpr uint32_t MODENUMBER = 0;
    PMACC_CONST_VECTOR(float_X, MODENUMBER + 1, LAGUERREMODES, 1.0);
    // This is just an example for a more complicated set of Laguerre modes
    //constexpr uint32_t MODENUMBER = 12;
    //PMACC_CONST_VECTOR(float_X, MODENUMBER + 1, LAGUERREMODES, -1.0, 0.0300519, 0.319461, -0.23783, 0.0954839, 0.0318653, -0.144547, 0.0249208, -0.111989, 0.0434385, -0.030038, -0.00896321, -0.0160788);

    struct GaussianBeamParam
    {
        /** unit: meter */
        static constexpr float_64 WAVE_LENGTH_SI = 0.8e-6;

        /** Convert the normalized laser strength parameter a0 to Volt per meter */
        static constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * ::picongpu::SI::ELECTRON_MASS_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI / ::picongpu::SI::ELECTRON_CHARGE_SI;

        /** unit: W / m^2 */
        // calculate: _A0 = 8.549297e-6 * sqrt( Intensity[W/m^2] ) * wavelength[m] (linearly polarized)

        /** unit: none */
        //static constexpr float_64 _A0  = 1.5;

        /** unit: Volt / meter */
        //static constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI;

        /** unit: Volt / meter */
        static constexpr float_64 AMPLITUDE_SI = 1.738e13;

        /** Pulse length: sigma of std. gauss for intensity (E^2)
         *  PULSE_LENGTH_SI = FWHM_of_Intensity   / [ 2*sqrt{ 2* ln(2) } ]
         *                                          [    2.354820045     ]
         *  Info:             FWHM_of_Intensity = FWHM_Illumination
         *                      = what a experimentalist calls "pulse duration"
         *
         *  unit: seconds (1 sigma) */
        static constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0;

        /** beam waist: distance from the axis where the pulse intensity (E^2)
         *              decreases to its 1/e^2-th part,
         *              at the focus position of the laser
         * W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) }
         *                             [   1.17741    ]
         *
         *  unit: meter */
        static constexpr float_64 W0_SI = 5.0e-6 / 1.17741;
        /** the distance to the laser focus in y-direction
         *  unit: meter */
        static constexpr float_64 FOCUS_POS_SI = 4.62e-5;

        /** The laser pulse will be initialized PULSE_INIT times of the PULSE_LENGTH
         *
         *  unit: none */
        static constexpr float_64 PULSE_INIT = 20.0;

        /** cell from top where the laser is initialized
         *
         * if `initPlaneY == 0` than the absorber are disabled.
         * if `initPlaneY > absorbercells negative Y` the negative absorber in y
         * direction is enabled
         *
         * valid ranges:
         *   - initPlaneY == 0
         *   - absorber cells negative Y < initPlaneY < cells in y direction of the top gpu
         */
        static constexpr uint32_t initPlaneY = 0;

        /** laser phase shift (no shift: 0.0)
         *
         * sin(omega*time + laser_phase): starts with phase=0 at center --> E-field=0 at center
         *
         * unit: rad, periodic in 2*pi
         */
        static constexpr float_X LASER_PHASE = 0.0;

        using LAGUERREMODES_t = defaults::LAGUERREMODES_t;
        static constexpr uint32_t MODENUMBER = defaults::MODENUMBER;

        /** Available polarisation types
         */
        enum PolarisationType
        {
            LINEAR_X = 1u,
            LINEAR_Z = 2u,
            CIRCULAR = 4u,
        };
        /** Polarization selection
         */
        static constexpr PolarisationType Polarisation = CIRCULAR;
    };

Gaussian Beam with Pulse Front Tilt

template<typename T_Params>
struct PulseFrontTilt : public picongpu::fields::laserProfiles::pulseFrontTilt::Unitless<T_Params>

Gaussian Beam laser profile with titled pulse front.

Template Parameters
  • T_Params: class parameter to configure the Gaussian Beam with pulse front titlt, see members of pulseFrontTilt::defaults::PulseFrontTiltParam for required members

    struct PulseFrontTiltParam
    {
        /** unit: meter */
        static constexpr float_64 WAVE_LENGTH_SI = 0.8e-6;

        /** Convert the normalized laser strength parameter a0 to Volt per meter */
        static constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * ::picongpu::SI::ELECTRON_MASS_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI / ::picongpu::SI::ELECTRON_CHARGE_SI;

        /** unit: W / m^2 */
        // calculate: _A0 = 8.549297e-6 * sqrt( Intensity[W/m^2] ) * wavelength[m] (linearly polarized)

        /** unit: none */
        //static constexpr float_64 _A0  = 1.5;

        /** unit: Volt / meter */
        //static constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI;

        /** unit: Volt / meter */
        static constexpr float_64 AMPLITUDE_SI = 1.738e13;

        /** Pulse length: sigma of std. gauss for intensity (E^2)
         *  PULSE_LENGTH_SI = FWHM_of_Intensity   / [ 2*sqrt{ 2* ln(2) } ]
         *                                          [    2.354820045     ]
         *  Info:             FWHM_of_Intensity = FWHM_Illumination
         *                      = what a experimentalist calls "pulse duration"
         *
         *  unit: seconds (1 sigma) */
        static constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0;

        /** beam waist: distance from the axis where the pulse intensity (E^2)
         *              decreases to its 1/e^2-th part,
         *              at the focus position of the laser
         * W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) }
         *                             [   1.17741    ]
         *
         *  unit: meter */
        static constexpr float_64 W0_SI = 5.0e-6 / 1.17741;

        /** the distance to the laser focus in y-direction
         *  unit: meter */
        static constexpr float_64 FOCUS_POS_SI = 4.62e-5;

        /** the tilt angle between laser propagation in y-direction and laser axis in
        *  x-direction (0 degree == no tilt)
        *  unit: degree */
        static constexpr float_64 TILT_X_SI = 0.0;

        /** The laser pulse will be initialized PULSE_INIT times of the PULSE_LENGTH
         *
         *  unit: none */
        static constexpr float_64 PULSE_INIT = 20.0;

        /** cell from top where the laser is initialized
         *
         * if `initPlaneY == 0` than the absorber are disabled.
         * if `initPlaneY > absorbercells negative Y` the negative absorber in y
         * direction is enabled
         *
         * valid ranges:
         *   - initPlaneY == 0
         *   - absorber cells negative Y < initPlaneY < cells in y direction of the top gpu
         */
        static constexpr uint32_t initPlaneY = 0;

        /** laser phase shift (no shift: 0.0)
         *
         * sin(omega*time + laser_phase): starts with phase=0 at center --> E-field=0 at center
         *
         * unit: rad, periodic in 2*pi
         */
        static constexpr float_X LASER_PHASE = 0.0;

        //! Available polarisation types
        enum PolarisationType
        {
            LINEAR_X = 1u,
            LINEAR_Z = 2u,
            CIRCULAR = 4u,
        };

        /** Polarization selection
         */
        static constexpr PolarisationType Polarisation = LINEAR_X;
    };

Wavepacket

template<typename T_Params>
struct Wavepacket : public picongpu::fields::laserProfiles::wavepacket::Unitless<T_Params>

Wavepacket with Gaussian spatial and temporal envelope.

Template Parameters
  • T_Params: class parameter to configure the Wavepacket profile, see members of wavepacket::defaults::WavepacketParam for required members

    struct WavepacketParam
    {
        /** unit: meter */
        static constexpr float_64 WAVE_LENGTH_SI = 0.8e-6;

        /** Convert the normalized laser strength parameter a0 to Volt per meter */
        static constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * ::picongpu::SI::ELECTRON_MASS_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI / ::picongpu::SI::ELECTRON_CHARGE_SI;

        /** unit: W / m^2 */
        // calculate: _A0 = 8.549297e-6 * sqrt( Intensity[W/m^2] ) * wavelength[m] (linearly polarized)

        /** unit: none */
        //static constexpr float_64 _A0  = 1.5;

        /** unit: Volt / meter */
        //static constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI;

        /** unit: Volt / meter */
        static constexpr float_64 AMPLITUDE_SI = 1.738e13;

        /** Stretch temporal profile by a constant plateau between the up and downramp
         *  unit: seconds */
        static constexpr float_64 LASER_NOFOCUS_CONSTANT_SI = 7.0 * WAVE_LENGTH_SI / ::picongpu::SI::SPEED_OF_LIGHT_SI;

        /** Pulse length: sigma of std. gauss for intensity (E^2)
         *  PULSE_LENGTH_SI = FWHM_of_Intensity   / [ 2*sqrt{ 2* ln(2) } ]
         *                                          [    2.354820045     ]
         *  Info:             FWHM_of_Intensity = FWHM_Illumination
         *                      = what a experimentalist calls "pulse duration"
         *
         *  unit: seconds (1 sigma) */
        static constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0;

        /** beam waist: distance from the axis where the pulse intensity (E^2)
         *              decreases to its 1/e^2-th part,
         *              at the focus position of the laser
         * W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) }
         *                             [   1.17741    ]
         *
         *  unit: meter */
        static constexpr float_64 W0_X_SI = 4.246e-6;
        static constexpr float_64 W0_Z_SI = W0_X_SI;

        /** The laser pulse will be initialized PULSE_INIT times of the PULSE_LENGTH
         *
         *  unit: none */
        static constexpr float_64 PULSE_INIT = 20.0;

        /** cell from top where the laser is initialized
         *
         * if `initPlaneY == 0` than the absorber are disabled.
         * if `initPlaneY > absorbercells negative Y` the negative absorber in y
         * direction is enabled
         *
         * valid ranges:
         *   - initPlaneY == 0
         *   - absorber cells negative Y < initPlaneY < cells in y direction of the top gpu
         */
        static constexpr uint32_t initPlaneY = 0;

        /** laser phase shift (no shift: 0.0)
         *
         * sin(omega*time + laser_phase): starts with phase=0 at center --> E-field=0 at center
         *
         * unit: rad, periodic in 2*pi
         */
        static constexpr float_X LASER_PHASE = 0.0;

        /** Available polarisation types
         */
        enum PolarisationType
        {
            LINEAR_X = 1u,
            LINEAR_Z = 2u,
            CIRCULAR = 4u,
        };
        /** Polarization selection
         */
        static constexpr PolarisationType Polarisation = LINEAR_X;
    };
} // namespace defaults

Wavepacket with Exponential Ramp and Prepulse

template<typename T_Params>
struct ExpRampWithPrepulse : public picongpu::fields::laserProfiles::expRampWithPrepulse::Unitless<T_Params>

Wavepacket with spatial Gaussian envelope and adjustable temporal shape.

Allows defining a prepulse and two regions of exponential preramp with independent slopes. The definition works by specifying three (t, intensity)- points, where time is counted from the very beginning in SI and the intensity (yes, intensity, not amplitude) is given in multiples of the main peak.

Be careful - problematic for few cycle pulses. Thought the rest is cloned from laserWavepacket, the correctionFactor is not included (this made a correction to the laser phase, which is necessary for very short pulses, since otherwise a test particle is, after the laser pulse has passed, not returned to immobility, as it should). Since the analytical solution is only implemented for the Gaussian regime, and we have mostly exponential regimes here, it was not retained here.

A Gaussian peak (optionally lengthened by a plateau) is preceded by two pieces of exponential preramps, defined by 3 (time, intensity)- -points.

The first two points get connected by an exponential, the 2nd and 3rd point are connected by another exponential, which is then extrapolated to the peak. The Gaussian is added everywhere, but typically contributes significantly only near the peak. It is advisable to set the third point far enough from the plateau (approx 3*FWHM), then the contribution from the Gaussian is negligible there, and the intensity can be set as measured from the laser profile.

Optionally a Gaussian prepulse can be added, given by the parameters of the relative intensity and time point. The time of the prepulse and the three preramp points are given in SI, the intensities are given as multiples of the peak intensity.

Template Parameters
  • T_Params: class parameter to configure the Gaussian Beam profile, see members of expRampWithPrepulse::defaults::ExpRampWithPrepulseParam for required members

    struct ExpRampWithPrepulseParam
    {
        // Intensities of prepulse and exponential preramp
        static constexpr float_X INT_RATIO_PREPULSE = 0.;
        static constexpr float_X INT_RATIO_POINT_1 = 1.e-8;
        static constexpr float_X INT_RATIO_POINT_2 = 1.e-4;
        static constexpr float_X INT_RATIO_POINT_3 = 1.e-4;

        // time-positions of prepulse and preramps points
        static constexpr float_64 TIME_PREPULSE_SI = -950.0e-15;
        static constexpr float_64 TIME_PEAKPULSE_SI = 0.0e-15;
        static constexpr float_64 TIME_POINT_1_SI = -1000.0e-15;
        static constexpr float_64 TIME_POINT_2_SI = -300.0e-15;
        static constexpr float_64 TIME_POINT_3_SI = -100.0e-15;

        /** unit: meter */
        static constexpr float_64 WAVE_LENGTH_SI = 0.8e-6;

        /** UNITCONV */
        static constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * ::picongpu::SI::ELECTRON_MASS_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI / ::picongpu::SI::ELECTRON_CHARGE_SI;

        /** unit: W / m^2 */
        // calculate: _A0 = 8.549297e-6 * sqrt( Intensity[W/m^2] ) * wavelength[m] (linearly polarized)

        /** unit: none */
        static constexpr float_64 _A0  = 20.;

        /** unit: Volt /meter */
        static constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI;

        /** unit: Volt /meter */
        //constexpr float_64 AMPLITUDE_SI = 1.738e13;

        /** Stretch temporal profile by a constant plateau between the up and downramp
         *  unit: seconds */
        static constexpr float_64 LASER_NOFOCUS_CONSTANT_SI = 0.0 * WAVE_LENGTH_SI / ::picongpu::SI::SPEED_OF_LIGHT_SI;

        /** Pulse length: sigma of std. gauss for intensity (E^2)
         *  PULSE_LENGTH_SI = FWHM_of_Intensity   / [ 2*sqrt{ 2* ln(2) } ]
         *                                          [    2.354820045     ]
         *  Info:             FWHM_of_Intensity = FWHM_Illumination
         *                      = what a experimentalist calls "pulse duration"
         *  unit: seconds (1 sigma) */
        static constexpr float_64 PULSE_LENGTH_SI = 3.0e-14 / 2.35482; // half of the time in which E falls to half its initial value (then I falls to half its value in 15fs, approx 6 wavelengths). Those are 4.8 wavelenghts.

        /** beam waist: distance from the axis where the pulse intensity (E^2)
         *              decreases to its 1/e^2-th part,
         *              WO_X_SI is this distance in x-direction
         *              W0_Z_SI is this distance in z-direction
         *              if both values are equal, the laser has a circular shape in x-z
         * W0_SI = FWHM_of_Intensity / sqrt{ 2* ln(2) }
         *                             [   1.17741    ]
         *  unit: meter */
        static constexpr float_64 W0_X_SI = 2.5 * WAVE_LENGTH_SI;
        static constexpr float_64 W0_Z_SI = W0_X_SI;

        /** The laser pulse will be initialized half of PULSE_INIT times of the PULSE_LENGTH before plateau
         *  and half at the end of the plateau
         *  unit: none */
        static constexpr float_64 RAMP_INIT = 16.0;

        /** cell from top where the laser is initialized
         *
         * if `initPlaneY == 0` than the absorber are disabled.
         * if `initPlaneY > absorbercells negative Y` the negative absorber in y
         * direction is enabled
         *
         * valid ranges:
         *   - initPlaneY == 0
         *   - absorber cells negative Y < initPlaneY < cells in y direction of the top gpu
         */
        static constexpr uint32_t initPlaneY = 0;

        /** laser phase shift (no shift: 0.0)
         *
         * sin(omega*time + laser_phase): starts with phase=0 at center --> E-field=0 at center
         *
         * unit: rad, periodic in 2*pi
         */
        static constexpr float_X LASER_PHASE = 0.0;

        /** Available polarisation types
         */
        enum PolarisationType
        {
            LINEAR_X = 1u,
            LINEAR_Z = 2u,
            CIRCULAR = 4u,
        };

        /** Polarization selection
         */
        static constexpr PolarisationType Polarisation = LINEAR_X;
    };
} // namespace defaults

Wavepacket with Polynomial Profile

template<typename T_Params>
struct Polynom : public picongpu::fields::laserProfiles::polynom::Unitless<T_Params>

Wavepacket with a polynomial temporal intensity shape.

Based on a wavepacket with Gaussian spatial envelope.

Template Parameters
  • T_Params: class parameter to configure the polynomial laser profile, see members of polynom::defaults::PolynomParam for required members

    struct PolynomParam
    {
        /** unit: meter */
        static constexpr float_64 WAVE_LENGTH_SI = 0.8e-6;

        /** Convert the normalized laser strength parameter a0 to Volt per meter */
        static constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * ::picongpu::SI::ELECTRON_MASS_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI / ::picongpu::SI::ELECTRON_CHARGE_SI;

        /** unit: W / m^2 */
        // calculate: _A0 = 8.549297e-6 * sqrt( Intensity[W/m^2] ) * wavelength[m] (linearly polarized)

        /** unit: none */
        //static constexpr float_64 _A0  = 1.5;

        /** unit: Volt / meter */
        //static constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI;

        /** unit: Volt / meter */
        static constexpr float_64 AMPLITUDE_SI = 1.738e13;

        /** Pulse length: sigma of std. gauss for intensity (E^2)
         *  PULSE_LENGTH_SI = FWHM_of_Intensity   / [ 2*sqrt{ 2* ln(2) } ]
         *                                          [    2.354820045     ]
         *  Info:             FWHM_of_Intensity = FWHM_Illumination
         *                      = what a experimentalist calls "pulse duration"
         *  unit: seconds (1 sigma) */
        static constexpr float_64 PULSE_LENGTH_SI = 4.0e-15;

        /** beam waist: distance from the axis where the pulse intensity (E^2)
         *              decreases to its 1/e^2-th part,
         *              at the focus position of the laser
         *  unit: meter
         */
        static constexpr float_64 W0_X_SI = 4.246e-6; // waist in x-direction
        static constexpr float_64 W0_Z_SI = W0_X_SI; // waist in z-direction

        /** cell from top where the laser is initialized
         *
         * if `initPlaneY == 0` than the absorber are disabled.
         * if `initPlaneY > absorbercells negative Y` the negative absorber in y
         * direction is enabled
         *
         * valid ranges:
         *   - initPlaneY == 0
         *   - absorber cells negative Y < initPlaneY < cells in y direction of the top gpu
         */
        static constexpr uint32_t initPlaneY = 0;

        /** laser phase shift (no shift: 0.0)
         *
         * sin(omega*time + laser_phase): starts with phase=0 at center --> E-field=0 at center
         *
         * unit: rad, periodic in 2*pi
         */
        static constexpr float_X LASER_PHASE = 0.0;

        /** Available polarization types
         */
        enum PolarisationType
        {
            LINEAR_X = 1u,
            LINEAR_Z = 2u,
            CIRCULAR = 4u,
        };
        /** Polarization selection
         */
        static constexpr PolarisationType Polarisation = LINEAR_X;
    };
} // namespace defaults
} // namespace gaussianBeam

    /** Wavepacket with a polynomial temporal intensity shape.
     *
     * Based on a wavepacket with Gaussian spatial envelope.
     *
     * @tparam T_Params class parameter to configure the polynomial laser profile,
     *                  see members of polynom::defaults::PolynomParam for
     *                  required members

Plane Wave

template<typename T_Params>
struct PlaneWave : public picongpu::fields::laserProfiles::planeWave::Unitless<T_Params>

Plane wave laser profile.

Defines a plane wave with temporally Gaussian envelope.

Template Parameters
  • T_Params: class parameter to configure the plane wave profile, see members of planeWave::defaults::PlaneWaveParam for required members

    struct PlaneWaveParam
    {
        /** unit: meter */
        static constexpr float_64 WAVE_LENGTH_SI = 0.8e-6;

        /** Convert the normalized laser strength parameter a0 to Volt per meter */
        static constexpr float_64 UNITCONV_A0_to_Amplitude_SI = -2.0 * PI / WAVE_LENGTH_SI * ::picongpu::SI::ELECTRON_MASS_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI * ::picongpu::SI::SPEED_OF_LIGHT_SI / ::picongpu::SI::ELECTRON_CHARGE_SI;

        /** unit: W / m^2 */
        // calculate: _A0 = 8.549297e-6 * sqrt( Intensity[W/m^2] ) * wavelength[m] (linearly polarized)

        /** unit: none */
        static constexpr float_64 _A0  = 1.5;

        /** unit: Volt / meter */
        static constexpr float_64 AMPLITUDE_SI = _A0 * UNITCONV_A0_to_Amplitude_SI;

        /** unit: Volt / meter */
        //static constexpr float_64 AMPLITUDE_SI = 1.738e13;

        /** Stretch temporal profile by a constant plateau between the up and downramp
         *  unit: seconds */
        static constexpr float_64 LASER_NOFOCUS_CONSTANT_SI = 13.34e-15;

        /** Pulse length: sigma of std. gauss for intensity (E^2)
         *  PULSE_LENGTH_SI = FWHM_of_Intensity   / [ 2*sqrt{ 2* ln(2) } ]
         *                                          [    2.354820045     ]
         *  Info:             FWHM_of_Intensity = FWHM_Illumination
         *                      = what a experimentalist calls "pulse duration"
         *  unit: seconds (1 sigma) */
        static constexpr float_64 PULSE_LENGTH_SI = 10.615e-15 / 4.0;

        /** cell from top where the laser is initialized
         *
         * if `initPlaneY == 0` than the absorber are disabled.
         * if `initPlaneY > absorbercells negative Y` the negative absorber in y
         * direction is enabled
         *
         * valid ranges:
         *   - initPlaneY == 0
         *   - absorber cells negative Y < initPlaneY < cells in y direction of the top gpu
         */
        static constexpr uint32_t initPlaneY = 0;

        /** The laser pulse will be initialized half of PULSE_INIT times of the PULSE_LENGTH before and after the plateau
         *  unit: none */
        static constexpr float_64 RAMP_INIT = 20.6146;

        /** laser phase shift (no shift: 0.0)
         *
         * sin(omega*time + laser_phase): starts with phase=0 at center --> E-field=0 at center
         *
         * unit: rad, periodic in 2*pi
         */
        static constexpr float_X LASER_PHASE = 0.0;

        /** Available polarization types
         */
        enum PolarisationType
        {
            LINEAR_X = 1u,
            LINEAR_Z = 2u,
            CIRCULAR = 4u,
        };
        /** Polarization selection
         */
        static constexpr PolarisationType Polarisation = LINEAR_X;
    };
} // namespace defaults

None

template<typename T_Params>
struct None : public picongpu::fields::laserProfiles::none::Unitless<T_Params>

Empty laser profile.

Does not define a laser profile but provides some hard-coded constants that are accessed directly in some places.

Template Parameters
  • T_Params: class parameter to configure the “no laser” profile, see members of none::defaults::NoneParam for required members