"""This module contains I/O operations for loading files that were created for theAtmospheric Radiation Measurement program supported by the Department of EnergyOffice of Science."""importcopyimportglobimportjsonimportreimporttarfileimporttempfileimporturllibimportwarningsfromdatetimeimportdatetime,timezonefromosimportPathLikefrompathlibimportPath,PosixPathimportnumpyasnpimportxarrayasxrfromcftimeimportnum2datefromnetCDF4importDatasetimportactimportact.utilsasutilsfromact.configimportDEFAULT_DATASTREAM_NAMEfromact.utils.io_utilsimportcleanup_files,is_gunzip_file,unpack_gzip,unpack_tar
[docs]defread_arm_netcdf(filenames,concat_dim=None,return_None=False,combine='by_coords',decode_times=True,use_cftime=True,use_base_time=False,combine_attrs='override',cleanup_qc=False,keep_variables=None,**kwargs,):""" Returns `xarray.Dataset` with stored data and metadata from a user-defined query of ARM-standard netCDF files from a single datastream. Has some procedures to ensure time is correctly fomatted in returned Dataset. Parameters ---------- filenames : str, pathlib.PosixPath, list of str, list of pathlib.PosixPath Name of file(s) to read. concat_dim : str Dimension to concatenate files along. return_None : boolean Catch IOError exception when file not found and return None. Default is False. combine : str String used by xarray.open_mfdataset() to determine how to combine data files into one Dataset. See Xarray documentation for options. decode_times : boolean Standard Xarray option to decode time values from int/float to python datetime values. Appears the default is to do this anyway but need this option to allow correct usage of use_base_time. use_cftime : boolean Option to use cftime library to parse the time units string and correctly establish the time values with a units string containing timezone offset. This is used because the Pandas units string parser does not correctly recognize time zone offset. Code will automatically detect cftime object and convert to datetime64 in returned Dataset. use_base_time : boolean Option to use ARM time variables base_time and time_offset. Useful when the time variable is not included (older files) or when the units attribute is incorrectly formatted. Will use the values of base_time and time_offset as seconds since epoch and create datetime64 values for time coordinate. If set will change decode_times and use_cftime to False. combine_attrs : str String indicating how to combine attrs of the datasets being merged cleanup_qc : boolean Call clean.cleanup() method to convert to standardized ancillary quality control variables. This will not allow any keyword options, so if non-default behavior is desired will need to call clean.cleanup() method on the dataset after reading the data. keep_variables : str or list of str Variable names to read from data file. Works by creating a list of variable names to exclude from reading and passing into open_mfdataset() via drop_variables keyword. Still allows use of drop_variables keyword for variables not listed in first file to read. **kwargs : keywords Keywords to pass through to xarray.open_mfdataset(). Returns ------- ds : xarray.Dataset (or None) ACT Xarray dataset (or None if no data file(s) found). Examples -------- This example will load the example sounding data used for unit testing. .. code-block :: python import act ds = act.io.arm.read_arm_netcdf(act.tests.sample_files.EXAMPLE_SONDE_WILDCARD) print(ds) """ds=Nonefilenames,cleanup_temp_directory=check_if_tar_gz_file(filenames)file_dates=[]file_times=[]# If requested to use base_time and time_offset, set keywords to correct attribute values# to pass into xarray open_mfdataset(). Need to turn off decode_times and use_cftime# or else will try to convert base_time and time_offset. Depending on values of attributes# may cause a failure.ifuse_base_time:decode_times=Falseuse_cftime=Falseifdecode_timesanduse_cftime:decode_times=xr.coders.CFDatetimeCoder(use_cftime=True)# Add function keywords to kwargs dictionary for passing into open_mfdataset.kwargs['combine']=combinekwargs['concat_dim']=concat_dimkwargs['decode_times']=decode_timesiflen(filenames)>1andnotisinstance(filenames,str):kwargs['combine_attrs']=combine_attrs# Check if keep_variables is set. If so determine correct drop_variablesifkeep_variablesisnotNone:drop_variables=Noneif'drop_variables'inkwargs.keys():drop_variables=kwargs['drop_variables']kwargs['drop_variables']=keep_variables_to_drop_variables(filenames,keep_variables,drop_variables=drop_variables)# Create an exception tuple to use with try statements. Doing it this way# so we can add the FileNotFoundError if requested. Can add more error# handling in the future.except_tuple=(ValueError,)ifreturn_None:except_tuple=except_tuple+(FileNotFoundError,OSError)try:# Read data file with Xarray functionds=xr.open_mfdataset(filenames,**kwargs)exceptexcept_tupleasexception:# If requested return None for File not found erroriftype(exception).__name__=='FileNotFoundError':returnNone# If requested return None for File not found erroriftype(exception).__name__=='OSError'andexception.args[0]=='no files to open':returnNone# Look at error message and see if could be nested error message. If so# update combine keyword and try again. This should allow reading files# without a time variable but base_time and time_offset variables.if(kwargs['combine']!='nested'andtype(exception).__name__=='ValueError'andexception.args[0]=='Could not find any dimension coordinates ''to use to order the datasets for concatenation'):kwargs['combine']='nested'ds=xr.open_mfdataset(filenames,**kwargs)else:# When all else fails raise the orginal exceptionraiseexception# If requested use base_time and time_offset to derive time. Assumes that the units# of both are in seconds and that the value is number of seconds since epoch.ifuse_base_time:time=num2date(ds['base_time'].values+ds['time_offset'].values,ds['base_time'].attrs['units'])time=time.astype('datetime64[ns]')# Need to use a new Dataset creation to correctly index time for use with# .group and .resample methods in Xarray Datasets.temp_ds=xr.Dataset({'time':(ds['time'].dims,time,ds['time'].attrs)})ds['time']=temp_ds['time']deltemp_dsforatt_namein['units','ancillary_variables']:try:delds['time'].attrs[att_name]exceptKeyError:pass# Xarray has issues reading a CF formatted time units string if it contains# timezone offset without a [+|-] preceeding timezone offset.# https://github.com/pydata/xarray/issues/3644# To ensure the times are read in correctly need to set use_cftime=True.# This will read in time as cftime object. But Xarray uses numpy datetime64# natively. This will convert the cftime time values to numpy datetime64.desired_time_precision='datetime64[ns]'forvar_namein['time','time_offset']:try:if'time'inds.dimsandtype(ds[var_name].values[0]).__module__.startswith('cftime.'):# If we just convert time to datetime64 the group, sel, and other Xarray# methods will not work correctly because time is not indexed. Need to# use the formation of a Dataset to correctly set the time indexing.temp_ds=xr.Dataset({var_name:(ds[var_name].dims,ds[var_name].values.astype(desired_time_precision),ds[var_name].attrs,)})ds[var_name]=temp_ds[var_name]deltemp_ds# If time_offset is in file try to convert base_time as wellifvar_name=='time_offset':ds['base_time'].values=ds['base_time'].values.astype(desired_time_precision)ds['base_time']=ds['base_time'].astype(desired_time_precision)exceptKeyError:pass# Check if "time" variable is not in the netCDF file. If so try to use# base_time and time_offset to make time variable. Basically a fix for incorrectly# formatted files. May require using decode_times=False to initially read the data.if'time'inds.dimsandnotnp.issubdtype(ds['time'].dtype,np.datetime64):try:ds['time']=ds['time_offset']except(KeyError,ValueError):pass# Adding support for wildcardsifisinstance(filenames,str):filenames=glob.glob(filenames)elifisinstance(filenames,PosixPath):filenames=[filenames]# Get file dates and times that were read in to the datasetfilenames.sort()forfinfilenames:f=Path(f).namepts=re.match(r'(^[a-zA-Z0-9]+)\.([0-9a-z]{2})\.([\d]{8})\.([\d]{6})\.([a-z]{2,3}$)',f)# If Not ARM format, read in first time for infoifptsisnotNone:pts=pts.groups()file_dates.append(pts[2])file_times.append(pts[3])else:iflen(ds['time'].shape)>0:dummy=ds['time'].values[0]else:dummy=ds['time'].valuesfile_dates.append(utils.numpy_to_arm_date(dummy))file_times.append(utils.numpy_to_arm_date(dummy,returnTime=True))# Add attributesds.attrs['_file_dates']=file_datesds.attrs['_file_times']=file_timesis_arm_file_flag=check_arm_standards(ds)# Ensure that we have _datastream set whether or no there's# a datastream attribute already.ifis_arm_file_flag==0:ds.attrs['_datastream']=DEFAULT_DATASTREAM_NAMEelse:ds.attrs['_datastream']=ds.attrs['datastream']ds.attrs['_arm_standards_flag']=is_arm_file_flagifcleanup_qc:ds.clean.cleanup()ifcleanup_temp_directory:cleanup_files(files=filenames)returnds
defkeep_variables_to_drop_variables(filenames,keep_variables,drop_variables=None):""" Returns a list of variable names to exclude from reading by passing into `Xarray.open_dataset` drop_variables keyword. This can greatly help reduce loading time and disk space use of the Dataset. When passed a netCDF file name, will open the file using the netCDF4 library to get list of variable names. There is less overhead reading the varible names using netCDF4 library than Xarray. If more than one filename is provided or string is used for shell syntax globbing, will use the first file in the list. Parameters ---------- filenames : str, pathlib.PosixPath or list of str Name of file(s) to read. keep_variables : str or list of str Variable names desired to keep. Do not need to list associated dimention names. These will be automatically kept as well. drop_variables : str or list of str Variable names to explicitly add to returned list. May be helpful if a variable exists in a file that is not in the first file in the list. Returns ------- drop_vars : list of str Variable names to exclude from returned Dataset by using drop_variables keyword when calling Xarray.open_dataset(). Examples -------- .. code-block :: python import act filename = '/data/datastream/hou/houkasacrcfrM1.a1/houkasacrcfrM1.a1.20220404.*.nc' drop_vars = act.io.arm.keep_variables_to_drop_variables( filename, ['lat','lon','alt','crosspolar_differential_phase'], drop_variables='variable_name_that_only_exists_in_last_file_of_the_day') """read_variables=[]return_variables=[]ifisinstance(keep_variables,str):keep_variables=[keep_variables]ifisinstance(drop_variables,str):drop_variables=[drop_variables]# If filenames is a list subset to first file name.ifisinstance(filenames,(list,tuple)):filename=filenames[0]# If filenames is a string, check if it needs to be expanded in shell# first. Then use first returned file name. Else use the string filename.elifisinstance(filenames,str):filename=glob.glob(filenames)iflen(filename)==0:returnreturn_variableselse:filename.sort()filename=filename[0]# Use netCDF4 library to extract the variable and dimension names.rootgrp=Dataset(filename,'r')read_variables=list(rootgrp.variables)# Loop over the variables to exclude needed coordinate dimention names.dims_to_keep=[]forvar_nameinkeep_variables:try:dims_to_keep.extend(list(rootgrp[var_name].dimensions))exceptIndexError:passrootgrp.close()# Remove names not matching keep_varibles excluding the associated coordinate dimentionsreturn_variables=set(read_variables)-set(keep_variables)-set(dims_to_keep)# Add drop_variables to listifdrop_variablesisnotNone:return_variables=set(return_variables)|set(drop_variables)returnlist(return_variables)
[docs]defcheck_arm_standards(ds):""" Checks to see if an xarray dataset conforms to ARM standards. Parameters ---------- ds : Xarray Dataset The dataset to check. Returns ------- flag : int The flag corresponding to whether or not the file conforms to ARM standards. Bit packed, so 0 for no, 1 for yes """the_flag=1<<0if'datastream'notinds.attrs.keys():the_flag=0# Check if the historical global attribute name is# used instead of updated name of 'datastream'. If so# correct the global attributes and flip flag.if'zeb_platform'inds.attrs.keys():ds.attrs['datastream']=copy.copy(ds.attrs['zeb_platform'])delds.attrs['zeb_platform']the_flag=1<<0returnthe_flag
[docs]defcreate_ds_from_arm_dod(proc,set_dims,version='',fill_value=-9999.0,scalar_fill_dim=None,local_file=False):""" Queries the ARM DOD api and builds a dataset based on the ARM DOD and the dimension sizes that are passed in. Parameters ---------- proc : string Process to create the dataset off of. This is normally in the format of inst.level. i.e. vdis.b1 or kazrge.a1. If local file is true, this points to the path of the .dod file. set_dims : dict Dictionary of dims from the DOD and the corresponding sizes. Time is required. Code will try and pull from DOD, unless set through this variable Note: names need to match exactly what is in the dod i.e. {'drop_diameter': 50, 'time': 1440} version : string Version number of the ingest to use. If not set, defaults to latest version fill_value : float Fill value for non-dimension variables. Dimensions cannot have duplicate values and are incrementally set (0, 1, 2) scalar_fill_dim : str Depending on how the dataset is set up, sometimes the scalar values are dimensioned to the main dimension. i.e. a lat/lon is set to have a dimension of time. This is a way to set it up similarly. local_file: bool If true, the DOD will be loaded from a file whose name is proc. If false, the DOD will be pulled from PCM. Returns ------- ds : xarray.Dataset ACT Xarray dataset populated with all variables and attributes. Examples -------- .. code-block :: python dims = {'time': 1440, 'drop_diameter': 50} ds = act.io.arm.create_ds_from_arm_dod( 'vdis.b1', dims, version='1.2', scalar_fill_dim='time') """# Set base url to get DOD informationiflocal_fileisFalse:base_url='https://pcm.arm.gov/pcm/api/dods/'# Get data from DOD apiwithurllib.request.urlopen(base_url+proc)asurl:data=json.loads(url.read().decode())else:withopen(proc)asfile:data=json.loads(file.read())# Check version numbers and alert if requested version in not availablekeys=list(data['versions'].keys())ifversionnotinkeys:warnings.warn(' '.join(['Version:',version,'not available or not specified. Using Version:',keys[-1]]),UserWarning,)version=keys[-1]# Create empty xarray datasetds=xr.Dataset()# Get the global attributes and add to datasetatts={}foraindata['versions'][version]['atts']:ifa['name']=='string':continueifa['value']isNone:a['value']=''atts[a['name']]=a['value']ds.attrs=atts# Get variable information and create dataarrays that are# then added to the dataset# If not passed in through set_dims, will look to the DOD# if not set in the DOD, then will raise errorvariables=data['versions'][version]['vars']dod_dims=data['versions'][version]['dims']fordindod_dims:ifd['name']notinlist(set_dims.keys()):ifd['length']>0:set_dims[d['name']]=d['length']else:raiseValueError('Dimension length not set in DOD for '+d['name']+', nor passed in through set_dim')forvinvariables:dims=v['dims']dim_shape=[]# Using provided dimension data, fill array accordingly for easy overwriteiflen(dims)==0:ifscalar_fill_dimisNone:data_na=fill_valueelse:data_na=np.full(set_dims[scalar_fill_dim],fill_value)v['dims']=scalar_fill_dimelse:fordindims:dim_shape.append(set_dims[d])iflen(dim_shape)==1andv['name']==dims[0]:data_na=np.arange(dim_shape[0])else:data_na=np.full(dim_shape,fill_value)# Get attribute information. Had to do some things to get to print to netcdfatts={}str_flag=Falseforainv['atts']:ifa['name']=='string':str_flag=Truecontinueifa['value']isNone:continueifstr_flaganda['name']=='units':continueatts[a['name']]=a['value']da=xr.DataArray(data=data_na,dims=v['dims'],name=v['name'],attrs=atts)ds[v['name']]=dareturnds
[docs]@xr.register_dataset_accessor('write')classWriteDataset:""" Class for cleaning up Dataset before writing to file. """def__init__(self,xarray_ds):self._ds=xarray_ds
[docs]defwrite_netcdf(self,cleanup_global_atts=True,cleanup_qc_atts=True,join_char='__',make_copy=True,cf_compliant=False,delete_global_attrs=['qc_standards_version','qc_method','qc_comment'],FillValue=True,cf_convention='CF-1.8',encoding={},**kwargs,):""" This is a wrapper around Dataset.to_netcdf to clean up the Dataset before writing to disk. Some things are added to global attributes during ACT reading process, and QC variables attributes are modified during QC cleanup process. This will modify before writing to disk to better match Climate & Forecast standards. Parameters ---------- cleanup_global_atts : boolean Option to cleanup global attributes by removing any global attribute that starts with an underscore. cleanup_qc_atts : boolean Option to convert attributes that would be written as string array to be a single character string. CF 1.7 does not allow string attribures. Will use a single space a delimeter between values and join_char to replace white space between words. join_char : str The character sting to use for replacing white spaces between words when converting a list of strings to single character string attributes. Main use is with the flag_meanings attribute. make_copy : boolean Make a copy before modifying Dataset to write. For large Datasets this may add processing time and memory. If modifying the Dataset is OK try setting to False. cf_compliant : boolean Option to output file with additional attributes to make file Climate & Forecast complient. May require runing .clean.cleanup() method on the dataset to fix other issues first. This does the best it can but it may not be truely complient. You should read the CF documents and try to make complient before writing to file. delete_global_attrs : list Optional global attributes to be deleted. Defaults to some standard QC attributes that are not needed. Can add more or set to None to not remove the attributes. FillValue : boolean Xarray assumes all float type variables had the missing value indicator converted to NaN upon reading. to_netcdf() will then write a _FillValue attribute set to NaN. Set FillValue to False to supress adding the _FillValue=NaN variable attribute to the written file. Set to True to allow to_netcdf() to add the attribute. If the Dataset variable already has a _FillValue attribute or a _FillValue key is provided in the encoding dictionary those will not be changed and a _FillValue will be written to NetCDF file. cf_convention : str The Climate and Forecast convention string to add to Conventions attribute. encoding : dict The encoding dictionary used with to_netcdf() method. **kwargs : keywords Keywords to pass through to Dataset.to_netcdf() Examples -------- .. code-block :: python ds.write.write_netcdf(path='output.nc') """ifmake_copy:ds=copy.deepcopy(self._ds)else:ds=self._dsifcleanup_global_atts:forattrinlist(ds.attrs):ifattr.startswith('_'):delds.attrs[attr]ifcleanup_qc_atts:check_atts=['flag_meanings','flag_assessments']forvar_nameinlist(ds.data_vars):if'standard_name'notinds[var_name].attrs.keys():continueifds[var_name].attrs['standard_name']!="quality_flag":continueforattr_nameincheck_atts:try:att_values=ds[var_name].attrs[attr_name]ifisinstance(att_values,(list,tuple)):att_values=[att_value.replace(' ',join_char)foratt_valueinatt_values]ds[var_name].attrs[attr_name]=' '.join(att_values)exceptKeyError:pass# Xarray makes an assumption that float type variables were read in and converted# missing value indicator to NaN. .to_netcdf() will then automatically assign# _FillValue attribute set to NaN when writing. If requested will set _FillValue# key in encoding to None which will supress to_netcdf() from adding a _FillValue.# If _FillValue attribute or _FillValue key in encoding is already set, will not# override and the _FillValue will be written to the file.ifnotFillValue:all_var_names=list(ds.coords.keys())+list(ds.data_vars)forvar_nameinall_var_names:if'_FillValue'inds[var_name].attrs:continueifvar_namenotinencoding.keys():encoding[var_name]={'_FillValue':None}elif'_FillValue'notinencoding[var_name].keys():encoding[var_name]['_FillValue']=Noneifdelete_global_attrsisnotNone:forattrindelete_global_attrs:try:delds.attrs[attr]exceptKeyError:passforvar_nameinlist(ds.keys()):if'string'inlist(ds[var_name].attrs.keys()):att=ds[var_name].attrs['string']ds[var_name].attrs[var_name+'_string']=attdelds[var_name].attrs['string']# If requested update global attributes and variables attributes for required# CF attributes.ifcf_compliant:# Get variable names and standard name for each variablevar_names=list(ds.keys())standard_names=[]forvar_nameinvar_names:try:standard_names.append(ds[var_name].attrs['standard_name'])exceptKeyError:standard_names.append(None)# Check if time varible has axis and standard_name attributecoord_name='time'try:ds[coord_name].attrs['axis']exceptKeyError:try:ds[coord_name].attrs['axis']='T'exceptKeyError:passtry:ds[coord_name].attrs['standard_name']exceptKeyError:try:ds[coord_name].attrs['standard_name']='time'exceptKeyError:pass# Try to determine type of dataset by coordinate dimention named time# and other factorstry:ds.attrs['FeatureType']exceptKeyError:dim_names=list(ds.dims)FeatureType=Noneifdim_names==['time']:FeatureType='timeSeries'eliflen(dim_names)==2and'time'indim_namesand'bound'indim_names:FeatureType='timeSeries'eliflen(dim_names)>=2and'time'indim_names:forvar_nameinvar_names:dims=list(ds[var_name].dims)iflen(dims)==2and'time'indims:prof_dim=list(set(dims)-{'time'})[0]ifds[prof_dim].values.size>2:FeatureType='timeSeriesProfile'breakifFeatureTypeisnotNone:ds.attrs['FeatureType']=FeatureType# Add axis and positive attributes to variables with standard_name# equal to 'altitude'alt_variables=[var_names[ii]forii,sninenumerate(standard_names)ifsn=='altitude']forvar_nameinalt_variables:try:ds[var_name].attrs['axis']exceptKeyError:ds[var_name].attrs['axis']='Z'try:ds[var_name].attrs['positive']exceptKeyError:ds[var_name].attrs['positive']='up'# Check if the Conventions global attribute lists the CF conventiontry:Conventions=ds.attrs['Conventions']Conventions=Conventions.split()cf_listed=FalseforiiinConventions:ifii.startswith('CF-'):cf_listed=Truebreakifnotcf_listed:Conventions.append(cf_convention)ds.attrs['Conventions']=' '.join(Conventions)exceptKeyError:ds.attrs['Conventions']=str(cf_convention)# Reorder global attributes to ensure history is lasttry:history=copy.copy(ds.attrs['history'])delds.attrs['history']ds.attrs['history']=historyexceptKeyError:passcurrent_time=datetime.now(timezone.utc).replace(microsecond=0)history_value=(f'Written to file by ACT-{act.__version__} 'f'with write_netcdf() at {current_time} UTC')if'history'inlist(ds.attrs.keys()):ds.attrs['history']+=f" ; {history_value}"else:ds.attrs['history']=history_valueds.to_netcdf(encoding=encoding,**kwargs)
[docs]defcheck_if_tar_gz_file(filenames):""" Unpacks gunzip and/or TAR file contents and returns Xarray Dataset ... Parameters ---------- filenames : str, pathlib.Path Filenames to check if gunzip and/or tar files. Returns ------- filenames : Paths to extracted files from gunzip or TAR files """cleanup=Falseifisinstance(filenames,(str,PathLike)):try:ifis_gunzip_file(filenames)ortarfile.is_tarfile(str(filenames)):tmpdirname=tempfile.mkdtemp()cleanup=Trueifis_gunzip_file(filenames):filenames=unpack_gzip(filenames,write_directory=tmpdirname)iftarfile.is_tarfile(str(filenames)):filenames=unpack_tar(filenames,write_directory=tmpdirname,randomize=False)exceptException:passreturnfilenames,cleanup
[docs]defread_arm_mmcr(filenames):""" Reads in ARM MMCR files and splits up the variables into specific mode variables based on what's in the files. MMCR files have the modes interleaved and are not readable using xarray so some modifications are needed ahead of time. Parameters ---------- filenames : str, pathlib.PosixPath or list of str Name of file(s) to read. Returns ------- ds : xarray.Dataset (or None) ACT Xarray dataset (or None if no data file(s) found). """# Sort the files to make sure they concatenate rightfilenames.sort()# Run through each file and read it in using netCDF4, then# read it in with xarraymulti_ds=[]forfinfilenames:nc=Dataset(f,'a',diskless=True)# Change heights name to range to read appropriately to xarrayif'heights'innc.dimensions:nc.renameDimension('heights','range')ifncisnotNone:ds=xr.open_dataset(xr.backends.NetCDF4DataStore(nc))multi_ds.append(ds)# Concatenate datasets togetheriflen(multi_ds)>1:ds=xr.concat(multi_ds,dim='time')else:ds=multi_ds[0]# Get mdoes and ranges with time/height modesmodes=ds['mode'].valuesmode_vars=[]forvinds:if'range'inds[v].dimsand'time'inds[v].dimsandlen(ds[v].dims)==2:mode_vars.append(v)# For each mode, run extract data variables if available# saves as individual variables in the file.forminmodes:iflen(ds['ModeDescription'].shape)>1:mode_desc=ds['ModeDescription'].values[0,m]ifnp.isnan(ds['heights'].values[0,m,:]).all():continuerange_data=ds['heights'].values[0,m,:]else:mode_desc=ds['ModeDescription'].values[m]ifnp.isnan(ds['heights'].values[m,:]).all():continuerange_data=ds['heights'].values[m,:]mode_desc=str(mode_desc).split('_')[-1][0:-1]mode_desc=str(mode_desc).split('\'')[0]idx=np.where(ds['ModeNum'].values==m)[0]idy=np.where(~np.isnan(range_data))[0]forvinmode_vars:new_var_name=v+'_'+mode_desctime_name='time_'+mode_descrange_name='range_'+mode_descdata=ds[v].values[idx,:]data=data[:,idy]attrs=ds[v].attrsda=xr.DataArray(data=data,coords={time_name:ds['time'].values[idx],range_name:range_data[idy]},dims=[time_name,range_name],attrs=attrs,)ds[new_var_name]=dareturnds