HydroGrid.py 12.8 KB
Newer Older
1 2 3 4 5 6 7 8
#
# Reads ans manages the model grid
#
import numpy as np
from netCDF4 import Dataset
# Doc : https://spacetelescope.github.io/spherical_geometry/api/spherical_geometry.polygon.SphericalPolygon.html
from spherical_geometry import polygon
import RPPtools as RPP
9
import os 
10
import sys
11
#
12 13
import getargs
config = getargs.SetupConfig()
14
EarthRadius=config.getfloat("OverAll", "EarthRadius", fallback=6370000.0)
15
FloodplainsFile=config.get("OverAll", "FloodplainsFile", fallback=None)
16 17 18 19 20 21 22
#
import getargs
log_master, log_world = getargs.getLogger(__name__)
INFO, DEBUG, ERROR = log_master.info, log_master.debug, log_world.error
INFO_ALL, DEBUG_ALL = log_world.info, log_world.debug
#
def getbox(ncdf, corners) :
23 24 25 26 27 28 29 30 31 32 33 34
    #
    # Check that the domain fits in the data available.
    lonrange=[getattrcontaining(ncdf, "nav_lon", "min"), getattrcontaining(ncdf, "nav_lon", "max")]
    latrange=[getattrcontaining(ncdf, "nav_lat", "min"), getattrcontaining(ncdf, "nav_lat", "max")]
    if len(lonrange[0]) <= 0 :
        lonrange=[np.min(ncdf.variables["nav_lon"][:,:]),np.max(ncdf.variables["nav_lon"][:,:])]
    if len(latrange[0]) <= 0 :
        latrange=[np.min(ncdf.variables["nav_lat"][:,:]),np.max(ncdf.variables["nav_lat"][:,:])]
    if np.min(corners[0][:]) < np.min(lonrange) or \
       np.max(corners[0][:]) > np.max(lonrange) or \
       np.min(corners[1][:]) < np.min(latrange) or \
       np.max(corners[1][:]) > np.max(latrange) :
POLCHER Jan's avatar
POLCHER Jan committed
35
        ERROR("The atmospheric domain does not fit inside the area covered by HydroFile. Please check the file.")
36 37
        sys.exit()
    #
38 39 40 41 42 43
    # Add a few points to tbe box to make sure to cover everything
    halo_pts=2
    ii=[np.argmin(abs(ncdf.variables["nav_lon"][0,:]-np.min(corners[0][:]))),
        np.argmin(abs(ncdf.variables["nav_lon"][0,:]-np.max(corners[0][:])))]
    jj=[np.argmin(abs(ncdf.variables["nav_lat"][:,0]-np.min(corners[1][:]))),
        np.argmin(abs(ncdf.variables["nav_lat"][:,0]-np.max(corners[1][:])))]
44

45
    return min(ii)-halo_pts, max(ii)+halo_pts, min(jj)-halo_pts, max(jj)+halo_pts
46 47 48
#
# Get the corners for a regular lat/lon grid
#
49
def centers(lon, lat) :
50
    jjm,iim = lon.shape
51 52
    #
    index = [[j,i] for i in range(iim) for j in range(jjm)]
53 54
    centersll = [[lon[j,i], lat[j,i]] for j,i in index]
    return index, centersll
55

56 57 58
def corners(subcentersll, hdlon, hdlat) :
    cornersll = [RPP.boxit(np.array([lo, la]), hdlon, hdlat, 2) for lo,la in subcentersll]
    Llats = [[0.0000001 if p[1] == 0 else np.sign(p[1])*89.99 if np.abs(p[1]) == 90 else p[1] for p in boxll] for boxll in cornersll]
59
    Llons = [[0.0000001 if p[0] == 0 else np.sign(p[0])*89.99 if np.abs(p[0]) == 90 else p[0] for p in boxll] for boxll in cornersll]
60 61
    cornerspoly = [ polygon.SphericalPolygon.from_lonlat(lons, lats, center=cent) for lons, lats, cent in zip(Llons, Llats, subcentersll)]
    #radiusll = [RPP.maxradius(np.array(cent), np.array(lons), np.array(lats)) for cent, lons, lats in zip(subcentersll,Llons,Llats)]
