pyPamtra

pyPamtra Class

class pyPamtra.core.pyPamtra[source]

Bases: object

Class for pamtra calculations. Initalisation fills dictonary ‘set’, ‘nmlSet’, ‘dimensions’, and ‘units’ with default values.

addCloudShape()[source]

adds cloud base and cloud top to the data to an existing pamtra profile clouds are detected with rh >= 0.95 and T >= 253.15 K if multiple cloud layers are found, only the uppermost and the lowest one are saved.

addIntegratedValues()[source]
addPIA(direction, applyAtmo=True, applyHydros=True)[source]

adds two-way path integrated attenuation to result dictionary

Parameters:
  • direction ({'bottom-up', 'top-down'}) – Assumed direction

  • applyAtmo (bool, optional) – use attenuation of the atmosphere (default: true)

  • applyHydros (bool, optional) – use attenuation of the hydrometeors (default: true)

Notes

Adds “PIA” to self.r

addProfile(profile, axis=0)[source]

Add additional dictionary “profiles” with profile information to axis “axis”. Number of height bins must be equal.

Parameters:
  • profile (Pamtra profile dict) – Pamtra profile to append

  • axis (int) –

addPseudoAdiabaticLWC()[source]

adds liquid water content to the data to an existing pamtra profile using a simple cloud model with pseudo adiabatic lapse rate clouds are detected with rh >= 0.95 and T >= 253.15 K existing cwc_q values are overriden!

addSpectralBroadening(EDR, wind_uv, beamwidth_deg, integration_time, frequency, kolmogorov=0.5)[source]

Fills array self.p[“airturb”] according to input data. Array are broadcasted to self._shape3D Only broadening by horizontal wind and turbuilence is considered as of now. Populates self.p[“airturb”]

EDRarrayE

Eddy dissipation rate (SI units)

wind_uvarray

horizontal wind field (SI units)

beamwidth_degfloat

full-width half-radiated one-way (deg)

integration_timefloat

radar integration time in s

frequencyfloat

frequency used to determien lower integration bound

kolmogorovfloat,optional

Kolmogorov constant, default 0.5

averageResultTbs(translatorDict)[source]

average several frequencies of the passive observations to account for channel width. Replaces self.r[“tb”]. As of know this works only in passive mode. Backups of teh original data are kept in self.r[“not_averaged_tb”] and self.set[“not_averaged_freqs”]

Parameters:

translatorDict (dict) –

with list of old and new frequencies. Note that the “new” frequency is only for naming of the channel. e.g.
translatorDict = {

90.0: [90.0], 120.15: [117.35, 120.15],

}

createFullProfile(timestamp, lat, lon, wind10u, wind10v, obs_height, hgt_lev, press_lev, temp_lev, relhum_lev, hydro_q, hydro_n, hydro_reff, radar_prop, sfc_type, sfc_model, sfc_refl, sfc_salinity, sfc_slf, sfc_sif)[source]

create complete Pamtra Profile

No Extras, no missing values are guessed. the Data is only reshaped

createProfile(**kwargs)[source]

Function to create Pamtra Profiles.

Variables ending on _lev mean level values, variables without are layer values (vectors one entry shorter!)!

Everything is needed in SI units, relhum is in %

The following variables are mandatory: hgt_lev, (temp_lev or temp), (press_lev or press) and (relhum_lev OR relhum)

The following variables are optional and guessed if not provided: “timestamp”,”lat”,”lon”,”wind10u”,”wind10v”,”hgt_lev”,”hydro_q”,”hydro_n”,”hydro_reff”,”obs_height”,”sfc_type”,”sfc_model”,”sfc_refl”,”sfc_salinity”

hydro_q, hydro_reff and hydro_n can also provided as hydro_q+no001, hydro_q+no002 etc etc

#

dumpForRT4(runFile, path, prefix='sca_')[source]

Write a files to run the passive model with the pure RT4 code by Frank Evans. runFile contains filenames of the files in path

This method takes profile (0,0). It ignores everything beyond that.

filterProfiles(condition)[source]

discard profiles, which do not fullfill “condition”

Parameter

conditionarray of bool

2D boolean array

Note: If the initial field is a 2D grid, the field is flattend after application

Applicate between CreateProfile and runPamtra

loadResultsFromNumpy(fname)[source]

