# coding=utf-8 """Compute the MOC for oceanic basins""" import numpy as np import six from bscearth.utils.log import Log import iris from iris.coords import DimCoord, AuxCoord from iris.cube import CubeList import iris.analysis from earthdiagnostics.constants import Basins from earthdiagnostics.diagnostic import Diagnostic, DiagnosticBasinListOption from earthdiagnostics.modelingrealm import ModelingRealms from earthdiagnostics.utils import Utils, TempFile class Moc(Diagnostic): """ Compute the MOC for oceanic basins :original author: Virginie Guemas :contributor: Javier Vegas-Regidor :created: March 2012 :last modified: June 2016 :param data_manager: data management object :type data_manager: DataManager :param startdate: startdate :type startdate: str :param member: member number :type member: int :param chunk: chunk's number :type chunk: int """ alias = 'moc' "Diagnostic alias for the configuration file" vsftmyz = 'vsftmyz' def __init__(self, data_manager, startdate, member, chunk, basins): Diagnostic.__init__(self, data_manager) self.startdate = startdate self.member = member self.chunk = chunk self.basins = basins def __str__(self): basins = [] basins.extend(self.basins.keys()) return 'MOC Startdate: {0.startdate} Member: {0.member} ' \ 'Chunk: {0.chunk} Basins: {1}'.format(self, basins) def __hash__(self): return hash(str(self)) def __eq__(self, other): if self._different_type(other): return False return self.startdate == other.startdate and self.member == other.member and self.chunk == other.chunk @classmethod def generate_jobs(cls, diags, options): """ Create a job for each chunk to compute the diagnostic :param diags: Diagnostics manager class :type diags: Diags :param options: None :type options: list[str] :return: """ basins = Basins() options_available = ( DiagnosticBasinListOption( 'basins', 'glob' ), ) options = cls.process_options(options, options_available) basins = options['basins'] if not basins: Log.error('Basins not recognized') return () try: e1v = iris.load_cube('mesh_hgr.nc', 'e1v') except iris.exceptions.ConstraintMismatchError: e1v = iris.load_cube('mesh_hgr.nc', 'e1v_0') try: e3v = iris.load_cube('mesh_hgr.nc', 'e3t_0') except iris.exceptions.ConstraintMismatchError: e3v = iris.load_cube('mesh_hgr.nc', 'e3t_0') e1v = iris.util.squeeze(e1v).data e3v = iris.util.squeeze(e3v).data if len(e3v.shape) == 1: e3v = np.expand_dims(e3v.data, 1) e3v = np.expand_dims(e3v, 2) else: e3v = e3v.data mesh = - e1v * e3v masks = {} basins.sort() for basin in basins: masks[basin.name] = Utils.get_mask(basin) * mesh / 1e6 job_list = list() for startdate, member, chunk in diags.config.experiment.get_chunk_list(): job_list.append(Moc(diags.data_manager, startdate, member, chunk, masks)) return job_list def request_data(self): """Request data required by the diagnostic""" self.variable_file = self.request_chunk(ModelingRealms.ocean, 'vo', self.startdate, self.member, self.chunk) def declare_data_generated(self): """Declare data to be generated by the diagnostic""" self.results = self.declare_chunk(ModelingRealms.ocean, Moc.vsftmyz, self.startdate, self.member, self.chunk) def compute(self): """Run the diagnostic""" data = iris.load_cube(self.variable_file.local_file) Log.debug(str(data)) try: data.coord('i') except iris.exceptions.CoordinateNotFoundError: dims = len(data.shape) data.add_dim_coord(iris.coords.DimCoord(np.arange(data.shape[dims - 1]), var_name='i'), dims - 1) try: data.coord('j') except iris.exceptions.CoordinateNotFoundError: dims = len(data.shape) data.add_dim_coord(iris.coords.DimCoord(np.arange(data.shape[dims - 2]), var_name='j'), dims - 2) for coord_name in ('model_level_number', 'Vertical V levels', 'lev'): if data.coords(coord_name): coord = data.coord(coord_name) coord.standard_name = 'depth' coord.long_name = 'depth' break moc_results = CubeList() for map_slice in data.slices_over('time'): # Force data loading map_slice.data Log.debug(str(map_slice)) for basin, mask in six.iteritems(self.basins): moc = map_slice.collapsed(('i', 'depth'), iris.analysis.SUM, weights=mask) moc.add_aux_coord( AuxCoord([basin], var_name='region', standard_name='region') ) moc_results.append(moc) results = moc_results.merge_cube() temp = TempFile.get() results.var_name = Moc.vsftmyz results.remove_coord('i') results.remove_coord('depth') results.remove_coord('longitude') results.units = 'Sverdrup' iris.save(results, temp) self.results.set_local_file(temp)