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

OGGM / oggm / 17036427554

18 Aug 2025 09:19AM UTC coverage: 84.387% (-0.03%) from 84.412%
17036427554

Pull #1795

github

web-flow
Merge f58c0836a into 0793d4813
Pull Request #1795: add new observation handling

208 of 217 new or added lines in 7 files covered. (95.85%)

28 existing lines in 2 files now uncovered.

12404 of 14699 relevant lines covered (84.39%)

3.78 hits per line

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

83.9
/oggm/graphics.py
1
"""Useful plotting functions"""
2
import os
9✔
3
import functools
9✔
4
import logging
9✔
5
from collections import OrderedDict
9✔
6
import itertools
9✔
7
import textwrap
9✔
8
import xarray as xr
9✔
9

10
import matplotlib
9✔
11
from matplotlib import colors
9✔
12
import matplotlib.pyplot as plt
9✔
13
import numpy as np
9✔
14
import shapely.geometry as shpg
9✔
15

16
try:
9✔
17
    import salem
9✔
18
except ImportError:
19
    pass
20

21
OGGM_CMAPS = dict()
9✔
22

23
from oggm.core.flowline import FileModel
9✔
24
from oggm import cfg, utils
9✔
25
from oggm.core import gis
9✔
26

27
# Module logger
28
log = logging.getLogger(__name__)
9✔
29

30

31
def set_oggm_cmaps():
9✔
32
    # Set global colormaps
33
    global OGGM_CMAPS
34
    OGGM_CMAPS['terrain'] = matplotlib.colormaps['terrain']
9✔
35
    OGGM_CMAPS['section_thickness'] = matplotlib.colormaps['YlGnBu']
9✔
36
    OGGM_CMAPS['glacier_thickness'] = matplotlib.colormaps['viridis']
9✔
37
    OGGM_CMAPS['ice_velocity'] = matplotlib.colormaps['Reds']
9✔
38

39

40
set_oggm_cmaps()
9✔
41

42

43
def truncate_colormap(cmap, minval=0.0, maxval=1.0, n=256):
9✔
44
    """Remove extreme colors from a colormap."""
45
    new_cmap = colors.LinearSegmentedColormap.from_list(
1✔
46
        'trunc({n},{a:.2f},{b:.2f})'.format(n=cmap.name, a=minval, b=maxval),
47
        cmap(np.linspace(minval, maxval, n)))
48
    return new_cmap
1✔
49

50

51
def gencolor_generator(n, cmap='Set1'):
9✔
52
    """ Color generator intended to work with qualitative color scales."""
53
    # don't use more than 9 discrete colors
54
    n_colors = min(n, 9)
1✔
55
    cmap = matplotlib.colormaps[cmap]
1✔
56
    colors = cmap(range(n_colors))
1✔
57
    for i in range(n):
1✔
58
        yield colors[i % n_colors]
1✔
59

60

61
def gencolor(n, cmap='Set1'):
9✔
62

63
    if isinstance(cmap, str):
1✔
64
        return gencolor_generator(n, cmap=cmap)
1✔
65
    else:
66
        return itertools.cycle(cmap)
×
67

68

69
def surf_to_nan(surf_h, thick):
9✔
70

71
    t1 = thick[:-2]
1✔
72
    t2 = thick[1:-1]
1✔
73
    t3 = thick[2:]
1✔
74
    pnan = ((t1 == 0) & (t2 == 0)) & ((t2 == 0) & (t3 == 0))
1✔
75
    surf_h[np.where(pnan)[0] + 1] = np.nan
1✔
76
    return surf_h
1✔
77

78

79
def _plot_map(plotfunc):
9✔
80
    """
81
    Decorator for common salem.Map plotting logic
82
    """
83
    commondoc = """
9✔
84

85
    Parameters
86
    ----------
87
    gdirs : [] or GlacierDirectory, required
88
        A single GlacierDirectory or a list of gdirs to plot.
89
    ax : matplotlib axes object, optional
90
        If None, uses own axis
91
    smap : Salem Map object, optional
92
        If None, makes a map from the first gdir in the list
93
    add_scalebar : Boolean, optional, default=True
94
        Adds scale bar to the plot
95
    add_colorbar : Boolean, optional, default=True
96
        Adds colorbar to axis
97
    horizontal_colorbar : Boolean, optional, default=False
98
        Horizontal colorbar instead
99
    title : str, optional
100
        If left to None, the plot decides whether it writes a title or not. Set
101
        to '' for no title.
102
    title_comment : str, optional
103
        add something to the default title. Set to none to remove default
104
    lonlat_contours_kwargs: dict, optional
105
        pass kwargs to salem.Map.set_lonlat_contours
106
    cbar_ax: ax, optional
107
        ax where to plot the colorbar
108
    autosave : bool, optional
109
        set to True to override to a default savefig filename (useful
110
        for multiprocessing)
111
    figsize : tuple, optional
112
        size of the figure
113
    savefig : str, optional
114
        save the figure to a file instead of displaying it
115
    savefig_kwargs : dict, optional
116
        the kwargs to plt.savefig
117
    extend_plot_limit : bool, optional
118
        set to True to extend the plotting limits for all provided gdirs grids
119
    """
120

121
    # Build on the original docstring
