### =========================================================================== ### ### Build coordnates mask ### ### =========================================================================== # creates file coordinates_mask.nc used by mosaix from NEMO netcdf files # coordinates.nc and bathymetry.nc (on the same grid as coordinates.nc) ## ## MOSAIX is under CeCILL_V2 licence. See "Licence_CeCILL_V2-en.txt" ## file for an english version of the licence and ## "Licence_CeCILL_V2-fr.txt" for a french version. ## ## Permission is hereby granted, free of charge, to any person or ## organization obtaining a copy of the software and accompanying ## documentation covered by this license (the "Software") to use, ## reproduce, display, distribute, execute, and transmit the ## Software, and to prepare derivative works of the Software, and to ## permit third-parties to whom the Software is furnished to do so, ## all subject to the following: ## ## Warning, to install, configure, run, use any of MOSAIX software or ## to read the associated documentation you'll need at least one (1) ## brain in a reasonably working order. Lack of this implement will ## void any warranties (either express or implied). Authors assumes ## no responsability for errors, omissions, data loss, or any other ## consequences caused directly or indirectly by the usage of his ## software by incorrectly or partially configured ## ## import numpy as np import xarray as xr import nemo import datetime, os, platform, argparse ## SVN information __Author__ = "$Author$" __Date__ = "$Date$" __Revision__ = "$Revision$" __Id__ = "$Id$" __HeadURL__ = "$HeadURL: $" ### ### ===== Handling command line parameters ================================================== # Creating a parser parser = argparse.ArgumentParser ( description = """Read coordinates and bathymetry to build masks, grid bounds and areas for MOSAIX""", epilog='-------- End of the help message --------', formatter_class=argparse.ArgumentDefaultsHelpFormatter) # Adding arguments parser.add_argument ('--model' , help='oce model name', type=str, default='eORCA1.2', choices=('paleORCA', 'ORCA2.3', 'ORCA2.4', 'eORCA1.2', 'eORCA025', 'eORCA025.1') ) parser.add_argument ('--bathy' , help='bathymetry file name', type=str, default='bathy_meter.nc' ) parser.add_argument ('--coord' , help='coordinates file name', type=str, default='coordinates.nc' ) parser.add_argument ('--fout' , help='output file name (given as an input to MOSAIX)', type=str, default='coordinates_mask.nc' ) parser.add_argument ('--ocePerio', help='periodicity of ocean grid, if not given is imposed according to --model', type=int, default=0, choices=range(1, 7) ) parser.add_argument ('--nemo4Ubug',help='reproduce NEMO lbc_lnk bug for U grid and periodicity 4', type=bool, default=False, choices=(True, False) ) parser.add_argument ('--straits' , help='modify grid metrics for selected strait points (depends on model name)', type=bool, default=False, choices=(True, False) ) parser.add_argument ('--coordPerio', help='impose periodicity of coordinates and areas through lbc', type=bool, default=False, choices=(True, False) ) parser.add_argument ('--maskbathy', help='use the same formula as in NEMO to compute mask_T from the bathymetry file. Does not always give the same result as NEMO ???', type=bool, default=False, choices=(True, False) ) # Parse command line myargs = parser.parse_args() # ocean model name model = myargs.model # bathymetry input file n_bathy = myargs.bathy # coordinates input file n_coord = myargs.coord # coordinates and mask output file n_out = myargs.fout # reproduce periodicity lbclnk bug for U and nperio = 4 nemoUbug = myargs.nemo4Ubug # change metrics (and bathymetry) for straits straits = myargs.straits # periodicity imposed on coordinates, metrics and areas coordperio = myargs.coordPerio # function used to compute mask_T from bathymetry maskbathynemo = myargs.maskbathy # type of grid periodicity nperio = myargs.ocePerio # check of periodicity with type of grid if model in ('eORCA1.2', 'paleORCA', 'ORCA025', 'eORCA025', 'eORCA025.1') : if nperio != 0 : if nperio != 6 : print(f'Warning ! model = {model} and ocePerio = {nperio} instead of 6 !') else : nperio = 6 if model in ('ORCA2.3', 'ORCA2.4') : if nperio != 0 : if nperio != 4 : print(f'Attention ! model = {model} and ocePerio = {nperio} instead of 4 !') else : nperio = 4 print(f' model = {model}\n ocePerio = {nperio}\n bathy = {n_bathy}\n coord = {n_coord}\n fout = {n_out}') print(f' Switchs : nemo4Ubug = {nemoUbug}, straits = {straits}, coordPerio = {coordperio}, maskbathy = {maskbathynemo}') ## #!!! bathymetry and coordinates files ## # open input files while removing the time dimension f_coord = xr.open_dataset (n_coord, decode_times=False).squeeze() f_bathy = xr.open_dataset (n_bathy, decode_times=False).squeeze() # Suppress time if necessary try : del f_coord['time'] print ('time successfully removed') except : pass print ('failed to suppress time') # rename latitude, longitude and grid variables in bathymetry Bathymetry = f_bathy['Bathymetry'].copy() nav_lon_grid_T = f_bathy['nav_lon'].data nav_lat_grid_T = f_bathy['nav_lat'].data Bathymetry = xr.DataArray (Bathymetry, coords = { "nav_lat_grid_T": (["y", "x"], nav_lat_grid_T), "nav_lon_grid_T": (["y", "x"], nav_lon_grid_T) } ) Bathymetry = Bathymetry.rename ({'y':'y_grid_T', 'x':'x_grid_T'}) # modify Bathymetry in straits for select grids (might change computed mask_T) if straits : # Open straits for usual grids if model == 'ORCA2.3' : # orca_r2: Gibraltar strait open Bathymetry[101,139] = 284. # orca_r2: Bab el Mandeb strait open Bathymetry[87,159] = 137. # impose periodicity of bathymetry Bathymetry = nemo.lbc (Bathymetry, nperio=nperio, cd_type='T') # suppress ocean points at southernmost position only if nperio in {3, 4, 5, 6} if nperio in (3, 4, 5, 6) : Bathymetry[0,:] = 0.0 ## #!!! Create masks from bathymetry ## # Creation of mask_T following choosed option maskbathy if maskbathynemo : # Use same formula as domzgr. mask_T = xr.where (Bathymetry - 1. + 0.1 >= 0.0, 1, 0).astype (dtype='f4') else : mask_T = xr.where (Bathymetry > 0.0, 1, 0).astype (dtype='f4') # Creation of U, V, W, F masks from mask_T mask_U = mask_T * mask_T.shift (x_grid_T=-1) mask_V = mask_T * mask_T.shift (y_grid_T=-1) mask_F = mask_T * mask_T.shift (y_grid_T=-1) * mask_T.shift (x_grid_T=-1) * mask_T.shift (y_grid_T=-1, x_grid_T=-1) mask_W = mask_T # loop on TUVFW to modify coordinates and attributes of mask_[TUVFW] and maskutil_[TUVFW] for cd_type in ['T', 'U', 'V', 'F', 'W'] : MaskName = 'mask_' + cd_type UtilName = 'maskutil_' + cd_type # impose periodicity of chosen grid model on masks locals()[MaskName] = nemo.lbc (locals()[MaskName], nperio=nperio, cd_type=cd_type, nemo_4U_bug=nemoUbug).astype (dtype='f4') # rename masks coordinates if cd_type != 'T' : locals()[MaskName] = locals()[MaskName].rename \ ( {'y_grid_T' : 'y_grid_'+cd_type, 'x_grid_T' : 'x_grid_'+cd_type, 'nav_lat_grid_T': 'nav_lat_grid_'+cd_type, 'nav_lon_grid_T': 'nav_lon_grid_'+cd_type} ) # create masks without duplicate points : maskutil_[TUVWF] locals()[UtilName] = nemo.lbc_mask (locals()[MaskName].copy(), nperio=nperio, cd_type=cd_type) #set name attribute of mask dataset locals()[MaskName].name = MaskName locals()[UtilName].name = UtilName # remove _FillVallue from mskutil locals()[MaskName].encoding['_FillValue'] = None locals()[UtilName].encoding['_FillValue'] = None # add masks attributes locals()[MaskName].attrs['cell_measures'] = 'area: area_grid_'+cd_type locals()[UtilName].attrs['cell_measures'] = 'area: area_grid_'+cd_type ## #!!! create grid coordinates from NEMO coordinates.nc file ## angle = { 'lon' : 'glam', 'lat' : 'gphi' } gridv = { 'T' : 't', 'U' : 'u', 'V' : 'v', 'F' : 'f', 'W' : 't' } for cd_type in ['T', 'U', 'V', 'F', 'W'] : for dir_type in ['lon', 'lat'] : coord_name = 'nav_' + dir_type + '_grid_' + cd_type dir_name = angle[dir_type]+gridv[cd_type] # impose or not periodicity on coordinates read from file if coordperio : locals()[coord_name] = nemo.lbc (f_coord[dir_name].copy(), nperio=nperio, cd_type=cd_type, nemo_4U_bug=nemoUbug) else : locals()[coord_name] = f_coord[dir_name].copy() locals()[coord_name] = locals()[coord_name].rename( {'y':'y_grid_'+cd_type, 'x':'x_grid_'+cd_type} ) # remove _FillValue and missing_value locals()[coord_name].encoding['_FillValue'] = None locals()[coord_name].encoding['missing_value'] = None # define name attribute locals()[coord_name].name = coord_name # add coordinates attributes locals()[coord_name].attrs['bounds'] ='bounds_' + dir_type + '_grid_' + cd_type locals()['nav_lon_grid_'+cd_type].attrs['standard_name'] = 'longitude' locals()['nav_lon_grid_'+cd_type].attrs['long_name'] = 'Longitude' locals()['nav_lon_grid_'+cd_type].attrs['units'] = 'degrees_east' locals()['nav_lat_grid_'+cd_type].attrs['standard_name'] = 'latitude' locals()['nav_lat_grid_'+cd_type].attrs['long_name'] = 'Latitude' locals()['nav_lat_grid_'+cd_type].attrs['units'] = 'degrees_north' ## #!!! compute areas of cells at coordinate points ## # create areas variables for grid cells from NEMO coordinates.nc file # create new variables e1 e2 to keep f_coord the same for cd_type in ['t', 'u', 'v', 'f'] : for axis in ['1', '2'] : coordName = 'e' + axis + cd_type locals()[coordName]=f_coord[coordName].copy() # remove zero values from areas # need to be define for the extended grid south of -80S # some point are undefined but you need to have e1 and e2 .NE. 0 locals()[coordName]=xr.where(locals()[coordName] == 0.0, 1.0e2, locals()[coordName]) # Correct areas for straits if straits : # ORCA R2 configuration if model == 'ORCA2.3' : # Gibraltar : e2u reduced to 20 km e2u[101,138:140] = 20.e3 # Bab el Mandeb : e2u reduced to 30 km # e1v reduced to 18 km e1v[87,159] = 18.e3 e2u[87,159] = 30.e3 # Danish Straits: e2u reduced to 10 km e2u[115,144:146] = 10.e3 # ORCA R1 configuration if model == 'eORCA1.2' : # Gibraltar : e2u reduced to 20 km e2u[240,281:283] = 20.e3 # Bhosporus : e2u reduced to 10 km e2u[247,313:315] = 10.e3 # Lombok : e1v reduced to 13 km e1v[163:165,43] = 13.e3 # Sumba : e1v reduced to 8 km e1v[163:165,47] = 8.e3 # Ombai : e1v reduced to 13 km e1v[163:165,52] = 13.e3 # Timor Passage : e1v reduced to 20 km #e1v[163:165,55] = 20.e3 # W Halmahera : e1v reduced to 30 km e1v[180:182,54] = 30.e3 # E Halmahera : e1v reduced to 50 km e1v[180:182,57] = 50.e3 # ORCA R05 configuration if model == 'ORCA.05' : # Reduced e2u at the Gibraltar Strait e2u[326,562:564] = 20.e3 # Reduced e2u at the Bosphore Strait e2u[342,626:628] = 10.e3 # Reduced e2u at the Sumba Strait e2u[231,92:94] = 40.e3 # Reduced e2u at the Ombai Strait e2u[231,102] = 15.e3 # Reduced e2u at the Palk Strait e2u[269,14] = 10.e3 # Reduced e1v at the Lombok Strait e1v[231:233,86] = 10.e3 # Reduced e1v at the Bab el Mandeb e1v[275,661] = 25.e3 # compute cells areas for cd_type in ['T', 'U', 'V', 'F', 'W'] : areaName = 'area_grid_' + cd_type if coordperio : locals()[areaName] = nemo.lbc (locals()['e1'+gridv[cd_type]]*locals()['e2'+gridv[cd_type]], nperio=nperio, cd_type=cd_type, nemo_4U_bug=nemoUbug) else : locals()[areaName] = locals()['e1'+gridv[cd_type]]*locals()['e2'+gridv[cd_type]] # rename indices locals()[areaName] = locals()[areaName].rename ({'y':'y_grid_'+cd_type, 'x':'x_grid_'+cd_type}) # add attributes locals()[areaName].name = areaName locals()[areaName].attrs['standard_name'] = 'cell_area' locals()[areaName].attrs['units'] = 'm2' # remove fill values locals()[areaName].encoding['_FillValue'] = None ## #!!! compute grid bounds !!! ## #--------------------------------------------------------------- # function to generate grid bounds from NEMO coordinates.nc file def set_bounds (cdgrd) : ''' Constructs lon/lat bounds Bounds are numerated counter clockwise, from bottom left See NEMO file OPA_SRC/IOM/iom.F90, ROUTINE set_grid_bounds, for more details ''' # Define offset of coordinate representing bottom-left corner if cdgrd in ['T', 'W'] : icnr = -1 ; jcnr = -1 corner_lon = f_coord['glamf'].copy() ; corner_lat = f_coord['gphif'].copy() center_lon = f_coord['glamt'].copy() ; center_lat = f_coord['gphit'].copy() if cdgrd == 'U' : icnr = 0 ; jcnr = -1 corner_lon = f_coord['glamv'].copy() ; corner_lat = f_coord['gphiv'].copy() center_lon = f_coord['glamu'].copy() ; center_lat = f_coord['gphiu'].copy() if cdgrd == 'V' : icnr = -1 ; jcnr = 0 corner_lon = f_coord['glamu'].copy() ; corner_lat = f_coord['gphiu'].copy() center_lon = f_coord['glamv'].copy() ; center_lat = f_coord['gphiv'].copy() if cdgrd == 'F' : icnr = -1 ; jcnr = -1 corner_lon = f_coord['glamt'].copy() ; corner_lat = f_coord['gphit'].copy() center_lon = f_coord['glamf'].copy() ; center_lat = f_coord['gphif'].copy() jpj, jpi = corner_lon.shape ; nvertex = 4 dims = ['y_grid_' + cdgrd, 'x_grid_' + cdgrd, 'nvertex_grid_' + cdgrd] bounds_lon = xr.DataArray (np.zeros ((jpj, jpi, nvertex)), dims=dims) bounds_lat = xr.DataArray (np.zeros ((jpj, jpi, nvertex)), dims=dims) idx = [(jcnr,icnr), (jcnr,icnr+1), (jcnr+1,icnr+1), (jcnr+1,icnr)] # Compute cell vertices that can be defined, # and complete with periodicity for nn in range (nvertex) : tmp = np.roll (corner_lon, shift=tuple(-1*np.array(idx[nn])), axis=(-2,-1)) bounds_lon[1:jpj,1:jpi,nn] = tmp[1:jpj,1:jpi] tmp = np.roll (corner_lat, shift=tuple(-1*np.array(idx[nn])), axis=(-2,-1)) bounds_lat[1:jpj,1:jpi,nn] = tmp[1:jpj,1:jpi] bounds_lon[:,:,nn] = nemo.lbc (bounds_lon[:,:,nn], nperio=nperio, cd_type=cdgrd, nemo_4U_bug=nemoUbug) bounds_lat[:,:,nn] = nemo.lbc (bounds_lat[:,:,nn], nperio=nperio, cd_type=cdgrd, nemo_4U_bug=nemoUbug) # Zero-size cells at closed boundaries if cell points provided, # otherwise they are closed cells with unrealistic bounds if not (nperio == 1 or nperio == 4 or nperio == 6) : for nn in range (nvertex) : bounds_lon[:,0,nn] = center_lon[:,0] # (West or jpni = 1), closed E-W bounds_lat[:,0,nn] = center_lat[:,0] bounds_lon[:,jpi,nn] = center_lon[:,jpi] # (East or jpni = 1), closed E-W bounds_lat[:,jpi,nn] = center_lat[:,jpi] if nperio != 2 : for nn in range (nvertex) : bounds_lon[0,:,nn] = center_lon[0,:] # (South or jpnj = 1), not symmetric bounds_lat[0,:,nn] = center_lat[0,:] if nperio < 3 : for nn in range (nvertex) : bounds_lon[jpj,:,nn] = center_lon[jpj,:] # (North or jpnj = 1), no north fold bounds_lat[jpj,:,nn] = center_lat[jpj,:] # Rotate cells at the north fold if nperio >= 3 : # Working array for location of northfold z_fld = nemo.lbc (np.ones ((jpj, jpi)), nperio=nperio, cd_type=cdgrd, psgn=-1., nemo_4U_bug=nemoUbug) z_fld = np.repeat((z_fld == -1.0)[...,np.newaxis],4,axis=2) # circular shift of 2 indices in bounds third index bounds_lon_tmp = np.roll (bounds_lon, shift=-2, axis=2) bounds_lat_tmp = np.roll (bounds_lat, shift=-2, axis=2) bounds_lon[:,:,:] = np.where (z_fld, bounds_lon_tmp[:,:,:] , bounds_lon[:,:,:] ) bounds_lat[:,:,:] = np.where (z_fld, bounds_lat_tmp[:,:,:] , bounds_lat[:,:,:] ) # Invert cells at the symmetric equator if nperio == 2 : bounds_lon_tmp = np.roll (bounds_lon, shift=-2, axis=2) bounds_lat_tmp = np.roll (bounds_lat, shift=-2, axis=2) bounds_lon[0,:,:] = bounds_lon[0,:,:] bounds_lat[0,:,:] = bounds_lat[0,:,:] #bounds_lon.attrs['coordinates'] = 'nav_lat_grid_' + cdgrd + ' nav_lon_grid_' + cdgrd #bounds_lat.attrs['coordinates'] = 'nav_lat_grid_' + cdgrd + ' nav_lon_grid_' + cdgrd bounds_lon.attrs['units'] = 'degrees_east' bounds_lat.attrs['units'] = 'degrees_north' bounds_lon.name = 'bounds_lon_grid_' + cdgrd bounds_lat.name = 'bounds_lat_grid_' + cdgrd # remove _FillValue bounds_lon.encoding['_FillValue'] = None bounds_lat.encoding['_FillValue'] = None return bounds_lon, bounds_lat #------------------------------------------------------------ bounds_lon_grid_T, bounds_lat_grid_T = set_bounds ('T') bounds_lon_grid_U, bounds_lat_grid_U = set_bounds ('U') bounds_lon_grid_V, bounds_lat_grid_V = set_bounds ('V') bounds_lon_grid_W, bounds_lat_grid_W = set_bounds ('W') bounds_lon_grid_F, bounds_lat_grid_F = set_bounds ('F') # build xarray dataset to be saved ds = xr.Dataset ({ 'mask_T' : mask_T , 'mask_U' : mask_U , 'mask_V' : mask_V , 'mask_W' : mask_W , 'mask_F' : mask_F , 'area_grid_T' : area_grid_T, 'area_grid_U' : area_grid_U, 'area_grid_V' : area_grid_V, 'area_grid_W' : area_grid_W, 'area_grid_F' : area_grid_F, 'maskutil_T' : maskutil_T , 'maskutil_U' : maskutil_U , 'maskutil_V' : maskutil_V , 'maskutil_W' : maskutil_W , 'maskutil_F' : maskutil_F , 'bounds_lon_grid_T': bounds_lon_grid_T, 'bounds_lat_grid_T': bounds_lat_grid_T, 'bounds_lon_grid_U': bounds_lon_grid_U, 'bounds_lat_grid_U': bounds_lat_grid_U, 'bounds_lon_grid_V': bounds_lon_grid_V, 'bounds_lat_grid_V': bounds_lat_grid_V, 'bounds_lon_grid_W': bounds_lon_grid_W, 'bounds_lat_grid_W': bounds_lat_grid_W, 'bounds_lon_grid_F': bounds_lon_grid_F, 'bounds_lat_grid_F': bounds_lat_grid_F, }) #replace nav_lon nav_lat with variables obtained from NEMO coordinates.nc file #by construction nav_lon nav_lat come from the bathymetry for cd_type in ['T', 'U', 'V', 'F', 'W'] : for dir_type in ['lon', 'lat'] : coord_name = 'nav_' + dir_type + '_grid_' + cd_type ds.coords[coord_name]=locals()[coord_name] ds.attrs['name'] = 'coordinates_mask' ds.attrs['description'] = 'coordinates and mask for MOSAIX' ds.attrs['title'] = 'coordinates_mask' ds.attrs['source'] = 'IPSL Earth system model' ds.attrs['group'] = 'ICMC IPSL Climate Modelling Center' ds.attrs['Institution'] = 'IPSL https.//www.ipsl.fr' ds.attrs['Model'] = model ds.attrs['timeStamp'] = '{:%Y-%b-%d %H:%M:%S}'.format (datetime.datetime.now ()) ds.attrs['history'] = 'Build from ' + n_coord + ' and ' + n_bathy ds.attrs['directory'] = os.getcwd () try: ds.attrs['user'] = os.getlogin () except: ds.attrs['user'] = 'NoUser' ds.attrs['HOSTNAME'] = platform.node () ds.attrs['Python'] = 'Python version: ' + platform.python_version () ds.attrs['xarray'] = 'xarray version: ' + xr.__version__ ds.attrs['OS'] = platform.system () ds.attrs['release'] = platform.release () ds.attrs['hardware'] = platform.machine () ds.attrs['SVN_Author'] = "$Author$" ds.attrs['SVN_Date'] = "$Date$" ds.attrs['SVN_Revision'] = "$Revision$" ds.attrs['SVN_Id'] = "$Id$" ds.attrs['SVN_HeadURL'] = "$HeadURL: http://forge.ipsl.jussieu.fr/igcmg/svn/TOOLS/MOSAIX/Build_coordinates_mask.py $" # save to output file ds.to_netcdf (n_out)