htrdrPy.data module

class htrdrPy.data.Data(radius, nTheta=None, nPhi=None, mass=None, gravity=None, name='')[source]

Bases: object

htrdrPy.Data is a class aiming to handle the optical and physical properties of the system and to create the input files for htrdr.

Examples

The first step is the creation of a instance of htrdrPy.Data:

>>> d = htrdrPy.Data(radius=1e6, nTheta=30, nPhi=50, name="Planet")

The next step is to provide the physical and radiative properties of the atmosphere and ground. Different methods exist depending on the case considered. In the following, we consider a 1D set of data forming an horizontally homogeneous planet.

>>> nLevel = 50
>>> nCoeff = 4
>>> nWavelength = 20
>>> weights = np.array(nWavelength * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelength, nCoeff)
>>> altitudes = np.linspace(0, 5e5, nLevel)
>>> temperatures = np.linspace(300, 500, nLevel)
>>> scatt = np.linspace(1e-8, 1e-2,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> absor = np.linspace(1e-5, 1e-1,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> asymm = np.linspace(0, 1,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:] = bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> surfTemp= 300
>>> surfAlb= np.ones(nWavelength) * 0.5
>>> d.makeAtmosphereFrom1D({
...         "nLevel": nLevel,
...         "nCoeff": nCoeff,
...         "nWavelength": nWavelength,
...         "weights": weights,
...         "altitude (m)": altitudes,
...         "temperature (K)": temperatures,
...         "scattering (m-1)": scatt,
...         "absorption (m-1)": absor,
...         "assymetry": asymm,
...         "wavelength": wavelengths,
...         "band low": bandsLow,
...         "band up": bandsUp
...         })
Mesh generator. Ntheta = 30, Nphi = 50, Nz = 50, r_min = 1000000.0, r_max = 1500000.0
Generating points...
Generating nodes...
Hexahedron & Octahedron generation completed. N_Hexahedron = 4802, N_Octahedron = 64827
Generating Tetrahedrons...
Assigning data to nodes ...
>>> d.makeGroundFrom1D(surfTemp, {
...         "kind": "lambertian",
...         "albedo": surfAlb,
...         "bands": np.array([bandsLow, bandsUp]).T
...         })
Mesh generator. Ntheta = 30, Nphi = 50, R = 1000000.0
Generating points...
Generating nodes...
triangles,rectangles generation completed. N_triangles = 98, N_rectangles = 1323
Generating triangles...
generating bin...

Once the data have been provided, the method htrdrPy.Data.writeInputs will generate the input file in an input_{name} folder.

>>> d.writeInputs()
generating surface mesh bin file...
Nnodes = 5586 Ncells = 2744 dim_node = 3 dim_cell = 3
bin generation completed.
generating surface properties bin file...
bin generation completed.
generating atmosphere mesh bin file...
Nnodes = 547428 Ncells = 338541 dim_node = 3 dim_cell = 4
bin generation completed.
generating gas temperature bin file...
bin generation completed.
generating gas properties bin file...
100%|_________________________________________________________________________________________________________________________________________________________| 20/20 [00:11<00:00,  1.69it/s]
bin generation completed.
generating haze properties bin file...
100%|_________________________________________________________________________________________________________________________________________________________| 20/20 [00:04<00:00,  4.95it/s]
bin generation completed.
generating haze phase function bin file...
bin generation completed.
Parameters:
  • radius (float) – Radius of the planet [m].

  • nTheta (int, optional) – Number of latitude points in the range [0°,360°] to be used. Not required if the tables provided are 2D or 3D already.

  • nPhi (int, optional) – Number of latitude points in the range [-180°,180°]. Not required if the tables provided are 3D already.

  • gravity (float, optional,) – Gravity of the planet [m/s2]. Requested only for LMDZ-PCM input/output files to extract the altitude grid from the geopotential. Alternatively, the user can provide the planet mass.

  • mass (float, optional,) – Mass of the planet [kg]. Requested only for LMDZ-PCM input/output files to extract the altitude grid from the geopotential. Alternatively, the user can provide the planet gravity.

  • name (str, optional, default uses a counter of the number of istances of htrdrPy.Data) – Name for the dataset, which will be used to name the input folders.

addParticle(dim, data, sameGrid=True, nTheta=None, nPhi=None)[source]

Add a particle (haze, cloud, …) in the atmosphere.

Parameters:
  • dim (int) – Dimension of the original data. If set to 0, the case is considered as a plan-parallel.

  • data (dict) –

    Dictionnary with the following items: - “nAngle”: int, optional

    Number of angles in the phase functions (not required if using the builtin Henyey-Greenstein phase function).

    • ”scattering (m-1)”: numpy.ndarray

      Array of scattering coefficients (shape=(nWavelength, nLevel, [nLat, nLon]), [m-1]).

    • ”absorption (m-1)”: numpy.ndarray

      Array of absorption coefficients (shape=(nWavelength, nLevel, [nLat, nLon]), [m-1]).

    • ”angles (°)”: numpy.ndarray, optional

      Array of angle values for the phase function (shape=(nAngle), float [°]). Only required when providing the phaseFunc.

    • ”phaseFunc”: numpy.ndarray, optional

      Array of discrete phase functions (shape = (nWavelength, nLevel, [nLat, nLon], nAngle)). Alternatively, the user can provide an array of asymmetry parameter.

    • ”asymmetry”: numpy.ndarray, optional

      Array of asymmetry parameter (shape=(nWavelength, nLevel, [nLat, nLon])). Alternatively, the user can provide the discrete phase function.

    • ”nWavelength”: int, optional

      Number of wavelengths. Required only if sameGrid is False.

    • ”wavelength”: numpy.ndarray, optional

      Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the band center. The value is used of the phase function. Required only if sameGrid is False.

    • ”band low”: numpy.ndarray, optional

      Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band. Required only if sameGrid is False.

    • ”band up”: numpy.ndarray, optional

      Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band. Required only if sameGrid is False.

    • ”nLevel”int, optional

      Number of levels. Required only if sameGrid is False.

    • ”nLat”int, optional

      Number of latitudes. Required only if sameGrid is False.

    • ”nLon”int, optional

      Number of longitudes. Required only if sameGrid is False.

    • ”altitude (m)”: numpy.ndarray, optional

      Array of altitudes (shape=(nLevel, nLat, nLon), [m]). Required only if sameGrid is False.

    • ”latitude (°)”numpy.ndarray, optional

      List of latitudes (shape=(nLat), [°]). Required only if sameGrid is False.

    • ”longitude (°)”numpy.ndarray, optional

      List of longitudes (shape=(nLon), [°]). Required only if sameGrid is False.

  • sameGrid (bool, default True) – Wheteher to use the same grid as for the gas. If set ti True, many keys of the data parameter are not required. This is meant to save computation time when haze and gas share the same grid.

  • nTheta (int, optional) – Number of latitude points in the range [0°,360°] to be used. Not requested if sameGrid is True, or dim is 2 or 3.

  • nPhi (int, optional, not requested if meshes loaded from file) – Number of latitude points in the range [-180°,180°]. Not requested if sameGrid is True, or dim is 3.

cleanInputs()[source]

Remove the inputs_{name} folder.

makeAtmosphereFrom1D(data)[source]

Generate a spherical atmosphere from single column data. The data are provided at the levels, i.e. at the interface between layers.

Parameters:

data (dict) –

Dictionnary with the following items:

  • ”nLevel”int

    Number of levels.

  • ”nCoeff”: int

    Number of quadrature points for the k-coeff.

  • ”nWavelength”: int

    Number of wavelengths.

  • ”weights”: numpy.ndarray

    The weigths of the nCoeff quadrature points (shape=(nCoeff)).

  • ”altitude (m)”: numpy.ndarray

    Array of altitudes (shape=(nLevel), [m]).

  • ”temperature (K)”: numpy.ndarray

    Array of temperatures (shape=(nLevel), [K]).

  • ”scattering (m-1)”: numpy.ndarray

    Array of scattering coefficients (shape=(nWavelength, nLevel), [m-1]).

  • ”absorption (m-1)”: numpy.ndarray

    Array of absorption coefficients (shape=(nWavelength, nLevel, nCoeff), [m-1]).

  • ”wavelength”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the band center. The value is used of the phase function.

  • ”band low”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

  • ”band up”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

Notes

The nTheta and nPhi parameters are used to define the resolution of the atmosphere mesh and are therefore mandatory.

Examples

>>> d = htrdrPy.Data(radius=1e6, nTheta=30, nPhi=50, name="Planet")
>>> nLevel = 50
>>> nCoeff = 4
>>> nWavelength = 20
>>> weights = np.array(nWavelength * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelength, nCoeff)
>>> altitudes = np.linspace(0, 5e5, nLevel)
>>> temperatures = np.linspace(300, 500, nLevel)
>>> scatt = np.linspace(1e-8, 1e-2,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> absor = np.linspace(1e-5, 1e-1,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> asymm = np.linspace(0, 1,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:] = bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> d.makeAtmosphereFrom1D({
...         "nLevel": nLevel,
...         "nCoeff": nCoeff,
...         "nWavelength": nWavelength,
...         "weights": weights,
...         "altitude (m)": altitudes,
...         "temperature (K)": temperatures,
...         "scattering (m-1)": scatt,
...         "absorption (m-1)": absor,
...         "asymmetry": asymm,
...         "wavelength": wavelengths,
...         "band low": bandsLow,
...         "band up": bandsUp
...         })
Mesh generator. Ntheta = 30, Nphi = 50, Nz = 50, r_min = 1000000.0, r_max = 1500000.0
Generating points...
Generating nodes...
Hexahedron & Octahedron generation completed. N_Hexahedron = 4802, N_Octahedron = 64827
Generating Tetrahedrons...
Assigning data to nodes ...
makeAtmosphereFrom1D_PP(data)[source]

Generate a plan parallel atmosphere from single column data. The data are provided at the levels, i.e. at the interface between layers.

Parameters:

data (dict) –

Dictionnary with the following items:

  • ”nLevel”int

    Number of levels.

  • ”nCoeff”: int

    Number of quadrature points for the k-coeff.

  • ”nWavelength”: int

    Number of wavelengths.

  • ”weights”: numpy.ndarray

    The weigths of the nCoeff quadrature points (shape=(nCoeff)).

  • ”altitude (m)”: numpy.ndarray

    Array of altitudes (shape=(nLevel), [m]).

  • ”temperature (K)”: numpy.ndarray

    Array of temperatures (shape=(nLevel), [K]).

  • ”scattering (m-1)”: numpy.ndarray

    Array of scattering coefficients (shape=(nWavelength, nLevel), [m-1]).

  • ”absorption (m-1)”: numpy.ndarray

    Array of absorption coefficients (shape=(nWavelength, nLevel, nCoeff), [m-1]).

  • ”wavelength”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the band center. The value is used of the phase function.

  • ”band low”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

  • ”band up”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

Notes

In plan-parallel mode, the radius parameter passed at the initialisation of the instance is used as the horizontal expansion of the ground. Make sure to use a large enough value. Also, the nTheta and nPhi parameters are not used.

Examples

>>> d = htrdrPy.Data(radius=1e6, nTheta=30, nPhi=50, name="Planet")
>>> nLevel = 50
>>> nCoeff = 4
>>> nWavelength = 20
>>> weights = np.array(nWavelength * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelength, nCoeff)
>>> altitudes = np.linspace(0, 5e5, nLevel)
>>> temperatures = np.linspace(300, 500, nLevel)
>>> scatt = np.linspace(1e-8, 1e-2,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> absor = np.linspace(1e-5, 1e-1,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> asymm = np.linspace(0, 1,
...     nLevel*nCoeff*nWavelength).reshape((nWavelength, nLevel, nCoeff))
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:] = bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> d.makeAtmosphereFrom1D_PP({
...         "nLevel": nLevel,
...         "nCoeff": nCoeff,
...         "nWavelength": nWavelength,
...         "weights": weights,
...         "altitude (m)": altitudes,
...         "temperature (K)": temperatures,
...         "scattering (m-1)": scatt,
...         "absorption (m-1)": absor,
...         "asymmetry": asymm,
...         "wavelength": wavelengths,
...         "band low": bandsLow,
...         "band up": bandsUp
...         })
Assigning data to nodes ...
makeAtmosphereFrom2D(data)[source]

Generate a spherical atmosphere from single column data. The data are provided at the levels, i.e. at the interface between layers.

Parameters:

data (dict) –

Dictionnary with the following items:

  • ”nLevel”int

    Number of levels.

  • ”nLat”int

    Number of latitudes.

  • ”nCoeff”: int

    Number of quadrature points for the k-coeff.

  • ”nWavelength”: int

    Number of wavelengths.

  • ”weights”: numpy.ndarray

    The weigths of the nCoeff quadrature points (shape=(nCoeff)).

  • ”altitude (m)”: numpy.ndarray

    Array of altitudes (shape=(nLevel, nLat), [m]).

  • ”latitude (°)”numpy.ndarray

    List of latitudes (shape=(nLat), [°]).

  • ”temperature (K)”: numpy.ndarray

    Array of temperatures (shape=(nLevel, nLat), [K]).

  • ”scattering (m-1)”: numpy.ndarray

    Array of scattering coefficients (shape=(nWavelength, nLevel, nLat), [m-1]).

  • ”absorption (m-1)”: numpy.ndarray

    Array of absorption coefficients (shape=(nWavelength, nLevel, nLat, nCoeff), [m-1]).

  • ”wavelength”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the band center. The value is used of the phase function.

  • ”band low”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

  • ”band up”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

Notes

The nPhi parameter provided at the initialisation of the instance is used to define the longitudinal resolution of the atmospheric mesh and is therefore manatory. The nTheta parameter is obtained from the length of the latitude provided.

Examples

>>> d = htrdrPy.Data(radius=1e6, nPhi=50, name="Planet")
>>> nLevel = 50
>>> nLat = 30
>>> nCoeff = 4
>>> nWavelength = 20
>>> weights = np.array(nWavelength * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelength, nCoeff)
>>> altitudes = np.tile(np.linspace(0, 5e5, nLevel), (nLat, 1)).T
>>> latitudes = np.linspace(-90, 90, nLat)
>>> temperatures = np.linspace(300, 500, nLevel*nLat).reshape(nLevel,nLat)
>>> scatt = np.linspace(1e-8, 1e-2,
...     nLevel*nLat*nCoeff*nWavelength).reshape((nWavelength, nLevel,
...                                                 nLat, nCoeff))
>>> absor = np.linspace(1e-5, 1e-1,
...     nLevel*nLat*nCoeff*nWavelength).reshape((nWavelength, nLevel,
...                                                 nLat, nCoeff))
>>> asymm = np.linspace(0, 1,
...     nLevel*nLat*nCoeff*nWavelength).reshape((nWavelength, nLevel,
...                                                 nLat, nCoeff))
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:] = bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> d.makeAtmosphereFrom2D({
...         "nLevel": nLevel,
...         "nLat": nLat,
...         "nCoeff": nCoeff,
...         "nWavelength": nWavelength,
...         "weights": weights,
...         "altitude (m)": altitudes,
...         "latitude (°)": latitudes,
...         "temperature (K)": temperatures,
...         "scattering (m-1)": scatt,
...         "absorption (m-1)": absor,
...         "asymmetry": asymm,
...         "wavelength": wavelengths,
...         "band low": bandsLow,
...         "band up": bandsUp
...         })
makeAtmosphereFrom3D(data)[source]

Generate a spherical atmosphere from single column data. The data are provided at the levels, i.e. at the interface between layers.

Parameters:

data (dict) –

Dictionnary with the following items:

  • ”nLevel”int

    Number of levels.

  • ”nLat”int

    Number of latitudes.

  • ”nLon”int

    Number of longitudes.

  • ”nCoeff”: int

    Number of quadrature points for the k-coeff.

  • ”nWavelength”: int

    Number of wavelengths.

  • ”weights”: numpy.ndarray

    The weigths of the nCoeff quadrature points (shape=(nCoeff)).

  • ”altitude (m)”: numpy.ndarray

    Array of altitudes (shape=(nLevel, nLat, nLon), [m]).

  • ”latitude (°)”numpy.ndarray

    List of latitudes (shape=(nLat), [°]).

  • ”longitude (°)”numpy.ndarray

    List of longitudes (shape=(nLon), [°]).

  • ”temperature (K)”: numpy.ndarray

    Array of temperatures (shape=(nLevel, nLat, nLon), [K]).

  • ”scattering (m-1)”: numpy.ndarray

    Array of scattering coefficients (shape=(nWavelength, nLevel, nLat, nLon), [m-1]).

  • ”absorption (m-1)”: numpy.ndarray

    Array of absorption coefficients (shape=(nWavelength, nLevel, nLat, nLon, nCoeff), [m-1]).

  • ”wavelength”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the band center. The value is used of the phase function.

  • ”band low”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

  • ”band up”: numpy.ndarray

    Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

Notes

The nTheta and nPhi parameters provided at the initialisation of the instance are not used and instead are defined from the length of the latitude (°) and longitude (°), respectively.

Examples

>>> d = htrdrPy.Data(radius=1e6, name="Planet")
>>> nLevel = 50
>>> nLat = 30
>>> nLon = 50
>>> nCoeff = 4
>>> nWavelength = 20
>>> weights = np.array(nWavelength * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelength, nCoeff)
>>> altitudes = np.moveaxis(np.tile(np.linspace(0, 5e5, nLevel),
...                                 (nLat, nLon, 1)), -1, 0)
>>> latitudes = np.linspace(-90, 90, nLat)
>>> longitudes = np.linspace(-180, 180, nLon)
>>> temperatures = np.linspace(300, 500, nLevel*nLat*nLon).reshape(nLevel,nLat, nLon)
>>> scatt = np.linspace(1e-8, 1e-2,
...     nLevel*nLat*nLon*nCoeff*nWavelength).reshape((nWavelength, nLevel,
...                                                 nLat, nLon, nCoeff))
>>> absor = np.linspace(1e-5, 1e-1,
...     nLevel*nLat*nLon*nCoeff*nWavelength).reshape((nWavelength, nLevel,
...                                                 nLat, nLon, nCoeff))
>>> asymm = np.linspace(0, 1,
...     nLevel*nLat*nLon*nCoeff*nWavelength).reshape((nWavelength, nLevel,
...                                                 nLat, nLon, nCoeff))
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:] = bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> d.makeAtmosphereFrom3D({
...         "nLevel": nLevel,
...         "nLat": nLat,
...         "nLon": nLon,
...         "nCoeff": nCoeff,
...         "nWavelength": nWavelength,
...         "weights": weights,
...         "altitude (m)": altitudes,
...         "latitude (°)": latitudes,
...         "longitude (°)": longitudes,
...         "temperature (K)": temperatures,
...         "scattering (m-1)": scatt,
...         "absorption (m-1)": absor,
...         "asymmetry": asymm,
...         "wavelength": wavelengths,
...         "band low": bandsLow,
...         "band up": bandsUp
...         })
makeGroundFrom1D(surfaceTemperature, brdf)[source]

Generates a spherical ground considering uniform temperature and optical properties.

Parameters:
  • surfaceTemperature (float) – Temperature of the surface [K].

  • brdf (dict) –

    Surface reflexion properties with the following items:

    • ”kind”{“lambertian”, “specular”})

      Kind of brdf function to use.

    • ”albedo”numpy.ndarray

      Wavelength dependent surface albedos (shape=(nWavelength)) .

    • ”wavelengths”numpy.ndarray, optional

      Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.

    • ”bands”numpy.ndarray, optional

      Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.

Notes

The nTheta and nPhi parameters provided at the initialisation of the instance are used to define the resolution of the ground mesh and are therefore manatory.

Examples

>>> d = htrdrPy.Data(radius=1e6, nTheta=30, nPhi=50, name="Planet")
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:], bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> surfTemp= 300
>>> surfAlb= np.ones(nWavelength) * 0.5
>>> d.makeGroundFrom1D(surfTemp, {
...         "kind": "lambertian",
...         "albedo": surfAlb,
...         "bands": np.array([bandsLow, bandsUp]).T
...         })
Mesh generator. Ntheta = 30, Nphi = 50, R = 1000000.0
Generating points...
Generating noeds...
triangles,rectangles generation completed. N_triangles = 98, N_rectangles = 1323
Generating triangles...
generating bin...
makeGroundFrom1D_PP(surfaceTemperature, brdf)[source]

Generates a plan-parallel ground considering uniform temperature and optical properties.

Parameters:
  • surfaceTemperature (float) – Temperature of the surface [K].

  • brdf (dict) –

    Surface reflexion properties with the following items:

    • ”kind”{“lambertian”, “specular”})

      Kind of brdf function to use.

    • ”albedo”numpy.ndarray

      Wavelength dependent surface albedos (shape=(nWavelength)) .

    • ”wavelengths”numpy.ndarray, optional

      Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.

    • ”bands”numpy.ndarray, optional

      Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.

Notes

In plan-parallel mode, the radius parameter passed at the initialisation of the instance is used as the horizontal expansion of the ground. Make sure to use a large enough value. Also, the nTheta and nPhi parameters are not used.

Examples

>>> d = htrdrPy.Data(radius=1e6, name="Planet")
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:], bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> surfTemp= 300
>>> surfAlb= np.ones(nWavelength) * 0.5
>>> d.makeGroundFrom1D_PP(surfTemp, {
...         "kind": "lambertian",
...         "albedo": surfAlb,
...         "bands": np.array([bandsLow, bandsUp]).T
...         })
generating bin...
makeGroundFrom2D(surfaceTemperature, brdf)[source]

Generates a spherical ground considering temperature and optical properties varying along the latitude.

Parameters:
  • surfaceTemperature (numpy.ndarray) – Temperature of the surface (shape=(nLat), [K]).

  • brdf (dict) –

    Surface reflexion properties with the following items:

    • ”kind”{“lambertian”, “specular”})

      Kind of brdf function to use.

    • ”albedo”numpy.ndarray

      Wavelength dependent surface albedos (shape=(nWavelength, nLat)).

    • ”latitude”numpy.ndarray

      List of latitudes (shape=(nLat), [°]).

    • ”wavelengths”numpy.ndarray, optional

      Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.

    • ”bands”numpy.ndarray, optional

      Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.