load complete pamtra object (profile,results,settings from (a) file(s)

Parameters:

fname (str) – filename or directory name

nmlSet

settings which are required for the nml file. keeping the order is important for fortran

p

test

readClassicPamtraProfile(inputFile, n_moments=1)[source]

read classical pamtra profile from file Input:

inputFile: str filename with path

readNmlFile(inputFile)[source]

read classical Pamtra Namelist File from inputFile

Parameter

nmlFile: str

filename with path

readPamtraProfile(inputFile)[source]

read lay or lev pamtra profile from file. Descriptot file must be defined before

Parameter

inputFile: str

filename with path

rescaleHeights(new_hgt_lev, new_hgt=[])[source]

Rescale Pamtra profile to new height grid

Parameters:
  • new_hgt_lev (Array) – Array to replace hgt_lev

  • new_hgt (Array, optional) – array to replace hgt_lev (default [], i.e. interpoated from new_hgt_lev)

runPamtra(freqs, checkData=True, failOnError=True)[source]

run Pamtra from python. Populates result dictionary ‘r’.

Parameters:
  • freqs (float of list of floats) – Frequencies in GHz.They must be sorted from small to large.

  • checkData (bool, optional) – Check input data for consitency (default True)

  • failOnError – raise Python Error in case of Fortran Error

runParallelPamtra(freqs, pp_local_workers='auto', pp_deltaF=1, pp_deltaX=0, pp_deltaY=0, checkData=True, timeout=None)[source]

run Pamtra parallel from python. Populates result dictionary ‘r’.

Parameters:
  • freqs (float of list of floats) – Frequencies in GHz.They must be sorted from small to large.

  • checkData (bool, optional) – Check input data for consitency (default True)

  • pp_local_workers (int or 'auto', optional) – number of parallel threads (default ‘auto’ correpons to number of cpus.)

  • pp_deltaF (int , optional) – Size of each thread in frequency domain. 0 means infinity. (default 1)

  • pp_deltaX (int , optional) – Size of each thread in X domain. 0 means infinity. (default 0)

  • pp_deltaY (int , optional) – Size of each thread in Y domain. 0 means infinity. (default 0)

  • timeout (int or None, optional) – Timeout for each thread in seconds (default None)

runPicklePamtra(freqs, picklePath='pyPamJobs', pp_deltaF=1, pp_deltaX=0, pp_deltaY=0, checkData=True, timeout=None, maxWait=3600)[source]

Special variant of runParallelPamtra writing Pickles to picklePath which are processed by another job.

runPicklePamtraSFTP(freqs, host, user, localPicklePath='pyPamJobs', remotePicklePath='pyPamJobs', pp_deltaF=1, pp_deltaX=0, pp_deltaY=0, checkData=True, timeout=None, maxWait=3600)[source]

Special variant of runParallelPamtra writing Pickles to picklePath which are send by SFTP.

tileProfiles(rep2D)[source]

repeat the profiles of Pamtra

Parameters:

rep2D (tuple(int,int)) – E.g. rep2D=(1,100) changes a Pamtra shape from (10,2) to (10,200).

writeNmlFile(nmlFile)[source]

write classical Pamtra Namelist File from nmlSet

Parameter

nmlFile: str

filename with path

writePamtraProfile(profileFile)[source]

write lay or lev (depending on file extension) pamtra profile to file.

Parameter

profileFile: str

filename with path

writeResultsToNetCDF(fname, profileVars='all', wpNames=[], ncForm='NETCDF3_CLASSIC', xarrayCompatibleOutput=False, ncCompression=False, lfracCompatibility=False)[source]

write the results to a netcdf file

Parameters:
  • fname (str) – filename with path

  • profileVars (list of str or 'all', optional) – list of variables of the profile to be saved. “all” saves all implmented ones (default all)

  • ncForm (str, optional) – netcdf file format, possible values are NETCDF3_CLASSIC, NETCDF3_64BIT, NETCDF4_CLASSIC, and NETCDF4 for the python-netcdf4 package (netcdf3 gives netcdf4 files for newer ubuntu versions!!”) NETCDF3 takes the “old” Scientific.IO.NetCDF module, which is a bit more convinient (default NETCDF3_CLASSIC)

  • wpNames (list of str, optional) – integrated values to be saved (default [])

  • xarrayCompatibleOutput (bool, optional, default(False)) –

    if True, make passive-only nc file xarray dataset compatible:

    set passive_polarisation to ‘V’ and ‘H’ rename dimension ‘outlevels’ to ‘outlevel’ and provide outlevel index vector

    Otherwise passive_polarisation is filled with ‘VH’ and outlevel dimension is ‘outlevels’ (default)

  • ncCompression (bool, optional) – If True, use netCDF4 compression. (NETCDF4 and NETCDF4_CLASSIC only). Otherwise save with default settings (default, no compression)

  • lfracCompatibility (bool, optional) – If True create ‘lfrac’ array filled with sfc_type. lfrac is deprecated since commit git:601f739a1c8419f84231ac8d5d4d2534e361c014 If False, do not create lfrac variable in netcdf (default)

writeResultsToNumpy(fname, seperateFiles=False)[source]

write the complete state of the session (profile,results,settings to npy pickles.

Parameters:
  • fname (str) – filename or - if seperateFiles - directory name

  • seperateFiles (bool, optional) – Write every variable to separate file. Required for very large data sets.

descriptorFile Class

class pyPamtra.descriptorFile.pamDescriptorFile(parent)[source]

Bases: object

class with the descriptor file content in data. In case you want to use 4D data, use data4D instead, the coreesponding column in data is automatically removed. Inializes with an empty array.

add4D(key, arr)[source]

Add 4D information for one hydrometeor property. I.e the 1D value with dimesnions (hydrometeors) in the df array is replaced by an array with (x,y,z,hydrometeor) dimensions.

Parameters:
  • key (str) – Name of hydrometeor property to be altered.

  • arr (array) – New 4D values

addFullSpectra()[source]

This function creates empty arrays which have 5 dimensions (x,y,z,hydrometeor,bin). They are not returned but added to the df object. Note that the hydrometeors have to be added to the df array before calling this routine with dummy values. Information not avaialable here (e.g. scattering) are still taken from the descriptor file. Note that all arrays of dataFullSpec mus be filled and you still need to define the hydrometeors if the df object! However, only name, liq_ice, nbins+1, scat_name, and vel_size_mod is used even though you have to define the other variables as well.

Make sure that the array dimensions are not changed, otherwise Pamtra might crash.

The following objects are added to the df object: :returns: * dataFullSpec[“rho_ds”] – particle density at the middle of the bin in [kg/m3]

  • dataFullSpec[“d_ds”] – particle diameter at the middle of the bin where the scattering properties shall be calculated

  • dataFullSpec[“d_bound_ds”] – particle diameter at the edges of the bins. 5th dimension bin is one longer.

  • dataFullSpec[“n_ds”] – particle number concentration in [#/m^3]. I.e. it is NOT normalized!

  • dataFullSpec[“mass_ds”] – particle mass at the middle of the bin in [kg]

  • dataFullSpec[“area_ds”] – particle cross section area at the middle of the bin in [m^2]

  • dataFullSpec[“as_ratio”] – particle aspect ratio in the middle of the bin [no unit]

  • dataFullSpec[“canting”] – particle canting the middle of the bin in [deg]

  • dataFullSpec[“fallvelocity”] – particle fall velocity in [m/s]. Only used if any value > 0 m/s, otherwise the fall velocity relation provided in the descriptor file is used.

Notes

Do not forget to set nmlSet[‘hydro_fullspec’] to True.

addHydrometeor(hydroTuple)[source]

Add hydrometeor to df array. See Descriptor File Structure for details.

Parameters:

hydroTuple (tuple) – tuple describing the hydrometeor. Consists of “hydro_name”, “as_ratio”, “liq_ice”, “rho_ms”, “a_ms”, “b_ms”, “alpha_as”, “beta_as”, “moment_in”, “nbin”, “dist_name”, “p_1”, “p_2”, “p_3”, “p_4”, “d_1”, “d_2”, “scat_name”, “vel_size_mod”,”canting”

readFile(fname)[source]

Read ASCII Pamtra descriptor file

Parameters:

fname (str) – filename

remove4D(key, val)[source]

Remove 4D information for one hydrometeor property.

Parameters:
  • key (str) – Name of hydrometeor property to be altered.

  • arr (array) – New 1 D values.

removeHydrometeor(hydroName)[source]

Remove hydrometeor from df array.

Parameters:

hydroName (str) – Name of hydrometeor to be removed.

writeFile(fname)[source]

Write ASCII Pamtra descriptor file

Parameters:

fname (str) – filename

pyPamtra.descriptorFile.riming_dependent_area_size(M, monomer)[source]

Return area size parameter a_A and b_A for given monomer type from the normalized rime mass M as defined in Seifert et al. (2019), M = mass of rime / mass of size equivalent graupel assuming graupel density of 700kg/m³ a_A and b_A are interpolated from Table 3 in Maherndl et al [in prep.], for M>0.8155 values for M=0.8155 are used Input: M (unitless) monomer type (column, dendrite, needle, plate, rosette, mean)

Output a_A, b_A

pyPamtra.descriptorFile.riming_dependent_mass_size(M, monomer)[source]

Return mass size parameter a_m and b_m for given monomer type from the normalized rime mass M as defined in Seifert et al. (2019), M = mass of rime / mass of size equivalent graupel assuming graupel density of 700kg/m³ a_m and b_m are interpolated from Table 2 in Maherndl et al [in prep.], for M>0.8155 values for M=0.8155 are used Input: M (unitless) monomer type (column, dendrite, needle, plate, rosette, mean)

Output a_m, b_m

pyPamtra.descriptorFile.riming_dependent_ssrga_name(M)[source]

Return string scattering name for ss-rayleigh-gans in form “ss-rayleigh-gans_kappa_beta_gamma_zeta” Calculate ssrga parameter kappa, beta, gamma, zeta1 and alpha_eff from the normalized rime mass M as defined in Seifert et al. (2019), M = mass of rime / mass of size equivalent graupel assuming graupel density of 700kg/m³

Input: M (unitless)

Output scatering name string

pyPamtra.descriptorFile.ssrga_parameter(M)[source]

Calculate ssrga parameter kappa, beta, gamma, zeta1 and alpha_eff from the normalized rime mass M as defined in Seifert et al. (2019), M = mass of rime / mass of size equivalent graupel assuming graupel density of 700kg/m³

Input: M (unitless)

Output kappa, beta, gamma, zeta1, alpha_eff

Import Routines

Collection of importer modules for the different input types.

They all return a pyPamtra object that holds the profile information in the .p object.

In addition ncToDict is provide to easily read necdf files into a dictionary.

pyPamtra.importer.createUsStandardProfile(pam=<pyPamtra.core.pyPamtra object>, **kwargs)[source]

Function to create clear sky US Standard Atmosphere.

hgt_lev is the only mandatory variables humidity will be set to zero if not provided, all other variables are guessed by “createProfile”

values provided in kwargs will be passed to “createProfile”, however, press_lev and temp_lev will overwrite us standard if provided

pyPamtra.importer.ncToDict(ncFilePath, keys='all', joinDimension='time', offsetKeys={}, ncLib='netCDF4', tmpDir='/tmp/', skipFiles=[])[source]

Load keys of netcdf file into dictionary that is returned. By using wildcards in ncFilePath, multiple files are possible. They are joint using joinDimension. Offsets of e.g.,, time vector can be corrected by offsetKeys with value to be corrected as key as correction key in value. gzip compressed netcdf with extension .gz are possible.

pyPamtra.importer.readCosmoDe1MomDataset(fnames, kind, descriptorFile, forecastIndex=1, colIndex=0, tmpDir='/tmp/', fnameInTar='', concatenateAxis=1, debug=False, verbosity=0, df_kind='default', constantFields=None, maxLevel=0, subGrid=None)[source]

read COMSO one moment data set

Parameters:
  • fnames ({str}) – file name, wildCards allowed! can be either nc, nc and gz, or nc.gz files of name fnameInTar within tar file

  • kind ({str}) – kind of COSMO file, right now only collum and synsat netcdf files are implemented

  • descriptorFile ({str}) – path and name of descriptor file

  • forecastIndex ({int}, optional) – index of forecast to be read (the default is 1)

  • colIndex ({int or list of int}, optional) – which collum should be taken? list allowed! (the default is 0)

  • tmpDir ({str}, optional) – [description] (the default is “/tmp/”)

  • fnameInTar ({str}, optional) – [description] (the default is “”)

  • concatenateAxis ({int}, optional) – concatenation along which axis (the default is 1)

  • debug ({bool}, optional) – stop and load debugger on exception (the default is False)

  • verbosity ({int}, optional) – verbosity of the module (the default is 0) (the default is 0)

  • df_kind ({str}, optional) – kind of data set. This defines the name of variables and dimensions. (the default is “default”)

  • constantFields ({str}, optional) – path and name of file with constant fields (the default is None)

  • maxLevel ({int}, optional) – Number of maximum levels to consider. (the default is 0)

  • subGrid ({[int,int,int,int]}, optional) – array with indices [lon_start,lon_end,lat_start,lat_end] ((1,1) in model corresponds to (0,0) in python!) (the default is None)

Returns:

pam

Return type:

pyPamtra Object

Raises:
  • IOError

  • TypeError

pyPamtra.importer.readCosmoDe2MomDataset(fnamesA, descriptorFile, fnamesN=None, kind='new', forecastIndex=1, tmpDir='/tmp/', fnameInTar='', debug=False, verbosity=0, df_kind='default', constantFields=None, maxLevel=0, subGrid=None)[source]

import COSMO 2-moment dataset

Parameters:
  • fnamesA ({str}) – path and names of files with atmospheric variables. wildCards allowed! can be either nc file, nc,gz file or nc.gz file of name fnameInTar within tar file

  • descriptorFile ({str}) – path and name of descriptor file

  • fnamesN ({str}, optional) – path and names of files with total numbers. wildCards allowed! can be either nc file, nc,gz file or nc.gz file of name fnameInTar within tar file (the default is None)

  • kind ({str}, optional) – can be either new (default, new NARVAL calculations) or old (old NARVAL calculations) (the default is ‘new’)

  • forecastIndex ({int}, optional) – index of forecast to be read (the default is 1)

  • tmpDir ({str}, optional) – [description] (the default is “/tmp/”)

  • fnameInTar ({str}, optional) – [description] (the default is “”)

  • debug ({bool}, optional) – stop and load debugger on exception (the default is False)

  • verbosity ({int}, optional) – verbosity of the module (the default is 0) (the default is 0)

  • df_kind ({str}, optional) – kind of data set. This defines the name of variables and dimensions. (the default is “default”)

  • constantFields ({str}, optional) – path and name of file with constant fields (the default is None)

  • maxLevel ({int}, optional) – Number of maximum levels to consider. (the default is 0)

  • subGrid ({[int,int,int,int]}, optional) – array with indices [lon_start,lon_end,lat_start,lat_end] ((1,1) in model corresponds to (0,0) in python!) (the default is None)

Returns:

pam

Return type:

pyPamtra Object

Raises:

IOError – on file error

pyPamtra.importer.readCosmoDe2MomDatasetOnFlightTrack(fnameA, descriptorFile, tmpDir='/tmp/', debug=False, verbosity=0, maxLevel=0)[source]

import COSMO 2-moment dataset extracted on a HALO flight track

Parameters:
  • fnameA ({str}) – path and names of file with atmospheric variables, wildCards allowed! can be either nc file, nc,gz file or nc.gz file of name fnameInTar within tar file

  • descriptorFile ({str}) – path and name of descriptor file

  • tmpDir ({str}, optional) – directory to unpack temporary files (the default is “/tmp/”)

  • debug ({bool}, optional) – stop and load debugger on exception (the default is False)

  • verbosity ({int}, optional) – verbosity of the module (the default is 0)

  • maxLevel ({int}, optional) – Number of maximum levels to consider. (the default is 0)

Returns:

pam

Return type:

pyPamtra Object

pyPamtra.importer.readCosmoReAn2km(constantFields, fname, descriptorFile, forecastIndex=1, tmpDir='/tmp/', fnameInTar='', debug=False, verbosity=0, maxLevel=51, subGrid=None)[source]

reads COSMO 2 km Re-Analyses

Parameters:
  • constantFields ({str}) – path and file name of constant file

  • fname ({str}) – path and file name of data file

  • descriptorFile ({str}) – path and name of descriptor file

  • forecastIndex ({int}, optional) – index of forecast to be read (the default is 1)

  • tmpDir ({str}, optional) – temporary directory to unpack data (the default is “/tmp/”, which [default_description])

  • fnameInTar ({str}, optional) – file name within tar file (the default is “”)

  • debug ({bool}, optional) – switch on debugging mode(the default is False)

  • verbosity ({int}, optional) – [description] (the default is 0)

  • maxLevel ({int}, optional) – index of maximum level (the default is 51)

  • subGrid ({[int,int,int,int]}, optional) – read subgrid of whole field (the default is None)

Return type:

pyPamtra object

pyPamtra.importer.readCosmoReAn6km(constantFields, fname, descriptorFile, forecastIndex=1, tmpDir='/tmp/', fnameInTar='', debug=False, verbosity=0, df_kind='default', maxLevel=41, subGrid=None)[source]
pyPamtra.importer.readECMWF(fname, constantFile, descriptorFile, landseamask, debug=False, verbosity=0, grid=[0, -1, 0, -1])[source]

read ECMWF IFS data

reads the data produced for the NAWDEX campaign by Heini Wernli at ETH Zuerich

Q specific humidity kg/kg LWC specific liquid water content kg/kg IWC specific ice water content kg/kg RWC specific rain water content kg/kg SWC specific snow water content kg/kg T temperature K - U u wind component m/s V v wind component m/s PS surface pressure Pa SLP mean sea level pressure Pa P pressure Pa SKT skin temperature K u_wind_10m 10 meter wind u component m/s v_wind_10m 10 meter wind v component m/s lsm land-sea-mask 0-1 (0=sea, 1=land)

Parameters:
  • fname ({str}) – path and name of data file

  • constantFile ({str}) – path and name of constant fields

  • descriptorFile ({str}) – path and name of descriptor file

  • landseamask ({str}) – path and name of landseamask file

  • debug ({bool}, optional) – stop and load debugger on exception (the default is False)

  • verbosity ({int}, optional) – verbosity of the module (the default is 0) (the default is 0)

  • grid ({[int,int,int,int]}, optional) – read only subgrid of whole model field

Returns:

pam

Return type:

pyPamtra Object

pyPamtra.importer.readHIRHAM(dataFile, additionalFile, topoFile, descriptorFile, grid=[0, 200, 0, 218], timestep=0, debug=False, verbosity=0)[source]

read HIRHAM output into PAMTRA

it reads the output produce by Ines Haberstedt at the AWI in Potsdam

Parameters:
  • dataFile ({str}) – path and name of data file

  • additionalFile ({str}) – path and name of additional file

  • topoFile ({str}) – path and name of file with topography information

  • descriptorFile ({str}) – path and name of descriptor file

  • grid ({[int,int,int,int]}, optional) – subgrid to read (the default is [0,200,0,218] the full grid)

  • timestep ({int}, optional) – whcih if the 8 available time steps should be read (the default is 0, max is 7)

  • debug ({bool}, optional) – stop and load debugger on exception (the default is False)

  • verbosity ({int}, optional) – verbosity of the module (the default is 0) (the default is 0)

Returns:

pam

Return type:

pyPamtra Object

pyPamtra.importer.readICON1momNWP(fname, descriptorFile, debug=False, verbosity=0, grid=[0, 463, 0, 6653])[source]

import ICON LEM 1-moment dataset output cross section over a site (Meteogram)

WARNING: This importer has been designed upon ICON 1 mom NWP output generated by Dr. Helene Bresson. Will most likely not work on other files.

fname = filename of the ICON output descriptorFile = pyPamtra descriptorFile object or filename of a proper formatted descriptorFile debug = flag causing stop and load of debugger upon raised exception verbosity = pyPamtra.pyVerbose verbosity level

22/03/2020 - Mario Mech - mario.mech@uni-koeln.de

pyPamtra.importer.readICON2mom(fname, descriptorFile, fextpar=None, finit=None, forcing='ICON', time=0, kind=None, debug=False, verbosity=0)[source]

import ICON LEM 2-moment dataset output in it’s pure version

fname = filename of the Icon output descriptorFile = pyPamtra descriptorFile object or filename of a proper formatted descriptorFile processed = flag to specify whether the files have been processed or not. (default = True) debug = flag causing stop and load of debugger upon raised exception verbosity = pyPamtra.pyVerbose verbosity level

Pamtra assumes that the apart from height, the first given dimension is latitude, here coordinates do not change across the dimension, but time does, thus generating a meteogram

09/11/2020 - Mario Mech - mario.mech@uni-koeln.de based on readIcon2momOnFlightTrack

pyPamtra.importer.readIcon1momMeteogram(fname, descriptorFile, debug=False, verbosity=0, timeidx=None, hydro_content=[1.0, 1.0, 1.0, 1.0, 1.0])[source]

import ICON LEM 1-moment dataset output cross section over a site (Meteogram)

WARNING: This importer has been designed upon Icon meteogram output generated by Dr. Vera Schemann Might not work on other files.

fname = filename of the Icon output descriptorFile = pyPamtra descriptorFile object or filename of a proper formatted descriptorFile debug = flag causing stop and load of debugger upon raised exception verbosity = pyPamtra.pyVerbose verbosity level timeidx = indexes of timestep to be computed. Used to slice just a few profiles out of the .nc file for rapid testing

Pamtra assumes that the apart from height, the first given dimension is latitude, here coordinates do not change across the dimension, but time does, thus generating a meteogram

03/03/2018 - Davide Ori - dori@uni-koeln.de

pyPamtra.importer.readIcon2momMeteogram(fname, descriptorFile, debug=False, verbosity=0, timeidx=None, hydro_content=[1.0, 1.0, 1.0, 1.0, 1.0, 1.0])[source]

import ICON LEM 2-moment dataset output cross section over a site (Meteogram)

WARNING: This importer has been designed upon Icon meteogram output generated by Dr. Vera Schemann Might not work on other files.

fname = filename of the Icon output descriptorFile = pyPamtra descriptorFile object or filename of a proper formatted descriptorFile debug = flag causing stop and load of debugger upon raised exception verbosity = pyPamtra.pyVerbose verbosity level timeidx = indexes of timestep to be computed. Used to slice just a few profiles out of the .nc file for rapid testing

Pamtra assumes that the apart from height, the first given dimension is latitude, here coordinates do not change across the dimension, but time does, thus generating a meteogram

03/03/2018 - Davide Ori - dori@uni-koeln.de

pyPamtra.importer.readIcon2momOnFlightTrack(fname, descriptorFile, time=0, kind='processed', constant_file=None, debug=False, verbosity=0)[source]

import ICON LEM 2-moment dataset output along a flight track

WARNING: This importer has been originally designed upon Icon output processed by Dr. Vera Schemann It has been later on extended to read extracted flight tracks without further processing. Might not work on other files.

fname = filename of the Icon output descriptorFile = pyPamtra descriptorFile object or filename of a proper formatted descriptorFile processed = flag to specify whether the files have been processed or not. (default = True) debug = flag causing stop and load of debugger upon raised exception verbosity = pyPamtra.pyVerbose verbosity level

Pamtra assumes that the apart from height, the first given dimension is latitude, here coordinates do not change across the dimension, but time does, thus generating a meteogram

16/05/2018 - Mario Mech - mario.mech@uni-koeln.de based on readIcon2momMeteogram by Davide Ori

pyPamtra.importer.readIconNwp1MomDataset(fname_fg, descriptorFile, debug=False, verbosity=0, constantFields=None, maxLevel=0, diagnostic=False)[source]

import ICON SRM 1-moment dataset

It is Icon with cosmo physics

fname_fg{str}

path and name of file with atmospheric variables.

descriptorFile{str}

path and name of descriptor file

debug{bool}, optional

stop and load debugger on exception (the default is False)

verbosity{int}, optional

verbosity of the module (the default is 0) (the default is 0)

constantFields{str}, optional

path and name of file with constant fields (the default is None)

maxLevel{int}, optional

Number of maximum levels to consider. (the default is 0)

diagnostic{bool}

Switch to use diagnostic qc and qi instead of prognostic ones is set to True (the default is False)

pam : pyPamtra Object

IOError

Remarks

This importer is adjusted to read regridded ICON-SRM simulations as they where run by Daniel Klocke within the HErZ-NARVALII framework. The Output consists of a _fg_ and a _cloud_ netcdf file. From the _fg_ file, which contains most atmospheric variables, we need data on grid 2. From the _cloud_ file, which contains qg, we need data on grid 1. Also there is the extpar_narval_nestTropAtl_R1250m.nc file wich contains time constant parameter.

Regridding

Use cdo to regrid the unstructured ICON data on a latlon grid. Nearest neighbor interpolation is used for regridding. The interpolation is defined in a Grid file. Example (filename “test_grid”):

` gridtype   = lonlat gridsize   = 8000 xname      = lon xlongname  = longitude xunits     = degrees_east yname      = lat ylongname  = latitude yunits     = degree_north xsize      = 100 ysize      = 80 xfirst     = -59.0 xinc       = 0.05 yfirst     = 10.0 yinc       = 0.05 `

Convert unstructured data to latlon in one cdo command:

cdo remapnn,test_grid -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc /dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016_test_grid.nc

Here, “narval_nestTropAtl_R1250m.nc” is the description of the unstructured data. Current location: /data/hamp/models/ICON/HErZ-NARVALII/GRIDS/narval_nestTropAtl_R1250m.nc

The corresponding command for _cloud_ data is:

cdo remapnn,test_grid -selgrid,1 -setgrid,narval_nestTropAtl_R1250m.nc /dei4_NARVALII_2016081700_cloud_DOM02_ML_0016.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016_test_grid.nc

Optimized regridding

To speedup the regridding of multiple ICON output one can save the interpolation weights in an intermediate netcdf file. All ICON output must be given on the same unstructured grid.

Generate the weights:

cdo gennn,test_grid -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc test_grid_weights.nc

The same weights can be used for _fg_ and _cloud_ (unconfirmed). Apply weights:

cdo remap,test_grid,test_grid_weights.nc -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016_test_grid.nc cdo remap,test_grid,test_grid_weights.nc -selgrid,1 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016_test_grid.nc cdo remapnn,test_grid,test_grid_weights.nc -setgrid,narval_nestTropAtl_R1250m extpar_narval_nestTropAtl_R1250m.nc extpar_narval_nestTropAtl_R1250m_test_grid.nc

pyPamtra.importer.readIconNwp1MomDataset_cells(fname_fg, descriptorFile, debug=False, verbosity=0, constantFields=None, maxLevel=0)[source]

import ICON SRM 1-moment dataset where certain cell were sellected

Such data gan be generated using $(cdo -selgridcell) It is ICON with cosmo physics

Use PAMTRAs internal “lat” index for cell. “lon”-dimension is 1.

fname_fg{str}

path and name of file with atmospheric variables. fname_fg has to include “_fg_”. A corresponding “_cloud_” has do exist. No wildCards allowed. Must be nc file.

descriptorFile{str}

path and name of descriptor file

debug{bool}, optional

stop and load debugger on exception (the default is False)

verbosity{int}, optional

verbosity of the module (the default is 0) (the default is 0)

constantFields{str}, optional

path and name of file with constant fields (the default is None)

maxLevel{int}, optional

Number of maximum levels to consider. (the default is 0)

pam : pyPamtra Object

IOError

Remarks

This importer is adjusted to read regridded ICON-SRM simulations as they where run by Daniel Klocke within the HErZ-NARVALII framework. The Output consists of a _fg_ and a _cloud_ netcdf file. From the _fg_ file, which contains most atmospheric variables, we need data on grid 2. From the _cloud_ file, which contains qg, we need data on grid 1. Also there is the extpar_narval_nestTropAtl_R1250m.nc file wich contains time constant parameter.

Regridding

Use cdo to regrid the unstructured ICON data on a latlon grid. Nearest neighbor interpolation is used for regridding. The interpolation is defined in a Grid file. Example (filename “test_grid”):

` gridtype   = lonlat gridsize   = 8000 xname      = lon xlongname  = longitude xunits     = degrees_east yname      = lat ylongname  = latitude yunits     = degree_north xsize      = 100 ysize      = 80 xfirst     = -59.0 xinc       = 0.05 yfirst     = 10.0 yinc       = 0.05 `

Convert unstructured data to latlon in one cdo command:

cdo remapnn,test_grid -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc /dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016_test_grid.nc

Here, “narval_nestTropAtl_R1250m.nc” is the description of the unstructured data. Current location: /data/hamp/models/ICON/HErZ-NARVALII/GRIDS/narval_nestTropAtl_R1250m.nc

The corresponding command for _cloud_ data is:

cdo remapnn,test_grid -selgrid,1 -setgrid,narval_nestTropAtl_R1250m.nc /dei4_NARVALII_2016081700_cloud_DOM02_ML_0016.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016_test_grid.nc

Optimized regridding

To speedup the regridding of multiple ICON output one can save the interpolation weights in an intermediate netcdf file. All ICON output must be given on the same unstructured grid.

Generate the weights:

cdo gennn,test_grid -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc test_grid_weights.nc

The same weights can be used for _fg_ and _cloud_ (unconfirmed). Apply weights:

cdo remap,test_grid,test_grid_weights.nc -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016_test_grid.nc cdo remap,test_grid,test_grid_weights.nc -selgrid,1 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016_test_grid.nc cdo remapnn,test_grid,test_grid_weights.nc -setgrid,narval_nestTropAtl_R1250m extpar_narval_nestTropAtl_R1250m.nc extpar_narval_nestTropAtl_R1250m_test_grid.nc

pyPamtra.importer.readIconNwp2MomDataset(fname_fg, descriptorFile, debug=False, verbosity=0, constantFields=None, maxLevel=0)[source]

import ICON SRM 2-moment dataset

It is ICON with cosmo physics

fname_fg{str}

path and name of file with atmospheric variables.

descriptorFile{str}

path and name of descriptor file

debug{bool}, optional

stop and load debugger on exception (the default is False)

verbosity{int}, optional

verbosity of the module (the default is 0) (the default is 0)

constantFields{str}, optional

path and name of file with constant fields (the default is None)

maxLevel{int}, optional

Number of maximum levels to consider. (the default is 0)

pam : pyPamtra Object

IOError

Remarks

This importer is adjusted to read regridded ICON-SRM simulations as they where run by Daniel Klocke within the HErZ-NARVALII framework. The Output consists of a _fg_ and a _cloud_ netcdf file. From the _fg_ file, which contains most atmospheric variables, we need data on grid 2. From the _cloud_ file, which contains qg, we need data on grid 1. Also there is the extpar_narval_nestTropAtl_R1250m.nc file wich contains time constant parameter.

Regridding

Use cdo to regrid the unstructured ICON data on a latlon grid. Nearest neighbor interpolation is used for regridding. The interpolation is defined in a Grid file. Example (filename “test_grid”):

` gridtype   = lonlat gridsize   = 8000 xname      = lon xlongname  = longitude xunits     = degrees_east yname      = lat ylongname  = latitude yunits     = degree_north xsize      = 100 ysize      = 80 xfirst     = -59.0 xinc       = 0.05 yfirst     = 10.0 yinc       = 0.05 `

Convert unstructured data to latlon in one cdo command:

cdo remapnn,test_grid -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc /dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016_test_grid.nc

Here, “narval_nestTropAtl_R1250m.nc” is the description of the unstructured data. Current location: /data/hamp/models/ICON/HErZ-NARVALII/GRIDS/narval_nestTropAtl_R1250m.nc

The corresponding command for _cloud_ data is:

cdo remapnn,test_grid -selgrid,1 -setgrid,narval_nestTropAtl_R1250m.nc /dei4_NARVALII_2016081700_cloud_DOM02_ML_0016.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016_test_grid.nc

Optimized regridding

To speedup the regridding of multiple ICON output one can save the interpolation weights in an intermediate netcdf file. All ICON output must be given on the same unstructured grid.

Generate the weights:

cdo gennn,test_grid -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc test_grid_weights.nc

The same weights can be used for _fg_ and _cloud_ (unconfirmed). Apply weights:

cdo remap,test_grid,test_grid_weights.nc -selgrid,2 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016.nc dei4_NARVALII_2016081700_fg_DOM02_ML_0016_test_grid.nc cdo remap,test_grid,test_grid_weights.nc -selgrid,1 -setgrid,narval_nestTropAtl_R1250m.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016.nc dei4_NARVALII_2016081700_cloud_DOM02_ML_0016_test_grid.nc cdo remapnn,test_grid,test_grid_weights.nc -setgrid,narval_nestTropAtl_R1250m extpar_narval_nestTropAtl_R1250m.nc extpar_narval_nestTropAtl_R1250m_test_grid.nc

pyPamtra.importer.readMesoNH(fnameBase, fnameExt, dataDir='.', debug=False, verbosity=0, dimX=160, dimY=160, dimZ=25, subGrid=None)[source]

read MesoNH output as produced by Jean-Pierre Chaboureau for the SIMGEO study

Parameters:
  • fnameBase ({str}) – base name of file name

  • fnameExt ({str}) – extension of file name

  • dataDir ({str}, optional) – path to data directory (the default is “.”)

  • debug ({bool}, optional) – stop and load debugger on exception (the default is False)

  • verbosity ({int}, optional) – verbosity of the module (the default is 0) (the default is 0)

  • dimX ({int}, optional) – dimension in x (the default is 160)

  • dimY ({int}, optional) – dimension in y (the default is 160)

  • dimZ ({int}, optional) – dimension in z (the default is 25)

  • subGrid ({[int,int,int,int]}, optional) – only read subgrid (the default is None)

Returns:

pam

Return type:

pyPamtra Object

pyPamtra.importer.readWrfDataset(fname, kind)[source]

read WRF data set

Parameters:
  • fname ({str}) – path and name of file atmospheric variables

  • kind ({str}) – kind of data file (up to now only gce is implemented)

Returns:

pam

Return type:

pyPamtra Object

Raises:

TypeError

Plot Routines

pyPamtra.plot.niceColors(length, cmap='hsv')[source]
pyPamtra.plot.plotTB(data, polarisation='H', angle=180.0, outlevel='space', frequency=0, xgrid='index', ygrid='index', levels=['min', 'max'], cmap=None, title='', fig=None)[source]

plot brightness temperatures from pamtra

input:

data: netcdf file OR pyPamtra object polarisatuion: H or V outlevel: space or ground frequency: if real number frequency, if int, index of frequency xgrid: index or? ygrid index or?

pyPamtra.plot.plotTBLine(data, polarisation='H', angle=180.0, outlevel='space', frequency=0, cmap=None, title='', fig=None, xIndex='all', yIndex=slice(0, 1, None), xVar=None, freqIndices='all', legendLoc='auto', relative=False)[source]

plot brightness temperatures from pamtra

input:

data: netcdf file OR pyPamtra object polarisatuion: H or V outlevel: space or ground frequency: if real number frequency, if int, index of frequency xgrid: index or? ygrid index or? relative: relative to 0,0?

MeteoSI Routines

pyPamtra.meteoSI.T_virt_q(T, q)[source]

Calculate the virtual temperature from air temperature and specific humidity.

Input: T is in K q is in kg/kg

Output: T_virt in K

pyPamtra.meteoSI.T_virt_rh(T, rh, p)[source]

Calculate the virtual temperature from air temperature, pressure, and relative humidity.

Input: T is in K rh is in Pa/Pa p is in Pa

Output: T_virt in K

pyPamtra.meteoSI.a2rh(a, T)[source]

Calculate the relative from absolute humidity and air temperature.

Input T is in K a in kg/kg Output rh in Pa/Pa

Source: Kraus: Chapter 8.1.2

pyPamtra.meteoSI.adiab(i, T, P, z)[source]

Adiabtic liquid water content assuming pseudoadiabatic lapse rate throughout the whole cloud layer. Thus the assumed temperature profile is differnt from the measured one Input: i no of levels T is in K p is in Pa z is in m Output: LWC in translated to Python from adiab.pro by mx.

pyPamtra.meteoSI.detect_liq_cloud(z, t, rh)[source]
pyPamtra.meteoSI.e2q(e, p)[source]

Calculate the specific humidity from water vapour pressure and air pressure.

Input: e is in Pa p is in Pa

Output q in kg/kg

pyPamtra.meteoSI.e_sat_gg_water(T)[source]

Calculates the saturation pressure over water after Goff and Gratch (1946). It is the most accurate that you can get for a temperture range from -90°C to +80°C. Source: Smithsonian Tables 1984, after Goff and Gratch 1946 http://cires.colorado.edu/~voemel/vp.html http://hurri.kean.edu/~yoh/calculations/satvap/satvap.html

Input: T in Kelvin. Output: e_sat_gg_water in Pa.

pyPamtra.meteoSI.mod_ad(T_cloud, p_cloud, z_cloud, fak)[source]
pyPamtra.meteoSI.moist_rho_q(p, T, q, *qm)[source]

Input p is in Pa T is in K q is in kg/kg Optional, several possible: qm is in kg/kg other species which contribute to the air mass! (ice, snow, cloud etc.)

Output: density of moist air [kg/m^3]

Example: moist_rho_q(p,T,q,q_ice,q_snow,q_rain,q_cloud,q_graupel,q_hail)

pyPamtra.meteoSI.moist_rho_rh(p, T, rh, *qm)[source]

Input: p is in Pa T is in K rh is in Pa/Pa Optional, several possible: qm is in kg/kg other species which contribute to the air mass! (ice, snow, cloud etc.)

Output: density of moist air [kg/m^3]

Example: moist_rho_rh(p,T,rh,q_ice,q_snow,q_rain,q_cloud,q_graupel,q_hail)

pyPamtra.meteoSI.pseudoAdiabLapseRate(T, Ws)[source]

Pseudoadiabatic lapse rate Input: T [K] thermodynamic temperature Ws [1] mixing ratio of saturation Output: PSEUDO [K/m] pseudoadiabatic lapse rate Constants: Grav [m/s2] : constant of acceleration

CP [J/(kg K)] : specific heat cap. at const. press Rair [J/(kg K)] : gas constant of dry air Rvapor [J/(kg K)] : gas constant of water vapor

Source: Rogers and Yau, 1989: A Short Course in Cloud Physics (III.Ed.), Pergamon Press, 293p. Page 32 translated to Python from pseudo1.pro by mx

pyPamtra.meteoSI.q2e(q, p)[source]

Calculate water vapour pressure from the specific humidity and air pressure.

Input: q in kg/kg p is in Pa

Output e is in Pa

pyPamtra.meteoSI.q2rh(q, T, p)[source]

Calculate relative humidity from specific humidity

Input: T is in K p is in Pa q in kg/kg

Output: rh is in Pa/Pa

pyPamtra.meteoSI.rh2a(rh, T)[source]

Calculate the absolute humidity from relative humidity, air temperature, and pressure.

Input T is in K rh is in Pa/Pa p is in Pa Output a in kg/m^3

Source: Kraus: Chapter 8.1.2

pyPamtra.meteoSI.rh2q(rh, T, p)[source]

Calculate the specific humidity from relative humidity, air temperature, and pressure.

Input: T is in K rh is in Pa/Pa p is in Pa

Output q in kg/kg

pyPamtra.meteoSI.rh_to_iwv(relhum_lev, temp_lev, press_lev, hgt_lev)[source]

Calculate the integrated water vapour

Input: T is in K rh is in Pa/Pa p is in Pa z is in m

Output iwv in kg/m^2

pyPamtra.meteoSI.vaphet(T)[source]

Compute specific heat of vaporisation Input : T [K] thermodynamic temperature Output : VAPHET [J/kg] specific heat of vaporisation Source: Liljequist, G.H. und K. Cehak, 1984: Allgemeine

Meteorologie (III.Ed.). Vieweg, 396p. Page 95

translated to Python from vaphet.pro by mx

Namelist Routines

Tools

pyPamtra.tools.formatExceptionInfo(maxTBlevel=5)[source]
class pyPamtra.tools.sftp2Cluster(machinename, username)[source]

Bases: object

mv(oldFile, newFile)[source]
put(filename, data)[source]
rm(filename, silent=True)[source]

Indices and tables