122
    plotfunc.__doc__ = '\n'.join((plotfunc.__doc__, commondoc))
9✔
123

124
    @functools.wraps(plotfunc)
9✔
125
    def newplotfunc(gdirs, ax=None, smap=None, add_colorbar=True, title=None,
9✔
126
                    title_comment=None, horizontal_colorbar=False,
127
                    lonlat_contours_kwargs=None, cbar_ax=None, autosave=False,
128
                    add_scalebar=True, figsize=None, savefig=None,
129
                    savefig_kwargs=None, extend_plot_limit=False,
130
                    **kwargs):
131

132
        dofig = False
3✔
133
        if ax is None:
3✔
134
            fig = plt.figure(figsize=figsize)
×
135
            ax = fig.add_subplot(111)
×
136
            dofig = True
×
137

138
        # Cast to list
139
        gdirs = utils.tolist(gdirs)
3✔
140

141
        if smap is None:
3✔
142
            if extend_plot_limit:
3✔
143
                grid_combined = utils.combine_grids(gdirs)
2✔
144
                mp = salem.Map(grid_combined, countries=False,
2✔
145
                               nx=grid_combined.nx)
146
            else:
147
                mp = salem.Map(gdirs[0].grid, countries=False,
1✔
148
                               nx=gdirs[0].grid.nx)
149
        else:
150
            mp = smap
3✔
151

152
        if lonlat_contours_kwargs is not None:
3✔
153
            mp.set_lonlat_contours(**lonlat_contours_kwargs)
×
154

155
        if add_scalebar:
3✔
156
            mp.set_scale_bar()
3✔
157
        out = plotfunc(gdirs, ax=ax, smap=mp, **kwargs)
3✔
158

159
        if add_colorbar and 'cbar_label' in out:
3✔
160
            cbprim = out.get('cbar_primitive', mp)
3✔
161
            if cbar_ax:
3✔
162
                cb = cbprim.colorbarbase(cbar_ax)
×
163
            else:
164
                if horizontal_colorbar:
3✔
165
                    cb = cbprim.append_colorbar(ax, "bottom", size="5%",
×
166
                                                pad=0.4)
167
                else:
168
                    cb = cbprim.append_colorbar(ax, "right", size="5%",
3✔
169
                                                pad=0.2)
170
            cb.set_label(out['cbar_label'])
3✔
171

172
        if title is None:
3✔
173
            if 'title' not in out:
3✔
174
                # Make a default one
175
                title = ''
3✔
176
                if len(gdirs) == 1:
3✔
177
                    gdir = gdirs[0]
1✔
178
                    title = gdir.rgi_id
1✔
179
                    if gdir.name is not None and gdir.name != '':
1✔
180
                        title += ': ' + gdir.name
1✔
181
                out['title'] = title
3✔
182

183
            if title_comment is None:
3✔
184
                title_comment = out.get('title_comment', '')
3✔
185

186
            out['title'] += title_comment
3✔
187
            ax.set_title(out['title'])
3✔
188
        else:
189
            ax.set_title(title)
×
190

191
        if dofig:
3✔
192
            plt.tight_layout()
×
193

194
        if autosave:
3✔
195
            savefig = os.path.join(cfg.PATHS['working_dir'], 'plots')
×
196
            utils.mkdir(savefig)
×
197
            savefig = os.path.join(savefig, plotfunc.__name__ + '_' +
×
198
                                   gdirs[0].rgi_id + '.png')
199

200
        if savefig is not None:
3✔
201
            plt.savefig(savefig, **savefig_kwargs)
×
202
            plt.close()
×
203

204
    return newplotfunc
9✔
205

206

207
def plot_googlemap(gdirs, ax=None, figsize=None, key=None):
9✔
208
    """Plots the glacier(s) over a googlemap."""
209

UNCOV
210
    if key is None:
×
UNCOV
211
        try:
×
UNCOV
212
            key = os.environ['STATIC_MAP_API_KEY']
×
213
        except KeyError:
×
214
            raise ValueError('You need to provide a Google API key'
×
215
                             ' or set the STATIC_MAP_API_KEY environment'
216
                             ' variable.')
217

UNCOV
218
    dofig = False
×
UNCOV
219
    if ax is None:
×
220
        fig = plt.figure(figsize=figsize)
×
221
        ax = fig.add_subplot(111)
×
222
        dofig = True
×
223

UNCOV
224
    gdirs = utils.tolist(gdirs)
×
225

UNCOV
226
    xx, yy = [], []
×
UNCOV
227
    for gdir in gdirs:
×
UNCOV
228
        xx.extend(gdir.extent_ll[0])
×
UNCOV
229
        yy.extend(gdir.extent_ll[1])
×
230

UNCOV
231
    gm = salem.GoogleVisibleMap(xx, yy, key=key, use_cache=False)
×
232

UNCOV
233
    img = gm.get_vardata()
×
UNCOV
234
    cmap = salem.Map(gm.grid, countries=False, nx=gm.grid.nx)
×
UNCOV
235
    cmap.set_rgb(img)
×
236

UNCOV
237
    for gdir in gdirs:
×
UNCOV
238
        cmap.set_shapefile(gdir.read_shapefile('outlines'))
×
239

UNCOV
240
    cmap.plot(ax)
