00001
00002 """
00003 Created on Mon Mar 12 17:41:54 2012
00004
00005 @author: sat kumar tomer
00006 @email: satkumartomer@gmail.com
00007 @website: www.ambhas.com
00008 """
00009
00010 from __future__ import division
00011 import numpy as np
00012 import xlrd
00013 from Scientific.IO import NetCDF as nc
00014 import datetime
00015 import matplotlib.pyplot as plt
00016 from BIP.Bayes.lhs import lhs
00017 from scipy import stats
00018
00019
00020 class RICHARDS_1D():
00021 """
00022 This is the main class of the RICHARDS_1D.
00023 This simulates the flow in unsaturated porus media
00024
00025 This will read the input data,
00026 do the processing
00027 and then write the output files
00028
00029 """
00030
00031 def __init__(self,input_file):
00032 """
00033 Input:
00034 input_file: the file which contains all the information
00035 including forcing and parameters.
00036 """
00037 self.input_file = input_file
00038
00039
00040 self._read_input()
00041
00042
00043 self.initialize()
00044
00045
00046 for t in range(self.max_t):
00047 self.t = t
00048
00049
00050 self._get_forcing()
00051
00052
00053 self._unsat()
00054
00055 self.nc_file.close()
00056
00057
00058
00059 def _read_input(self):
00060 """
00061 This checks if all the required input sheets are present in the xls file,
00062 read the data from input file, which can be used later in other functions
00063 """
00064
00065
00066 input_sheets = ['ind', 'forcing', 'initial_condition', 'units', 'temporal_info',
00067 'spatial_info', 'soil_hyd_par', 'output_par']
00068
00069
00070 self._check_sheets(input_sheets, self.input_file)
00071
00072
00073 self._read_ind()
00074
00075
00076 self._read_spatial()
00077
00078
00079 self._read_temporal()
00080
00081
00082 self._read_units()
00083
00084
00085 self._read_initial_condition()
00086
00087
00088 self._read_shp()
00089
00090
00091 self._read_forcing()
00092
00093
00094 self._read_ofile_name()
00095
00096
00097 output_message = 'Input data reading completed sucessfully'
00098 self._colored_output(output_message, 32)
00099
00100 def _check_sheets(self, check_sheets, check_file):
00101 """
00102 This functions check if all the sheets needed to model are present
00103 in check_file
00104
00105 """
00106
00107 foo = xlrd.open_workbook(check_file)
00108 check_sheet_names = foo.sheet_names()
00109
00110 for check_sheet in check_sheets:
00111 if check_sheet not in check_sheet_names:
00112 output_message = 'The sheet ' + check_sheet + ' is missing'
00113 self._colored_output(output_message,31)
00114
00115
00116 def _read_ind(self):
00117 """
00118 Read the ind sheet
00119 legend stores the information about the indices of other properties,
00120 which would be used by all other properties reading functions
00121 """
00122 book = xlrd.open_workbook(self.input_file)
00123 sheet = book.sheet_by_name('ind')
00124
00125 ind = {}
00126 for i in range(sheet.nrows-1):
00127 ind[str(sheet.cell_value(i+1,0))] = int(sheet.cell_value(i+1,1))
00128
00129 self.ind = ind
00130
00131 def _read_spatial(self):
00132 """
00133 Read the spatial info
00134 """
00135 book = xlrd.open_workbook(self.input_file)
00136 sheet = book.sheet_by_name('spatial_info')
00137
00138 j = self.ind['spatial_info']
00139 no_layer = int(sheet.cell_value(j,1))
00140 dz = sheet.cell_value(j,2)
00141
00142 self.no_layer = no_layer
00143 self.dz = dz
00144
00145 def _read_temporal(self):
00146 """
00147 Read the temporal info
00148 """
00149 book = xlrd.open_workbook(self.input_file)
00150 sheet = book.sheet_by_name('temporal_info')
00151
00152 j = self.ind['temporal_info']
00153 dt = sheet.cell_value(j,1)
00154 final_time = sheet.cell_value(j,2)
00155
00156 self.dt_flux = dt
00157 self.final_time = final_time
00158
00159
00160 def _read_units(self):
00161 """
00162 read the units of the forcing data
00163 """
00164 book = xlrd.open_workbook(self.input_file)
00165 sheet = book.sheet_by_name('units')
00166
00167 j = self.ind['units']
00168 forcing_units = {}
00169 for i in range(sheet.ncols-1):
00170 forcing_units[str(sheet.cell_value(0,i+1))] = str(sheet.cell_value(j,i+1))
00171 self.forcing_units = forcing_units
00172
00173 def _read_initial_condition(self):
00174 """
00175 read initial condition
00176 """
00177
00178 j = self.ind['initial_condition']
00179
00180 book = xlrd.open_workbook(self.input_file)
00181 sheet = book.sheet_by_name('initial_condition')
00182 theta_0 = sheet.cell_value(j,1)
00183 self.theta = np.tile(theta_0,self.no_layer)
00184
00185 def _read_shp(self):
00186 """
00187 read the soil hydraulic parameters
00188 """
00189
00190 j = self.ind['soil_hyd_par']
00191
00192 book = xlrd.open_workbook(self.input_file)
00193 sheet = book.sheet_by_name('soil_hyd_par')
00194 soil_par = {}
00195 soil_par['thetar'] = sheet.cell_value(j,1)
00196 soil_par['thetas'] = sheet.cell_value(j,2)
00197 soil_par['alpha'] = sheet.cell_value(j,3)
00198 soil_par['n'] = sheet.cell_value(j,4)
00199 soil_par['Ks'] = sheet.cell_value(j,5)
00200 soil_par['l'] = sheet.cell_value(j,6)
00201 soil_par['evap_0'] = sheet.cell_value(j,7)
00202 soil_par['evap_1'] = sheet.cell_value(j,8)
00203 soil_par['m'] = 1-1/soil_par['n']
00204 self.soil_par = soil_par
00205
00206
00207 def _read_forcing(self):
00208 """
00209 read the forcing data from xls file
00210 """
00211 book = xlrd.open_workbook(self.input_file)
00212 sheet = book.sheet_by_name('forcing')
00213
00214 data_len = sheet.nrows-1
00215 year = np.zeros(data_len)
00216 doy = np.zeros(data_len)
00217 rain = np.zeros(data_len)
00218 pet = np.zeros(data_len)
00219
00220 for i in xrange(data_len):
00221 year[i] = sheet.cell_value(i+1,0)
00222 doy[i] = sheet.cell_value(i+1,1)
00223 rain[i] = sheet.cell_value(i+1,2)
00224 pet[i] = sheet.cell_value(i+1,3)
00225
00226
00227 self.year = year
00228 self.doy = doy
00229
00230
00231 if self.forcing_units['rain'] == 'mm':
00232 self.rain = rain/1000.0
00233 elif self.forcing_units['rain'] == 'm':
00234 self.rain = rain
00235 else:
00236 raise ValueError("The units of rain should be either 'mm' or 'm' ")
00237
00238 if self.forcing_units['pet'] == 'mm':
00239 self.pet = pet/1000.0
00240 elif self.forcing_units['pet'] == 'm':
00241 self.pet = pet
00242 else:
00243 raise ValueError("The units of PET should be either 'mm' or 'm' ")
00244
00245
00246
00247 def _read_ofile_name(self):
00248 """
00249 read the forcing data from xls file
00250 """
00251 book = xlrd.open_workbook(self.input_file)
00252 sheet = book.sheet_by_name('output_par')
00253 j = self.ind['output_par']
00254 self.ofile_name = str(sheet.cell_value(j,1))
00255
00256
00257 def _colored_output(self, output_message, color):
00258 """
00259 This functions print the output_message in the color
00260 Input:
00261 output_messgae: the text you want to print
00262 color: the color in which you want to print text, it could be one of:
00263 30: Gray
00264 31: Red
00265 32: Green
00266 33: Yellow
00267 34: Blue
00268 35: Magneta
00269 66: Cyan
00270 37: White
00271 38: Crimson
00272 41: Highlighted Red
00273 42: Highlighted Green
00274 43: Highlighted Brown
00275 44: Highlighted Blue
00276 45: Highlighted Magenta
00277 46: Highlighted Cyan
00278 47: Highlighted Gray
00279 48: Highlighted Crimson
00280 Output:
00281 This returns None, but print the output in python shell
00282 """
00283
00284 print(("\033[31m" +output_message+ "\033[0m").replace('31',str(color)))
00285
00286 def _get_forcing(self):
00287 """
00288 this will give the forcing at time t
00289 forcing are given in terms of L/T
00290 """
00291 self.rain_cur = self.rain[self.t]/self.dt_flux
00292 self.pet_cur = self.pet[self.t]/self.dt_flux
00293
00294 self.cur_year = self.year[self.t]
00295 self.cur_doy = self.doy[self.t]
00296
00297
00298 def smcf(self, theta, thetar, thetas, alpha, m, n):
00299 """
00300 smcf: calculate the smc
00301 """
00302 Se = (theta-thetar)/(thetas-thetar)
00303 Se[Se<=0] = 0.0001
00304 Se[Se>=1] = 0.9999
00305 try:
00306 smc = alpha*(thetas-thetar)*m*n*pow(Se,1/m+1)*pow(pow(Se,-1/m)-1,m)
00307 except:
00308 print thetas,thetar,m,n,Se,
00309 return smc
00310
00311 def theta2psi(self,theta, thetar, thetas, m, n, alpha):
00312 """
00313 theta2psi: given the theta calculate the psi
00314 """
00315 Se = (theta-thetar)/(thetas-thetar)
00316 Se[Se<=0] = 0.0001
00317 Se[Se>=1] = 0.9999
00318 psi = -(1/alpha)*pow(pow(Se,-1/m)-1,1/n)
00319 psi[psi<-1e6] = -1e6
00320 return psi
00321
00322 def psi2theta(self,psi, thetar, thetas, alpha, m, n):
00323 """
00324 psi2theta: given the theta calculate the pressure head
00325 """
00326 if (psi>=0):
00327 theta = thetas
00328 elif psi<-1e6:
00329 theta = 1.01*thetar
00330 else:
00331 theta = thetar+(thetas-thetar)*pow(1+pow(abs(alpha*psi),n),-m)
00332 return theta
00333
00334 def theta2kr(self,theta, thetar, thetas, m, l, Ks):
00335 """
00336 theta2kr: given the theta, calculate the kr
00337 """
00338 Se = (theta-thetar)/(thetas-thetar)
00339 Se[Se<0] = 0.00001
00340 Se[Se>1] = 0.99999
00341 kr = Ks*(pow(Se,l))*pow(1-pow(1-pow(Se,1/m),m),2)
00342 kr[Se<0] = 0
00343 kr[Se>1] = Ks
00344
00345 return kr
00346
00347 def initialize(self):
00348 """
00349 this initializes all the required variables
00350 and open the netcdf file for writting
00351 """
00352 max_t = int(self.final_time/self.dt_flux)
00353
00354 self.max_t = max_t
00355 self.iter_dt = 1
00356
00357
00358 file = nc.NetCDFFile(self.ofile_name, 'w')
00359 setattr(file, 'title', 'output of the model ambhas.richards')
00360 now = datetime.datetime.now()
00361 setattr(file, 'description', 'The model was run at %s'%(now.ctime()))
00362 file.createDimension('depth', self.no_layer)
00363 file.createDimension('time', self.max_t+1)
00364
00365
00366 varDims = 'depth',
00367 depth = file.createVariable('depth', 'd', varDims)
00368 depth.units = 'm'
00369 depth[:] = np.tile(self.dz,self.no_layer).cumsum()-self.dz/2
00370
00371
00372 varDims = 'time',
00373 self.nc_year = file.createVariable('year', 'd', varDims)
00374 self.nc_doy = file.createVariable('doy', 'd', varDims)
00375
00376
00377 varDims = 'depth','time'
00378 self.nc_sm = file.createVariable('sm','d', varDims)
00379 self.nc_sm.units = 'v/v'
00380 self.nc_sm[:,0] = self.theta
00381
00382
00383 varDims = 'time',
00384 self.nc_aet = file.createVariable('aet','d',varDims)
00385 self.nc_aet.units = 'mm'
00386 self.nc_recharge = file.createVariable('recharge','d',varDims)
00387 self.nc_recharge.units = 'mm'
00388
00389
00390 setattr(file, 'thetar', self.soil_par['thetar'])
00391 setattr(file, 'thetas', self.soil_par['thetas'])
00392 setattr(file, 'alpha', self.soil_par['alpha'])
00393 setattr(file, 'n', self.soil_par['n'])
00394 setattr(file, 'Ks', self.soil_par['Ks'])
00395 setattr(file, 'l', self.soil_par['l'])
00396
00397 self.nc_file = file
00398
00399
00400 def _unsat(self):
00401 """
00402 top boundary: atmoshpheric
00403 bottom boundary: gravity drainage
00404 """
00405
00406 thetar = self.soil_par['thetar']
00407 thetas = self.soil_par['thetas']
00408 alpha = self.soil_par['alpha']
00409 n = self.soil_par['n']
00410 m = self.soil_par['m']
00411 l = self.soil_par['l']
00412 Ks = self.soil_par['Ks']
00413 nz = self.no_layer
00414
00415 theta = 1.0*self.theta
00416
00417
00418
00419 iter_dt = max(24,int(np.ceil(self.rain_cur*self.dt_flux*1000/0.15)))
00420 self.iter_dt = int(max(iter_dt,0.75*self.iter_dt))
00421
00422
00423
00424
00425
00426 recharge_day = 0
00427 aet_day = 0
00428
00429
00430 for i in range(self.iter_dt):
00431 dt = self.dt_flux/self.iter_dt
00432
00433
00434 smi = (self.theta[0]-self.soil_par['evap_0'])/(self.soil_par['evap_1']-self.soil_par['evap_0'])
00435 if smi<0: smi=0
00436 if smi>1: smi=1
00437 aet = smi*self.pet_cur
00438 Bvalue = self.rain_cur-aet
00439
00440 K = self.theta2kr(theta,thetar,thetas,m,l,Ks)
00441 smc = self.smcf(theta,thetar,thetas,alpha,m,n)
00442 psi = self.theta2psi(theta,thetar,thetas,m,n,alpha)
00443
00444
00445 Kmid = np.empty(nz+1)
00446 Kmid[0] = 0
00447 for i in range(1,nz):
00448 Kmid[i] = 0.5*(K[i]+K[i-1])
00449 Kmid[nz] = K[nz-1]
00450
00451
00452 A = np.empty(nz)
00453 B = np.empty(nz)
00454 C = np.empty(nz)
00455 D = np.empty(nz)
00456 dz = self.dz
00457 dz2 = dz**2
00458
00459 for i in range(nz):
00460 A[i] = -(Kmid[i]/dz2)
00461 B[i] = smc[i]/dt+(Kmid[i+1]+Kmid[i])/dz2
00462 C[i] = A[i]
00463 D[i] = smc[i]*psi[i]/dt-(Kmid[i+1]-Kmid[i])/dz
00464
00465 i = 0
00466 A[0] = 0
00467 B[0] = smc[i]/dt+(Kmid[1])/dz2
00468 D[0] = smc[i]*psi[i]/dt+(Bvalue-Kmid[1])/dz
00469
00470
00471 B[nz-1] = smc[nz-1]/dt+(Kmid[nz])/dz2
00472 C[nz-1] = 0
00473 D[nz-1] = smc[nz-1]*psi[nz-1]/dt-(Kmid[nz]-Kmid[nz-1])/dz
00474
00475
00476 beta = np.empty(nz)
00477 gamma = np.empty(nz)
00478 u = np.empty(nz)
00479 beta[0] = B[0]
00480 gamma[0] = D[0]/beta[0]
00481 for i in range(1,nz):
00482 beta[i] = B[i]-(A[i]*C[i-1])/(beta[i-1])
00483 gamma[i] = (D[i]-A[i]*gamma[i-1])/(beta[i])
00484
00485 u[nz-1] = gamma[nz-1]
00486 for i in range(nz-2,-1,-1):
00487 u[i] = gamma[i]-(C[i]*u[i+1])/beta[i]
00488
00489
00490 J = np.empty(nz+1)
00491 for i in range(1,nz):
00492 J[i] = Kmid[i]*(1-(u[i]-u[i-1])/dz)
00493 J[0] = Bvalue
00494 J[nz] = Kmid[nz]
00495
00496 J[nz] = J[nz]
00497
00498 flux = np.diff(J)*dt/dz
00499 theta = theta - flux
00500
00501 if theta[0]>thetas:
00502 theta[theta>thetas] = 0.99*thetas
00503
00504 aet_day += aet*dt
00505 recharge_day += J[nz]*dt
00506
00507 self.theta = theta
00508
00509
00510 self.nc_year[self.t] = (self.cur_year)
00511 self.nc_doy[self.t] = (self.cur_doy)
00512 self.nc_sm[:,self.t+1] = theta
00513 self.nc_recharge[self.t] = recharge_day
00514 self.nc_aet[self.t] = aet_day
00515
00516
00517 if self.t == int(0.25*self.max_t):
00518 output_message = '25 % completed'
00519 self._colored_output(output_message, 32)
00520
00521 elif self.t == int(0.5*self.max_t):
00522 output_message = '50 % completed'
00523 self._colored_output(output_message, 32)
00524
00525 elif self.t == int(0.75*self.max_t):
00526 output_message = '75 % completed'
00527 self._colored_output(output_message, 32)
00528
00529 elif self.t == self.max_t-1:
00530 output_message = '100 % completed'
00531 self._colored_output(output_message, 32)
00532
00533
00534
00535
00536 class RICHARDS_1D_ENKF(RICHARDS_1D):
00537 """
00538 This is the main class of the Ensemble Kalman Filter (EnKF)
00539 coupled with the one dimensional unsaturated model based on the
00540 RICHARDS equation. The model is given in the class RICHARDS_1D.
00541
00542 This will read the input data,
00543 do the processing
00544 and then write the output files
00545
00546 """
00547
00548 def __init__(self,input_file):
00549 """
00550 Input:
00551 input_file: the file which contains all the information
00552 including forcing and parameters.
00553 """
00554 self.input_file = input_file
00555 self.n_ens = 10
00556
00557 self._read_input()
00558
00559
00560 self.initialize()
00561
00562
00563 for t in range(self.max_t):
00564 self.t = t
00565
00566
00567 self._get_forcing()
00568
00569
00570 self._perturb_soil_par_ens()
00571
00572
00573 for ens in range(self.n_ens):
00574 self.ens = ens
00575
00576 self._unsat_ens()
00577
00578
00579 self._enkf_par_depth()
00580
00581
00582 self._write_output()
00583
00584
00585 self.nc_file.close()
00586
00587
00588 def _enkf(self):
00589 """
00590 ensemble kalman filter
00591 """
00592
00593 x = self.theta_ens
00594 x_bar = np.tile(x.mean(axis=0),(10,1))
00595 x_x_bar = x-x_bar
00596 cov_xx = np.dot(x_x_bar.T,x_x_bar)
00597
00598
00599
00600 e = np.zeros((self.n_ens, self.no_layer))
00601 e[:,0] = self.theta_ens[:,0] - self.meas_ssm[self.t]
00602 e = e + 0.02*np.random.normal(size=(self.n_ens,self.no_layer))
00603 cov_ee = np.dot(e.T, e)
00604
00605
00606 K = np.dot(cov_xx,np.linalg.pinv(cov_xx+cov_ee))
00607
00608
00609 d = np.zeros(self.no_layer)
00610 usm = np.zeros((self.n_ens,self.no_layer))
00611
00612 for ens in range(self.n_ens):
00613 d[0] = self.meas_ssm[self.t] - x[ens,0]
00614 usm[ens,:] = (x[ens,:] + np.dot(K,d))
00615
00616 self.theta_ens = usm
00617
00618 def _enkf_par(self):
00619 """
00620 ensemble kalman filter
00621 """
00622
00623
00624 x = self.theta_ens + 0.0*np.random.normal(size=self.no_layer)
00625 thetar = self.soil_par_ens['thetar']
00626 thetas = self.soil_par_ens['thetas']
00627 alpha = self.soil_par_ens['alpha']
00628 n = self.soil_par_ens['n']
00629 Ks = self.soil_par_ens['Ks']
00630 l = self.soil_par_ens['l']
00631 soil_par = (np.vstack([thetar, thetas, alpha, n, Ks, l])).T
00632 X = np.hstack([x, soil_par])
00633
00634
00635 X_bar = np.tile(X.mean(axis=0),(10,1))
00636 X_X_bar = X-X_bar
00637 cov_XX = np.dot(X_X_bar.T,X_X_bar) + 1e-6*np.eye(self.no_layer+6)
00638 cov_XX = 0.5*(cov_XX + cov_XX.T)
00639
00640
00641
00642 e = np.zeros((self.n_ens, self.no_layer+6))
00643 e[:,0] = self.meas_ssm[self.t-1] - self.theta_ens[:,0].mean()
00644
00645 v = 0.025*np.random.normal(size=(self.n_ens,self.no_layer+6))
00646 v = v-np.tile(v.mean(axis=0),(self.n_ens,1))
00647 ev = e + v
00648 cov_ee = np.dot(e.T, e) + 1e-6*np.eye(self.no_layer+6)
00649 cov_ee = 0.5*(cov_ee + cov_ee.T)
00650
00651
00652 K = np.dot(cov_XX, np.linalg.pinv(cov_XX+cov_ee))
00653
00654
00655 v = np.random.normal(size=(self.n_ens,))
00656 v = v-v.mean()
00657 e[:,0] = e[:,0]+0.0005*v
00658 K = 0.5*(K + K.T)
00659 usm_par = X + np.dot(K,e.T).T
00660 self.usm_par = usm_par
00661
00662
00663
00664
00665 theta_ens = usm_par[:,:self.no_layer]
00666 theta_ens[theta_ens<0] = 0
00667 theta_ens[theta_ens>1] = 1
00668 v = 0.01*np.random.normal(size=(theta_ens.shape))
00669 v = v-np.tile(v.mean(axis=0),(self.n_ens,1))
00670 self.theta_ens = theta_ens + v
00671
00672 thetar = usm_par[:,self.no_layer+0]
00673 thetas = usm_par[:,self.no_layer+1]
00674 alpha = usm_par[:,self.no_layer+2]
00675 n = usm_par[:,self.no_layer+3]
00676 Ks = usm_par[:,self.no_layer+4]
00677 l = usm_par[:,self.no_layer+5]
00678
00679
00680
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692
00693
00694 thetar[thetar>self.thetar_max] = self.thetar_max
00695 thetas[thetas>self.thetas_max] = self.thetas_max
00696 alpha[alpha>self.alpha_max] = self.alpha_max
00697 n[n>self.n_max] = self.n_max
00698 Ks[Ks>self.Ks_max] = self.Ks_max
00699 l[l>self.l_max] = self.l_max
00700
00701 thetar[thetar<self.thetar_min] = self.thetar_min
00702 thetas[thetas<self.thetas_min] = self.thetas_min
00703 alpha[alpha<self.alpha_min] = self.alpha_min
00704 n[n<self.n_min] = self.n_min
00705
00706 l[l<self.l_min] = self.l_min
00707
00708
00709 self.soil_par_ens['thetar'] = thetar
00710 self.soil_par_ens['thetas'] = thetas
00711 self.soil_par_ens['alpha'] = alpha
00712 self.soil_par_ens['n'] = n
00713 self.soil_par_ens['Ks'] = Ks
00714 self.soil_par_ens['l'] = l
00715
00716 self.K = K
00717 self.cov_ee = cov_ee
00718 self.cov_XX = cov_XX
00719
00720 def _enkf_par_depth(self):
00721 """
00722 ensemble kalman filter
00723 """
00724
00725
00726 x = self.theta_ens + 0.0*np.random.normal(size=self.no_layer)
00727 thetar = self.soil_par_ens['thetar']
00728 thetas = self.soil_par_ens['thetas']
00729 alpha = self.soil_par_ens['alpha']
00730 n = self.soil_par_ens['n']
00731 Ks = self.soil_par_ens['Ks']
00732 l = self.soil_par_ens['l']
00733 soil_par = (np.vstack([thetar, thetas, alpha, n, Ks, l])).T
00734 X = np.hstack([x, soil_par])
00735
00736
00737 X_bar = np.tile(X.mean(axis=0),(10,1))
00738 X_X_bar = X-X_bar
00739 cov_XX = np.dot(X_X_bar.T,X_X_bar) + 1e-6*np.eye(self.no_layer+6)
00740 cov_XX = 0.5*(cov_XX + cov_XX.T)
00741
00742
00743
00744 e = np.zeros((self.n_ens, self.no_layer+6))
00745 ev = np.zeros((self.n_ens, self.no_layer+6))
00746 obs = np.zeros(self.no_layer)
00747 obs = self.a+self.b*self.meas_ssm[self.t-1]
00748
00749
00750 for i in range(20):
00751 e[:,i] = obs[i] - self.theta_ens[:,i].mean()
00752 self.e = e
00753 for i in range(self.no_layer+6):
00754 if i<self.no_layer:
00755 ev[:,i] = e[:,i] + 0.005*(i**2+1)*np.random.normal(size=(self.n_ens))
00756 else:
00757 ev[:,i] = e[:,i] + 0.000025*np.random.normal(size=(self.n_ens))
00758
00759
00760
00761
00762 cov_ee = np.dot(ev.T, ev) + 1e-6*np.eye(self.no_layer+6)
00763 cov_ee = 0.5*(cov_ee + cov_ee.T)
00764
00765
00766 K = np.dot(cov_XX, np.linalg.pinv(cov_XX+cov_ee))
00767
00768
00769 v = np.random.normal(size=(self.n_ens,))
00770 v = v-v.mean()
00771 e[:,0] = e[:,0]+0.0005*v
00772 K = 0.5*(K + K.T)
00773 usm_par = X + np.dot(K,e.T).T
00774 self.usm_par = usm_par
00775
00776
00777
00778
00779 theta_ens = usm_par[:,:self.no_layer]
00780 theta_ens[theta_ens<0] = 0
00781 theta_ens[theta_ens>1] = 1
00782 v = 0.01*np.random.normal(size=(theta_ens.shape))
00783 v = v-np.tile(v.mean(axis=0),(self.n_ens,1))
00784 self.theta_ens = theta_ens + v
00785
00786 thetar = usm_par[:,self.no_layer+0]
00787 thetas = usm_par[:,self.no_layer+1]
00788 alpha = usm_par[:,self.no_layer+2]
00789 n = usm_par[:,self.no_layer+3]
00790 Ks = usm_par[:,self.no_layer+4]
00791 l = usm_par[:,self.no_layer+5]
00792
00793
00794
00795
00796
00797
00798
00799
00800
00801
00802
00803
00804
00805
00806
00807
00808 thetar[thetar>self.thetar_max] = self.thetar_max
00809 thetas[thetas>self.thetas_max] = self.thetas_max
00810 alpha[alpha>self.alpha_max] = self.alpha_max
00811 n[n>self.n_max] = self.n_max
00812 Ks[Ks>self.Ks_max] = self.Ks_max
00813 l[l>self.l_max] = self.l_max
00814
00815 thetar[thetar<self.thetar_min] = self.thetar_min
00816 thetas[thetas<self.thetas_min] = self.thetas_min
00817 alpha[alpha<self.alpha_min] = self.alpha_min
00818 n[n<self.n_min] = self.n_min
00819 Ks[Ks<self.Ks_min] = self.Ks_min
00820 l[l<self.l_min] = self.l_min
00821
00822
00823 self.soil_par_ens['thetar'] = thetar
00824 self.soil_par_ens['thetas'] = thetas
00825 self.soil_par_ens['alpha'] = alpha
00826 self.soil_par_ens['n'] = n
00827 self.soil_par_ens['Ks'] = Ks
00828 self.soil_par_ens['l'] = l
00829
00830 self.K = K
00831 self.cov_ee = cov_ee
00832 self.cov_XX = cov_XX
00833
00834 def initialize(self):
00835 """
00836 this initializes all the required variables
00837 and open the netcdf file for writting
00838 also generates the initial ensemble of the soil hydraulic parameters
00839 """
00840 max_t = int(self.final_time/self.dt_flux)
00841
00842 self.max_t = max_t
00843 self.iter_dt = 1
00844
00845
00846 file = nc.NetCDFFile(self.ofile_name, 'w')
00847 setattr(file, 'title', 'output of the model ambhas.richards')
00848 now = datetime.datetime.now()
00849 setattr(file, 'description', 'The model was run at %s'%(now.ctime()))
00850 file.createDimension('depth', self.no_layer)
00851 file.createDimension('time', self.max_t+1)
00852 file.createDimension('ensemble', self.n_ens)
00853
00854
00855 varDims = 'depth',
00856 depth = file.createVariable('depth', 'd', varDims)
00857 depth.units = 'm'
00858 depth[:] = np.tile(self.dz,self.no_layer).cumsum()-self.dz/2
00859
00860
00861 varDims = 'time',
00862 self.nc_year = file.createVariable('year', 'd', varDims)
00863 self.nc_doy = file.createVariable('doy', 'd', varDims)
00864
00865
00866 varDims = 'ensemble', 'depth', 'time'
00867 self.nc_sm = file.createVariable('sm','d', varDims)
00868 self.nc_sm.units = 'v/v'
00869 self.nc_sm[:,:,0] = self.theta_ens
00870
00871
00872 varDims = 'time',
00873 self.nc_aet = file.createVariable('aet','d',varDims)
00874 self.nc_aet.units = 'mm'
00875 self.nc_recharge = file.createVariable('recharge','d',varDims)
00876 self.nc_recharge.units = 'mm'
00877
00878
00879 varDims = 'ensemble','time'
00880 self.nc_thetar = file.createVariable('thetar','d',varDims)
00881 self.nc_thetar.units = 'v/v'
00882 self.nc_thetas = file.createVariable('thetas','d',varDims)
00883 self.nc_thetas.units = 'v/v'
00884 self.nc_alpha = file.createVariable('alpha','d',varDims)
00885 self.nc_alpha.units = '1/m'
00886 self.nc_n = file.createVariable('n','d',varDims)
00887 self.nc_n.units = '-'
00888 self.nc_Ks = file.createVariable('Ks','d',varDims)
00889 self.nc_Ks.units = 'm/s'
00890 self.nc_l = file.createVariable('l','d',varDims)
00891 self.nc_l.units = '-'
00892
00893 self.nc_file = file
00894
00895
00896 self._generate_soil_par_ens()
00897
00898 def _read_input(self):
00899 """
00900 This checks if all the required input sheets are present in the xls file,
00901 read the data from input file, which can be used later in other functions
00902 """
00903
00904
00905 input_sheets = ['ind', 'forcing', 'initial_condition', 'units', 'temporal_info',
00906 'spatial_info', 'soil_hyd_par_ens', 'output_par']
00907
00908
00909 self._check_sheets(input_sheets, self.input_file)
00910
00911
00912 self._read_ind()
00913
00914
00915 self._read_spatial()
00916
00917
00918 self._read_temporal()
00919
00920
00921 self._read_units()
00922
00923
00924 self._read_initial_condition()
00925
00926
00927 self._read_shp_ens()
00928
00929
00930 self._read_forcing()
00931
00932
00933 self._read_ofile_name()
00934
00935
00936 self._read_measured()
00937
00938
00939 self._read_ab()
00940
00941
00942 output_message = 'Input data reading completed sucessfully'
00943 self._colored_output(output_message, 32)
00944
00945 def _read_ab(self):
00946 """
00947 read the intercept and slope of the relationship of the surface soil
00948 moisture with the profile soil moisture
00949 """
00950 book = xlrd.open_workbook(self.input_file)
00951 sheet = book.sheet_by_name('ab')
00952
00953 data_len = sheet.nrows-1
00954 a = np.zeros(data_len)
00955 b = np.zeros(data_len)
00956
00957 for i in xrange(data_len):
00958 a[i] = sheet.cell_value(i+1,1)
00959 b[i] = sheet.cell_value(i+1,2)
00960
00961 self.a = a
00962 self.b = b
00963
00964
00965 def _read_measured(self):
00966 """
00967 read the measured surface soil moisture (ssm) data
00968 """
00969 book = xlrd.open_workbook(self.input_file)
00970 sheet = book.sheet_by_name('forcing')
00971
00972 data_len = sheet.nrows-1
00973 meas_ssm = np.zeros(data_len)
00974 j = self.ind['meas_sm']
00975
00976 for i in xrange(data_len):
00977 meas_ssm[i] = sheet.cell_value(i+1,3+j)
00978
00979 self.meas_ssm = meas_ssm
00980
00981 def _read_initial_condition(self):
00982 """
00983 read initial condition
00984 """
00985
00986 j = self.ind['initial_condition']
00987
00988 book = xlrd.open_workbook(self.input_file)
00989 sheet = book.sheet_by_name('initial_condition')
00990 theta_0 = sheet.cell_value(j,1)
00991 self.theta_ens = theta_0 + 0.05*np.random.normal(size=(self.n_ens,self.no_layer))
00992
00993 def _read_shp_ens(self):
00994 """
00995 read the information about the ensemble of the soil hydraulic parameters
00996 the information being read is the min, max, best estimate and its uncertainty
00997
00998 """
00999
01000 j = self.ind['soil_hyd_par_ens']
01001
01002 book = xlrd.open_workbook(self.input_file)
01003 sheet = book.sheet_by_name('soil_hyd_par_ens')
01004 self.thetar_min, self.thetar_max = sheet.cell_value(j+1,1), sheet.cell_value(j+1,7)
01005 self.thetas_min, self.thetas_max = sheet.cell_value(j+1,2), sheet.cell_value(j+1,8)
01006 self.alpha_min, self.alpha_max = sheet.cell_value(j+1,3), sheet.cell_value(j+1,9)
01007 self.n_min, self.n_max = sheet.cell_value(j+1,4), sheet.cell_value(j+1,10)
01008 self.Ks_min, self.Ks_max = sheet.cell_value(j+1,5), sheet.cell_value(j+1,11)
01009 self.l_min, self.l_max = sheet.cell_value(j+1,6), sheet.cell_value(j+1,12)
01010
01011 shp_ens = {}
01012 shp_ens['thetar'] = sheet.cell_value(j+1,13), sheet.cell_value(j+1,19)
01013 shp_ens['thetas'] = sheet.cell_value(j+1,14), sheet.cell_value(j+1,20)
01014 shp_ens['alpha'] = sheet.cell_value(j+1,15), sheet.cell_value(j+1,21)
01015 shp_ens['n'] = sheet.cell_value(j+1,16), sheet.cell_value(j+1,22)
01016 shp_ens['Ks'] = sheet.cell_value(j+1,17), sheet.cell_value(j+1,23)
01017 shp_ens['l'] = sheet.cell_value(j+1,18), sheet.cell_value(j+1,24)
01018
01019 self.shp_ens = shp_ens
01020
01021
01022 def _generate_soil_par_ens(self):
01023 """
01024 this uses the LHS to generate the ensemble of the parameters
01025
01026 this also computes the perturbation needed to perturb the parameters
01027 which is done in another function
01028 """
01029
01030
01031 v = np.random.normal(size=(self.n_ens,6))
01032 v = v-np.tile(v.mean(axis=0),(self.n_ens,1))
01033 thetar = self.shp_ens['thetar'][0] + v[:,0]*self.shp_ens['thetar'][1]
01034 thetas = self.shp_ens['thetas'][0] + v[:,1]*self.shp_ens['thetas'][1]
01035 alpha = self.shp_ens['alpha'][0] + v[:,2]*self.shp_ens['alpha'][1]
01036 n = self.shp_ens['n'][0] + v[:,3]*self.shp_ens['n'][1]
01037 Ks = self.shp_ens['Ks'][0] + v[:,4]*self.shp_ens['Ks'][1]
01038 l = self.shp_ens['l'][0] + v[:,5]*self.shp_ens['l'][1]
01039
01040
01041 thetar[thetar>self.thetar_max] = self.thetar_max
01042 thetas[thetas>self.thetas_max] = self.thetas_max
01043 alpha[alpha>self.alpha_max] = self.alpha_max
01044 n[n>self.n_max] = self.n_max
01045 Ks[Ks>self.Ks_max] = self.Ks_max
01046 l[l>self.l_max] = self.l_max
01047
01048 thetar[thetar<self.thetar_min] = self.thetar_min
01049 thetas[thetas<self.thetas_min] = self.thetas_min
01050 alpha[alpha<self.alpha_min] = self.alpha_min
01051 n[n<self.n_min] = self.n_min
01052 Ks[Ks<self.Ks_min] = self.Ks_min
01053 l[l<self.l_min] = self.l_min
01054
01055 soil_par_ens = {}
01056 soil_par_ens['thetar'] = thetar
01057 soil_par_ens['thetas'] = thetas
01058 soil_par_ens['alpha'] = alpha
01059 soil_par_ens['n'] = n
01060 soil_par_ens['Ks'] = Ks
01061 soil_par_ens['l'] = l
01062
01063 self.soil_par_ens = soil_par_ens
01064
01065
01066 soil_pert = {}
01067 soil_pert['thetar'] = (self.thetar_max - self.thetar_min)*0.1/100.0
01068 soil_pert['thetas'] = (self.thetas_max - self.thetas_min)*0.1/100.0
01069 soil_pert['alpha'] = (self.alpha_max - self.alpha_min)*0.1/100.0
01070 soil_pert['n'] = (self.n_max - self.n_min)*0.1/100.0
01071 soil_pert['Ks'] = (self.Ks_max - self.Ks_min)*0.1/100.0
01072 soil_pert['l'] = (self.l_max - self.l_min)*0.1/100.0
01073
01074 self.soil_pert = soil_pert
01075
01076 def _perturb_soil_par_ens(self):
01077 """
01078 this functions perturb the soil hydraulic parameters
01079 using the gaussian random variables
01080 """
01081 v = np.random.normal(size=(self.n_ens,6))
01082 v = v-np.tile(v.mean(axis=0),(self.n_ens,1))
01083 thetar = self.soil_par_ens['thetar']+self.soil_pert['thetar']*v[:,0]
01084 thetas = self.soil_par_ens['thetas']+self.soil_pert['thetas']*v[:,1]
01085 alpha = self.soil_par_ens['alpha']+self.soil_pert['alpha']*v[:,2]
01086 n = self.soil_par_ens['n']+self.soil_pert['n']*v[:,3]
01087 Ks = self.soil_par_ens['Ks']+self.soil_pert['Ks']*(v[:,4]-v[:,4].mean())
01088 l = self.soil_par_ens['l']+self.soil_pert['l']*v[:,5]
01089
01090
01091
01092 thetar[thetar>self.thetar_max] = self.thetar_max
01093 thetas[thetas>self.thetas_max] = self.thetas_max
01094 alpha[alpha>self.alpha_max] = self.alpha_max
01095 n[n>self.n_max] = self.n_max
01096 Ks[Ks>self.Ks_max] = self.Ks_max
01097 l[l>self.l_max] = self.l_max
01098
01099 thetar[thetar<self.thetar_min] = self.thetar_min
01100 thetas[thetas<self.thetas_min] = self.thetas_min
01101 alpha[alpha<self.alpha_min] = self.alpha_min
01102 n[n<self.n_min] = self.n_min
01103 Ks[Ks<self.Ks_min] = self.Ks_min
01104 l[l<self.l_min] = self.l_min
01105
01106 soil_par_ens = {}
01107 soil_par_ens['thetar'] = thetar
01108 soil_par_ens['thetas'] = thetas
01109 soil_par_ens['alpha'] = alpha
01110 soil_par_ens['n'] = n
01111 soil_par_ens['Ks'] = Ks
01112 soil_par_ens['l'] = l
01113
01114 self.soil_par_ens = soil_par_ens
01115
01116
01117 def _unsat_ens(self):
01118 """
01119 top boundary: atmoshpheric
01120 bottom boundary: gravity drainage
01121 """
01122 ens = self.ens
01123 thetar = self.soil_par_ens['thetar'][ens]
01124 thetas = self.soil_par_ens['thetas'][ens]
01125 alpha = self.soil_par_ens['alpha'][ens]
01126 n = self.soil_par_ens['n'][ens]
01127 l = self.soil_par_ens['l'][ens]
01128 Ks = self.soil_par_ens['Ks'][ens]
01129 m = 1-1/n
01130 nz = self.no_layer
01131 evap_0 = thetar+0.25*(thetas-thetar)
01132 evap_1 = thetar+0.75*(thetas-thetar)
01133
01134
01135 theta = 1.0*self.theta_ens[ens]
01136
01137
01138
01139 iter_dt = max(24,int(np.ceil(self.rain_cur*self.dt_flux*1000/0.15)))
01140 self.iter_dt = int(max(iter_dt,0.75*self.iter_dt))
01141
01142
01143
01144
01145
01146 recharge_day = 0
01147 aet_day = 0
01148
01149
01150 for i in range(self.iter_dt):
01151 dt = self.dt_flux/self.iter_dt
01152
01153
01154 smi = (theta[0]-evap_0)/(evap_1-evap_0)
01155 if smi<0: smi=0
01156 if smi>1: smi=1
01157 aet = smi*self.pet_cur
01158 Bvalue = self.rain_cur-aet
01159
01160 K = self.theta2kr(theta,thetar,thetas,m,l,Ks)
01161 smc = self.smcf(theta,thetar,thetas,alpha,m,n)
01162 psi = self.theta2psi(theta,thetar,thetas,m,n,alpha)
01163
01164
01165
01166 Kmid = np.empty(nz+1)
01167 Kmid[0] = 0
01168 for i in range(1,nz):
01169 Kmid[i] = 0.5*(K[i]+K[i-1])
01170 Kmid[nz] = K[nz-1]
01171
01172
01173 A = np.empty(nz)
01174 B = np.empty(nz)
01175 C = np.empty(nz)
01176 D = np.empty(nz)
01177 dz = self.dz
01178 dz2 = dz**2
01179
01180 for i in range(nz):
01181 A[i] = -(Kmid[i]/dz2)
01182 B[i] = smc[i]/dt+(Kmid[i+1]+Kmid[i])/dz2
01183 C[i] = A[i]
01184 D[i] = smc[i]*psi[i]/dt-(Kmid[i+1]-Kmid[i])/dz
01185
01186 i = 0
01187 A[0] = 0
01188 B[0] = smc[i]/dt+(Kmid[1])/dz2
01189 D[0] = smc[i]*psi[i]/dt+(Bvalue-Kmid[1])/dz
01190
01191
01192 B[nz-1] = smc[nz-1]/dt+(Kmid[nz])/dz2
01193 C[nz-1] = 0
01194 D[nz-1] = smc[nz-1]*psi[nz-1]/dt-(Kmid[nz]-Kmid[nz-1])/dz
01195
01196
01197 beta = np.empty(nz)
01198 gamma = np.empty(nz)
01199 u = np.empty(nz)
01200 beta[0] = B[0]
01201 gamma[0] = D[0]/beta[0]
01202 for i in range(1,nz):
01203 beta[i] = B[i]-(A[i]*C[i-1])/(beta[i-1])
01204 gamma[i] = (D[i]-A[i]*gamma[i-1])/(beta[i])
01205
01206 u[nz-1] = gamma[nz-1]
01207 for i in range(nz-2,-1,-1):
01208 u[i] = gamma[i]-(C[i]*u[i+1])/beta[i]
01209
01210
01211 J = np.empty(nz+1)
01212 for i in range(1,nz):
01213 J[i] = Kmid[i]*(1-(u[i]-u[i-1])/dz)
01214 J[0] = Bvalue
01215 J[nz] = Kmid[nz]
01216
01217 J[nz] = J[nz]
01218
01219 flux = np.diff(J)*dt/dz
01220 theta = theta - flux
01221
01222 theta[theta>thetas] = 0.99*thetas
01223 theta[theta<thetar] = 1.01*thetar
01224
01225 aet_day += aet*dt
01226 recharge_day += J[nz]*dt
01227
01228
01229 self.theta_ens[ens] = theta
01230
01231 def _write_output(self):
01232 """
01233 this functions writes the output at each time step
01234 """
01235
01236 self.nc_year[self.t] = (self.cur_year)
01237 self.nc_doy[self.t] = (self.cur_doy)
01238 self.nc_sm[:,:,self.t+1] = self.theta_ens
01239
01240
01241 self.nc_thetar[:,self.t] = self.soil_par_ens['thetar']
01242 self.nc_thetas[:,self.t] = self.soil_par_ens['thetas']
01243 self.nc_alpha[:,self.t] = self.soil_par_ens['alpha']
01244 self.nc_n[:,self.t] = self.soil_par_ens['n']
01245 self.nc_Ks[:,self.t] = self.soil_par_ens['Ks']
01246 self.nc_l[:,self.t] = self.soil_par_ens['l']
01247
01248
01249
01250 class RICHARDS_1D_GLUE(RICHARDS_1D):
01251 """
01252 This is the main class of the Ensemble Kalman Filter (EnKF)
01253 coupled with the one dimensional unsaturated model based on the
01254 RICHARDS equation. The model is given in the class RICHARDS_1D.
01255
01256 This will read the input data,
01257 do the processing
01258 and then write the output files
01259
01260 """
01261
01262 def __init__(self,input_file):
01263 """
01264 Input:
01265 input_file: the file which contains all the information
01266 including forcing and parameters.
01267 """
01268 self.input_file = input_file
01269 self.n_ens = 1000
01270
01271
01272 self._read_input()
01273
01274
01275 self.initialize()
01276
01277
01278 for ens in range(self.n_ens):
01279 self.ens = ens
01280
01281 self._shp_cur()
01282
01283
01284
01285 self._read_initial_condition()
01286
01287
01288 for t in range(self.max_t):
01289 self.t = t
01290
01291
01292 self._get_forcing()
01293
01294
01295 self._unsat()
01296
01297 output_message = '%d out of %d ensemble completed'%(ens,self.n_ens)
01298 self._colored_output(output_message, 41)
01299
01300 self.nc_file.close()
01301
01302 def _read_input(self):
01303 """
01304 This checks if all the required input sheets are present in the xls file,
01305 read the data from input file, which can be used later in other functions
01306 """
01307
01308
01309 input_sheets = ['ind', 'forcing', 'initial_condition', 'units', 'temporal_info',
01310 'spatial_info', 'soil_hyd_par_ens', 'output_par']
01311
01312
01313 self._check_sheets(input_sheets, self.input_file)
01314
01315
01316 self._read_ind()
01317
01318
01319 self._read_spatial()
01320
01321
01322 self._read_temporal()
01323
01324
01325 self._read_units()
01326
01327
01328 self._read_initial_condition()
01329
01330
01331 self._read_shp_ens()
01332
01333
01334 self._read_forcing()
01335
01336
01337 self._read_ofile_name()
01338
01339
01340
01341
01342 output_message = 'Input data reading completed sucessfully'
01343 self._colored_output(output_message, 32)
01344
01345
01346 def _read_shp_ens(self):
01347 """
01348 read the information about the ensemble of the soil hydraulic parameters
01349 the information being read is the min and max
01350
01351 this also uses LHS to generate the ensemble of the parameters
01352 """
01353
01354 j = self.ind['soil_hyd_par_ens']
01355
01356 book = xlrd.open_workbook(self.input_file)
01357 sheet = book.sheet_by_name('soil_hyd_par_ens')
01358 thetar_min, thetar_max = sheet.cell_value(j+1,1), sheet.cell_value(j+1,7)
01359 thetas_min, thetas_max = sheet.cell_value(j+1,2), sheet.cell_value(j+1,8)
01360 alpha_min, alpha_max = sheet.cell_value(j+1,3), sheet.cell_value(j+1,9)
01361 n_min, n_max = sheet.cell_value(j+1,4), sheet.cell_value(j+1,10)
01362 Ks_min, Ks_max = sheet.cell_value(j+1,5), sheet.cell_value(j+1,11)
01363 l_min, l_max = sheet.cell_value(j+1,6), sheet.cell_value(j+1,12)
01364
01365 v = lhs(stats.uniform,[],(6,self.n_ens))
01366
01367 shp_ens = {}
01368 shp_ens['thetar'] = thetar_min + (thetar_max-thetar_min)*v[0,:]
01369 shp_ens['thetas'] = thetas_min + (thetas_max-thetas_min)*v[1,:]
01370 shp_ens['alpha'] = alpha_min + (alpha_max-alpha_min)*v[2,:]
01371 shp_ens['n'] = n_min + (n_max-n_min)*v[3,:]
01372 shp_ens['Ks'] = Ks_min + (Ks_max-Ks_min)*v[4,:]
01373 shp_ens['l'] = l_min + (l_max-l_min)*v[5,:]
01374
01375 self.shp_ens = shp_ens
01376
01377 def _shp_cur(self):
01378 """
01379 read the current soil hydraulic parameters
01380 """
01381 soil_par = {}
01382 soil_par['thetar'] = self.shp_ens['thetar'][self.ens]
01383 soil_par['thetas'] = self.shp_ens['thetas'][self.ens]
01384 soil_par['alpha'] = self.shp_ens['alpha'][self.ens]
01385 soil_par['n'] = self.shp_ens['n'][self.ens]
01386 soil_par['Ks'] = self.shp_ens['Ks'][self.ens]
01387 soil_par['l'] = self.shp_ens['l'][self.ens]
01388 soil_par['m'] = 1-1/soil_par['n']
01389
01390
01391 soil_par['evap_1'] = self.psi2theta(-0.33, soil_par['thetar'], soil_par['thetas'],
01392 soil_par['alpha'], soil_par['m'], soil_par['n'])
01393
01394 soil_par['evap_0'] = self.psi2theta(-15, soil_par['thetar'], soil_par['thetas'],
01395 soil_par['alpha'], soil_par['m'], soil_par['n'])
01396
01397 self.soil_par = soil_par
01398
01399 self.nc_thetar[self.ens] = soil_par['thetar']
01400 self.nc_thetas[self.ens] = soil_par['thetas']
01401 self.nc_alpha[self.ens] = soil_par['alpha']
01402 self.nc_n[self.ens] = soil_par['n']
01403 self.nc_Ks[self.ens] = soil_par['Ks']
01404 self.nc_l[self.ens] = soil_par['l']
01405
01406
01407 def initialize(self):
01408 """
01409 this initializes all the required variables
01410 and open the netcdf file for writting
01411 """
01412 max_t = int(self.final_time/self.dt_flux)
01413
01414 self.max_t = max_t
01415 self.iter_dt = 1
01416
01417
01418 file = nc.NetCDFFile(self.ofile_name, 'w')
01419 setattr(file, 'title', 'output of the model ambhas.richards_glue')
01420 now = datetime.datetime.now()
01421 setattr(file, 'description', 'The model was run at %s'%(now.ctime()))
01422 file.createDimension('depth', self.no_layer)
01423 file.createDimension('time', self.max_t+1)
01424 file.createDimension('ensemble', self.n_ens)
01425
01426
01427 varDims = 'depth',
01428 depth = file.createVariable('depth', 'd', varDims)
01429 depth.units = 'm'
01430 depth[:] = np.tile(self.dz,self.no_layer).cumsum()-self.dz/2
01431
01432
01433 varDims = 'time',
01434 self.nc_year = file.createVariable('year', 'd', varDims)
01435 self.nc_doy = file.createVariable('doy', 'd', varDims)
01436
01437
01438 varDims = 'ensemble', 'depth', 'time'
01439 self.nc_sm = file.createVariable('sm','d', varDims)
01440 self.nc_sm.units = 'v/v'
01441 self.nc_sm[:,:,0] = self.theta
01442
01443
01444 varDims = 'ensemble','time'
01445 self.nc_aet = file.createVariable('aet','d',varDims)
01446 self.nc_aet.units = 'mm'
01447 self.nc_recharge = file.createVariable('recharge','d',varDims)
01448 self.nc_recharge.units = 'mm'
01449
01450
01451 varDims = 'ensemble',
01452 self.nc_thetar = file.createVariable('thetar','d',varDims)
01453 self.nc_thetas = file.createVariable('thetas','d',varDims)
01454 self.nc_alpha = file.createVariable('alpha','d',varDims)
01455 self.nc_n = file.createVariable('n','d',varDims)
01456 self.nc_Ks = file.createVariable('Ks','d',varDims)
01457 self.nc_l = file.createVariable('l','d',varDims)
01458
01459 self.nc_file = file
01460
01461 def _read_initial_condition(self):
01462 """
01463 read initial condition
01464 """
01465
01466 j = self.ind['initial_condition']
01467
01468 book = xlrd.open_workbook(self.input_file)
01469 sheet = book.sheet_by_name('initial_condition')
01470
01471 data_len = sheet.nrows-1
01472 theta = np.zeros(data_len)
01473
01474 for i in xrange(data_len):
01475 theta[i] = sheet.cell_value(i+1,j-1)
01476
01477 self.theta = theta
01478
01479
01480 def _unsat(self):
01481 """
01482 top boundary: atmoshpheric
01483 bottom boundary: gravity drainage
01484 """
01485
01486 thetar = self.soil_par['thetar']
01487 thetas = self.soil_par['thetas']
01488 alpha = self.soil_par['alpha']
01489 n = self.soil_par['n']
01490 m = self.soil_par['m']
01491 l = self.soil_par['l']
01492 Ks = self.soil_par['Ks']
01493 nz = self.no_layer
01494
01495 theta = 1.0*self.theta
01496
01497
01498
01499 iter_dt = max(24,int(np.ceil(self.rain_cur*self.dt_flux*1000/0.15)))
01500 self.iter_dt = int(max(iter_dt,0.75*self.iter_dt))
01501
01502
01503
01504
01505
01506 recharge_day = 0
01507 aet_day = 0
01508
01509
01510 for i in range(self.iter_dt):
01511 dt = self.dt_flux/self.iter_dt
01512
01513
01514 smi = (self.theta[0]-self.soil_par['evap_0'])/(self.soil_par['evap_1']-self.soil_par['evap_0'])
01515 if smi<0: smi=0
01516 if smi>1: smi=1
01517 aet = smi*self.pet_cur
01518 Bvalue = self.rain_cur-aet
01519
01520 K = self.theta2kr(theta,thetar,thetas,m,l,Ks)
01521 smc = self.smcf(theta,thetar,thetas,alpha,m,n)
01522 psi = self.theta2psi(theta,thetar,thetas,m,n,alpha)
01523
01524
01525 Kmid = np.empty(nz+1)
01526 Kmid[0] = 0
01527 for i in range(1,nz):
01528 Kmid[i] = 0.5*(K[i]+K[i-1])
01529 Kmid[nz] = K[nz-1]
01530
01531
01532 A = np.empty(nz)
01533 B = np.empty(nz)
01534 C = np.empty(nz)
01535 D = np.empty(nz)
01536 dz = self.dz
01537 dz2 = dz**2
01538
01539 for i in range(nz):
01540 A[i] = -(Kmid[i]/dz2)
01541 B[i] = smc[i]/dt+(Kmid[i+1]+Kmid[i])/dz2
01542 C[i] = A[i]
01543 D[i] = smc[i]*psi[i]/dt-(Kmid[i+1]-Kmid[i])/dz
01544
01545 i = 0
01546 A[0] = 0
01547 B[0] = smc[i]/dt+(Kmid[1])/dz2
01548 D[0] = smc[i]*psi[i]/dt+(Bvalue-Kmid[1])/dz
01549
01550
01551 B[nz-1] = smc[nz-1]/dt+(Kmid[nz])/dz2
01552 C[nz-1] = 0
01553 D[nz-1] = smc[nz-1]*psi[nz-1]/dt-(Kmid[nz]-Kmid[nz-1])/dz
01554
01555
01556 beta = np.empty(nz)
01557 gamma = np.empty(nz)
01558 u = np.empty(nz)
01559 beta[0] = B[0]
01560 gamma[0] = D[0]/beta[0]
01561 for i in range(1,nz):
01562 beta[i] = B[i]-(A[i]*C[i-1])/(beta[i-1])
01563 gamma[i] = (D[i]-A[i]*gamma[i-1])/(beta[i])
01564
01565 u[nz-1] = gamma[nz-1]
01566 for i in range(nz-2,-1,-1):
01567 u[i] = gamma[i]-(C[i]*u[i+1])/beta[i]
01568
01569
01570 J = np.empty(nz+1)
01571 for i in range(1,nz):
01572 J[i] = Kmid[i]*(1-(u[i]-u[i-1])/dz)
01573 J[0] = Bvalue
01574 J[nz] = Kmid[nz]
01575
01576 J[nz] = J[nz]
01577
01578 flux = np.diff(J)*dt/dz
01579 theta = theta - flux
01580
01581 theta[theta>thetas] = 0.99*thetas
01582 theta[theta<thetar] = 1.01*thetar
01583
01584 aet_day += aet*dt
01585 recharge_day += J[nz]*dt
01586
01587 self.theta = theta
01588
01589
01590 self.nc_year[self.t] = (self.cur_year)
01591 self.nc_doy[self.t] = (self.cur_doy)
01592 self.nc_sm[self.ens,:,self.t+1] = theta
01593 self.nc_recharge[self.ens,self.t] = recharge_day
01594 self.nc_aet[self.ens,self.t] = aet_day
01595
01596
01597 if self.t == int(0.25*self.max_t):
01598 output_message = '25 % completed'
01599 self._colored_output(output_message, 32)
01600
01601 elif self.t == int(0.5*self.max_t):
01602 output_message = '50 % completed'
01603 self._colored_output(output_message, 32)
01604
01605 elif self.t == int(0.75*self.max_t):
01606 output_message = '75 % completed'
01607 self._colored_output(output_message, 32)
01608
01609 elif self.t == self.max_t-1:
01610 output_message = '100 % completed'
01611 self._colored_output(output_message, 32)
01612
01613
01614 if __name__=='__main__':
01615
01616
01617
01618
01619
01620
01621
01622
01623
01624
01625
01626
01627
01628
01629 maddur_glue = RICHARDS_1D_GLUE('/home/tomer/richards/input/maddur_glue_57.xls')
01630
01631