62
    #
63
    return cornersll, cornerspoly
64
#
65
def gather(x, index, default = 0) :
66 67
    y=[]
    for ia in index :
68
      try:
69 70 71 72
        if (ia[:,0] == [-1, -1]).all():
            y.append([default])
        else:
            y.append(list(x[ia[0,i],ia[1,i]] for i in range(ia.shape[1]) ))
73 74 75
      except:
        print("Error - the hydrological file may not recover all the domain", ia)
        sys.exit()
76 77 78 79 80 81 82 83 84
    return y
#
def getattrcontaining(nc, varname, substr) :
    att=[]
    for s in nc.variables[varname].ncattrs() :
        if s.lower().find(substr) >= 0 :
            att.append(nc.variables[varname].getncattr(s))
    return att
#
85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109
def calctopoindex(nf, istr, iend, jstr, jend, index, hydrogrid) :
    # The routing code in Fortran (i=1, j=2)
    # 8 --- 1 --- 2
    # |     |     |
    # 7 --- X --- 3
    # |     |     |
    # 6 --- 5 --- 4
    #
    inc = np.array([[0,1,1,1,0,-1,-1,-1],[-1,-1,0,1,1,1,0,-1]])
    elevation_precision=0.1
    #
    jlen,ilen = nf.variables["orog"].shape
    missing = nf.variables["orog"].missing_value
    land = np.where(nf.variables["orog"][jstr:jend,istr:iend] < missing)
    topoind=np.zeros((jend-jstr,iend-istr))
    topoind[:,:] = np.nan
    #
    for j,i in zip(land[:][0],land[:][1]) :
        trip = int(nf.variables["trip"][jstr+j,istr+i]-1)
        if trip < 8 :
            inf=(istr+i+inc[0,trip])%ilen
            jnf=(jstr+j+inc[1,trip])%jlen
            dz=max(elevation_precision,nf.variables["orog"][jstr+j,istr+i]-nf.variables["orog"][jnf,inf])
            dl=nf.variables["disto"][jstr+j,istr+i]-nf.variables["disto"][jnf,inf]
        else :
POLCHER Jan's avatar
POLCHER Jan committed
110 111 112 113 114 115
            if trip > 97 :
                # Going into the ocean which is at elevation 0
                dz = max(elevation_precision,nf.variables["orog"][jstr+j,istr+i])
            else :
                # Endorehic basin so we use thte lowest dz possible.
                dz = elevation_precision
116
            #
POLCHER Jan's avatar
POLCHER Jan committed
117 118 119 120 121 122
            if nf.variables["disto"][jstr+j,istr+i] > 0 :
                dl = nf.variables["disto"][jstr+j,istr+i]
            else :
                # Some other reasonable assumption
                rx,ry = hydrogrid.resolution(i,j)
                dl = (rx+ry)*0.5*0.5
123 124 125 126 127 128 129
        topoind[j,i]=np.sqrt(dl**3./(dz*1.0e6))
        #
    minitopo = np.nanmin(topoind)
    if minitopo <= np.finfo(float).eps :
        ERROR("Topoind close to zero encoutered.")
        sys.exit()
    #
130
    return gather(topoind, index, default = 10), "Topographic index which serves for the residence time", minitopo
131
#
132
class HydroGrid :
133
    def __init__(self, lolacorners, wfile) :
134
        #
135 136 137
        self.source=config.get("OverAll", "HydroFile")
        INFO("Opening in HydroGrid : "+self.source)
        self.ncfile=Dataset(self.source,'r')
138 139 140 141 142 143 144 145 146 147 148 149 150 151
        #
        # Test if floodplains is correct
        if FloodplainsFile is not None:
           try:
              nfp = Dataset(FloodplainsFile, "r")
              njgfp, nigfp = nfp.variables["floodplains"].shape
              njg, nig = self.ncfile.variables["nav_lon"].shape
              if ((njg != njgfp) or (nig != nigfp)):
                 ERROR("Invalid FloodplainsFile format. Hydrogrod: {0} / FPfile: {1}".format((njg, nig),(njgfp, nigfp)))
                 sys.exit()
           except:
              ERROR("Invalid FloodplainsFile :{0}".format(FloodplainsFile))
              sys.exit()
        #
