-
Notifications
You must be signed in to change notification settings - Fork 271
/
construction_data_containers.py
2919 lines (2666 loc) · 106 KB
/
construction_data_containers.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 fileinput
import io
import os
import warnings
import zipfile
from collections.abc import Sized
from functools import wraps
from re import finditer
from tempfile import NamedTemporaryFile, TemporaryFile
import numpy as np
from more_itertools import always_iterable
from tqdm import tqdm
from yt.config import ytcfg
from yt.data_objects.field_data import YTFieldData
from yt.data_objects.selection_objects.data_selection_objects import (
YTSelectionContainer1D,
YTSelectionContainer2D,
YTSelectionContainer3D,
)
from yt.fields.field_exceptions import NeedsGridType, NeedsOriginalGrid
from yt.frontends.sph.data_structures import ParticleDataset
from yt.funcs import get_memory_usage, iter_fields, mylog, only_on_root
from yt.geometry import particle_deposit as particle_deposit
from yt.geometry.coordinates.cartesian_coordinates import all_data
from yt.loaders import load_uniform_grid
from yt.units.unit_object import Unit
from yt.units.yt_array import YTArray, uconcatenate
from yt.utilities.exceptions import (
YTNoAPIKey,
YTNotInsideNotebook,
YTParticleDepositionNotImplemented,
YTTooManyVertices,
)
from yt.utilities.grid_data_format.writer import write_to_gdf
from yt.utilities.lib.cyoctree import CyOctree
from yt.utilities.lib.interpolators import ghost_zone_interpolate
from yt.utilities.lib.marching_cubes import march_cubes_grid, march_cubes_grid_flux
from yt.utilities.lib.misc_utilities import fill_region, fill_region_float
from yt.utilities.lib.pixelization_routines import (
interpolate_sph_grid_gather,
interpolate_sph_positions_gather,
normalization_1d_utility,
normalization_3d_utility,
pixelize_sph_kernel_arbitrary_grid,
)
from yt.utilities.lib.quad_tree import QuadTree
from yt.utilities.minimal_representation import MinimalProjectionData
from yt.utilities.parallel_tools.parallel_analysis_interface import (
communication_system,
parallel_objects,
parallel_root_only,
)
class YTStreamline(YTSelectionContainer1D):
"""
This is a streamline, which is a set of points defined as
being parallel to some vector field.
This object is typically accessed through the Streamlines.path
function. The resulting arrays have their dimensionality
reduced to one, and an ordered list of points at an (x,y)
tuple along `axis` are available, as is the `t` field, which
corresponds to a unitless measurement along the ray from start
to end.
Parameters
----------
positions : array-like
List of streamline positions
length : float
The magnitude of the distance; dts will be divided by this
fields : list of strings, optional
If you want the object to pre-retrieve a set of fields, supply them
here. This is not necessary.
ds : dataset object
Passed in to access the index
kwargs : dict of items
Any additional values are passed as field parameters that can be
accessed by generated fields.
Examples
--------
>>> from yt.visualization.api import Streamlines
>>> streamlines = Streamlines(ds, [0.5] * 3)
>>> streamlines.integrate_through_volume()
>>> stream = streamlines.path(0)
>>> fig, ax = plt.subplots()
>>> ax.set_yscale("log")
>>> ax.plot(stream["t"], stream[("gas", "density")], "-x")
"""
_type_name = "streamline"
_con_args = ("positions",)
sort_by = "t"
def __init__(self, positions, length=1.0, fields=None, ds=None, **kwargs):
YTSelectionContainer1D.__init__(self, ds, fields, **kwargs)
self.positions = positions
self.dts = np.empty_like(positions[:, 0])
self.dts[:-1] = np.sqrt(
np.sum((self.positions[1:] - self.positions[:-1]) ** 2, axis=1)
)
self.dts[-1] = self.dts[-1]
self.length = length
self.dts /= length
self.ts = np.add.accumulate(self.dts)
self._set_center(self.positions[0])
self.set_field_parameter("center", self.positions[0])
self._dts, self._ts = {}, {}
# self._refresh_data()
def _get_list_of_grids(self):
# Get the value of the line at each LeftEdge and RightEdge
LE = self.ds.grid_left_edge
RE = self.ds.grid_right_edge
# Check left faces first
min_streampoint = np.min(self.positions, axis=0)
max_streampoint = np.max(self.positions, axis=0)
p = np.all((min_streampoint <= RE) & (max_streampoint > LE), axis=1)
self._grids = self.index.grids[p]
def _get_data_from_grid(self, grid, field):
# No child masking here; it happens inside the mask cut
mask = self._get_cut_mask(grid)
if field == "dts":
return self._dts[grid.id]
if field == "t":
return self._ts[grid.id]
return grid[field].flat[mask]
def _get_cut_mask(self, grid):
points_in_grid = np.all(self.positions > grid.LeftEdge, axis=1) & np.all(
self.positions <= grid.RightEdge, axis=1
)
pids = np.where(points_in_grid)[0]
mask = np.zeros(points_in_grid.sum(), dtype="int64")
dts = np.zeros(points_in_grid.sum(), dtype="float64")
ts = np.zeros(points_in_grid.sum(), dtype="float64")
for mi, (i, pos) in enumerate(zip(pids, self.positions[points_in_grid])):
if not points_in_grid[i]:
continue
ci = ((pos - grid.LeftEdge) / grid.dds).astype("int")
if grid.child_mask[ci[0], ci[1], ci[2]] == 0:
continue
for j in range(3):
ci[j] = min(ci[j], grid.ActiveDimensions[j] - 1)
mask[mi] = np.ravel_multi_index(ci, grid.ActiveDimensions)
dts[mi] = self.dts[i]
ts[mi] = self.ts[i]
self._dts[grid.id] = dts
self._ts[grid.id] = ts
return mask
class YTProj(YTSelectionContainer2D):
_key_fields = YTSelectionContainer2D._key_fields + ["weight_field"]
_con_args = ("axis", "field", "weight_field")
_container_fields = ("px", "py", "pdx", "pdy", "weight_field")
def __init__(
self,
field,
axis,
weight_field=None,
center=None,
ds=None,
data_source=None,
style=None,
method="integrate",
field_parameters=None,
max_level=None,
):
super().__init__(axis, ds, field_parameters)
# Style is deprecated, but if it is set, then it trumps method
# keyword. TODO: Remove this keyword and this check at some point in
# the future.
if style is not None:
method = style
if method == "sum":
self.method = "integrate"
self._sum_only = True
else:
self.method = method
self._sum_only = False
if self.method == "mip":
self.func = np.max
elif self.method == "integrate":
self.func = np.sum # for the future
else:
raise NotImplementedError(self.method)
self._set_center(center)
self._projected_units = {}
if data_source is None:
data_source = self.ds.all_data()
if max_level is not None:
data_source.max_level = max_level
for k, v in data_source.field_parameters.items():
if k not in self.field_parameters or self._is_default_field_parameter(k):
self.set_field_parameter(k, v)
self.data_source = data_source
if weight_field is None:
self.weight_field = weight_field
else:
self.weight_field = self._determine_fields(weight_field)[0]
for f in self._determine_fields(field):
nodal_flag = self.ds._get_field_info(f).nodal_flag
if any(nodal_flag):
raise RuntimeError(
"Nodal fields are currently not supported for projections."
)
@property
def blocks(self):
return self.data_source.blocks
@property
def field(self):
return [k for k in self.field_data.keys() if k not in self._container_fields]
def get_data(self, fields=None):
fields = self._determine_fields(fields)
# We need a new tree for every single set of fields we add
if len(fields) == 0:
return
if isinstance(self.ds, ParticleDataset):
return
tree = self._get_tree(len(fields))
# This only needs to be done if we are in parallel; otherwise, we can
# safely build the mesh as we go.
if communication_system.communicators[-1].size > 1:
for chunk in self.data_source.chunks([], "io", local_only=False):
self._initialize_chunk(chunk, tree)
_units_initialized = False
with self.data_source._field_parameter_state(self.field_parameters):
for chunk in parallel_objects(
self.data_source.chunks([], "io", local_only=True)
):
if not _units_initialized:
self._initialize_projected_units(fields, chunk)
_units_initialized = True
self._handle_chunk(chunk, fields, tree)
# if there's less than nprocs chunks, units won't be initialized
# on all processors, so sync with _projected_units on rank 0
projected_units = self.comm.mpi_bcast(self._projected_units)
self._projected_units = projected_units
# Note that this will briefly double RAM usage
if self.method == "mip":
merge_style = -1
op = "max"
elif self.method == "integrate":
merge_style = 1
op = "sum"
else:
raise NotImplementedError
# TODO: Add the combine operation
xax = self.ds.coordinates.x_axis[self.axis]
yax = self.ds.coordinates.y_axis[self.axis]
ox = self.ds.domain_left_edge[xax].v
oy = self.ds.domain_left_edge[yax].v
px, py, pdx, pdy, nvals, nwvals = tree.get_all(False, merge_style)
nvals = self.comm.mpi_allreduce(nvals, op=op)
nwvals = self.comm.mpi_allreduce(nwvals, op=op)
np.multiply(px, self.ds.domain_width[xax], px)
np.add(px, ox, px)
np.multiply(pdx, self.ds.domain_width[xax], pdx)
np.multiply(py, self.ds.domain_width[yax], py)
np.add(py, oy, py)
np.multiply(pdy, self.ds.domain_width[yax], pdy)
if self.weight_field is not None:
# If there are 0s remaining in the weight vals
# this will not throw an error, but silently
# return nans for vals where dividing by 0
# Leave as NaNs to be auto-masked by Matplotlib
with np.errstate(invalid="ignore"):
np.divide(nvals, nwvals[:, None], nvals)
# We now convert to half-widths and center-points
data = {}
code_length = self.ds.domain_width.units
data["px"] = self.ds.arr(px, code_length)
data["py"] = self.ds.arr(py, code_length)
data["weight_field"] = nwvals
data["pdx"] = self.ds.arr(pdx, code_length)
data["pdy"] = self.ds.arr(pdy, code_length)
data["fields"] = nvals
# Now we run the finalizer, which is ignored if we don't need it
field_data = np.hsplit(data.pop("fields"), len(fields))
for fi, field in enumerate(fields):
mylog.debug("Setting field %s", field)
input_units = self._projected_units[field]
self[field] = self.ds.arr(field_data[fi].ravel(), input_units)
for i in list(data.keys()):
self[i] = data.pop(i)
mylog.info("Projection completed")
self.tree = tree
def to_pw(self, fields=None, center="c", width=None, origin="center-window"):
r"""Create a :class:`~yt.visualization.plot_window.PWViewerMPL` from this
object.
This is a bare-bones mechanism of creating a plot window from this
object, which can then be moved around, zoomed, and on and on. All
behavior of the plot window is relegated to that routine.
"""
pw = self._get_pw(fields, center, width, origin, "Projection")
return pw
def plot(self, fields=None):
if hasattr(self.data_source, "left_edge") and hasattr(
self.data_source, "right_edge"
):
left_edge = self.data_source.left_edge
right_edge = self.data_source.right_edge
center = (left_edge + right_edge) / 2.0
width = right_edge - left_edge
xax = self.ds.coordinates.x_axis[self.axis]
yax = self.ds.coordinates.y_axis[self.axis]
lx, rx = left_edge[xax], right_edge[xax]
ly, ry = left_edge[yax], right_edge[yax]
width = (rx - lx), (ry - ly)
else:
width = self.ds.domain_width
center = self.ds.domain_center
pw = self._get_pw(fields, center, width, "native", "Projection")
try:
pw.show()
except YTNotInsideNotebook:
pass
return pw
def _initialize_projected_units(self, fields, chunk):
for field in self.data_source._determine_fields(fields):
if field in self._projected_units:
continue
finfo = self.ds._get_field_info(*field)
if finfo.units is None:
# First time calling a units="auto" field, infer units and cache
# for future field accesses.
finfo.units = str(chunk[field].units)
field_unit = Unit(finfo.output_units, registry=self.ds.unit_registry)
if self.method == "mip" or self._sum_only:
path_length_unit = Unit(registry=self.ds.unit_registry)
else:
ax_name = self.ds.coordinates.axis_name[self.axis]
path_element_name = ("index", f"path_element_{ax_name}")
path_length_unit = self.ds.field_info[path_element_name].units
path_length_unit = Unit(
path_length_unit, registry=self.ds.unit_registry
)
# Only convert to appropriate unit system for path
# elements that aren't angles
if not path_length_unit.is_dimensionless:
path_length_unit = path_length_unit.get_base_equivalent(
unit_system=self.ds.unit_system
)
if self.weight_field is None:
self._projected_units[field] = field_unit * path_length_unit
else:
self._projected_units[field] = field_unit
class YTParticleProj(YTProj):
"""
A projection operation optimized for SPH particle data.
"""
_type_name = "particle_proj"
def __init__(
self,
field,
axis,
weight_field=None,
center=None,
ds=None,
data_source=None,
style=None,
method="integrate",
field_parameters=None,
max_level=None,
):
super().__init__(
field,
axis,
weight_field,
center,
ds,
data_source,
style,
method,
field_parameters,
max_level,
)
def _handle_chunk(self, chunk, fields, tree):
raise NotImplementedError("Particle projections have not yet been implemented")
class YTQuadTreeProj(YTProj):
"""
This is a data object corresponding to a line integral through the
simulation domain.
This object is typically accessed through the `proj` object that
hangs off of index objects. YTQuadTreeProj is a projection of a
`field` along an `axis`. The field can have an associated
`weight_field`, in which case the values are multiplied by a weight
before being summed, and then divided by the sum of that weight; the
two fundamental modes of operating are direct line integral (no
weighting) and average along a line of sight (weighting.) What makes
`proj` different from the standard projection mechanism is that it
utilizes a quadtree data structure, rather than the old mechanism for
projections. It will not run in parallel, but serial runs should be
substantially faster. Note also that lines of sight are integrated at
every projected finest-level cell.
Parameters
----------
field : string
This is the field which will be "projected" along the axis. If
multiple are specified (in a list) they will all be projected in
the first pass.
axis : int
The axis along which to slice. Can be 0, 1, or 2 for x, y, z.
weight_field : string
If supplied, the field being projected will be multiplied by this
weight value before being integrated, and at the conclusion of the
projection the resultant values will be divided by the projected
`weight_field`.
center : array_like, optional
The 'center' supplied to fields that use it. Note that this does
not have to have `coord` as one value. Strictly optional.
data_source : `yt.data_objects.data_containers.YTSelectionContainer`, optional
If specified, this will be the data source used for selecting
regions to project.
method : string, optional
The method of projection to be performed.
"integrate" : integration along the axis
"mip" : maximum intensity projection
"sum" : same as "integrate", except that we don't multiply by the path length
WARNING: The "sum" option should only be used for uniform resolution grid
datasets, as other datasets may result in unphysical images.
style : string, optional
The same as the method keyword. Deprecated as of version 3.0.2.
Please use method keyword instead.
field_parameters : dict of items
Values to be passed as field parameters that can be
accessed by generated fields.
Examples
--------
>>> ds = load("RedshiftOutput0005")
>>> prj = ds.proj(("gas", "density"), 0)
>>> print(proj[("gas", "density")])
"""
_type_name = "quad_proj"
def __init__(
self,
field,
axis,
weight_field=None,
center=None,
ds=None,
data_source=None,
style=None,
method="integrate",
field_parameters=None,
max_level=None,
):
super().__init__(
field,
axis,
weight_field,
center,
ds,
data_source,
style,
method,
field_parameters,
max_level,
)
if not self.deserialize(field):
self.get_data(field)
self.serialize()
@property
def _mrep(self):
return MinimalProjectionData(self)
def deserialize(self, fields):
if not ytcfg.get("yt", "serialize"):
return False
for field in fields:
self[field] = None
deserialized_successfully = False
store_file = self.ds.parameter_filename + ".yt"
if os.path.isfile(store_file):
deserialized_successfully = self._mrep.restore(store_file, self.ds)
if deserialized_successfully:
mylog.info("Using previous projection data from %s", store_file)
for field, field_data in self._mrep.field_data.items():
self[field] = field_data
if not deserialized_successfully:
for field in fields:
del self[field]
return deserialized_successfully
def serialize(self):
if not ytcfg.get("yt", "serialize"):
return
self._mrep.store(self.ds.parameter_filename + ".yt")
def _get_tree(self, nvals):
xax = self.ds.coordinates.x_axis[self.axis]
yax = self.ds.coordinates.y_axis[self.axis]
xd = self.ds.domain_dimensions[xax]
yd = self.ds.domain_dimensions[yax]
bounds = (
self.ds.domain_left_edge[xax],
self.ds.domain_right_edge[xax],
self.ds.domain_left_edge[yax],
self.ds.domain_right_edge[yax],
)
return QuadTree(
np.array([xd, yd], dtype="int64"), nvals, bounds, method=self.method
)
def _initialize_chunk(self, chunk, tree):
icoords = chunk.icoords
xax = self.ds.coordinates.x_axis[self.axis]
yax = self.ds.coordinates.y_axis[self.axis]
i1 = icoords[:, xax]
i2 = icoords[:, yax]
ilevel = chunk.ires * self.ds.ires_factor
tree.initialize_chunk(i1, i2, ilevel)
def _handle_chunk(self, chunk, fields, tree):
mylog.debug(
"Adding chunk (%s) to tree (%0.3e GB RAM)",
chunk.ires.size,
get_memory_usage() / 1024.0,
)
if self.method == "mip" or self._sum_only:
dl = self.ds.quan(1.0, "")
else:
ax_name = self.ds.coordinates.axis_name[self.axis]
dl = chunk["index", f"path_element_{ax_name}"]
# This is done for cases where our path element does not have a CGS
# equivalent. Once "preferred units" have been implemented, this
# will not be necessary at all, as the final conversion will occur
# at the display layer.
if not dl.units.is_dimensionless:
dl.convert_to_units(self.ds.unit_system["length"])
v = np.empty((chunk.ires.size, len(fields)), dtype="float64")
for i, field in enumerate(fields):
d = chunk[field] * dl
v[:, i] = d
if self.weight_field is not None:
w = chunk[self.weight_field]
np.multiply(v, w[:, None], v)
np.multiply(w, dl, w)
else:
w = np.ones(chunk.ires.size, dtype="float64")
icoords = chunk.icoords
xax = self.ds.coordinates.x_axis[self.axis]
yax = self.ds.coordinates.y_axis[self.axis]
i1 = icoords[:, xax]
i2 = icoords[:, yax]
ilevel = chunk.ires * self.ds.ires_factor
tree.add_chunk_to_tree(i1, i2, ilevel, v, w)
class YTCoveringGrid(YTSelectionContainer3D):
"""A 3D region with all data extracted to a single, specified
resolution. Left edge should align with a cell boundary, but
defaults to the closest cell boundary.
Parameters
----------
level : int
The resolution level data to which data will be gridded. Level
0 is the root grid dx for that dataset.
left_edge : array_like
The left edge of the region to be extracted. Specify units by supplying
a YTArray, otherwise code length units are assumed.
dims : array_like
Number of cells along each axis of resulting covering_grid
fields : array_like, optional
A list of fields that you'd like pre-generated for your object
num_ghost_zones : integer, optional
The number of padding ghost zones used when accessing fields.
Examples
--------
>>> cube = ds.covering_grid(2, left_edge=[0.0, 0.0, 0.0], dims=[128, 128, 128])
"""
_spatial = True
_type_name = "covering_grid"
_con_args = ("level", "left_edge", "ActiveDimensions")
_container_fields = (
("index", "dx"),
("index", "dy"),
("index", "dz"),
("index", "x"),
("index", "y"),
("index", "z"),
)
_base_grid = None
def __init__(
self,
level,
left_edge,
dims,
fields=None,
ds=None,
num_ghost_zones=0,
use_pbar=True,
field_parameters=None,
):
if field_parameters is None:
center = None
else:
center = field_parameters.get("center", None)
YTSelectionContainer3D.__init__(self, center, ds, field_parameters)
self.level = level
self.left_edge = self._sanitize_edge(left_edge)
self.ActiveDimensions = self._sanitize_dims(dims)
rdx = self.ds.domain_dimensions * self.ds.relative_refinement(0, level)
# normalize dims as a non-zero dim array
dims = np.array(list(always_iterable(dims)))
rdx[np.where(dims - 2 * num_ghost_zones <= 1)] = 1 # issue 602
self.base_dds = self.ds.domain_width / self.ds.domain_dimensions
self.dds = self.ds.domain_width / rdx.astype("float64")
self.right_edge = self.left_edge + self.ActiveDimensions * self.dds
self._num_ghost_zones = num_ghost_zones
self._use_pbar = use_pbar
self.global_startindex = (
np.rint((self.left_edge - self.ds.domain_left_edge) / self.dds).astype(
"int64"
)
+ self.ds.domain_offset
)
self._setup_data_source()
self.get_data(fields)
def get_global_startindex(self):
r"""Get the global start index of the covering grid."""
return self.global_startindex
def to_xarray(self, fields=None):
r"""Export this fixed-resolution object to an xarray Dataset
This function will take a regularized grid and optionally a list of
fields and return an xarray Dataset object. If xarray is not
importable, this will raise ImportError.
Parameters
----------
fields : list of strings or tuple field names, default None
If this is supplied, it is the list of fields to be exported into
the data frame. If not supplied, whatever fields presently exist
will be used.
Returns
-------
arr : Dataset
The data contained in the object.
Examples
--------
>>> dd = ds.r[::256j, ::256j, ::256j]
>>> xf1 = dd.to_xarray([("gas", "density"), ("gas", "temperature")])
>>> dd[("gas", "velocity_magnitude")]
>>> xf2 = dd.to_xarray()
"""
import xarray as xr
data = {}
coords = {}
for f in fields or self.field_data.keys():
data[f] = {
"dims": (
"x",
"y",
"z",
),
"data": self[f],
"attrs": {"units": str(self[f].uq)},
}
# We have our data, so now we generate both our coordinates and our metadata.
LE = self.LeftEdge + self.dds / 2.0
RE = self.RightEdge - self.dds / 2.0
N = self.ActiveDimensions
u = str(LE.uq)
for i, ax in enumerate("xyz"):
coords[ax] = {
"dims": (ax,),
"data": np.mgrid[LE[i] : RE[i] : N[i] * 1j],
"attrs": {"units": u},
}
return xr.Dataset.from_dict({"data_vars": data, "coords": coords})
@property
def icoords(self):
ic = np.indices(self.ActiveDimensions).astype("int64")
return np.column_stack(
[i.ravel() + gi for i, gi in zip(ic, self.get_global_startindex())]
)
@property
def fwidth(self):
fw = np.ones((self.ActiveDimensions.prod(), 3), dtype="float64")
fw *= self.dds
return fw
@property
def fcoords(self):
LE = self.LeftEdge + self.dds / 2.0
RE = self.RightEdge - self.dds / 2.0
N = self.ActiveDimensions
fc = np.mgrid[
LE[0] : RE[0] : N[0] * 1j,
LE[1] : RE[1] : N[1] * 1j,
LE[2] : RE[2] : N[2] * 1j,
]
return np.column_stack([f.ravel() for f in fc])
@property
def ires(self):
tr = np.ones(self.ActiveDimensions.prod(), dtype="int64")
tr *= self.level
return tr
def set_field_parameter(self, name, val):
super().set_field_parameter(name, val)
if self._data_source is not None:
self._data_source.set_field_parameter(name, val)
def _sanitize_dims(self, dims):
if not isinstance(dims, Sized):
dims = [dims] * len(self.ds.domain_left_edge)
if len(dims) != len(self.ds.domain_left_edge):
raise RuntimeError(
"Length of dims must match the dimensionality of the dataset"
)
return np.array(dims, dtype="int32")
def _sanitize_edge(self, edge):
if not isinstance(edge, Sized):
edge = [edge] * len(self.ds.domain_left_edge)
if len(edge) != len(self.ds.domain_left_edge):
raise RuntimeError(
"Length of edges must match the dimensionality of the dataset"
)
if hasattr(edge, "units"):
if edge.units.registry is self.ds.unit_registry:
return edge
edge_units = edge.units.copy()
edge_units.registry = self.ds.unit_registry
else:
edge_units = "code_length"
return self.ds.arr(edge, edge_units, dtype="float64")
def _reshape_vals(self, arr):
if len(arr.shape) == 3:
return arr
return arr.reshape(self.ActiveDimensions, order="C")
@property
def shape(self):
return tuple(self.ActiveDimensions.tolist())
def _setup_data_source(self):
self._data_source = self.ds.region(self.center, self.left_edge, self.right_edge)
self._data_source.min_level = 0
self._data_source.max_level = self.level
# This triggers "special" behavior in the RegionSelector to ensure we
# select *cells* whose bounding boxes overlap with our region, not just
# their cell centers.
self._data_source.loose_selection = True
def get_data(self, fields=None):
if fields is None:
return
fields = self._determine_fields(fields)
fields_to_get = [f for f in fields if f not in self.field_data]
fields_to_get = self._identify_dependencies(fields_to_get)
if len(fields_to_get) == 0:
return
try:
fill, gen, part, alias = self._split_fields(fields_to_get)
except NeedsGridType as e:
if self._num_ghost_zones == 0:
raise RuntimeError(
"Attempting to access a field that needs ghost zones, but "
"num_ghost_zones = %s. You should create the covering grid "
"with nonzero num_ghost_zones." % self._num_ghost_zones
) from e
else:
raise
# checking if we have a sph particles
if len(part) == 0:
is_sph_field = False
else:
is_sph_field = self.ds.field_info[part[0]].is_sph_field
if len(part) > 0 and len(alias) == 0:
if is_sph_field:
self._fill_sph_particles(fields)
for field in fields:
if field in gen:
gen.remove(field)
else:
self._fill_particles(part)
if len(fill) > 0:
self._fill_fields(fill)
for a, f in sorted(alias.items()):
if f.sampling_type == "particle" and not is_sph_field:
self[a] = self._data_source[f]
else:
self[a] = f(self)
self.field_data[a].convert_to_units(f.output_units)
if len(gen) > 0:
part_gen = []
cell_gen = []
for field in gen:
finfo = self.ds.field_info[field]
if finfo.sampling_type == "particle":
part_gen.append(field)
else:
cell_gen.append(field)
self._generate_fields(cell_gen)
for p in part_gen:
self[p] = self._data_source[p]
def _split_fields(self, fields_to_get):
fill, gen = self.index._split_fields(fields_to_get)
particles = []
alias = {}
for field in gen:
finfo = self.ds._get_field_info(*field)
if finfo._function.__name__ == "_TranslationFunc":
alias[field] = finfo
continue
try:
finfo.check_available(self)
except NeedsOriginalGrid:
fill.append(field)
for field in fill:
finfo = self.ds._get_field_info(*field)
if finfo.sampling_type == "particle":
particles.append(field)
gen = [f for f in gen if f not in fill and f not in alias]
fill = [f for f in fill if f not in particles]
return fill, gen, particles, alias
def _fill_particles(self, part):
for p in part:
self[p] = self._data_source[p]
def _fill_sph_particles(self, fields):
# checks that we have the field and gets information
fields = [f for f in fields if f not in self.field_data]
if len(fields) == 0:
return
smoothing_style = getattr(self.ds, "sph_smoothing_style", "scatter")
normalize = getattr(self.ds, "use_sph_normalization", True)
bounds, size = self._get_grid_bounds_size()
period = self.ds.coordinates.period.copy()
if hasattr(period, "in_units"):
period = period.in_units("code_length").d
# TODO maybe there is a better way of handling this
is_periodic = int(any(self.ds.periodicity))
if smoothing_style == "scatter":
for field in fields:
fi = self.ds._get_field_info(field)
ptype = fi.name[0]
if ptype not in self.ds._sph_ptypes:
raise KeyError(f"{ptype} is not a SPH particle type!")
buff = np.zeros(size, dtype="float64")
if normalize:
buff_den = np.zeros(size, dtype="float64")
pbar = tqdm(desc=f"Interpolating SPH field {field}")
for chunk in self._data_source.chunks([field], "io"):
px = chunk[(ptype, "particle_position_x")].in_base("code").d
py = chunk[(ptype, "particle_position_y")].in_base("code").d
pz = chunk[(ptype, "particle_position_z")].in_base("code").d
hsml = chunk[(ptype, "smoothing_length")].in_base("code").d
mass = chunk[(ptype, "particle_mass")].in_base("code").d
dens = chunk[(ptype, "density")].in_base("code").d
field_quantity = chunk[field].d
pixelize_sph_kernel_arbitrary_grid(
buff,
px,
py,
pz,
hsml,
mass,
dens,
field_quantity,
bounds,
pbar=pbar,
check_period=is_periodic,
period=period,
)
if normalize:
pixelize_sph_kernel_arbitrary_grid(
buff_den,
px,
py,
pz,
hsml,
mass,
dens,
np.ones(dens.shape[0]),
bounds,
pbar=pbar,
check_period=is_periodic,
period=period,
)
if normalize:
normalization_3d_utility(buff, buff_den)
self[field] = self.ds.arr(buff, fi.units)
pbar.close()
if smoothing_style == "gather":
num_neighbors = getattr(self.ds, "num_neighbors", 32)
for field in fields:
buff = np.zeros(size, dtype="float64")
fields_to_get = [
"particle_position",
"density",
"particle_mass",
"smoothing_length",
field[1],
]
all_fields = all_data(self.ds, field[0], fields_to_get, kdtree=True)
fi = self.ds._get_field_info(field)
interpolate_sph_grid_gather(
buff,
all_fields["particle_position"],
bounds,
all_fields["smoothing_length"],
all_fields["particle_mass"],
all_fields["density"],
all_fields[field[1]].in_units(fi.units),
self.ds.index.kdtree,
use_normalization=normalize,
num_neigh=num_neighbors,
)
self[field] = self.ds.arr(buff, fi.units)
def _fill_fields(self, fields):
fields = [f for f in fields if f not in self.field_data]
if len(fields) == 0:
return
output_fields = [
np.zeros(self.ActiveDimensions, dtype="float64") for field in fields
]
domain_dims = self.ds.domain_dimensions.astype(
"int64"
) * self.ds.relative_refinement(0, self.level)
refine_by = self.ds.refine_by
if not isinstance(self.ds.refine_by, Sized):
refine_by = [refine_by, refine_by, refine_by]
refine_by = np.array(refine_by, dtype="i8")
for chunk in parallel_objects(self._data_source.chunks(fields, "io")):
input_fields = [chunk[field] for field in fields]
# NOTE: This usage of "refine_by" is actually *okay*, because it's
# being used with respect to iref, which is *already* scaled!