Diagnostics
We describe here the modules used to analyze the output of the model with diagnostics. Please read the User guide to learn how to use them.
Diagnostic base classes
Abstract base classes defining the diagnostics of the model and used to analyze its outputs.
Description of the classes
Diagnostic
: General base class.FieldDiagnostic
: General base class for diagnostics returning model’s fields.
Warning
These are abstract base class, they must be subclassed to create new diagnostics!
- class qgs.diagnostics.base.Diagnostic(model_params, dimensional)[source]
Bases:
ABC
General base class to create diagnostics.
- Parameters
- set_data(time, data)[source]
Provide the model data to the diagnostic.
- Parameters
time (ndarray) – The time (in nondimensional timeunits) corresponding to the data. Its length should match the length of the last axis of the provided data.
data (ndarray) – The model output data that the user want to convert using the diagnostic. Should be a 2D array of shape (
ndim
, number_of_timesteps).
- class qgs.diagnostics.base.FieldDiagnostic(model_params, dimensional)[source]
Bases:
Diagnostic
General base class for field diagnostic on the model’s domain. Should provide a spatial gridded representation of the fields.
- Parameters
- animate(output='animate', style='image', ax=None, figsize=(16, 9), contour_labels=True, color_bar=True, show_time=True, stride=1, plot_kwargs=None, oro_kwargs=None, anim_kwargs=None, show=True)[source]
Return the output of the plot method animated over time.
- Parameters
output (str, optional) –
Define the kind of animation being created. Can be:
animate: Create and show a
ipywidgets.widgets.interactive
widget. Works only in Jupyter notebooks.show: Create and show an animation with the
matplotlib.animation
module. Works only in IPython or Python.
style (str, optional) –
The style of the plot. Can be:
image: show the fields as images with a given colormap specified in the plot_kwargs argument.
contour: show the fields as contour superimposed on the image of the orographic height (see the oro_kwargs below).
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
contour_labels (bool, optional) – If style is set to contour, specify if the contours must be labelled with their value or not. Default to True.
color_bar (bool, optional) – Specify if a color bar must be drawn beside the plot or not. Default to True.
show_time (bool, optional) – Show the timestamp of the field on the plot or not. Default to True.
stride (int, optional) – Specify the time step of the animation. Works only with output set to animate.
plot_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.oro_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method plotting the image of the orography if style is set to contour.anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation. Works only with output set to show.show (bool, optional) – Whether to plot or not the animation.
- Returns
The animation object or the callable to update the widget, depending on the value of the output and show parameters.
- Return type
FuncAnimation or DisplayHandle or callable
- property grid_shape
Return the shape of the grid of points covering the model’s domain.
- movie(output='html', filename='', style='image', ax=None, figsize=(16, 9), contour_labels=True, color_bar=True, show_time=True, plot_kwargs=None, oro_kwargs=None, anim_kwargs=None)[source]
Create and return a movie of the output of the plot method animated over time.
- Parameters
output (str, optional) –
Define the kind of movie being created. Can be:
jshtml: Generate an interactive HTML representation of the animation.
html5: Generate the movie as HTML5 code.
html: Output the movie as a HTML video tag.
ihtml: Output the interactive movie as a HTML video tag.
save: Save the movie in MP4 format (H264 codec).
Default to html.
filename (str, optional) – Filename (and path) where to save the movie. Needed if output is set to save.
style (str, optional) –
The style of the plot. Can be:
image: show the fields as images with a given colormap specified in the plot_kwargs argument.
contour: show the fields as contour superimposed on the image of the orographic height (see the oro_kwargs below).
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
contour_labels (bool, optional) – If style is set to contour, specify if the contours must be labelled with their value or not. Default to True.
color_bar (bool, optional) – Specify if a color bar must be drawn beside the plot or not. Default to True.
show_time (bool, optional) – Show the timestamp of the field on the plot or not. Default to True.
plot_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.oro_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method plotting the image of the orography if style is set to contour.anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation.
- Returns
The animation object or the HTML code or tag.
- Return type
FuncAnimation or HTML code or HTML tag
- plot(time_index, style='image', ax=None, figsize=(16, 9), contour_labels=True, color_bar=True, show_time=True, plot_kwargs=None, oro_kwargs=None)[source]
Plot the field of the provided data at the given time index.
- Parameters
time_index (int) – The time index of the data.
style (str, optional) –
The style of the plot. Can be:
image: show the fields as images with a given colormap specified in the plot_kwargs argument.
contour: show the fields as contour superimposed on the image of the orographic height (if it exists, see the oro_kwargs below).
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
contour_labels (bool, optional) – If style is set to contour, specify if the contours must be labelled with their value or not. Default to True.
color_bar (bool, optional) – Specify if a color bar must be drawn beside the plot or not. Default to True.
show_time (bool, optional) – Show the timestamp of the field on the plot or not. Default to True.
plot_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.oro_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method plotting the image of the orography if style is set to contour.
- Returns
An axes where the data were plotted.
- Return type
- plot_grid_point(i, j, ax=None, figsize=(16, 9), plot_kwargs=None)[source]
Plot the time serie of the field at a given grid point.
- plot_kwargs
Arguments to pass to the
matplotlib.axes.Axes.plot()
method.- Type
dict, optional
- Returns
An axes where the data were plotted.
- Return type
- class qgs.diagnostics.base.FieldPointDiagnostic(model_params, dimensional)[source]
Bases:
Diagnostic
General base class to give field values over time at a given point of the domain.
Warning
Not yet implemented.
- class qgs.diagnostics.base.ProfileDiagnostic(model_params, dimensional)[source]
Bases:
Diagnostic
General base class for profile diagnostic on the model’s domain. Should provide a 1D representation of the fields averages or section.
- Parameters
- animate(output='animate', ax=None, figsize=(16, 9), show_time=True, stride=1, plot_kwargs=None, anim_kwargs=None, show=True)[source]
Return the output of the plot method animated over time.
- Parameters
output (str, optional) –
Define the kind of animation being created. Can be:
animate: Create and show a
ipywidgets.widgets.interactive
widget. Works only in Jupyter notebooks.show: Create and show an animation with the
matplotlib.animation
module. Works only in IPython or Python.
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
show_time (bool, optional) – Show the timestamp of the field on the plot or not. Default to True.
stride (int, optional) – Specify the time step of the animation. Works only with output set to animate.
plot_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation. Works only with output set to show.show (bool, optional) – Whether to plot or not the animation.
- Returns
The animation object or the callable to update the widget, depending on the value of the output and show parameters.
- Return type
FuncAnimation or DisplayHandle or callable
- movie(output='html', filename='', ax=None, figsize=(16, 9), plot_kwargs=None, anim_kwargs=None)[source]
Create and return a movie of the output of the plot method animated over time. Show a red dot moving and depicting the current value of the model’s selected variables.
- Parameters
output (str, optional) –
Define the kind of movie being created. Can be:
jshtml: Generate an interactive HTML representation of the animation.
html5: Generate the movie as HTML5 code.
html: Output the movie as a HTML video tag.
ihtml: Output the interactive movie as a HTML video tag.
save: Save the movie in MP4 format (H264 codec).
Default to html.
filename (str, optional) – Filename (and path) where to save the movie. Needed if output is set to save.
ax (Axes, optional) – An axes on which to plot. If not provided, create a new one.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
plot_kwargs (dict, optional) – Arguments to pass to the background plot.
anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation.
- Returns
The animation object or the HTML code or tag.
- Return type
FuncAnimation or HTML code or HTML tag
- plot(time_index=0, ax=None, figsize=(16, 9), show_time=True, plot_kwargs=None, **kwargs)[source]
Plot the multiple profile diagnostic provided.
- Parameters
time_index (int) – The time index of the data. Not used in this subclass.
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
show_time (bool, optional) – Show the timestamp of the field on the plot or not. Default to True.
plot_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.
- Returns
An axes where the data were plotted.
- Return type
Diagnostic streamfunction classes
Classes defining streamfunction fields diagnostics.
Description of the classes
AtmosphericStreamfunctionDiagnostic
: General base class for atmospheric streamfunction fields diagnostic.LowerLayerAtmosphericStreamfunctionDiagnostic
: Diagnostic giving the lower layer atmospheric streamfunction fields \(\psi^3_{\rm a}\).UpperLayerAtmosphericStreamfunctionDiagnostic
: Diagnostic giving the upper layer atmospheric streamfunction fields \(\psi^1_{\rm a}\).MiddleAtmosphericStreamfunctionDiagnostic
: Diagnostic giving the middle atmospheric streamfunction fields \(\psi_{\rm a}\).OceanicLayerStreamfunctionDiagnostic
: Diagnostic giving the oceanic layer streamfunction fields \(\psi_{\rm o}\).
- class qgs.diagnostics.streamfunctions.AtmosphericStreamfunctionDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
FieldDiagnostic
General base class for atmospheric streamfunction fields diagnostic. Provide a spatial gridded representation of the fields. This is an abstract base class, it must be subclassed to create new diagnostics!
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool) – Indicate if the output diagnostic must be dimensionalized or not.
- class qgs.diagnostics.streamfunctions.LowerLayerAtmosphericStreamfunctionDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
AtmosphericStreamfunctionDiagnostic
Diagnostic giving the lower layer atmospheric streamfunction fields \(\psi^3_{\rm a}\). Computed as \(\psi^3_{\rm a} = \psi_{\rm a} - \theta_{\rm a}\) where \(\psi_{\rm a}\) and \(\theta_{\rm a}\) are respectively the barotropic and baroclinic streamfunctions. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.streamfunctions.MiddleAtmosphericStreamfunctionDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True, geopotential=False)[source]
Bases:
AtmosphericStreamfunctionDiagnostic
Diagnostic giving the middle atmospheric streamfunction fields \(\psi_{\rm a}\) at 500hPa, i.e. the barotropic streamfunction of the system. See also Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
geopotential (bool, optional) – Dimensionalize the field in geopotential height (in meter). Default to False.
- class qgs.diagnostics.streamfunctions.OceanicLayerStreamfunctionDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True, conserved=True)[source]
Bases:
OceanicStreamfunctionDiagnostic
Diagnostic giving the oceanic layer streamfunction fields \(\psi_{\rm o}\).
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
conserved (bool, optional) – Whether to plot the conserved oceanic fields or not. Default to True.
- class qgs.diagnostics.streamfunctions.OceanicStreamfunctionDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
FieldDiagnostic
General base class for oceanic streamfunction fields diagnostic. Provide a spatial gridded representation of the fields.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool) – Indicate if the output diagnostic must be dimensionalized or not.
- class qgs.diagnostics.streamfunctions.UpperLayerAtmosphericStreamfunctionDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
AtmosphericStreamfunctionDiagnostic
Diagnostic giving the upper layer atmospheric streamfunction fields \(\psi^1_{\rm a}\). Computed as \(\psi^1_{\rm a} = \psi_{\rm a} + \theta_{\rm a}\) where \(\psi_{\rm a}\) and \(\theta_{\rm a}\) are respectively the barotropic and baroclinic streamfunctions. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
Diagnostic temperature classes
Classes defining temperature fields diagnostics.
Description of the classes
AtmosphericTemperatureDiagnostic
: General base class for atmospheric temperature fields diagnostics.MiddleAtmosphericTemperatureDiagnostic
: Diagnostic giving the middle atmospheric anomaly fields \(T_{\rm a}\).MiddleAtmosphericTemperatureAnomalyDiagnostic
: Diagnostic giving the middle atmospheric temperature anomaly fields \(\delta T_{\rm a}\).OceanicTemperatureDiagnostic
: General base class for oceanic temperature fields diagnostic.OceanicLayerTemperatureDiagnostic
: Diagnostic giving the oceanic layer temperature fields \(T_{\rm o}\).OceanicLayerTemperatureAnomalyDiagnostic
: Diagnostic giving the oceanic layer temperature anomaly fields \(\delta T_{\rm o}\).GroundTemperatureDiagnostic
: Diagnostic giving the ground layer temperature fields \(T_{\rm g}\).GroundTemperatureAnomalyDiagnostic
: Diagnostic giving the ground layer temperature anomaly fields \(\delta T_{\rm g}\).AtmosphericTemperatureMeridionalGradientDiagnostic
: General base class for meridional gradient of atmospheric temperature fields diagnostics.MiddleAtmosphericTemperatureMeridionalGradientDiagnostic
: Diagnostic giving the meridional gradient of the middle atmospheric temperature fields \(\partial_y T_{\rm a}\).
- class qgs.diagnostics.temperatures.AtmosphericTemperatureDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
FieldDiagnostic
General base class for atmospheric temperature fields diagnostic. Provide a spatial gridded representation of the fields. This is an abstract base class, it must be subclassed to create new diagnostics!
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool) – Indicate if the output diagnostic must be dimensionalized or not.
- class qgs.diagnostics.temperatures.AtmosphericTemperatureMeridionalGradientDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
DifferentialFieldDiagnostic
General base class for atmospheric temperature fields meridional gradient diagnostic. Provide a spatial gridded representation of the fields. This is an abstract base class, it must be subclassed to create new diagnostics!
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool) – Indicate if the output diagnostic must be dimensionalized or not.
- class qgs.diagnostics.temperatures.GroundTemperatureAnomalyDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
FieldDiagnostic
Diagnostic giving the ground temperature anomaly fields \(\delta T_{\rm g}\).
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.temperatures.GroundTemperatureDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
GroundTemperatureAnomalyDiagnostic
Diagnostic giving the ground temperature fields \(T_{\rm g} = T_{{\rm g}, 0} + \delta T_{\rm g}\), where \(T_{{\rm g}, 0}\) is the reference temperature
T0
or the 0-th order dynamic temperature.- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.temperatures.MiddleAtmosphericTemperatureAnomalyDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
AtmosphericTemperatureDiagnostic
Diagnostic giving the middle atmospheric temperature anomaly fields \(\delta T_{\rm a}\) at 500hPa. It is identified with the baroclinic streamfunction \(\theta_{\rm a}\) of the system. See also Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.temperatures.MiddleAtmosphericTemperatureDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
AtmosphericTemperatureDiagnostic
Diagnostic giving the middle atmospheric temperature fields \(T_{\rm a} = T_{{\rm a}, 0} + \delta T_{\rm a}\) at 500hPa, where \(T_{{\rm a}, 0}\) is the reference temperature
T0
or the 0-th order dynamic temperature.- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
Notes
Only works if the heat exchange scheme is activated, i.e. does not work with the Newton cooling scheme.
- class qgs.diagnostics.temperatures.MiddleAtmosphericTemperatureMeridionalGradientDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
AtmosphericTemperatureMeridionalGradientDiagnostic
Diagnostic giving the meridional gradient of the middle atmospheric temperature fields \(\partial_y T_{\rm a}\) at 500hPa. It is identified with the meridional gradient of the baroclinic streamfunction \(\partial_y \theta_{\rm a}\) of the system. See also Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.temperatures.OceanicLayerTemperatureAnomalyDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
OceanicTemperatureDiagnostic
Diagnostic giving the oceanic layer temperature anomaly fields \(\delta T_{\rm o}\).
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.temperatures.OceanicLayerTemperatureDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
OceanicTemperatureDiagnostic
Diagnostic giving the oceanic layer temperature fields \(T_{\rm o} = T_{{\rm o}, 0} + \delta T_{\rm o}\), where \(T_{{\rm o}, 0}\) is the reference temperature
T0
or the 0-th order dynamic temperature.- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.temperatures.OceanicTemperatureDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Bases:
FieldDiagnostic
General base class for atmospheric temperature fields diagnostic. Provide a spatial gridded representation of the fields. This is an abstract base class, it must be subclassed to create new diagnostics!
- Parameters
Diagnostic variables classes
Classes defining multiple scalar diagnostics (variables) of the model and used to analyze its outputs.
Description of the classes
VariablesDiagnostic
: General class to get and show the scalar variables of the models.GeopotentialHeightDifferenceDiagnostic
: Class to compute and show the geopotential height difference between points of the model’s domain.
- class qgs.diagnostics.variables.GeopotentialHeightDifferenceDiagnostic(points_list, model_params, dimensional)[source]
Bases:
VariablesDiagnostic
Class to compute and show the geopotential height difference between points of the model’s domain.
- Parameters
- class qgs.diagnostics.variables.VariablesDiagnostic(variable_list, model_params, dimensional)[source]
Bases:
Diagnostic
General class to create multiple scalar diagnostics based on the variables of the model.
- Parameters
- animate(variables='all', output='animate', style='2Dscatter', background=None, ax=None, figsize=(16, 9), show_time=True, stride=1, plot_kwargs=None, anim_kwargs=None, show=True)[source]
Return the output of the plot method animated over time. Show a red dot moving and depicting the current value of the model’s selected variables.
- Parameters
variables (str or list(int)) – List of the model variables to consider as diagnostics. Default to all, i.e. select all the variables of the model.
output (str, optional) –
Define the kind of animation being created. Can be:
animate: Create and show a
ipywidgets.widgets.interactive
widget. Works only in Jupyter notebooks.show: Create and show an animation with the
matplotlib.animation
module. Works only in IPython or Python.
style (str, optional) –
The style of the plot. Can be:
timeserie: Plot all the selected variables as a function of time.
moving-timeserie: Plot all the selected variables as a function of time. Draw the lines as the time evolves.
2Dscatter: Plot the first two selected variables on a 2D scatter plot.
3Dscatter: Plot the first three selected variables on a 3D scatter plot.
background (VariablesDiagnostic, optional) – The variables diagnostic data used as background for the evolving red dot. If None, use the current VariablesDiagnostic instance. Default to None.
ax (Axes, optional) – An axes on which to plot. If not provided, create a new one.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
show_time (bool, optional) – Show the timestamp on the plot or not. Only valid for scatter plots. Default to True.
stride (int, optional) – Specify the time step of the animation. Works only with output set to animate.
plot_kwargs (dict, optional) – Arguments to pass to the background plot.
anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation. Works only with output set to show.show (bool, optional) – Whether to plot or not the animation.
- Returns
The animation object or the callable to update the widget, depending on the value of the output and show parameters.
- Return type
FuncAnimation or DisplayHandle or callable
- movie(variables='all', output='html', filename='', style='2Dscatter', background=None, ax=None, figsize=(16, 9), plot_kwargs=None, anim_kwargs=None)[source]
Create and return a movie of the output of the plot method animated over time. Show a red dot moving and depicting the current value of the model’s selected variables.
- Parameters
variables (str or list(int)) – List of the model variables to consider as diagnostics. Default to all, i.e. select all the variables of the model.
output (str, optional) –
Define the kind of movie being created. Can be:
jshtml: Generate an interactive HTML representation of the animation.
html5: Generate the movie as HTML5 code.
html: Output the movie as a HTML video tag.
ihtml: Output the interactive movie as a HTML video tag.
save: Save the movie in MP4 format (H264 codec).
Default to html.
filename (str, optional) – Filename (and path) where to save the movie. Needed if output is set to save.
style (str, optional) –
The style of the plot. Can be:
timeserie: Plot all the selected variables as a function of time.
moving-timeserie: Plot all the selected variables as a function of time. Draw the lines as the time evolves.
2Dscatter: Plot the first two selected variables on a 2D scatter plot.
3Dscatter: Plot the first three selected variables on a 3D scatter plot.
background (VariablesDiagnostic, optional) – The variables diagnostic data used as background for the evolving red dot. If None, use the current VariablesDiagnostic instance. Default to None.
ax (Axes, optional) – An axes on which to plot. If not provided, create a new one.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
plot_kwargs (dict, optional) – Arguments to pass to the background plot.
anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation.
- Returns
The animation object or the HTML code or tag.
- Return type
FuncAnimation or HTML code or HTML tag
- plot(time_index=0, variables='all', style='timeserie', ax=None, figsize=(16, 9), plot_kwargs=None, **kwargs)[source]
Plot the multiple scalar diagnostic provided.
- Parameters
time_index (int) – The time index of the data. Not used in this subclass.
variables (str or list(int)) – List of the model variables to consider as diagnostics. Default to all, i.e. select all the variables of the model.
style (str, optional) –
The style of the plot. Can be:
timeserie: Plot all the selected variables as a function of time.
2Dscatter: Plot the first two selected variables on a 2D scatter plot.
3Dscatter: Plot the first three selected variables on a 3D scatter plot.
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
plot_kwargs (dict, optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.
- Returns
An axes where the data were plotted.
- Return type
Multidiagnostic class
This class is used to analyze and plot simultaneously several diagnostic together.
- class qgs.diagnostics.multi.FieldsDiagnosticsList(diagnostics_list=None)[source]
General base class for plotting multiple diagnostics on a single axe. The diagnostics must be provided as a list. Assumes that the first diagnostic in the list set the parameters that are not specified.
- Parameters
diagnostics_list (list) – List of initialized auxialiary diagnostics to plot with the main diagnostic.
- animate(output='animate', style='image', ax=None, figsize=(16, 9), contour_labels=True, color_bar=True, show_time=True, stride=1, plot_kwargs=None, oro_kwargs=None, anim_kwargs=None, show=True)[source]
Return the output of the plot method animated over time. Almost all the parameters (except animate, ax, figsize, stride, anim_kwargs and show) and fig should be lists corresponding to diagnostics in the list. If a single parameter is provided, it applies to all the diagnostics.
- Parameters
output (str, optional) –
Define the kind of animation being created. Can be:
animate: Create and show a
ipywidgets.widgets.interactive
widget. Works only in Jupyter notebooks.show: Create and show an animation with the
matplotlib.animation
module. Works only in IPython or Python.
The style of the plot. Can be:
image: show the fields as images with a given colormap specified in the plot_kwargs argument.
contour: show the fields as contour superimposed on the image of the orographic height (see the oro_kwargs below).
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
contour_labels (list(bool), optional) – If style is set to contour, specify if the contours must be labelled with their value or not. Default to True.
color_bar (list(bool), optional) – Specify if a color bar must be drawn beside the plot or not. Default to True.
show_time (list(bool), optional) – Show the timestamp of the field on the plot or not. Default to True.
stride (int, optional) – Specify the time step of the animation. Works only with output set to animate.
plot_kwargs (list(dict), optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.oro_kwargs (list(dict), optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method plotting the image of the orography if style is set to contour.anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation. Works only with output set to show.show (bool, optional) – Whether to plot or not the animation.
- Returns
The animation object or the callable to update the widget, depending on the value of the output and show parameters.
- Return type
FuncAnimation or DisplayHandle or callable
- append_diagnostic(diagnostic)[source]
Method to add an auxiliary diagnostic to the list.
- Parameters
diagnostic (Diagnostic) – The diagnostic to add to the list.
- movie(output='html', filename='', style='image', ax=None, figsize=(16, 9), contour_labels=True, color_bar=True, show_time=True, plot_kwargs=None, oro_kwargs=None, anim_kwargs=None)[source]
Create and return a movie of the output of the plot method animated over time. Almost all the parameters (except output, filename, ax, figsize, and anim_kwargs) and fig should be lists corresponding to diagnostics in the list. If a single parameter is provided, it applies to all the diagnostics.
- Parameters
output (str, optional) –
Define the kind of movie being created. Can be:
jshtml: Generate an interactive HTML representation of the animation.
html5: Generate the movie as HTML5 code.
html: Output the movie as a HTML video tag.
ihtml: Output the interactive movie as a HTML video tag.
save: Save the movie in MP4 format (H264 codec).
Default to html.
filename (str, optional) – Filename (and path) where to save the movie. Needed if output is set to save.
style (lits(str), optional) –
The style of the plot. Can be:
image: show the fields as images with a given colormap specified in the plot_kwargs argument.
contour: show the fields as contour superimposed on the image of the orographic height (see the oro_kwargs below).
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
contour_labels (list(bool), optional) – If style is set to contour, specify if the contours must be labelled with their value or not. Default to True.
color_bar (list(bool), optional) – Specify if a color bar must be drawn beside the plot or not. Default to True.
show_time (list(bool), optional) – Show the timestamp of the field on the plot or not. Default to True.
plot_kwargs (list(dict), optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.oro_kwargs (list(dict), optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method plotting the image of the orography if style is set to contour.anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation.
- Returns
The animation object or the HTML code or tag.
- Return type
FuncAnimation or HTML code or HTML tag
- plot(time_index, style='image', ax=None, figsize=(16, 9), contour_labels=True, color_bar=True, show_time=True, plot_kwargs=None, oro_kwargs=None)[source]
Plot the field of the provided diagnostics at the given time index. Almost all the parameters (except ax and figsize) and fig should be lists corresponding to diagnostics in the list. If a single parameter is provided, it applies to all the diagnostics.
- Parameters
The style of the plot. Can be:
image: show the fields as images with a given colormap specified in the plot_kwargs argument.
contour: show the fields as contour superimposed on the image of the orographic height (if it exists, see the oro_kwargs below).
ax (Axes, optional) – An axes on which to plot the fields.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple.
contour_labels (list(bool), optional) – If style is set to contour, specify if the contours must be labelled with their value or not. Default to True.
color_bar (list(bool), optional) – Specify if a color bar must be drawn beside the plot or not. Default to True.
show_time (list(bool), optional) – Show the timestamp of the field on the plot or not. Default to True.
plot_kwargs (list(dict), optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method if style is set to image, or to thematplotlib.axes.Axes.contour()
method if style is set to contour.oro_kwargs (list(dict), optional) – Arguments to pass to the
matplotlib.axes.Axes.imshow()
method plotting the image of the orography if style is set to contour.
- Returns
An axes where the data were plotted.
- Return type
- set_data(time, data, index=None)[source]
Provide the model data to the index-th diagnostic.
- Parameters
time (ndarray) – The time (in nondimensional timeunits) corresponding to the data. Its length should match the length of the last axis of the provided data.
data (ndarray) – The model output data that the user want to convert using the diagnostic. Should be a 2D array of shape (
ndim
, number_of_timesteps).index (int or None) – The index of the diagnostic in the list to provide the data to.
- class qgs.diagnostics.multi.MultiDiagnostic(nrows, ncols)[source]
class analyze and plot simultaneously several diagnostic together. The diagnostics information are plotted in arrays of fixed dimensions, defining the total number of diagnostics that the object can hold.
- Parameters
- add_diagnostic(diagnostic, position=None, diagnostic_kwargs=None, plot_kwargs=None)[source]
Method to add a diagnostic to the list.
- Parameters
diagnostic (Diagnostic) – Diagnostic to add.
position (tuple(int), optional) – 2-tuple specifying the position of the diagnostic in the plotting array. Find a free spot in the array if not specified. If the plotting array is full, a position must be provided to overwrite a already defined diagnostic.
diagnostic_kwargs (dict, optional) – Dictionary of arguments to pass to the plot, animate and movie method of the provided diagnostic.
plot_kwargs (dict, optional) – Specific plot_kwargs argument to pass to the plot, animate and movie method of the provided diagnostic. If provided, overwrite the one possibly present in the diagnostic_kwargs above.
- animate(output='animate', figure=None, figsize=(16, 9), stride=1, anim_kwargs=None)[source]
Return the output of the plot method animated over time.
- Parameters
output (str, optional) –
Define the kind of animation being created. Can be:
animate: Create and show a
ipywidgets.widgets.interactive
widget. Works only in Jupyter notebooks.show: Create and show an animation with the
matplotlib.animation
module. Works only in IPython or Python.
figure (Figure, optional) – The Matplotlib figure instance where the data are plotted. Update the figure attribute of the object. If not provided, use the figure attribute of the object.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple. Used only if a new figure must be created.
stride (int, optional) – Specify the time step of the animation. Works only with output set to animate.
anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation. Works only with output set to show.
- Returns
The animation object.
- Return type
FuncAnimation or DisplayHandle
- property diagnostic_positions
Position occupied by each diagnostic in the plotting array.
- property diagnostics_list
The list of stored diagnostics.
- Type
- movie(output='html', filename='', figure=None, figsize=(16, 9), anim_kwargs=None)[source]
Create and return a movie of the output of the plot method animated over time.
- Parameters
output (str, optional) –
Define the kind of movie being created. Can be:
jshtml: Generate an interactive HTML representation of the animation.
html5: Generate the movie as HTML5 code.
html: Output the movie as a HTML video tag.
ihtml: Output the interactive movie as a HTML video tag.
save: Save the movie in MP4 format (H264 codec).
Default to html.
filename (str, optional) – Filename (and path) where to save the movie. Needed if output is set to save.
figure (Figure, optional) – The Matplotlib figure instance where the data are plotted. Update the figure attribute of the object. If not provided, use the figure attribute of the object.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple. Used only if a new figure must be created.
anim_kwargs (dict, optional) – Arguments to pass to the
matplotlib.animation.FuncAnimation
instantiation method. Specify the parameters of the animation. Works only with output set to show.
- Returns
The animation object or the HTML code or tag.
- Return type
FuncAnimation or HTML code or HTML tag
- plot(time_index, figure=None, figsize=(16, 9), tight_layout=False)[source]
Plot the fields of the provided data at the given time index.
- Parameters
time_index (int) – The time index of the data.
figure (Figure, optional) – The Matplotlib figure instance where the data are plotted. Update the figure attribute of the object. If not provided, use the figure attribute of the object.
figsize (tuple(float), optional) – The size of the figure in inches as a 2-tuple. Used only if a new figure must be created.
tight_layout (bool, optional) – Enforce a tight layout of the diagnostics axes.
Differential diagnostic base class
Abstract base classes defining diagnostics on differnentiated grids.
Description of the classes
DifferentialFieldDiagnostic
: Base class for diagnostics returning model’s fields based on differential grids.
Warning
These are abstract base class, they must be subclassed to create new diagnostics!
- class qgs.diagnostics.differential.DifferentialFieldDiagnostic(model_params, dimensional)[source]
General base class for differential fields diagnostic. This is an abstract base class, it must be subclassed to create new diagnostics!
- Parameters
Diagnostic wind classes
Classes defining wind fields diagnostics.
Description of the classes
AtmosphericWindDiagnostic
: General base class for atmospheric wind diagnostic.LowerLayerAtmosphericVWindDiagnostic
: Diagnostic giving the lower layer atmospheric V wind fields \(\partial_x \psi^3_{\rm a}\).LowerLayerAtmosphericUWindDiagnostic
: Diagnostic giving the lower layer atmospheric U wind fields \(- \partial_y \psi^3_{\rm a}\).MiddleAtmosphericVWindDiagnostic
: Diagnostic giving the middle atmospheric V wind fields \(\partial_x \psi_{\rm a}\).MiddleAtmosphericUWindDiagnostic
: Diagnostic giving the middle atmospheric U wind fields \(- \partial_y \psi_{\rm a}\).UpperLayerAtmosphericVWindDiagnostic
: Diagnostic giving the upper layer atmospheric V wind fields \(\partial_x \psi^1_{\rm a}\).UpperLayerAtmosphericUWindDiagnostic
: Diagnostic giving the upper layer atmospheric U wind fields \(- \partial_y \psi^1_{\rm a}\).LowerLayerAtmosphericWindIntensityDiagnostic
: Diagnostic giving the lower layer atmospheric horizontal wind intensity fields.MiddleAtmosphericWindIntensityDiagnostic
: Diagnostic giving the middle atmospheric horizontal wind intensity fields.UpperLayerAtmosphericWindIntensityDiagnostic
: Diagnostic giving the upper layer atmospheric horizontal wind intensity fields.MiddleLayerVerticalVelocity
: Diagnostic giving the middle atmospheric layer vertical wind intensity fields.
- class qgs.diagnostics.wind.AtmosphericWindDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
General base class for atmospheric wind fields diagnostic. Provide a spatial gridded representation of the fields. This is an abstract base class, it must be subclassed to create new diagnostics!
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool) – Indicate if the output diagnostic must be dimensionalized or not.
- class qgs.diagnostics.wind.LowerLayerAtmosphericUWindDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the lower layer atmospheric U wind fields \(- \partial_y \psi^3_{\rm a}\). Computed as \(- \partial_y \psi^3_{\rm a} = - \partial_y \psi_{\rm a} + \partial_y \theta_{\rm a}\) where \(\psi_{\rm a}\) and \(\theta_{\rm a}\) are respectively the barotropic and baroclinic streamfunctions. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.LowerLayerAtmosphericVWindDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the lower layer atmospheric V wind fields \(\partial_x \psi^3_{\rm a}\). Computed as \(\partial_x \psi^3_{\rm a} = \partial_x \psi_{\rm a} - \partial_x \theta_{\rm a}\) where \(\psi_{\rm a}\) and \(\theta_{\rm a}\) are respectively the barotropic and baroclinic streamfunctions. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.LowerLayerAtmosphericWindIntensityDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the lower layer atmospheric horizontal wind intensity fields.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.MiddleAtmosphericUWindDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the middle atmospheric U wind fields \(- \partial_y \psi_{\rm a}\) where \(\psi_{\rm a}\) is the barotropic streamfunction. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.MiddleAtmosphericVWindDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the middle atmospheric V wind fields \(\partial_x \psi_{\rm a}\) where \(\psi_{\rm a}\) is the barotropic streamfunction. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.MiddleAtmosphericWindIntensityDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the middle atmospheric horizontal wind intensity fields.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.MiddleLayerVerticalVelocity(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the middle atmospheric layer vertical wind intensity fields \(\omega\). See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- set_data(time, data)[source]
Provide the model data to the diagnostic.
- Parameters
time (ndarray) – The time (in nondimensional timeunits) corresponding to the data. Its length should match the length of the last axis of the provided data.
data (ndarray) – The model output data that the user want to convert using the diagnostic. Should be a 2D array of shape (
ndim
, number_of_timesteps).
- class qgs.diagnostics.wind.UpperLayerAtmosphericUWindDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the upper layer atmospheric U wind fields \(- \partial_y \psi^1_{\rm a}\). Computed as \(- \partial_y \psi^1_{\rm a} = - \partial_y \psi_{\rm a} - \partial_y \theta_{\rm a}\) where \(\psi_{\rm a}\) and \(\theta_{\rm a}\) are respectively the barotropic and baroclinic streamfunctions. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.UpperLayerAtmosphericVWindDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the upper layer atmospheric V wind fields \(\partial_x \psi^1_{\rm a}\). Computed as \(\partial_x \psi^1_{\rm a} = \partial_x \psi_{\rm a} + \partial_x \theta_{\rm a}\) where \(\psi_{\rm a}\) and \(\theta_{\rm a}\) are respectively the barotropic and baroclinic streamfunctions. See also the Atmospheric component and Mid-layer equations and the thermal wind relation sections.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
- class qgs.diagnostics.wind.UpperLayerAtmosphericWindIntensityDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True)[source]
Diagnostic giving the lower layer atmospheric horizontal wind intensity fields.
- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
Diagnostic eddy classes
Classes defining eddy related fields diagnostics.
Description of the classes
MiddleAtmosphericEddyHeatFluxDiagnostic
: Diagnostic giving the middle atmospheric eddy heat flux field.MiddleAtmosphericEddyHeatFluxProfileDiagnostic
: Diagnostic giving the middle atmospheric eddy heat flux zonally averaged profile.
- class qgs.diagnostics.eddy.MiddleAtmosphericEddyHeatFluxDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True, temp_mean_state=None, vwind_mean_state=None, heat_capacity=None)[source]
Diagnostic giving the middle atmospheric eddy heat flux field. Computed as \(v'_{\rm a} \, T'_{\rm a}\) and scaled with the atmospheric specific heat capicity if available (through the heat_capacity argument or the
gamma
parameter).- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
temp_mean_state (MiddleAtmosphericTemperatureDiagnostic, optional) – A temperature diagnostic with a long trajectory as data to compute the mean temperature field. If not provided, compute the mean with the data stored in the object.
vwind_mean_state (MiddleAtmosphericVWindDiagnostic, optional) – A \(v\) wind diagnostic with a long trajectory as data to compute the mean wind field. If not provided, compute the mean with the data stored in the object.
heat_capacity (float, optional) – The air specific heat capacity. If not provided, uses the one of
gamma
if available or or let the heat flux in K m s^{-1}.
- class qgs.diagnostics.eddy.MiddleAtmosphericEddyHeatFluxProfileDiagnostic(model_params, delta_x=None, delta_y=None, dimensional=True, temp_mean_state=None, vwind_mean_state=None, heat_capacity=None)[source]
Diagnostic giving the middle atmospheric eddy heat flux zonally averaged profile. Computed as \(\Phi_{\rm e} = \overline{v'_{\rm a} \, T'_{\rm a}} = \frac{n}{2\pi} \, \int_0^{2\pi/n} \Phi_{\rm e} \, \mathrm{d} x\) where \(v'_{\rm a} \, T'_{\rm a}\) is the eddy heat flux scaled with the atmospheric specific heat capicity if available (through the heat_capacity argument or the
gamma
parameter).- Parameters
model_params (QgParams) – An instance of the model parameters.
delta_x (float, optional) – Spatial step in the zonal direction x for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
delta_y (float, optional) – Spatial step in the meridional direction y for the gridded representation of the field. If not provided, take an optimal guess based on the provided model’s parameters.
dimensional (bool, optional) – Indicate if the output diagnostic must be dimensionalized or not. Default to True.
temp_mean_state (MiddleAtmosphericTemperatureDiagnostic, optional) – A temperature diagnostic with a long trajectory as data to compute the mean temperature field. If not provided, compute the mean with the data stored in the object.
vwind_mean_state (MiddleAtmosphericVWindDiagnostic, optional) – A \(v\) wind diagnostic with a long trajectory as data to compute the mean wind field. If not provided, compute the mean with the data stored in the object.
heat_capacity (float, optional) – The air specific heat capacity. If not provided, uses the one of
gamma
if available or or let the heat flux in K m s^{-1}.