×
UNCOV
241
    title = ''
×
UNCOV
242
    if len(gdirs) == 1:
×
UNCOV
243
        title = gdir.rgi_id
×
UNCOV
244
        if gdir.name is not None and gdir.name != '':
×
UNCOV
245
            title += ': ' + gdir.name
×
UNCOV
246
    ax.set_title(title)
×
247

UNCOV
248
    if dofig:
×
249
        plt.tight_layout()
×
250

251

252
@_plot_map
9✔
253
def plot_raster(gdirs, var_name=None, cmap='viridis', ax=None, smap=None):
9✔
254
    """Plot any raster from the gridded_data file."""
255

256
    # Files
257
    gdir = gdirs[0]
1✔
258

259
    with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
1✔
260
        var = nc.variables[var_name]
1✔
261
        data = var[:]
1✔
262
        description = var.long_name
1✔
263
        description += ' [{}]'.format(var.units)
1✔
264

265
    smap.set_data(data)
1✔
266

267
    smap.set_cmap(cmap)
1✔
268

269
    for gdir in gdirs:
1✔
270
        crs = gdir.grid.center_grid
1✔
271

272
        try:
1✔
273
            geom = gdir.read_pickle('geometries')
1✔
274
            # Plot boundaries
275
            poly_pix = geom['polygon_pix']
1✔
276
            smap.set_geometry(poly_pix, crs=crs, fc='none',
1✔
277
                              alpha=0.3, zorder=2, linewidth=.2)
278
            poly_pix = utils.tolist(poly_pix)
1✔
279
            for _poly in poly_pix:
1✔
280
                for l in _poly.interiors:
1✔
281
                    smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
1✔
282
        except FileNotFoundError:
×
283
            smap.set_shapefile(gdir.read_shapefile('outlines'))
×
284

285
    smap.plot(ax)
1✔
286

287
    return dict(cbar_label='\n'.join(textwrap.wrap(description, 30)))
1✔
288

289

290
@_plot_map
9✔
291
def plot_domain(gdirs, ax=None, smap=None, use_netcdf=False):
9✔
292
    """Plot the glacier directory.
293

294
    Parameters
295
    ----------
296
    gdirs
297
    ax
298
    smap
299
    use_netcdf : bool
300
        use output of glacier_masks instead of geotiff DEM
301
    """
302

303
    # Files
304
    gdir = gdirs[0]
1✔
305
    if use_netcdf:
1✔
306
        with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
×
307
            topo = nc.variables['topo'][:]
×
308
    else:
309
        topo = gis.read_geotiff_dem(gdir)
1✔
310
    try:
1✔
311
        smap.set_data(topo)
1✔
312
    except ValueError:
×
313
        pass
×
314

315
    cm = truncate_colormap(OGGM_CMAPS['terrain'], minval=0.25, maxval=1.0)
1✔
316
    smap.set_plot_params(cmap=cm)
1✔
317

318
    for gdir in gdirs:
1✔
319
        crs = gdir.grid.center_grid
1✔
320

321
        try:
1✔
322
            geom = gdir.read_pickle('geometries')
1✔
323

324
            # Plot boundaries
325
            poly_pix = geom['polygon_pix']
1✔
326
            smap.set_geometry(poly_pix, crs=crs, fc='white',
1✔
327
                              alpha=0.3, zorder=2, linewidth=.2)
328
            poly_pix = utils.tolist(poly_pix)
1✔
329
            for _poly in poly_pix:
1✔
330
                for l in _poly.interiors:
1✔
331
                    smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
1✔
332
        except FileNotFoundError:
×
333
            smap.set_shapefile(gdir.read_shapefile('outlines'))
×
334

335
    smap.plot(ax)
1✔
336

337
    return dict(cbar_label='Alt. [m]')
1✔
338

339

340
@_plot_map
9✔
341
def plot_centerlines(gdirs, ax=None, smap=None, use_flowlines=False,
9✔
342
                     add_downstream=False, lines_cmap='Set1',
343
                     add_line_index=False, use_model_flowlines=False):
344
    """Plots the centerlines of a glacier directory."""
345

346
    if add_downstream and not use_flowlines:
1✔
347
        raise ValueError('Downstream lines can be plotted with flowlines only')
×
348

349
    # Files
350
    filename = 'centerlines'
1✔
351
    if use_model_flowlines:
1✔
352
        filename = 'model_flowlines'
×
353
    elif use_flowlines:
1✔
354
        filename = 'inversion_flowlines'
1✔
355

356
    gdir = gdirs[0]
1✔
357
    with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
1✔
358
        topo = nc.variables['topo'][:]
1✔
359

360
    cm = truncate_colormap(OGGM_CMAPS['terrain'], minval=0.25, maxval=1.0)
1✔
361
    smap.set_plot_params(cmap=cm)
1✔
362
    smap.set_data(topo)
1✔
363
    for gdir in gdirs:
1✔
364
        crs = gdir.grid.center_grid
1✔
365
        geom = gdir.read_pickle('geometries')
1✔
366

367
        # Plot boundaries
368
        poly_pix = geom['polygon_pix']
1✔
369

370
        smap.set_geometry(poly_pix, crs=crs, fc='white',
1✔
371
                          alpha=0.3, zorder=2, linewidth=.2)
372
        poly_pix = utils.tolist(poly_pix)
1✔
373
        for _poly in poly_pix:
1✔
374
            for l in _poly.interiors:
1✔
375
                smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
1✔
376

377
        # plot Centerlines
378
        cls = gdir.read_pickle(filename)
1✔
379

380
        # Go in reverse order for red always being the longest
381
        cls = cls[::-1]
1✔
382
        nl = len(cls)
1✔
383
        color = gencolor(len(cls) + 1, cmap=lines_cmap)
1✔
384
        for i, (l, c) in enumerate(zip(cls, color)):
1✔
385
            if add_downstream and not gdir.is_tidewater and l is cls[0]:
1✔
386
                line = gdir.read_pickle('downstream_line')['full_line']
1✔
387
            else:
388
                line = l.line
1✔
389

390
            smap.set_geometry(line, crs=crs, color=c,
1✔
391
                              linewidth=2.5, zorder=50)
392

393
            text = '{}'.format(nl - i - 1) if add_line_index else None
1✔
394
            smap.set_geometry(l.head, crs=gdir.grid, marker='o',
1✔
395
                              markersize=60, alpha=0.8, color=c, zorder=99,
396
                              text=text)
397

398
            for j in l.inflow_points:
1✔
399
                smap.set_geometry(j, crs=crs, marker='o',
1✔
400
                                  markersize=40, edgecolor='k', alpha=0.8,
401
                                  zorder=99, facecolor='none')
402

403
    smap.plot(ax)
1✔
404
    return dict(cbar_label='Alt. [m]')
1✔
405

406

407
@_plot_map
9✔
408
def plot_catchment_areas(gdirs, ax=None, smap=None, lines_cmap='Set1',
9✔
409
                         mask_cmap='Set2'):
410
    """Plots the catchments out of a glacier directory.
411
    """
412

413
    gdir = gdirs[0]
1✔
414
    if len(gdirs) > 1:
1✔
415
        raise NotImplementedError('Cannot plot a list of gdirs (yet)')
416

417
    with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
1✔
418
        topo = nc.variables['topo'][:]
1✔
419
        mask = nc.variables['glacier_mask'][:] * np.nan
1✔
420

421
    smap.set_topography(topo)
1✔
422

423
    crs = gdir.grid.center_grid
1✔
424
    geom = gdir.read_pickle('geometries')
1✔
425

426
    # Plot boundaries
427
    poly_pix = geom['polygon_pix']
1✔
428
    smap.set_geometry(poly_pix, crs=crs, fc='none', zorder=2,
1✔
429
                      linewidth=.2)
430
    for l in poly_pix.interiors:
1✔
431
        smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
1✔
432

433
    # plot Centerlines
434
    cls = gdir.read_pickle('centerlines')[::-1]
1✔
435
    color = gencolor(len(cls) + 1, cmap=lines_cmap)
1✔
436
    for l, c in zip(cls, color):
1✔
437
        smap.set_geometry(l.line, crs=crs, color=c,
1✔
438
                          linewidth=2.5, zorder=50)
439

440
    # catchment areas
441
    cis = gdir.read_pickle('geometries')['catchment_indices']
1✔
442
    for j, ci in enumerate(cis[::-1]):
1✔
443
        mask[tuple(ci.T)] = j+1
1✔
444

445
    smap.set_cmap(mask_cmap)
1✔
446
    smap.set_data(mask)
1✔
447
    smap.plot(ax)
1✔
448

449
    return {}
1✔
450

451

452
@_plot_map
9✔
453
def plot_catchment_width(gdirs, ax=None, smap=None, corrected=False,
9✔
454
                         add_intersects=False, add_touches=False,
455
                         lines_cmap='Set1'):
456
    """Plots the catchment widths out of a glacier directory.
457
    """
458

459
    gdir = gdirs[0]
1✔
460
    with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
1✔
461
        topo = nc.variables['topo'][:]
1✔
462
    # Dirty optim
463
    try:
1✔
464
        smap.set_topography(topo)
1✔
465
    except ValueError:
1✔
466
        pass
1✔
467

468
    # Maybe plot touches
469
    xis, yis, cis = [], [], []
1✔
470
    ogrid = smap.grid
1✔
471

472
    for gdir in gdirs:
1✔
473
        crs = gdir.grid.center_grid
1✔
474
        geom = gdir.read_pickle('geometries')
1✔
475

476
        # Plot boundaries
477
        poly_pix = geom['polygon_pix']
1✔
478
        smap.set_geometry(poly_pix, crs=crs, fc='none', zorder=2,
1✔
479
                          linewidth=.2)
480
        for l in poly_pix.interiors:
1✔
481
            smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
1✔
482

483
        # Plot intersects
484
        if add_intersects and gdir.has_file('intersects'):
1✔
485
            gdf = gdir.read_shapefile('intersects')
1✔
486
            smap.set_shapefile(gdf, color='k', linewidth=3.5, zorder=3)
1✔
487

488
        # plot Centerlines
489
        cls = gdir.read_pickle('inversion_flowlines')[::-1]
1✔
490
        color = gencolor(len(cls) + 1, cmap=lines_cmap)
1✔
491
        for l, c in zip(cls, color):
