cmorizer.py 24.1 KB
Newer Older
# coding=utf-8
import glob
import shutil
import uuid

import netCDF4
import os
from datetime import datetime

import pygrib
from autosubmit.config.log import Log
from autosubmit.date.chunk_date_lib import parse_date, chunk_end_date, previous_day, date2str, add_months

from earthdiagnostics.variable import Variable
from earthdiagnostics.utils import TempFile, Utils


class Cmorizer(object):
Javier Vegas-Regidor's avatar
Javier Vegas-Regidor committed
    """
    Class to manage CMORization

    :param data_manager: experiment's data manager
    :type data_manager: DataManager
    :param startdate: startdate to cmorize
    :type startdate: str
    :param member: member to cmorize
    :type member: int

    """

    def __init__(self, data_manager, startdate, member):
        self.data_manager = data_manager
        self.startdate = startdate
        self.member = member
        self.config = data_manager.config
        self.experiment = self.config.experiment
        self.cmor = self.config.cmor
        self.member_str = self.experiment.get_member_str(member)
        self.original_files_path = os.path.join(self.config.data_dir, self.experiment.expid, 'original_files',
                                                self.startdate, self.member_str, 'outputs')

    def cmorize_ocean(self):
Javier Vegas-Regidor's avatar
Javier Vegas-Regidor committed
        """
        CMORizes ocean files from MMO files
        :return:
        """
        if not self.cmor.ocean:
            return
        self._unpack_ocean_files('MMO')
        self._unpack_ocean_files('diags')

    def _unpack_ocean_files(self, prefix):
        tar_folder = os.path.join(self.original_files_path, '{0}*'.format(prefix))
        tar_files = glob.glob(tar_folder)
        tar_files.sort()
        errors = list()
        count = 1
        for tarfile in tar_files:
            Log.info('Unpacking oceanic file {0}/{1}'.format(count, len(tar_files)))
            self._unpack_tar(tarfile)
            Log.result('Oceanic file {0}/{1} finished'.format(count, len(tar_files)))
            count += 1
        return errors

    def _unpack_tar(self, tarfile):
        Log.info('Unpacking {0}', tarfile)

        scratch_dir = os.path.join(self.config.scratch_dir, 'CMOR')

        if os.path.exists(scratch_dir):
            shutil.rmtree(scratch_dir)
        os.makedirs(scratch_dir)
        Utils.untar((tarfile,), scratch_dir)
        errors = Utils.unzip(glob.glob(os.path.join(scratch_dir, '*.gz')))

        if os.path.basename(tarfile).startswith('MMA'):
            temp = TempFile.get()
            for filename in glob.glob(os.path.join(scratch_dir, 'MMA_*_SH_*.nc')):
                Utils.cdo.sp2gpl(options='-O', input=filename, output=temp)
                shutil.move(temp, filename)

            sh_files = glob.glob(os.path.join(scratch_dir, 'MMA_*_SH_*.nc'))
            Utils.cdo.mergetime(input=sh_files, output=os.path.join(scratch_dir, 'sh.nc'))

            gg_files = glob.glob(os.path.join(scratch_dir, 'MMA_*_GG_*.nc'))
            Utils.cdo.mergetime(input=gg_files, output=os.path.join(scratch_dir, 'gg.nc'))

            for filename in sh_files + gg_files:
                os.remove(filename)

            Utils.nco.ncks(input=os.path.join(scratch_dir, 'sh.nc'),
                           output=os.path.join(scratch_dir, 'gg.nc'), options='-A')
            os.remove(os.path.join(scratch_dir, 'sh.nc'))

            tar_startdate = tarfile[0:-4].split('_')[5].split('-')
            new_name = 'MMA_1m_{0[0]}_{0[1]}.nc'.format(tar_startdate)
            shutil.move(os.path.join(scratch_dir, 'gg.nc'), os.path.join(scratch_dir, new_name))

        for filename in glob.glob(os.path.join(scratch_dir, '*.nc')):
            self._cmorize_nc_file(filename)
        return errors

    def cmorize_atmos(self):