152
        istr, iend, jstr, jend = getbox(self.ncfile, lolacorners)
153 154
        self.hdlon = np.mean(np.abs(np.diff(self.ncfile.variables["nav_lon"][0,:])))
        self.hdlat = np.mean(np.abs(np.diff(self.ncfile.variables["nav_lat"][:,0])))
155
        self.hd = max(self.hdlon,self.hdlat)
156 157 158 159 160
        # In case the box is on the edge of the routing.nc file
        if istr<0:
            istr = 0  
        if iend >=self.ncfile.variables["nav_lon"][0,:].shape[0]: 
            iend = self.ncfile.variables["nav_lon"][0,:].shape[0]
161 162 163 164
        if jstr<0:
            jstr = 0
        if iend >=self.ncfile.variables["nav_lat"][0,:].shape[0]:
            iend = self.ncfile.variables["nav_lat"][0,:].shape[0]
165
        self.box=[istr, iend, jstr, jend]   
166 167 168
        self.lon=np.copy(self.ncfile.variables["nav_lon"][jstr:jend,istr:iend])
        self.lat=np.copy(self.ncfile.variables["nav_lat"][jstr:jend,istr:iend])
        self.jjm,self.iim = self.lon.shape
169 170 171 172
        DEBUG("# Dimensions :"+str(self.iim)+" -- "+str(self.jjm))
        if self.iim > 0 and self.jjm > 0 :
            DEBUG("# Range Lon :"+str(np.min(self.lon))+" -- "+str(np.max(self.lon)))
            DEBUG("# Range Lat :"+str(np.min(self.lat))+" -- "+str(np.max(self.lat)))
173
        #
174 175
        if wfile is None or not os.path.exists(wfile):
            self.index, self.centers = centers(self.lon, self.lat) 
176 177

    def select(self, c, r) :
178
        indices = [i for i in range(len(self.centers)) if (RPP.loladist(np.array(c),np.array(self.centers[i])) <= r+1.1*self.hd)]
179 180 181
        subcenters = [self.centers[i] for i in indices]
        polyll, polylist = corners(subcenters, self.hdlon, self.hdlat)
        return indices, polylist
182

183 184 185 186 187 188 189 190 191 192 193 194 195
    def resolution(self,i,j) :
        if i+1 >= self.iim or j+1 >= self.jjm :
            inx=i-1; jnx=j-1
        else :
            inx=i+1; jnx=j+1
        dlon = np.radians(np.abs(self.lon[j,i]-self.lon[j,inx]))
        alon = np.cos(np.radians(self.lat[j,i]))**2*(np.sin(dlon/2))**2
        clon = 2*np.arctan2(np.sqrt(alon), np.sqrt(1-alon))
        dlat = np.radians(np.abs(self.lat[j,i]-self.lat[jnx,i]))
        alat = (np.sin(dlat/2))**2
        clat = 2*np.arctan2(np.sqrt(alat), np.sqrt(1-alat))
        return EarthRadius*clon, EarthRadius*clat
    
196
class HydroData :
197
    def __init__(self, nf, box, index, part, hydrogrid) :
198
        istr, iend, jstr, jend = box[:]
199 200 201
        #
        # Flow direction
        #
202
        self.trip=gather(nf.variables["trip"][jstr:jend,istr:iend].astype(np.float32), index, default = 98)
203 204
        self.tripdesc=nf.variables["trip"].long_name
        #
205 206
        # ID of basin
        #
207
        self.basins=gather(nf.variables["basins"][jstr:jend,istr:iend], index, 999)
208 209 210 211 212
        self.basinsdesc=nf.variables["basins"].long_name
        att = getattrcontaining(nf, "basins", "max")
        if len(att) > 0 :
            self.basinsmax=att[0]
        else :
213
            INFO("We need to scan full file to find the largest basin ID")
214
            # This variable seems not to be used further