1✔
492
            smap.set_geometry(l.line, crs=crs, color=c,
1✔
493
                              linewidth=2.5, zorder=50)
494
            if corrected:
1✔
495
                for wi, cur, (n1, n2) in zip(l.widths, l.line.coords,
1✔
496
                                             l.normals):
497
                    _l = shpg.LineString([shpg.Point(cur + wi / 2. * n1),
1✔
498
                                          shpg.Point(cur + wi / 2. * n2)])
499

500
                    smap.set_geometry(_l, crs=crs, color=c,
1✔
501
                                      linewidth=0.6, zorder=50)
502
            else:
503
                for wl, wi in zip(l.geometrical_widths, l.widths):
1✔
504
                    col = c if np.isfinite(wi) else 'grey'
1✔
505
                    for w in wl.geoms:
1✔
506
                        smap.set_geometry(w, crs=crs, color=col,
1✔
507
                                          linewidth=0.6, zorder=50)
508

509
            if add_touches:
1✔
510
                pok = np.where(l.is_rectangular)
1✔
511
                if np.size(pok[0]) != 0:
1✔
512
                    xi, yi = l.line.xy
1✔
513
                    xi, yi = ogrid.transform(np.asarray(xi)[pok],
1✔
514
                                             np.asarray(yi)[pok], crs=crs)
515
                    xis.append(xi)
1✔
516
                    yis.append(yi)
1✔
517
                    cis.append(c)
1✔
518

519
    smap.plot(ax)
1✔
520
    for xi, yi, c in zip(xis, yis, cis):
1✔
521
        ax.scatter(xi, yi, color=c, s=20, zorder=51)
1✔
522

523
    return {}
1✔
524

525

526
@_plot_map
9✔
527
def plot_inversion(gdirs, ax=None, smap=None, linewidth=3, vmax=None,
9✔
528
                   plot_var='thick', cbar_label='Section thickness (m)',
529
                   color_map='YlGnBu'):
530
    """Plots the result of the inversion out of a glacier directory.
531
       Default is thickness (m). Change plot_var to u_surface or u_integrated
532
       for velocity (m/yr)."""
533

534
    gdir = gdirs[0]
3✔
535
    with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
3✔
536
        topo = nc.variables['topo'][:]
3✔
537

538
    # Dirty optim
539
    try:
3✔
540
        smap.set_topography(topo)
3✔
541
    except ValueError:
2✔
542
        pass
2✔
543

544
    toplot_var = np.array([])
3✔
545
    toplot_lines = []
3✔
546
    toplot_crs = []
3✔
547
    vol = []
3✔
548
    for gdir in gdirs:
3✔
549
        crs = gdir.grid.center_grid
3✔
550
        geom = gdir.read_pickle('geometries')
3✔
551
        inv = gdir.read_pickle('inversion_output')
3✔
552
        # Plot boundaries
553
        poly_pix = geom['polygon_pix']
3✔
554
        smap.set_geometry(poly_pix, crs=crs, fc='none', zorder=2,
3✔
555
                          linewidth=.2)
556
        for l in poly_pix.interiors:
3✔
557
            smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
3✔
558

559
        # Plot Centerlines
560
        cls = gdir.read_pickle('inversion_flowlines')
3✔
561
        for l, c in zip(cls, inv):
3✔
562

563
            smap.set_geometry(l.line, crs=crs, color='gray',
3✔
564
                              linewidth=1.2, zorder=50)
565

566
            toplot_var = np.append(toplot_var, c[plot_var])
3✔
567
            for wi, cur, (n1, n2) in zip(l.widths, l.line.coords, l.normals):
3✔
568
                line = shpg.LineString([shpg.Point(cur + wi / 2. * n1),
3✔
569
                                        shpg.Point(cur + wi / 2. * n2)])
570
                toplot_lines.append(line)
3✔
571
                toplot_crs.append(crs)
3✔
572
            vol.extend(c['volume'])
3✔
573

574
    dl = salem.DataLevels(cmap=matplotlib.colormaps[color_map],
3✔
575
                          data=toplot_var, vmin=0, vmax=vmax)
576
    colors = dl.to_rgb()
3✔
577
    for l, c, crs in zip(toplot_lines, colors, toplot_crs):
3✔
578
        smap.set_geometry(l, crs=crs, color=c,
3✔
579
                          linewidth=linewidth, zorder=50)
580

581
    smap.plot(ax)
3✔
582
    out = dict(cbar_label=cbar_label,
3✔
583
                cbar_primitive=dl)
584

585
    if plot_var == 'thick':
3✔
586
        out['title_comment'] = ' ({:.2f} km3)'.format(np.nansum(vol) * 1e-9)
3✔
587

588
    return out
3✔
589

590

591
@_plot_map
9✔
592
def plot_distributed_thickness(gdirs, ax=None, smap=None, varname_suffix=''):
9✔
593
    """Plots the result of the inversion out of a glacier directory.
594

595
    Method: 'alt' or 'interp'
596
    """
597

598
    gdir = gdirs[0]
1✔
599

600
    with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
1✔
601
        topo = nc.variables['topo'][:]
1✔
602

603
    try:
1✔
604
        smap.set_topography(topo)
1✔
605
    except ValueError:
×
606
        pass
×
607

608
    for gdir in gdirs:
