ESMValGroup / ESMValCore

ESMValCore: A community tool for pre-processing data from Earth system models in CMIP and running analysis scripts.
https://www.esmvaltool.org
Apache License 2.0
40 stars 37 forks source link

Dataset problem: Missing latitude and longitude bounds in FGOALS-g3 tos and siconc #614

Open npgillett opened 4 years ago

npgillett commented 4 years ago

Coordinate variables latitude and longitude in FGOALS-g3 tos data do not have attribute bounds:

acrnrng@lxwrk2:~/cmip6_data/CMIP/CAS/FGOALS-g3/historical/r1i1p1f1/Omon/tos/gn/v20191107$ ncdump -h tos_Omon_FGOALS-g3_historical_r1i1p1f1_gn_185001-201412.nc
netcdf tos_Omon_FGOALS-g3_historical_r1i1p1f1_gn_185001-201412 {
dimensions:
        time = UNLIMITED ; // (1980 currently)
        j = 218 ;
        i = 360 ;
        bnds = 2 ;
variables:
        double time(time) ;
                time:bounds = "time_bnds" ;
                time:units = "days since 0001-01-01 00:00:00" ;
                time:calendar = "365_day" ;
                time:axis = "T" ;
                time:long_name = "time" ;
                time:standard_name = "time" ;
        double time_bnds(time, bnds) ;
        int j(j) ;
                j:units = "1" ;
                j:long_name = "cell index along second dimension" ;
        int i(i) ;
                i:units = "1" ;
                i:long_name = "cell index along first dimension" ;
        double latitude(j, i) ;
                latitude:standard_name = "latitude" ;
                latitude:long_name = "latitude" ;
                latitude:units = "degrees_north" ;
                latitude:missing_value = 1.e+20 ;
                latitude:_FillValue = 1.e+20 ;
        double longitude(j, i) ;
                longitude:standard_name = "longitude" ;
                longitude:long_name = "longitude" ;
                longitude:units = "degrees_east" ;
                longitude:missing_value = 1.e+20 ;
                longitude:_FillValue = 1.e+20 ;
        float tos(time, j, i) ;

This causes the following problem when esmvaltool is run:

  File "/home/rng/miniconda3/envs/esmvaltool2/lib/python3.7/site-packages/esmvalcore/preprocessor/_regrid_esmpy.py", line 130, in cube_to_empty_field
    circular = is_lon_circular(lon)
  File "/home/rng/miniconda3/envs/esmvaltool2/lib/python3.7/site-packages/esmvalcore/preprocessor/_regrid_esmpy.py", line 114, in is_lon_circular
    seam = (lon.bounds[1:-1, -1, (1, 2)]
TypeError: 'NoneType' object is not subscriptable
INFO    [7195] If you suspect this is a bug or need help, please open an issue on https://github.com/ESMValGroup/ESMValTool/issues and attach the run/recipe_*.yml and run/main_log_de
bug.txt files from the output directory.

main_log_debug.txt

I tried using this fix for BCC-CSM1-1: https://github.com/ESMValGroup/ESMValCore/blob/master/esmvalcore/cmor/_fixes/cmip5/bcc_csm1_1.py

But that didn't work because it just interpolates the bounds used forgrid_latitude and grid_longitude, and these bounds aren't set in the case of FGOALS-g3.

How should I proceed? Is it possible to guess the bounds in some other way, or is this a case where we just have to ask the modelling centre to fix this and republish the data? siconc also seems to have the same problem.

Thanks!

Relative path to data: CAS/FGOALS-g3/historical/r1i1p1f1/Omon/tos/gn/v20191107/tos_Omon_FGOALS-g3_historical_r1i1p1f1_gn_185001-201412.nc

npgillett commented 4 years ago

I worked out a way round this - by just guessing the bounds for grid_latitude and grid_longitude:


class OceanFixGrid(Fix):

    """Fixes for tos, siconc in FGOALS-g3."""

    def fix_data(self, cube):
        """Fix data.
        Calculate missing latitude/longitude boundaries using interpolation.
        Based on a similar fix for BCC-CSM2-MR
        Parameters
        ----------
        cube: iris.cube.Cube
            Input cube to fix.
        Returns
        -------
        iris.cube.Cube
        """

        rlat = cube.coord('grid_latitude').points
        rlon = cube.coord('grid_longitude').points

        #Guess coordinate bounds in rlat, rlon (following BCC-CSM2-MR-1).
        rlat_idx_bnds = np.zeros((len(rlat),2))
        rlat_idx_bnds[:,0]=np.arange(len(rlat))-0.5
        rlat_idx_bnds[:,1]=np.arange(len(rlat))+0.5
        rlat_idx_bnds[0,0]=0.
        rlat_idx_bnds[len(rlat)-1,1]=len(rlat)
        rlon_idx_bnds = np.zeros((len(rlon),2))
        rlon_idx_bnds[:,0]=np.arange(len(rlon))-0.5
        rlon_idx_bnds[:,1]=np.arange(len(rlon))+0.5

        # Calculate latitude/longitude vertices by interpolation
        lat_vertices = []
        lon_vertices = []
        for (i, j) in [(0, 0), (0, 1), (1, 1), (1, 0)]:
            (rlat_v, rlon_v) = np.meshgrid(rlat_idx_bnds[:, i],
                                           rlon_idx_bnds[:, j],
                                           indexing='ij')
            lat_vertices.append(
                map_coordinates(cube.coord('latitude').points,
                                [rlat_v, rlon_v],
                                mode='nearest'))
            lon_vertices.append(
                map_coordinates(cube.coord('longitude').points,
                                [rlat_v, rlon_v],
                                mode='wrap'))
        lat_vertices = np.array(lat_vertices)
        lon_vertices = np.array(lon_vertices)
        lat_vertices = np.moveaxis(lat_vertices, 0, -1)
        lon_vertices = np.moveaxis(lon_vertices, 0, -1)

        # Copy vertices to cube
        cube.coord('latitude').bounds = lat_vertices
        cube.coord('longitude').bounds = lon_vertices

        return cube

Once I have checked over my CMIP6 fixes again, I plan to put in a pull request for them.