/Users/huttone/Devel/sedflux-new/sedflux/trunk/ew/plume/plumeconc.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 /*
00022  * PlumeConc     Calculate the plume concentrations for PLUME
00023  *
00024  *     Author:          M.D. Morehead
00025  *     Original:     May 1998
00026  *
00027  *
00028  */
00029 #include "plumeinput.h"
00030 #include "plumevars.h"
00031 #include <math.h>
00032 #include <utils/utils.h>
00033 
00034 int plumeconc(Plume_enviro *env , Plume_grid *grid , Plume_options *opt )
00035 {
00036    int ii, jj, nn;
00037    double aa, bb, ucor, *ccor, plugwidth;
00038    double dl, utest, uu, XX, v1, v2;
00039    Plume_river river = *(env->river);
00040    Plume_ocean ocean = *(env->ocean);
00041    Plume_sediment *sedload = env->sed;
00042 //   double erfc(double);
00043 
00044 //---     
00045 // Set the minimum velocity cutoff used in the following calculations
00046 // The value chosen reduces mass balance errors at the edges of the plume
00047 //---
00048    utest  = 0.05*river.u0;
00049    dl     = 0.5*(double)(grid->dx+grid->dy); // Average grid length
00050 
00051 // Allocate memory for the concentration correction 
00052    if( (ccor = (double *) calloc(env->n_grains,sizeof(double))) == NULL )
00053    {
00054        fprintf(stderr," PlumeArray ERROR: memory allocation failed \n");
00055        fprintf(stderr,"                failed on xcl, ycl, dcl \n");
00056        eh_exit(1);
00057    }
00058 
00059    for ( nn=0 ; nn<env->n_grains ; nn++ )
00060    {
00061 
00062       if ( river.Cs[nn] > .001 )
00063       {
00064 
00065          for( ii=0 ; ii<grid->lx ; ii++ )
00066          {
00067 
00068 #ifdef DBG
00069             fprintf( stderr ,
00070                      "\r  PlumeConc  ii = %d, lx = %d, ly = %d, nn = %d, ngrains = %d    " , 
00071                      ii+1 , grid->lx , grid->ly , nn+1 , env->n_grains );
00072 #endif
00073    
00074             //---
00075             // Mass and Velocity Corrections for flow out of the edges
00076             // of a Fjord
00077             //---
00078             if( opt->fjrd && grid->xval[ii] > 0 )
00079             {
00080                aa = sqrt( river.b0/(sqpi*C1*grid->xval[ii]) );
00081                bb = 1/(sqtwo*C1*grid->xval[ii]);
00082    
00083                ucor = 2.*river.u0*aa*aa*erfc(bb*0.5*(grid->ymax-grid->ymin))/((grid->ymax-grid->ymin)*bb);
00084                ccor[nn] = 2.*river.Cs[nn]*aa*aa*erfc(bb*0.5*(grid->ymax-grid->ymin))/((grid->ymax-grid->ymin)*bb);
00085             }
00086             else
00087             {
00088                ucor     = 0;
00089                ccor[nn] = 0;
00090             }
00091    
00092             for ( jj=0 ; jj<grid->ly ; jj++ )
00093             {
00094                // 'zone of flow establishment'
00095                if( grid->dist[ii][jj][2] < plg*river.b0 )
00096                {
00097    
00098                   plugwidth  = -grid->dist[ii][jj][2]/(2.*plg) + river.b0/2.;
00099    
00100                   if( grid->dist[ii][jj][3] < plugwidth )
00101                   {
00102                      grid->ccnc[ii][jj][nn] = river.Cs[nn]+ccor[nn];
00103                      grid->ualb[ii][jj]     = river.u0+ucor;
00104                   }
00105                   else
00106                   {
00107                      v1 = grid->dist[ii][jj][3]
00108                         + 0.5*sqpi*C1*grid->dist[ii][jj][2]
00109                         - river.b0/2.;
00110                      v2 = mx( (sqtwo*C1*grid->dist[ii][jj][2]), 0.01);
00111    
00112                      grid->ccnc[ii][jj][nn] = river.Cs[nn]*exp(-sq(v1/v2))+ccor[nn];
00113                      grid->ualb[ii][jj]     = river.u0*exp(-sq(v1/v2))+ucor;
00114                   }
00115                }
00116                else // 'zone of established flow'
00117                {
00118                   v1 = river.b0/(sqpi*C1*grid->dist[ii][jj][2]);
00119                   v2 = grid->dist[ii][jj][3]/(sqtwo*C1*grid->dist[ii][jj][2]);
00120    
00121                   grid->ccnc[ii][jj][nn] = river.Cs[nn]*sqrt(v1)*exp(-sq(v2))
00122                                          + ccor[nn];
00123                   grid->ualb[ii][jj] = river.u0*sqrt(v1)*exp(-sq(v2))+ucor;
00124                }
00125    
00126                // scale surface concentration by: t = x/u
00127                if ( opt->fjrd )
00128                   uu = (river.u0 + grid->dist[ii][jj][4] + 7.*grid->ualb[ii][jj])/9.;
00129                else
00130                   uu = (river.u0 + grid->dist[ii][jj][4] + 3.*grid->ualb[ii][jj])/5.;
00131                XX = sqrt(sq(grid->dist[ii][jj][2]) + sq(grid->dist[ii][jj][3]));
00132    
00133                // Calculate Non-conservative Concentration
00134                if ( grid->ualb[ii][jj] > utest )
00135                   grid->ncnc[ii][jj][nn] = grid->ccnc[ii][jj][nn]
00136                                          * exp(-sedload[nn].lambda*XX/uu)
00137                                          + ocean.Cw;
00138                else
00139                   grid->ncnc[ii][jj][nn] = ocean.Cw;
00140    
00141                //---
00142                // Calculate Deposit Thickness, scale time by local u and local x
00143                //  C do 1/rho => kg/m^3 m m^3/kg => m (~/dt)
00144                //  Vol/Area = do
00145                //  dt/day => dTOs*u/l  => (#ofs/day)*1/(#ofs/dt) => dt/day
00146                //  m/dt * dt/day => m (~/day)
00147                //---
00148                if (    grid->ncnc[ii][jj][nn] > ocean.Cw
00149                     && grid->ualb[ii][jj] > utest )
00150                {
00151                   grid->deps[ii][jj][nn] = grid->ncnc[ii][jj][nn] 
00152                                *( exp(sedload[nn].lambda*dl/grid->ualb[ii][jj]) - 1. )
00153                                *( river.d0*dTOs*grid->ualb[ii][jj] )
00154                                /( sedload[nn].rho*dl );
00155                }
00156                else
00157                   grid->deps[ii][jj][nn] = 0.0;
00158    
00159                // Calculate the concentration of a conservative tracer (like Salinity)
00160                if( opt->o1 )
00161                   grid->sln[ii][jj] = (ocean.Sw-ocean.So)
00162                                     * (1-grid->ccnc[ii][jj][0]/(river.Cs[0]-ocean.Cw))
00163                                     + ocean.So;
00164             }
00165          }
00166       }
00167    }
00168 
00169    free(ccor);
00170 
00171 #ifdef DBG
00172    fprintf(stderr,"\n");
00173 #endif
00174 
00175    return 0;
00176 
00177 } // end of PlumeConc
00178 
00179 /* 
00180  * estimate sediment/velocity lost out the edges of fjord runs
00181  * 
00182  * let     x = local distance from the river mouth
00183  *      ym = 1/2 width of the fjord (assumes river is at the center)
00184  * 
00185  * define     a = sqrt(bo/(sqrt(pi)*C1*x))
00186  *      b = 1/(sqrt(2)*C1*x)
00187  * 
00188  * then     u = uc*a*(2/sqrt(pi))*exp(-(b*y)^2)
00189  *      
00190  * with     uc = uo*a*exp(-(b*0)^2)          local centerline velocity
00191  * 
00192  * with     erfc(ym) = 1-erf(ym) = 2/sqrt(pi) Integral(ym->inf)[*exp(-y^2)];
00193  * 
00194  * then     erfc(u(ym)) = (1/b)*uo*a^2*erfc(b*ym)
00195  * 
00196  * then multiply by 2 (account for +- sides) and scale to width and
00197  * grid spacing dy/(ymax-ymin)
00198  * 
00199  */

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