1✔
609
        grids_file = gdir.get_filepath('gridded_data')
1✔
610
        with utils.ncDataset(grids_file) as nc:
1✔
611
            import warnings
1✔
612
            with warnings.catch_warnings():
1✔
613
                # https://github.com/Unidata/netcdf4-python/issues/766
614
                warnings.filterwarnings("ignore", category=RuntimeWarning)
1✔
615
                vn = 'distributed_thickness' + varname_suffix
1✔
616
                thick = nc.variables[vn][:]
1✔
617
                mask = nc.variables['glacier_mask'][:]
1✔
618

619
        thick = np.where(mask, thick, np.nan)
1✔
620

621
        crs = gdir.grid.center_grid
1✔
622

623
        # Plot boundaries
624
        # Try to read geometries.pkl as the glacier boundary,
625
        # if it can't be found, we use the shapefile to instead.
626
        try:
1✔
627
            geom = gdir.read_pickle('geometries')
1✔
628
            poly_pix = geom['polygon_pix']
1✔
629
            smap.set_geometry(poly_pix, crs=crs, fc='none', zorder=2, linewidth=.2)
1✔
630
            for l in poly_pix.interiors:
1✔
631
                smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
1✔
632
        except FileNotFoundError:
1✔
633
            smap.set_shapefile(gdir.read_shapefile('outlines'), fc='none')
1✔
634
        smap.set_data(thick, crs=crs, overplot=True)
1✔
635

636
    smap.set_plot_params(cmap=OGGM_CMAPS['glacier_thickness'])
1✔
637
    smap.plot(ax)
1✔
638

639
    return dict(cbar_label='Glacier thickness [m]')
1✔
640

641

642
@_plot_map
9✔
643
def plot_modeloutput_map(gdirs, ax=None, smap=None, model=None,
9✔
644
                         vmax=None, linewidth=3, filesuffix='',
645
                         modelyr=None, plotting_var='thickness'):
646
    """Plots the result of the model output.
647

648
    Parameters
649
    ----------
650
    gdirs
651
    ax
652
    smap
653
    model
654
    vmax
655
    linewidth
656
    filesuffix
657
    modelyr
658
    plotting_var : str
659
        Defines which variable should be plotted. Options are 'thickness'
660
        (default) and 'velocity'. If you want to plot velocity the flowline
661
        diagnostics of the run are needed (set
662
        cfg.PARAMS['store_fl_diagnostics'] = True, before the
663
        actual simulation) and be aware that there is no velocity available for
664
        the first year of the simulation.
665

666
    Returns
667
    -------
668

669
    """
670

671
    gdir = gdirs[0]
3✔
672
    with utils.ncDataset(gdir.get_filepath('gridded_data')) as nc:
3✔
673
        topo = nc.variables['topo'][:]
3✔
674

675
    # Dirty optim
676
    try:
3✔
677
        smap.set_topography(topo)
3✔
678
    except ValueError:
2✔
679
        pass
2✔
680

681
    toplot_var = np.array([])
3✔
682
    toplot_lines = []
3✔
683
    toplot_crs = []
3✔
684

685
    if model is None:
3✔
686
        models = []
2✔
687
        for gdir in gdirs:
2✔
688
            model = FileModel(gdir.get_filepath('model_geometry',
2✔
689
                                                filesuffix=filesuffix))
690
            model.run_until(modelyr)
2✔
691
            models.append(model)
2✔
692
    else:
693
        models = utils.tolist(model)
1✔
694

695
    if modelyr is None:
3✔
696
        modelyr = models[0].yr
1✔
697

698
    for gdir, model in zip(gdirs, models):
3✔
699
        geom = gdir.read_pickle('geometries')
3✔
700
        poly_pix = geom['polygon_pix']
3✔
701

702
        crs = gdir.grid.center_grid
3✔
703
        smap.set_geometry(poly_pix, crs=crs, fc='none', zorder=2, linewidth=.2)
3✔
704

705
        poly_pix = utils.tolist(poly_pix)
3✔
706
        for _poly in poly_pix:
3✔
707
            for l in _poly.interiors:
3✔
708
                smap.set_geometry(l, crs=crs, color='black', linewidth=0.5)
3✔
709

710
        if plotting_var == 'velocity':
3✔
711
            f_fl_diag = gdir.get_filepath('fl_diagnostics',
×
712
                                          filesuffix=filesuffix)
713

714
        # plot Centerlines
715
        cls = model.fls
3✔
716
        for fl_id, l in enumerate(cls):
3✔
717
            smap.set_geometry(l.line, crs=crs, color='gray',
3✔
718
                              linewidth=1.2, zorder=50)
719
            if plotting_var == 'thickness':
3✔
720
                toplot_var = np.append(toplot_var, l.thick)
3✔
721
            elif plotting_var == 'velocity':
×
722
                with xr.open_dataset(f_fl_diag, group=f'fl_{fl_id}') as ds:
×
723
                    toplot_var = np.append(toplot_var,
×
724
                                           ds.sel(dict(time=modelyr)).ice_velocity_myr)
725
            widths = l.widths.copy()
3✔
726
            widths = np.where(l.thick > 0, widths, 0.)
3✔
727
            for wi, cur, (n1, n2) in zip(widths, l.line.coords, l.normals):
