-
Notifications
You must be signed in to change notification settings - Fork 8
Expand file tree
/
Copy pathtest_circle.py
More file actions
355 lines (300 loc) Β· 12.1 KB
/
test_circle.py
File metadata and controls
355 lines (300 loc) Β· 12.1 KB
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
import textwrap
from copy import deepcopy
from datetime import datetime
import numpy as np
import pandas as pd
import pytest
import xugrid as xu
import imod
from imod.logging import LoggerType, LogLevel
from imod.mf6.validation_settings import ValidationSettings
from imod.mf6.write_context import WriteContext
def test_simulation_write_and_run(circle_model, tmp_path):
imod.logging.configure(
LoggerType.PYTHON,
log_level=LogLevel.DEBUG,
add_default_file_handler=True,
add_default_stream_handler=True,
)
simulation = circle_model
with pytest.raises(
RuntimeError, match="Simulation circle has not been written yet."
):
simulation.run()
modeldir = tmp_path / "circle"
mask = deepcopy(simulation["GWF_1"].domain)
mask.values[:, 133] = -1
simulation.mask_all_models(mask)
simulation.write(modeldir, binary=False, use_absolute_paths=True)
simulation.run()
# Test opening with different time unit, see if time unit is correctly
# converted to seconds
head = imod.mf6.open_hds(
modeldir / "GWF_1/GWF_1.hds",
modeldir / "GWF_1/disv.disv.grb",
simulation_start_time="01-01-1999",
time_unit="s",
)
assert isinstance(head, xu.UgridDataArray)
assert np.all(
head["time"].values[0]
== np.array(datetime(1999, 1, 1, 0, 0, 1), dtype="datetime64[ns]")
)
assert head.dims == ("time", "layer", "mesh2d_nFaces")
assert head.shape == (2, 2, 216)
def test_gwfmodel_render(circle_model, tmp_path):
simulation = circle_model
globaltimes = simulation["time_discretization"]["time"].values
gwfmodel = simulation["GWF_1"]
write_context1 = WriteContext()
actual = gwfmodel._render("GWF_1", write_context1)
path = "GWF_1"
expected = textwrap.dedent(
f"""\
begin options
end options
begin packages
disv6 {path}/disv.disv disv
chd6 {path}/chd.chd chd
ic6 {path}/ic.ic ic
npf6 {path}/npf.npf npf
sto6 {path}/sto.sto sto
oc6 {path}/oc.oc oc
rch6 {path}/rch.rch rch
end packages
"""
)
assert actual == expected
validation_context = ValidationSettings(True)
write_context2 = WriteContext(tmp_path)
gwfmodel._write("GWF_1", globaltimes, write_context2, validation_context)
assert (tmp_path / "GWF_1" / "GWF_1.nam").is_file()
assert (tmp_path / "GWF_1").is_dir()
def test_simulation_write_and_run_evt(circle_model_evt, tmp_path):
simulation = circle_model_evt
with pytest.raises(
RuntimeError, match="Simulation circle has not been written yet."
):
simulation.run()
modeldir = tmp_path / "circle_evt"
simulation.write(modeldir, binary=True)
simulation.run()
head = imod.mf6.open_hds(
modeldir / "GWF_1/GWF_1.hds", modeldir / "GWF_1/disv.disv.grb"
)
assert isinstance(head, xu.UgridDataArray)
assert head.dims == ("time", "layer", "mesh2d_nFaces")
assert head.shape == (2, 2, 216)
def test_simulation_write_and_run_evt__no_segments(circle_model_evt, tmp_path):
simulation = circle_model_evt
evt = simulation["GWF_1"].pop("evt")
evt_ds = evt.dataset
evt_ds = evt_ds.drop_vars(["proportion_rate", "proportion_depth"])
evt_dict = {key: evt_ds[key] for key in evt_ds.keys()}
simulation["GWF_1"]["evt"] = imod.mf6.Evapotranspiration(**evt_dict)
with pytest.raises(
RuntimeError, match="Simulation circle has not been written yet."
):
simulation.run()
modeldir = tmp_path / "circle_evt"
simulation.write(modeldir, binary=True)
simulation.run()
head = imod.mf6.open_hds(
modeldir / "GWF_1/GWF_1.hds", modeldir / "GWF_1/disv.disv.grb"
)
assert isinstance(head, xu.UgridDataArray)
assert head.dims == ("time", "layer", "mesh2d_nFaces")
assert head.shape == (2, 2, 216)
def test_gwfmodel_render_evt(circle_model_evt, tmp_path):
simulation = circle_model_evt
globaltimes = simulation["time_discretization"]["time"].values
gwfmodel = simulation["GWF_1"]
write_context1 = WriteContext()
actual = gwfmodel._render("GWF_1", write_context1)
path = "GWF_1"
expected = textwrap.dedent(
f"""\
begin options
end options
begin packages
disv6 {path}/disv.disv disv
chd6 {path}/chd.chd chd
ic6 {path}/ic.ic ic
npf6 {path}/npf.npf npf
sto6 {path}/sto.sto sto
oc6 {path}/oc.oc oc
rch6 {path}/rch.rch rch
evt6 {path}/evt.evt evt
end packages
"""
)
assert actual == expected
validation_context = ValidationSettings(True)
write_context2 = WriteContext(tmp_path)
gwfmodel._write("GWF_1", globaltimes, write_context2, validation_context)
assert (tmp_path / "GWF_1" / "GWF_1.nam").is_file()
assert (tmp_path / "GWF_1").is_dir()
def test_simulation_write_and_run_transport(circle_model_transport, tmp_path):
simulation = circle_model_transport
with pytest.raises(
RuntimeError, match="Simulation circle has not been written yet."
):
simulation.run()
modeldir = tmp_path / "circle_transport"
simulation.write(modeldir)
simulation.run()
head = simulation.open_head()
concentration = simulation.open_concentration()
assert isinstance(head, xu.UgridDataArray)
assert head.dims == ("time", "layer", "mesh2d_nFaces")
assert head.shape == (2, 2, 216)
assert isinstance(concentration, xu.UgridDataArray)
assert concentration.dims == ("time", "layer", "mesh2d_nFaces")
assert concentration.shape == (2, 2, 216)
def test_simulation_write_and_run_transport_vsc(circle_model_transport_vsc, tmp_path):
simulation = circle_model_transport_vsc
with pytest.raises(
RuntimeError, match="Simulation circle has not been written yet."
):
circle_model_transport_vsc.run()
modeldir = tmp_path / "circle_transport"
simulation.write(modeldir)
simulation.run()
head = simulation.open_head()
concentration = simulation.open_concentration()
assert isinstance(head, xu.UgridDataArray)
assert head.dims == ("time", "layer", "mesh2d_nFaces")
assert head.shape == (2, 2, 216)
assert isinstance(concentration, xu.UgridDataArray)
assert concentration.dims == ("time", "layer", "mesh2d_nFaces")
assert concentration.shape == (2, 2, 216)
def run_simulation(simulation, modeldir):
idomain = simulation["GWF_1"]["disv"]["idomain"].compute()
time = simulation["time_discretization"].dataset.coords["time"].values
simulation.write(modeldir)
simulation.run()
head = (
simulation.open_head(simulation_start_time=time[0])
.compute()
.reindex_like(idomain)
)
concentration = (
simulation.open_concentration()
.compute(simulation_start_time=time[0])
.reindex_like(idomain)
)
return head, concentration
def test_simulation_clip_and_state_at_boundary(circle_model_transport, tmp_path):
# Arrange
simulation = circle_model_transport
idomain = simulation["GWF_1"]["disv"]["idomain"].compute()
head, concentration = run_simulation(simulation, tmp_path / "full")
states_for_boundary = {
"GWF_1": head.isel(time=-1, drop=True),
"transport": concentration.isel(time=-1, drop=True),
}
# Act
half_simulation = simulation.clip_box(
x_max=0.1, states_for_boundary=states_for_boundary
)
# Assert
# Test if model dims halved
idomain_half = half_simulation["GWF_1"]["disv"]["idomain"]
dim = idomain_half.grid.face_dimension
np.testing.assert_array_equal(idomain_half.sizes[dim] / idomain.sizes[dim], 0.5)
assert (
half_simulation["transport"]["cnc_clipped"]["concentration"].notnull().sum()
== 20
)
assert half_simulation["GWF_1"]["chd_clipped"]["head"].notnull().sum() == 20
# Test if model runs
head_half, concentration_half = run_simulation(half_simulation, tmp_path / "half")
assert head_half.shape == (2, 2, 108)
assert concentration_half.shape == (2, 2, 108)
def test_simulation_clip_and_constant_state_at_boundary__transient_chd(
circle_model_transport, tmp_path
):
# Arrange
simulation = circle_model_transport
idomain = simulation["GWF_1"]["disv"]["idomain"].compute()
time = simulation["time_discretization"].dataset.coords["time"].values
simulation["GWF_1"]["chd"].dataset["head"] = (
simulation["GWF_1"]["chd"].dataset["head"].expand_dims(time=time)
)
head, concentration = run_simulation(simulation, tmp_path / "full")
states_for_boundary = {
"GWF_1": head.isel(time=-1, drop=True),
"transport": concentration.isel(time=-1, drop=True),
}
# Act
half_simulation = simulation.clip_box(
x_max=0.1, states_for_boundary=states_for_boundary
)
# Assert
# Test if model dims halved
idomain_half = half_simulation["GWF_1"]["disv"]["idomain"]
dim = idomain_half.grid.face_dimension
np.testing.assert_array_equal(idomain_half.sizes[dim] / idomain.sizes[dim], 0.5)
# Test if the clipped boundary conditions have the correct dimension and values
conc_clipped = half_simulation["transport"]["cnc_clipped"]["concentration"]
assert conc_clipped.sizes["layer"] == 2
assert "time" not in conc_clipped.dims
assert conc_clipped.notnull().sum() == 20
head_clipped = half_simulation["GWF_1"]["chd_clipped"]["head"]
assert head_clipped.sizes["layer"] == 2
assert head_clipped.sizes["time"] == 2
assert head_clipped.notnull().sum() == 20 * 2
# Test if model runs
head_half, concentration_half = run_simulation(half_simulation, tmp_path / "half")
assert head_half.shape == (2, 2, 108)
assert concentration_half.shape == (2, 2, 108)
def test_simulation_clip_and_transient_state_at_boundary__transient_chd(
circle_model_transport, tmp_path
):
"""
Test with a transient chd boundary condition and transient states for
boundary conditions. The time steps of the chd package are outside of the
time domain of the states_for_boundary.
Extend the time discretization as it was before
https://github.com/Deltares/imod-python/pull/1029, to allow clipping longer
timeseries.
"""
# Arrange
simulation = circle_model_transport
simulation.create_time_discretization(
additional_times=pd.date_range(start="2000-01-01", end="2001-01-01", freq="W")
)
idomain = simulation["GWF_1"]["disv"]["idomain"].compute()
time = simulation["time_discretization"].dataset.coords["time"].values
simulation["GWF_1"]["chd"].dataset["head"] = (
simulation["GWF_1"]["chd"].dataset["head"].expand_dims(time=time[:5])
)
head, concentration = run_simulation(simulation, tmp_path / "full")
states_for_boundary = {
"GWF_1": head.isel(time=slice(11, -11)),
"transport": concentration.isel(time=slice(11, -11)),
}
# Act
half_simulation = simulation.clip_box(
x_max=0.1, states_for_boundary=states_for_boundary
)
# Assert
# Test if model dims halved
idomain_half = half_simulation["GWF_1"]["disv"]["idomain"]
dim = idomain_half.grid.face_dimension
np.testing.assert_array_equal(idomain_half.sizes[dim] / idomain.sizes[dim], 0.5)
# Test if the clipped boundary conditions have the correct dimension and values
# Transport data for just 30 time steps (from 12 up to 41)
conc_clipped = half_simulation["transport"]["cnc_clipped"]["concentration"]
assert conc_clipped.sizes["layer"] == 2
assert conc_clipped.sizes["time"] == 30
assert conc_clipped.notnull().sum() == 20 * 30
# Head data for chd for just 30 time steps (from 12 up to 41)
head_clipped = half_simulation["GWF_1"]["chd_clipped"]["head"]
assert head_clipped.sizes["layer"] == 2
assert head_clipped.sizes["time"] == 30
assert head_clipped.notnull().sum() == 20 * 30
# Test if model runs
head_half, concentration_half = run_simulation(half_simulation, tmp_path / "half")
assert head_half.shape == (52, 2, 108)
assert concentration_half.shape == (52, 2, 108)