/Users/huttone/Devel/sedflux-new/sedflux/trunk/ew/inflow/inflow.c

Go to the documentation of this file.
00001 //---
00002 //
00003 // This file is part of sedflux.
00004 //
00005 // sedflux is free software; you can redistribute it and/or modify
00006 // it under the terms of the GNU General Public License as published by
00007 // the Free Software Foundation; either version 2 of the License, or
00008 // (at your option) any later version.
00009 //
00010 // sedflux is distributed in the hope that it will be useful,
00011 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00012 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00013 // GNU General Public License for more details.
00014 //
00015 // You should have received a copy of the GNU General Public License
00016 // along with sedflux; if not, write to the Free Software
00017 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00018 //
00019 //---
00020 
00021 #include <stdlib.h>
00022 #include <stdio.h>
00023 #include <math.h>
00024 #include <limits.h>
00025 #include "inflow.h"
00026 #include <sed/sed_sedflux.h>
00027 #include <utils/utils.h>
00028 
00074 gboolean
00075 inflow( double day        , double x[]         , double slope[]  , double width[] ,
00076         int n_nodes       , double dx          , double x_dep    , double river_width ,
00077         double river_vel  , double river_depth , double dc       , double *gzF0 ,
00078         double *grain_dia , double *lambda     , double *rho_sed , double *rho_grain,
00079         int n_grains      , double rho_rw      , double rho_flow , Inflow_const_st* c ,
00080         double **deposited, double** eroded    , FILE *fpout )
00081 {
00082    double **save_data;
00083 
00084    // Node parameters.
00085    double q, q0, J, J0, conc, conc0;
00086    double h, h0, u, u0;
00087    double rho, rho0, rhoF, rhoF0, rhoS, rhoS0;
00088    double *Cgrain, *Cgrain0, *Jgrain, *Jgrain0;
00089    double RI, g0, E, A, A1, A2, A3;
00090    // Constants.
00091    double Ea, Eb, sua, sub, Cd, tanPhi, mu, rhoSW;
00092    double gravity;
00093    // Variables.
00094    double mass=0,massin=0,massout,masseroded=0;
00095    double gamma, beta, sineBeta, cosBeta;
00096    double denominatorExponent;
00097    double max_grain_dia, max_ws;
00098    double *Cgrain_init, *gzF;
00099    double *phe, erosion, *fDep, *fEro, uCritical;
00100    double muEffective, *ws;
00101    double depth=0, max_depth=0;
00102    int i, n;
00103    Inflow_phe_query_st phe_query;
00104 
00105    Ea     = c->e_a;
00106    Eb     = c->e_b;
00107    sua    = c->sua;
00108    sub    = c->sub;
00109    Cd     = c->c_drag;
00110    tanPhi = c->tan_phi;
00111    mu     = c->mu_water;
00112    rhoSW  = c->rho_sea_water;
00113 
00114    if ( fpout )
00115       save_data = eh_new_2( double , 13 , n_nodes );
00116 
00117    Jgrain      = eh_new( double , n_grains );
00118    Jgrain0     = eh_new( double , n_grains );
00119    Cgrain      = eh_new( double , n_grains );
00120    Cgrain_init = eh_new( double , n_grains );
00121    Cgrain0     = eh_new( double , n_grains );
00122    ws          = eh_new( double , n_grains );
00123    fEro        = eh_new( double , n_grains );
00124    fDep        = eh_new( double , n_grains );
00125    phe         = eh_new( double , n_grains );
00126    gzF         = eh_new( double , n_grains );
00127 
00128    gravity = sed_gravity();
00129    denominatorExponent = exp(1.0)-1.;
00130    
00131    // Initial conditions.
00132    for (n=0;n<n_grains;n++)
00133       gzF[n] = gzF0[n];
00134 //   h0 = river_depth*river_width/width[0];
00135    u0    = river_vel;
00136    q0    = dc;
00137    rho0  = rho_rw;
00138    rhoF0 = rho_flow;
00139    rhoS0 = eh_dbl_array_mean_weighted(gzF,n_grains,rho_grain);
00140    rhoS0 = eh_dbl_array_mean_weighted(gzF,n_grains,rho_grain);
00141    conc0 = (rhoF0-rho0)/(rhoS0-rho0);
00142    J0    = q0*conc0;
00143    for (n=0;n<n_grains;n++)
00144    {
00145       Cgrain0[n] = conc0*gzF[n];
00146       Cgrain_init[n] = Cgrain0[n];
00147       Jgrain0[n] = Cgrain0[n]*q0;
00148    }
00149 
00150    for (n=0;n<n_grains;n++)
00151       massin = conc0*q0*rho_grain[n];
00152 
00153 //   fprintf(stderr,"inflow : Mass of sediment entering the flow (kg) : %f\n",massin*day);
00154    
00155 /*
00156    for (i=0;i<n_nodes && conc0 > 1e-6 ;i++)
00157 */
00158    for (i=0;i<n_nodes ;i++)
00159    {
00160       // Find slope.
00161       beta = -slope[i];
00162 
00163       if ( fabs(beta) > .2)
00164          beta = .2;
00165 
00166       sineBeta = sin(beta);
00167       cosBeta = cos(beta);
00168       
00169       // Richardson number of the flow.
00170       g0 = gravity*(rhoS0-rho0)/rho0;
00171       RI = g0*cosBeta*J0/(u0*u0*u0*width[i]);
00172       
00173       // Entrainment coefficient.
00174       if ( sineBeta > 0.01 )
00175          E = Ea / (Eb+RI);
00176       else if ( sineBeta < 0.01 && sineBeta > 0 )
00177          E = 0.072*sineBeta;
00178       else
00179          E = 0.;
00180 
00181       // Gravity.
00182       if ( sineBeta > 0 )
00183       {
00184          A1 = g0*J0*sineBeta/u0/q0;
00185 /*
00186          Cd = 0.004;
00187 */
00188       }
00189       else
00190       {
00191          A1 = 0.;
00192          A1 = g0*J0*sineBeta/u0/q0;
00193 /*
00194          Cd = 0.004;
00195 */
00196       }
00197 
00198       depth += sineBeta;
00199       if ( depth > max_depth )
00200          max_depth = depth;
00201       if ( depth < max_depth )
00202       {
00203          if ( sineBeta > 0 )
00204             A1 *= 0.05;
00205          else
00206             A1 *= 0.05;
00207       }
00208       
00209       // Friction.
00210       A2 = -(E+Cd)*u0*u0*width[i]/q0;
00211       
00212       // Internal friction.
00213       gamma = tanPhi*(exp(J0/q0)-1.)/denominatorExponent;
00214       A3 = -.1*g0*J0*cosBeta*gamma/u0/q0;
00215 
00216 /*
00217       if ( A3 < -.02 ) A3 = -.02;
00218       if ( A2 < -.02 ) A2 = -.02;
00219 */
00220 
00221 /*      
00222       Cd = 0.004;
00223 */
00224       
00225       A = A1 + A2 + A3;
00226       
00227       u = u0 + A*dx;
00228 
00229 // Drift velocity      
00230       if ( u < 0.01 )
00231          u = 0.01;
00232 
00233       if ( u <= 0 )
00234       {
00235          fprintf(stderr,"inflow : Error.  Velocity dropped below zero m/s\n");
00236          if ( fpout )
00237          {
00238             for (n=0;n<13;n++)
00239                fwrite( save_data[n] , n_nodes , sizeof(double) , fpout );
00240             eh_free_2(save_data);
00241          }
00242          return FALSE;         
00243       }
00244 /*
00245       q = q0;
00246       h = q0/width[i]/u;
00247 */
00248 
00249       // New volume discharge.
00250       q = q0 + E*u*width[i]*dx;
00251       
00252       // New thickness.
00253       h = q/width[i]/u;
00254       
00255       /* New fluid density.
00256       */
00257       rho = rho0+E*u*width[i]*(rhoSW-rho0)/q0*dx;
00258       
00259       // Depth of erosion for a day.
00260       erosion = (Cd*rhoF0*u0*u0-sub)/sua;
00261 
00262       if ( erosion < 0 ) erosion = 0.;
00263 
00264       // Need to determine what phe is at each bin.
00265       // For the standalone version it will simply be an input and 
00266       // will be constant over the run.  For the version in SEDFLUX,
00267       // we'll have to write a routine to average the top 'erosion'
00268       // bins.
00269       phe_query.x           = x[i];
00270       phe_query.dx          = dx;
00271       phe_query.erode_depth = erosion;
00272       phe_query.phe         = phe;
00273 
00274       (*(c->get_phe))( &phe_query , c->get_phe_data );
00275 
00276       erosion = phe_query.erode_depth;
00277 
00278 // Find the maximum grain diameter that has a concentration greater than .01.
00279       for (n=0,max_grain_dia=0.;n<n_grains;n++)
00280          if ( grain_dia[n] > max_grain_dia )
00281             max_grain_dia = grain_dia[n];
00282 
00283       // Settling velocities.
00284       muEffective = mu*(1.+2.5*conc0);
00285 /*
00286        ws[n] = g0/18.*grain_dia[n]*grain_dia[n]/muEffective;
00287        ws[n] = (rhoS0-rhoF0)*gravity*grain_dia[n]*grain_dia[n]/18./(rhoF0*muEffective);
00288 */
00289          
00290       max_ws = g0/18.*max_grain_dia*max_grain_dia/muEffective;
00291 
00292       /* Critical velocity for deposition.
00293       uCritical = ws[n]/sqrt(.016/8.);
00294       uCritical = ws[n]/sqrt(Cd);
00295       */
00296 
00297       uCritical = max_ws/sqrt(Cd);
00298 
00299       for (n=0;n<n_grains;n++)
00300       {
00301          // Rate of erosion (m^2/s).
00302          fEro[n] = erosion*phe[n]*S_DAYS_PER_SECOND*width[i];
00303 
00304          // Deposition rate (m^2/s).
00305          // NOTE: we have divided the deposition rates by 10 to account
00306          // for the lambda obtained from plume data.
00307          if ( u >= uCritical )
00308             fDep[n] = 0.;
00309          else
00310             fDep[n] = -lambda[n]*Jgrain0[n]/u*(1-pow(u/uCritical,2.))/10.;
00311 /*
00312             fDep[n] = -lambda[n]*Jgrain0[n]/u*(1-pow(u/uCritical,1/(n+1)))/10.;
00313             fDep[n] = -lambda[n]*Jgrain0[n]/u*(1-sqrt(u/uCritical));
00314             fDep[n] = -lambda[n]*Jgrain0[n]/u*(1-sqrt(u*u/uCritical/uCritical));
00315             fDep[n] = -lambda[n]*Jgrain0[n]/u*(1-(u*u/uCritical/uCritical));
00316 */
00317 
00318 /*         
00319          if ( (x[i]-x[0]) <= x_dep )
00320             fDep[n] = 0.;
00321 */
00322          
00323          Jgrain[n] = Jgrain0[n] + (fDep[n]+fEro[n])*dx;
00324          
00325          // Deposit thickness.  Convert the deposit sediment to the appropriate
00326          // porosity.  The erosion part is given in meters of bottom sediment
00327          // eroded and so is already at the correct porosity.
00328 /*
00329 #if defined(INFLOW_STANDALONE)
00330          deposit[n][i] = -(fDep[n]*rho_grain[n]/rho_sed[n]+fEro[n])*day/width[i];
00331 #else
00332          deposit[n][i] = -(fDep[n]*rho_grain[n]/rho_sed[n])*day/width[i];
00333 #endif
00334 */
00335          deposited[n][i] = -(fDep[n]*rho_grain[n]/rho_sed[n])*day/width[i];
00336          eroded   [n][i] = -(fEro[n])*day/width[i];
00337 
00338          mass       += -fDep[n]*rho_grain[n]*dx;
00339          masseroded +=  fEro[n]*rho_grain[n]*dx;
00340       }
00341 
00342       // Get new grain size fractions in the flow.
00343       for (n=0,J=0;n<n_grains;n++)
00344          J += Jgrain[n];
00345       for (n=0;n<n_grains;n++)
00346          gzF[n] = Jgrain[n]/J;
00347       
00348       rhoS = eh_dbl_array_mean_weighted(gzF,n_grains,rho_grain);
00349       for (n=0,conc=0;n<n_grains;n++)
00350       {
00351          Cgrain[n] = Jgrain[n]/q;
00352          conc += Cgrain[n];
00353       }
00354 
00355       if ( fpout )
00356       {
00357          save_data[0][i] = conc0;
00358          save_data[1][i] = rhoF0;
00359          save_data[2][i] = rho0;
00360          save_data[3][i] = A/dx;
00361          save_data[4][i] = A1;
00362          save_data[5][i] = A2;
00363          save_data[6][i] = A3;
00364          save_data[7][i] = E;
00365          save_data[8][i] = h0;
00366          save_data[9][i] = u0;
00367          save_data[10][i] = RI;
00368          for (n=0;n<n_grains;n++)
00369          {
00370             save_data[11][i] += fEro[n];
00371             save_data[12][i] += fDep[n];
00372          }
00373       }
00374 
00375       rhoF = conc*(eh_dbl_array_mean_weighted(gzF,n_grains,rho_grain)-rho)+rho;
00376 
00377       // Save the node values.
00378       rho0 = rho;
00379       rhoS0 = rhoS;
00380       rhoF0 = rhoF;
00381       h0 = h;
00382       u0 = u;
00383       q0 = q;
00384       J0 = J;
00385       conc0 = conc;
00386       for (n=0;n<n_grains;n++)
00387       {
00388          Jgrain0[n] = Jgrain[n];
00389          Cgrain0[n] = Cgrain[n];
00390       }
00391       
00392    }
00393 
00394    massout = mass;
00395 //   fprintf(stderr,"inflow : Mass deposited by the flow (kg) : %f\n",mass*day);
00396 //   fprintf(stderr,"inflow : Mass eroded by the flow (kg) : %f\n",masseroded*day);
00397 //   fprintf(stderr,"inflow : Mass balance error of : %f\n",(massin+masseroded-massout)/massin);
00398 
00399 /* Deposit the sediment with the appropriate porosity.
00400    for (i=0;i<n_nodes;i++)
00401       for (n=0;n<n_grains;n++)
00402          deposit[n][i] *= (rho_grain[n]/rho_sed[n]);
00403 */
00404 
00405 
00406    if ( fpout )
00407    {
00408       for ( n=0 ; n<13 ; n++ )
00409          fwrite( save_data[n] , n_nodes , sizeof(double) , fpout );
00410       eh_free_2( save_data );
00411    }
00412    
00413    eh_free( Jgrain      );
00414    eh_free( Jgrain0     );
00415    eh_free( Cgrain      );
00416    eh_free( Cgrain0     );
00417    eh_free( ws          );
00418    eh_free( fDep        );
00419    eh_free( fEro        );
00420    eh_free( phe         );
00421    eh_free( gzF         );
00422    eh_free( Cgrain_init );
00423 
00424    return TRUE;
00425 }
00426  

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