3✔
728
                line = shpg.LineString([shpg.Point(cur + wi/2. * n1),
3✔
729
                                        shpg.Point(cur + wi/2. * n2)])
730
                toplot_lines.append(line)
3✔
731
                toplot_crs.append(crs)
3✔
732

733
    if plotting_var == 'thickness':
3✔
734
        cmap = OGGM_CMAPS['section_thickness']
3✔
735
        cbar_label = 'Section thickness [m]'
3✔
736
    elif plotting_var == 'velocity':
×
737
        cmap = OGGM_CMAPS['ice_velocity']
×
738
        cbar_label = 'Ice velocity [m yr-1]'
×
739
    dl = salem.DataLevels(cmap=cmap,
3✔
740
                          data=toplot_var, vmin=0, vmax=vmax)
741
    colors = dl.to_rgb()
3✔
742
    for l, c, crs in zip(toplot_lines, colors, toplot_crs):
3✔
743
        smap.set_geometry(l, crs=crs, color=c,
3✔
744
                          linewidth=linewidth, zorder=50)
745
    smap.plot(ax)
3✔
746
    return dict(cbar_label=cbar_label,
3✔
747
                cbar_primitive=dl,
748
                title_comment=' -- year: {:d}'.format(np.int64(model.yr)))
749

750

751
def plot_modeloutput_section(model=None, ax=None, title=''):
9✔
752
    """Plots the result of the model output along the flowline.
753

754
    Parameters
755
    ----------
756
    model: obj
757
        either a FlowlineModel or a list of model flowlines.
758
    fig
759
    title
760
    """
761

762
    try:
1✔
763
        fls = model.fls
1✔
764
    except AttributeError:
×
765
        fls = model
×
766

767
    if ax is None:
1✔
768
        fig = plt.figure(figsize=(12, 6))
×
769
        ax = fig.add_axes([0.07, 0.08, 0.7, 0.84])
×
770
    else:
771
        fig = plt.gcf()
1✔
772

773
    # Compute area histo
774
    area = np.array([])
1✔
775
    height = np.array([])
1✔
776
    bed = np.array([])
1✔
777
    for cls in fls:
1✔
778
        a = cls.widths_m * cls.dx_meter * 1e-6
1✔
779
        a = np.where(cls.thick > 0, a, 0)
1✔
780
        area = np.concatenate((area, a))
1✔
781
        height = np.concatenate((height, cls.surface_h))
1✔
782
        bed = np.concatenate((bed, cls.bed_h))
1✔
783
    ylim = [bed.min(), height.max()]
1✔
784

785
    # Plot histo
786
    posax = ax.get_position()
1✔
787
    posax = [posax.x0 + 2 * posax.width / 3.0,
1✔
788
             posax.y0,  posax.width / 3.0,
789
             posax.height]
790
    axh = fig.add_axes(posax, frameon=False)
1✔
791

792
    axh.hist(height, orientation='horizontal', range=ylim, bins=20,
1✔
793
             alpha=0.3, weights=area)
794
    axh.invert_xaxis()
1✔
795
    axh.xaxis.tick_top()
1✔
796
    axh.set_xlabel('Area incl. tributaries (km$^2$)')
1✔
797
    axh.xaxis.set_label_position('top')
1✔
798
    axh.set_ylim(ylim)
1✔
799
    axh.yaxis.set_ticks_position('right')
1✔
800
    axh.set_yticks([])
1✔
801
    axh.axhline(y=ylim[1], color='black', alpha=1)  # qick n dirty trick
1✔
802

803
    # plot Centerlines
804
    cls = fls[-1]
1✔
805
    x = np.arange(cls.nx) * cls.dx * cls.map_dx
1✔
806

807
    # Plot the bed
808
    ax.plot(x, cls.bed_h, color='k', linewidth=2.5, label='Bed (Parab.)')
1✔
809

810
    # Where trapezoid change color
811
    if hasattr(cls, '_do_trapeze') and cls._do_trapeze:
1✔
812
        bed_t = cls.bed_h * np.nan
1✔
813
        pt = cls.is_trapezoid & (~cls.is_rectangular)
1✔
814
        bed_t[pt] = cls.bed_h[pt]
1✔
815
        ax.plot(x, bed_t, color='rebeccapurple', linewidth=2.5,
1✔
816
                label='Bed (Trap.)')
817
        bed_t = cls.bed_h * np.nan
1✔
818
        bed_t[cls.is_rectangular] = cls.bed_h[cls.is_rectangular]
1✔
819
        ax.plot(x, bed_t, color='crimson', linewidth=2.5, label='Bed (Rect.)')
1✔
820

821
    # Plot glacier
822
    surfh = surf_to_nan(cls.surface_h, cls.thick)
1✔
823
    ax.plot(x, surfh, color='#003399', linewidth=2, label='Glacier')
1✔
824

825
    # Plot tributaries
826
    for i, l in zip(cls.inflow_indices, cls.inflows):
1✔
827
        if l.thick[-1] > 0:
1✔
828
            ax.plot(x[i], cls.surface_h[i], 's', markerfacecolor='#993399',
1✔
829
                    markeredgecolor='k',
830
                    label='Tributary (active)')
831
        else:
832
            ax.plot(x[i], cls.surface_h[i], 's', markerfacecolor='w',
×
833
                    markeredgecolor='k',
834
                    label='Tributary (inactive)')
835
    if getattr(model, 'do_calving', False):
1✔
836
        ax.hlines(model.water_level, x[0], x[-1], linestyles=':', color='C0')
×
837

838
    ax.set_ylim(ylim)
1✔
839

840
    ax.spines['top'].set_color('none')
1✔
841
    ax.xaxis.set_ticks_position('bottom')
1✔
842
    ax.set_xlabel('Distance along flowline (m)')
1✔
843
    ax.set_ylabel('Altitude (m)')
1✔
844

845
    # Title
846
    ax.set_title(title, loc='left')
1✔
847

848
    # Legend
849
    handles, labels = ax.get_legend_handles_labels()
1✔
850
    by_label = OrderedDict(zip(labels, handles))
1✔
851
    ax.legend(list(by_label.values()), list(by_label.keys()),
1✔
852
              bbox_to_anchor=(1.34, 1.0),
853
              frameon=False)
854

855

856
def plot_modeloutput_section_withtrib(model=None, fig=None, title=''):
9✔
857
    """Plots the result of the model output along the flowline.
858

859
    Parameters
860
    ----------
861
    model: obj
862
        either a FlowlineModel or a list of model flowlines.
863
    fig
864
    title
865
    """
866

867
    try:
1✔
868
        fls = model.fls
1✔
869
    except AttributeError:
×
870
        fls = model
×
871

872
    n_tribs = len(fls) - 1
1✔
873

874
    axs = []
1✔
875
    if n_tribs == 0:
1✔
876
        if fig is None:
×
877
            fig = plt.figure(figsize=(8, 5))
×
878
        axmaj = fig.add_subplot(111)
×
879
    elif n_tribs <= 3:
1✔
880
        if fig is None:
1✔
881
            fig = plt.figure(figsize=(14, 10))
×
882
        axmaj = plt.subplot2grid((2, 3), (1, 0), colspan=3)
1✔
883
        for i in np.arange(n_tribs):
1✔
884
            axs.append(plt.subplot2grid((2, 3), (0, i)))
1✔
885
    elif n_tribs <= 6:
×
886
        if fig is None:
×
887
            fig = plt.figure(figsize=(14, 10))
×
888
        axmaj = plt.subplot2grid((3, 3), (2, 0), colspan=3)
×
889
        for i in np.arange(n_tribs):
×
890
            j = 0
×
891
            if i >= 3:
×
892
                i -= 3
×
893
                j = 1
×
894
            axs.append(plt.subplot2grid((3, 3), (j, i)))
×
895
    else:
896
        raise NotImplementedError()
897

898
    for i, cls in enumerate(fls):
1✔
899
        if i == n_tribs:
1✔
900
            ax = axmaj
1✔
901
        else:
902
            ax = axs[i]
1✔
903

904
        x = np.arange(cls.nx) * cls.dx * cls.map_dx
1✔
905

906
        # Plot the bed
907
        ax.plot(x, cls.bed_h, color='k', linewidth=2.5, label='Bed (Parab.)')
1✔
908

909
        # Where trapezoid change color
910
        if hasattr(cls, '_do_trapeze') and cls._do_trapeze:
1✔
911
            bed_t = cls.bed_h * np.nan
1✔
912
            pt = cls.is_trapezoid & (~cls.is_rectangular)
1✔
913
            bed_t[pt] = cls.bed_h[pt]
1✔
914
            ax.plot(x, bed_t, color='rebeccapurple', linewidth=2.5,
1✔
915
                    label='Bed (Trap.)')
916
            bed_t = cls.bed_h * np.nan
1✔
917
            bed_t[cls.is_rectangular] = cls.bed_h[cls.is_rectangular]
1✔
918
            ax.plot(x, bed_t, color='crimson', linewidth=2.5,
1✔
919
                    label='Bed (Rect.)')
920

921
        # Plot glacier
922
        surfh = surf_to_nan(cls.surface_h, cls.thick)
1✔
923
        ax.plot(x, surfh, color='#003399', linewidth=2, label='Glacier')
1✔
924

925
        # Plot tributaries
926
        for i, l in zip(cls.inflow_indices, cls.inflows):
1✔
927
            if l.thick[-1] > 0:
1✔
928
                ax.plot(x[i], cls.surface_h[i], 's', color='#993399',
1✔
929
                        label='Tributary (active)')
930
            else:
931
                ax.plot(x[i], cls.surface_h[i], 's', color='none',
×
932
                        label='Tributary (inactive)')
933

934
        ax.spines['top'].set_color('none')
1✔
935
        ax.xaxis.set_ticks_position('bottom')
1✔
936
        ax.set_xlabel('Distance along flowline (m)')
1✔
937
        ax.set_ylabel('Altitude (m)')
1✔
938

939
    # Title
940
    plt.title(title, loc='left')
1✔
941

942
    # Legend
943
    handles, labels = ax.get_legend_handles_labels()
1✔
944
    by_label = OrderedDict(zip(labels, handles))
1✔
945
    ax.legend(list(by_label.values()), list(by_label.keys()),
1✔
946
              loc='best', frameon=False)
947
    fig.tight_layout()
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