:py:mod:`opendrift.models.openoil.openoil` ========================================== .. py:module:: opendrift.models.openoil.openoil .. autoapi-nested-parse:: OpenOil is a 3D oil drift module bundled within the OpenDrift framework. The oil weathering calculations is based on the NOAA-ERR-ERD OilLibrary package, which is installed as a dependency. The code for evaporation and emulsification in OpenOil is borrowed from the NOAA PyGnome code, and adapted to the OpenDrift architecture. Example of ship leaking oil along the coast of Northern Norway ############################################################## .. image:: https://dl.dropboxusercontent.com/s/ty6dmqf0oohewky/oilspill_tromsoe.gif?dl=0 Simulation of Deepwater Horizon (Macondo) accident, initiated from satellite images ################################################################################### .. image:: https://dl.dropboxusercontent.com/s/ghi7crtmwpyjgto/macondo_simulation.gif?dl=0 Satellite images provided by Prof. Chuanmin Hu, and ocean model output provided by Prof. Vassiliki Kourafalou Example oil budget for a simulation ################################### .. image:: https://dl.dropboxusercontent.com/s/pb0h6tlev9pnoh3/oil_budget_draugen.png?dl=0 Oil properties affecting the drift *********************************** The vertical (and thus indirectly also the horizontal) motion of oil (droplets) is affected by oil density and droplet diameters. When using the NOAA oil weathering model (``o = OpenOil(weathering_model='noaa')``), which is the default, the density is obtained from the NOAA database according to the oiltype selected when seeding. This value can not be overridden by the user, and it will also change during the simulation due to oil weathering processes (evaporation and emulsification). The droplet diameter may be given explicitly when seeding, e.g.: .. testcode:: o = OpenOil() o.set_config('environment:constant:x_wind', 0) o.set_config('environment:constant:y_wind', 0) o.set_config('environment:constant:x_sea_water_velocity', 0) o.set_config('environment:constant:y_sea_water_velocity', 0) o.seed_elements(4, 60, number=100, time=datetime.now(), diameter=1e-5) In this case, the diameter will not change during the simulation, which is useful e.g. for sensitivity tests. The same diameter will be used for all elements for this example, but an array of the same length as the number of elements may also be provided. If a constant droplet diameter is not given by the user, it will be chosen randomly within given config limits for a subsea spill ('blowout'), and modified after any later wave breaking event. Oil droplets seeded under sea surface (z<0) will be assigned initial diameters between the following limits, typical for a subsea blowout (Johansen, 2000):: .. code:: o.set_config('seed:droplet_diameter_min_subsea', 0.0005) # 0.5 mm o.set_config('seed:droplet_diameter_max_subsea', 0.005) # 5 mm Alternatively, the user can specify normal or lognormal initial subsea droplet size distributions, which are later modified by wave breaking events. In these cases the user must specify the mean and standard deviation of the distribution:: .. code:: o.set_config('seed:droplet_size_distribution','lognormal') o.set_config('seed:droplet_diameter_mu',0.001) # 1 mm o.set_config('seed:droplet_diameter_sigma',0.0008) # 0.8 mm Note that these config settings must be adjusted before the seeding call. After each wave breaking event, a new droplet diameter will be chosen based on the config setting for droplet size distribution. Module Contents --------------- Classes ~~~~~~~ .. autoapisummary:: opendrift.models.openoil.openoil.Oil opendrift.models.openoil.openoil.OpenOil Attributes ~~~~~~~~~~ .. autoapisummary:: opendrift.models.openoil.openoil.logger .. py:data:: logger .. py:class:: Oil(**kwargs) Bases: :py:obj:`opendrift.models.oceandrift.Lagrangian3DArray` Extending LagrangianArray with variables relevant for oil particles. Initialises a LagrangianArray with given properties. Args: Keyword arguments (kwargs) with names corresponding to the OrderedDict 'variables' of the class, and corresponding values. The values must be ndarrays of equal length, or scalars. All (or none) variables must be given, unless a default value is specified in the OrderedDict 'variables' An empty object may be created by giving no input. .. py:attribute:: variables .. py:class:: OpenOil(weathering_model='noaa', *args, **kwargs) Bases: :py:obj:`opendrift.models.oceandrift.OceanDrift` Open source oil trajectory model based on the OpenDrift framework. Developed at MET Norway based on oil weathering parameterisations found in open/published litterature. Under construction. Initialise OpenDriftSimulation Args: seed: integer or None. A given integer will yield identical random numbers drawn each simulation. Random numbers are e.g. used to distribute particles spatially when seeding, and may be used by modules (subclasses) for e.g. diffusion. Specifying a fixed value (default: 0) is useful for sensitivity tests. With seed = None, different random numbers will be drawn for subsequent runs, even with identical configuration/input. iomodule: name of module used to export data default: netcdf, see :py:mod:`opendrift.io` for more alternatives. `iomodule` is module/filename without preceeding `io_` loglevel: set to 0 (default) to retrieve all debug information. Provide a higher value (e.g. 20) to receive less output. Use the string 'custom' to configure logging from outside. logtime: if True, a time stamp is given for each logging line. logtime can also be given as a python time specifier (e.g. '%H:%M:%S') .. py:attribute:: ElementType .. py:attribute:: required_variables .. py:attribute:: status_colors .. py:attribute:: duplicate_oils :value: ['ALVHEIM BLEND, STATOIL', 'DRAUGEN, STATOIL', 'EKOFISK BLEND 2000', 'EKOFISK BLEND, STATOIL',... .. py:method:: update_surface_oilfilm_thickness() The mass of oil is summed within a grid of 100x100 cells covering the oil at a given time. Each oil particle within each cell is given a film thickness as the amount of oil divided by the cell area. .. py:method:: biodegradation() .. py:method:: biodegradation_half_time() Oil biodegradation with exponential decay .. py:method:: biodegradation_adcroft() Oil biodegradation function based on the article: Adcroft et al. (2010), Simulations of underwater plumes of dissolved oil in the Gulf of Mexico. .. py:method:: disperse() .. py:method:: oil_weathering() .. py:method:: prepare_run() .. py:method:: oil_weathering_noaa() Oil weathering scheme adopted from NOAA PyGNOME model: https://github.com/NOAA-ORR-ERD/PyGnome .. py:method:: disperse_noaa() .. py:method:: plot_droplet_spectrum() Plotting distribution of droplet radii, for debugging .. py:method:: evaporation_noaa() .. py:method:: emulsification_noaa() .. py:method:: update_terminal_velocity(Tprofiles=None, Sprofiles=None, z_index=None) Calculate terminal velocity for oil droplets according to: * Tkalich et al. (2002): Vertical mixing of oil droplets by breaking waves * Marine Pollution Bulletin 44, 1219-1229 If profiles of temperature and salt are passed into this function, they will be interpolated from the profiles. if not, T,S will be fetched from reader. .. py:method:: oil_wave_entrainment_rate() .. py:method:: prepare_vertical_mixing() Calculate entrainment probability before main loop .. py:method:: surface_wave_mixing(time_step_seconds) Mix surface oil into water column. .. py:method:: surface_stick() set surfaced particles to exactly zero depth to let them form a slick .. py:method:: get_wave_breaking_droplet_diameter() .. py:method:: get_wave_breaking_droplet_diameter_liz2017() .. py:method:: get_wave_breaking_droplet_diameter_johansen2015() .. py:method:: resurface_elements(minimum_depth=None) Oil elements reaching surface (or above) form slick, not droplet .. py:method:: advect_oil() .. py:method:: update() Update positions and properties of oil particles. .. py:method:: get_oil_budget() Get oil budget for the current simulation The oil budget consists of the following categories: * surface: the sum of variable mass_oil for all active elements where z = 0 * submerged: the sum of variable mass_oil for all active elements where z < 0 * stranded: the sum of variable mass_oil for all elements which are stranded * evaporated: the sum of variable mass_evaporated for all elements * dispersed: the sum of variable mass_dispersed for all elements The sum (total mass) shall equal the mass released. Note that the mass of oil is conserved, whereas the volume may change with changes in density and water uptake (emulsification). Therefore mass should be used for budgets, eventually converted to volume (by dividing on density) in the final step before presentation. Note that mass_oil is the mass of pure oil. The mass of oil emulsion (oil containing entrained water droplets) can be calculated as: .. code:: mass_emulsion = mass_oil / (1 - water_fraction) I.e. water_fraction = 0 means pure oil, water_fraction = 0.5 means mixture of 50% oil and 50% water, and water_fraction = 0.9 (which is maximum) means 10% oil and 90% water. .. py:method:: plot_oil_budget(filename=None, ax=None, show_watercontent_and_viscosity=True, show_wind_and_current=True) .. py:method:: get_oil_name() .. py:method:: cumulative_oil_entrainment_fraction() Returns the fraction of oil elements which has been entrained vs time .. py:method:: plot_oil_watercontent_and_viscosity(ax=None, show=True) .. py:method:: set_oiltype(oiltype) Sets the oil type by specifying the name, the first match will be chosen. See the `ADIOS database `_ for a list. OpenDrift provides a small set of extra oils. .. py:method:: set_oiltype_by_id(oiltypeid) Sets the oil type by specifying the ADIOS ID. See the `ADIOS database `_ for a list. OpenDrift provides a small set of extra oils. .. py:method:: set_oiltype_by_json(json) Sets the oil type by specifing a JSON dict. The format should be the same as the ADIOS database. See the `ADIOS database `_ for a list. .. py:method:: set_oiltype_from_file(path) Sets the oil type by specifing a JSON file. The format should be the same as the ADIOS database. See the `ADIOS database `_ for a list. >>> o = OpenOil() >>> o.set_oiltype_from_file('opendrift/models/openoil/adios/extra_oils/AD04001.json') .. py:method:: store_oil_seed_metadata(**kwargs) .. py:method:: seed_elements(*args, **kwargs) Seed elements with given position(s), time and properties. Arguments: lon: scalar or array central longitude(s). lat: scalar or array central latitude(s). radius: scalar or array radius in meters around each lon-lat pair, within which particles will be randomly seeded. number: integer, total number of particles to be seeded If number is None, the number of elements is the length of lon/lat or time if these are arrays. Otherwise the number of elements are obtained from the config-default. time: datenum or list The time at which particles are seeded/released. If time is a list with two elements, elements are seeded continously from start/first to end/last time. If time is a list with more than two elements, the number of elements is equal to len(time) and are seeded as a time series. radius_type: string If 'gaussian' (default), the radius is the standard deviation in x-y-directions. If 'uniform', elements are spread evenly and always inside a circle with the given radius. kwargs: keyword arguments containing properties/attributes and values corresponding to the actual particle type (ElementType). These are forwarded to the ElementType class. All properties for which there are no default value must be specified. .. py:method:: seed_cone(*args, **kwargs) Seed elements along a transect/cone between two points/times Arguments: lon: scalar or list with 2 elements [lon0, lon1] lat: scalar or list with 2 elements [lat0, lat] time: datetime or list with 2 elements [t0, t1] radius: scalar or list with 2 elements [r0, r1] Unit: meters number (int): The number of elements. If this is None, the number of elements is taken from configuration. Elements are seeded along a transect from (lon0, lat0) with uncertainty radius r0 at time t0, towards (lon1, lat1) with uncertainty radius r1 at time t1. If r0 != r1, the unceetainty radius is linearly changed along the transect, thus outlining a "cone". .. py:method:: seed_from_gml(gmlfile, num_elements=1000, *args, **kwargs) Read oil slick contours from GML file, and seed particles within. .. py:method:: seed_from_geotiff_thickness(filename, number=50000, *args, **kwargs) Seed from files as provided by Prof. Chuanmin Hu .. py:method:: _substance_name()