/Users/huttone/Devel/sedflux-new/sedflux/trunk/ew/hydrotrend/hydrorain.c

Go to the documentation of this file.
00001 /*
00002  * HydroRain.c          Calculates the daily rain fall and rain derived runoff
00003  *                      for each altitude bin.  Also calculates the rain
00004  *                      derived groundwater flow, evaporation and time lag.
00005  *
00006  *      Physically based precipitation/groundwater/evaporation routine based on:
00007  *      M. Sivapalan, J.K. Ruprecht, N.R. Viney, 1996.  Water and Salt Balance
00008  *              Modelling to predict the effects of Land-Use changes in forested
00009  *              catchments. 1. Small Catchment water balance model.
00010  *              Hydrological Processes, V.10, P.393-411.
00011  *
00012  *      Original:       M.D. Morehead   July 1996       precip.F
00013  *      1st C Version:  M.D. Morehead   1998
00014  *      Author2:        A.J. Kettner    October 2002
00015  *
00016  * Variable             Def.Location    Type    Units   Usage
00017  * --------             ------------    ----    -----   -----
00018  * Mannualin    HydroRain.c             double  m^3/a   annual mass input to rain routine
00019  * Mannualout   HydroRain.c             double  m^3/a   annual mass output from rain routine
00020  * Minput               HydroRain.c             double  m^3             daily mass input to rain routine
00021  * Mout                 HydroRain.c             double  m^3             daily mass input to rain routine
00022  * Mwrap                HydroRain.c             double  m^3             mass input to rain routine from previous year
00023  * Pground              HydroRain.c             double  m/day   ground precipitation
00024  * Qi1, Qi2             HydroRain.c             double  m^3/s   temporary variables for Qinfiltration
00025  * Qinfiltexcess HydroRain.c    double  m^3/s   river discharge from infiltration excess
00026  * Qinfiltration HydroRain.c    double  m^3/s   = Pground-Qrain[ii]
00027  * Qrainlag             HydroRain.c             double  m^3/s   temporary array for routing
00028  * Qsaturationexcess HydroRain.c double m^3/s   river discharge from saturation excess
00029  * Qsslag               HydroRain.c             double  m^3/s   temporary array for routing
00030  * Qtogw                HydroRain.c             double  m^3/s   Rain derived discharge to the GW pool
00031  * err                  various                 int             -               error flag, halts program
00032  * ii                   various                 int             -               temporary loop counter
00033  * infiltrate   HydroRain.c             double  m/day   infiltration capacity
00034  * jj                   various                 int             -               temporary loop counter
00035  * kk                   various                 int             -               temporary loop counter
00036  * ratio1               HydroRain.c             double  -               ratio of the fullness of the GW pool
00037  * ratio2               HydroRain.c             double  -               linear interpolation for infiltration
00038  * shldday[]    HydroRain.c             double  m^3/s   shoulder discharge array
00039  *
00040  */
00041 
00042 #include <stdlib.h>
00043 
00044 #include "hydroclimate.h"
00045 #include "hydroparams.h"
00046 #include "hydrotimeser.h"
00047 
00048 #ifdef DBG
00049 #include "hydroinout.h"
00050 #endif
00051 
00052 /*----------------------
00053  *  Start of HydroRain
00054  *----------------------*/
00055 int hydrorain()
00056 {
00057 
00058 /*-------------------
00059  *  Local Variables
00060  *-------------------*/
00061 #ifdef DBG
00062 FILE *fid;
00063 #endif
00064 
00065 double infiltrate, Pground[daysiy];
00066 double Qsaturationexcess[daysiy], Qinfiltexcess[daysiy];
00067 double Qinfiltration[daysiy], Qi1, Qi2, Qtogw, shldday[maxday];
00068 double  Mout, Mannualin, Mannualout;
00069 double Qrainlag[maxday], Qsslag[maxday];
00070 double ratio1, ratio2;
00071 int ii, jj, kk, err;
00072 
00073 /*------------------------
00074  *  Initialize Variables
00075  *------------------------*/
00076 err = 0;
00077 MPrain = 0.0;
00078 Mannualin = 0.0;
00079 Mannualout = 0.0;
00080 for( ii=0; ii<maxday; ii++ ) {
00081    Qrainlag[ii] = 0.0;
00082    Qsslag[ii]   = 0.0;
00083    shldday[ii]  = 0.0;
00084 }
00085 
00086 /*-------------------------
00087  *  Loop through the year
00088  *-------------------------*/
00089 for( ii=0; ii<daysiy; ii++ ) {
00090 
00091    /*----------------------------------------------------------------------
00092     *  Find the area being rained on
00093     *  Temperature must be above freezing
00094     *    All rain (P with T>=0.0) is handled here.
00095     *    All snow/ice (P with T<0.0) is handled in HydroSnow/HydroGlacier
00096     *----------------------------------------------------------------------*/
00097    rainarea[ii] = 0.0;
00098    for( kk=0; kk<nelevbins && kk<FLAindex[ii]; kk++ )
00099         rainarea[ii] += areabins[kk];
00100 
00101    /*--------------------------------------------------------------
00102     *  Get the remaining GW pool from the previous day
00103     *  the first day of the year is taken care of in HydroTrend.c
00104     *--------------------------------------------------------------*/
00105    if( ii > 0 )
00106       gwstore[ii] = gwstore[ii-1];
00107 
00108    /*---------------------------------------------------------------------------
00109     *  Calculate groundwater discharge (subsurface storm flow) and evaporation
00110     *  Qss (m^3/s), Egw (m/day), gwstore (m^3)
00111     *  Only drain down to gwmin
00112     *
00113     *  Having this first builds in a one delay lag to the GW discharge
00114     *  This also helps prevent the Glacier and Snow from totally filling
00115     *  the pool, and leaves some room for Rain.
00116     *
00117     *  The lines E=max(E,0.0) and Qss=max(Qss,0.0) prevent rare errors due
00118     *  to very small numbers when gwstore is close to gwmin.
00119     *---------------------------------------------------------------------------*/
00120    Egw[ii] = mn( alphagwe[ep]*pow((gwstore[ii]-gwmin[ep])/(gwmax[ep]-gwmin[ep]),betagwe[ep]), (gwstore[ii]-gwmin[ep])/totalarea[ep] );
00121    Egw[ii] = mx( Egw[ii], 0.0 );
00122    gwstore[ii] -= Egw[ii] * totalarea[ep];
00123 
00124    Qss[ii] = mn( alphass[ep]*pow((gwstore[ii]-gwmin[ep])/(gwmax[ep]-gwmin[ep]),betass[ep]), (gwstore[ii]-gwmin[ep])/dTOs );
00125    Qss[ii] = mx( Qss[ii], 0.0 );
00126    gwstore[ii] -= Qss[ii] * dTOs;
00127 
00128    /*-------------------------------------------------------------------------------
00129     *  Add the other sources to the Groundwater pool
00130     *  allocate overflow to Qexceedgw and add to river discharge in HydroSumFlow.c
00131     *-------------------------------------------------------------------------------*/
00132    if( gwstore[ii] < gwmax[ep] )
00133        Qtogw = mn( Qicetogw[ii]+Qnivaltogw[ii], (gwmax[ep]-gwstore[ii])/dTOs );
00134    else
00135        Qtogw = 0.0;
00136 
00137    gwstore[ii]  += Qtogw*dTOs;
00138    Qexceedgw[ii] = Qicetogw[ii] + Qnivaltogw[ii] - Qtogw;
00139 
00140    /*--------------------------------------
00141     *  Calculate commonly used quantities
00142     *--------------------------------------*/
00143    ratio1  = (gwstore[ii]-gwmin[ep])/(gwmax[ep]-gwmin[ep]);
00144    ratio2  = (Ko[ep]-pcr[ep])/(pmax[ep]-pcr[ep]);
00145 
00146    /*------------------------------------------------------------
00147     *  Calculate the Ground Precipitation (m/day)
00148     *  the remainder is evaporated from canopy interception
00149     *  keep track of evaporation (m/day) for mass balance check
00150     *------------------------------------------------------------*/
00151    Pground[ii] = mx( alphag[ep] + betag[ep]*Pdaily[ii], 0.0 );
00152    Ecanopy[ii] = Pdaily[ii] - Pground[ii];
00153 
00154    /*----------------------------------------------------
00155     *  Calculate Surface Runoff (m^3/s)
00156     *  Sum of saturation excess and infiltration excess
00157     *  remainder goes to groundwater pool
00158     *----------------------------------------------------*/
00159 
00160    /*-----------------------------------------------------      
00161     *  Calculate the Saturation Excess
00162     *  this is all the rain that falls on saturated soil
00163     *  and cannot go into the groundwater pool
00164     *-----------------------------------------------------*/
00165     if( gwstore[ii] > gwmin[ep] )
00166        Qsaturationexcess[ii] = (Pground[ii]*rainarea[ii]/dTOs)*alphac[ep]*pow(ratio1,betac[ep]);
00167     else
00168        Qsaturationexcess[ii] = 0.0;
00169 
00170    /*--------------------------------------------------------------------
00171     *  Calculate Theoretical Infiltration Rate based on Precip  (m/day)
00172     *    a) below pcr, 100% infiltration
00173     *    b) above pmax, max hydr. cond.
00174     *    c) linear fit in between
00175     *--------------------------------------------------------------------*/
00176    if( Pground[ii] < pcr[ep] )
00177       infiltrate = Pground[ii];
00178    else if( Pground[ii] > pmax[ep] )
00179       infiltrate = Ko[ep];
00180    else infiltrate = Pground[ii]*ratio2 + pcr[ep]*(1-ratio2);
00181 
00182    /*-------------------------------------------------------------------------
00183     *  Calculate infiltration "discharge" to the groundwater pool (m^3/s)
00184     *  1st) calculate the theoretical infiltration discharge
00185     *  2nd) calculate the discharge that fills the groundwater pool in a day
00186     *  take the minimum so as not to overfill the pool
00187     *-------------------------------------------------------------------------*/
00188    if( gwstore[ii] > gwmin[ep] )
00189       Qi1 = infiltrate*(1.0-ratio1)*rainarea[ii]/dTOs;
00190    else
00191       Qi1 = infiltrate*rainarea[ii]/dTOs;
00192    Qi2 = (gwmax[ep]-gwstore[ii])/dTOs;
00193 
00194    Qinfiltration[ii] = mn( Qi1, Qi2 );
00195    gwstore[ii] += Qinfiltration[ii]*dTOs;
00196 
00197    /*--------------------------------------------------
00198     *  Calculate the infiltration excess (m^3/s)
00199     *  this is the runoff that:
00200     *    a) is not due to saturation excess
00201     *    b) and can not get into the groundwater pool
00202     *--------------------------------------------------*/
00203    Qinfiltexcess[ii] = mx( 0.0, \
00204         Pground[ii]*rainarea[ii]/dTOs - Qsaturationexcess[ii] - Qinfiltration[ii]);
00205 
00206    /*-------------------------------------------------------
00207     *  Calculate Surface Runoff (m^3/s)
00208     *    Sum of saturation excess and infiltration excess
00209     *    remainder went to groundwater pool or evaporation
00210     *-------------------------------------------------------*/
00211    Qrain[ii] = Qsaturationexcess[ii] + Qinfiltexcess[ii];
00212 
00213    /*---------------------------------------------------------------------------
00214     *  Add the carryover from the previous year
00215     *
00216     *  NOTE: adding at this stage assumes the wrapped rain is due to the wave
00217     *  propagating down the river and that this water does not have additional
00218     *  contributions to the Evap or gw pools
00219     *---------------------------------------------------------------------------*/
00220         
00221    if( ii < wrapday )
00222       Qrain[ii] += Qrainwrap[ii];
00223 
00224    /*------------------------------------------------
00225     *  Add to the daily rain mass balance (m^3/day)
00226     *------------------------------------------------*/
00227    Minput = Pdaily[ii]*rainarea[ii];
00228    if( ii < maxday-daysiy )
00229       Minput += Qrainwrap[ii]*dTOs;
00230    Mout = Ecanopy[ii]*rainarea[ii] + (Qrain[ii] + Qinfiltration[ii])*dTOs;
00231 
00232    /*-----------------------------------------------
00233     *  Add to the annual rain mass balance (m^3/a)
00234     *-----------------------------------------------*/
00235    MPrain       += Pdaily[ii]*rainarea[ii];
00236    Mannualin    += Minput;
00237    Mannualout   += Mout;
00238 
00239    /*---------------------------------------------
00240     *  Check the daily rain mass balance (m^3/d)
00241     *---------------------------------------------*/
00242    if( Minput > 0. && fabs(Mout-Minput)/Minput > masscheck ) {
00243         fprintf( stderr, "ERROR: in HydroRain.c: \n" );
00244         fprintf( stderr, "   Daily Mass Balance Error: Mout != Minput \n\n" );
00245         fprintf( stderr, "   fabs(Mout-Minput)/Minput > masscheck \n" );
00246         fprintf( stderr, "   note: masscheck set in HydroParams.h \n" );
00247         fprintf( stderr, "   masscheck = %e (m^3) \n", masscheck );
00248         fprintf( stderr, "   fabs(Mout-Minput)/Minput = %e (m^3) \n\n", fabs(Mout-Minput)/Minput );
00249     fprintf( stderr, " \t Minput = Pdaily + Qrainwrap \n" );
00250     fprintf( stderr, " \t Minput \t\t = %e (m^3) \n", Minput );
00251     fprintf( stderr, " \t Pdaily[ii]*rainarea[ii] \t = %e (m^3) \n", Pdaily[ii]*rainarea[ii] );
00252         if( ii < maxday-daysiy )
00253        fprintf( stderr, " \t Qrainwrap[ii] \t\t = %e (m^3) \n\n", Qrainwrap[ii]*dTOs );
00254         else
00255         fprintf( stderr, " \t Qrainwrap[ii] \t\t = 0.0 (m^3) \n\n" );
00256     fprintf( stderr, " \t Mout = Ecanopy + Qrain + Qinfiltration \n" );
00257     fprintf( stderr, " \t Mout \t\t = %e (m^3) \n", Mout );
00258     fprintf( stderr, " \t Ecanopy[ii]*rainarea[ii] \t\t = %e (m^3) \n", Ecanopy[ii]*rainarea[ii] );
00259     fprintf( stderr, " \t Qrain[ii]*dTOs \t\t = %e (m^3) \n", Qrain[ii]*dTOs );
00260     fprintf( stderr, " \t Qinfiltration[ii]*dTOs \t = %e (m^3) \n\n", Qinfiltration[ii]*dTOs );
00261         fprintf( stderr," \t Day \t d \t %d \n", ii+1 );
00262         fprintf( stderr," \t Year \t d \t %d \n", yr+1 );
00263         fprintf( stderr," \t Epoch \t - \t %d \n", ep+1 );
00264     exit(-1);
00265    }
00266 
00267    /*------------------------------------------------
00268     *  Check for over/under-filled groundwater pool
00269     *  return error
00270     *------------------------------------------------*/
00271    if( gwstore[ii] < gwmin[ep]/1.001 || gwstore[ii] > gwmax[ep]*1.001 ) {
00272       fprintf(stderr,"WARNING in HydroRain.c \n");
00273       fprintf(stderr,"   gwstore out of range. \n");
00274       fprintf(stderr,"   criteria:  gwmin < gwstore < gwmax \n");
00275       fprintf(stderr,"      gwmax   (m^3) = %e \n", gwmax[ep]    );
00276       fprintf(stderr,"      gwstore (m^3) = %e \n", gwstore[ii]  );
00277       fprintf(stderr,"      gwmin   (m^3) = %e \n", gwmin[ep]    );
00278       fprintf(stderr,"      Qss[ii] (m^3) = %e \n", Qss[ii]*dTOs     );
00279       fprintf(stderr,"      Egw[ii] (m^3) = %e \n", Egw[ii]*totalarea[ep] );
00280       err = 1;
00281    }  
00282 }       /* end day loop */
00283 
00284 /*------------------------------------------------------------
00285  *  Use shoulder on Rain events if the basin is large enough
00286  *------------------------------------------------------------*/
00287    if( (int)( basinlength[ep]/(avgvel[ep]*dTOs) ) >= 3 || Rvol[ep] > 0.0)  {
00288 
00289    /*------------------------------------------------------------
00290     *  Create the shoulder events (Murray's version of routing)
00291     *    there is one left (preceeding) day scaled as:
00292     *      shoulderleft*event
00293     *    the main event is scaled down to:
00294     *      shouldermain*event
00295     *    there are 1 or more right (following days) scaled to:
00296     *      shoulderright[]*event
00297     *
00298     *    1.0 = Sum(shoulderleft+shouldermain+shoulderright[])
00299     *------------------------------------------------------------*/
00300       ii = 0;
00301       if( Qrain[ii] > 0.0 ) {
00302           shldday[ii] += shoulderleft*Qrain[ii];
00303           for( jj=0; jj<shouldern-2; jj++ )
00304               shldday[ii+jj+1] += shoulderright[jj]*Qrain[ii];
00305           Qrain[ii] = shouldermain*Qrain[ii];
00306       }
00307       for( ii=1; ii<daysiy; ii++ ) {
00308           if( Qrain[ii] > 0.0 ) {
00309               shldday[ii-1] += shoulderleft*Qrain[ii];
00310               for( jj=0; jj<shouldern-2; jj++ )
00311                   shldday[ii+jj+1] += shoulderright[jj]*Qrain[ii];
00312               Qrain[ii] = shouldermain*Qrain[ii];
00313           }
00314       }
00315 
00316    /*----------------------------------------------------------------
00317     *  Add the shoulder events and the main events to get the total
00318     *  rain derived discharge
00319     *----------------------------------------------------------------*/
00320       for( ii=0; ii<maxday; ii++ )
00321           Qrain[ii] += shldday[ii];
00322    } /* endif shoulder */
00323 
00324 /*---------------------------------------------------------------------
00325  *  (Mark's version of routing)
00326  *  Remove any lagged discharge from the initial rain discharge event
00327  *---------------------------------------------------------------------*/
00328 for( ii=0; ii<daysiy; ii++ ) {
00329    for( kk=0; kk<nelevbins && kk<FLAindex[ii]; kk++ ) {
00330       Qrainlag[ii+distbins[kk]] += Qrain[ii]*areabins[kk]/rainarea[ii];
00331       Qsslag[ii+distbins[kk]]   +=   Qss[ii]*areabins[kk]/rainarea[ii];
00332       Qrain[ii] -= Qrain[ii]*areabins[kk]/rainarea[ii];
00333       Qss[ii]   -=   Qss[ii]*areabins[kk]/rainarea[ii];
00334    }  
00335 }
00336 
00337 /*------------------------------------
00338  *  Add the lagged discharge back in
00339  *------------------------------------*/
00340 for( ii=0; ii<maxday; ii++ ) {
00341    Qrain[ii] += Qrainlag[ii];
00342    Qss[ii]   += Qsslag[ii];   
00343 }
00344 
00345 /*--------------------------------
00346  *  Add the GW wrapped discharge
00347  *--------------------------------*/
00348 for( ii=0; ii<wrapday; ii++ )
00349    Qss[ii]   += Qsswrap[ii];
00350 
00351 /*---------------------------------------
00352  *  Check the annual mass balance (m^3)
00353  *---------------------------------------*/
00354 if( fabs(Mannualout-Mannualin)/Mannualin > masscheck ) {
00355      fprintf( stderr, "ERROR: in HydroRain.c: \n" );
00356      fprintf( stderr, "   Annual Mass Balance Error: Mout != Minput \n\n" );
00357 
00358      fprintf( stderr, "   fabs(Mannualout-Mannualin)/Mannualin > masscheck \n" );
00359      fprintf( stderr, "   note: masscheck set in HydroParams.h \n" );
00360      fprintf( stderr, "   masscheck = %e (m^3) \n", masscheck );
00361      fprintf( stderr, "   fabs(Mannualout-Mannualin)/Mannualin = %e (m^3) \n\n", fabs(Mannualout-Mannualin)/Mannualin );
00362 
00363      fprintf( stderr, " \t Mannualin = Sum(Pdaily+Qrainwrap) \t\t = %e (m^3) \n", Mannualin );
00364      fprintf( stderr, " \t Mannualout = Sum(Ecanopy+Qrain+Qinfiltration) \t = %e (m^3) \n\n", Mannualout );
00365 
00366      fprintf( stderr," \t Year \t d \t %d \n", yr+1 );
00367      fprintf( stderr," \t Epoch \t - \t %d \n\n", ep+1 );
00368 
00369      exit(-1);
00370 }
00371 
00372 #ifdef DBG
00373   /*--------------------------------------------
00374    *  print out summary stats for the log file
00375    *--------------------------------------------*/
00376   fprintf( fidlog, "\n HydroRain Summary values: \n\n" );
00377   fprintf( fidlog, " \t Mannualin = Sum(Pdaily+Qrainwrap) \t\t = %e (m^3) \n", Mannualin );
00378   fprintf( fidlog, " \t Mannualout = Sum(Ecanopy+Qrain+Qinfiltration) \t = %e (m^3) \n\n", Mannualout );
00379 
00380   /*---------------------------------------------------------------------------
00381    *  print out precipitation for error checking
00382    *  the % at the begining of each text line is for a comment line in matlab
00383    *  this enables the file to be read directly into matlab for plotting
00384    *---------------------------------------------------------------------------*/
00385   if( tblstart[ep] <= yr && yr <= tblend[ep] ) {
00386     if( (fid = fopen("hydro.p1","a+")) == NULL) {
00387       printf("  HydroRain ERROR: Unable to open the precipitation file hydro.p1 \n");
00388       printf("     non-fatal error, continueing. \n\n");
00389     }
00390     else {
00391         fprintf( fid,"%%\n%%\n%% HydroRain output: \n%%\n" );
00392         fprintf( fid,"%%Daily predicted precipitation and rain runoff for epoch %d \n%%\n", ep+1 );
00393         fprintf( fid,"%%Year \t Day \t Qrain(m^3/s) \t Qsat(m^3/s) \t Qine(m^3/s) \t Qinf(m^3/s) \t Qss(m^3/s) \n" );
00394         fprintf( fid,"%%---- \t --- \t ------------ \t ------------ \t ----------- \t ----------- \t ----------- \n" );
00395         for(ii=0; ii<daysiy; ii++ )
00396            fprintf( fid,"%d \t %d \t %f \t %f \t %f \t %f \t %f \n", yr, ii+1, \
00397                 Qrain[ii], Qsaturationexcess[ii], Qinfiltexcess[ii],Qinfiltration[ii],Qss[ii] );
00398         fclose(fid);
00399     }
00400   }
00401 
00402   if( tblstart[ep] <= yr && yr <= tblend[ep] ) {
00403     if( (fid = fopen("hydro.p2","a+")) == NULL) {
00404       printf("  HydroRain ERROR: Unable to open the precipitation file hydro.p2 \n");
00405       printf("     non-fatal error, continueing. \n\n");
00406     }
00407     else {
00408         fprintf( fid,"%%\n%%\n%% HydroRain output: \n%%\n" );
00409         fprintf( fid,"%%Daily predicted precipitation and rain runoff for epoch %d \n%%\n", ep+1 );
00410         fprintf( fid,"%%Year \t Day \t Pg(m^3/s) \t Ecanopy(m^3/s)\t Egw(m^3/s) \t GWstore(m^3/s) \n" );
00411         fprintf( fid,"%%---- \t --- \t ------------- \t ------------ \t ------------  \t ------------\n" );
00412         for(ii=0; ii<daysiy; ii++ )
00413            fprintf( fid,"%d \t %d \t %f \t %f \t %f \t %f \n", yr, ii+1, \
00414                 Pground[ii]*rainarea[ii]/dTOs, Ecanopy[ii]*rainarea[ii]/dTOs, Egw[ii]*totalarea[ep]/dTOs, \
00415                 gwstore[ii]/dTOs );
00416         fclose(fid);
00417     }
00418   }
00419 #endif
00420 
00421 return(err);
00422 }       /* end of HydroRain.c */

Generated on Fri Jan 4 18:04:14 2008 for sedflux by  doxygen 1.5.2