Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Wet Bulb Sounding Computation fails to converge after 50 iterations, value is nan #1705

Closed
jsillin opened this issue Feb 11, 2021 · 7 comments · Fixed by #1811
Closed

Wet Bulb Sounding Computation fails to converge after 50 iterations, value is nan #1705

jsillin opened this issue Feb 11, 2021 · 7 comments · Fixed by #1811
Labels
Area: Calc Pertains to calculations Subarea: Thermo Pertains to thermodynamic calculations and their physical correctness Type: Bug Something is not working like it should
Milestone

Comments

@jsillin
Copy link

jsillin commented Feb 11, 2021

I have arrays with temperature, pressure, and relative humidity corresponding to each mandatory pressure level off the GFS from NOMADS, and am trying to produce profiles of wet bulb temperature with these data.

The particular arrays with which I am working at the moment are:

# dew point, computed using metpy's dew point calculation tool
[18.26819029083491 17.461145433226555 16.336553308346822 13.115820973342743 6.143452735189838 -0.2623966010678873 -3.3370926082535926 -6.550658968811679 -9.53219938820244 -14.238884651017894 -22.66240089860936 -30.84751911458423 -33.68444429930677 -33.758997581695766 -44.84378577284087 -43.464276398501696 -60.91091595765517 -59.34560239336541 -77.78116414272493 -77.3133221938806 -82.27649829610058 -83.63567550304444 -84.52241968178798 nan -90.08663879090643 -88.4437696852317 nan -89.92561301616831 nan nan nan nan nan nan] degree_Celsius

# temperature
<xarray.DataArray 'temperature' (lev: 34)>
array([ 21.38790283,  19.25753174,  17.29430542,  16.52999878,
        16.88999023,  15.28999634,  12.13189697,   8.73193359,
         5.1416748 ,   1.52350464,  -2.15553589,  -6.14998169,
       -11.08999634, -17.96802979, -25.26040649, -32.41000061,
       -41.09917908, -45.71267395, -48.46999512, -59.97120056,
       -74.6106781 , -79.78999329, -74.53684082, -63.49831848,
       -58.1831604 , -53.39454041, -48.2122467 , -39.66942444,
       -34.30999451, -33.35002747, -23.31001892, -14.55002441,
       -11.74678955, -25.58999634])
Coordinates:
    time     datetime64[ns] 2021-02-11T08:00:00
  * lev      (lev) float64 1e+03 975.0 950.0 925.0 900.0 ... 5.0 3.0 2.0 1.0 0.4
    crs      object Projection: latitude_longitude
    lat      float64 27.15
    lon      float64 282.0

# relative humidity
<xarray.DataArray 'rh' (lev: 34)>
array([82.40000153, 89.34000244, 94.1       , 80.24000244, 49.10000076,
       34.54000168, 33.82000046, 33.28000183, 33.78000183, 29.89999962,
       19.09999962, 12.18000011, 13.62000027, 23.77999992, 14.32000027,
       32.41999969,  9.96000023, 19.98000031,  2.02000008,  8.66000004,
       29.78000031, 53.14000092, 20.20000114,         nan,  0.80000001,
        0.60000002,         nan,  0.1       ,  0.        ,  0.        ,
        0.        ,  0.        ,  0.        ,         nan])
Coordinates:
    time     datetime64[ns] 2021-02-11T08:00:00
  * lev      (lev) float64 1e+03 975.0 950.0 925.0 900.0 ... 5.0 3.0 2.0 1.0 0.4
    crs      object Projection: latitude_longitude
    lat      float64 27.15
    lon      float64 282.0
Attributes:
    long_name:  ** (1000 975 950 925 900.. 7 5 3 2 1) relative humidity [%]

when I attempt to combine these into wet bulb using the mpcalc.wet_bulb_temperature method

sound_wb = mpcalc.wet_bulb_temperature(sound_pres.values*units.hPa,sound_temps.values*units.degC,sound_dp)

I get hit with

