/
model.py
3238 lines (2512 loc) · 116 KB
/
model.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
'''This file is part of AeoLiS.
AeoLiS is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
AeoLiS is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with AeoLiS. If not, see <http://www.gnu.org/licenses/>.
AeoLiS Copyright (C) 2015 Bas Hoonhout
bas.hoonhout@deltares.nl b.m.hoonhout@tudelft.nl
Deltares Delft University of Technology
Unit of Hydraulic Engineering Faculty of Civil Engineering and Geosciences
Boussinesqweg 1 Stevinweg 1
2629 HVDelft 2628CN Delft
The Netherlands The Netherlands
'''
from __future__ import absolute_import, division
import os
import imp
import sys
import time
import glob
import logging
import warnings
import operator
import numpy as np
import scipy.sparse
import pickle
import scipy.sparse.linalg
import matplotlib.pyplot as plt
from datetime import timedelta
from bmi.api import IBmi
from functools import reduce
# package modules
import aeolis.inout
import aeolis.bed
import aeolis.avalanching
import aeolis.wind
import aeolis.threshold
import aeolis.transport
import aeolis.hydro
import aeolis.netcdf
import aeolis.constants
import aeolis.erosion
import aeolis.vegetation
import aeolis.fences
import aeolis.gridparams
from aeolis.utils import *
class StreamFormatter(logging.Formatter):
def format(self, record):
if record.levelname == 'INFO':
return record.getMessage()
else:
return '%s: %s' % (record.levelname, record.getMessage())
# initialize logger
logger = logging.getLogger(__name__)
__version__ = ''
__root__ = os.path.dirname(__file__)
try:
__version__ = open(os.path.join(__root__, 'VERSION')).read().strip()
except:
logger.warning('WARNING: Unknown model version.')
class ModelState(dict):
'''Dictionary-like object to store model state
Model state variables are mutable by default, but can be set
immutable. In the latter case any actions that set the immutable
model state variable are ignored.
'''
def __init__(self, *args, **kwargs):
self.ismutable = set()
super(ModelState, self).__init__(*args, **kwargs)
def __setitem__(self, k, v):
if k not in self.keys() or k in self.ismutable:
super(ModelState, self).__setitem__(k, v)
self.set_mutable(k)
def set_mutable(self, k):
self.ismutable.add(k)
def set_immutable(self, k):
if k in self.ismutable:
self.ismutable.remove(k)
class AeoLiS(IBmi):
'''AeoLiS model class
AeoLiS is a process-based model for simulating supply-limited
aeolian sediment transport. This model class is compatible with
the Basic Model Interface (BMI) and provides basic model
operations, like initialization, time stepping, finalization and
data exchange. For higher level operations, like a progress
indicator and netCDF4 output is refered to the AeoLiS model
runner class, see :class:`~model.AeoLiSRunner`.
Examples
--------
>>> with AeoLiS(configfile='aeolis.txt') as model:
>>> while model.get_current_time() <= model.get_end_time():
>>> model.update()
>>> model = AeoLiS(configfile='aeolis.txt')
>>> model.initialize()
>>> zb = model.get_var('zb')
>>> model.set_var('zb', zb + 1)
>>> for i in range(10):
>>> model.update(60.) # step 60 seconds forward
>>> model.finalize()
'''
def __init__(self, configfile):
'''Initialize class
Parameters
----------
configfile : str
Model configuration file. See :func:`~inout.read_configfile()`.
'''
self.t = 0.
self.dt = 0.
self.configfile = ''
self.l = {} # previous spatial grids
self.s = ModelState() # spatial grids
self.p = {} # parameters
self.c = {} # counters
self.configfile = configfile
def __enter__(self):
self.initialize()
return self
def __exit__(self, *args):
self.finalize()
def initialize(self):
'''Initialize model
Read model configuration file and initialize parameters and
spatial grids dictionary and load bathymetry and bed
composition.
'''
# read configuration file
self.p = aeolis.inout.read_configfile(self.configfile)
aeolis.inout.check_configuration(self.p)
# set nx, ny and nfractions
if self.p['xgrid_file'].ndim == 2:
self.p['ny'], self.p['nx'] = self.p['xgrid_file'].shape
# change from number of points to number of cells
self.p['nx'] -= 1
self.p['ny'] -= 1
else:
self.p['nx'] = len(self.p['xgrid_file'])
self.p['nx'] -= 1
self.p['ny'] = 0
#self.p['nfractions'] = len(self.p['grain_dist'])
self.p['nfractions'] = len(self.p['grain_size'])
# initialize time
self.t = self.p['tstart']
# get model dimensions
nx = self.p['nx']
ny = self.p['ny']
nl = self.p['nlayers']
nf = self.p['nfractions']
# initialize spatial grids
for var, dims in self.dimensions().items():
self.s[var] = np.zeros(self._dims2shape(dims))
self.l[var] = self.s[var].copy()
# initialize grid parameters
self.s, self.p = aeolis.gridparams.initialize(self.s, self.p)
# initialize bed composition
self.s = aeolis.bed.initialize(self.s, self.p)
# initialize wind model
self.s = aeolis.wind.initialize(self.s, self.p)
#initialize vegetation model
self.s = aeolis.vegetation.initialize(self.s, self.p)
#initialize fence model
self.s = aeolis.fences.initialize(self.s, self.p)
# Create interpretation information
if self.p['visualization']:
aeolis.inout.visualize_grid(self.s, self.p)
aeolis.inout.visualize_timeseries(self.p, self.t)
def update(self, dt=-1):
'''Time stepping function
Takes a single step in time. Interpolates wind and
hydrodynamic time series to the current time, updates the soil
moisture, mixes the bed due to wave action, computes wind
velocity threshold and the equilibrium sediment transport
concentration. Subsequently runs one of the available
numerical schemes to compute the instantaneous sediment
concentration and pickup for the next time step and updates
the bed accordingly.
For explicit schemes the time step is maximized by the
Courant-Friedrichs-Lewy (CFL) condition. See
:func:`~model.AeoLiS.set_timestep()`.
Parameters
----------
dt : float, optional
Time step in seconds. The time step specified in the model
configuration file is used in case dt is smaller than
zero. For explicit numerical schemes the time step is
maximized by the CFL confition.
'''
self.p['_time'] = self.t
# store previous state
self.l = self.s.copy()
self.l['zb'] = self.s['zb'].copy()
self.l['dzbavg'] = self.s['dzbavg'].copy()
# interpolate wind time series
self.s = aeolis.wind.interpolate(self.s, self.p, self.t)
# Rotate gridparams, such that the grids is alligned horizontally
self.s = self.grid_rotate(self.p['alpha'])
if np.sum(self.s['uw']) != 0:
self.s = aeolis.wind.shear(self.s, self.p)
#compute sand fence shear
if self.p['process_fences']:
self.s = aeolis.fences.update_fences(self.s, self.p)
# compute vegetation shear
if self.p['process_vegetation']:
self.s = aeolis.vegetation.vegshear(self.s, self.p)
# determine optimal time step
self.dt_prev = self.dt
if not self.set_timestep(dt):
return
# interpolate hydrodynamic time series
self.s = aeolis.hydro.interpolate(self.s, self.p, self.t)
self.s = aeolis.hydro.update(self.s, self.p, self.dt, self.t)
# mix top layer
self.s = aeolis.bed.mixtoplayer(self.s, self.p)
# compute threshold
if np.sum(self.s['uw']) != 0:
self.s = aeolis.threshold.compute(self.s, self.p)
# compute saltation velocity and equilibrium transport
self.s = aeolis.transport.equilibrium(self.s, self.p)
# compute instantaneous transport
if self.p['scheme'] == 'euler_forward':
self.s.update(self.euler_forward())
elif self.p['scheme'] == 'euler_backward':
self.s.update(self.euler_backward())
elif self.p['scheme'] == 'crank_nicolson':
self.s.update(self.crank_nicolson())
else:
logger.log_and_raise('Unknown scheme [%s]' % self.p['scheme'], exc=ValueError)
# update bed
self.s = aeolis.bed.update(self.s, self.p)
# avalanching
self.s = aeolis.avalanching.angele_of_repose(self.s, self.p)
self.s = aeolis.avalanching.avalanche(self.s, self.p)
# reset original bed in marine zone (wet)
self.s = aeolis.bed.wet_bed_reset(self.s, self.p)
# calculate average bed level change over time
self.s = aeolis.bed.average_change(self.l, self.s, self.p)
# compute dune erosion
if self.p['process_dune_erosion']:
self.s = aeolis.erosion.run_ph12(self.s, self.p, self.t)
self.s = aeolis.avalanching.angele_of_repose(self.s, self.p) #Since the aeolian module is only run for winds above threshold, also run avalanching routine here
self.s = aeolis.avalanching.avalanche(self.s, self.p)
# grow vegetation
if self.p['process_vegetation']:
self.s = aeolis.vegetation.germinate(self.s, self.p)
self.s = aeolis.vegetation.grow(self.s, self.p)
# increment time
self.t += self.dt * self.p['accfac']
self._count('time')
# Rotate gridparams back to original grid orientation
self.s = self.grid_rotate(-self.p['alpha'])
# Visualization of the model results after the first time step as a check for interpretation
if self.c['time'] == 1 and self.p['visualization']:
aeolis.inout.visualize_spatial(self.s, self.p)
def finalize(self):
'''Finalize model'''
pass
def get_current_time(self):
'''
Returns
-------
float
Current simulation time
'''
return self.t
def get_end_time(self):
'''
Returns
-------
float
Final simulation time
'''
return self.p['tstop']
def get_start_time(self):
'''
Returns
-------
float
Initial simulation time
'''
return self.p['tstart']
def get_var(self, var):
'''Returns spatial grid or model configuration parameter
If the given variable name matches with a spatial grid, the
spatial grid is returned. If not, the given variable name is
matched with a model configuration parameter. If a match is
found, the parameter value is returned. Otherwise, nothing is
returned.
Parameters
----------
var : str
Name of spatial grid or model configuration parameter
Returns
-------
np.ndarray or int, float, str or list
Spatial grid or model configuration parameter
Examples
--------
>>> # returns bathymetry grid
... model.get_var('zb')
>>> # returns simulation duration
... model.get_var('tstop')
See Also
--------
model.AeoLiS.set_var
'''
if var in self.s:
if var in ['Ct', 'Cu']:
return self.s[var] / self.p['accfac']
else:
return self.s[var]
elif var in self.p:
return self.p[var]
else:
return None
def get_var_count(self):
'''
Returns
-------
int
Number of spatial grids
'''
return len(self.s)
def get_var_name(self, i):
'''Returns name of spatial grid by index (in alphabetical order)
Parameters
----------
i : int
Index of spatial grid
Returns
-------
str or -1
Name of spatial grid or -1 in case index exceeds the number of grids
'''
if len(self.s) > i:
return sorted(self.s.keys())[i]
else:
return -1
def get_var_rank(self, var):
'''Returns rank of spatial grid
Parameters
----------
var : str
Name of spatial grid
Returns
-------
int
Rank of spatial grid or -1 if not found
'''
if var in self.s:
return len(self.s[var].shape)
else:
return -1
def get_var_shape(self, var):
'''Returns shape of spatial grid
Parameters
----------
var : str
Name of spatial grid
Returns
-------
tuple or int
Dimensions of spatial grid or -1 if not found
'''
if var in self.s:
return self.s[var].shape
else:
return -1
def get_var_type(self, var):
'''Returns variable type of spatial grid
Parameters
----------
var : str
Name of spatial grid
Returns
-------
str or int
Variable type of spatial grid or -1 if not found
'''
if var in self.s:
return 'double'
else:
return -1
def inq_compound(self):
logger.log_and_raise('Method not yet implemented [inq_compound]', exc=NotImplementedError)
def inq_compound_field(self):
logger.log_and_raise('Method not yet implemented [inq_compound_field]', exc=NotImplementedError)
def set_var(self, var, val):
'''Sets spatial grid or model configuration parameter
If the given variable name matches with a spatial grid, the
spatial grid is set. If not, the given variable name is
matched with a model configuration parameter. If a match is
found, the parameter value is set. Otherwise, nothing is set.
Parameters
----------
var : str
Name of spatial grid or model configuration parameter
val : np.ndarray or int, float, str or list
Spatial grid or model configuration parameter
Examples
--------
>>> # set bathymetry grid
... model.set_var('zb', np.array([[0.,0., ... ,0.]]))
>>> # set simulation duration
... model.set_var('tstop', 3600.)
See Also
--------
model.AeoLiS.get_var
'''
if var in self.s:
self.s[var] = val
elif var in self.p:
self.p[var] = val
def set_var_index(self, i, val):
'''Set spatial grid by index (in alphabetical order)
Parameters
----------
i : int
Index of spatial grid
val : np.ndarray
Spatial grid
'''
var = self.get_var_name(i)
self.set_var(var, val)
def set_var_slice(self):
logger.log_and_raise('Method not yet implemented [set_var_slice]', exc=NotImplementedError)
def set_timestep(self, dt=-1.):
'''Determine optimal time step
If no time step is given the optimal time step is
determined. For explicit numerical schemes the time step is
based in the Courant-Frierichs-Lewy (CFL) condition. For
implicit numerical schemes the time step specified in the
model configuration file is used. Alternatively, a preferred
time step is given that is maximized by the CFL condition in
case of an explicit numerical scheme.
Returns True except when:
1. No time step could be determined, for example when there is
no wind and the numerical scheme is explicit. In this case the
time step is set arbitrarily to one second.
2. Or when the time step is smaller than -1. In this case the
time is updated with the absolute value of the time step, but
no model execution is performed. This funcionality can be used
to skip fast-forward in time.
Parameters
----------
df : float, optional
Preferred time step
Returns
-------
bool
False if determination of time step was unsuccessful, True otherwise
'''
if dt > 0.:
self.dt = dt
elif dt < -1:
self.dt = dt
self.t += np.abs(dt)
return False
else:
self.dt = self.p['dt']
if self.p['scheme'] == 'euler_forward':
if self.p['CFL'] > 0.:
dtref = np.max(np.abs(self.s['uws']) / self.s['ds']) + \
np.max(np.abs(self.s['uwn']) / self.s['dn'])
if dtref > 0.:
self.dt = np.minimum(self.dt, self.p['CFL'] / dtref)
else:
self.dt = np.minimum(self.dt, 1.)
return False
if self.p['max_bedlevel_change'] != 999. and np.max(self.s['dzb']) != 0. and self.dt_prev != 0.:
dt_zb = self.dt_prev * self.p['max_bedlevel_change'] / np.max(self.s['dzb'])
self.dt = np.minimum(self.dt, dt_zb)
self.p['dt_opt'] = self.dt
return True
def grid_rotate(self, angle):
s = self.s
p = self.p
s['x'], s['y'] = rotate(s['x'], s['y'], angle, origin=(np.mean(s['x']), np.mean(s['y'])))
s['taus'], s['taun'] = rotate(s['taus'], s['taun'], angle, origin=(0, 0))
s['taus0'], s['taun0'] = rotate(s['taus0'], s['taun0'], angle, origin=(0, 0))
s['ustars'], s['ustarn'] = rotate(s['ustars'], s['ustarn'], angle, origin=(0, 0))
s['ustars0'], s['ustarn0'] = rotate(s['ustars0'], s['ustarn0'], angle, origin=(0, 0))
s['uws'], s['uwn'] = rotate(s['uws'], s['uwn'], angle, origin=(0, 0))
for i in range(p['nfractions']):
s['qs'][:,:,i], s['qn'][:,:,i] = rotate(s['qs'][:,:,i], s['qn'][:,:,i], angle, origin=(0, 0))
self.s['udir'] += angle
return s
def euler_forward(self):
'''Convenience function for explicit solver based on Euler forward scheme
See Also
--------
model.AeoLiS.solve
'''
if self.p['solver'].lower() == 'trunk':
solve = self.solve(alpha=0., beta=1)
elif self.p['solver'].lower() == 'pieter':
solve = self.solve_pieter(alpha=0., beta=1)
elif self.p['solver'].lower() == 'steadystate':
solve = self.solve_steadystate()
elif self.p['solver'].lower() == 'steadystatepieter':
solve = self.solve_steadystatepieter()
return solve
def euler_backward(self):
'''Convenience function for implicit solver based on Euler backward scheme
See Also
--------
model.AeoLiS.solve
'''
if self.p['solver'].lower() == 'trunk':
solve = self.solve(alpha=1., beta=1)
elif self.p['solver'].lower() == 'pieter':
solve = self.solve_pieter(alpha=1., beta=1)
elif self.p['solver'].lower() == 'steadystate':
solve = self.solve_steadystate()
elif self.p['solver'].lower() == 'steadystatepieter':
solve = self.solve_steadystatepieter()
return solve
def crank_nicolson(self):
'''Convenience function for semi-implicit solver based on Crank-Nicolson scheme
See Also
--------
model.AeoLiS.solve
'''
if self.p['solver'].lower() == 'trunk':
solve = self.solve(alpha=.5, beta=1)
elif self.p['solver'].lower() == 'pieter':
solve = self.solve_pieter(alpha=.5, beta=1)
elif self.p['solver'].lower() == 'steadystate':
solve = self.solve_steadystate()
elif self.p['solver'].lower() == 'steadystatepieter':
solve = self.solve_steadystatepieter()
return solve
def solve_steadystate(self):
'''Implements the steady state solution
'''
# upwind scheme:
beta = 1.
l = self.l
s = self.s
p = self.p
Ct = s['Ct'].copy()
pickup = s['pickup'].copy()
# compute transport weights for all sediment fractions
w_init, w_air, w_bed = aeolis.transport.compute_weights(s, p)
if self.t == 0.:
# use initial guess for first time step
if p['grain_dist'] != None:
w = p['grain_dist'].reshape((1,1,-1))
w = w.repeat(p['ny']+1, axis=0)
w = w.repeat(p['nx']+1, axis=1)
else:
w = w_init.copy()
else:
w = w_init.copy()
# set model state properties that are added to warnings and errors
logprops = dict(minwind=s['uw'].min(),
maxdrop=(l['uw']-s['uw']).max(),
time=self.t,
dt=self.dt)
nf = p['nfractions']
us = np.zeros((p['ny']+1,p['nx']+1))
un = np.zeros((p['ny']+1,p['nx']+1))
us_plus = np.zeros((p['ny']+1,p['nx']+1))
un_plus = np.zeros((p['ny']+1,p['nx']+1))
us_min = np.zeros((p['ny']+1,p['nx']+1))
un_min = np.zeros((p['ny']+1,p['nx']+1))
Cs = np.zeros(us.shape)
Cn = np.zeros(un.shape)
Cs_plus = np.zeros(us.shape)
Cn_plus = np.zeros(un.shape)
Cs_min = np.zeros(us.shape)
Cn_min = np.zeros(un.shape)
for i in range(nf):
us[:,:] = s['us'][:,:,i]
un[:,:] = s['un'][:,:,i]
us_plus[:,1:] = s['us'][:,:-1,i]
un_plus[1:,:] = s['un'][:-1,:,i]
us_min[:,:-1] = s['us'][:,1:,i]
un_min[:-1,:] = s['un'][1:,:,i]
#boundary values
us[:,0] = s['us'][:,0,i]
un[0,:] = s['un'][0,:,i]
us_plus[:,0] = s['us'][:,0,i]
un_plus[0,:] = s['un'][0,:,i]
us_min[:,-1] = s['us'][:,-1,i]
un_min[-1,:] = s['un'][-1,:,i]
# define matrix coefficients to solve linear system of equations
Cs = s['dn'] * s['dsdni'] * us[:,:]
Cn = s['ds'] * s['dsdni'] * un[:,:]
Cs_plus = s['dn'] * s['dsdni'] * us_plus[:,:]
Cn_plus = s['ds'] * s['dsdni'] * un_plus[:,:]
Cs_min = s['dn'] * s['dsdni'] * us_min[:,:]
Cn_min = s['ds'] * s['dsdni'] * un_min[:,:]
Ti = 1 / p['T']
beta = abs(beta)
if beta >= 1.:
# define upwind direction
ixs = np.asarray(us[:,:] >= 0., dtype=float)
ixn = np.asarray(un[:,:] >= 0., dtype=float)
sgs = 2. * ixs - 1.
sgn = 2. * ixn - 1.
else:
# or centralizing weights
ixs = beta + np.zeros(us)
ixn = beta + np.zeros(un)
sgs = np.zeros(us)
sgn = np.zeros(un)
# initialize matrix diagonals
A0 = np.zeros(s['zb'].shape)
Apx = np.zeros(s['zb'].shape)
Ap1 = np.zeros(s['zb'].shape)
Ap2 = np.zeros(s['zb'].shape)
Amx = np.zeros(s['zb'].shape)
Am1 = np.zeros(s['zb'].shape)
Am2 = np.zeros(s['zb'].shape)
# populate matrix diagonals
A0 = sgs * Cs + sgn * Cn + Ti
Apx = Cn_min * (1. - ixn)
Ap1 = Cs_min * (1. - ixs)
Amx = -Cn_plus * ixn
Am1 = -Cs_plus * ixs
# add boundaries
A0[:,0] = 1.
Apx[:,0] = 0.
Amx[:,0] = 0.
Am2[:,0] = 0.
Am1[:,0] = 0.
A0[:,-1] = 1.
Apx[:,-1] = 0.
Ap1[:,-1] = 0.
Ap2[:,-1] = 0.
Amx[:,-1] = 0.
if p['boundary_offshore'] == 'flux':
Ap2[:,0] = 0.
Ap1[:,0] = 0.
elif p['boundary_offshore'] == 'constant':
Ap2[:,0] = 0.
Ap1[:,0] = 0.
elif p['boundary_offshore'] == 'uniform':
Ap2[:,0] = 0.
Ap1[:,0] = -1.
elif p['boundary_offshore'] == 'gradient':
Ap2[:,0] = s['ds'][:,1] / s['ds'][:,2]
Ap1[:,0] = -1. - s['ds'][:,1] / s['ds'][:,2]
elif p['boundary_offshore'] == 'circular':
logger.log_and_raise('Cross-shore cricular boundary condition not yet implemented', exc=NotImplementedError)
else:
logger.log_and_raise('Unknown offshore boundary condition [%s]' % self.p['boundary_offshore'], exc=ValueError)
if p['boundary_onshore'] == 'flux':
Am2[:,-1] = 0.
Am1[:,-1] = 0.
elif p['boundary_onshore'] == 'constant':
Am2[:,-1] = 0.
Am1[:,-1] = 0.
elif p['boundary_onshore'] == 'uniform':
Am2[:,-1] = 0.
Am1[:,-1] = -1.
elif p['boundary_onshore'] == 'gradient':
Am2[:,-1] = s['ds'][:,-2] / s['ds'][:,-3]
Am1[:,-1] = -1. - s['ds'][:,-2] / s['ds'][:,-3]
elif p['boundary_offshore'] == 'circular':
logger.log_and_raise('Cross-shore cricular boundary condition not yet implemented', exc=NotImplementedError)
else:
logger.log_and_raise('Unknown onshore boundary condition [%s]' % self.p['boundary_onshore'], exc=ValueError)
if p['boundary_lateral'] == 'constant':
A0[0,:] = 1.
Apx[0,:] = 0.
Ap1[0,:] = 0.
Amx[0,:] = 0.
Am1[0,:] = 0.
A0[-1,:] = 1.
Apx[-1,:] = 0.
Ap1[-1,:] = 0.
Amx[-1,:] = 0.
Am1[-1,:] = 0.
#logger.log_and_raise('Lateral constant boundary condition not yet implemented', exc=NotImplementedError)
elif p['boundary_lateral'] == 'uniform':
logger.log_and_raise('Lateral uniform boundary condition not yet implemented', exc=NotImplementedError)
elif p['boundary_lateral'] == 'gradient':
logger.log_and_raise('Lateral gradient boundary condition not yet implemented', exc=NotImplementedError)
elif p['boundary_lateral'] == 'circular':
pass
else:
logger.log_and_raise('Unknown lateral boundary condition [%s]' % self.p['boundary_lateral'], exc=ValueError)
# construct sparse matrix
if p['ny'] > 0:
j = p['nx']+1
A = scipy.sparse.diags((Apx.flatten()[:j],
Amx.flatten()[j:],
Am2.flatten()[2:],
Am1.flatten()[1:],
A0.flatten(),
Ap1.flatten()[:-1],
Ap2.flatten()[:-2],
Apx.flatten()[j:],
Amx.flatten()[:j]),
(-j*p['ny'],-j,-2,-1,0,1,2,j,j*p['ny']), format='csr')
else:
A = scipy.sparse.diags((Am2.flatten()[2:],
Am1.flatten()[1:],
A0.flatten(),
Ap1.flatten()[:-1],
Ap2.flatten()[:-2]),
(-2,-1,0,1,2), format='csr')
# solve transport for each fraction separately using latest
# available weights
# renormalize weights for all fractions equal or larger
# than the current one such that the sum of all weights is
# unity
w = aeolis.transport.renormalize_weights(w, i)
# iteratively find a solution of the linear system that
# does not violate the availability of sediment in the bed
for n in range(p['max_iter']):
self._count('matrixsolve')
# compute saturation levels
ix = s['Cu'] > 0.
S_i = np.zeros(s['Cu'].shape)
S_i[ix] = s['Ct'][ix] / s['Cu'][ix]
s['S'] = S_i.sum(axis=-1)
# create the right hand side of the linear system
y_i = np.zeros(s['zb'].shape)
y_i[:,1:-1] = (
(w[:,1:-1,i] * s['Cuf'][:,1:-1,i] * Ti) * (1. - s['S'][:,1:-1]) +
(w[:,1:-1,i] * s['Cu'][:,1:-1,i] * Ti) * s['S'][:,1:-1]
)
# add boundaries
if p['boundary_offshore'] == 'flux':
y_i[:,0] = p['offshore_flux'] * s['Cu0'][:,0,i]
if p['boundary_onshore'] == 'flux':
y_i[:,-1] = p['onshore_flux'] * s['Cu0'][:,-1,i]
if p['boundary_offshore'] == 'constant':
y_i[:,0] = p['constant_offshore_flux'] / s['u'][:,0,i]
if p['boundary_onshore'] == 'constant':
y_i[:,-1] = p['constant_onshore_flux'] / s['u'][:,-1,i]
# solve system with current weights
Ct_i = scipy.sparse.linalg.spsolve(A, y_i.flatten())
Ct_i = prevent_tiny_negatives(Ct_i, p['max_error'])
# check for negative values
if Ct_i.min() < 0.:
ix = Ct_i < 0.