cf.Field.creation_commands¶
-
Field.creation_commands(representative_data=False, namespace='cf', indent=0, string=True)[source]¶ Return the commands that would create the field construct.
See also
Parameters: - representative_data:
bool, optional Return one-line representations of
Datainstances, which are not executable code but prevent the data being converted in its entirety to a string representation.- namespace:
str, optional The namespace containing classes of the
cfpackage. This is prefixed to the class name in commands that instantiate instances ofcfobjects. By default, namespace is'cf', i.e. it is assumed thatcfwas imported asimport cf.- Parameter example:
If
cfwas imported asimport cf as cfpthen setnamespace='cfp'- Parameter example:
If
cfwas imported asfrom cf import *then setnamespace=''
- indent:
int, optional Indent each line by this many spaces. Ignored if string is False.
- string:
bool, optional If False then return each command as an element of a
list. By default the the commands are concatenated into a string, with a new line inserted between each command.
Returns: Examples:
>>> q = cf.example_field(0) >>> print(q) Field: specific_humidity (ncvar%q) ---------------------------------- Data : specific_humidity(latitude(5), longitude(8)) 1 Cell methods : area: mean Dimension coords: latitude(5) = [-75.0, ..., 75.0] degrees_north : longitude(8) = [22.5, ..., 337.5] degrees_east : time(1) = [2019-01-01 00:00:00] >>> print(q.creation_commands()) f = cf.Field() # f.set_properties({'Conventions': 'CF-1.7', 'project': 'research', 'standard_name': 'specific_humidity', 'units': '1'}) f.nc_set_variable('q') # # domain_axis c = cf.DomainAxis(size=5) c.nc_set_dimension('lat') f.set_construct(c, key='domainaxis0') # # domain_axis c = cf.DomainAxis(size=8) c.nc_set_dimension('lon') f.set_construct(c, key='domainaxis1') # # domain_axis c = cf.DomainAxis(size=1) f.set_construct(c, key='domainaxis2') # # field data data = cf.Data([[0.007, 0.034, 0.003, 0.014, 0.018, 0.037, 0.024, 0.029], [0.023, 0.036, 0.045, 0.062, 0.046, 0.073, 0.006, 0.066], [0.11, 0.131, 0.124, 0.146, 0.087, 0.103, 0.057, 0.011], [0.029, 0.059, 0.039, 0.07, 0.058, 0.072, 0.009, 0.017], [0.006, 0.036, 0.019, 0.035, 0.018, 0.037, 0.034, 0.013]], units='1', dtype='f8') f.set_data(data, axes=('domainaxis0', 'domainaxis1')) # # dimension_coordinate c = cf.DimensionCoordinate() c.set_properties({'units': 'degrees_north', 'standard_name': 'latitude'}) c.nc_set_variable('lat') data = cf.Data([-75.0, -45.0, 0.0, 45.0, 75.0], units='degrees_north', dtype='f8') c.set_data(data) b = cf.Bounds() b.set_properties({'units': 'degrees_north'}) b.nc_set_variable('lat_bnds') data = cf.Data([[-90.0, -60.0], [-60.0, -30.0], [-30.0, 30.0], [30.0, 60.0], [60.0, 90.0]], units='degrees_north', dtype='f8') b.set_data(data) c.set_bounds(b) f.set_construct(c, axes=('domainaxis0',), key='dimensioncoordinate0', copy=False) # # dimension_coordinate c = cf.DimensionCoordinate() c.set_properties({'units': 'degrees_east', 'standard_name': 'longitude'}) c.nc_set_variable('lon') data = cf.Data([22.5, 67.5, 112.5, 157.5, 202.5, 247.5, 292.5, 337.5], units='degrees_east', dtype='f8') c.set_data(data) b = cf.Bounds() b.set_properties({'units': 'degrees_east'}) b.nc_set_variable('lon_bnds') data = cf.Data([[0.0, 45.0], [45.0, 90.0], [90.0, 135.0], [135.0, 180.0], [180.0, 225.0], [225.0, 270.0], [270.0, 315.0], [315.0, 360.0]], units='degrees_east', dtype='f8') b.set_data(data) c.set_bounds(b) f.set_construct(c, axes=('domainaxis1',), key='dimensioncoordinate1', copy=False) # # dimension_coordinate c = cf.DimensionCoordinate() c.set_properties({'units': 'days since 2018-12-01', 'standard_name': 'time'}) c.nc_set_variable('time') data = cf.Data([31.0], units='days since 2018-12-01', dtype='f8') c.set_data(data) f.set_construct(c, axes=('domainaxis2',), key='dimensioncoordinate2', copy=False) # # cell_method c = cf.CellMethod() c.method = 'mean' c.axes = ('area',) f.set_construct(c)
>>> print(q.creation_commands(representative_data=True, namespace='', indent=4)) f = Field() # f.set_properties({'Conventions': 'CF-1.7', 'project': 'research', 'standard_name': 'specific_humidity', 'units': '1'}) f.nc_set_variable('q') # # domain_axis c = DomainAxis(size=5) c.nc_set_dimension('lat') f.set_construct(c, key='domainaxis0') # # domain_axis c = DomainAxis(size=8) c.nc_set_dimension('lon') f.set_construct(c, key='domainaxis1') # # domain_axis c = DomainAxis(size=1) f.set_construct(c, key='domainaxis2') # # field data data = <CF Data(5, 8): [[0.007, ..., 0.013]] 1> # Representative data f.set_data(data, axes=('domainaxis0', 'domainaxis1')) # # dimension_coordinate c = DimensionCoordinate() c.set_properties({'units': 'degrees_north', 'standard_name': 'latitude'}) c.nc_set_variable('lat') data = <CF Data(5): [-75.0, ..., 75.0] degrees_north> # Representative data c.set_data(data) b = Bounds() b.set_properties({'units': 'degrees_north'}) b.nc_set_variable('lat_bnds') data = <CF Data(5, 2): [[-90.0, ..., 90.0]] degrees_north> # Representative data b.set_data(data) c.set_bounds(b) f.set_construct(c, axes=('domainaxis0',), key='dimensioncoordinate0', copy=False) # # dimension_coordinate c = DimensionCoordinate() c.set_properties({'units': 'degrees_east', 'standard_name': 'longitude'}) c.nc_set_variable('lon') data = <CF Data(8): [22.5, ..., 337.5] degrees_east> # Representative data c.set_data(data) b = Bounds() b.set_properties({'units': 'degrees_east'}) b.nc_set_variable('lon_bnds') data = <CF Data(8, 2): [[0.0, ..., 360.0]] degrees_east> # Representative data b.set_data(data) c.set_bounds(b) f.set_construct(c, axes=('domainaxis1',), key='dimensioncoordinate1', copy=False) # # dimension_coordinate c = DimensionCoordinate() c.set_properties({'units': 'days since 2018-12-01', 'standard_name': 'time'}) c.nc_set_variable('time') data = <CF Data(1): [2019-01-01 00:00:00]> # Representative data c.set_data(data) f.set_construct(c, axes=('domainaxis2',), key='dimensioncoordinate2', copy=False) # # cell_method c = CellMethod() c.method = 'mean' c.axes = ('area',) f.set_construct(c)
- representative_data: