source: trunk/LMDZ.COMMON/libf/evolution/deftank/visu_evol_layering.py @ 3819

Last change on this file since 3819 was 3819, checked in by jbclement, 7 days ago

PEM:
Addition of 2 figures in the outputs of "visu_evol_layering.py": (i) the stratification over time colored using RGB ternary mix of H2O ice, CO2 ice and dust; (ii) the dust-to-ice ratio in the stratification over time.
JBC

  • Property svn:executable set to *
File size: 30.5 KB
Line 
1#!/usr/bin/env python3
2#######################################################################################################
3### Python script to output stratification data over time from "restartpem#.nc" files               ###
4### and to plot orbital parameters from "obl_ecc_lsp.asc"                                           ###
5#######################################################################################################
6
7import os
8import sys
9import numpy as np
10from glob import glob
11from netCDF4 import Dataset
12import matplotlib.pyplot as plt
13from mpl_toolkits.axes_grid1.inset_locator import inset_axes
14from matplotlib.colors import LinearSegmentedColormap, LogNorm
15from scipy.interpolate import interp1d
16
17
18def get_user_inputs():
19    """
20    Prompt the user for:
21      - folder_path: directory containing NetCDF files (default: "starts")
22      - base_name:   base filename (default: "restartpem")
23      - infofile:    name of the PEM info file (default: "info_PEM.txt")
24    Validates existence of folder and infofile before returning.
25    """
26    folder_path = input(
27        "Enter the folder path containing the NetCDF files "
28        "(press Enter for default [starts]): "
29    ).strip() or "starts"
30    while not os.path.isdir(folder_path):
31        print(f"  » \"{folder_path}\" does not exist or is not a directory.")
32        folder_path = input(
33            "Enter a valid folder path (press Enter for default [starts]): "
34        ).strip() or "starts"
35
36    base_name = input(
37        "Enter the base name of the NetCDF files "
38        "(press Enter for default [restartpem]): "
39    ).strip() or "restartpem"
40
41    infofile = input(
42        "Enter the name of the PEM info file "
43        "(press Enter for default [info_PEM.txt]): "
44    ).strip() or "info_PEM.txt"
45    while not os.path.isfile(infofile):
46        print(f"  » \"{infofile}\" does not exist or is not a file.")
47        infofile = input(
48            "Enter a valid PEM info filename (press Enter for default [info_PEM.txt]): "
49        ).strip() or "info_PEM.txt"
50
51    orbfile = input(
52        "Enter the name of the orbital parameters ASCII file "
53        "(press Enter for default [obl_ecc_lsp.asc]): "
54    ).strip() or "obl_ecc_lsp.asc"
55    while not os.path.isfile(orbfile):
56        print(f"  » \"{orbfile}\" does not exist or is not a file.")
57        orbfile = input(
58            "Enter a valid orbital parameters ASCII filename (press Enter for default [obl_ecc_lsp.asc]): "
59        ).strip() or "info_PEM.txt"
60
61    return folder_path, base_name, infofile, orbfile
62
63
64def list_netcdf_files(folder_path, base_name):
65    """
66    List and sort all NetCDF files matching the pattern {base_name}#.nc
67    in folder_path. Returns a sorted list of full file paths.
68    """
69    pattern = os.path.join(folder_path, f"{base_name}[0-9]*.nc")
70    all_files = glob(pattern)
71    if not all_files:
72        return []
73
74    def extract_index(pathname):
75        fname = os.path.basename(pathname)
76        idx_str = fname[len(base_name):-3]
77        return int(idx_str) if idx_str.isdigit() else float('inf')
78
79    sorted_files = sorted(all_files, key=extract_index)
80    return sorted_files
81
82
83def open_sample_dataset(file_path):
84    """
85    Open a single NetCDF file and extract:
86      - ngrid, nslope
87      - longitude, latitude
88    Returns (ngrid, nslope, longitude_array, latitude_array).
89    """
90    with Dataset(file_path, 'r') as ds:
91        ngrid = ds.dimensions['physical_points'].size
92        nslope = ds.dimensions['nslope'].size
93        longitude = ds.variables['longitude'][:].copy()
94        latitude = ds.variables['latitude'][:].copy()
95    return ngrid, nslope, longitude, latitude
96
97
98def collect_stratification_variables(files, base_name):
99    """
100    Scan all files to collect:
101      - variable names for each stratification property
102      - max number of strata (max_nb_str)
103      - global min base elevation and max top elevation
104    Returns:
105      - var_info: dict mapping each property_name -> sorted list of var names
106      - max_nb_str: int
107      - min_base_elev: float
108      - max_top_elev: float
109    """
110    max_nb_str = 0
111    min_base_elev = np.inf
112    max_top_elev = -np.inf
113
114    property_markers = {
115        'heights':   'stratif_slope',    # "..._top_elevation"
116        'co2_ice':   'h_co2ice',
117        'h2o_ice':   'h_h2oice',
118        'dust':      'h_dust',
119        'pore':      'h_pore',
120        'pore_ice':  'poreice_volfrac'
121    }
122    var_info = {prop: set() for prop in property_markers}
123
124    for file_path in files:
125        with Dataset(file_path, 'r') as ds:
126            if 'nb_str_max' in ds.dimensions:
127                max_nb_str = max(max_nb_str, ds.dimensions['nb_str_max'].size)
128
129            nslope = ds.dimensions['nslope'].size
130            for k in range(1, nslope + 1):
131                var_name = f"stratif_slope{k:02d}_top_elevation"
132                if var_name in ds.variables:
133                    arr = ds.variables[var_name][:]
134                    min_base_elev = min(min_base_elev, np.min(arr))
135                    max_top_elev = max(max_top_elev, np.max(arr))
136                    var_info['heights'].add(var_name)
137
138            for full_var in ds.variables:
139                for prop, marker in property_markers.items():
140                    if (marker in full_var) and prop != 'heights':
141                        var_info[prop].add(full_var)
142
143    for prop in var_info:
144        var_info[prop] = sorted(var_info[prop])
145
146    return var_info, max_nb_str, min_base_elev, max_top_elev
147
148
149def load_full_datasets(files):
150    """
151    Open all NetCDF files and return a list of Dataset objects.
152    (They should be closed by the caller after use.)
153    """
154    return [Dataset(fp, 'r') for fp in files]
155
156
157def extract_stratification_data(datasets, var_info, ngrid, nslope, max_nb_str):
158    """
159    Build:
160      - heights_data[t_idx][isl] = 2D array (ngrid, n_strata_current) of top_elevations.
161      - raw_prop_arrays[prop] = 4D array (ngrid, ntime, nslope, max_nb_str) of per-strata values.
162    Returns:
163      - heights_data: list (ntime) of lists (nslope) of 2D arrays
164      - raw_prop_arrays: dict mapping each property_name -> 4D array
165      - ntime: number of time steps (files)
166    """
167    ntime = len(datasets)
168
169    heights_data = [
170        [None for _ in range(nslope)]
171        for _ in range(ntime)
172    ]
173    for t_idx, ds in enumerate(datasets):
174        for var_name in var_info['heights']:
175            slope_idx = int(var_name.split("slope")[1].split("_")[0]) - 1
176            if 0 <= slope_idx < nslope:
177                raw = ds.variables[var_name][0, :, :]  # (n_strata, ngrid)
178                heights_data[t_idx][slope_idx] = raw.# (ngrid, n_strata)
179
180    raw_prop_arrays = {}
181    for prop in var_info:
182        if prop == 'heights':
183            continue
184        raw_prop_arrays[prop] = np.zeros((ngrid, ntime, nslope, max_nb_str), dtype=np.float32)
185
186    def slope_index_from_var(vname):
187        return int(vname.split("slope")[1].split("_")[0]) - 1
188
189    for prop in raw_prop_arrays:
190        slope_map = {}
191        for vname in var_info[prop]:
192            isl = slope_index_from_var(vname)
193            if 0 <= isl < nslope:
194                slope_map[isl] = vname
195
196        arr = raw_prop_arrays[prop]
197        for t_idx, ds in enumerate(datasets):
198            for isl, var_name in slope_map.items():
199                raw = ds.variables[var_name][0, :, :]  # (n_strata, ngrid)
200                n_strata_current = raw.shape[0]
201                arr[:, t_idx, isl, :n_strata_current] = raw.T
202
203    return heights_data, raw_prop_arrays, ntime
204
205
206def normalize_to_fractions(raw_prop_arrays):
207    """
208    Given raw_prop_arrays for 'co2_ice', 'h2o_ice', 'dust', 'pore' (in meters),
209    normalize each set of strata so that the sum of those four = 1 per cell.
210    Returns:
211      - frac_arrays: dict mapping same keys -> 4D arrays of fractions (0..1).
212    """
213    co2 = raw_prop_arrays['co2_ice']
214    h2o = raw_prop_arrays['h2o_ice']
215    dust = raw_prop_arrays['dust']
216    pore = raw_prop_arrays['pore']
217
218    total = co2 + h2o + dust + pore
219    mask = total > 0.0
220
221    frac_co2 = np.zeros_like(co2, dtype=np.float32)
222    frac_h2o = np.zeros_like(h2o, dtype=np.float32)
223    frac_dust = np.zeros_like(dust, dtype=np.float32)
224    frac_pore = np.zeros_like(pore, dtype=np.float32)
225
226    frac_co2[mask] = co2[mask] / total[mask]
227    frac_h2o[mask] = h2o[mask] / total[mask]
228    frac_dust[mask] = dust[mask] / total[mask]
229    frac_pore[mask] = pore[mask] / total[mask]
230
231    return {
232        'co2_ice': frac_co2,
233        'h2o_ice': frac_h2o,
234        'dust':     frac_dust,
235        'pore':     frac_pore
236    }
237
238
239def read_infofile(file_name):
240    """
241    Reads "info_PEM.txt". Expects:
242      - First line: parameters where the 3rd value is martian_to_earth conversion factor.
243      - Each subsequent line: floats where first value is simulation timestamp (in Mars years).
244    Returns:
245      - date_time: 1D numpy array of timestamps (Mars years)
246      - martian_to_earth: float conversion factor
247    """
248    date_time = []
249    with open(file_name, 'r') as fp:
250        first = fp.readline().split()
251        martian_to_earth = float(first[2])
252        for line in fp:
253            parts = line.strip().split()
254            if not parts:
255                continue
256            try:
257                date_time.append(float(parts[0]))
258            except ValueError:
259                continue
260    return np.array(date_time, dtype=np.float64), martian_to_earth
261
262
263def get_yes_no_input(prompt: str) -> bool:
264    """
265    Prompt the user with a yes/no question. Returns True for yes, False for no.
266    """
267    while True:
268        choice = input(f"{prompt} (y/n): ").strip().lower()
269        if choice in ['y', 'yes']:
270            return True
271        elif choice in ['n', 'no']:
272            return False
273        else:
274            print("Please respond with y or n.")
275
276
277def prompt_discretization_step(max_top_elev):
278    """
279    Prompt for a positive float dz such that 0 < dz <= max_top_elev.
280    """
281    while True:
282        entry = input(
283            "Enter the discretization step of the reference grid for the elevation [m]: "
284        ).strip()
285        try:
286            dz = float(entry)
287            if dz <= 0:
288                print("  » Discretization step must be strictly positive!")
289                continue
290            if dz > max_top_elev:
291                print(
292                    f"  » {dz:.3e} m is greater than the maximum top elevation "
293                    f"({max_top_elev:.3e} m). Please enter a smaller value."
294                )
295                continue
296            return dz
297        except ValueError:
298            print("  » Invalid numeric value. Please try again.")
299
300
301def interpolate_data_on_refgrid(
302    heights_data,
303    prop_arrays,
304    min_base_for_interp,
305    max_top_elev,
306    dz,
307    exclude_sub=False
308):
309    """
310    Build a reference elevation grid and interpolate strata fractions onto it.
311
312    Returns:
313      - ref_grid: 1D array of elevations (nz,)
314      - gridded_data: dict mapping each property_name to 4D array
315        (ngrid, ntime, nslope, nz) with interpolated fractions.
316      - top_index: 3D array (ngrid, ntime, nslope) of ints:
317        number of levels covered by the topmost stratum.
318    """
319    if exclude_sub and (dz > max_top_elev):
320        ref_grid = np.array([0.0, max_top_elev], dtype=np.float32)
321    else:
322        ref_grid = np.arange(min_base_for_interp, max_top_elev + dz/2, dz)
323    nz = len(ref_grid)
324    print(f"> Number of reference grid points = {nz}")
325
326    sample_prop = next(iter(prop_arrays.values()))
327    ngrid, ntime, nslope, max_nb_str = sample_prop.shape
328
329    gridded_data = {
330        prop: np.full((ngrid, ntime, nslope, nz), -1.0, dtype=np.float32)
331        for prop in prop_arrays
332    }
333    top_index = np.zeros((ngrid, ntime, nslope), dtype=np.int32)
334
335    for ig in range(ngrid):
336        for t_idx in range(ntime):
337            for isl in range(nslope):
338                h_mat = heights_data[t_idx][isl]
339                if h_mat is None:
340                    continue
341
342                raw_h = h_mat[ig, :]
343                h_all = np.full((max_nb_str,), np.nan, dtype=np.float32)
344                n_strata_current = raw_h.shape[0]
345                h_all[:n_strata_current] = raw_h
346
347                if exclude_sub:
348                    epsilon = 1e-6
349                    valid_mask = (h_all >= -epsilon)
350                else:
351                    valid_mask = (~np.isnan(h_all)) & (h_all != 0.0)
352
353                if not np.any(valid_mask):
354                    continue
355
356                h_valid = h_all[valid_mask]
357                top_h = np.max(h_valid)
358                i_zmax = np.searchsorted(ref_grid, top_h, side='right')
359                top_index[ig, t_idx, isl] = i_zmax
360                if i_zmax == 0:
361                    continue
362
363                for prop, arr in prop_arrays.items():
364                    prop_profile_all = arr[ig, t_idx, isl, :]
365                    prop_profile = prop_profile_all[valid_mask]
366                    if prop_profile.size == 0:
367                        continue
368
369                    f_interp = interp1d(
370                        h_valid,
371                        prop_profile,
372                        kind='next',
373                        bounds_error=False,
374                        fill_value=-1.0
375                    )
376                    gridded_data[prop][ig, t_idx, isl, :i_zmax] = f_interp(ref_grid[:i_zmax])
377
378    return ref_grid, gridded_data, top_index
379
380
381def plot_stratification_over_time(
382    gridded_data,
383    ref_grid,
384    top_index,
385    heights_data,
386    date_time,
387    exclude_sub=False,
388    output_folder="."
389):
390    """
391    For each grid point and slope, generate a 2×2 figure of:
392      - CO2 ice fraction
393      - H2O ice fraction
394      - Dust fraction
395      - Pore fraction
396    """
397    prop_names = ['co2_ice', 'h2o_ice', 'dust', 'pore']
398    titles = ["CO2 ice", "H2O ice", "Dust", "Pore"]
399    cmap = plt.get_cmap('turbo').copy()
400    cmap.set_under('white')
401    vmin, vmax = 0.0, 1.0
402
403    sample_prop = next(iter(gridded_data.values()))
404    ngrid, ntime, nslope, nz = sample_prop.shape
405
406    if exclude_sub:
407        positive_indices = np.where(ref_grid >= 0.0)[0]
408        sub_ref_grid = ref_grid[positive_indices]
409    else:
410        positive_indices = np.arange(nz)
411        sub_ref_grid = ref_grid
412
413    for ig in range(ngrid):
414        for isl in range(nslope):
415            fig, axes = plt.subplots(2, 2, figsize=(10, 8))
416            fig.suptitle(
417                f"Content variation over time for (Grid point {ig+1}, Slope {isl+1})",
418                fontsize=14,
419                fontweight='bold'
420            )
421
422            # Precompute valid stratum tops per time
423            valid_tops_per_time = []
424            for t_idx in range(ntime):
425                raw_h = heights_data[t_idx][isl][ig, :]
426                h_all = raw_h[~np.isnan(raw_h)]
427                if exclude_sub:
428                    h_all = h_all[h_all >= 0.0]
429                valid_tops_per_time.append(np.unique(h_all))
430
431            for idx, prop in enumerate(prop_names):
432                ax = axes.flat[idx]
433                data_3d = gridded_data[prop][ig, :, isl, :]
434                mat_full = data_3d.T
435                mat = mat_full[positive_indices, :].copy()
436                mat[mat < 0.0] = np.nan
437
438                # Mask above top stratum
439                for t_idx in range(ntime):
440                    i_zmax = top_index[ig, t_idx, isl]
441                    if i_zmax <= positive_indices[0]:
442                        mat[:, t_idx] = np.nan
443                    else:
444                        count_z = np.count_nonzero(positive_indices < i_zmax)
445                        mat[count_z:, t_idx] = np.nan
446
447                im = ax.pcolormesh(
448                    date_time,
449                    sub_ref_grid,
450                    mat,
451                    cmap=cmap,
452                    shading='auto',
453                    vmin=vmin,
454                    vmax=vmax
455                )
456                ax.set_title(titles[idx], fontsize=12)
457                ax.set_xlabel("Time (Mars years)")
458                ax.set_ylabel("Elevation (m)")
459
460            fig.subplots_adjust(right=0.88)
461            fig.tight_layout(rect=[0, 0, 0.88, 1.0])
462            cbar_ax = fig.add_axes([0.90, 0.15, 0.02, 0.7])
463            fig.colorbar(im, cax=cbar_ax, orientation='vertical', label="Content")
464
465            fname = os.path.join(
466                output_folder, f"layering_evolution_ig{ig+1}_is{isl+1}.png"
467            )
468            fig.savefig(fname, dpi=150)
469
470
471def plot_stratification_rgb_over_time(
472    gridded_data,
473    ref_grid,
474    top_index,
475    heights_data,
476    date_time,
477    exclude_sub=False,
478    output_folder="."
479):
480    """
481    Plot stratification over time colored using RGB ternary mix of H2O ice (blue), CO2 ice (violet), and dust (orange).
482    Includes a triangular legend showing the mix proportions.
483    """
484
485    # Define constant colors
486    violet = np.array([255, 0, 255], dtype=float) / 255
487    blue   = np.array([0, 0, 255], dtype=float) / 255
488    orange = np.array([255, 165, 0], dtype=float) / 255
489
490    # Prepare elevation mask
491    mask_elev = (ref_grid >= 0.0) if exclude_sub else np.ones_like(ref_grid, dtype=bool)
492    elev = ref_grid[mask_elev]
493
494    # Generate legend image once
495    res = 300
496    u = np.linspace(0, 1, res)
497    v = np.linspace(0, np.sqrt(3)/2, res)
498    X, Y = np.meshgrid(u, v)
499    V_bary = 2 * Y / np.sqrt(3)
500    U_bary = X - 0.5 * V_bary
501    W_bary = 1 - U_bary - V_bary
502    mask_triangle = (U_bary >= 0) & (V_bary >= 0) & (W_bary >= 0)
503
504    legend_rgb = (
505        U_bary[..., None] * violet
506        + V_bary[..., None] * orange
507        + W_bary[..., None] * blue
508    )
509    legend_rgb = np.clip(legend_rgb, 0.0, 1.0)
510    legend_rgba = np.zeros((res, res, 4))
511    legend_rgba[..., :3] = legend_rgb
512    legend_rgba[..., 3] = mask_triangle.astype(float)
513
514    # Loop over grid and slope
515    h2o = gridded_data['h2o_ice']
516    co2 = gridded_data['co2_ice']
517    dust = gridded_data['dust']
518    ngrid, ntime, nslope, nz = h2o.shape
519
520    for ig in range(ngrid):
521        for isl in range(nslope):
522            # Compute RGB stratification over time
523            rgb = np.ones((nz, ntime, 3), dtype=float)
524            for t in range(ntime):
525                mask_z = np.arange(nz) < top_index[ig, t, isl]
526                if not mask_z.any():
527                    continue
528                cH2O = np.clip(h2o[ig, t, isl, mask_z], 0, None)
529                cCO2 = np.clip(co2[ig, t, isl, mask_z], 0, None)
530                cDust = np.clip(dust[ig, t, isl, mask_z], 0, None)
531                total = cH2O + cCO2 + cDust
532                total[total == 0] = 1.0
533                fH2O = cH2O / total
534                fCO2 = cCO2 / total
535                fDust = cDust / total
536                mix = (
537                    np.outer(fH2O, blue)
538                    + np.outer(fCO2, violet)
539                    + np.outer(fDust, orange)
540                )
541                mix = np.clip(mix, 0.0, 1.0)
542                rgb[mask_z, t, :] = mix
543
544            display_rgb = rgb[mask_elev, :, :]
545
546            # Create figure with legend
547            fig, (ax_main, ax_leg) = plt.subplots(
548                1, 2, figsize=(12, 5), dpi=200,
549                gridspec_kw={'width_ratios': [5, 1]}
550            )
551
552            # Main stratification panel
553            ax_main.imshow(
554                display_rgb,
555                aspect='auto',
556                extent=[date_time[0], date_time[-1], elev.min(), elev.max()],
557                interpolation='nearest',
558                origin='lower'
559            )
560            ax_main.set_facecolor('white')
561            ax_main.set_title(f"Ternary mix over time (Grid point {ig+1}, Slope {isl+1})", fontweight='bold')
562            ax_main.set_xlabel("Time (Mars years)")
563            ax_main.set_ylabel("Elevation (m)")
564
565            # Legend panel
566            ax_leg.imshow(
567                legend_rgba,
568                extent=[0, 1, 0, np.sqrt(3)/2],
569                origin='lower',
570                interpolation='nearest'
571            )
572
573            # Draw triangle border
574            triangle = np.array([[0, 0], [1, 0], [0.5, np.sqrt(3)/2], [0, 0]])
575            ax_leg.plot(triangle[:, 0], triangle[:, 1], 'k-', linewidth=1)
576
577            # Dashed gridlines
578            ticks = np.linspace(0.25, 0.75, 3)
579            for f in ticks:
580                ax_leg.plot([1 - f, 0.5 * (1 - f)], [0, (1 - f)*np.sqrt(3)/2], '--', color='k', linewidth=0.5)
581                ax_leg.plot([f, f + 0.5 * (1 - f)], [0, (1 - f)*np.sqrt(3)/2], '--', color='k', linewidth=0.5)
582                y = (np.sqrt(3)/2) * f
583                ax_leg.plot([0.5 * f, 1 - 0.5 * f], [y, y], '--', color='k', linewidth=0.5)
584
585            # Legend labels
586            ax_leg.text(0, -0.05, 'H2O ice', ha='center', va='top', fontsize=8)
587            ax_leg.text(1, -0.05, 'CO2 ice', ha='center', va='top', fontsize=8)
588            ax_leg.text(0.5, np.sqrt(3)/2 + 0.05, 'Dust', ha='center', va='bottom', fontsize=8)
589            ax_leg.axis('off')
590
591            plt.tight_layout()
592
593            # Save figure
594            fname = os.path.join(output_folder, f"layering_rgb_evolution_ig{ig+1}_is{isl+1}.png")
595            fig.savefig(fname, dpi=150, bbox_inches='tight')
596
597
598def plot_dust_to_ice_ratio_over_time(
599    gridded_data,
600    ref_grid,
601    top_index,
602    heights_data,
603    date_time,
604    exclude_sub=False,
605    output_folder="."
606):
607    """
608    Plot the dust-to-ice ratio in the stratification over time,
609    using a blue-to-orange colormap:
610    - blue: ice-dominated (low dust-to-ice ratio)
611    - orange: dust-dominated (high dust-to-ice ratio)
612    """
613    h2o = gridded_data['h2o_ice']
614    dust = gridded_data['dust']
615    ngrid, ntime, nslope, nz = h2o.shape
616
617    # Elevation mask
618    if exclude_sub:
619        elevation_mask = (ref_grid >= 0.0)
620        elev = ref_grid[elevation_mask]
621    else:
622        elevation_mask = np.ones_like(ref_grid, dtype=bool)
623        elev = ref_grid
624
625    # Define custom blue-to-orange colormap
626    blue = np.array([0, 0, 255], dtype=float) / 255
627    orange = np.array([255, 165, 0], dtype=float) / 255
628    custom_cmap = LinearSegmentedColormap.from_list('BlueOrange', [blue, orange], N=256)
629
630    # Log‑ratio bounds and small epsilon to avoid log(0)
631    vmin, vmax = -2, 1
632    epsilon = 1e-6
633
634    # Loop over grids and slopes
635    for ig in range(ngrid):
636        for isl in range(nslope):
637            log_ratio_array = np.full((nz, ntime), np.nan, dtype=np.float32)
638
639            # Compute log10(dust/ice) profile at each time step
640            for t in range(ntime):
641                zmax = top_index[ig, t, isl]
642                if zmax <= 0:
643                    continue
644
645                h2o_profile = np.clip(h2o[ig, t, isl, :zmax], 0, None)
646                dust_profile = np.clip(dust[ig, t, isl, :zmax], 0, None)
647
648                with np.errstate(divide='ignore', invalid='ignore'):
649                    ratio_profile = np.where(
650                        h2o_profile > 0,
651                        dust_profile / h2o_profile,
652                        10**(vmax + 1)
653                    )
654                    log_ratio = np.log10(ratio_profile + epsilon)
655                    log_ratio = np.clip(log_ratio, vmin, vmax)
656
657                log_ratio_array[:zmax, t] = log_ratio
658
659            # Convert back to linear ratio and apply elevation mask
660            ratio_array = 10**log_ratio_array
661            ratio_display = ratio_array[elevation_mask, :]
662
663            # Plot
664            fig, ax = plt.subplots(figsize=(8, 6), dpi=150)
665            im = ax.imshow(
666                ratio_display,
667                aspect='auto',
668                extent=[date_time[0], date_time[-1], elev.min(), elev.max()],
669                origin='lower',
670                interpolation='nearest',
671                cmap='managua_r',
672                norm=LogNorm(vmin=10**vmin, vmax=10**vmax)
673            )
674
675            # Add colorbar with simplified ratio labels
676            cbar = fig.colorbar(im, ax=ax, orientation='vertical')
677            cbar.set_label('Dust / H₂O ice (ratio)')
678
679            # Define custom ticks and labels
680            ticks = [1e-2, 1e-1, 1, 1e1]
681            labels = ['1:100', '1:10', '1:1', '10:1']
682            cbar.set_ticks(ticks)
683            cbar.set_ticklabels(labels)
684
685            # Save figure
686            plt.tight_layout()
687            fname = os.path.join(
688                output_folder,
689                f"dust_to_ice_ratio_grid{ig+1}_slope{isl+1}.png"
690            )
691            fig.savefig(fname, dpi=150)
692
693
694def plot_strata_count_and_total_height(heights_data, date_time, output_folder="."):
695    """
696    For each grid point and slope, plot:
697      - Number of strata vs time
698      - Total deposit height vs time
699    """
700    ntime = len(heights_data)
701    nslope = len(heights_data[0])
702    ngrid = heights_data[0][0].shape[0]
703
704    for ig in range(ngrid):
705        for isl in range(nslope):
706            n_strata_t = np.zeros(ntime, dtype=int)
707            total_height_t = np.zeros(ntime, dtype=float)
708
709            for t_idx in range(ntime):
710                h_mat = heights_data[t_idx][isl]
711                raw_h = h_mat[ig, :]
712                valid_mask = (~np.isnan(raw_h)) & (raw_h != 0.0)
713                if np.any(valid_mask):
714                    h_valid = raw_h[valid_mask]
715                    n_strata_t[t_idx] = h_valid.size
716                    total_height_t[t_idx] = np.max(h_valid)
717
718            fig, axes = plt.subplots(2, 1, figsize=(8, 6), sharex=True)
719            fig.suptitle(
720                f"Strata count & total height over time for (Grid point {ig+1}, Slope {isl+1})",
721                fontsize=14,
722                fontweight='bold'
723            )
724
725            axes[0].plot(date_time, n_strata_t, marker='+', linestyle='-')
726            axes[0].set_ylabel("Number of strata")
727            axes[0].grid(True)
728
729            axes[1].plot(date_time, total_height_t, marker='+', linestyle='-')
730            axes[1].set_xlabel("Time (Mars years)")
731            axes[1].set_ylabel("Total height (m)")
732            axes[1].grid(True)
733
734            fig.tight_layout(rect=[0, 0, 1, 0.95])
735            fname = os.path.join(
736                output_folder, f"strata_count_height_ig{ig+1}_is{isl+1}.png"
737            )
738            fig.savefig(fname, dpi=150)
739
740
741def read_orbital_data(orb_file, martian_to_earth):
742    """
743    Read the .asc file containing obliquity, eccentricity and Ls p.
744    Columns:
745      0 = time in thousand Martian years
746      1 = obliquity (deg)
747      2 = eccentricity
748      3 = Ls p (deg)
749    Converts times to Earth years.
750    """
751    data = np.loadtxt(orb_file)
752    dates_mka = data[:, 0]
753    dates_yr = dates_mka * 1e3 / martian_to_earth
754    obliquity = data[:, 1]
755    eccentricity = data[:, 2]
756    lsp = data[:, 3]
757    return dates_yr, obliquity, eccentricity, lsp
758
759
760def plot_orbital_parameters(infofile, orb_file, date_time, output_folder="."):
761    """
762    Plot the evolution of obliquity, eccentricity and Ls p
763    versus simulated time.
764    """
765    # Read conversion factor from infofile
766    _, martian_to_earth = read_infofile(infofile)
767
768    # Read orbital data
769    dates_yr, obl, ecc, lsp = read_orbital_data(orb_file, martian_to_earth)
770
771    # Interpolate orbital parameters at simulation dates (date_time)
772    obl_interp = interp1d(dates_yr, obl, kind='linear', bounds_error=False, fill_value="extrapolate")(date_time)
773    ecc_interp = interp1d(dates_yr, ecc, kind='linear', bounds_error=False, fill_value="extrapolate")(date_time)
774    lsp_interp = interp1d(dates_yr, lsp, kind='linear', bounds_error=False, fill_value="extrapolate")(date_time)
775
776    # Plot
777    fig, axes = plt.subplots(3, 1, figsize=(8, 10), sharex=True)
778    fig.suptitle("Orbital Parameters vs Simulated Time", fontsize=14, fontweight='bold')
779
780    axes[0].plot(date_time, obl_interp, 'r+', linestyle='-')
781    axes[0].set_ylabel("Obliquity (°)")
782    axes[0].grid(True)
783
784    axes[1].plot(date_time, ecc_interp, 'b+', linestyle='-')
785    axes[1].set_ylabel("Eccentricity")
786    axes[1].grid(True)
787
788    axes[2].plot(date_time, lsp_interp, 'g+', linestyle='-')
789    axes[2].set_ylabel("Ls p (°)")
790    axes[2].set_xlabel("Time (Mars years)")
791    axes[2].grid(True)
792
793    plt.tight_layout(rect=[0, 0, 1, 0.96])
794    fname = os.path.join(output_folder, "orbital_parameters.png")
795    fig.savefig(fname, dpi=150)
796
797
798def main():
799    # 1) Get user inputs
800    folder_path, base_name, infofile, orbfile = get_user_inputs()
801
802    # 2) List and verify NetCDF files
803    files = list_netcdf_files(folder_path, base_name)
804    if not files:
805        print(f"No NetCDF files named \"{base_name}#.nc\" found in \"{folder_path}\".")
806        sys.exit(1)
807    print(f"> Found {len(files)} NetCDF file(s).")
808
809    # 3) Open one sample to get grid dimensions & coordinates
810    sample_file = files[0]
811    ngrid, nslope, longitude, latitude = open_sample_dataset(sample_file)
812    print(f"> ngrid  = {ngrid}, nslope = {nslope}")
813
814    # 4) Collect variable info + global min/max elevations
815    var_info, max_nb_str, min_base_elev, max_top_elev = collect_stratification_variables(files, base_name)
816    print(f"> max strata per slope = {max_nb_str}")
817    print(f"> min base elev = {min_base_elev:.3f} m, max top elev = {max_top_elev:.3f} m")
818
819    # 5) Load full datasets
820    datasets = load_full_datasets(files)
821
822    # 6) Extract stratification data
823    heights_data, raw_prop_arrays, ntime = extract_stratification_data(datasets, var_info, ngrid, nslope, max_nb_str)
824
825    # 7) Close datasets
826    for ds in datasets:
827        ds.close()
828
829    # 8) Normalize to fractions
830    frac_arrays = normalize_to_fractions(raw_prop_arrays)
831
832    # 9) Ask whether to include subsurface
833    show_subsurface = get_yes_no_input("Show subsurface layers?")
834    exclude_sub = not show_subsurface
835    if exclude_sub:
836        min_base_for_interp = 0.0
837        print("> Interpolating only elevations >= 0 m (surface strata).")
838    else:
839        min_base_for_interp = min_base_elev
840        print(f"> Interpolating full depth down to {min_base_elev:.3f} m.")
841
842    # 10) Prompt discretization step
843    dz = prompt_discretization_step(max_top_elev)
844
845    # 11) Build reference grid and interpolate
846    ref_grid, gridded_data, top_index = interpolate_data_on_refgrid(
847        heights_data, frac_arrays, min_base_for_interp, max_top_elev, dz, exclude_sub=exclude_sub
848    )
849
850    # 12) Read timestamps and conversion factor from infofile
851    date_time, martian_to_earth = read_infofile(infofile)
852    if date_time.size != ntime:
853        print(f"Warning: {date_time.size} timestamps vs {ntime} NetCDF files.")
854
855    # 13) Plot stratification data over time
856    plot_stratification_over_time(
857        gridded_data, ref_grid, top_index, heights_data, date_time,
858        exclude_sub=exclude_sub, output_folder="."
859    )
860    plot_stratification_rgb_over_time(
861        gridded_data, ref_grid, top_index, heights_data, date_time,
862        exclude_sub=exclude_sub, output_folder="."
863    )
864    plot_dust_to_ice_ratio_over_time(
865        gridded_data, ref_grid, top_index, heights_data, date_time,
866        exclude_sub=exclude_sub, output_folder="."
867    )
868    plot_strata_count_and_total_height(heights_data, date_time, output_folder=".")
869
870    # 14) Plot orbital parameters
871    plot_orbital_parameters(infofile, orbfile, date_time, output_folder=".")
872
873    # 15) Show all figures
874    plt.show()
875
876
877if __name__ == "__main__":
878    main()
879
Note: See TracBrowser for help on using the repository browser.