-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathfield.py
1900 lines (1704 loc) · 97 KB
/
field.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
import collections
import datetime
import math
from ctypes import c_float
from ctypes import c_int
from ctypes import POINTER
from ctypes import pointer
from ctypes import Structure
import dask.array as da
import numpy as np
import xarray as xr
from pathlib import Path
import parcels.tools.interpolation_utils as i_u
from .fieldfilebuffer import (NetcdfFileBuffer, DeferredNetcdfFileBuffer,
DaskFileBuffer, DeferredDaskFileBuffer)
from .grid import CGrid
from .grid import Grid
from .grid import GridCode
from parcels.tools.converters import Geographic
from parcels.tools.converters import GeographicPolar
from parcels.tools.converters import TimeConverter
from parcels.tools.converters import UnitConverter
from parcels.tools.converters import unitconverters_map
from parcels.tools.statuscodes import FieldOutOfBoundError
from parcels.tools.statuscodes import FieldOutOfBoundSurfaceError
from parcels.tools.statuscodes import FieldSamplingError
from parcels.tools.statuscodes import TimeExtrapolationError
from parcels.tools.loggers import logger
__all__ = ['Field', 'VectorField', 'SummedField', 'NestedField']
class Field(object):
"""Class that encapsulates access to field data.
:param name: Name of the field
:param data: 2D, 3D or 4D numpy array of field data.
1. If data shape is [xdim, ydim], [xdim, ydim, zdim], [xdim, ydim, tdim] or [xdim, ydim, zdim, tdim],
whichever is relevant for the dataset, use the flag transpose=True
2. If data shape is [ydim, xdim], [zdim, ydim, xdim], [tdim, ydim, xdim] or [tdim, zdim, ydim, xdim],
use the flag transpose=False
3. If data has any other shape, you first need to reorder it
:param lon: Longitude coordinates (numpy vector or array) of the field (only if grid is None)
:param lat: Latitude coordinates (numpy vector or array) of the field (only if grid is None)
:param depth: Depth coordinates (numpy vector or array) of the field (only if grid is None)
:param time: Time coordinates (numpy vector) of the field (only if grid is None)
:param mesh: String indicating the type of mesh coordinates and
units used during velocity interpolation: (only if grid is None)
1. spherical: Lat and lon in degree, with a
correction for zonal velocity U near the poles.
2. flat (default): No conversion, lat/lon are assumed to be in m.
:param timestamps: A numpy array containing the timestamps for each of the files in filenames, for loading
from netCDF files only. Default is None if the netCDF dimensions dictionary includes time.
:param grid: :class:`parcels.grid.Grid` object containing all the lon, lat depth, time
mesh and time_origin information. Can be constructed from any of the Grid objects
:param fieldtype: Type of Field to be used for UnitConverter when using SummedFields
(either 'U', 'V', 'Kh_zonal', 'Kh_meridional' or None)
:param transpose: Transpose data to required (lon, lat) layout
:param vmin: Minimum allowed value on the field. Data below this value are set to zero
:param vmax: Maximum allowed value on the field. Data above this value are set to zero
:param time_origin: Time origin (TimeConverter object) of the time axis (only if grid is None)
:param interp_method: Method for interpolation. Options are 'linear' (default), 'nearest',
'linear_invdist_land_tracer', 'cgrid_velocity', 'cgrid_tracer' and 'bgrid_velocity'
:param allow_time_extrapolation: boolean whether to allow for extrapolation in time
(i.e. beyond the last available time snapshot)
:param time_periodic: To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object).
The last value of the time series can be provided (which is the same as the initial one) or not (Default: False)
This flag overrides the allow_time_interpolation and sets it to False
:param chunkdims_name_map (opt.): gives a name map to the FieldFileBuffer that declared a mapping between chunksize name, NetCDF dimension and Parcels dimension;
required only if currently incompatible OCM field is loaded and chunking is used by 'chunksize' (which is the default)
For usage examples see the following tutorials:
* `Nested Fields <https://nbviewer.jupyter.org/github/OceanParcels/parcels/blob/master/parcels/examples/tutorial_NestedFields.ipynb>`_
* `Summed Fields <https://nbviewer.jupyter.org/github/OceanParcels/parcels/blob/master/parcels/examples/tutorial_SummedFields.ipynb>`_
"""
def __init__(self, name, data, lon=None, lat=None, depth=None, time=None, grid=None, mesh='flat', timestamps=None,
fieldtype=None, transpose=False, vmin=None, vmax=None, time_origin=None,
interp_method='linear', allow_time_extrapolation=None, time_periodic=False, gridindexingtype='nemo', **kwargs):
if not isinstance(name, tuple):
self.name = name
self.filebuffername = name
else:
self.name, self.filebuffername = name
self.data = data
time_origin = TimeConverter(0) if time_origin is None else time_origin
if grid:
if grid.defer_load and isinstance(data, np.ndarray):
raise ValueError('Cannot combine Grid from defer_loaded Field with np.ndarray data. please specify lon, lat, depth and time dimensions separately')
self.grid = grid
else:
self.grid = Grid.create_grid(lon, lat, depth, time, time_origin=time_origin, mesh=mesh)
self.igrid = -1
# self.lon, self.lat, self.depth and self.time are not used anymore in parcels.
# self.grid should be used instead.
# Those variables are still defined for backwards compatibility with users codes.
self.lon = self.grid.lon
self.lat = self.grid.lat
self.depth = self.grid.depth
self.fieldtype = self.name if fieldtype is None else fieldtype
if self.grid.mesh == 'flat' or (self.fieldtype not in unitconverters_map.keys()):
self.units = UnitConverter()
elif self.grid.mesh == 'spherical':
self.units = unitconverters_map[self.fieldtype]
else:
raise ValueError("Unsupported mesh type. Choose either: 'spherical' or 'flat'")
self.timestamps = timestamps
if type(interp_method) is dict:
if self.name in interp_method:
self.interp_method = interp_method[self.name]
else:
raise RuntimeError('interp_method is a dictionary but %s is not in it' % name)
else:
self.interp_method = interp_method
self.gridindexingtype = gridindexingtype
if self.interp_method in ['bgrid_velocity', 'bgrid_w_velocity', 'bgrid_tracer'] and \
self.grid.gtype in [GridCode.RectilinearSGrid, GridCode.CurvilinearSGrid]:
logger.warning_once('General s-levels are not supported in B-grid. RectilinearSGrid and CurvilinearSGrid can still be used to deal with shaved cells, but the levels must be horizontal.')
self.fieldset = None
if allow_time_extrapolation is None:
self.allow_time_extrapolation = True if len(self.grid.time) == 1 else False
else:
self.allow_time_extrapolation = allow_time_extrapolation
self.time_periodic = time_periodic
if self.time_periodic is not False and self.allow_time_extrapolation:
logger.warning_once("allow_time_extrapolation and time_periodic cannot be used together.\n \
allow_time_extrapolation is set to False")
self.allow_time_extrapolation = False
if self.time_periodic is True:
raise ValueError("Unsupported time_periodic=True. time_periodic must now be either False or the length of the period (either float in seconds or datetime.timedelta object.")
if self.time_periodic is not False:
if isinstance(self.time_periodic, datetime.timedelta):
self.time_periodic = self.time_periodic.total_seconds()
if not np.isclose(self.grid.time[-1] - self.grid.time[0], self.time_periodic):
if self.grid.time[-1] - self.grid.time[0] > self.time_periodic:
raise ValueError("Time series provided is longer than the time_periodic parameter")
self.grid._add_last_periodic_data_timestep = True
self.grid.time = np.append(self.grid.time, self.grid.time[0] + self.time_periodic)
self.grid.time_full = self.grid.time
self.vmin = vmin
self.vmax = vmax
if not self.grid.defer_load:
self.data = self.reshape(self.data, transpose)
# Hack around the fact that NaN and ridiculously large values
# propagate in SciPy's interpolators
lib = np if isinstance(self.data, np.ndarray) else da
self.data[lib.isnan(self.data)] = 0.
if self.vmin is not None:
self.data[self.data < self.vmin] = 0.
if self.vmax is not None:
self.data[self.data > self.vmax] = 0.
if self.grid._add_last_periodic_data_timestep:
self.data = lib.concatenate((self.data, self.data[:1, :]), axis=0)
self._scaling_factor = None
# Variable names in JIT code
self.dimensions = kwargs.pop('dimensions', None)
self.indices = kwargs.pop('indices', None)
self.dataFiles = kwargs.pop('dataFiles', None)
if self.grid._add_last_periodic_data_timestep and self.dataFiles is not None:
self.dataFiles = np.append(self.dataFiles, self.dataFiles[0])
self._field_fb_class = kwargs.pop('FieldFileBuffer', None)
self.netcdf_engine = kwargs.pop('netcdf_engine', 'netcdf4')
self.loaded_time_indices = []
self.creation_log = kwargs.pop('creation_log', '')
self.chunksize = kwargs.pop('chunksize', None)
self.netcdf_chunkdims_name_map = kwargs.pop('chunkdims_name_map', None)
self.grid.depth_field = kwargs.pop('depth_field', None)
if self.grid.depth_field == 'not_yet_set':
assert self.grid.z4d, 'Providing the depth dimensions from another field data is only available for 4d S grids'
# data_full_zdim is the vertical dimension of the complete field data, ignoring the indices.
# (data_full_zdim = grid.zdim if no indices are used, for A- and C-grids and for some B-grids). It is used for the B-grid,
# since some datasets do not provide the deeper level of data (which is ignored by the interpolation).
self.data_full_zdim = kwargs.pop('data_full_zdim', None)
self.data_chunks = []
self.c_data_chunks = []
self.nchunks = []
self.chunk_set = False
self.filebuffers = [None] * 3
if len(kwargs) > 0:
raise SyntaxError('Field received an unexpected keyword argument "%s"' % list(kwargs.keys())[0])
@classmethod
def get_dim_filenames(cls, filenames, dim):
if isinstance(filenames, str) or not isinstance(filenames, collections.abc.Iterable):
return [filenames]
elif isinstance(filenames, dict):
assert dim in filenames.keys(), \
'filename dimension keys must be lon, lat, depth or data'
filename = filenames[dim]
if isinstance(filename, str):
return [filename]
else:
return filename
else:
return filenames
@staticmethod
def collect_timeslices(timestamps, data_filenames, _grid_fb_class, dimensions, indices, netcdf_engine):
if timestamps is not None:
dataFiles = []
for findex in range(len(data_filenames)):
for f in [data_filenames[findex], ] * len(timestamps[findex]):
dataFiles.append(f)
timeslices = np.array([stamp for file in timestamps for stamp in file])
time = timeslices
else:
timeslices = []
dataFiles = []
for fname in data_filenames:
with _grid_fb_class(fname, dimensions, indices, netcdf_engine=netcdf_engine) as filebuffer:
ftime = filebuffer.time
timeslices.append(ftime)
dataFiles.append([fname] * len(ftime))
timeslices = np.array(timeslices)
time = np.concatenate(timeslices)
dataFiles = np.concatenate(np.array(dataFiles))
if time.size == 1 and time[0] is None:
time[0] = 0
time_origin = TimeConverter(time[0])
time = time_origin.reltime(time)
if not np.all((time[1:] - time[:-1]) > 0):
id_not_ordered = np.where(time[1:] < time[:-1])[0][0]
raise AssertionError(
'Please make sure your netCDF files are ordered in time. First pair of non-ordered files: %s, %s'
% (dataFiles[id_not_ordered], dataFiles[id_not_ordered + 1]))
return time, time_origin, timeslices, dataFiles
@classmethod
def from_netcdf(cls, filenames, variable, dimensions, indices=None, grid=None,
mesh='spherical', timestamps=None, allow_time_extrapolation=None, time_periodic=False,
deferred_load=True, **kwargs):
"""Create field from netCDF file
:param filenames: list of filenames to read for the field. filenames can be a list [files] or
a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data)
In the latetr case, time values are in filenames[data]
:param variable: Tuple mapping field name to variable name in the NetCDF file.
:param dimensions: Dictionary mapping variable names for the relevant dimensions in the NetCDF file
:param indices: dictionary mapping indices for each dimension to read from file.
This can be used for reading in only a subregion of the NetCDF file.
Note that negative indices are not allowed.
:param mesh: String indicating the type of mesh coordinates and
units used during velocity interpolation:
1. spherical (default): Lat and lon in degree, with a
correction for zonal velocity U near the poles.
2. flat: No conversion, lat/lon are assumed to be in m.
:param timestamps: A numpy array of datetime64 objects containing the timestamps for each of the files in filenames.
Default is None if dimensions includes time.
:param allow_time_extrapolation: boolean whether to allow for extrapolation in time
(i.e. beyond the last available time snapshot)
Default is False if dimensions includes time, else True
:param time_periodic: boolean whether to loop periodically over the time component of the FieldSet
This flag overrides the allow_time_interpolation and sets it to False
:param deferred_load: boolean whether to only pre-load data (in deferred mode) or
fully load them (default: True). It is advised to deferred load the data, since in
that case Parcels deals with a better memory management during particle set execution.
deferred_load=False is however sometimes necessary for plotting the fields.
:param gridindexingtype: The type of gridindexing. Either 'nemo' (default) or 'mitgcm' are supported.
See also the Grid indexing documentation on oceanparcels.org
:param chunksize: size of the chunks in dask loading
For usage examples see the following tutorial:
* `Timestamps <https://nbviewer.jupyter.org/github/OceanParcels/parcels/blob/master/parcels/examples/tutorial_timestamps.ipynb>`_
"""
# Ensure the timestamps array is compatible with the user-provided datafiles.
if timestamps is not None:
if isinstance(filenames, list):
assert len(filenames) == len(timestamps), 'Outer dimension of timestamps should correspond to number of files.'
elif isinstance(filenames, dict):
for k in filenames.keys():
if k not in ['lat', 'lon', 'depth', 'time']:
assert(len(filenames[k]) == len(timestamps)), 'Outer dimension of timestamps should correspond to number of files.'
else:
raise TypeError("Filenames type is inconsistent with manual timestamp provision."
+ "Should be dict or list")
if isinstance(variable, str): # for backward compatibility with Parcels < 2.0.0
variable = (variable, variable)
assert len(variable) == 2, 'The variable tuple must have length 2. Use FieldSet.from_netcdf() for multiple variables'
data_filenames = cls.get_dim_filenames(filenames, 'data')
lonlat_filename = cls.get_dim_filenames(filenames, 'lon')
if isinstance(filenames, dict):
assert len(lonlat_filename) == 1
if lonlat_filename != cls.get_dim_filenames(filenames, 'lat'):
raise NotImplementedError('longitude and latitude dimensions are currently processed together from one single file')
lonlat_filename = lonlat_filename[0]
if 'depth' in dimensions:
depth_filename = cls.get_dim_filenames(filenames, 'depth')
if isinstance(filenames, dict) and len(depth_filename) != 1:
raise NotImplementedError('Vertically adaptive meshes not implemented for from_netcdf()')
depth_filename = depth_filename[0]
netcdf_engine = kwargs.pop('netcdf_engine', 'netcdf4')
indices = {} if indices is None else indices.copy()
for ind in indices:
if len(indices[ind]) == 0:
raise RuntimeError('Indices for %s can not be empty' % ind)
assert np.min(indices[ind]) >= 0, \
('Negative indices are currently not allowed in Parcels. '
+ 'This is related to the non-increasing dimension it could generate '
+ 'if the domain goes from lon[-4] to lon[6] for example. '
+ 'Please raise an issue on https://github.com/OceanParcels/parcels/issues '
+ 'if you would need such feature implemented.')
interp_method = kwargs.pop('interp_method', 'linear')
if type(interp_method) is dict:
if variable[0] in interp_method:
interp_method = interp_method[variable[0]]
else:
raise RuntimeError('interp_method is a dictionary but %s is not in it' % variable[0])
_grid_fb_class = NetcdfFileBuffer
with _grid_fb_class(lonlat_filename, dimensions, indices, netcdf_engine) as filebuffer:
lon, lat = filebuffer.lonlat
indices = filebuffer.indices
# Check if parcels_mesh has been explicitly set in file
if 'parcels_mesh' in filebuffer.dataset.attrs:
mesh = filebuffer.dataset.attrs['parcels_mesh']
if 'depth' in dimensions:
with _grid_fb_class(depth_filename, dimensions, indices, netcdf_engine, interp_method=interp_method) as filebuffer:
filebuffer.name = filebuffer.parse_name(variable[1])
if dimensions['depth'] == 'not_yet_set':
depth = filebuffer.depth_dimensions
kwargs['depth_field'] = 'not_yet_set'
else:
depth = filebuffer.depth
data_full_zdim = filebuffer.data_full_zdim
else:
indices['depth'] = [0]
depth = np.zeros(1)
data_full_zdim = 1
kwargs['data_full_zdim'] = data_full_zdim
if len(data_filenames) > 1 and 'time' not in dimensions and timestamps is None:
raise RuntimeError('Multiple files given but no time dimension specified')
if grid is None:
# Concatenate time variable to determine overall dimension
# across multiple files
time, time_origin, timeslices, dataFiles = cls.collect_timeslices(timestamps, data_filenames,
_grid_fb_class, dimensions,
indices, netcdf_engine)
grid = Grid.create_grid(lon, lat, depth, time, time_origin=time_origin, mesh=mesh)
grid.timeslices = timeslices
kwargs['dataFiles'] = dataFiles
elif grid is not None and ('dataFiles' not in kwargs or kwargs['dataFiles'] is None):
# ==== means: the field has a shared grid, but may have different data files, so we need to collect the
# ==== correct file time series again.
_, _, _, dataFiles = cls.collect_timeslices(timestamps, data_filenames, _grid_fb_class,
dimensions, indices, netcdf_engine)
kwargs['dataFiles'] = dataFiles
chunksize = kwargs.get('chunksize', None)
grid.chunksize = chunksize
if 'time' in indices:
logger.warning_once('time dimension in indices is not necessary anymore. It is then ignored.')
if 'full_load' in kwargs: # for backward compatibility with Parcels < v2.0.0
deferred_load = not kwargs['full_load']
if grid.time.size <= 3 or deferred_load is False:
deferred_load = False
if chunksize not in [False, None]:
if deferred_load:
_field_fb_class = DeferredDaskFileBuffer
else:
_field_fb_class = DaskFileBuffer
elif deferred_load:
_field_fb_class = DeferredNetcdfFileBuffer
else:
_field_fb_class = NetcdfFileBuffer
kwargs['FieldFileBuffer'] = _field_fb_class
if not deferred_load:
# Pre-allocate data before reading files into buffer
data_list = []
ti = 0
for tslice, fname in zip(grid.timeslices, data_filenames):
with _field_fb_class(fname, dimensions, indices, netcdf_engine,
interp_method=interp_method, data_full_zdim=data_full_zdim,
chunksize=chunksize) as filebuffer:
# If Field.from_netcdf is called directly, it may not have a 'data' dimension
# In that case, assume that 'name' is the data dimension
filebuffer.name = filebuffer.parse_name(variable[1])
buffer_data = filebuffer.data
if len(buffer_data.shape) == 2:
data_list.append(buffer_data.reshape(sum(((len(tslice), 1), buffer_data.shape), ())))
elif len(buffer_data.shape) == 3:
if len(filebuffer.indices['depth']) > 1:
data_list.append(buffer_data.reshape(sum(((1,), buffer_data.shape), ())))
else:
if type(tslice) not in [list, np.ndarray, da.Array, xr.DataArray]:
tslice = [tslice]
data_list.append(buffer_data.reshape(sum(((len(tslice), 1), buffer_data.shape[1:]), ())))
else:
data_list.append(buffer_data)
if type(tslice) not in [list, np.ndarray, da.Array, xr.DataArray]:
tslice = [tslice]
ti += len(tslice)
lib = np if isinstance(data_list[0], np.ndarray) else da
data = lib.concatenate(data_list, axis=0)
else:
grid.defer_load = True
grid.ti = -1
data = DeferredArray()
data.compute_shape(grid.xdim, grid.ydim, grid.zdim, grid.tdim, len(grid.timeslices))
if allow_time_extrapolation is None:
allow_time_extrapolation = False if 'time' in dimensions else True
kwargs['dimensions'] = dimensions.copy()
kwargs['indices'] = indices
kwargs['time_periodic'] = time_periodic
kwargs['netcdf_engine'] = netcdf_engine
return cls(variable, data, grid=grid, timestamps=timestamps,
allow_time_extrapolation=allow_time_extrapolation, interp_method=interp_method, **kwargs)
@classmethod
def from_xarray(cls, da, name, dimensions, mesh='spherical', allow_time_extrapolation=None,
time_periodic=False, **kwargs):
"""Create field from xarray Variable
:param da: Xarray DataArray
:param name: Name of the Field
:param dimensions: Dictionary mapping variable names for the relevant dimensions in the DataArray
:param mesh: String indicating the type of mesh coordinates and
units used during velocity interpolation:
1. spherical (default): Lat and lon in degree, with a
correction for zonal velocity U near the poles.
2. flat: No conversion, lat/lon are assumed to be in m.
:param allow_time_extrapolation: boolean whether to allow for extrapolation in time
(i.e. beyond the last available time snapshot)
Default is False if dimensions includes time, else True
:param time_periodic: boolean whether to loop periodically over the time component of the FieldSet
This flag overrides the allow_time_interpolation and sets it to False
"""
data = da.data
interp_method = kwargs.pop('interp_method', 'linear')
time = da[dimensions['time']].values if 'time' in dimensions else np.array([0])
depth = da[dimensions['depth']].values if 'depth' in dimensions else np.array([0])
lon = da[dimensions['lon']].values
lat = da[dimensions['lat']].values
time_origin = TimeConverter(time[0])
time = time_origin.reltime(time)
grid = Grid.create_grid(lon, lat, depth, time, time_origin=time_origin, mesh=mesh)
return cls(name, data, grid=grid, allow_time_extrapolation=allow_time_extrapolation,
interp_method=interp_method, **kwargs)
def reshape(self, data, transpose=False):
# Ensure that field data is the right data type
if not isinstance(data, (np.ndarray, da.core.Array)):
data = np.array(data)
if not data.dtype == np.float32:
logger.warning_once("Casting field data to np.float32")
data = data.astype(np.float32)
lib = np if isinstance(data, np.ndarray) else da
if transpose:
data = lib.transpose(data)
if self.grid.lat_flipped:
data = lib.flip(data, axis=-2)
if self.grid.xdim == 1 or self.grid.ydim == 1:
data = lib.squeeze(data) # First remove all length-1 dimensions in data, so that we can add them below
if self.grid.xdim == 1 and len(data.shape) < 4:
if lib == da:
raise NotImplementedError('Length-one dimensions with field chunking not implemented, as dask does not have an `expand_dims` method. Use chunksize=None')
data = lib.expand_dims(data, axis=-1)
if self.grid.ydim == 1 and len(data.shape) < 4:
if lib == da:
raise NotImplementedError('Length-one dimensions with field chunking not implemented, as dask does not have an `expand_dims` method. Use chunksize=None')
data = lib.expand_dims(data, axis=-2)
if self.grid.tdim == 1:
if len(data.shape) < 4:
data = data.reshape(sum(((1,), data.shape), ()))
if self.grid.zdim == 1:
if len(data.shape) == 4:
data = data.reshape(sum(((data.shape[0],), data.shape[2:]), ()))
if len(data.shape) == 4:
errormessage = ('Field %s expecting a data shape of [tdim, zdim, ydim, xdim]. '
'Flag transpose=True could help to reorder the data.' % self.name)
assert data.shape[0] == self.grid.tdim, errormessage
assert data.shape[2] == self.grid.ydim - 2 * self.grid.meridional_halo, errormessage
assert data.shape[3] == self.grid.xdim - 2 * self.grid.zonal_halo, errormessage
if self.gridindexingtype == 'pop':
assert data.shape[1] == self.grid.zdim or data.shape[1] == self.grid.zdim-1, errormessage
else:
assert data.shape[1] == self.grid.zdim, errormessage
else:
assert (data.shape == (self.grid.tdim,
self.grid.ydim - 2 * self.grid.meridional_halo,
self.grid.xdim - 2 * self.grid.zonal_halo)), \
('Field %s expecting a data shape of [tdim, ydim, xdim]. '
'Flag transpose=True could help to reorder the data.' % self.name)
if self.grid.meridional_halo > 0 or self.grid.zonal_halo > 0:
data = self.add_periodic_halo(zonal=self.grid.zonal_halo > 0, meridional=self.grid.meridional_halo > 0, halosize=max(self.grid.meridional_halo, self.grid.zonal_halo), data=data)
return data
def set_scaling_factor(self, factor):
"""Scales the field data by some constant factor.
:param factor: scaling factor
For usage examples see the following tutorial:
* `Unit converters <https://nbviewer.jupyter.org/github/OceanParcels/parcels/blob/master/parcels/examples/tutorial_unitconverters.ipynb>`_
"""
if self._scaling_factor:
raise NotImplementedError(('Scaling factor for field %s already defined.' % self.name))
self._scaling_factor = factor
if not self.grid.defer_load:
self.data *= factor
def set_depth_from_field(self, field):
"""Define the depth dimensions from another (time-varying) field
See `this tutorial <https://nbviewer.jupyter.org/github/OceanParcels/parcels/blob/master/parcels/examples/tutorial_timevaryingdepthdimensions.ipynb>`_
for a detailed explanation on how to set up time-evolving depth dimensions
"""
self.grid.depth_field = field
if self.grid != field.grid:
field.grid.depth_field = field
def __getitem__(self, key):
# TODO: ideally, we'd like to use isinstance(key, ParticleAssessor) here, but that results in cyclic imports between Field and ParticleSet. Could/should be fixed in #913?
if hasattr(key, 'set_index'):
return self.eval(key.time, key.depth, key.lat, key.lon, key)
else:
return self.eval(*key)
def calc_cell_edge_sizes(self):
"""Method to calculate cell sizes based on numpy.gradient method
Currently only works for Rectilinear Grids"""
if not self.grid.cell_edge_sizes:
if self.grid.gtype in (GridCode.RectilinearZGrid, GridCode.RectilinearSGrid):
self.grid.cell_edge_sizes['x'] = np.zeros((self.grid.ydim, self.grid.xdim), dtype=np.float32)
self.grid.cell_edge_sizes['y'] = np.zeros((self.grid.ydim, self.grid.xdim), dtype=np.float32)
x_conv = GeographicPolar() if self.grid.mesh == 'spherical' else UnitConverter()
y_conv = Geographic() if self.grid.mesh == 'spherical' else UnitConverter()
for y, (lat, dy) in enumerate(zip(self.grid.lat, np.gradient(self.grid.lat))):
for x, (lon, dx) in enumerate(zip(self.grid.lon, np.gradient(self.grid.lon))):
self.grid.cell_edge_sizes['x'][y, x] = x_conv.to_source(dx, lon, lat, self.grid.depth[0])
self.grid.cell_edge_sizes['y'][y, x] = y_conv.to_source(dy, lon, lat, self.grid.depth[0])
self.cell_edge_sizes = self.grid.cell_edge_sizes
else:
logger.error(('Field.cell_edge_sizes() not implemented for ', self.grid.gtype, 'grids.',
'You can provide Field.grid.cell_edge_sizes yourself',
'by in e.g. NEMO using the e1u fields etc from the mesh_mask.nc file'))
exit(-1)
def cell_areas(self):
"""Method to calculate cell sizes based on cell_edge_sizes
Currently only works for Rectilinear Grids"""
if not self.grid.cell_edge_sizes:
self.calc_cell_edge_sizes()
return self.grid.cell_edge_sizes['x'] * self.grid.cell_edge_sizes['y']
def search_indices_vertical_z(self, z):
grid = self.grid
z = np.float32(z)
if grid.depth[-1] > grid.depth[0]:
if z < grid.depth[0]:
# Since MOM5 is indexed at cell bottom, allow z at depth[0] - dz where dz = (depth[1] - depth[0])
if self.gridindexingtype == "mom5" and z > 2*grid.depth[0] - grid.depth[1]:
return (-1, z / grid.depth[0])
else:
raise FieldOutOfBoundSurfaceError(0, 0, z, field=self)
elif z > grid.depth[-1]:
raise FieldOutOfBoundError(0, 0, z, field=self)
depth_indices = grid.depth <= z
if z >= grid.depth[-1]:
zi = len(grid.depth) - 2
else:
zi = depth_indices.argmin() - 1 if z >= grid.depth[0] else 0
else:
if z > grid.depth[0]:
raise FieldOutOfBoundSurfaceError(0, 0, z, field=self)
elif z < grid.depth[-1]:
raise FieldOutOfBoundError(0, 0, z, field=self)
depth_indices = grid.depth >= z
if z <= grid.depth[-1]:
zi = len(grid.depth) - 2
else:
zi = depth_indices.argmin() - 1 if z <= grid.depth[0] else 0
zeta = (z-grid.depth[zi]) / (grid.depth[zi+1]-grid.depth[zi])
return (zi, zeta)
def search_indices_vertical_s(self, x, y, z, xi, yi, xsi, eta, ti, time):
grid = self.grid
if self.interp_method in ['bgrid_velocity', 'bgrid_w_velocity', 'bgrid_tracer']:
xsi = 1
eta = 1
if time < grid.time[ti]:
ti -= 1
if grid.z4d:
if ti == len(grid.time)-1:
depth_vector = (1-xsi)*(1-eta) * grid.depth[-1, :, yi, xi] + \
xsi*(1-eta) * grid.depth[-1, :, yi, xi+1] + \
xsi*eta * grid.depth[-1, :, yi+1, xi+1] + \
(1-xsi)*eta * grid.depth[-1, :, yi+1, xi]
else:
dv2 = (1-xsi)*(1-eta) * grid.depth[ti:ti+2, :, yi, xi] + \
xsi*(1-eta) * grid.depth[ti:ti+2, :, yi, xi+1] + \
xsi*eta * grid.depth[ti:ti+2, :, yi+1, xi+1] + \
(1-xsi)*eta * grid.depth[ti:ti+2, :, yi+1, xi]
tt = (time-grid.time[ti]) / (grid.time[ti+1]-grid.time[ti])
assert tt >= 0 and tt <= 1, 'Vertical s grid is being wrongly interpolated in time'
depth_vector = dv2[0, :] * (1-tt) + dv2[1, :] * tt
else:
depth_vector = (1-xsi)*(1-eta) * grid.depth[:, yi, xi] + \
xsi*(1-eta) * grid.depth[:, yi, xi+1] + \
xsi*eta * grid.depth[:, yi+1, xi+1] + \
(1-xsi)*eta * grid.depth[:, yi+1, xi]
z = np.float32(z)
if depth_vector[-1] > depth_vector[0]:
depth_indices = depth_vector <= z
if z >= depth_vector[-1]:
zi = len(depth_vector) - 2
else:
zi = depth_indices.argmin() - 1 if z >= depth_vector[0] else 0
if z < depth_vector[zi]:
raise FieldOutOfBoundSurfaceError(0, 0, z, field=self)
elif z > depth_vector[zi+1]:
raise FieldOutOfBoundError(x, y, z, field=self)
else:
depth_indices = depth_vector >= z
if z <= depth_vector[-1]:
zi = len(depth_vector) - 2
else:
zi = depth_indices.argmin() - 1 if z <= depth_vector[0] else 0
if z > depth_vector[zi]:
raise FieldOutOfBoundSurfaceError(0, 0, z, field=self)
elif z < depth_vector[zi+1]:
raise FieldOutOfBoundError(x, y, z, field=self)
zeta = (z - depth_vector[zi]) / (depth_vector[zi+1]-depth_vector[zi])
return (zi, zeta)
def reconnect_bnd_indices(self, xi, yi, xdim, ydim, sphere_mesh):
if xi < 0:
if sphere_mesh:
xi = xdim-2
else:
xi = 0
if xi > xdim-2:
if sphere_mesh:
xi = 0
else:
xi = xdim-2
if yi < 0:
yi = 0
if yi > ydim-2:
yi = ydim-2
if sphere_mesh:
xi = xdim - xi
return xi, yi
def search_indices_rectilinear(self, x, y, z, ti=-1, time=-1, particle=None, search2D=False):
grid = self.grid
if grid.xdim > 1 and (not grid.zonal_periodic):
if x < grid.lonlat_minmax[0] or x > grid.lonlat_minmax[1]:
raise FieldOutOfBoundError(x, y, z, field=self)
if grid.ydim > 1 and (y < grid.lonlat_minmax[2] or y > grid.lonlat_minmax[3]):
raise FieldOutOfBoundError(x, y, z, field=self)
if grid.xdim > 1:
if grid.mesh != 'spherical':
lon_index = grid.lon < x
if lon_index.all():
xi = len(grid.lon) - 2
else:
xi = lon_index.argmin() - 1 if lon_index.any() else 0
xsi = (x-grid.lon[xi]) / (grid.lon[xi+1]-grid.lon[xi])
if xsi < 0:
xi -= 1
xsi = (x-grid.lon[xi]) / (grid.lon[xi+1]-grid.lon[xi])
elif xsi > 1:
xi += 1
xsi = (x-grid.lon[xi]) / (grid.lon[xi+1]-grid.lon[xi])
else:
lon_fixed = grid.lon.copy()
indices = lon_fixed >= lon_fixed[0]
if not indices.all():
lon_fixed[indices.argmin():] += 360
if x < lon_fixed[0]:
lon_fixed -= 360
lon_index = lon_fixed < x
if lon_index.all():
xi = len(lon_fixed) - 2
else:
xi = lon_index.argmin() - 1 if lon_index.any() else 0
xsi = (x-lon_fixed[xi]) / (lon_fixed[xi+1]-lon_fixed[xi])
if xsi < 0:
xi -= 1
xsi = (x-lon_fixed[xi]) / (lon_fixed[xi+1]-lon_fixed[xi])
elif xsi > 1:
xi += 1
xsi = (x-lon_fixed[xi]) / (lon_fixed[xi+1]-lon_fixed[xi])
else:
xi, xsi = -1, 0
if grid.ydim > 1:
lat_index = grid.lat < y
if lat_index.all():
yi = len(grid.lat) - 2
else:
yi = lat_index.argmin() - 1 if lat_index.any() else 0
eta = (y-grid.lat[yi]) / (grid.lat[yi+1]-grid.lat[yi])
if eta < 0:
yi -= 1
eta = (y-grid.lat[yi]) / (grid.lat[yi+1]-grid.lat[yi])
elif eta > 1:
yi += 1
eta = (y-grid.lat[yi]) / (grid.lat[yi+1]-grid.lat[yi])
else:
yi, eta = -1, 0
if grid.zdim > 1 and not search2D:
if grid.gtype == GridCode.RectilinearZGrid:
# Never passes here, because in this case, we work with scipy
try:
(zi, zeta) = self.search_indices_vertical_z(z)
except FieldOutOfBoundError:
raise FieldOutOfBoundError(x, y, z, field=self)
except FieldOutOfBoundSurfaceError:
raise FieldOutOfBoundSurfaceError(x, y, z, field=self)
elif grid.gtype == GridCode.RectilinearSGrid:
(zi, zeta) = self.search_indices_vertical_s(x, y, z, xi, yi, xsi, eta, ti, time)
else:
zi, zeta = -1, 0
if not ((0 <= xsi <= 1) and (0 <= eta <= 1) and (0 <= zeta <= 1)):
raise FieldSamplingError(x, y, z, field=self)
if particle:
particle.xi[self.igrid] = xi
particle.yi[self.igrid] = yi
particle.zi[self.igrid] = zi
return (xsi, eta, zeta, xi, yi, zi)
def search_indices_curvilinear(self, x, y, z, ti=-1, time=-1, particle=None, search2D=False):
if particle:
xi = particle.xi[self.igrid]
yi = particle.yi[self.igrid]
else:
xi = int(self.grid.xdim / 2) - 1
yi = int(self.grid.ydim / 2) - 1
xsi = eta = -1
grid = self.grid
invA = np.array([[1, 0, 0, 0],
[-1, 1, 0, 0],
[-1, 0, 0, 1],
[1, -1, 1, -1]])
maxIterSearch = 1e6
it = 0
tol = 1.e-10
if not grid.zonal_periodic:
if x < grid.lonlat_minmax[0] or x > grid.lonlat_minmax[1]:
if grid.lon[0, 0] < grid.lon[0, -1]:
raise FieldOutOfBoundError(x, y, z, field=self)
elif x < grid.lon[0, 0] and x > grid.lon[0, -1]: # This prevents from crashing in [160, -160]
raise FieldOutOfBoundError(x, y, z, field=self)
if y < grid.lonlat_minmax[2] or y > grid.lonlat_minmax[3]:
raise FieldOutOfBoundError(x, y, z, field=self)
while xsi < -tol or xsi > 1+tol or eta < -tol or eta > 1+tol:
px = np.array([grid.lon[yi, xi], grid.lon[yi, xi+1], grid.lon[yi+1, xi+1], grid.lon[yi+1, xi]])
if grid.mesh == 'spherical':
px[0] = px[0]+360 if px[0] < x-225 else px[0]
px[0] = px[0]-360 if px[0] > x+225 else px[0]
px[1:] = np.where(px[1:] - px[0] > 180, px[1:]-360, px[1:])
px[1:] = np.where(-px[1:] + px[0] > 180, px[1:]+360, px[1:])
py = np.array([grid.lat[yi, xi], grid.lat[yi, xi+1], grid.lat[yi+1, xi+1], grid.lat[yi+1, xi]])
a = np.dot(invA, px)
b = np.dot(invA, py)
aa = a[3]*b[2] - a[2]*b[3]
bb = a[3]*b[0] - a[0]*b[3] + a[1]*b[2] - a[2]*b[1] + x*b[3] - y*a[3]
cc = a[1]*b[0] - a[0]*b[1] + x*b[1] - y*a[1]
if abs(aa) < 1e-12: # Rectilinear cell, or quasi
eta = -cc / bb
else:
det2 = bb*bb-4*aa*cc
if det2 > 0: # so, if det is nan we keep the xsi, eta from previous iter
det = np.sqrt(det2)
eta = (-bb+det)/(2*aa)
if abs(a[1]+a[3]*eta) < 1e-12: # this happens when recti cell rotated of 90deg
xsi = ((y-py[0])/(py[1]-py[0]) + (y-py[3])/(py[2]-py[3])) * .5
else:
xsi = (x-a[0]-a[2]*eta) / (a[1]+a[3]*eta)
if xsi < 0 and eta < 0 and xi == 0 and yi == 0:
raise FieldOutOfBoundError(x, y, 0, field=self)
if xsi > 1 and eta > 1 and xi == grid.xdim-1 and yi == grid.ydim-1:
raise FieldOutOfBoundError(x, y, 0, field=self)
if xsi < -tol:
xi -= 1
elif xsi > 1+tol:
xi += 1
if eta < -tol:
yi -= 1
elif eta > 1+tol:
yi += 1
(xi, yi) = self.reconnect_bnd_indices(xi, yi, grid.xdim, grid.ydim, grid.mesh)
it += 1
if it > maxIterSearch:
print('Correct cell not found after %d iterations' % maxIterSearch)
raise FieldOutOfBoundError(x, y, 0, field=self)
xsi = max(0., xsi)
eta = max(0., eta)
xsi = min(1., xsi)
eta = min(1., eta)
if grid.zdim > 1 and not search2D:
if grid.gtype == GridCode.CurvilinearZGrid:
try:
(zi, zeta) = self.search_indices_vertical_z(z)
except FieldOutOfBoundError:
raise FieldOutOfBoundError(x, y, z, field=self)
elif grid.gtype == GridCode.CurvilinearSGrid:
(zi, zeta) = self.search_indices_vertical_s(x, y, z, xi, yi, xsi, eta, ti, time)
else:
zi = -1
zeta = 0
if not ((0 <= xsi <= 1) and (0 <= eta <= 1) and (0 <= zeta <= 1)):
raise FieldSamplingError(x, y, z, field=self)
if particle:
particle.xi[self.igrid] = xi
particle.yi[self.igrid] = yi
particle.zi[self.igrid] = zi
return (xsi, eta, zeta, xi, yi, zi)
def search_indices(self, x, y, z, ti=-1, time=-1, particle=None, search2D=False):
if self.grid.gtype in [GridCode.RectilinearSGrid, GridCode.RectilinearZGrid]:
return self.search_indices_rectilinear(x, y, z, ti, time, particle=particle, search2D=search2D)
else:
return self.search_indices_curvilinear(x, y, z, ti, time, particle=particle, search2D=search2D)
def interpolator2D(self, ti, z, y, x, particle=None):
(xsi, eta, _, xi, yi, _) = self.search_indices(x, y, z, particle=particle)
if self.interp_method == 'nearest':
xii = xi if xsi <= .5 else xi+1
yii = yi if eta <= .5 else yi+1
return self.data[ti, yii, xii]
elif self.interp_method in ['linear', 'bgrid_velocity']:
val = (1-xsi)*(1-eta) * self.data[ti, yi, xi] + \
xsi*(1-eta) * self.data[ti, yi, xi+1] + \
xsi*eta * self.data[ti, yi+1, xi+1] + \
(1-xsi)*eta * self.data[ti, yi+1, xi]
return val
elif self.interp_method == 'linear_invdist_land_tracer':
land = np.isclose(self.data[ti, yi:yi+2, xi:xi+2], 0.)
nb_land = np.sum(land)
if nb_land == 4:
return 0
elif nb_land > 0:
val = 0
w_sum = 0
for j in range(2):
for i in range(2):
distance = pow((eta - j), 2) + pow((xsi - i), 2)
if np.isclose(distance, 0):
if land[j][i] == 1: # index search led us directly onto land
return 0
else:
return self.data[ti, yi+j, xi+i]
elif land[i][j] == 0:
val += self.data[ti, yi+j, xi+i] / distance
w_sum += 1 / distance
return val / w_sum
else:
val = (1 - xsi) * (1 - eta) * self.data[ti, yi, xi] + \
xsi * (1 - eta) * self.data[ti, yi, xi + 1] + \
xsi * eta * self.data[ti, yi + 1, xi + 1] + \
(1 - xsi) * eta * self.data[ti, yi + 1, xi]
return val
elif self.interp_method in ['cgrid_tracer', 'bgrid_tracer']:
return self.data[ti, yi+1, xi+1]
elif self.interp_method in ['cgrid_velocity']:
raise RuntimeError("%s is a scalar field. cgrid_velocity interpolation method should be used for vector fields (e.g. FieldSet.UV)" % self.name)
else:
raise RuntimeError(self.interp_method+" is not implemented for 2D grids")
def interpolator3D(self, ti, z, y, x, time, particle=None):
(xsi, eta, zeta, xi, yi, zi) = self.search_indices(x, y, z, ti, time, particle=particle)
if self.interp_method == 'nearest':
xii = xi if xsi <= .5 else xi+1
yii = yi if eta <= .5 else yi+1
zii = zi if zeta <= .5 else zi+1
return self.data[ti, zii, yii, xii]
elif self.interp_method == 'cgrid_velocity':
# evaluating W velocity in c_grid
if self.gridindexingtype == 'nemo':
f0 = self.data[ti, zi, yi+1, xi+1]
f1 = self.data[ti, zi+1, yi+1, xi+1]
elif self.gridindexingtype == 'mitgcm':
f0 = self.data[ti, zi, yi, xi]
f1 = self.data[ti, zi+1, yi, xi]
return (1-zeta) * f0 + zeta * f1
elif self.interp_method == 'linear_invdist_land_tracer':
land = np.isclose(self.data[ti, zi:zi+2, yi:yi+2, xi:xi+2], 0.)
nb_land = np.sum(land)
if nb_land == 8:
return 0
elif nb_land > 0:
val = 0
w_sum = 0
for k in range(2):
for j in range(2):
for i in range(2):
distance = pow((zeta - k), 2) + pow((eta - j), 2) + pow((xsi - i), 2)
if np.isclose(distance, 0):
if land[k][j][i] == 1: # index search led us directly onto land
return 0
else:
return self.data[ti, zi+i, yi+j, xi+k]
elif land[k][j][i] == 0:
val += self.data[ti, zi+k, yi+j, xi+i] / distance
w_sum += 1 / distance
return val / w_sum
else:
data = self.data[ti, zi, :, :]
f0 = (1 - xsi) * (1 - eta) * data[yi, xi] + \
xsi * (1 - eta) * data[yi, xi + 1] + \
xsi * eta * data[yi + 1, xi + 1] + \
(1 - xsi) * eta * data[yi + 1, xi]
data = self.data[ti, zi + 1, :, :]
f1 = (1 - xsi) * (1 - eta) * data[yi, xi] + \
xsi * (1 - eta) * data[yi, xi + 1] + \
xsi * eta * data[yi + 1, xi + 1] + \
(1 - xsi) * eta * data[yi + 1, xi]
return (1 - zeta) * f0 + zeta * f1
elif self.interp_method in ['linear', 'bgrid_velocity', 'bgrid_w_velocity']:
if self.interp_method == 'bgrid_velocity':
if self.gridindexingtype == 'mom5':
zeta = 1.
else:
zeta = 0.
elif self.interp_method == 'bgrid_w_velocity':
eta = 1.
xsi = 1.
data = self.data[ti, zi, :, :]
f0 = (1-xsi)*(1-eta) * data[yi, xi] + \
xsi*(1-eta) * data[yi, xi+1] + \
xsi*eta * data[yi+1, xi+1] + \
(1-xsi)*eta * data[yi+1, xi]
if self.gridindexingtype == 'pop' and zi >= self.grid.zdim-2:
# Since POP is indexed at cell top, allow linear interpolation of W to zero in lowest cell
return (1-zeta) * f0
data = self.data[ti, zi+1, :, :]
f1 = (1-xsi)*(1-eta) * data[yi, xi] + \
xsi*(1-eta) * data[yi, xi+1] + \
xsi*eta * data[yi+1, xi+1] + \
(1-xsi)*eta * data[yi+1, xi]
if self.interp_method == 'bgrid_w_velocity' and self.gridindexingtype == 'mom5' and zi == -1:
# Since MOM5 is indexed at cell bottom, allow linear interpolation of W to zero in uppermost cell
return zeta * f1
else: