htrdrPy.data module
- class htrdrPy.data.Data(radius, nTheta=None, nPhi=None, mass=None, gravity=None, name='')[source]
Bases:
objecthtrdrPy.Datais 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 >>> nWavelengths = 20 >>> weights = np.array(nWavelengths * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelengths, nCoeff) >>> altitudes = np.linspace(0, 5e5, nLevel) >>> temperatures = np.linspace(300, 500, nLevel) >>> scatt = np.linspace(1e-8, 1e-2, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> absor = np.linspace(1e-5, 1e-1, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> asymm = np.linspace(0, 1, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> wavelengths = np.linspace(2e7, 9e7, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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(nWavelengths) * 0.5 >>> d.makeAtmosphereFrom1D({ ... "nLevel": nLevel, ... "nCoeff": nCoeff, ... "nWavelengths": nWavelengths, ... "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.writeInputswill 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, not requested if meshes loaded from file) – Number of latitude points in the range [0°,360°] to be used.
nPhi (int, optional, not requested if meshes loaded from file) – Number of latitude points in the range [-180°,180°].
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.
- 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.
- ”nWavelengths”: int
Number of wavelengths.
- ”nAngle”: int, optional
Number of angles in the phase functions (not required if using the builtin Henyey-Greenstein phase function).
- ”weights”:
numpy.ndarray The weigths of the nCoeff quadrature points (shape=(nWavelengths, nCoeff)).
- ”weights”:
- ”altitude (m)”:
numpy.ndarray Array of altitudes (shape=(nLevel), [m]).
- ”altitude (m)”:
- ”temperature (K)”:
numpy.ndarray Array of temperatures (shape=(nLevel), [K]).
- ”temperature (K)”:
- ”scattering (m-1)”:
numpy.ndarray Array of scattering coefficients (shape=(nWavelengths, nLevel, nCoeff), [m-1]).
- ”scattering (m-1)”:
- ”absorption (m-1)”:
numpy.ndarray Array of absorption coefficients (shape=(nWavelengths, nLevel, nCoeff), [m-1]).
- ”absorption (m-1)”:
- ”angles (°)”:
numpy.ndarray, optional Array of angle values for the phase function (shape=(nAngle), float [°]). Only required when providing the
phaseFunc.
- ”angles (°)”:
- ”phaseFunc”:
numpy.ndarray, optional Array of discrete phase functions (shape = (nWavelengths, nLevel, nAngle, nCoeff)). Alternatively, the user can provide an array of asymmetry parameter.
- ”phaseFunc”:
- ”asymmetry”:
numpy.ndarray, optional Array of asymmetry parameter (shape=(nWavelengths, nLevel, nCoeff)). Alternatively, the user can provide the discrete phase function.
- ”asymmetry”:
- ”wavelength”:
numpy.ndarray, optional Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the band center. The value is used of the phase function.
- ”wavelength”:
- ”band low”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band low”:
- ”band up”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band up”:
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 >>> nWavelengths = 20 >>> weights = np.array(nWavelengths * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelengths, nCoeff) >>> altitudes = np.linspace(0, 5e5, nLevel) >>> temperatures = np.linspace(300, 500, nLevel) >>> scatt = np.linspace(1e-8, 1e-2, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> absor = np.linspace(1e-5, 1e-1, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> asymm = np.linspace(0, 1, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> wavelengths = np.linspace(2e7, 9e7, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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, ... "nWavelengths": nWavelengths, ... "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.
- ”nWavelengths”: int
Number of wavelengths.
- ”nAngle”: int, optional
Number of angles in the phase functions (not required if using the builtin Henyey-Greenstein phase function).
- ”weights”:
numpy.ndarray The weigths of the nCoeff quadrature points (shape=(nWavelengths, nCoeff)).
- ”weights”:
- ”altitude (m)”:
numpy.ndarray Array of altitudes (shape=(nLevel), [m]).
- ”altitude (m)”:
- ”temperature (K)”:
numpy.ndarray Array of temperatures (shape=(nLevel), [K]).
- ”temperature (K)”:
- ”scattering (m-1)”:
numpy.ndarray Array of scattering coefficients (shape=(nWavelengths, nLevel, nCoeff), [m-1]).
- ”scattering (m-1)”:
- ”absorption (m-1)”:
numpy.ndarray Array of absorption coefficients (shape=(nWavelengths, nLevel, nCoeff), [m-1]).
- ”absorption (m-1)”:
- ”angles (°)”:
numpy.ndarray, optional Array of angle values for the phase function (shape=(nAngle), float [°]). Only required when providing the
phaseFunc.
- ”angles (°)”:
- ”phaseFunc”:
numpy.ndarray, optional Array of discrete phase functions (shape = (nWavelengths, nLevel, nAngle, nCoeff)). Alternatively, the user can provide an array of asymmetry parameter.
- ”phaseFunc”:
- ”asymmetry”:
numpy.ndarray, optional Array of asymmetry parameter (shape=(nWavelengths, nLevel, nCoeff)). Alternatively, the user can provide the discrete phase function.
- ”asymmetry”:
- ”wavelength”:
numpy.ndarray, optional Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the band center. The value is used of the phase function.
- ”wavelength”:
- ”band low”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band low”:
- ”band up”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band up”:
Notes
In plan-parallel mode, the
radiusparameter 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 >>> nWavelengths = 20 >>> weights = np.array(nWavelengths * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelengths, nCoeff) >>> altitudes = np.linspace(0, 5e5, nLevel) >>> temperatures = np.linspace(300, 500, nLevel) >>> scatt = np.linspace(1e-8, 1e-2, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> absor = np.linspace(1e-5, 1e-1, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> asymm = np.linspace(0, 1, ... nLevel*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, nCoeff)) >>> wavelengths = np.linspace(2e7, 9e7, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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, ... "nWavelengths": nWavelengths, ... "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.
- ”nWavelengths”: int
Number of wavelengths.
- ”nAngle”: int, optional
Number of angles in the phase functions (not required if using the builtin Henyey-Greenstein phase function).
- ”weights”:
numpy.ndarray The weigths of the nCoeff quadrature points (shape=(nWavelengths, nCoeff)).
- ”weights”:
- ”altitude (m)”:
numpy.ndarray Array of altitudes (shape=(nLevel, nLat), [m]).
- ”altitude (m)”:
- ”latitude (°)”
numpy.ndarray List of latitudes (shape=(nLat), [°]).
- ”latitude (°)”
- ”temperature (K)”:
numpy.ndarray Array of temperatures (shape=(nLevel, nLat), [K]).
- ”temperature (K)”:
- ”scattering (m-1)”:
numpy.ndarray Array of scattering coefficients (shape=(nWavelengths, nLevel, nLat, nCoeff), [m-1]).
- ”scattering (m-1)”:
- ”absorption (m-1)”:
numpy.ndarray Array of absorption coefficients (shape=(nWavelengths, nLevel, nLat, nCoeff), [m-1]).
- ”absorption (m-1)”:
- ”angles (°)”:
numpy.ndarray, optional Array of angle values for the phase function (shape=(nAngle), float [°]). Only required when providing the
phaseFunc.
- ”angles (°)”:
- ”phaseFunc”:
numpy.ndarray, optional Array of discrete phase functions (shape = (nWavelengths, nLevel, nLat, nAngle, nCoeff)). Alternatively, the user can provide an array of asymmetry parameter.
- ”phaseFunc”:
- ”asymmetry”:
numpy.ndarray, optional Array of asymmetry parameter (shape=(nWavelengths, nLevel, nLat, nCoeff)). Alternatively, the user can provide the discrete phase function.
- ”asymmetry”:
- ”wavelength”:
numpy.ndarray, optional Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the band center. The value is used of the phase function.
- ”wavelength”:
- ”band low”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band low”:
- ”band up”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band up”:
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
latitudeprovided.Examples
>>> d = htrdrPy.Data(radius=1e6, nPhi=50, name="Planet") >>> nLevel = 50 >>> nLat = 30 >>> nCoeff = 4 >>> nWavelengths = 20 >>> weights = np.array(nWavelengths * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelengths, 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*nWavelengths).reshape((nWavelengths, nLevel, ... nLat, nCoeff)) >>> absor = np.linspace(1e-5, 1e-1, ... nLevel*nLat*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, ... nLat, nCoeff)) >>> asymm = np.linspace(0, 1, ... nLevel*nLat*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, ... nLat, nCoeff)) >>> wavelengths = np.linspace(2e7, 9e7, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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, ... "nWavelengths": nWavelengths, ... "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.
- ”nWavelengths”: int
Number of wavelengths.
- ”nAngle”: int, optional
Number of angles in the phase functions (not required if using the builtin Henyey-Greenstein phase function).
- ”weights”:
numpy.ndarray The weigths of the nCoeff quadrature points (shape=(nWavelengths, nCoeff)).
- ”weights”:
- ”altitude (m)”:
numpy.ndarray Array of altitudes (shape=(nLevel, nLat, nLon), [m]).
- ”altitude (m)”:
- ”latitude (°)”
numpy.ndarray List of latitudes (shape=(nLat), [°]).
- ”latitude (°)”
- ”longitude (°)”
numpy.ndarray List of longitudes (shape=(nLon), [°]).
- ”longitude (°)”
- ”temperature (K)”:
numpy.ndarray Array of temperatures (shape=(nLevel, nLat, nLon), [K]).
- ”temperature (K)”:
- ”scattering (m-1)”:
numpy.ndarray Array of scattering coefficients (shape=(nWavelengths, nLevel, nLat, nLon, nCoeff), [m-1]).
- ”scattering (m-1)”:
- ”absorption (m-1)”:
numpy.ndarray Array of absorption coefficients (shape=(nWavelengths, nLevel, nLat, nLon, nCoeff), [m-1]).
- ”absorption (m-1)”:
- ”angles (°)”:
numpy.ndarray, optional Array of angle values for the phase function (shape=(nAngle), float [°]). Only required when providing the
phaseFunc.
- ”angles (°)”:
- ”phaseFunc”:
numpy.ndarray, optional Array of discrete phase functions (shape = (nWavelengths, nLevel, nLat, nLon, nAngle, nCoeff)). Alternatively, the user can provide an array of asymmetry parameter.
- ”phaseFunc”:
- ”asymmetry”:
numpy.ndarray, optional Array of asymmetry parameter (shape=(nWavelengths, nLevel, nLat, nLon, nCoeff)). Alternatively, the user can provide the discrete phase function.
- ”asymmetry”:
- ”wavelength”:
numpy.ndarray, optional Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the band center. The value is used of the phase function.
- ”wavelength”:
- ”band low”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band low”:
- ”band up”:
numpy.ndarray Array of wavelength (shape = (nWavelengths), [m]). The values corresponds to the lower boundary of the band.
- ”band up”:
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 (°)andlongitude (°), respectively.Examples
>>> d = htrdrPy.Data(radius=1e6, name="Planet") >>> nLevel = 50 >>> nLat = 30 >>> nLon = 50 >>> nCoeff = 4 >>> nWavelengths = 20 >>> weights = np.array(nWavelengths * [0.2, 0.3, 0.3, 0.2]).reshape(nWavelengths, 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*nWavelengths).reshape((nWavelengths, nLevel, ... nLat, nLon, nCoeff)) >>> absor = np.linspace(1e-5, 1e-1, ... nLevel*nLat*nLon*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, ... nLat, nLon, nCoeff)) >>> asymm = np.linspace(0, 1, ... nLevel*nLat*nLon*nCoeff*nWavelengths).reshape((nWavelengths, nLevel, ... nLat, nLon, nCoeff)) >>> wavelengths = np.linspace(2e7, 9e7, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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, ... "nWavelengths": nWavelengths, ... "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 ... })
- makeFromLMDZ(LMDZinput, LMDZouput, weights=None, keys={'altitude': 'altitude', 'assym': 'gv', 'extinction': 'kv', 'geopotential': 'pphi', 'geopotentialSurf': 'pphis', 'latitude': 'lat', 'longitude': 'lon', 'pressure': 'p', 'ssa': 'wv', 'temp': 'temp', 'tsurf': 'tsurf', 'wavelength': 'wavelength_vi'}, wavelength=None, time=-1, hg=False, nAngle=181, phaseFunc=<function Data.<lambda>>)[source]
Generate an heterogeneous sphere from LMDZ output and input files.
- Parameters:
LMDZinput (str) – Path to the netCDF file containing surface information, to be read.
LMDZouput (str) – Path to the netCDF file containing atmosphere information, to be read.
weights (
numpy.ndarray, optional) – Array containing the weights of the Gaussian quadrature points for correlated-k data (shape=(nWeight)). If not provided, assumes no correlated-k are used. It is not required if wavelengths are separted (see wavelengths).keys ({str : str}, default htrdrPy.keysLMDZtitan_vi) –
Dictionnary of correspondance between keys used in this method and keys from the input/output file. It must contain the following keys:
- ’tsurf’default ‘tsurf’
Key for the surface temperature.
- ’temp’default ‘temp’
Key for the atmospheric temperature.
- ’wavelength’default ‘wavelength_vi’
Key for the wavelengths.
- ’ssa’: default ‘wv’
Key for the single scatering albedo.
- ’extinction’default ‘kv’
Key for the extinction coefficeint.
- ’assym’default ‘gv’
Key for the assymetry parameter.
- ’geopotential’default ‘pphi’
Key for the geopotential.
- ’altitude’default ‘altitude’
Key for the altitude.
- ’latitude’default ‘lat’
Key for the latitude.
- ’longitude’default ‘lon’
Key for the longitude.
- ’pressure’default ‘p’
Key for the pressure.
wavelengths ({str : dict}, optional) –
Dictionnary containing the wavelength bands data. It must be composed of 1 sub-dictionnary for every band.
- ’index’dict
Dictionnary containing the data for a specific band (‘index’ is the index refering to the band in the GCM output file, e.g. ‘v_23’ for the last visible band in Titan GCM). The sub-dictionnary must contain the following items:
- ’wavelength’float
Central wavelength of the band (in m).
- ’low’float
Lower bound of the band (in m).
- ’up’float
Upper bound of the band (in m).
- ’weights’
numpy.ndarray, optional The weights associated to the different k-correlated coefficents (shape=(nWeights)). If no k-coeff are used, omit the key.
- ’weights’
time (int, default -1) – Time index to use, default is last.
hg (bool, default False) – Whether or not to use the Henyey-Greenstein phase function built in htrdr
nAngle (int, optional, default 181) – Number of angles to use in the discrete phase functions. (If
hg=TruenAngle is omitted).phaseFunc (func, optional, default = 1 + g cos(theta)) – Function to calculate the discrete phase function. The function must take in first argument the asymetry parameter and in second argument the angle (in rad). If
hg=TruephaseFunc is omitted.
Notes
If the surface temperature is not provided in the output file, it will be read from the input file.
- 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)) .
- ”albedo”
- ”wavelengths”
numpy.ndarray, optional Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.
- ”wavelengths”
- ”bands”
numpy.ndarray, optional Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.
- ”bands”
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, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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(nWavelengths) * 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)) .
- ”albedo”
- ”wavelengths”
numpy.ndarray, optional Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.
- ”wavelengths”
- ”bands”
numpy.ndarray, optional Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.
- ”bands”
Notes
In plan-parallel mode, the
radiusparameter 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, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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(nWavelengths) * 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)).
- ”albedo”
- ”latitude”
numpy.ndarray List of latitudes (shape=(nLat), [°]).
- ”latitude”
- ”wavelengths”
numpy.ndarray, optional Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.
- ”wavelengths”
- ”bands”
numpy.ndarray, optional Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.
- ”bands”
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
latitudeprovided.Examples
>>> d = htrdrPy.Data(radius=1e6, nPhi=50, name="Planet") >>> nLat = 30 >>> wavelengths = np.linspace(2e7, 9e7, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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((nWavelengths, 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)).
- ”albedo”
- ”latitude”
numpy.ndarray List of latitudes (shape=(nLat), [°]).
- ”latitude”
- ”longitude”
numpy.ndarray List of longitudes (shape=(nLon), [°]).
- ”longitude”
- ”wavelengths”
numpy.ndarray, optional Wavelengths where the albedo is defined (shape=(nWavelength), [m]). Alternatively, the user can specify the bands with the “bands” keyword.
- ”wavelengths”
- ”bands”
numpy.ndarray, optional Wavelengths bands where the albedo is defined (shape=(nWavelength,2), [m]). The values corresponds to the bands limits.
- ”bands”
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
latitudeandlongitude, respectively.Examples
>>> d = htrdrPy.Data(radius=1e6, name="Planet") >>> nLat = 30 >>> nLon = 50 >>> wavelengths = np.linspace(2e7, 9e7, nWavelengths) >>> bandsLow = np.zeros(nWavelengths) >>> bandsUp = np.zeros(nWavelengths) >>> 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((nWavelengths, 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...
- 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
octreeFileis 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.