215 216 217
            ##self.basinsmax = part.domainmax(np.ma.max(np.ma.masked_where(self.basins < 1.e10, self.basins)))
            ##self.basinsmax = part.domainmax(np.max(np.array(self.basins)[np.array(self.basins) < missing]))
            self.basinsmax = 999
218
        #
219
        # Distance to ocean
220
        #
221
        self.disto=gather(nf.variables["disto"][jstr:jend,istr:iend].astype(np.float32), index, default = 500)
222 223
        self.distodesc=nf.variables["disto"].long_name
        #
224 225
        # Flow accumulation
        #
226
        self.fac=gather(nf.variables["fac"][jstr:jend,istr:iend].astype(np.float32), index, default = 0)
227
        self.facdesc=nf.variables["fac"].long_name
228
        #
229 230 231
        # If Orography is present in the file then we can compute the topographic index.
        # Else the topographic index needs to be present in the file.
        #
POLCHER Jan's avatar
POLCHER Jan committed
232
        if "orog" in nf.variables.keys() and "disto" in nf.variables.keys() :
233
            self.orog = gather(nf.variables["orog"][jstr:jend,istr:iend].astype(np.float32), index, default = 0)
234
            self.orogdesc=nf.variables["orog"].long_name
235 236
            self.topoind, self.topoinddesc, mintopo = calctopoindex(nf, istr, iend, jstr, jend, index, hydrogrid)
            self.topoindmin = part.domainmin(mintopo)
237
        elif "topoind" in nf.variables.keys() and "orog" in nf.variables.keys() :
238
            self.orog = gather(np.zeros((jend-jstr,iend-istr)).astype(np.float32), index)
239
            self.topoind=gather(nf.variables["topoind"][jstr:jend,istr:iend].astype(np.float32), index, default = 10)
240 241 242 243 244 245 246
            self.topoinddesc=nf.variables["topoind"].long_name
            att = getattrcontaining(nf, "topoind", "min")
            if len(att) > 0 :
                self.topoindmin=att[0]
            else :
                INFO("We need to scan full file to find minimum topoind over domain")
                self.topoindmin=np.min(np.where(nf.variables["topoind"][:,:] <  1.e15))
POLCHER Jan's avatar
POLCHER Jan committed
247 248 249
        else :
            ERROR("Not enough information in the HydroFile in order to compute topographic residence time.")
            sys.exit()
250
        #
251
        #
252 253 254 255 256
        if FloodplainsFile is not None:
            nfp = Dataset(FloodplainsFile, "r")
            self.floodplains = gather(nfp.variables["floodplains"][jstr:jend,istr:iend].astype(np.float32), index, 0)
            #self.floodplainsdesc=nfp.variables["floodplains"].long_name
        elif "floodplains" in nf.variables.keys():
257
            self.floodplains = gather(nf.variables["floodplains"][jstr:jend,istr:iend].astype(np.float32), index, 0)
258 259 260
            self.floodplainsdesc=nf.variables["floodplains"].long_name
        else:
            self.floodplains = gather(np.zeros((jend-jstr,iend-istr)).astype(np.float32), index)
261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288

class HydroParameter :
    # Class to set, compute and manage the time constants of the routing scheme.
    def __init__(self, hydrogrid) :
        if hydrogrid.hd >= 0.5 :
            # Case of the Fekete/Vörösmarty hydrological data set
            self.stream_tcst = 0.24
            self.fast_tcst = 3.0
            self.slow_tcst = 25.0
            self.flood_tcst = 4.0
            self.swamp_tcst = 0.2
        elif hydrogrid.hd >= 0.016 :
            # Case for MERIT
            self.stream_tcst = 0.7
            self.fast_tcst = 1.0
            self.slow_tcst = 10.0
            self.flood_tcst = 4.0
            self.swamp_tcst = 0.2
        elif hydrogrid.hd >= 0.008 :
            # Case for HydroSEHDS
            self.stream_tcst = 0.01
            self.fast_tcst = 0.5
            self.slow_tcst = 7.0
            self.flood_tcst = 4.0
            self.swamp_tcst = 0.2
        else :
            print("For the resolution ",hydrogrid.hd," deg. we do not yet have a parameter set.")
            sys.exit