Notes

The nPhi parameter provided at the initialisation of the instance is used to define the longitudinal resolution of the ground mesh and is therefore manatory. The nTheta parameter is obtained from the length of the latitude provided.

Examples

>>> d = htrdrPy.Data(radius=1e6, nPhi=50, name="Planet")
>>> nLat = 30
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:], bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> latitudes = np.linspace(-90, 90, nLat)
>>> surfTemp = 300 * np.ones_like(latitudes)
>>> surfAlb= np.ones((nWavelength, nLat)) * 0.5
>>> d.makeGroundFrom2D(surfTemp, {
...         "kind": "lambertian",
...         "albedo": surfAlb,
...         "latitude": latitudes,
...         "bands": np.array([bandsLow, bandsUp]).T
...         })
Mesh generator. Ntheta = 30, Nphi = 50, R = 1000000.0
Generating points...
Generating noeds...
triangles,rectangles generation completed. N_triangles = 98, N_rectangles = 1323
Generating triangles...
generating bin...
makeGroundFrom3D(surfaceTemperature, brdf)[source]

Generates a spherical ground considering temperature and optical properties varying along the latitude and longitude.

Parameters:
  • surfaceTemperature (numpy.ndarray) – Temperature of the surface (shape=(nLat, nLon), [K]).

  • brdf (dict) –

    Surface reflexion properties with the following items:

    • ”kind”{“lambertian”, “specular”})

      Kind of brdf function to use.

    • ”albedo”numpy.ndarray

      Wavelength dependent surface albedos (shape=(nWavelength, nLat, nLon)).

    • ”latitude”numpy.ndarray

      List of latitudes (shape=(nLat), [°]).

    • ”longitude”numpy.ndarray

      List of longitudes (shape=(nLon), [°]).

    • ”wavelengths”numpy.ndarray, optional

      Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.

    • ”bands”numpy.ndarray, optional

      Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.

