• Home
  • Features
  • Pricing
  • Docs
  • Announcements
  • Sign In

OGGM / oggm / 4861042853

pending completion
4861042853

push

github

GitHub
add to whats-new (#1569)

11181 of 13110 relevant lines covered (85.29%)

3.41 hits per line

Source File
Press 'n' to go to next uncovered line, 'b' for previous

89.86
/oggm/core/massbalance.py
1
"""Mass balance models - next generation"""
2

3
# Built ins
4
import logging
9✔
5
import os
9✔
6
# External libs
7
import cftime
9✔
8
import numpy as np
9✔
9
import xarray as xr
9✔
10
import pandas as pd
9✔
11
from scipy.interpolate import interp1d
9✔
12
from scipy import optimize
9✔
13
# Locals
14
import oggm.cfg as cfg
9✔
15
from oggm.cfg import SEC_IN_YEAR, SEC_IN_MONTH
9✔
16
from oggm.utils import (SuperclassMeta, get_geodetic_mb_dataframe,
9✔
17
                        floatyear_to_date, date_to_floatyear,
18
                        monthly_timeseries, ncDataset, get_temp_bias_dataframe,
19
                        clip_min, clip_max, clip_array, clip_scalar,
20
                        weighted_average_1d, lazy_property)
21
from oggm.exceptions import (InvalidWorkflowError, InvalidParamsError,
9✔
22
                             MassBalanceCalibrationError)
23
from oggm import entity_task
9✔
24

25
# Module logger
26
log = logging.getLogger(__name__)
9✔
27

28
# Climate relevant global params - not optimised
29
MB_GLOBAL_PARAMS = ['temp_default_gradient',
9✔
30
                    'temp_all_solid',
31
                    'temp_all_liq',
32
                    'temp_melt']
33

34

35
class MassBalanceModel(object, metaclass=SuperclassMeta):
9✔
36
    """Interface and common logic for all mass balance models used in OGGM.
37

38
    All mass balance models should implement this interface.
39

40
    Attributes
41
    ----------
42
    valid_bounds : [float, float]
43
        The altitudinal bounds where the MassBalanceModel is valid. This is
44
        necessary for automated ELA search.
45
    hemisphere : str, {'nh', 'sh'}
46
        Used for certain methods - if the hydrological year is requested.
47
    rho : float, default: ``cfg.PARAMS['ice_density']``
48
        Density of ice
49
    """
50

51
    def __init__(self):
9✔
52
        """ Initialize."""
53
        self.valid_bounds = None
9✔
54
        self.hemisphere = None
9✔
55
        self.rho = cfg.PARAMS['ice_density']
9✔
56

57
    def __repr__(self):
58
        """String Representation of the mass balance model"""
59
        summary = ['<oggm.MassBalanceModel>']
60
        summary += ['  Class: ' + self.__class__.__name__]
61
        summary += ['  Attributes:']
62
        # Add all scalar attributes
63
        for k, v in self.__dict__.items():
64
            if np.isscalar(v) and not k.startswith('_'):
65
                nbform = '    - {}: {}'
66
                summary += [nbform.format(k, v)]
67
        return '\n'.join(summary) + '\n'
68

69
    def get_monthly_mb(self, heights, year=None, fl_id=None, fls=None):
9✔
70
        """Monthly mass balance at given altitude(s) for a moment in time.
71

72
        Units: [m s-1], or meters of ice per second
73

74
        Note: `year` is optional because some simpler models have no time
75
        component.
76

77
        Parameters
78
        ----------
79
        heights: ndarray
80
            the atitudes at which the mass balance will be computed
81
        year: float, optional
82
            the time (in the "floating year" convention)
83
        fl_id: float, optional
84
            the index of the flowline in the fls array (might be ignored
85
            by some MB models)
86
        fls: list of flowline instances, optional
87
            the flowlines array, in case the MB model implementation needs
88
            to know details about the glacier geometry at the moment the
89
            MB model is called
90

91
        Returns
92
        -------
93
        the mass balance (same dim as `heights`) (units: [m s-1])
94
        """
95
        raise NotImplementedError()
96

97
    def get_annual_mb(self, heights, year=None, fl_id=None, fls=None):
9✔
98
        """Like `self.get_monthly_mb()`, but for annual MB.
99

100
        For some simpler mass balance models ``get_monthly_mb()` and
101
        `get_annual_mb()`` can be equivalent.
102

103
        Units: [m s-1], or meters of ice per second
104

105
        Note: `year` is optional because some simpler models have no time
106
        component.
107

108
        Parameters
109
        ----------
110
        heights: ndarray
111
            the altitudes at which the mass balance will be computed
112
        year: float, optional
113
            the time (in the "floating year" convention)
114
        fl_id: float, optional
115
            the index of the flowline in the fls array (might be ignored
116
            by some MB models)
117
        fls: list of flowline instances, optional
118
            the flowlines array, in case the MB model implementation needs
119
            to know details about the glacier geometry at the moment the
120
            MB model is called
121

122
        Returns
123
        -------
124
        the mass balance (same dim as `heights`) (units: [m s-1])
125
        """
126
        raise NotImplementedError()
127

128
    def get_specific_mb(self, heights=None, widths=None, fls=None,
9✔
129
                        year=None):
130
        """Specific mb for this year and a specific glacier geometry.
131

132
         Units: [mm w.e. yr-1], or millimeter water equivalent per year
133

134
        Parameters
135
        ----------
136
        heights: ndarray
137
            the altitudes at which the mass balance will be computed.
138
            Overridden by ``fls`` if provided
139
        widths: ndarray
140
            the widths of the flowline (necessary for the weighted average).
141
            Overridden by ``fls`` if provided
142
        fls: list of flowline instances, optional
143
            Another way to get heights and widths - overrides them if
144
            provided.
145
        year: float, optional
146
            the time (in the "floating year" convention)
147

148
        Returns
149
        -------
150
        the specific mass balance (units: mm w.e. yr-1)
151
        """
152

153
        if len(np.atleast_1d(year)) > 1:
6✔
154
            out = [self.get_specific_mb(heights=heights, widths=widths,
5✔
155
                                        fls=fls, year=yr)
156
                   for yr in year]
157
            return np.asarray(out)
5✔
158

159
        if fls is not None:
6✔
160
            mbs = []
6✔
161
            widths = []
6✔
162
            for i, fl in enumerate(fls):
6✔
163
                _widths = fl.widths
6✔
164
                try:
6✔
165
                    # For rect and parabola don't compute spec mb
166
                    _widths = np.where(fl.thick > 0, _widths, 0)
6✔
167
                except AttributeError:
5✔
168
                    pass
5✔
169
                widths = np.append(widths, _widths)
6✔
170
                mbs = np.append(mbs, self.get_annual_mb(fl.surface_h,
6✔
171
                                                        fls=fls, fl_id=i,
172
                                                        year=year))
173
        else:
174
            mbs = self.get_annual_mb(heights, year=year)
3✔
175

176
        return weighted_average_1d(mbs, widths) * SEC_IN_YEAR * self.rho
6✔
177

178
    def get_ela(self, year=None, **kwargs):
9✔
179
        """Compute the equilibrium line altitude for a given year.
180

181
        Parameters
182
        ----------
183
        year: float, optional
184
            the time (in the "floating year" convention)
185
        **kwargs: any other keyword argument accepted by self.get_annual_mb
186
        Returns
187
        -------
188
        the equilibrium line altitude (ELA, units: m)
189
        """
190

191
        if len(np.atleast_1d(year)) > 1:
2✔
192
            return np.asarray([self.get_ela(year=yr, **kwargs) for yr in year])
1✔
193

194
        if self.valid_bounds is None:
2✔
195
            raise ValueError('attribute `valid_bounds` needs to be '
×
196
                             'set for the ELA computation.')
197

198
        # Check for invalid ELAs
199
        b0, b1 = self.valid_bounds
2✔
200
        if (np.any(~np.isfinite(
2✔
201
                self.get_annual_mb([b0, b1], year=year, **kwargs))) or
202
                (self.get_annual_mb([b0], year=year, **kwargs)[0] > 0) or
203
                (self.get_annual_mb([b1], year=year, **kwargs)[0] < 0)):
204
            return np.NaN
1✔
205

206
        def to_minimize(x):
2✔
207
            return (self.get_annual_mb([x], year=year, **kwargs)[0] *
2✔
208
                    SEC_IN_YEAR * self.rho)
209
        return optimize.brentq(to_minimize, *self.valid_bounds, xtol=0.1)
2✔
210

211
    def is_year_valid(self, year):
9✔
212
        """Checks if a given date year be simulated by this model.
213

214
        Parameters
215
        ----------
216
        year: float, optional
217
            the time (in the "floating year" convention)
218

219
        Returns
220
        -------
221
        True if this year can be simulated by the model
222
        """
223
        raise NotImplementedError()
224

225

226
class ScalarMassBalance(MassBalanceModel):
9✔
227
    """Constant mass balance, everywhere."""
228

229
    def __init__(self, mb=0.):
9✔
230
        """ Initialize.
231

232
        Parameters
233
        ----------
234
        mb : float
235
            Fix the mass balance to a certain value (unit: [mm w.e. yr-1])
236
        """
237
        super(ScalarMassBalance, self).__init__()
2✔
238
        self.hemisphere = 'nh'
2✔
239
        self.valid_bounds = [-2e4, 2e4]  # in m
2✔
240
        self._mb = mb
2✔
241

242
    def get_monthly_mb(self, heights, **kwargs):
9✔
243
        mb = np.asarray(heights) * 0 + self._mb
×
244
        return mb / SEC_IN_YEAR / self.rho
×
245

246
    def get_annual_mb(self, heights, **kwargs):
9✔
247
        mb = np.asarray(heights) * 0 + self._mb
2✔
248
        return mb / SEC_IN_YEAR / self.rho
2✔
249

250
    def is_year_valid(self, year):
9✔
251
        return True
×
252

253

254
class LinearMassBalance(MassBalanceModel):
9✔
255
    """Constant mass balance as a linear function of altitude.
256

257
    Attributes
258
    ----------
259
    ela_h: float
260
        the equilibrium line altitude (units: [m])
261
    grad: float
262
        the mass balance gradient (unit: [mm w.e. yr-1 m-1])
263
    max_mb: float
264
        Cap the mass balance to a certain value (unit: [mm w.e. yr-1])
265
    temp_bias
266
    """
267

268
    def __init__(self, ela_h, grad=3., max_mb=None):
9✔
269
        """ Initialize.
270

271
        Parameters
272
        ----------
273
        ela_h: float
274
            Equilibrium line altitude (units: [m])
275
        grad: float
276
            Mass balance gradient (unit: [mm w.e. yr-1 m-1])
277
        max_mb: float
278
            Cap the mass balance to a certain value (unit: [mm w.e. yr-1])
279
        """
280
        super(LinearMassBalance, self).__init__()
8✔
281
        self.hemisphere = 'nh'
8✔
282
        self.valid_bounds = [-1e4, 2e4]  # in m
8✔
283
        self.orig_ela_h = ela_h
8✔
284
        self.ela_h = ela_h
8✔
285
        self.grad = grad
8✔
286
        self.max_mb = max_mb
8✔
287
        self._temp_bias = 0
8✔
288

289
    @property
9✔
290
    def temp_bias(self):
9✔
291
        """Change the ELA following a simple rule: + 1K -> ELA + 150 m
292

293
        A "temperature bias" doesn't makes much sense in the linear MB
294
        context, but we implemented a simple empirical rule:
295
        + 1K -> ELA + 150 m
296
        """
297
        return self._temp_bias
×
298

299
    @temp_bias.setter
9✔
300
    def temp_bias(self, value):
9✔
301
        self.ela_h = self.orig_ela_h + value * 150
1✔
302
        self._temp_bias = value
1✔
303

304
    def get_monthly_mb(self, heights, **kwargs):
9✔
305
        mb = (np.asarray(heights) - self.ela_h) * self.grad
6✔
306
        if self.max_mb is not None:
6✔
307
            clip_max(mb, self.max_mb, out=mb)
×
308
        return mb / SEC_IN_YEAR / self.rho
6✔
309

310
    def get_annual_mb(self, heights, **kwargs):
9✔
311
        return self.get_monthly_mb(heights, **kwargs)
6✔
312

313
    def is_year_valid(self, year):
9✔
314
        return True
×
315

316

317
class MonthlyTIModel(MassBalanceModel):
9✔
318
    """Monthly temperature index model.
319
    """
320
    def __init__(self, gdir,
9✔
321
                 filename='climate_historical',
322
                 input_filesuffix='',
323
                 fl_id=None,
324
                 melt_f=None,
325
                 temp_bias=None,
326
                 prcp_fac=None,
327
                 bias=0,
328
                 ys=None,
329
                 ye=None,
330
                 repeat=False,
331
                 check_calib_params=True,
332
                 ):
333
        """Initialize.
334

335
        Parameters
336
        ----------
337
        gdir : GlacierDirectory
338
            the glacier directory
339
        filename : str, optional
340
            set to a different BASENAME if you want to use alternative climate
341
            data. Default is 'climate_historical'
342
        input_filesuffix : str, optional
343
            append a suffix to the filename (useful for GCM runs).
344
        fl_id : int, optional
345
            if this flowline has been calibrated alone and has specific
346
            model parameters.
347
        melt_f : float, optional
348
            set to the value of the melt factor you want to use,
349
            here the unit is kg m-2 day-1 K-1
350
            (the default is to use the calibrated value).
351
        temp_bias : float, optional
352
            set to the value of the temperature bias you want to use
353
            (the default is to use the calibrated value).
354
        prcp_fac : float, optional
355
            set to the value of the precipitation factor you want to use
356
            (the default is to use the calibrated value).
357
        bias : float, optional
358
            set to the alternative value of the calibration bias [mm we yr-1]
359
            you want to use (the default is to use the calibrated value)
360
            Note that this bias is *substracted* from the computed MB. Indeed:
361
            BIAS = MODEL_MB - REFERENCE_MB.
362
        ys : int
363
            The start of the climate period where the MB model is valid
364
            (default: the period with available data)
365
        ye : int
366
            The end of the climate period where the MB model is valid
367
            (default: the period with available data)
368
        repeat : bool
369
            Whether the climate period given by [ys, ye] should be repeated
370
            indefinitely in a circular way
371
        check_calib_params : bool
372
            OGGM will try hard not to use wrongly calibrated parameters
373
            by checking the global parameters used during calibration and
374
            the ones you are using at run time. If they don't match, it will
375
            raise an error. Set to "False" to suppress this check.
376
        """
377

378
        super(MonthlyTIModel, self).__init__()
5✔
379
        self.valid_bounds = [-1e4, 2e4]  # in m
5✔
380
        self.fl_id = fl_id  # which flowline are we the model of?
5✔
381
        self.gdir = gdir
5✔
382

383
        if melt_f is None:
5✔
384
            melt_f = self.calib_params['melt_f']
5✔
385

386
        if temp_bias is None:
5✔
387
            temp_bias = self.calib_params['temp_bias']
5✔
388

389
        if prcp_fac is None:
5✔
390
            prcp_fac = self.calib_params['prcp_fac']
5✔
391

392
        # Check the climate related params to the GlacierDir to make sure
393
        if check_calib_params:
5✔
394
            mb_calib = self.calib_params['mb_global_params']
5✔
395
            for k, v in mb_calib.items():
5✔
396
                if v != cfg.PARAMS[k]:
5✔
397
                    msg = ('You seem to use different mass balance parameters '
×
398
                           'than used for the calibration: '
399
                           f"you use cfg.PARAMS['{k}']={cfg.PARAMS[k]} while "
400
                           f"it was calibrated with cfg.PARAMS['{k}']={v}. "
401
                           'Set `check_calib_params=False` to ignore this '
402
                           'warning.')
403
                    raise InvalidWorkflowError(msg)
×
404
            src = self.calib_params['baseline_climate_source']
5✔
405
            src_calib = gdir.get_climate_info()['baseline_climate_source']
5✔
406
            if src != src_calib:
5✔
407
                msg = (f'You seem to have calibrated with the {src} '
×
408
                       f"climate data while this gdir was calibrated with "
409
                       f"{src_calib}. Set `check_calib_params=False` to "
410
                       f"ignore this warning.")
411
                raise InvalidWorkflowError(msg)
×
412

413
        self.melt_f = melt_f
5✔
414
        self.bias = bias
5✔
415

416
        # Global parameters
417
        self.t_solid = cfg.PARAMS['temp_all_solid']
5✔
418
        self.t_liq = cfg.PARAMS['temp_all_liq']
5✔
419
        self.t_melt = cfg.PARAMS['temp_melt']
5✔
420

421
        # check if valid prcp_fac is used
422
        if prcp_fac <= 0:
5✔
423
            raise InvalidParamsError('prcp_fac has to be above zero!')
×
424
        default_grad = cfg.PARAMS['temp_default_gradient']
5✔
425

426
        # Public attrs
427
        self.hemisphere = gdir.hemisphere
5✔
428
        self.repeat = repeat
5✔
429

430
        # Private attrs
431
        # to allow prcp_fac to be changed after instantiation
432
        # prescribe the prcp_fac as it is instantiated
433
        self._prcp_fac = prcp_fac
5✔
434
        # same for temp bias
435
        self._temp_bias = temp_bias
5✔
436

437
        # Read climate file
438
        fpath = gdir.get_filepath(filename, filesuffix=input_filesuffix)
5✔
439
        with ncDataset(fpath, mode='r') as nc:
5✔
440
            # time
441
            time = nc.variables['time']
5✔
442
            time = cftime.num2date(time[:], time.units, calendar=time.calendar)
5✔
443
            ny, r = divmod(len(time), 12)
5✔
444
            if r != 0:
5✔
445
                raise ValueError('Climate data should be N full years')
×
446

447
            # We check for calendar years
448
            if (time[0].month != 1) or (time[-1].month != 12):
5✔
449
                raise InvalidWorkflowError('We now work exclusively with '
×
450
                                           'calendar years.')
451

452
            # Quick trick because we know the size of our array
453
            years = np.repeat(np.arange(time[-1].year - ny + 1,
5✔
454
                                        time[-1].year + 1), 12)
455
            pok = slice(None)  # take all is default (optim)
5✔
456
            if ys is not None:
5✔
457
                pok = years >= ys
1✔
458
            if ye is not None:
5✔
459
                try:
1✔
460
                    pok = pok & (years <= ye)
1✔
461
                except TypeError:
×
462
                    pok = years <= ye
×
463

464
            self.years = years[pok]
5✔
465
            self.months = np.tile(np.arange(1, 13), ny)[pok]
5✔
466

467
            # Read timeseries and correct it
468
            self.temp = nc.variables['temp'][pok].astype(np.float64) + self._temp_bias
5✔
469
            self.prcp = nc.variables['prcp'][pok].astype(np.float64) * self._prcp_fac
5✔
470

471
            grad = self.prcp * 0 + default_grad
5✔
472
            self.grad = grad
5✔
473
            self.ref_hgt = nc.ref_hgt
5✔
474
            self.climate_source = nc.climate_source
5✔
475
            self.ys = self.years[0]
5✔
476
            self.ye = self.years[-1]
5✔
477

478
    def __repr__(self):
479
        """String Representation of the mass balance model"""
480
        summary = ['<oggm.MassBalanceModel>']
481
        summary += ['  Class: ' + self.__class__.__name__]
482
        summary += ['  Attributes:']
483
        # Add all scalar attributes
484
        done = []
485
        for k in ['hemisphere', 'climate_source', 'melt_f', 'prcp_fac', 'temp_bias', 'bias']:
486
            done.append(k)
487
            v = self.__getattribute__(k)
488
            if k == 'climate_source':
489
                if v.endswith('.nc'):
490
                    v = os.path.basename(v)
491
            nofloat = ['hemisphere', 'climate_source']
492
            nbform = '    - {}: {}' if k in nofloat else '    - {}: {:.02f}'
493
            summary += [nbform.format(k, v)]
494
        for k, v in self.__dict__.items():
495
            if np.isscalar(v) and not k.startswith('_') and k not in done:
496
                nbform = '    - {}: {}'
497
                summary += [nbform.format(k, v)]
498
        return '\n'.join(summary) + '\n'
499

500
    @property
9✔
501
    def monthly_melt_f(self):
9✔
502
        return self.melt_f * 365 / 12
5✔
503

504
    # adds the possibility of changing prcp_fac
505
    # after instantiation with properly changing the prcp time series
506
    @property
9✔
507
    def prcp_fac(self):
9✔
508
        """Precipitation factor (default: cfg.PARAMS['prcp_fac'])
509

510
        Called factor to make clear that it is a multiplicative factor in
511
        contrast to the additive temperature bias
512
        """
513
        return self._prcp_fac
1✔
514

515
    @prcp_fac.setter
9✔
516
    def prcp_fac(self, new_prcp_fac):
9✔
517
        # just to check that no invalid prcp_factors are used
518
        if np.any(np.asarray(new_prcp_fac) <= 0):
3✔
519
            raise InvalidParamsError('prcp_fac has to be above zero!')
1✔
520

521
        if len(np.atleast_1d(new_prcp_fac)) == 12:
3✔
522
            # OK so that's monthly stuff
523
            new_prcp_fac = np.tile(new_prcp_fac, len(self.prcp) // 12)
1✔
524

525
        self.prcp *= new_prcp_fac / self._prcp_fac
3✔
526
        self._prcp_fac = new_prcp_fac
3✔
527

528
    @property
9✔
529
    def temp_bias(self):
9✔
530
        """Add a temperature bias to the time series"""
531
        return self._temp_bias
3✔
532

533
    @temp_bias.setter
9✔
534
    def temp_bias(self, new_temp_bias):
9✔
535

536
        if len(np.atleast_1d(new_temp_bias)) == 12:
4✔
537
            # OK so that's monthly stuff
538
            new_temp_bias = np.tile(new_temp_bias, len(self.temp) // 12)
1✔
539

540
        self.temp += new_temp_bias - self._temp_bias
4✔
541
        self._temp_bias = new_temp_bias
4✔
542

543
    @lazy_property
9✔
544
    def calib_params(self):
9✔
545
        if self.fl_id is None:
5✔
546
            return self.gdir.read_json('mb_calib')
5✔
547
        else:
548
            try:
×
549
                return self.gdir.read_json('mb_calib',
×
550
                                           filesuffix=f'_fl{self.fl_id}')
551
            except FileNotFoundError:
×
552
                return self.gdir.read_json('mb_calib')
×
553

554
    def is_year_valid(self, year):
9✔
555
        return self.ys <= year <= self.ye
5✔
556

557
    def get_monthly_climate(self, heights, year=None):
9✔
558
        """Monthly climate information at given heights.
559

560
        Note that prcp is corrected with the precipitation factor and that
561
        all other model biases (temp and prcp) are applied.
562

563
        Returns
564
        -------
565
        (temp, tempformelt, prcp, prcpsol)
566
        """
567

568
        y, m = floatyear_to_date(year)
2✔
569
        if self.repeat:
2✔
570
            y = self.ys + (y - self.ys) % (self.ye - self.ys + 1)
×
571
        if not self.is_year_valid(y):
2✔
572
            raise ValueError('year {} out of the valid time bounds: '
×
573
                             '[{}, {}]'.format(y, self.ys, self.ye))
574
        pok = np.where((self.years == y) & (self.months == m))[0][0]
2✔
575

576
        # Read already (temperature bias and precipitation factor corrected!)
577
        itemp = self.temp[pok]
2✔
578
        iprcp = self.prcp[pok]
2✔
579
        igrad = self.grad[pok]
2✔
580

581
        # For each height pixel:
582
        # Compute temp and tempformelt (temperature above melting threshold)
583
        npix = len(heights)
2✔
584
        temp = np.ones(npix) * itemp + igrad * (heights - self.ref_hgt)
2✔
585
        tempformelt = temp - self.t_melt
2✔
586
        clip_min(tempformelt, 0, out=tempformelt)
2✔
587

588
        # Compute solid precipitation from total precipitation
589
        prcp = np.ones(npix) * iprcp
2✔
590
        fac = 1 - (temp - self.t_solid) / (self.t_liq - self.t_solid)
2✔
591
        prcpsol = prcp * clip_array(fac, 0, 1)
2✔
592

593
        return temp, tempformelt, prcp, prcpsol
2✔
594

595
    def _get_2d_annual_climate(self, heights, year):
9✔
596
        # Avoid code duplication with a getter routine
597
        year = np.floor(year)
5✔
598
        if self.repeat:
5✔
599
            year = self.ys + (year - self.ys) % (self.ye - self.ys + 1)
1✔
600
        if not self.is_year_valid(year):
5✔
601
            raise ValueError('year {} out of the valid time bounds: '
×
602
                             '[{}, {}]'.format(year, self.ys, self.ye))
603
        pok = np.where(self.years == year)[0]
5✔
604
        if len(pok) < 1:
5✔
605
            raise ValueError('Year {} not in record'.format(int(year)))
×
606

607
        # Read already (temperature bias and precipitation factor corrected!)
608
        itemp = self.temp[pok]
5✔
609
        iprcp = self.prcp[pok]
5✔
610
        igrad = self.grad[pok]
5✔
611

612
        # For each height pixel:
613
        # Compute temp and tempformelt (temperature above melting threshold)
614
        heights = np.asarray(heights)
5✔
615
        npix = len(heights)
5✔
616
        grad_temp = np.atleast_2d(igrad).repeat(npix, 0)
5✔
617
        grad_temp *= (heights.repeat(12).reshape(grad_temp.shape) -
5✔
618
                      self.ref_hgt)
619
        temp2d = np.atleast_2d(itemp).repeat(npix, 0) + grad_temp
5✔
620
        temp2dformelt = temp2d - self.t_melt
5✔
621
        clip_min(temp2dformelt, 0, out=temp2dformelt)
5✔
622

623
        # Compute solid precipitation from total precipitation
624
        prcp = np.atleast_2d(iprcp).repeat(npix, 0)
5✔
625
        fac = 1 - (temp2d - self.t_solid) / (self.t_liq - self.t_solid)
5✔
626
        prcpsol = prcp * clip_array(fac, 0, 1)
5✔
627

628
        return temp2d, temp2dformelt, prcp, prcpsol
5✔
629

630
    def get_annual_climate(self, heights, year=None):
9✔
631
        """Annual climate information at given heights.
632

633
        Note that prcp is corrected with the precipitation factor and that
634
        all other model biases (temp and prcp) are applied.
635

636
        Returns
637
        -------
638
        (temp, tempformelt, prcp, prcpsol)
639
        """
640
        t, tmelt, prcp, prcpsol = self._get_2d_annual_climate(heights, year)
×
641
        return (t.mean(axis=1), tmelt.sum(axis=1),
×
642
                prcp.sum(axis=1), prcpsol.sum(axis=1))
643

644
    def get_monthly_mb(self, heights, year=None, add_climate=False, **kwargs):
9✔
645

646
        t, tmelt, prcp, prcpsol = self.get_monthly_climate(heights, year=year)
2✔
647
        mb_month = prcpsol - self.monthly_melt_f * tmelt
2✔
648
        mb_month -= self.bias * SEC_IN_MONTH / SEC_IN_YEAR
2✔
649
        if add_climate:
2✔
650
            return (mb_month / SEC_IN_MONTH / self.rho, t, tmelt,
1✔
651
                    prcp, prcpsol)
652
        return mb_month / SEC_IN_MONTH / self.rho
2✔
653

654
    def get_annual_mb(self, heights, year=None, add_climate=False, **kwargs):
9✔
655

656
        t, tmelt, prcp, prcpsol = self._get_2d_annual_climate(heights, year)
5✔
657
        mb_annual = np.sum(prcpsol - self.monthly_melt_f * tmelt, axis=1)
5✔
658
        mb_annual = (mb_annual - self.bias) / SEC_IN_YEAR / self.rho
5✔
659
        if add_climate:
5✔
660
            return (mb_annual, t.mean(axis=1), tmelt.sum(axis=1),
1✔
661
                    prcp.sum(axis=1), prcpsol.sum(axis=1))
662
        return mb_annual
5✔
663

664

665
class ConstantMassBalance(MassBalanceModel):
9✔
666
    """Constant mass balance during a chosen period.
667

668
    This is useful for equilibrium experiments. Note that is is the "correct"
669
    way to represent the average mass balance over a given period.
670
    See: https://oggm.org/2021/08/05/mean-forcing/
671

672
    Attributes
673
    ----------
674
    y0 : int
675
        the center year of the period
676
    halfsize : int
677
        the halfsize of the period
678
    years : ndarray
679
        the years of the period
680
    """
681

682
    def __init__(self, gdir, mb_model_class=MonthlyTIModel,
9✔
683
                 y0=None, halfsize=15,
684
                 **kwargs):
685
        """Initialize
686

687
        Parameters
688
        ----------
689
        gdir : GlacierDirectory
690
            the glacier directory
691
        mb_model_class : MassBalanceModel class
692
            the MassBalanceModel to use for the constant climate
693
        y0 : int, required
694
            the year at the center of the period of interest.
695
        halfsize : int, optional
696
            the half-size of the time window (window size = 2 * halfsize + 1)
697
        **kwargs:
698
            keyword arguments to pass to the mb_model_class
699
        """
700

701
        super().__init__()
4✔
702
        self.mbmod = mb_model_class(gdir,
4✔
703
                                    **kwargs)
704

705
        if y0 is None:
4✔
706
            raise InvalidParamsError('Please set `y0` explicitly')
×
707

708
        # This is a quick'n dirty optimisation
709
        try:
4✔
710
            fls = gdir.read_pickle('model_flowlines')
4✔
711
            h = []
4✔
712
            for fl in fls:
4✔
713
                # We use bed because of overdeepenings
714
                h = np.append(h, fl.bed_h)
4✔
715
                h = np.append(h, fl.surface_h)
4✔
716
            zminmax = np.round([np.min(h)-50, np.max(h)+2000])
4✔
717
        except FileNotFoundError:
1✔
718
            # in case we don't have them
719
            with ncDataset(gdir.get_filepath('gridded_data')) as nc:
1✔
720
                if np.isfinite(nc.min_h_dem):
1✔
721
                    # a bug sometimes led to non-finite
722
                    zminmax = [nc.min_h_dem-250, nc.max_h_dem+1500]
1✔
723
                else:
724
                    zminmax = [nc.min_h_glacier-1250, nc.max_h_glacier+1500]
×
725
        self.hbins = np.arange(*zminmax, step=10)
4✔
726
        self.valid_bounds = self.hbins[[0, -1]]
4✔
727
        self.y0 = y0
4✔
728
        self.halfsize = halfsize
4✔
729
        self.years = np.arange(y0-halfsize, y0+halfsize+1)
4✔
730
        self.hemisphere = gdir.hemisphere
4✔
731

732
    @property
9✔
733
    def temp_bias(self):
9✔
734
        """Temperature bias to add to the original series."""
735
        return self.mbmod.temp_bias
3✔
736

737
    @temp_bias.setter
9✔
738
    def temp_bias(self, value):
9✔
739
        for attr_name in ['_lazy_interp_yr', '_lazy_interp_m']:
4✔
740
            if hasattr(self, attr_name):
4✔
741
                delattr(self, attr_name)
2✔
742
        self.mbmod.temp_bias = value
4✔
743

744
    @property
9✔
745
    def prcp_fac(self):
9✔
746
        """Precipitation factor to apply to the original series."""
747
        return self.mbmod.prcp_fac
1✔
748

749
    @prcp_fac.setter
9✔
750
    def prcp_fac(self, value):
9✔
751
        for attr_name in ['_lazy_interp_yr', '_lazy_interp_m']:
1✔
752
            if hasattr(self, attr_name):
1✔
753
                delattr(self, attr_name)
1✔
754
        self.mbmod.prcp_fac = value
1✔
755

756
    @property
9✔
757
    def bias(self):
9✔
758
        """Residual bias to apply to the original series."""
759
        return self.mbmod.bias
2✔
760

761
    @bias.setter
9✔
762
    def bias(self, value):
9✔
763
        self.mbmod.bias = value
1✔
764

765
    @lazy_property
9✔
766
    def interp_yr(self):
9✔
767
        # annual MB
768
        mb_on_h = self.hbins*0.
4✔
769
        for yr in self.years:
4✔
770
            mb_on_h += self.mbmod.get_annual_mb(self.hbins, year=yr)
4✔
771
        return interp1d(self.hbins, mb_on_h / len(self.years))
4✔
772

773
    @lazy_property
9✔
774
    def interp_m(self):
9✔
775
        # monthly MB
776
        months = np.arange(12)+1
2✔
777
        interp_m = []
2✔
778
        for m in months:
2✔
779
            mb_on_h = self.hbins*0.
2✔
780
            for yr in self.years:
2✔
781
                yr = date_to_floatyear(yr, m)
2✔
782
                mb_on_h += self.mbmod.get_monthly_mb(self.hbins, year=yr)
2✔
783
            interp_m.append(interp1d(self.hbins, mb_on_h / len(self.years)))
2✔
784
        return interp_m
2✔
785

786
    def is_year_valid(self, year):
9✔
787
        return True
×
788

789
    def get_monthly_climate(self, heights, year=None):
9✔
790
        """Average climate information at given heights.
791

792
        Note that prcp is corrected with the precipitation factor and that
793
        all other biases (precipitation, temp) are applied
794

795
        Returns
796
        -------
797
        (temp, tempformelt, prcp, prcpsol)
798
        """
799
        _, m = floatyear_to_date(year)
2✔
800
        yrs = [date_to_floatyear(y, m) for y in self.years]
2✔
801
        heights = np.atleast_1d(heights)
2✔
802
        nh = len(heights)
2✔
803
        shape = (len(yrs), nh)
2✔
804
        temp = np.zeros(shape)
2✔
805
        tempformelt = np.zeros(shape)
2✔
806
        prcp = np.zeros(shape)
2✔
807
        prcpsol = np.zeros(shape)
2✔
808
        for i, yr in enumerate(yrs):
2✔
809
            t, tm, p, ps = self.mbmod.get_monthly_climate(heights, year=yr)
2✔
810
            temp[i, :] = t
2✔
811
            tempformelt[i, :] = tm
2✔
812
            prcp[i, :] = p
2✔
813
            prcpsol[i, :] = ps
2✔
814
        return (np.mean(temp, axis=0),
2✔
815
                np.mean(tempformelt, axis=0),
816
                np.mean(prcp, axis=0),
817
                np.mean(prcpsol, axis=0))
818

819
    def get_annual_climate(self, heights, year=None):
9✔
820
        """Average climate information at given heights.
821

822
        Note that prcp is corrected with the precipitation factor and that
823
        all other biases (precipitation, temp) are applied
824

825
        Returns
826
        -------
827
        (temp, tempformelt, prcp, prcpsol)
828
        """
829
        yrs = monthly_timeseries(self.years[0], self.years[-1],
2✔
830
                                 include_last_year=True)
831
        heights = np.atleast_1d(heights)
2✔
832
        nh = len(heights)
2✔
833
        shape = (len(yrs), nh)
2✔
834
        temp = np.zeros(shape)
2✔
835
        tempformelt = np.zeros(shape)
2✔
836
        prcp = np.zeros(shape)
2✔
837
        prcpsol = np.zeros(shape)
2✔
838
        for i, yr in enumerate(yrs):
2✔
839
            t, tm, p, ps = self.mbmod.get_monthly_climate(heights, year=yr)
2✔
840
            temp[i, :] = t
2✔
841
            tempformelt[i, :] = tm
2✔
842
            prcp[i, :] = p
2✔
843
            prcpsol[i, :] = ps
2✔
844
        # Note that we do not weight for number of days per month:
845
        # this is consistent with OGGM's calendar
846
        return (np.mean(temp, axis=0),
2✔
847
                np.mean(tempformelt, axis=0) * 12,
848
                np.mean(prcp, axis=0) * 12,
849
                np.mean(prcpsol, axis=0) * 12)
850

851
    def get_monthly_mb(self, heights, year=None, add_climate=False, **kwargs):
9✔
852
        yr, m = floatyear_to_date(year)
2✔
853
        if add_climate:
2✔
854
            t, tmelt, prcp, prcpsol = self.get_monthly_climate(heights, year=year)
2✔
855
            return self.interp_m[m-1](heights), t, tmelt, prcp, prcpsol
2✔
856
        return self.interp_m[m-1](heights)
1✔
857

858
    def get_annual_mb(self, heights, year=None, add_climate=False, **kwargs):
9✔
859
        mb = self.interp_yr(heights)
4✔
860
        if add_climate:
4✔
861
            t, tmelt, prcp, prcpsol = self.get_annual_climate(heights)
1✔
862
            return mb, t, tmelt, prcp, prcpsol
1✔
863
        return mb
4✔
864

865

866
class RandomMassBalance(MassBalanceModel):
9✔
867
    """Random shuffle of all MB years within a given time period.
868

869
    This is useful for finding a possible past glacier state or for sensitivity
870
    experiments.
871

872
    Note that this is going to be sensitive to extreme years in certain
873
    periods, but it is by far more physically reasonable than other
874
    approaches based on gaussian assumptions.
875
    """
876

877
    def __init__(self, gdir, mb_model_class=MonthlyTIModel,
9✔
878
                 y0=None, halfsize=15, seed=None,
879
                 all_years=False, unique_samples=False, prescribe_years=None,
880
                 **kwargs):
881
        """Initialize.
882

883
        Parameters
884
        ----------
885
        gdir : GlacierDirectory
886
            the glacier directory
887
        mb_model_class : MassBalanceModel class
888
            the MassBalanceModel to use for the random shuffle
889
        y0 : int, required
890
            the year at the center of the period of interest.
891
        halfsize : int, optional
892
            the half-size of the time window (window size = 2 * halfsize + 1)
893
        seed : int, optional
894
            Random seed used to initialize the pseudo-random number generator.
895
        all_years : bool
896
            if True, overwrites ``y0`` and ``halfsize`` to use all available
897
            years.
898
        unique_samples: bool
899
            if true, chosen random mass balance years will only be available
900
            once per random climate period-length
901
            if false, every model year will be chosen from the random climate
902
            period with the same probability
903
        prescribe_years : pandas Series
904
            instead of random samples, take a series of (i, y) pairs where
905
            (i) is the simulation year index and (y) is the year to pick in the
906
            original timeseries. Overrides `y0`, `halfsize`, `all_years`,
907
            `unique_samples` and `seed`.
908
        **kwargs:
909
            keyword arguments to pass to the mb_model_class
910
        """
911

912
        super().__init__()
3✔
913
        self.valid_bounds = [-1e4, 2e4]  # in m
3✔
914
        self.mbmod = mb_model_class(gdir, **kwargs)
3✔
915

916
        # Climate period
917
        self.prescribe_years = prescribe_years
3✔
918

919
        if self.prescribe_years is None:
3✔
920
            # Normal stuff
921
            self.rng = np.random.RandomState(seed)
3✔
922
            if all_years:
3✔
923
                self.years = self.mbmod.years
×
924
            else:
925
                if y0 is None:
3✔
926
                    raise InvalidParamsError('Please set `y0` explicitly')
×
927
                self.years = np.arange(y0 - halfsize, y0 + halfsize + 1)
3✔
928
        else:
929
            self.rng = None
1✔
930
            self.years = self.prescribe_years.index
1✔
931

932
        self.yr_range = (self.years[0], self.years[-1] + 1)
3✔
933
        self.ny = len(self.years)
3✔
934
        self.hemisphere = gdir.hemisphere
3✔
935

936
        self._state_yr = dict()
3✔
937

938
        # Sampling without replacement
939
        self.unique_samples = unique_samples
3✔
940
        if self.unique_samples:
3✔
941
            self.sampling_years = self.years
1✔
942

943
    @property
9✔
944
    def temp_bias(self):
9✔
945
        """Temperature bias to add to the original series."""
946
        return self.mbmod.temp_bias
1✔
947

948
    @temp_bias.setter
9✔
949
    def temp_bias(self, value):
9✔
950
        """Temperature bias to add to the original series."""
951
        for attr_name in ['_lazy_interp_yr', '_lazy_interp_m']:
1✔
952
            if hasattr(self, attr_name):
1✔
953
                delattr(self, attr_name)
×
954
        self.mbmod.temp_bias = value
1✔
955

956
    @property
9✔
957
    def prcp_fac(self):
9✔
958
        """Precipitation factor to apply to the original series."""
959
        return self.mbmod.prcp_fac
1✔
960

961
    @prcp_fac.setter
9✔
962
    def prcp_fac(self, value):
9✔
963
        """Precipitation factor to apply to the original series."""
964
        for attr_name in ['_lazy_interp_yr', '_lazy_interp_m']:
1✔
965
            if hasattr(self, attr_name):
1✔
966
                delattr(self, attr_name)
×
967
        self.mbmod.prcp_fac = value
1✔
968

969
    @property
9✔
970
    def bias(self):
9✔
971
        """Residual bias to apply to the original series."""
972
        return self.mbmod.bias
1✔
973

974
    @bias.setter
9✔
975
    def bias(self, value):
9✔
976
        """Residual bias to apply to the original series."""
977
        self.mbmod.bias = value
1✔
978

979
    def is_year_valid(self, year):
9✔
980
        return True
×
981

982
    def get_state_yr(self, year=None):
9✔
983
        """For a given year, get the random year associated to it."""
984
        year = int(year)
3✔
985
        if year not in self._state_yr:
3✔
986
            if self.prescribe_years is not None:
3✔
987
                self._state_yr[year] = self.prescribe_years.loc[year]
1✔
988
            else:
989
                if self.unique_samples:
3✔
990
                    # --- Sampling without replacement ---
991
                    if self.sampling_years.size == 0:
1✔
992
                        # refill sample pool when all years were picked once
993
                        self.sampling_years = self.years
×
994
                    # choose one year which was not used in the current period
995
                    _sample = self.rng.choice(self.sampling_years)
1✔
996
                    # write chosen year to dictionary
997
                    self._state_yr[year] = _sample
1✔
998
                    # update sample pool: remove the chosen year from it
999
                    self.sampling_years = np.delete(
1✔
1000
                        self.sampling_years,
1001
                        np.where(self.sampling_years == _sample))
1002
                else:
1003
                    # --- Sampling with replacement ---
1004
                    self._state_yr[year] = self.rng.randint(*self.yr_range)
3✔
1005
        return self._state_yr[year]
3✔
1006

1007
    def get_monthly_mb(self, heights, year=None, **kwargs):
9✔
1008
        ryr, m = floatyear_to_date(year)
1✔
1009
        ryr = date_to_floatyear(self.get_state_yr(ryr), m)
1✔
1010
        return self.mbmod.get_monthly_mb(heights, year=ryr, **kwargs)
1✔
1011

1012
    def get_annual_mb(self, heights, year=None, **kwargs):
9✔
1013
        ryr = self.get_state_yr(int(year))
3✔
1014
        return self.mbmod.get_annual_mb(heights, year=ryr, **kwargs)
3✔
1015

1016

1017
class UncertainMassBalance(MassBalanceModel):
9✔
1018
    """Adding uncertainty to a mass balance model.
1019

1020
    There are three variables for which you can add uncertainty:
1021
    - temperature (additive bias)
1022
    - precipitation (multiplicative factor)
1023
    - residual (a bias in units of MB)
1024
    """
1025

1026
    def __init__(self, basis_model,
9✔
1027
                 rdn_temp_bias_seed=None, rdn_temp_bias_sigma=0.1,
1028
                 rdn_prcp_fac_seed=None, rdn_prcp_fac_sigma=0.1,
1029
                 rdn_bias_seed=None, rdn_bias_sigma=100):
1030
        """Initialize.
1031

1032
        Parameters
1033
        ----------
1034
        basis_model : MassBalanceModel
1035
            the model to which you want to add the uncertainty to
1036
        rdn_temp_bias_seed : int
1037
            the seed of the random number generator
1038
        rdn_temp_bias_sigma : float
1039
            the standard deviation of the random temperature error
1040
        rdn_prcp_fac_seed : int
1041
            the seed of the random number generator
1042
        rdn_prcp_fac_sigma : float
1043
            the standard deviation of the random precipitation error
1044
            (to be consistent this should be renamed prcp_fac as well)
1045
        rdn_bias_seed : int
1046
            the seed of the random number generator
1047
        rdn_bias_sigma : float
1048
            the standard deviation of the random MB error
1049
        """
1050
        super(UncertainMassBalance, self).__init__()
1✔
1051
        # the aim here is to change temp_bias and prcp_fac so
1052
        self.mbmod = basis_model
1✔
1053
        self.hemisphere = basis_model.hemisphere
1✔
1054
        self.valid_bounds = self.mbmod.valid_bounds
1✔
1055
        self.is_year_valid = self.mbmod.is_year_valid
1✔
1056
        self.rng_temp = np.random.RandomState(rdn_temp_bias_seed)
1✔
1057
        self.rng_prcp = np.random.RandomState(rdn_prcp_fac_seed)
1✔
1058
        self.rng_bias = np.random.RandomState(rdn_bias_seed)
1✔
1059
        self._temp_sigma = rdn_temp_bias_sigma
1✔
1060
        self._prcp_sigma = rdn_prcp_fac_sigma
1✔
1061
        self._bias_sigma = rdn_bias_sigma
1✔
1062
        self._state_temp = dict()
1✔
1063
        self._state_prcp = dict()
1✔
1064
        self._state_bias = dict()
1✔
1065

1066
    @property
9✔
1067
    def temp_bias(self):
9✔
1068
        """Temperature bias to add to the original series."""
1069
        return self.mbmod.temp_bias
×
1070

1071
    @temp_bias.setter
9✔
1072
    def temp_bias(self, value):
9✔
1073
        """Temperature bias to add to the original series."""
1074
        for attr_name in ['_lazy_interp_yr', '_lazy_interp_m']:
×
1075
            if hasattr(self, attr_name):
×
1076
                delattr(self, attr_name)
×
1077
        self.mbmod.temp_bias = value
×
1078

1079
    @property
9✔
1080
    def prcp_fac(self):
9✔
1081
        """Precipitation factor to apply to the original series."""
1082
        return self.mbmod.prcp_fac
×
1083

1084
    @prcp_fac.setter
9✔
1085
    def prcp_fac(self, value):
9✔
1086
        """Precipitation factor to apply to the original series."""
1087
        self.mbmod.prcp_fac = value
×
1088

1089
    def _get_state_temp(self, year):
9✔
1090
        year = int(year)
1✔
1091
        if year not in self._state_temp:
1✔
1092
            self._state_temp[year] = self.rng_temp.randn() * self._temp_sigma
1✔
1093
        return self._state_temp[year]
1✔
1094

1095
    def _get_state_prcp(self, year):
9✔
1096
        year = int(year)
1✔
1097
        if year not in self._state_prcp:
1✔
1098
            self._state_prcp[year] = self.rng_prcp.randn() * self._prcp_sigma
1✔
1099
        return self._state_prcp[year]
1✔
1100

1101
    def _get_state_bias(self, year):
9✔
1102
        year = int(year)
1✔
1103
        if year not in self._state_bias:
1✔
1104
            self._state_bias[year] = self.rng_bias.randn() * self._bias_sigma
1✔
1105
        return self._state_bias[year]
1✔
1106

1107
    def get_monthly_mb(self, heights, year=None, **kwargs):
9✔
1108
        raise NotImplementedError()
1109

1110
    def get_annual_mb(self, heights, year=None, fl_id=None, **kwargs):
9✔
1111

1112
        # Keep the original biases and add a random error
1113
        _t = self.mbmod.temp_bias
1✔
1114
        _p = self.mbmod.prcp_fac
1✔
1115
        _b = self.mbmod.bias
1✔
1116
        self.mbmod.temp_bias = self._get_state_temp(year) + _t
1✔
1117
        self.mbmod.prcp_fac = self._get_state_prcp(year) + _p
1✔
1118
        self.mbmod.bias = self._get_state_bias(year) + _b
1✔
1119
        try:
1✔
1120
            out = self.mbmod.get_annual_mb(heights, year=year, fl_id=fl_id)
1✔
1121
        except BaseException:
×
1122
            self.mbmod.temp_bias = _t
×
1123
            self.mbmod.prcp_fac = _p
×
1124
            self.mbmod.bias = _b
×
1125
            raise
×
1126
        # Back to normal
1127
        self.mbmod.temp_bias = _t
1✔
1128
        self.mbmod.prcp_fac = _p
1✔
1129
        self.mbmod.bias = _b
1✔
1130
        return out
1✔
1131

1132

1133
class MultipleFlowlineMassBalance(MassBalanceModel):
9✔
1134
    """Handle mass balance at the glacier level instead of flowline level.
1135

1136
    Convenience class doing not much more than wrapping a list of mass balance
1137
    models, one for each flowline.
1138

1139
    This is useful for real-case studies, where each flowline might have
1140
    different model parameters.
1141

1142
    Attributes
1143
    ----------
1144
    fls : list
1145
        list of flowline objects
1146

1147
    """
1148

1149
    def __init__(self, gdir, fls=None, mb_model_class=MonthlyTIModel,
9✔
1150
                 use_inversion_flowlines=False,
1151
                 input_filesuffix='',
1152
                 **kwargs):
1153
        """Initialize.
1154

1155
        Parameters
1156
        ----------
1157
        gdir : GlacierDirectory
1158
            the glacier directory
1159
        fls : list, optional
1160
            list of flowline objects to use (defaults to 'model_flowlines')
1161
        mb_model_class : MassBalanceModel class
1162
            the MassBalanceModel to use (default is MonthlyTIModel,
1163
            alternatives are e.g. ConstantMassBalance...)
1164
        use_inversion_flowlines: bool, optional
1165
            use 'inversion_flowlines' instead of 'model_flowlines'
1166
        kwargs : kwargs to pass to mb_model_class
1167
        """
1168

1169
        # Read in the flowlines
1170
        if use_inversion_flowlines:
4✔
1171
            fls = gdir.read_pickle('inversion_flowlines')
4✔
1172

1173
        if fls is None:
4✔
1174
            try:
4✔
1175
                fls = gdir.read_pickle('model_flowlines')
4✔
1176
            except FileNotFoundError:
1✔
1177
                raise InvalidWorkflowError('Need a valid `model_flowlines` '
1✔
1178
                                           'file. If you explicitly want to '
1179
                                           'use `inversion_flowlines`, set '
1180
                                           'use_inversion_flowlines=True.')
1181

1182
        self.fls = fls
4✔
1183

1184
        # Initialise the mb models
1185
        self.flowline_mb_models = []
4✔
1186
        for fl in self.fls:
4✔
1187
            # Merged glaciers will need different climate files, use filesuffix
1188
            if (fl.rgi_id is not None) and (fl.rgi_id != gdir.rgi_id):
4✔
1189
                rgi_filesuffix = '_' + fl.rgi_id + input_filesuffix
×
1190
            else:
1191
                rgi_filesuffix = input_filesuffix
4✔
1192

1193
            self.flowline_mb_models.append(
4✔
1194
                mb_model_class(gdir, input_filesuffix=rgi_filesuffix,
1195
                               **kwargs))
1196

1197
        self.valid_bounds = self.flowline_mb_models[-1].valid_bounds
4✔
1198
        self.hemisphere = gdir.hemisphere
4✔
1199

1200
    @property
9✔
1201
    def temp_bias(self):
9✔
1202
        """Temperature bias to add to the original series."""
1203
        return self.flowline_mb_models[0].temp_bias
3✔
1204

1205
    @temp_bias.setter
9✔
1206
    def temp_bias(self, value):
9✔
1207
        """Temperature bias to add to the original series."""
1208
        for mbmod in self.flowline_mb_models:
4✔
1209
            mbmod.temp_bias = value
4✔
1210

1211
    @property
9✔
1212
    def prcp_fac(self):
9✔
1213
        """Precipitation factor to apply to the original series."""
1214
        return self.flowline_mb_models[0].prcp_fac
1✔
1215

1216
    @prcp_fac.setter
9✔
1217
    def prcp_fac(self, value):
9✔
1218
        """Precipitation factor to apply to the original series."""
1219
        for mbmod in self.flowline_mb_models:
1✔
1220
            mbmod.prcp_fac = value
1✔
1221

1222
    @property
9✔
1223
    def bias(self):
9✔
1224
        """Residual bias to apply to the original series."""
1225
        return self.flowline_mb_models[0].bias
3✔
1226

1227
    @bias.setter
9✔
1228
    def bias(self, value):
9✔
1229
        """Residual bias to apply to the original series."""
1230
        for mbmod in self.flowline_mb_models:
1✔
1231
            mbmod.bias = value
1✔
1232

1233
    def is_year_valid(self, year):
9✔
1234
        return self.flowline_mb_models[0].is_year_valid(year)
×
1235

1236
    def get_monthly_mb(self, heights, year=None, fl_id=None, **kwargs):
9✔
1237

1238
        if fl_id is None:
2✔
1239
            raise ValueError('`fl_id` is required for '
×
1240
                             'MultipleFlowlineMassBalance!')
1241

1242
        return self.flowline_mb_models[fl_id].get_monthly_mb(heights,
2✔
1243
                                                             year=year,
1244
                                                             **kwargs)
1245

1246
    def get_annual_mb(self, heights, year=None, fl_id=None, **kwargs):
9✔
1247

1248
        if fl_id is None:
4✔
1249
            raise ValueError('`fl_id` is required for '
×
1250
                             'MultipleFlowlineMassBalance!')
1251

1252
        return self.flowline_mb_models[fl_id].get_annual_mb(heights,
4✔
1253
                                                            year=year,
1254
                                                            **kwargs)
1255

1256
    def get_annual_mb_on_flowlines(self, fls=None, year=None):
9✔
1257
        """Get the MB on all points of the glacier at once.
1258

1259
        Parameters
1260
        ----------
1261
        fls: list, optional
1262
            the list of flowlines to get the mass balance from. Defaults
1263
            to self.fls
1264
        year: float, optional
1265
            the time (in the "floating year" convention)
1266
        Returns
1267
        -------
1268
        Tuple of (heights, widths, mass_balance) 1D arrays
1269
        """
1270

1271
        if fls is None:
2✔
1272
            fls = self.fls
2✔
1273

1274
        heights = []
2✔
1275
        widths = []
2✔
1276
        mbs = []
2✔
1277
        for i, fl in enumerate(fls):
2✔
1278
            h = fl.surface_h
2✔
1279
            heights = np.append(heights, h)
2✔
1280
            widths = np.append(widths, fl.widths)
2✔
1281
            mbs = np.append(mbs, self.get_annual_mb(h, year=year, fl_id=i))
2✔
1282

1283
        return heights, widths, mbs
2✔
1284

1285
    def get_specific_mb(self, heights=None, widths=None, fls=None,
9✔
1286
                        year=None):
1287

1288
        if heights is not None or widths is not None:
4✔
1289
            raise ValueError('`heights` and `widths` kwargs do not work with '
×
1290
                             'MultipleFlowlineMassBalance!')
1291

1292
        if fls is None:
4✔
1293
            fls = self.fls
4✔
1294

1295
        if len(np.atleast_1d(year)) > 1:
4✔
1296
            out = [self.get_specific_mb(fls=fls, year=yr) for yr in year]
4✔
1297
            return np.asarray(out)
4✔
1298

1299
        mbs = []
4✔
1300
        widths = []
4✔
1301
        for i, (fl, mb_mod) in enumerate(zip(fls, self.flowline_mb_models)):
4✔
1302
            _widths = fl.widths
4✔
1303
            try:
4✔
1304
                # For rect and parabola don't compute spec mb
1305
                _widths = np.where(fl.thick > 0, _widths, 0)
4✔
1306
            except AttributeError:
4✔
1307
                pass
4✔
1308
            widths = np.append(widths, _widths)
4✔
1309
            mb = mb_mod.get_annual_mb(fl.surface_h, year=year, fls=fls, fl_id=i)
4✔
1310
            mbs = np.append(mbs, mb * SEC_IN_YEAR * mb_mod.rho)
4✔
1311

1312
        return weighted_average_1d(mbs, widths)
4✔
1313

1314
    def get_ela(self, year=None, **kwargs):
9✔
1315

1316
        # ELA here is not without ambiguity.
1317
        # We compute a mean weighted by area.
1318

1319
        if len(np.atleast_1d(year)) > 1:
2✔
1320
            return np.asarray([self.get_ela(year=yr) for yr in year])
1✔
1321

1322
        elas = []
2✔
1323
        areas = []
2✔
1324
        for fl_id, (fl, mb_mod) in enumerate(zip(self.fls,
2✔
1325
                                                 self.flowline_mb_models)):
1326
            elas = np.append(elas, mb_mod.get_ela(year=year, fl_id=fl_id,
2✔
1327
                                                  fls=self.fls))
1328
            areas = np.append(areas, np.sum(fl.widths))
2✔
1329

1330
        return weighted_average_1d(elas, areas)
2✔
1331

1332

1333
def calving_mb(gdir):
9✔
1334
    """Calving mass-loss in specific MB equivalent.
1335

1336
    This is necessary to calibrate the mass balance.
1337
    """
1338

1339
    if not gdir.is_tidewater:
5✔
1340
        return 0.
5✔
1341

1342
    # Ok. Just take the calving rate from cfg and change its units
1343
    # Original units: km3 a-1, to change to mm a-1 (units of specific MB)
1344
    rho = cfg.PARAMS['ice_density']
3✔
1345
    return gdir.inversion_calving_rate * 1e9 * rho / gdir.rgi_area_m2
3✔
1346

1347

1348
def decide_winter_precip_factor(gdir):
9✔
1349
    """Utility function to decide on a precip factor based on winter precip."""
1350

1351
    # We have to decide on a precip factor
1352
    if 'W5E5' not in cfg.PARAMS['baseline_climate']:
1✔
1353
        raise InvalidWorkflowError('prcp_fac from_winter_prcp is only '
×
1354
                                   'compatible with the W5E5 climate '
1355
                                   'dataset!')
1356

1357
    # get non-corrected winter daily mean prcp (kg m-2 day-1)
1358
    # it is easier to get this directly from the raw climate files
1359
    fp = gdir.get_filepath('climate_historical')
1✔
1360
    with xr.open_dataset(fp).prcp as ds_pr:
1✔
1361
        # just select winter months
1362
        if gdir.hemisphere == 'nh':
1✔
1363
            m_winter = [10, 11, 12, 1, 2, 3, 4]
1✔
1364
        else:
1365
            m_winter = [4, 5, 6, 7, 8, 9, 10]
×
1366

1367
        ds_pr_winter = ds_pr.where(ds_pr['time.month'].isin(m_winter), drop=True)
1✔
1368

1369
        # select the correct 41 year time period
1370
        ds_pr_winter = ds_pr_winter.sel(time=slice('1979-01-01', '2019-12-01'))
1✔
1371

1372
        # check if we have the full time period: 41 years * 7 months
1373
        text = ('the climate period has to go from 1979-01 to 2019-12,',
1✔
1374
                'use W5E5 or GSWP3_W5E5 as baseline climate and',
1375
                'repeat the climate processing')
1376
        assert len(ds_pr_winter.time) == 41 * 7, text
1✔
1377
        w_prcp = float((ds_pr_winter / ds_pr_winter.time.dt.daysinmonth).mean())
1✔
1378

1379
    # from MB sandbox calibration to winter MB
1380
    # using t_melt=-1, cte lapse rate, monthly resolution
1381
    a, b = cfg.PARAMS['winter_prcp_fac_ab']
1✔
1382
    prcp_fac = a * np.log(w_prcp) + b
1✔
1383
    # don't allow extremely low/high prcp. factors!!!
1384
    return clip_scalar(prcp_fac,
1✔
1385
                       cfg.PARAMS['prcp_fac_min'],
1386
                       cfg.PARAMS['prcp_fac_max'])
1387

1388

1389
@entity_task(log, writes=['mb_calib'])
9✔
1390
def mb_calibration_from_wgms_mb(gdir, **kwargs):
9✔
1391
    """Calibrate for in-situ, annual MB.
1392

1393
    This only works for glaciers which have WGMS data!
1394

1395
    For now this just calls mb_calibration_from_scalar_mb internally,
1396
    but could be cleverer than that if someone wishes to implement it.
1397

1398
    Parameters
1399
    ----------
1400
    **kwargs : any kwarg accepted by mb_calibration_from_scalar_mb
1401
    """
1402

1403
    # Note that this currently does not work for hydro years (WGMS uses hydro)
1404
    # A way to go would be to teach the mb models to use calendar years
1405
    # internally but still output annual MB in hydro convention.
1406
    mbdf = gdir.get_ref_mb_data()
1✔
1407
    # Keep only valid values
1408
    mbdf = mbdf.loc[~mbdf['ANNUAL_BALANCE'].isnull()]
1✔
1409
    return mb_calibration_from_scalar_mb(gdir,
1✔
1410
                                         ref_mb=mbdf['ANNUAL_BALANCE'].mean(),
1411
                                         ref_mb_years=mbdf.index.values,
1412
                                         **kwargs)
1413

1414

1415
@entity_task(log, writes=['mb_calib'])
9✔
1416
def mb_calibration_from_geodetic_mb(gdir, *,
9✔
1417
                                    ref_period=None,
1418
                                    write_to_gdir=True,
1419
                                    overwrite_gdir=False,
1420
                                    override_missing=None,
1421
                                    informed_threestep=False,
1422
                                    calibrate_param1='melt_f',
1423
                                    calibrate_param2=None,
1424
                                    calibrate_param3=None,
1425
                                    mb_model_class=MonthlyTIModel,
1426
                                    ):
1427
    """Calibrate for geodetic MB data from Hugonnet et al., 2021.
1428

1429
    The data table can be obtained with utils.get_geodetic_mb_dataframe().
1430
    It is equivalent to the original data from Hugonnet, but has some outlier
1431
    values filtered. See `this notebook` for more details.
1432

1433
    The problem of calibrating many unknown parameters on geodetic data is
1434
    currently unsolved. This is OGGM's current take, based on trial and
1435
    error and based on ideas from the literature.
1436

1437
    Parameters
1438
    ----------
1439
    gdir : :py:class:`oggm.GlacierDirectory`
1440
        the glacier directory to calibrate
1441
    ref_period : str, default: PARAMS['geodetic_mb_period']
1442
        one of '2000-01-01_2010-01-01', '2010-01-01_2020-01-01',
1443
        '2000-01-01_2020-01-01'. If `ref_mb` is set, this should still match
1444
        the same format but can be any date.
1445
    write_to_gdir : bool
1446
        whether to write the results of the calibration to the glacier
1447
        directory. If True (the default), this will be saved as `mb_calib.json`
1448
        and be used by the MassBalanceModel class as parameters in subsequent
1449
        tasks.
1450
    overwrite_gdir : bool
1451
        if a `mb_calib.json` exists, this task won't overwrite it per default.
1452
        Set this to True to enforce overwriting (i.e. with consequences for the
1453
        future workflow).
1454
    override_missing : scalar
1455
        if the reference geodetic data is not available, use this value instead
1456
        (mostly for testing with exotic datasets, but could be used to open
1457
        the door to using other datasets).
1458
    informed_threestep : bool
1459
        the magic method Fabi found out one day before release.
1460
        Overrides the calibrate_param order below.
1461
    calibrate_param1 : str
1462
        in the three-step calibration, the name of the first parameter
1463
        to calibrate (one of 'melt_f', 'temp_bias', 'prcp_fac').
1464
    calibrate_param2 : str
1465
        in the three-step calibration, the name of the second parameter
1466
        to calibrate (one of 'melt_f', 'temp_bias', 'prcp_fac'). If not
1467
        set and the algorithm cannot match observations, it will raise an
1468
        error.
1469
    calibrate_param3 : str
1470
        in the three-step calibration, the name of the third parameter
1471
        to calibrate (one of 'melt_f', 'temp_bias', 'prcp_fac'). If not
1472
        set and the algorithm cannot match observations, it will raise an
1473
        error.
1474
    mb_model_class : MassBalanceModel class
1475
        the MassBalanceModel to use for the calibration. Needs to use the
1476
        same parameters as MonthlyTIModel (the default): melt_f,
1477
        temp_bias, prcp_fac.
1478

1479
    Returns
1480
    -------
1481
    the calibrated parameters as dict
1482
    """
1483
    if not ref_period:
4✔
1484
        ref_period = cfg.PARAMS['geodetic_mb_period']
4✔
1485

1486
    # Get the reference data
1487
    ref_mb_err = np.NaN
4✔
1488
    try:
4✔
1489
        ref_mb_df = get_geodetic_mb_dataframe().loc[gdir.rgi_id]
4✔
1490
        ref_mb_df = ref_mb_df.loc[ref_mb_df['period'] == ref_period]
4✔
1491
        # dmdtda: in meters water-equivalent per year -> we convert to kg m-2 yr-1
1492
        ref_mb = float(ref_mb_df['dmdtda']) * 1000
4✔
1493
        ref_mb_err = float(ref_mb_df['err_dmdtda']) * 1000
4✔
1494
    except KeyError:
×
1495
        if override_missing is None:
×
1496
            raise
×
1497
        ref_mb = override_missing
×
1498

1499
    temp_bias = 0
4✔
1500
    if cfg.PARAMS['use_temp_bias_from_file']:
4✔
1501
        climinfo = gdir.get_climate_info()
1✔
1502
        if 'w5e5' not in climinfo['baseline_climate_source'].lower():
1✔
1503
            raise InvalidWorkflowError('use_temp_bias_from_file currently '
×
1504
                                       'only available for W5E5 data.')
1505
        bias_df = get_temp_bias_dataframe()
1✔
1506
        ref_lon = climinfo['baseline_climate_ref_pix_lon']
1✔
1507
        ref_lat = climinfo['baseline_climate_ref_pix_lat']
1✔
1508
        # Take nearest
1509
        dis = ((bias_df.lon_val - ref_lon)**2 + (bias_df.lat_val - ref_lat)**2)**0.5
1✔
1510
        sel_df = bias_df.iloc[np.argmin(dis)]
1✔
1511
        temp_bias = sel_df['median_temp_bias_w_err_grouped']
1✔
1512
        assert np.isfinite(temp_bias), 'Temp bias not finite?'
1✔
1513

1514
    if informed_threestep:
4✔
1515
        if not cfg.PARAMS['use_temp_bias_from_file']:
1✔
1516
            raise InvalidParamsError('With `informed_threestep` you need to '
×
1517
                                     'set `use_temp_bias_from_file`.')
1518
        if not cfg.PARAMS['use_winter_prcp_fac']:
1✔
1519
            raise InvalidParamsError('With `informed_threestep` you need to '
×
1520
                                     'set `use_winter_prcp_fac`.')
1521

1522
        # Some magic heuristics - we just decide to calibrate
1523
        # precip -> melt_f -> temp but informed by previous data.
1524

1525
        # Temp bias was decided anyway, we keep as previous value and
1526
        # allow it to vary as last resort
1527

1528
        # We use the precip factor but allow it to vary between 0.8, 1.2 of
1529
        # the previous value (uncertainty).
1530
        prcp_fac = decide_winter_precip_factor(gdir)
1✔
1531
        mi, ma = cfg.PARAMS['prcp_fac_min'], cfg.PARAMS['prcp_fac_max']
1✔
1532
        prcp_fac_min = clip_scalar(prcp_fac * 0.8, mi, ma)
1✔
1533
        prcp_fac_max = clip_scalar(prcp_fac * 1.2, mi, ma)
1✔
1534

1535
        return mb_calibration_from_scalar_mb(gdir,
1✔
1536
                                             ref_mb=ref_mb,
1537
                                             ref_mb_err=ref_mb_err,
1538
                                             ref_period=ref_period,
1539
                                             write_to_gdir=write_to_gdir,
1540
                                             overwrite_gdir=overwrite_gdir,
1541
                                             calibrate_param1='prcp_fac',
1542
                                             calibrate_param2='melt_f',
1543
                                             calibrate_param3='temp_bias',
1544
                                             prcp_fac=prcp_fac,
1545
                                             prcp_fac_min=prcp_fac_min,
1546
                                             prcp_fac_max=prcp_fac_max,
1547
                                             temp_bias=temp_bias,
1548
                                             mb_model_class=mb_model_class,
1549
                                             )
1550

1551
    else:
1552
        return mb_calibration_from_scalar_mb(gdir,
4✔
1553
                                             ref_mb=ref_mb,
1554
                                             ref_mb_err=ref_mb_err,
1555
                                             ref_period=ref_period,
1556
                                             write_to_gdir=write_to_gdir,
1557
                                             overwrite_gdir=overwrite_gdir,
1558
                                             calibrate_param1=calibrate_param1,
1559
                                             calibrate_param2=calibrate_param2,
1560
                                             calibrate_param3=calibrate_param3,
1561
                                             temp_bias=temp_bias,
1562
                                             mb_model_class=mb_model_class,
1563
                                             )
1564

1565

1566
@entity_task(log, writes=['mb_calib'])
9✔
1567
def mb_calibration_from_scalar_mb(gdir, *,
9✔
1568
                                  ref_mb=None,
1569
                                  ref_mb_err=None,
1570
                                  ref_period=None,
1571
                                  ref_mb_years=None,
1572
                                  write_to_gdir=True,
1573
                                  overwrite_gdir=False,
1574
                                  calibrate_param1='melt_f',
1575
                                  calibrate_param2=None,
1576
                                  calibrate_param3=None,
1577
                                  melt_f=None,
1578
                                  melt_f_min=None,
1579
                                  melt_f_max=None,
1580
                                  prcp_fac=None,
1581
                                  prcp_fac_min=None,
1582
                                  prcp_fac_max=None,
1583
                                  temp_bias=None,
1584
                                  temp_bias_min=None,
1585
                                  temp_bias_max=None,
1586
                                  mb_model_class=MonthlyTIModel,
1587
                                  ):
1588
    """Determine the mass balance parameters from a scalar mass-balance value.
1589

1590
    This calibrates the mass balance parameters using a reference average
1591
    MB data over a given period (e.g. average in-situ SMB or geodetic MB).
1592
    This flexible calibration allows to calibrate three parameters one after
1593
    another. The first parameter is varied between two chosen values (a range)
1594
    until the ref MB value is matched. If this fails, the second parameter
1595
    can be changed, etc.
1596

1597
    This can be used for example to apply the "three-step calibration"
1598
    introduced by Huss & Hock 2015, but you can choose any order of
1599
    calibration.
1600

1601
    This task can be called by other, "higher level" tasks, for example
1602
    :py:func:`oggm.core.massbalance.mb_calibration_from_geodetic_mb` or
1603
    :py:func:`oggm.core.massbalance.mb_calibration_from_wgms_mb`.
1604

1605
    Note that this does not compute the apparent mass balance at
1606
    the same time - users need to run `apparent_mb_from_any_mb after`
1607
    calibration.
1608

1609
    Parameters
1610
    ----------
1611
    gdir : :py:class:`oggm.GlacierDirectory`
1612
        the glacier directory to calibrate
1613
    ref_mb : float, required
1614
        the reference mass balance to match (units: kg m-2 yr-1)
1615
        It is required here - if you want to use available observations,
1616
        use :py:func:`oggm.core.massbalance.mb_calibration_from_geodetic_mb`
1617
        or :py:func:`oggm.core.massbalance.mb_calibration_from_wgms_mb`.
1618
    ref_mb_err : float, optional
1619
        currently only used for logging - it is not used in the calibration.
1620
    ref_period : str, optional
1621
        date format - for example '2000-01-01_2010-01-01'. If this is not
1622
        set, ref_mb_years needs to be set.
1623
    ref_mb_years : tuple of length 2 (range) or list of years.
1624
        convenience kwarg to override ref_period. If a tuple of length 2 is
1625
        given, all years between this range (excluding the last one) are used.
1626
        If a list  of years is given, all these will be used (useful for
1627
        data with gaps)
1628
    write_to_gdir : bool
1629
        whether to write the results of the calibration to the glacier
1630
        directory. If True (the default), this will be saved as `mb_calib.json`
1631
        and be used by the MassBalanceModel class as parameters in subsequent
1632
        tasks.
1633
    overwrite_gdir : bool
1634
        if a `mb_calib.json` exists, this task won't overwrite it per default.
1635
        Set this to True to enforce overwriting (i.e. with consequences for the
1636
        future workflow).
1637
    mb_model_class : MassBalanceModel class
1638
        the MassBalanceModel to use for the calibration. Needs to use the
1639
        same parameters as MonthlyTIModel (the default): melt_f,
1640
        temp_bias, prcp_fac.
1641
    calibrate_param1 : str
1642
        in the three-step calibration, the name of the first parameter
1643
        to calibrate (one of 'melt_f', 'temp_bias', 'prcp_fac').
1644
    calibrate_param2 : str
1645
        in the three-step calibration, the name of the second parameter
1646
        to calibrate (one of 'melt_f', 'temp_bias', 'prcp_fac'). If not
1647
        set and the algorithm cannot match observations, it will raise an
1648
        error.
1649
    calibrate_param3 : str
1650
        in the three-step calibration, the name of the third parameter
1651
        to calibrate (one of 'melt_f', 'temp_bias', 'prcp_fac'). If not
1652
        set and the algorithm cannot match observations, it will raise an
1653
        error.
1654
    melt_f: float
1655
        the default value to use as melt factor (or the starting value when
1656
        optimizing MB). Defaults to cfg.PARAMS['melt_f'].
1657
    melt_f_min: float
1658
        the minimum accepted value for the melt factor during optimisation.
1659
        Defaults to cfg.PARAMS['melt_f_min'].
1660
    melt_f_max: float
1661
        the maximum accepted value for the melt factor during optimisation.
1662
        Defaults to cfg.PARAMS['melt_f_max'].
1663
    prcp_fac: float
1664
        the default value to use as precipitation scaling factor
1665
        (or the starting value when optimizing MB). Defaults to the method
1666
        chosen in `params.cfg` (winter prcp or global factor).
1667
    prcp_fac_min: float
1668
        the minimum accepted value for the precipitation scaling factor during
1669
        optimisation. Defaults to cfg.PARAMS['prcp_fac_min'].
1670
    prcp_fac_max: float
1671
        the maximum accepted value for the precipitation scaling factor during
1672
        optimisation. Defaults to cfg.PARAMS['prcp_fac_max'].
1673
    temp_bias: float
1674
        the default value to use as temperature bias (or the starting value when
1675
        optimizing MB). Defaults to 0.
1676
    temp_bias_min: float
1677
        the minimum accepted value for the temperature bias during optimisation.
1678
        Defaults to cfg.PARAMS['temp_bias_min'].
1679
    temp_bias_max: float
1680
        the maximum accepted value for the temperature bias during optimisation.
1681
        Defaults to cfg.PARAMS['temp_bias_max'].
1682
    """
1683

1684
    # Param constraints
1685
    if melt_f_min is None:
5✔
1686
        melt_f_min = cfg.PARAMS['melt_f_min']
5✔
1687
    if melt_f_max is None:
5✔
1688
        melt_f_max = cfg.PARAMS['melt_f_max']
5✔
1689
    if prcp_fac_min is None:
5✔
1690
        prcp_fac_min = cfg.PARAMS['prcp_fac_min']
5✔
1691
    if prcp_fac_max is None:
5✔
1692
        prcp_fac_max = cfg.PARAMS['prcp_fac_max']
5✔
1693
    if temp_bias_min is None:
5✔
1694
        temp_bias_min = cfg.PARAMS['temp_bias_min']
5✔
1695
    if temp_bias_max is None:
5✔
1696
        temp_bias_max = cfg.PARAMS['temp_bias_max']
5✔
1697
    if ref_mb_years is not None and ref_period is not None:
5✔
1698
        raise InvalidParamsError('Cannot set `ref_mb_years` and `ref_period` '
×
1699
                                 'at the same time.')
1700

1701
    fls = gdir.read_pickle('inversion_flowlines')
5✔
1702

1703
    # Let's go
1704
    # Climate period
1705
    if ref_mb_years is not None:
5✔
1706
        if len(ref_mb_years) > 2:
3✔
1707
            years = np.asarray(ref_mb_years)
1✔
1708
            ref_period = 'custom'
1✔
1709
        else:
1710
            years = np.arange(*ref_mb_years)
2✔
1711
            ref_period = f'{ref_mb_years[0]}-01-01_{ref_mb_years[1]}-01-01'
2✔
1712
    elif ref_period is not None:
5✔
1713
        y0, y1 = ref_period.split('_')
5✔
1714
        y0 = int(y0.split('-')[0])
5✔
1715
        y1 = int(y1.split('-')[0])
5✔
1716
        years = np.arange(y0, y1)
5✔
1717
    else:
1718
        raise InvalidParamsError('One of `ref_mb_years` or `ref_period` '
×
1719
                                 'is required for calibration.')
1720

1721
    # Do we have a calving glacier?
1722
    cmb = calving_mb(gdir)
5✔
1723
    if cmb != 0:
5✔
1724
        raise NotImplementedError('Calving with geodetic MB is not implemented '
1725
                                  'yet, but it should actually work. Well keep '
1726
                                  'you posted!')
1727

1728
    # Ok, regardless on how we want to calibrate, we start with defaults
1729
    if melt_f is None:
5✔
1730
        melt_f = cfg.PARAMS['melt_f']
5✔
1731

1732
    if prcp_fac is None:
5✔
1733
        if cfg.PARAMS['use_winter_prcp_fac']:
5✔
1734
            # Some sanity check
1735
            if cfg.PARAMS['prcp_fac'] is not None:
×
1736
                raise InvalidWorkflowError("Set PARAMS['prcp_fac'] to None "
×
1737
                                           "if using PARAMS['winter_prcp_factor'].")
1738
            prcp_fac = decide_winter_precip_factor(gdir)
×
1739
        else:
1740
            prcp_fac = cfg.PARAMS['prcp_fac']
5✔
1741
            if prcp_fac is None:
5✔
1742
                raise InvalidWorkflowError("Set either PARAMS['use_winter_prcp_fac'] "
×
1743
                                           "or PARAMS['winter_prcp_factor'].")
1744

1745
    if temp_bias is None:
5✔
1746
        temp_bias = 0
4✔
1747

1748
    # Create the MB model we will calibrate
1749
    mb_mod = mb_model_class(gdir,
5✔
1750
                            melt_f=melt_f,
1751
                            temp_bias=temp_bias,
1752
                            prcp_fac=prcp_fac,
1753
                            check_calib_params=False)
1754

1755
    # Check that the years are available
1756
    for y in years:
5✔
1757
        if not mb_mod.is_year_valid(y):
5✔
1758
            raise ValueError(f'year {y} out of the valid time bounds: '
×
1759
                             f'[{mb_mod.ys}, {mb_mod.ye}]')
1760

1761
    if calibrate_param1 == 'melt_f':
5✔
1762
        min_range, max_range = melt_f_min, melt_f_max
5✔
1763
    elif calibrate_param1 == 'prcp_fac':
3✔
1764
        min_range, max_range = prcp_fac_min, prcp_fac_max
2✔
1765
    elif calibrate_param1 == 'temp_bias':
2✔
1766
        min_range, max_range = temp_bias_min, temp_bias_max
2✔
1767
    else:
1768
        raise InvalidParamsError("calibrate_param1 must be one of "
×
1769
                                 "['melt_f', 'prcp_fac', 'temp_bias']")
1770

1771
    def to_minimize(x, model_attr):
5✔
1772
        # Set the new attr value
1773
        setattr(mb_mod, model_attr, x)
5✔
1774
        out = mb_mod.get_specific_mb(fls=fls, year=years).mean()
5✔
1775
        return np.mean(out - ref_mb)
5✔
1776

1777
    try:
5✔
1778
        optim_param1 = optimize.brentq(to_minimize,
5✔
1779
                                       min_range, max_range,
1780
                                       args=(calibrate_param1,)
1781
                                       )
1782
    except ValueError:
2✔
1783
        if not calibrate_param2:
2✔
1784
            raise RuntimeError(f'{gdir.rgi_id}: ref mb not matched. '
1✔
1785
                               f'Try to set calibrate_param2.')
1786

1787
        # Check which direction we need to go
1788
        diff_1 = to_minimize(min_range, calibrate_param1)
2✔
1789
        diff_2 = to_minimize(max_range, calibrate_param1)
2✔
1790
        optim_param1 = min_range if abs(diff_1) < abs(diff_2) else max_range
2✔
1791
        setattr(mb_mod, calibrate_param1, optim_param1)
2✔
1792

1793
        # Second step
1794
        if calibrate_param2 == 'melt_f':
2✔
1795
            min_range, max_range = melt_f_min, melt_f_max
1✔
1796
        elif calibrate_param2 == 'prcp_fac':
1✔
1797
            min_range, max_range = prcp_fac_min, prcp_fac_max
×
1798
        elif calibrate_param2 == 'temp_bias':
1✔
1799
            min_range, max_range = temp_bias_min, temp_bias_max
1✔
1800
        else:
1801
            raise InvalidParamsError("calibrate_param2 must be one of "
×
1802
                                     "['melt_f', 'prcp_fac', 'temp_bias']")
1803
        try:
2✔
1804
            optim_param2 = optimize.brentq(to_minimize,
2✔
1805
                                           min_range, max_range,
1806
                                           args=(calibrate_param2,)
1807
                                           )
1808
        except ValueError:
1✔
1809
            # Third step
1810
            if not calibrate_param3:
1✔
1811
                raise RuntimeError(f'{gdir.rgi_id}: ref mb not matched. '
1✔
1812
                                   f'Try to set calibrate_param3.')
1813

1814
            # Check which direction we need to go
1815
            diff_1 = to_minimize(min_range, calibrate_param2)
1✔
1816
            diff_2 = to_minimize(max_range, calibrate_param2)
1✔
1817
            optim_param2 = min_range if abs(diff_1) < abs(diff_2) else max_range
1✔
1818
            setattr(mb_mod, calibrate_param2, optim_param2)
1✔
1819

1820
            # Third step
1821
            if calibrate_param3 == 'melt_f':
1✔
1822
                min_range, max_range = melt_f_min, melt_f_max
1✔
1823
            elif calibrate_param3 == 'prcp_fac':
×
1824
                min_range, max_range = prcp_fac_min, prcp_fac_max
×
1825
            elif calibrate_param3 == 'temp_bias':
×
1826
                min_range, max_range = temp_bias_min, temp_bias_max
×
1827
            else:
1828
                raise InvalidParamsError("calibrate_param3 must be one of "
×
1829
                                         "['melt_f', 'prcp_fac', 'temp_bias']")
1830
            try:
1✔
1831
                optim_param3 = optimize.brentq(to_minimize,
1✔
1832
                                               min_range, max_range,
1833
                                               args=(calibrate_param3,)
1834
                                               )
1835
            except ValueError:
1✔
1836
                raise RuntimeError(f'{gdir.rgi_id}: we tried very hard but we '
1✔
1837
                                   f'could not find a combination of '
1838
                                   f'parameters that works for this ref mb.')
1839

1840
            if calibrate_param3 == 'melt_f':
1✔
1841
                melt_f = optim_param3
1✔
1842
            elif calibrate_param3 == 'prcp_fac':
×
1843
                prcp_fac = optim_param3
×
1844
            elif calibrate_param3 == 'temp_bias':
×
1845
                temp_bias = optim_param3
×
1846

1847
        if calibrate_param2 == 'melt_f':
2✔
1848
            melt_f = optim_param2
1✔
1849
        elif calibrate_param2 == 'prcp_fac':
1✔
1850
            prcp_fac = optim_param2
×
1851
        elif calibrate_param2 == 'temp_bias':
1✔
1852
            temp_bias = optim_param2
1✔
1853

1854
    if calibrate_param1 == 'melt_f':
5✔
1855
        melt_f = optim_param1
5✔
1856
    elif calibrate_param1 == 'prcp_fac':
3✔
1857
        prcp_fac = optim_param1
2✔
1858
    elif calibrate_param1 == 'temp_bias':
2✔
1859
        temp_bias = optim_param1
2✔
1860

1861
    # Store parameters
1862
    df = gdir.read_json('mb_calib', allow_empty=True)
5✔
1863
    df['rgi_id'] = gdir.rgi_id
5✔
1864
    df['bias'] = 0
5✔
1865
    df['melt_f'] = melt_f
5✔
1866
    df['prcp_fac'] = prcp_fac
5✔
1867
    df['temp_bias'] = temp_bias
5✔
1868
    # What did we try to match?
1869
    df['reference_mb'] = ref_mb
5✔
1870
    df['reference_mb_err'] = ref_mb_err
5✔
1871
    df['reference_period'] = ref_period
5✔
1872

1873
    # Add the climate related params to the GlacierDir to make sure
1874
    # other tools cannot fool around without re-calibration
1875
    df['mb_global_params'] = {k: cfg.PARAMS[k] for k in MB_GLOBAL_PARAMS}
5✔
1876
    df['baseline_climate_source'] = gdir.get_climate_info()['baseline_climate_source']
5✔
1877
    # Write
1878
    if write_to_gdir:
5✔
1879
        if gdir.has_file('mb_calib') and not overwrite_gdir:
5✔
1880
            raise InvalidWorkflowError('`mb_calib.json` already exists for '
×
1881
                                       'this repository. Set `overwrite_gdir` '
1882
                                       'to True if you want to overwrite '
1883
                                       'a previous calibration.')
1884
        gdir.write_json(df, 'mb_calib')
5✔
1885
    return df
5✔
1886

1887

1888
def _check_terminus_mass_flux(gdir, fls):
9✔
1889
    # Check that we have done this correctly
1890
    rho = cfg.PARAMS['ice_density']
5✔
1891
    cmb = calving_mb(gdir)
5✔
1892

1893
    # This variable is in "sensible" units normalized by width
1894
    flux = fls[-1].flux_out
5✔
1895
    aflux = flux * (gdir.grid.dx ** 2) / rho * 1e-9  # km3 ice per year
5✔
1896

1897
    # If not marine and a bit far from zero, warning
1898
    if cmb == 0 and not np.allclose(flux, 0, atol=0.01):
5✔
1899
        log.info('(%s) flux should be zero, but is: '
×
1900
                 '%.4f km3 ice yr-1', gdir.rgi_id, aflux)
1901

1902
    # If not marine and quite far from zero, error
1903
    if cmb == 0 and not np.allclose(flux, 0, atol=1):
5✔
1904
        msg = ('({}) flux should be zero, but is: {:.4f} km3 ice yr-1'
×
1905
               .format(gdir.rgi_id, aflux))
1906
        raise MassBalanceCalibrationError(msg)
×
1907

1908

1909
@entity_task(log, writes=['inversion_flowlines', 'linear_mb_params'])
9✔
1910
def apparent_mb_from_linear_mb(gdir, mb_gradient=3., ela_h=None):
9✔
1911
    """Compute apparent mb from a linear mass balance assumption (for testing).
1912

1913
    This is for testing currently, but could be used as alternative method
1914
    for the inversion quite easily.
1915

1916
    Parameters
1917
    ----------
1918
    gdir : :py:class:`oggm.GlacierDirectory`
1919
        the glacier directory to process
1920
    """
1921

1922
    # Do we have a calving glacier?
1923
    cmb = calving_mb(gdir)
3✔
1924
    is_calving = cmb != 0.
3✔
1925

1926
    # Get the height and widths along the fls
1927
    h, w = gdir.get_inversion_flowline_hw()
3✔
1928

1929
    # Now find the ELA till the integrated mb is zero
1930
    from oggm.core.massbalance import LinearMassBalance
3✔
1931

1932
    def to_minimize(ela_h):
3✔
1933
        mbmod = LinearMassBalance(ela_h, grad=mb_gradient)
3✔
1934
        smb = mbmod.get_specific_mb(heights=h, widths=w)
3✔
1935
        return smb - cmb
3✔
1936

1937
    if ela_h is None:
3✔
1938
        ela_h = optimize.brentq(to_minimize, -1e5, 1e5)
3✔
1939

1940
    # For each flowline compute the apparent MB
1941
    rho = cfg.PARAMS['ice_density']
3✔
1942
    fls = gdir.read_pickle('inversion_flowlines')
3✔
1943
    # Reset flux
1944
    for fl in fls:
3✔
1945
        fl.flux = np.zeros(len(fl.surface_h))
3✔
1946
    # Flowlines in order to be sure
1947
    mbmod = LinearMassBalance(ela_h, grad=mb_gradient)
3✔
1948
    for fl in fls:
3✔
1949
        mbz = mbmod.get_annual_mb(fl.surface_h) * cfg.SEC_IN_YEAR * rho
3✔
1950
        fl.set_apparent_mb(mbz, is_calving=is_calving)
3✔
1951

1952
    # Check and write
1953
    _check_terminus_mass_flux(gdir, fls)
3✔
1954
    gdir.write_pickle(fls, 'inversion_flowlines')
3✔
1955
    gdir.write_pickle({'ela_h': ela_h, 'grad': mb_gradient},
3✔
1956
                      'linear_mb_params')
1957

1958

1959
@entity_task(log, writes=['inversion_flowlines'])
9✔
1960
def apparent_mb_from_any_mb(gdir, mb_model=None,
9✔
1961
                            mb_model_class=MonthlyTIModel,
1962
                            mb_years=None):
1963
    """Compute apparent mb from an arbitrary mass balance profile.
1964

1965
    This searches for a mass balance residual to add to the mass balance
1966
    profile so that the average specific MB is zero.
1967

1968
    Parameters
1969
    ----------
1970
    gdir : :py:class:`oggm.GlacierDirectory`
1971
        the glacier directory to process
1972
    mb_model : :py:class:`oggm.core.massbalance.MassBalanceModel`
1973
        the mass balance model to use - if None, will use the
1974
        one given by mb_model_class
1975
    mb_model_class : MassBalanceModel class
1976
        the MassBalanceModel class to use, default is MonthlyTIModel
1977
    mb_years : array, or tuple of length 2 (range)
1978
        the array of years from which you want to average the MB for (for
1979
        mb_model only). If an array of length 2 is given, all years
1980
        between this range (excluding the last one) are used.
1981
        Default is to pick all years from the reference
1982
        geodetic MB period, i.e. PARAMS['geodetic_mb_period'].
1983
        It does not matter much for the final result, but it should be a
1984
        period long enough to have a representative MB gradient.
1985
    """
1986

1987
    # Do we have a calving glacier?
1988
    cmb = calving_mb(gdir)
5✔
1989
    is_calving = cmb != 0
5✔
1990

1991
    # For each flowline compute the apparent MB
1992
    fls = gdir.read_pickle('inversion_flowlines')
5✔
1993

1994
    if mb_model is None:
5✔
1995
        mb_model = mb_model_class(gdir)
5✔
1996

1997
    if mb_years is None:
5✔
1998
        mb_years = cfg.PARAMS['geodetic_mb_period']
5✔
1999
        y0, y1 = mb_years.split('_')
5✔
2000
        y0 = int(y0.split('-')[0])
5✔
2001
        y1 = int(y1.split('-')[0])
5✔
2002
        mb_years = np.arange(y0, y1, 1)
5✔
2003

2004
    if len(mb_years) == 2:
5✔
2005
        # Range
2006
        mb_years = np.arange(*mb_years, 1)
5✔
2007

2008
    # Unchanged SMB
2009
    o_smb = np.mean(mb_model.get_specific_mb(fls=fls, year=mb_years))
5✔
2010

2011
    def to_minimize(residual_to_opt):
5✔
2012
        return o_smb + residual_to_opt - cmb
5✔
2013

2014
    residual = optimize.brentq(to_minimize, -1e5, 1e5)
5✔
2015

2016
    # Reset flux
2017
    for fl in fls:
5✔
2018
        fl.reset_flux()
5✔
2019

2020
    # Flowlines in order to be sure
2021
    rho = cfg.PARAMS['ice_density']
5✔
2022
    for fl_id, fl in enumerate(fls):
5✔
2023
        mbz = 0
5✔
2024
        for yr in mb_years:
5✔
2025
            mbz += mb_model.get_annual_mb(fl.surface_h, year=yr,
5✔
2026
                                          fls=fls, fl_id=fl_id)
2027
        mbz = mbz / len(mb_years)
5✔
2028
        fl.set_apparent_mb(mbz * cfg.SEC_IN_YEAR * rho + residual,
5✔
2029
                           is_calving=is_calving)
2030
        if fl_id < len(fls) and fl.flux_out < -1e3:
5✔
2031
            log.warning('({}) a tributary has a strongly negative flux. '
3✔
2032
                        'Inversion works but is physically quite '
2033
                        'questionable.'.format(gdir.rgi_id))
2034

2035
    # Check and write
2036
    _check_terminus_mass_flux(gdir, fls)
5✔
2037
    gdir.add_to_diagnostics('apparent_mb_from_any_mb_residual', residual)
5✔
2038
    gdir.write_pickle(fls, 'inversion_flowlines')
5✔
2039

2040

2041
@entity_task(log)
9✔
2042
def fixed_geometry_mass_balance(gdir, ys=None, ye=None, years=None,
9✔
2043
                                monthly_step=False,
2044
                                use_inversion_flowlines=True,
2045
                                climate_filename='climate_historical',
2046
                                climate_input_filesuffix='',
2047
                                temperature_bias=None,
2048
                                precipitation_factor=None,
2049
                                mb_model_class=MonthlyTIModel):
2050
    """Computes the mass balance with climate input from e.g. CRU or a GCM.
2051

2052
    Parameters
2053
    ----------
2054
    gdir : :py:class:`oggm.GlacierDirectory`
2055
        the glacier directory to process
2056
    ys : int
2057
        start year of the model run (default: from the climate file)
2058
        date)
2059
    ye : int
2060
        end year of the model run (default: from the climate file)
2061
    years : array of ints
2062
        override ys and ye with the years of your choice
2063
    monthly_step : bool
2064
        whether to store the diagnostic data at a monthly time step or not
2065
        (default is yearly)
2066
    use_inversion_flowlines : bool
2067
        whether to use the inversion flowlines or the model flowlines
2068
    climate_filename : str
2069
        name of the climate file, e.g. 'climate_historical' (default) or
2070
        'gcm_data'
2071
    climate_input_filesuffix: str
2072
        filesuffix for the input climate file
2073
    temperature_bias : float
2074
        add a bias to the temperature timeseries
2075
    precipitation_factor: float
2076
        multiply a factor to the precipitation time series
2077
        default is None and means that the precipitation factor from the
2078
        calibration is applied which is cfg.PARAMS['prcp_fac']
2079
    mb_model_class : MassBalanceModel class
2080
        the MassBalanceModel class to use, default is MonthlyTIModel
2081
    """
2082

2083
    if monthly_step:
3✔
2084
        raise NotImplementedError('monthly_step not implemented yet')
2085

2086
    mbmod = MultipleFlowlineMassBalance(gdir, mb_model_class=mb_model_class,
3✔
2087
                                        filename=climate_filename,
2088
                                        use_inversion_flowlines=use_inversion_flowlines,
2089
                                        input_filesuffix=climate_input_filesuffix)
2090

2091
    if temperature_bias is not None:
3✔
2092
        mbmod.temp_bias = temperature_bias
×
2093
    if precipitation_factor is not None:
3✔
2094
        mbmod.prcp_fac = precipitation_factor
×
2095

2096
    if years is None:
3✔
2097
        if ys is None:
3✔
2098
            ys = mbmod.flowline_mb_models[0].ys
3✔
2099
        if ye is None:
3✔
2100
            ye = mbmod.flowline_mb_models[0].ye
3✔
2101
        years = np.arange(ys, ye + 1)
3✔
2102

2103
    odf = pd.Series(data=mbmod.get_specific_mb(year=years),
3✔
2104
                    index=years)
2105
    return odf
3✔
2106

2107

2108
@entity_task(log)
9✔
2109
def compute_ela(gdir, ys=None, ye=None, years=None, climate_filename='climate_historical',
9✔
2110
                temperature_bias=None, precipitation_factor=None, climate_input_filesuffix='',
2111
                mb_model_class=MonthlyTIModel):
2112

2113
    """Computes the ELA of a glacier for a for given years and climate.
2114

2115
    Parameters
2116
    ----------
2117
    gdir : :py:class:`oggm.GlacierDirectory`
2118
        the glacier directory to process
2119
    ys : int
2120
        start year
2121
    ye : int
2122
        end year
2123
    years : array of ints
2124
        override ys and ye with the years of your choice
2125
    climate_filename : str
2126
        name of the climate file, e.g. 'climate_historical' (default) or
2127
        'gcm_data'
2128
    climate_input_filesuffix : str
2129
        filesuffix for the input climate file
2130
    temperature_bias : float
2131
        add a bias to the temperature timeseries
2132
    precipitation_factor: float
2133
        multiply a factor to the precipitation time series
2134
        default is None and means that the precipitation factor from the
2135
        calibration is applied which is cfg.PARAMS['prcp_fac']
2136
    mb_model_class : MassBalanceModel class
2137
        the MassBalanceModel class to use, default is MonthlyTIModel
2138
    """
2139

2140
    mbmod = mb_model_class(gdir, filename=climate_filename,
1✔
2141
                           input_filesuffix=climate_input_filesuffix)
2142

2143
    if temperature_bias is not None:
1✔
2144
        mbmod.temp_bias = temperature_bias
1✔
2145
    if precipitation_factor is not None:
1✔
2146
        mbmod.prcp_fac = precipitation_factor
×
2147

2148
    mbmod.valid_bounds = [-10000, 20000]
1✔
2149

2150
    if years is None:
1✔
2151
        years = np.arange(ys, ye+1)
1✔
2152

2153
    ela = []
1✔
2154
    for yr in years:
1✔
2155
        ela = np.append(ela, mbmod.get_ela(year=yr))
1✔
2156

2157
    odf = pd.Series(data=ela, index=years)
1✔
2158
    return odf
1✔
STATUS · Troubleshooting · Open an Issue · Sales · Support · CAREERS · ENTERPRISE · START FREE · SCHEDULE DEMO
ANNOUNCEMENTS · TWITTER · TOS & SLA · Supported CI Services · What's a CI service? · Automated Testing

© 2025 Coveralls, Inc