Traceback (most recent call last):
  File "gfs_hrly_ptype_soundings_wetbulb.py", line 327, in <module>
    sound_wb = mpcalc.wet_bulb_temperature(sound_pres.values*units.hPa,sound_temps.values*units.degC,sound_dp)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\xarray.py", line 677, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\units.py", line 320, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\calc\thermo.py", line 2410, in wet_bulb_temperature
    lcl_pressure, lcl_temperature = lcl(press, temp, dewp)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\xarray.py", line 677, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\units.py", line 320, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\calc\thermo.py", line 357, in lcl
    lcl_p = so.fixed_point(_lcl_iter, pressure.m, args=(pressure.m, w, temperature),
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\scipy\optimize\minpack.py", line 937, in fixed_point
    return _fixed_point_helper(func, x0, args, xtol, maxiter, use_accel)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\scipy\optimize\minpack.py", line 891, in _fixed_point_helper
    raise RuntimeError(msg)
RuntimeError: Failed to converge after 50 iterations, value is nan

The expected output should be an array of wet bulb temperature values.

I am running metpy 0.12.2 and python 3.8.5 on windows with conda

I was chatting with @jthielen about this over on the discussion board earlier today, and perhaps that thread may have more information that would be helpful #1702

@jsillin jsillin added the Type: Bug Something is not working like it should label Feb 11, 2021
@jthielen
Copy link
Collaborator

jthielen commented Feb 12, 2021

I was able to reproduce this on latest master branch with the following copy/pasteable example derived from @jsillin's example:

import numpy as np
from metpy.units import units
import metpy.calc as mpcalc

pressure = [
    1000.0, 975.0, 950.0, 925.0, 900.0,
    850.0, 800.0, 750.0, 700.0, 650.0,
    600.0, 550.0, 500.0, 450.0, 400.0,
    350.0, 300.0, 250.0, 200.0, 150.0,
    100.0, 70.0, 50.0, 40.0, 30.0,
    20.0, 15.0, 10.0, 7.0, 5.0,
    3.0, 2.0, 1.0, 0.4
] * units.hPa
dewpoint = [
    18.26819029083491, 17.461145433226555, 16.336553308346822, 13.115820973342743, 6.143452735189838,
    -0.2623966010678873, -3.3370926082535926, -6.550658968811679, -9.53219938820244, -14.238884651017894,
    -22.66240089860936, -30.84751911458423, -33.68444429930677, -33.758997581695766, -44.84378577284087,
    -43.464276398501696, -60.91091595765517, -59.34560239336541, -77.78116414272493, -77.3133221938806,
    -82.27649829610058, -83.63567550304444, -84.52241968178798, np.nan, -90.08663879090643,
    -88.4437696852317, np.nan, -89.92561301616831, np.nan, np.nan,
    np.nan, np.nan, np.nan, np.nan
] * units.degC
temperature = [
    21.38790283,  19.25753174,  17.29430542,  16.52999878,
    16.88999023,  15.28999634,  12.13189697,   8.73193359,
    5.1416748 ,   1.52350464,  -2.15553589,  -6.14998169,
    -11.08999634, -17.96802979, -25.26040649, -32.41000061,
    -41.09917908, -45.71267395, -48.46999512, -59.97120056,
    -74.6106781 , -79.78999329, -74.53684082, -63.49831848,
    -58.1831604 , -53.39454041, -48.2122467 , -39.66942444,
    -34.30999451, -33.35002747, -23.31001892, -14.55002441,
    -11.74678955, -25.58999634
] * units.degC

mpcalc.wet_bulb_temperature(pressure, temperature, dewpoint)

I'm pretty sure this is a case of metpy.calc.lcl, and so in turn metpy.calc.wet_bulb_temperature, not safely handling NaNs.

For a temporary workaround, the NaN values can be removed, and then those positions replaced with NaNs in the output with the following

nan_mask = np.isnan(pressure) | np.isnan(temperature) | np.isnan(dewpoint)
idx = np.arange(pressure.size)[~nan_mask]
wetbulb_valid_only = mpcalc.wet_bulb_temperature(pressure[idx], temperature[idx], dewpoint[idx])
wetbulb_full = np.full(pressure.size, np.nan) * wetbulb_valid_only.units
wetbulb_full[idx] = wetbulb_valid_only

returning

[19.238071735308814 18.033294139060633 16.65179610640866 14.341431051131467 10.713278865013098 7.2703265785039 4.63234087372236 1.8324379627895773 -1.0154897545814394 -4.337334561885717 -8.23596994210175 -11.902727397896111 -15.669313544076992 -20.78875735056887 -27.342629368898884 -33.42946313024462 -41.8026422159221 -46.17279371976847 -48.99179569697857 -60.13602538549741 -74.63516192059605 -79.80028104362006 -74.59295016865613 nan -59.20059644897026 -55.76040402608365 nan -49.692666440433335 nan nan nan nan nan nan] degree_Celsius

(there's probably a cleaner way to do this, but this is the first workaround that came to mind)

@jthielen jthielen added Area: Calc Pertains to calculations Subarea: Thermo Pertains to thermodynamic calculations and their physical correctness labels Feb 12, 2021
@jthielen jthielen added this to the 1.0.1 milestone Feb 12, 2021
@jsillin
Copy link
Author

jsillin commented Mar 4, 2021

hi @jthielen I've been using this workaround for a bit but started to have some troubles today. I'm feeding it many many sounding arrays and lots will work until all of the sudden I get hit with:

[5.52528414196141 4.229539910567291 2.895846874131037 1.5887699256220666 0.20390642244782042 -2.654507913641146 -5.632318081430397 -8.712560046968134 -12.172700146679722 -14.720381827125612 -13.122452644198267 -17.154515813521805 -19.689468383789066 -25.915908813476566 -32.74919738769529 -40.69527893066405 -49.62999877929687 -66.24284775499198 -69.81460466032047 -84.71954949433804 -81.78304949954776 -84.90852851671951 -87.67998463101769 nan -88.34171272692159 -88.00866993036921 nan -91.40986729544244 -93.52735138498159 nan nan nan nan nan] degree_Celsius
<xarray.DataArray 'temperature' (lev: 34)>
array([ 15.87064819,  14.44068604,  13.0194519 ,  11.60307007,
        10.10307007,   7.009198  ,   3.78922729,   0.46185303,
        -2.76993408,  -3.07687378,  -8.50082397, -14.13804321,
       -19.68946838, -25.91590881, -32.74919739, -40.69527893,
       -49.62999878, -55.8655304 , -59.51889038, -53.16999817,
       -59.6422821 , -60.88999634, -61.15001526, -59.4051178 ,
       -59.03000183, -56.10873108, -51.39000244, -48.59403381,
       -45.96583252, -39.56842651, -33.40492249, -22.87559509,
       -11.52073364, -27.73213196])
Coordinates:
    time     datetime64[ns] 2021-03-04T12:00:00
  * lev      (lev) float64 1e+03 975.0 950.0 925.0 900.0 ... 5.0 3.0 2.0 1.0 0.4
    crs      object Projection: latitude_longitude
    lat      float64 40.25
    lon      float64 254.3
Traceback (most recent call last):
  File "get_soundingmaps.py", line 349, in <module>
    smap.plot_soundings(fig,ax1,prs_temps,prs_relh,sfc_pressure,centerlat,centerlon,domainsize,model,cape=cape_bool,wetbulb=wetbulb_bool)
  File "C:\Users\Jack\Documents\GitHub\plotting\dev\soundingmaps.py", line 110, in plot_soundings
    sound_wb = spt.wetbulb_with_nan(sound_pres,sound_temps.data*units.degC,sound_dp)
  File "C:\Users\Jack\Documents\GitHub\plotting\dev\supplementary_tools.py", line 10, in wetbulb_with_nan
    wetbulb_valid_only = mpcalc.wet_bulb_temperature(pressure[idx], temperature[idx], dewpoint[idx])
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\xarray.py", line 677, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\units.py", line 320, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\calc\thermo.py", line 2413, in wet_bulb_temperature
    ret[...] = moist_adiabat_temperatures[-1].magnitude
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\pint\quantity.py", line 1834, in __getitem__
    return type(self)(self._magnitude[key], self._units)
IndexError: index -1 is out of bounds for axis 0 with size 0

where the first two arrays are the inputs dewpoint and temperature

@jsillin
Copy link
Author

jsillin commented Mar 4, 2021

update: this also happens when I just use the normal mpcalc.wet_bulb_temperature function with data that contains no NaN values...

Works 20 or 30 times in a for loop then all of the sudden doesn't.

Input data for a time when things worked:

[-3.4900620834804332 -4.745168817082063 -5.996857243612559 -7.3056621650349225 -11.01318205991729 -17.000809564988824 -18.66687795344326 -27.285804163612678 -21.239492164330507 -20.652501068519875 -25.68935160925538 -32.444251243020844 -33.75374378319033 -39.11298035000362 -46.61876182607188 -55.01624807787972 -65.65506180078174] degree_Celsius
<xarray.DataArray 'temperature' (lev: 17)>
array([ 11.94066772,  10.51946411,   9.10308228,   7.6230835 ,
         8.44921265,   8.10922852,   4.2418457 ,   0.95004883,
        -3.83686523,  -8.80083618, -14.2180481 , -19.72946777,
       -25.37590332, -32.12919006, -39.19527893, -47.16999817,
       -55.94552917])
Coordinates:
    time     datetime64[ns] 2021-03-04T12:00:00
  * lev      (lev) float64 975.0 950.0 925.0 900.0 ... 400.0 350.0 300.0 250.0
    crs      object Projection: latitude_longitude
    lat      float64 40.25
    lon      float64 257.2

input for the one that failed:

[4.229539910567291 2.895846874131037 1.5887699256220666 0.20390642244782042 -2.654507913641146 -5.632318081430397 -8.712560046968134 -12.172700146679722 -14.720381827125612 -13.122452644198267 -17.154515813521805 -19.689468383789066 -25.915908813476566 -32.74919738769529 -40.69527893066405 -49.62999877929687 -66.24284775499198] degree_Celsius
<xarray.DataArray 'temperature' (lev: 17)>
array([ 14.44068604,  13.0194519 ,  11.60307007,  10.10307007,
         7.009198  ,   3.78922729,   0.46185303,  -2.76993408,
        -3.07687378,  -8.50082397, -14.13804321, -19.68946838,
       -25.91590881, -32.74919739, -40.69527893, -49.62999878,
       -55.8655304 ])
Coordinates:
    time     datetime64[ns] 2021-03-04T12:00:00
  * lev      (lev) float64 975.0 950.0 925.0 900.0 ... 400.0 350.0 300.0 250.0
    crs      object Projection: latitude_longitude
    lat      float64 40.25
    lon      float64 254.3
Traceback (most recent call last):
  File "get_soundingmaps.py", line 350, in <module>
    smap.plot_soundings(fig,ax1,prs_temps,prs_relh,sfc_pressure,centerlat,centerlon,domainsize,model,cape=cape_bool,wetbulb=wetbulb_bool)
  File "C:\Users\Jack\Documents\GitHub\plotting\dev\soundingmaps.py", line 110, in plot_soundings
    sound_wb = mpcalc.wet_bulb_temperature(sound_pres,sound_temps.data*units.degC,sound_dp)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\xarray.py", line 677, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\units.py", line 320, in wrapper
    return func(*args, **kwargs)
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\metpy\calc\thermo.py", line 2413, in wet_bulb_temperature
    ret[...] = moist_adiabat_temperatures[-1].magnitude
  File "C:\Users\Jack\anaconda3\envs\metpy2\lib\site-packages\pint\quantity.py", line 1834, in __getitem__
    return type(self)(self._magnitude[key], self._units)
IndexError: index -1 is out of bounds for axis 0 with size 0

@jthielen
Copy link
Collaborator

jthielen commented Mar 4, 2021

This most recent error looks awfully similar to #1332, which was fixed in MetPy 1.0 (which I noticed you aren't using based on the code in the traceback). Does it only/always occur when you have perfectly saturated layer(s) in your profile? Is it feasible for you to try these same profiles in MetPy 1.0?

@jsillin
Copy link
Author

jsillin commented Mar 5, 2021

hmm I have made a new conda environment with 1.0 and so far not so good:

C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\coding\times.py:113: SerializationWarning: Ambiguous reference date string: 1-1-1 00:00:0.0. The first value is assumed to be the year hence will be padded with zeros to remove the ambiguity (the padded reference date string is: 0001-1-1 00:00:0.0). To remove this message, remove the ambiguity by padding your reference date strings with zeros.
  warnings.warn(warning_msg, SerializationWarning)
Traceback (most recent call last):
  File "get_soundingmaps.py", line 265, in <module>
    data = ds.metpy.parse_cf()
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\metpy\xarray.py", line 712, in parse_cf
    subset = xr.merge([self.parse_cf(single_varname, coordinates=coordinates)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\metpy\xarray.py", line 712, in <listcomp>
    subset = xr.merge([self.parse_cf(single_varname, coordinates=coordinates)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\metpy\xarray.py", line 750, in parse_cf
    if crs is None and not check_axis(var, 'latitude', 'longitude'):
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\metpy\xarray.py", line 1045, in check_axis
    and str(var.metpy.units)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\metpy\xarray.py", line 131, in units
    if isinstance(self._data_array.data, units.Quantity):
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\dataarray.py", line 625, in data
    return self.variable.data
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\variable.py", line 374, in data
    return self.values
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\variable.py", line 554, in values
    return _as_array_or_item(self._data)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\variable.py", line 287, in _as_array_or_item
    data = np.asarray(data)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\numpy\core\_asarray.py", line 83, in asarray
    return array(a, dtype, copy=False, order=order)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\indexing.py", line 693, in __array__
    self._ensure_cached()
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\indexing.py", line 690, in _ensure_cached
    self.array = NumpyIndexingAdapter(np.asarray(self.array))
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\numpy\core\_asarray.py", line 83, in asarray
    return array(a, dtype, copy=False, order=order)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\indexing.py", line 663, in __array__
    return np.asarray(self.array, dtype=dtype)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\numpy\core\_asarray.py", line 83, in asarray
    return array(a, dtype, copy=False, order=order)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\indexing.py", line 568, in __array__
    return np.asarray(array[self.key], dtype=None)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\numpy\core\_asarray.py", line 83, in asarray
    return array(a, dtype, copy=False, order=order)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\coding\variables.py", line 70, in __array__
    return self.func(self.array)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\coding\variables.py", line 138, in _apply_mask
    data = np.asarray(data, dtype=dtype)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\numpy\core\_asarray.py", line 83, in asarray
    return array(a, dtype, copy=False, order=order)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\indexing.py", line 568, in __array__
    return np.asarray(array[self.key], dtype=None)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\backends\netCDF4_.py", line 85, in __getitem__
    return indexing.explicit_indexing_adapter(
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\core\indexing.py", line 853, in explicit_indexing_adapter
    result = raw_indexing_method(raw_key.tuple)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\backends\netCDF4_.py", line 98, in _getitem
    array = getitem(original_array, key)
  File "C:\Users\Jack\anaconda3\envs\metpynew\lib\site-packages\xarray\backends\common.py", line 53, in robust_getitem
    return array[key]
  File "src\netCDF4\_netCDF4.pyx", line 4420, in netCDF4._netCDF4.Variable.__getitem__
  File "src\netCDF4\_netCDF4.pyx", line 5363, in netCDF4._netCDF4.Variable._get
  File "src\netCDF4\_netCDF4.pyx", line 1950, in netCDF4._netCDF4._ensure_nc_success
RuntimeError: NetCDF: file not found

takes >20min for that to pop up (usually <1min to read in data)

@jthielen
Copy link
Collaborator

jthielen commented Mar 5, 2021

Ooh, now that's a strange one. I can't know for sure without playing around with the exact data and process your using, but I think there are two things going on:

  1. NOMADS (or wherever you're getting the data from) is being flaky (due to that the warning about weird time units, and that I've seen that NetCDF error when it can't read properly from the remote source)
  2. Based on your traceback, it looks like we've discovered a (pretty serious) bug in parse_cf: if your dataset is not CF-compliant (i.e., missing the grid_mapping attribute), it will load all the data into memory in the process of performing unit checks for lat/lon coordinates. This is likely why it is taking so long. I'll open a separate issue for this to hopefully be resolved soon in 1.0.1.

In the meantime, you could try adding the CRS manually with .metpy.assign_crs(grid_mapping_name='latitude_longitude') instead of using .metpy.parse_cf (assuming you're on a lat/lon grid, otherwise, adjust the CF arguments accordingly).

@jsillin
Copy link
Author

jsillin commented Mar 5, 2021

The exact data in question is NOMADS 18z GFS and is working fine on several other scripts using 0.12.1.

I think your idea of a bug in parse_cf is a good one. Will give the manual CRS a try. Thanks for the suggestion!

dopplershift added a commit to dopplershift/MetPy that referenced this issue Apr 14, 2021
This required fixing lcl and moist_lapse to work with nans, or at least
work better.
dopplershift added a commit to dopplershift/MetPy that referenced this issue Apr 14, 2021
This required fixing lcl and moist_lapse to work with nans, or at least
work better.
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
Area: Calc Pertains to calculations Subarea: Thermo Pertains to thermodynamic calculations and their physical correctness Type: Bug Something is not working like it should
Projects
None yet
Development

Successfully merging a pull request may close this issue.

2 participants