Notes

The nTheta and nPhi parameters provided at the initialisation of the instance are not used and instead are defined from the length of the latitude and longitude, respectively.

Examples

>>> d = htrdrPy.Data(radius=1e6, name="Planet")
>>> nLat = 30
>>> nLon = 50
>>> wavelengths = np.linspace(2e7, 9e7, nWavelength)
>>> bandsLow = np.zeros(nWavelength)
>>> bandsUp = np.zeros(nWavelength)
>>> bandsLow[1:], bandsUp[:-1] = 0.5 * (wavelengths[1:] + wavelengths[:-1])
>>> bandsLow[0] = 1.5 * wavelengths[0] - 0.5 * wavelengths[1]
>>> bandsUp[-1] =  1.5 * wavelengths[-1] - 0.5 * wavelengths[-2]
>>> latitudes = np.linspace(-90, 90, nLat)
>>> longitudes = np.linspace(-180, 180, nLon)
>>> surfTemp = 300 * np.ones((nLat, nLon))
>>> surfAlb= np.ones((nWavelength, nLat, nLon)) * 0.5
>>> d.makeGroundFrom3D(surfTemp, {
...         "kind": "lambertian",
...         "albedo": surfAlb,
...         "latitude": latitudes,
...         "longitude": longitudes,
...         "bands": np.array([bandsLow, bandsUp]).T
...         })
Mesh generator. Ntheta = 30, Nphi = 50, R = 1000000.0
Generating points...
Generating noeds...
triangles,rectangles generation completed. N_triangles = 98, N_rectangles = 1323
Generating triangles...
generating bin...
makeMixture(data, dim)[source]