Javier Vegas-Regidor's avatar
Javier Vegas-Regidor committed
        """
        CMORizes atmospheric data, from grib or MMA files
        :return:
        """
        if not self.cmor.atmosphere:
            return

        if self.cmor.use_grib:
            grb_path = os.path.join(self.original_files_path, '*.grb')
            gribfiles = glob.glob(grb_path)
        else:
            gribfiles = tuple()

        if len(gribfiles) == 0:
            tar_files = glob.glob(
                os.path.join(self.original_files_path, 'MMA*'))
            tar_files.sort()
            count = 1
            for tarfile in tar_files:
                Log.info('Unpacking atmospheric file {0}/{1}'.format(count, len(tar_files)))
                self._unpack_tar(tarfile)
                Log.result('Atmospheric file {0}/{1} finished'.format(count, len(tar_files)))
                count += 1
    def _cmorize_grib(self):
        count = 1
        atmos_timestep = None
        chunk_start = parse_date(self.startdate)

        while os.path.exists(os.path.join(self.original_files_path, self._get_grib_filename('GG', chunk_start))) or \
                os.path.exists(os.path.join(self.original_files_path, self._get_grib_filename('SH', chunk_start))):

            chunk_end = chunk_end_date(chunk_start, self.experiment.chunk_size, 'month', 'standard')
            chunk_end = previous_day(chunk_end, 'standard')
            Log.info('CMORizing chunk {0}-{1}', date2str(chunk_start), date2str(chunk_end))
            for grid in ('SH', 'GG'):
                Log.info('Processing {0} variables', grid)

                if not os.path.exists(os.path.join(self.original_files_path,
                                                   self._get_grib_filename(grid, chunk_start))):
                    continue

                for month in range(0, self.experiment.chunk_size):
                    current_month = add_months(chunk_start, month, 'standard')
                    original_gribfile = os.path.join(self.original_files_path,
                                                     self._get_grib_filename(grid, current_month))
                    Log.info('Processing month {1}', grid, date2str(current_month))
                    gribfile = os.path.join(self.config.scratch_dir, os.path.basename(original_gribfile))
                    if not os.path.isfile(gribfile):
                        Log.info('Copying file...', grid, date2str(current_month))
                        shutil.copy(original_gribfile, gribfile)

                    if atmos_timestep is None:
                        atmos_timestep = self._get_atmos_timestep(gribfile)

                    prev_gribfile = os.path.join(self.config.scratch_dir,
                                                 self._get_grib_filename(grid,
                                                                         add_months(current_month, -1, 'standard')))
                    if os.path.exists(prev_gribfile):
                        self._merge_grib_files(current_month, prev_gribfile, gribfile)
                        full_file = 'ICM'
                    else:
                        full_file = gribfile

                    Log.info('Unpacking... ')
                    # remap on regular Gauss grid
                    if grid == 'SH':
                        Utils.cdo.splitparam(input='-sp2gpl {0}'.format(full_file), output=gribfile + '_',
                                             options='-f nc4')
                    else:
                        Utils.cdo.splitparam(input=full_file, output=gribfile + '_', options='-R -f nc4')
                        # total precipitation (remove negative values)
                        Utils.cdo.setcode(228, input='-setmisstoc,0 -setvrange,0,Inf -add '
                                                     '{0}_{{142,143}}.128.nc'.format(gribfile),
                                          output='{0}_228.128.nc'.format(gribfile))
                    Utils.remove_file('ICM')
                    next_gribfile = os.path.join(self.original_files_path,
                                                 self._get_grib_filename(grid,
                                                                         add_months(current_month, 1, 'standard')))

                    if not os.path.exists(next_gribfile):
                        os.remove(gribfile)

                    cdo_reftime = parse_date(self.startdate).strftime('%Y-%m-%d,00:00')

                    self._ungrib_vars(cdo_reftime, gribfile, current_month.month, '{0}hr'.format(atmos_timestep))
                    self._ungrib_vars(cdo_reftime, gribfile, current_month.month, '1d')
                    self._ungrib_vars(cdo_reftime, gribfile, current_month.month, '1m')

                    for splited_file in glob.glob('{0}_*.128.nc'.format(gribfile)):
                        os.remove(splited_file)

                    Log.result('Month {0}, {1} variables finished', date2str(current_month), grid)
                    count += 1

                self._merge_and_cmorize_atmos(chunk_start, chunk_end, grid, '1m')
                self._merge_and_cmorize_atmos(chunk_start, chunk_end, grid, '1d')
                self._merge_and_cmorize_atmos(chunk_start, chunk_end, grid,
                                              '{0}hr'.format(atmos_timestep))
            chunk_start = chunk_end_date(chunk_start, self.experiment.chunk_size, 'month', 'standard')

    def _get_grib_filename(self, grid, month):
        return 'ICM{0}{1}+{2}.grb'.format(grid, self.experiment.expid, date2str(month)[:-2])

    def _get_atmos_timestep(self, gribfile):
        Log.info('Getting timestep...')
        grib_handler = pygrib.open(gribfile)
        mes1 = grib_handler.message(1)
        mes2 = grib_handler.readline()
        while mes2.analDate == mes1.analDate:
            mes2 = grib_handler.readline()
        atmos_timestep = mes2.analDate - mes1.analDate
        atmos_timestep = int(atmos_timestep.total_seconds() / 3600)
        self.experiment.atmos_timestep = atmos_timestep
        grib_handler.close()
        return atmos_timestep

    def _cmorize_nc_file(self, filename):
        Log.info('Processing file {0}', filename)
        temp = TempFile.get()
        handler = Utils.openCdf(filename)
        variables = handler.variables.keys()
        handler.close()
        if not self.cmor.any_required(variables):
            os.remove(filename)
            return
        Utils.execute_shell_command(["nccopy", "-4", "-d4", "-s", filename, temp])
        shutil.move(temp, filename)
        file_parts = os.path.basename(filename).split('_')
        if self.experiment.expid in [file_parts[1], file_parts[2]]:
            frequency = 'm'
        elif self.experiment.expid == file_parts[0]:
            try:
                parse_date(file_parts[1])
                frequency = 'm'
            except ValueError:
                frequency = file_parts[1][1].lower()
        else:
            frequency = file_parts[1][1].lower()
        variables = dict()
        variables['time_counter'] = 'time'
        variables['time_counter_bnds'] = 'time_bnds'
        variables['tbnds'] = 'bnds'
        variables['nav_lat'] = 'lat'
        variables['nav_lon'] = 'lon'
        variables['x'] = 'i'
        variables['y'] = 'j'
        Utils.rename_variables(filename, variables, False, True)
        handler = Utils.openCdf(filename)
        self._add_common_attributes(frequency, handler)
        self._update_time_variables(handler)
        handler.sync()
        temp = TempFile.get()
        Log.info('Splitting file {0}', filename)
        for variable in handler.variables.keys():
            if variable in ('lon', 'lat', 'time', 'time_bnds', 'leadtime', 'lev', 'icethi',
                            'deptht', 'depthu', 'depthw', 'depthv', 'time_centered', 'time_centered_bounds',
                            'deptht_bounds', 'depthu_bounds', 'depthv_bounds', 'depthw_bounds',
                            'time_counter_bounds', 'ncatice',
                            'nav_lat_grid_V', 'nav_lat_grid_U', 'nav_lat_grid_T',
                            'nav_lon_grid_V', 'nav_lon_grid_U', 'nav_lon_grid_T',
                            'depth', 'depth_2', 'depth_3', 'depth_4',
                            'mlev', 'hyai', 'hybi', 'hyam', 'hybm'):
                continue
            self.extract_variable(filename, handler, frequency, temp, variable)
        Log.result('File {0} cmorized!', filename)
        handler.close()
        os.remove(filename)

    def extract_variable(self, file_path, handler, frequency, temp, variable):
        """
        Extracts a variable from a file and creates the CMOR file

        :param file_path: path to the file
        :type file_path: str
        :param handler: netCDF4 handler for the file
        :type handler: netCDF$.Dataset
        :param frequency: variable's frequency
        :type frequency: str
        :param temp: temporal file to use
        :type temp: str
        :param variable: variable's name
        :type variable: str
        """
        file_parts = os.path.basename(file_path).split('_')
        var_cmor = Variable.get_variable(variable)
        if var_cmor is None:
            return
        if not self.cmor.cmorize(var_cmor):
            return
        if frequency == 'd':
            frequency = 'day'
        elif frequency == 'm':
            frequency = 'mon'
        elif frequency == 'h':
            frequency = '6hr'
        else:
            raise Exception('Frequency {0} not supported'.format(frequency))
        Utils.nco.ncks(input=file_path, output=temp, options='-v {0}'.format(variable))
        if var_cmor.domain == 'ocean':
            Utils.rename_variables(temp, {'deptht': 'lev', 'depthu': 'lev', 'depthw': 'lev', 'depthv': 'lev',
                                          'depth': 'lev'}, False, True)
        elif var_cmor.domain in ('land', 'landIce'):
            Utils.rename_variables(temp, {'depth': 'sdepth', 'depth_2': 'sdepth', 'depth_3': 'sdepth',
                                          'depth_4': 'sdepth'}, False, True)
        elif var_cmor.domain == 'atmos':
            Utils.rename_variables(temp, {'depth': 'plev'}, False, True)

        handler_cmor = Utils.openCdf(temp)
        Utils.copy_variable(handler, handler_cmor, 'lon', False)
        Utils.copy_variable(handler, handler_cmor, 'lat', False)
        if 'time' in handler_cmor.dimensions.keys():
            Utils.copy_variable(handler, handler_cmor, 'leadtime', False)
        handler_cmor.close()

        if var_cmor.basin is None:
            region = None
        else:
            region = var_cmor.basin.fullname

        if file_parts[0] == self.experiment.expid or file_parts[0].startswith('ORCA') or \
                file_parts[0] in ('MMA', 'MMO'):
            # Model output
            date_str = '{0}-{1}'.format(file_parts[2][0:6], file_parts[3][0:6])
        elif file_parts[1] == self.experiment.expid:
            # Files generated by the old version of the diagnostics
            date_str = '{0}-{1}'.format(file_parts[4][0:6], file_parts[5][0:6])
        else:
            Log.error('Variable {0} can not be cmorized. Original filename does not match a recognized pattern',
                      var_cmor.short_name)
            return

        self.data_manager.send_file(temp, var_cmor.domain, var_cmor.short_name, self.startdate, self.member,
                                    frequency=frequency, rename_var=variable, date_str=date_str, region=region,

    @staticmethod
    def _merge_grib_files(current_month, prev_gribfile, gribfile):
        Log.info('Merging data from different files...')
        fd = open('rules_files', 'w')
        fd.write('if (dataDate >= {0.year}{0.month:02}01) {{ write ; }}\n'.format(current_month))
        fd.close()
        # get first timestep for each month from previous file (if possible)
        if os.path.exists('ICM'):
            os.remove('ICM')
        Utils.execute_shell_command('grib_filter -o ICM rules_files '
                                    '{0} {1}'.format(os.path.basename(prev_gribfile),
                                                     os.path.basename(gribfile)))
        os.remove('rules_files')
        Utils.remove_file(prev_gribfile)

    def _ungrib_vars(self, cdo_reftime, gribfile, month, frequency):
        Log.info('Preparing {0} variables'.format(frequency))
        var_codes = self.config.cmor.get_variables(frequency)
        for var_code in var_codes:
            if not os.path.exists('{0}_{1}.128.nc'.format(gribfile, var_code)):
                continue
            new_units = None

            cdo_operator = '-selmon,{0}'.format(month)
            if frequency in ('month', 'monthly', 'mon', '1m'):
                if var_code == 201:
                    cdo_operator = "-monmean -daymax {0}".format(cdo_operator)
                elif var_code == 202:
                    cdo_operator = "-monmean -daymax {0}".format(cdo_operator)
                else:
                    cdo_operator = "-monmean {0} ".format(cdo_operator)
            elif frequency in ('day', 'daily', '1d'):
                if var_code == 201:
                    cdo_operator = "-daymax {0} ".format(cdo_operator)
                elif var_code == 202:
                    cdo_operator = "-daymin {0} ".format(cdo_operator)
                else:
                    cdo_operator = "-daymean {0} ".format(cdo_operator)

            if var_code in (144, 146, 147, 169, 175, 176, 177, 179, 180, 181, 182, 201, 202, 205, 212, 228):
                cdo_operator = '{0} -shifttime,-{1}hours'.format(cdo_operator, self.experiment.atmos_timestep)

            if var_code == 129:
                # geopotential
                new_units = "m"
                cdo_operator = "-divc,9.81 {0}".format(cdo_operator)
            elif var_code in (146, 147, 169, 175, 176, 177, 179, 212):
                # radiation
                new_units = "W m-2"
                cdo_operator = "-divc,{0} {1}".format(self.experiment.atmos_timestep * 3600, cdo_operator)
            elif var_code in (180, 181):
                # momentum flux
                new_units = "N m-2"
                cdo_operator = "-divc,{0} {1}".format(self.experiment.atmos_timestep * 3600, cdo_operator)
            elif var_code in (144, 182, 205, 228):
                # precipitation/evaporation/runoff
                new_units = "kg m-2 s-1"
                cdo_operator = "-mulc,1000 -divc,{0}".format(self.experiment.atmos_timestep * 3600)

            levels = self.config.cmor.get_levels(frequency, var_code)
            if levels:
                cdo_operator = "{0} -sellevel,{1}".format(cdo_operator, levels)

            Utils.execute_shell_command('cdo -t ecmwf setreftime,{0} '
                                        '{1} {2}_{3}.128.nc '
                                        '{2}_{3}_{4}.nc'.format(cdo_reftime, cdo_operator,
                                                                gribfile, var_code, frequency))
            h_var_file = '{0}_{1}_{2}.nc'.format(gribfile, var_code, frequency)

            handler = Utils.openCdf(h_var_file)
            if new_units:
                for var in handler.variables.values():
                    if 'code' in var.ncattrs() and var.code == var_code:
                        var.units = new_units
                        break

            var_name = None
            for key in handler.variables.keys():
                if key + '_2' in handler.variables and key not in handler.dimensions:
                    var_name = key
            handler.close()

            if var_name is not None:
                Utils.nco.ncks(input='{0}_{1}_1m.nc'.format(gribfile, var_code),
                               output='{0}_{1}_1m.nc'.format(gribfile, var_code),
                               options='-O -v {0}'.format(var_name))

    def _merge_and_cmorize_atmos(self, chunk_start, chunk_end, grid, frequency):
        merged_file = 'MMA_{0}_{1}_{2}_{3}.nc'.format(frequency, date2str(chunk_start), date2str(chunk_end), grid)
        files = glob.glob(os.path.join(self.config.scratch_dir,
                                       '{0}_*_{1}.nc'.format(self._get_grib_filename(grid, chunk_start), frequency)))
        for first_file in files:
            shutil.move(first_file, merged_file)
            current_month = add_months(chunk_start, 1, 'standard')
            while current_month < chunk_end:
                month_file = first_file.replace('+{0}.grb'.format(date2str(chunk_start)[:-2]),
                                                '+{0}.grb'.format(date2str(current_month)[:-2]))
                Utils.concat_variables(month_file, merged_file, True)
                current_month = add_months(current_month, 1, 'standard')

            self._cmorize_nc_file(merged_file)

    def _update_time_variables(self, handler):
        time_var = handler.variables['time']
        times = Utils.get_datetime_from_netcdf(handler)
        if type(times[0]) is not datetime:
            for x in range(0, times.shape[0]):
                # noinspection PyProtectedMember
                times[x] = times[x]._to_real_datetime()
        time_var[:] = netCDF4.date2num(times, 'days since 1850-01-01', 'standard')
        if 'axis_nbounds' in handler.dimensions:
            handler.renameDimension('axis_nbounds', 'bnds')

        if 'time_counter_bounds' in handler.variables:
            handler.renameVariable('time_counter_bounds', 'time_bnds')
            handler.sync()
        if 'time_bnds' in handler.variables:
            time_bounds_var = handler.variables['time_bnds']
            time_var.bounds = "time_bnds"
            try:
                units = time_bounds_var.units
            except AttributeError:
                # we suppose that if not units are present, they match the ones in the time variable
                units = time_var.units
            time_bounds_var.units = units

            time_bounds = Utils.get_datetime_from_netcdf(handler, 'time_bnds')
            if type(time_bounds[0, 0]) is not datetime:
                for x in range(0, time_bounds.shape[0]):
                    for y in range(0, time_bounds.shape[1]):
                        # noinspection PyProtectedMember
                        time_bounds[x, y] = time_bounds[x, y]._to_real_datetime()
            time_bounds_var[:] = netCDF4.date2num(time_bounds, 'days since 1850-01-01', 'standard')
        time_var.units = 'days since 1850-01-01'
        time_var.time_origin = "1850-01-01"
        time_var.calendar = 'standard'
        time_var.long_name = "Verification time of the forecast"
        time_var.standard_name = "time"
        time_var.axis = "T"
        if 'leadtime' in handler.variables:
            var = handler.variables['leadtime']
        else:
            var = handler.createVariable('leadtime', float, 'time')
        var.units = "days"
        var.long_name = "Time elapsed since the start of the forecast"
        var.standard_name = "forecast_period"
        leadtime = (Utils.get_datetime_from_netcdf(handler) - parse_date(self.startdate))
        for lt in range(0, leadtime.shape[0]):
            var[lt] = leadtime[lt].days

    def _add_common_attributes(self, frequency, handler):
        cmor = self.config.cmor
        experiment = self.config.experiment
        handler.associated_experiment = cmor.associated_experiment
        handler.batch = '{0}{1}'.format(experiment.institute, datetime.now().strftime('%Y-%m-%d(T%H:%M:%SZ)'))
        handler.contact = 'Pierre-Antoine Bretonnière, pierre-antoine.bretonniere@bsc.es , ' \
                          'Javier Vegas-Regidor, javier.vegas@bsc.es '
        handler.Conventions = 'CF-1.6'
        handler.creation_date = datetime.now().strftime('%Y-%m-%d(T%H:%M:%SZ)')
        handler.experiment_id = experiment.experiment_name
        handler.forecast_reference_time = parse_date(self.startdate).strftime('%Y-%m-%d(T%H:%M:%SZ)')
        if frequency == 'd':
            handler.frequency = 'day'
        elif frequency == 'm':
            handler.frequency = 'mon'
        handler.institute_id = experiment.institute
        handler.institution = experiment.institute
        handler.initialization_method = cmor.initialization_method
        handler.initialization_description = cmor.initialization_description
        handler.physics_version = cmor.physics_version
        handler.physics_description = cmor.physics_description
        handler.model_id = experiment.model
        handler.associated_model = cmor.associated_model
        handler.project_id = 'SPECS'
        handler.realization = str(self.member + 1)
        handler.source = cmor.source
        handler.startdate = 'S{0}'.format(self.startdate)
        handler.tracking_id = str(uuid.uuid1())
        handler.title = "{0} model output prepared for SPECS {1}".format(experiment.model, experiment.experiment_name)