Skip to content

Commit

Permalink
Merge pull request #22 from PlasmaFAIR/fix-arrays-without-grids
Browse files Browse the repository at this point in the history
Fix reading arrays without grids
  • Loading branch information
ZedThree authored Oct 17, 2024
2 parents 10e8a9d + 4ac7c36 commit 55515d3
Show file tree
Hide file tree
Showing 6 changed files with 159 additions and 2 deletions.
10 changes: 8 additions & 2 deletions src/sdf_xarray/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -275,12 +275,18 @@ def _grid_species_name(grid_name: str) -> str:
continue

if isinstance(value, Constant) or value.grid is None:
# No grid, so just a scalar
data_attrs = {}
if value.units is not None:
data_attrs["units"] = value.units

data_vars[key] = Variable((), value.data, data_attrs)
# We don't have a grid, either because it's just a
# scalar, or because it's an array over something
# else. We have no more information, so just make up
# some (hopefully) unique dimension names
shape = getattr(value.data, "shape", ())
dims = [f"dim_{key}_{n}" for n, _ in enumerate(shape)]

data_vars[key] = Variable(dims, value.data, attrs=data_attrs)
continue

if value.is_point_data:
Expand Down
Binary file added tests/example_array_no_grids/0000.sdf
Binary file not shown.
Binary file added tests/example_array_no_grids/0001.sdf
Binary file not shown.
9 changes: 9 additions & 0 deletions tests/example_array_no_grids/README.md
Original file line number Diff line number Diff line change
@@ -0,0 +1,9 @@
# Test files

These files were generated using EPOCH commit aafed39 and the
`input.deck` here, which is modified from a run by Joel Adams. It is
not intended to be physical. Run with:

echo "." | mpirun -np 2 $EPOCH1D

where `$EPOCH1D` is the path to the `epoch1d` executable
119 changes: 119 additions & 0 deletions tests/example_array_no_grids/input.deck
Original file line number Diff line number Diff line change
@@ -0,0 +1,119 @@
begin:constant
nel = 1e+23
intens = 2e+22
omega = 2.0 * pi * c / (1.0e-6)
den_crit = critical(omega)
scale = 1e-06
den_max = 5.0 * den_crit
den_maxpoint = 4e-05
den_contrast = 1.0
amax = 1.0
end:constant

begin:control
nx = 8
nparticles = nx * 2
nsteps = 1
t_end = 2e-15
x_min = -1e-05
x_max = 1e-05
dt_multiplier = 0.8
end:control

begin:qed
use_qed = F
qed_start_time = 0
produce_photons = F
photon_energy_min = 50 * kev
produce_pairs = F
photon_dynamics = F
end:qed

begin:collisions
use_collisions = T
coulomb_log = auto
collide = all
end:collisions

begin:boundaries
bc_x_min = simple_laser
bc_x_max = simple_laser
end:boundaries

begin:species
name = Electron
fraction = 0.5
dump = F
temperature = 0
number_density = if (x gt 0e-6, nel * 1.0e6, 0)
number_density_min = 1
identify:electron
end:species

begin:species
name = Ion
fraction = 0.5
dump = F
number_density = number_density(Electron)
temperature = temperature_x(Electron)
number_density_min = 1
identify:proton
end:species

begin:species
name = Photon
nparticles = 0
dump = F
identify:photon
end:species

begin:output_global
force_final_to_be_restartable = T
end:output_global

begin:output
name = normal
use_offset_grid = F
dt_snapshot = 5e-15
particles = never
px = never
py = never
pz = never
vx = never
vy = never
vz = never
charge = never
mass = never
particle_weight = never
species_id = never
grid = always
ex = never
ey = never
ez = never
bx = never
by = never
bz = never
jx = never
jy = never
jz = never
average_particle_energy = always + species
mass_density = never + species
charge_density = never
number_density = never + species
temperature = never + species
distribution_functions = always
particle_probes = never
absorption = always
total_energy_sum = always + species
end:output

begin:laser
boundary = x_min
intensity = intens * 1.0e4
omega = omega
polarisation = 0.0
phase = 0.0
t_profile = gauss(time, 40.0e-15, 30.0e-15)
t_start = 0.0
t_end = end
end:laser
23 changes: 23 additions & 0 deletions tests/test_basic.py
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,7 @@
EXAMPLE_MISMATCHED_FILES_DIR = (
pathlib.Path(__file__).parent / "example_mismatched_files"
)
EXAMPLE_ARRAYS_DIR = pathlib.Path(__file__).parent / "example_array_no_grids"


def test_basic():
Expand Down Expand Up @@ -101,3 +102,25 @@ def test_erroring_on_mismatched_jobid_files():
compat="override",
preprocess=SDFPreprocess(),
)


def test_arrays_with_no_grids():
with xr.open_dataset(EXAMPLE_ARRAYS_DIR / "0001.sdf") as df:
laser_phase = "laser_x_min_phase"
assert laser_phase in df
assert df[laser_phase].shape == (1,)

random_states = "Random States"
assert random_states in df
assert df[random_states].shape == (8,)


def test_arrays_with_no_grids_multifile():
df = xr.open_mfdataset(EXAMPLE_ARRAYS_DIR.glob("*.sdf"), preprocess=SDFPreprocess())
laser_phase = "laser_x_min_phase"
assert laser_phase in df
assert df[laser_phase].shape == (2, 1)

random_states = "Random States"
assert random_states in df
assert df[random_states].shape == (2, 8)

0 comments on commit 55515d3

Please sign in to comment.