Generates the atmosphere from a mixture of gas and particles.

Parameters:
  • data (dict) –

    Dictionnary with the following items: - “nLevel” : int

    Number of levels. Required only if sameGrid is False.

    • ”nCoeff”: int

      Number of quadrature points for the k-coeff.

    • ”nLat”int, optional

      Number of latitudes. Required only if dim is 2 or 3.

    • ”nLon”int, optional

      Number of longitudes. Required only if dim is 3.

    • ”nAngle”: int, optional

      Number of angles in the phase functions (not required if using the builtin Henyey-Greenstein phase function).

    • ”altitude (m)”: numpy.ndarray

      Array of altitudes (shape=(nLevel, nLat, nLon), [m]).

    • ”latitude (°)”numpy.ndarray, optional

      List of latitudes (shape=(nLat), [°]). Required only if dim is 2 or 3.

    • ”longitude (°)”numpy.ndarray, optional

      List of longitudes (shape=(nLon), [°]). Required only if dim is 3.

    • ”nWavelength”: int

      Number of wavelengths.

    • ”wavelength”: numpy.ndarray

      Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the band center. The value is used of the phase function.

    • ”band low”: numpy.ndarray

      Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

    • ”band up”: numpy.ndarray

      Array of wavelength (shape = (nWavelength), [m]). The values corresponds to the lower boundary of the band.

    • ”scattering (m-1)”: numpy.ndarray

      Array of scattering coefficients (shape=(nWavelength, nLevel, [nLat, nLon]), [m-1]).

    • ”absorption (m-1)”: numpy.ndarray

      Array of absorption coefficients (shape=(nWavelength, nLevel, [nLat, nLon], nCoeff), [m-1]).

    • ”angles (°)”: numpy.ndarray, optional

      Array of angle values for the phase function (shape=(nAngle), float [°]). Only required when providing the phaseFunc.

    • ”phaseFunc”: numpy.ndarray, optional

      Array of discrete phase functions (shape = (nWavelength, nLevel, [nLat, nLon], nAngle)). Alternatively, the user can provide an array of asymmetry parameter.

    • ”asymmetry”: numpy.ndarray, optional

      Array of asymmetry parameter (shape=(nWavelength, nLevel, [nLat, nLon])). Alternatively, the user can provide the discrete phase function.

  • dim (int) – Dimension of the original data. If set to 0, the case is considered as a plan-parallel.

writeInputAtmosphere()[source]

Write the atmosphere binary input files for htrdr-planets

writeInputGround()[source]

Write the ground binary input files for htrdr-planets

writeInputs(octree_def=512, opthick=1, nthOctree=8, procOctree='master', octreeFile='')[source]

Write the binary input files for htrdr-planets and precalculate octrees if the octreeFile is passed.

Parameters:
  • octree_def (int or str, default 512) – Maximal defintion of the octree grid.

  • opthick (float or str, default 1) – Optical thickness threshold to assess the merge of cells.

  • nthOctree (int, default 8) – Number of threads to use for octree computation.

  • procOctree ({'all', 'master'}, default 'master') – Which process must realize the octree calculation (useless if no storage). Put to ‘all’ if the processes do not share the disk space.

  • octreeFile (str, optional) – Filename to use for the stroring of octrees. If not provided, octrees are not stored on disk. If a file is provided, a “blank run” of htrdr-planets is launched to precalculate octrees. The file is located in the outputs_{name} folder.

Examples

>>> d = htrdrPy.Data(radius=1e6, nPhi=50, name="Planet")
...
>>> d.writeInputs(octreeFile="octree.bin")

This will generate the input files and start a very quick run of htrdr-planets to generate an octree file that will be stored in outputs_Planet/octree.bin

>>> d.writeInputs(octree_def=1024)

This will generate the input files and keep the information on the octree_def for the future calculation of the octrees. The octrees are not precalculated since they will only be stored in memory.

writeVTKfiles()[source]

Write the VTK and obj files for view in paraview or other visualisation software handling obj and VTK files.

writeVTKfilesAtmosphere()[source]

Write the atmosphere VTK and obj files.

writeVTKfilesGround()[source]

Write the ground VTK and obj files.