/Users/huttone/Devel/sedflux-new/sedflux/trunk/ew/squall/squall.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 //#define DISABLE_WATCH_POINTS
00022 //#define DISABLE_CHECKS
00023 
00024 #include <stdio.h>
00025 #include <string.h>
00026 #include <glib.h>
00027 
00028 #include "squall.h"
00029 
00030 typedef struct
00031 {
00032    gsize count;
00033    double lower_edge;
00034    double upper_edge;
00035    gpointer data;
00036 }
00037 Eh_histogram;
00038 
00039 typedef struct
00040 {
00041    double f;
00042    gssize ind;
00043 }
00044 Hist_data;
00045 
00046 Eh_histogram **eh_create_histogram( double dx , gssize n );
00047 void eh_destroy_histogram( Eh_histogram **h );
00048 Eh_histogram **make_depth_histogram( double *x , int n_x , double dz , gssize n_z );
00049 
00050 double get_weibull_deposition_rate( double x , double alpha , double beta );
00051 double get_threshold_depth( double wave_length , double wave_height ,
00052                             double wave_period , double grain_size );
00053 
00063 gboolean squall( Sed_cube p , double time_step_in_years )
00064 {
00065    int i, j, k, n, n_grains = sed_sediment_env_n_types();
00066    int *zone;
00067    int i_w, i_s, i_c, i_l, i_m, i_0;
00068    double h_w, h_c;
00069    double c_v = SQUALL_DEFAULT_C_V;
00070    double *e       = eh_new0( double , sed_cube_n_y(p) );
00071    double *g       = eh_new0( double , sed_cube_n_y(p) );
00072    double *f       = eh_new0( double , n_grains );
00073    double *is_moveable = eh_new0( double , n_grains );
00074    double depth, h, dep, dep_fraction;
00075    double distance_to_h_w, total;
00076    double *fraction_total;
00077    double *threshold_depth = eh_new( double , n_grains );
00078    double **dep_thickness  = eh_new_2( double , sed_cube_n_y(p) , n_grains );
00079    double *profile_depth;
00080    double *total_erosion   = eh_new( double , n_grains );
00081    double dz;
00082    int n_z_bins;
00083    gssize dep_column;
00084    Eh_histogram **hist;
00085    gboolean back_barrier_is_on = FALSE;
00086    Sed_cell add_cell;
00087    Sed_cell top_cell           = sed_cell_new( n_grains );
00088    Sed_cell lag_cell           = sed_cell_new( n_grains );
00089    Sed_cell removed_cell       = sed_cell_new( n_grains );
00090    Sed_cell erosion_cell       = sed_cell_new( n_grains );
00091    Sed_cell* dep_cell = sed_cell_list_new( sed_cube_n_y(p) , n_grains );
00092    Sed_cell bb_cell, sf_cell;
00093    double mass_in, mass_out;
00094 
00095    h_c = sed_cube_wave_height( p );
00096    h_w = get_deep_water_wave_base( h_c );
00097 
00098 //eh_make_note( h_w=50 );
00099 
00100 h_c = 5;
00101 
00102    //---
00103    // calculate the locations of the various domains for erosion and
00104    // deposition.  there will be a maximum of 6 zones:
00105    //  (1) continental
00106    //  (2) lagoon
00107    //  (3) barrier island
00108    //  (4) shoreface
00109    //  (5) breaking
00110    //  (6) deep water
00111    // however, zones (2) and (3) will only be present if a barrier island
00112    // has formed.  this is not always the case.  thus, there may only be
00113    // zones (1), (4), (5), and (6).
00114    //---
00115    zone = get_zone_boundaries( p , h_w , h_c );
00116 
00117    i_m = zone[0];
00118    i_l = zone[1];
00119    i_c = zone[2];
00120    i_s = zone[3];
00121    i_w = zone[4];
00122 
00123 //   eh_make_note( h_w = h_shelf );
00124 //   eh_make_note( i_w = zone[5] );
00125 
00126    distance_to_h_w = ( i_w - i_s )*sed_cube_y_res( p );
00127 
00128 /*
00129    h_c = -sed_get_depth_from_profile( p , i_c );
00130    h_w =  sed_get_depth_from_profile( p , i_w );
00131 */
00132 /*
00133    for ( i=i_c ; i<i_w ; i++ )
00134    {
00135       if ( sed_get_depth_from_profile( p , i )>h_w )
00136       {
00137          eh_watch_int( i );
00138          eh_watch_int( i_w );
00139          eh_watch_dbl( h_w );
00140          eh_watch_int( sed_get_profile_river_mouth( p ) );
00141          eh_watch_dbl( sed_get_depth_from_profile( p , i ) );
00142       }
00143       if ( -sed_get_depth_from_profile( p , i )>h_c )
00144       {
00145          eh_watch_int( i );
00146          eh_watch_int( i_c );
00147          eh_watch_dbl( h_c );
00148          eh_watch_int( sed_get_profile_river_mouth( p ) );
00149          eh_watch_dbl( sed_get_depth_from_profile( p , i ) );
00150       }
00151    }
00152 */
00153    eh_free( zone );
00154 
00155    //---
00156    // phase 1:
00157    //
00158    // calculate erosion rates.
00159    //---
00160    if ( back_barrier_is_on )
00161       i = i_c;
00162    else
00163       i = i_s;
00164    for ( ; i<i_w ; i++ )
00165    {
00166 //      e[i] = get_erosion_rate_from_profile( p , i , i_c , i_w )
00167 //      e[i] = get_erosion_rate_from_profile( p , i , h_c , h_w )
00168 //           * time_step_in_years;
00169       h    = sed_cube_water_depth( p , 0 , i );
00170       g[i] = pow( (h-h_w) / (-h_c-h_w) , 3. );
00171 /*
00172 eh_require( g[i]>=0 );
00173 eh_require( g[i]<=1 );
00174 if ( g[i]<0 || g[i]>1 )
00175 {
00176    eh_watch_dbl(  h   );
00177    eh_watch_dbl(  h_w );
00178    eh_watch_dbl( -h_c );
00179 }
00180 */
00181 
00182       e[i] = (h_w/10.)*SQUALL_DEFAULT_C_E*g[i]*time_step_in_years;
00183 //      e[i] = SQUALL_DEFAULT_C_E*g[i]*time_step_in_years;
00184       if ( e[i]<0 )
00185          e[i] = 0;
00186    }
00187 /*
00188    if ( back_barrier_is_on )
00189       i = i_c;
00190    else
00191       i = i_s;
00192    for ( ; i<i_shelf ; i++ )
00193    {
00194       e_shelf[i] = ( pow((h-h_shelf)/(-h_c-h_shelf),3.) - g[i] )
00195                  * SQUALL_DEFAULT_C_E*time_step_in_years;
00196       e_shelf[i] = 0.;
00197    }
00198 eh_message( "e_shelf set to zero." );
00199 */
00200 
00201    //---
00202    // remove the eroded sediment and mix it together.
00203    //---
00204    for ( i=i_c ; i<i_w ; i++ )
00205    {
00206       get_moveable_grains( sed_cube_water_depth( p , 0 , i ) ,
00207                            sed_cube_wave_height( p )        ,
00208                            sed_cube_wave_period( p )        ,
00209                            NULL                                 ,
00210                            is_moveable );
00211 
00212       sed_column_extract_top( sed_cube_col(p,i) , e[i] , lag_cell );
00213 
00214       sed_cell_separate_fraction( lag_cell    ,
00215                                   is_moveable ,
00216                                   removed_cell );
00217 
00218       sed_column_add_cell( sed_cube_col(p,i) , lag_cell );
00219       sed_cell_add( erosion_cell , removed_cell );
00220    }
00221    sed_cell_set_facies( erosion_cell , S_FACIES_WAVE );
00222 
00223    //---
00224    // deposit sediment in the backbarrier (phase 2) and shoreface (phase 3).
00225    // the total eroded sediment is divided between the backbarrier and the
00226    // shoreface according to c_v.
00227    //---
00228    if ( i_m<=0 || !back_barrier_is_on )
00229       c_v = 0;
00230    bb_cell = sed_cell_dup( erosion_cell );
00231    sf_cell = sed_cell_dup( erosion_cell );
00232    sed_cell_resize( bb_cell ,      c_v*sed_cell_size(erosion_cell) );
00233    sed_cell_resize( sf_cell , (1.-c_v)*sed_cell_size(erosion_cell) );
00234 
00235    //---
00236    // phase 2:
00237    //
00238    // deposit sediment in the backbarrier.  begin on the barrier island, x_c
00239    // and more landward to x_m.
00240    //----
00241 
00242    if ( i_m > 0 && back_barrier_is_on )
00243    {
00244       double* grain_size = sed_sediment_property( NULL , &sed_type_grain_size_in_meters );
00245 //      for ( i=i_c ; i>i_m && sed_get_cell_thickness(bb_cell)>1e-5 ; i-- )
00246       for ( i=i_c-1 ; i>i_m && sed_cell_size(bb_cell)>1e-5 ; i-- )
00247       {
00248          depth = sed_cube_water_depth( p , 0 , i );
00249    
00250          for ( n=0 ; n<n_grains ; n++ )
00251          {
00252             memset(f,0,sizeof(double)*n_grains);
00253             f[n] = 1;
00254    
00255             h = get_travel_dist( grain_size[n] ,
00256                                  depth         ,
00257                                  sed_cube_y_res(p) );
00258             dep = sed_cell_size( bb_cell )
00259                 * sed_cell_fraction( bb_cell , n )
00260                 * sed_cube_y_res( p )
00261                 / h;
00262 
00263             sed_cell_move( bb_cell , dep_cell[i] , f , dep );
00264          }
00265       }
00266       eh_free( grain_size );
00267    }
00268 
00269    //---
00270    // phase 3:
00271    // 
00272    // deposit sediment on the shoreface.  any sediment that was able to be
00273    // deposited in the backbarrier during phase 2 is added to the 
00274    // shoreface sediment.
00275    //---
00276    sed_cell_add( sf_cell , bb_cell );
00277 
00278 //   for ( i=i_c+1 ; i<p->size && sed_get_cell_thickness(sf_cell)>1e-5 ; i++ )
00279 //   for ( i=i_c ; i<p->size && sed_get_cell_thickness(sf_cell)>1e-5 ; i++ )
00280    if ( back_barrier_is_on )
00281       i_0 = i_c;
00282    else
00283       i_0 = i_s;
00284 
00285    for ( mass_in=0,n=0 ; n<n_grains ; n++ )
00286    {
00287       total_erosion[n] = sed_cell_size( sf_cell )
00288                        * sed_cell_fraction( sf_cell , n );
00289       mass_in += total_erosion[n];
00290    }
00291    mass_out = 0;
00292 
00293 {
00294    double* grain_size = sed_sediment_property( NULL , &sed_type_grain_size_in_meters );
00295    for ( k=0 ; k<1 ; k++ )
00296    {
00297 
00298       //---
00299       // divide the profile up into depth bins.  the histogram includes an 
00300       // index to the column at that depth, and the fraction of the bin that
00301       // the column occupies.
00302       //---
00303       dz = sed_cube_z_res(p);
00304       profile_depth = eh_new( double , sed_cube_n_y(p) );
00305       for ( i=0 ; i<sed_cube_n_y(p) ; i++ )
00306          profile_depth[i] = sed_cube_water_depth(p,0,i);
00307       n_z_bins = (int)(profile_depth[sed_cube_n_y(p)-1]/dz);
00308       hist = make_depth_histogram( profile_depth , sed_cube_n_y(p) , dz , n_z_bins );
00309       eh_free( profile_depth );
00310    
00311       for ( n=0 ; n<n_grains ; n++ )
00312       {
00313          threshold_depth[n] = get_threshold_depth(
00314                                  sed_cube_wave_length(p) ,
00315                                  sed_cube_wave_height(p) ,
00316                                  sed_cube_wave_period(p) ,
00317                                  grain_size[n] );
00318       }
00319    
00320       for ( i=0 ; i<sed_cube_n_y(p) ; i++ )
00321          for ( n=0 ; n<n_grains ; n++ )
00322             dep_thickness[i][n] = 0;
00323    
00324       fraction_total = eh_new( double , n_grains );
00325       for ( n=0 ; n<n_grains ; n++ )
00326          fraction_total[n] = 0.;
00327    
00328       for ( i=0 ; i<n_z_bins ; i++ )
00329          for ( n=0 ; n<n_grains ; n++ )
00330             if ( hist[i]->count > 0 )
00331                fraction_total[n] += get_weibull_deposition_rate(
00332                                        i*dz ,
00333                                        2 ,
00334                                        threshold_depth[n] )*dz;
00335    
00336       for ( i=0 ; i<n_z_bins && sed_cell_size(sf_cell)>1e-5 ; i++ )
00337       {
00338          depth = i*dz;
00339    
00340          for ( n=0 ; n<n_grains ; n++ )
00341          {
00342             memset(f,0,sizeof(double)*n_grains);
00343             f[n] = 1;
00344    
00345             if ( depth>0 && hist[i]->count>0 )
00346             {
00347                dep_fraction = -exp( -pow((depth+dz)/threshold_depth[n],2) )
00348                             +  exp( -pow(depth/threshold_depth[n],2) );
00349    /*
00350                dep_fraction = get_weibull_deposition_rate(
00351                                  depth ,
00352                                  (n==0)?1:2 ,
00353                                  threshold_depth[n] )*dz;
00354    */
00355             }
00356             else
00357                dep_fraction = 0;
00358    
00359    //         if ( dep > depth+h_c )
00360    //            dep = depth+h_c;
00361    
00362             for ( total=0,j=0 ; j<hist[i]->count ; j++ )
00363                total += ((Hist_data*)(hist[i]->data))[j].f;
00364    
00365             for ( j=0 ; j<hist[i]->count ; j++ )
00366             {
00367    //eh_require( ((Hist_data*)(hist[i]->data))[j].f<=1.0001 );
00368    //eh_require( ((Hist_data*)(hist[i]->data))[j].f>-0.0001 );
00369                dep_column  = ((Hist_data*)(hist[i]->data))[j].ind;
00370                dep        = dep_fraction
00371                           / fraction_total[n]
00372                           * total_erosion[n]
00373                           * ((Hist_data*)(hist[i]->data))[j].f
00374                           / total
00375                           / 1.;
00376    
00377    //            sed_move_sediment_to_cell( sf_cell , dep_cell[dep_column] , f ,
00378    //                                       dep     , n_grains );
00379    
00380                dep_thickness[dep_column][n] += dep;
00381                mass_out += dep;
00382             }
00383          }
00384       }
00385    
00386       eh_destroy_histogram( hist );
00387    
00388       //---
00389       // add the deposited sediment back to the profile.
00390       //---
00391       add_cell = sed_cell_new( n_grains );
00392       sed_cell_set_age( add_cell , sed_cube_age_in_years(p) );
00393       sed_cell_set_facies( add_cell , S_FACIES_WAVE );
00394       mass_out = 0;
00395    //   for ( i=0 ; i<p->size ; i++ )
00396       for ( i=1 ; i<sed_cube_n_y(p) ; i++ )
00397       {
00398    if ( sed_cube_water_depth( p , 0 , i ) > 0 )
00399    {
00400          sed_cell_resize( add_cell , 0. );
00401 //         sed_add_vector_to_cell( add_cell , dep_thickness[i] , n_grains );
00402          sed_cell_add_amount( add_cell , dep_thickness[i] );
00403          if ( sed_cell_size(add_cell) > sed_cube_water_depth(p,0,i) )
00404             sed_cell_resize( add_cell , sed_cube_water_depth(p,0,i) );
00405    if (   sed_cube_water_depth(p,0,i)
00406         - sed_cell_size(add_cell)
00407         - sed_cube_water_depth(p,0,i-1) < 0. )
00408       sed_cell_resize( add_cell ,
00409                          sed_cube_water_depth(p,0,i)
00410                        - sed_cube_water_depth(p,0,i-1)
00411                        - .0 );
00412          if ( sed_cell_size( add_cell ) > 0 )
00413          {
00414             mass_out += sed_cell_size( add_cell );
00415             sed_column_add_cell( sed_cube_col(p,i) , add_cell );
00416          }
00417    }
00418       }
00419    
00420       sed_cell_destroy( add_cell );
00421       eh_free( fraction_total );
00422    
00423    }
00424 
00425    eh_free( grain_size );
00426 }
00427 
00428    eh_free_2( dep_thickness );
00429 
00430 //eh_watch_dbl( mass_in );
00431 //eh_watch_dbl( mass_out );
00432 
00433    sed_cell_destroy( sf_cell            );
00434    sed_cell_destroy( bb_cell            );
00435    sed_cell_destroy( top_cell           );
00436    sed_cell_destroy( lag_cell           );
00437    sed_cell_destroy( removed_cell       );
00438    sed_cell_destroy( erosion_cell       );
00439 
00440    sed_cell_list_destroy( dep_cell );
00441 
00442    eh_free( is_moveable );
00443    eh_free( total_erosion );
00444    eh_free( e       );
00445    eh_free( g       );
00446    eh_free( f );
00447    eh_free( threshold_depth );
00448 
00449    return TRUE;
00450 }
00451 
00452 Eh_histogram **eh_create_histogram( double dx , gssize n )
00453 {
00454    gsize i;
00455    Eh_histogram **h;
00456 
00457    h = eh_new( Eh_histogram* , n+1 );
00458    for ( i=0 ; i<n ; i++ )
00459    {
00460       h[i] = eh_new( Eh_histogram , 1 );
00461       h[i]->count = 0;
00462       h[i]->lower_edge = i*dx;
00463       h[i]->upper_edge = (i+1)*dx;
00464       h[i]->data = NULL;
00465    }
00466    h[n] = NULL;
00467 
00468    return h;
00469 }
00470 
00471 void eh_destroy_histogram( Eh_histogram **h )
00472 {
00473    gssize i;
00474 
00475    for ( i=0 ; h[i] ; i++ )
00476    {
00477       eh_free( h[i]->data );
00478       eh_free( h[i] );
00479    }
00480    eh_free( h );
00481 }
00482 
00483 void add_to_hist( Eh_histogram **hist , gssize bin , double f , gssize ind );
00484 
00485 Eh_histogram **make_depth_histogram( double *x , int n_x , double dz , gssize n_z )
00486 {
00487    double dx = 50;
00488    gssize i, j;
00489    gssize lower_z_bin, upper_z_bin;
00490    double lower_edge, upper_edge;
00491    double lower_x, upper_x;
00492    double f;
00493    double slope;
00494    Eh_histogram **hist;
00495 
00496    eh_require( n_z>0 );
00497 
00498    hist = eh_create_histogram( dz , n_z );
00499 
00500    for ( i=1 ; i<n_x-1 ; i++ )
00501    {
00502 
00503       for ( j=0 ; j<1 ; j++ )
00504       {
00505 
00506          lower_x = x[i];
00507          upper_x = x[i+1];
00508 
00509          if ( j==0 )
00510          {
00511             lower_x = .5*(x[i] + x[i-1]);
00512             upper_x = x[i];
00513          }
00514          else
00515          {
00516             lower_x = x[i];
00517             upper_x = .5*(x[i+1] + x[i]);
00518          }
00519 
00520          lower_x = .5*(x[i] + x[i-1]);
00521          upper_x = .5*(x[i+1] + x[i]);
00522 
00523          if ( upper_x<lower_x )
00524             swap_dbl( upper_x , lower_x );
00525 
00526          lower_z_bin = floor(lower_x/dz);
00527          upper_z_bin = floor(upper_x/dz);
00528 
00529          if ( fabs(upper_x-upper_z_bin*dz) < 1e-5 )
00530             upper_z_bin--;
00531          if ( upper_z_bin<lower_z_bin )
00532             swap_int( upper_z_bin , lower_z_bin );
00533 
00534          if ( upper_z_bin<n_z && lower_z_bin>=0 )
00535          {
00536 
00537             lower_edge = lower_z_bin*dz;
00538             upper_edge = (upper_z_bin+1)*dz;
00539 
00540             if ( lower_z_bin == upper_z_bin )
00541             {
00542                f = sqrt( dx*dx + pow(upper_x-lower_x,2.) );
00543                f = dx;
00544                f = 1;
00545                f = (upper_x-lower_x)/dz;
00546                add_to_hist( hist , lower_z_bin , fabs(f) , i );
00547             }
00548             else
00549             {
00550                slope = atan(dx/(upper_x-lower_x));
00551 
00552                f = (lower_edge+dz-lower_x)/sin(slope);
00553                f = (lower_edge+dz-lower_x)/tan(slope);
00554                f = 1;
00555                f = ((lower_edge+dz)-lower_x)/dz;
00556 
00557                add_to_hist( hist , lower_z_bin , fabs(f) , i );
00558                for ( j=lower_z_bin+1 ; j<upper_z_bin ; j++ )
00559                {
00560                   f = dz/sin(slope);
00561                   f = dz/tan(slope);
00562                   f = 1;
00563                   add_to_hist( hist , j , f , i );
00564                }
00565                f = (upper_x-(upper_edge-dz))/sin(slope);
00566                f = (upper_x-(upper_edge-dz))/tan(slope);
00567                f = 1;
00568                f = (upper_x - (upper_edge-dz))/dz;
00569                add_to_hist( hist , upper_z_bin , fabs(f) , i );
00570 //               add_to_hist( hist , upper_z_bin , f/2. , i+1 );
00571             }
00572          }
00573       }
00574 
00575    }
00576 
00577    return hist;
00578 }
00579 
00580 void add_to_hist( Eh_histogram **hist , gssize bin , double f , gssize ind )
00581 {
00582    Eh_histogram *this_bin = hist[bin];
00583    
00584 if ( this_bin->count > 1000 )
00585 {
00586    eh_watch_int( bin );
00587    eh_watch_int( this_bin->count );
00588    eh_watch_dbl( this_bin->lower_edge );
00589    eh_watch_dbl( this_bin->upper_edge );
00590 }
00591    this_bin->count = this_bin->count + 1;
00592 if ( this_bin->count > 1000 )
00593 {
00594    eh_watch_int( bin );
00595    eh_watch_int( this_bin->count );
00596    eh_watch_dbl( this_bin->lower_edge );
00597    eh_watch_dbl( this_bin->upper_edge );
00598 }
00599    this_bin->data  = g_renew( Hist_data , this_bin->data , this_bin->count );
00600    ((Hist_data*)(this_bin->data))[this_bin->count-1].f   = f;
00601    ((Hist_data*)(this_bin->data))[this_bin->count-1].ind = ind;
00602 }
00603 
00604 
00605 double get_weibull_deposition_rate( double x , double alpha , double beta )
00606 {
00607    if ( fabs(alpha-1.) < 1e-5 )
00608       return exp( -x/beta )/beta;
00609    else
00610       return alpha*pow(x,alpha-1)/pow(beta,alpha)*exp(-pow(x/beta,alpha));
00611 }
00612 
00622 double get_erosion_rate_from_profile( Sed_cube p , int i ,
00623                                       double h_c  , double h_w )
00624 {
00625 //   double alpha_sf    = get_shoreface_slope( p , i_s+1 , i_w );
00626    double alpha_sf    = 1;
00627    double h           = sed_cube_water_depth( p , 0 , i   );
00628    double h_wave      = sed_cube_wave_height( p );
00629    double h_wave_fair = h_wave;
00630 
00631    return get_erosion_rate( h , h_w , h_c , alpha_sf , h_wave , h_wave_fair );
00632 }
00633 
00642 double get_shoreface_slope( Sed_cube p , int i_c , int i_w )
00643 {
00644    double rise = sed_cube_water_depth( p , 0 , i_w )
00645                - sed_cube_water_depth( p , 0 , i_c );
00646    double run  = sed_cube_col_y( p,i_w ) - sed_cube_col_y( p,i_c );
00647 
00648    return rise/run;
00649 }
00650 
00662 double get_erosion_rate( double h        , double h_w    , double h_c ,
00663                          double alpha_sf , double h_wave , double h_wave_fair )
00664 {
00665    double c_e = SQUALL_DEFAULT_C_E;
00666    double c_d = get_coastal_dissipation( alpha_sf );
00667    double g   = get_erosion_efficiency( h , h_w , h_c );
00668 c_d = 1;
00669 /*
00670 eh_watch_dbl( c_e );
00671 eh_watch_dbl( c_d );
00672 eh_watch_dbl( c_w );
00673 eh_watch_dbl( g   );
00674 */
00675 //   return c_e*c_d*c_w*g;
00676    return c_e*g;
00677 }
00678 
00687 double get_erosion_efficiency( double h , double h_w , double h_c )
00688 {
00689    double m = SQUALL_DEFAULT_M;
00690    double g = pow( (h-h_w)/(h_c-h_w) , m );
00691 //   double g = pow( (h_w-h)/(h_w-h_c) , m );
00692 
00693    eh_require( g>=0 );
00694    eh_require( g<=1 );
00695 
00696    if ( g<0 || g>1 )
00697    {
00698       eh_watch_dbl( g   );
00699       eh_watch_dbl( h   );
00700       eh_watch_dbl( h_w );
00701       eh_watch_dbl( h_c );
00702    }
00703 
00704    return g;
00705 }
00706 
00714 double get_coastal_dissipation( double alpha_sf )
00715 {
00716    double alpha_ref = SQUALL_DEFAULT_ALPHA_REF;
00717    return 1+(alpha_sf-alpha_ref)/alpha_ref;
00718 }
00719 
00727 double get_coastal_wave_energy( double wave_height_actual ,
00728                                 double wave_height_fair )
00729 {
00730    return wave_height_actual/wave_height_fair;
00731 }
00732 
00741 double get_wave_base( double wave_length )
00742 {
00743    return wave_length/2.;
00744 }
00745 
00755 double get_breaking_wave_base( double wave_height )
00756 {
00757    return get_wave_base( wave_height*7. );
00758 }
00759 
00769 double get_deep_water_wave_base( double wave_height )
00770 {
00771    return get_wave_base(wave_height*25.);
00772 }
00773 
00792 int *get_zone_boundaries( Sed_cube p , double h_w , double h_c )
00793 {
00794    int i;
00795    int *zone = eh_new( int , 5 );
00796 
00797    for ( i=sed_cube_n_y(p)-1 ; i>0 && sed_cube_water_depth(p,0,i)>h_w ; i-- );
00798    zone[4] = i;
00799 
00800    for ( i=zone[4] ; i>0 && sed_cube_water_depth(p,0,i)>0 ; i-- );
00801    zone[3] = i;
00802 
00803    for ( i=zone[3] ;
00804             i>0
00805          && - sed_cube_water_depth(p,0,i)<h_c
00806          && (  sed_cube_col_y(p,zone[3])
00807              - sed_cube_col_y(p,i) < 1500. ) ;
00808          i-- );
00809    zone[2] = i;
00810 
00811    for ( i=zone[2] ; i>0 && sed_cube_water_depth(p,0,i)<0 ; i-- );
00812    zone[1] = i;
00813 
00814    for ( i=zone[1] ; i>0 && sed_cube_water_depth(p,0,i)>0 ; i-- );
00815    zone[0] = i;
00816 
00817    return zone;
00818 }
00819 
00832 double get_travel_dist( double grain_size_in_m , double depth , double dx )
00833 {
00834    double a   = SQUALL_DEFAULT_A;
00835    double z_l = SQUALL_DEFAULT_Z_L;
00836    double h_star = get_non_dim_travel_dist( grain_size_in_m );
00837    double z = -depth / z_l;
00838    double h;
00839 
00840    if ( dx!=50 )
00841       h_star = dx / ( 1 - pow( 1-50/h_star , dx/50 ) );
00842 
00843    h = h_star*( 1.+exp(a*z) );
00844 
00845    return h;
00846 }
00847 
00857 double get_non_dim_travel_dist( double grain_size_in_m )
00858 {
00859    double grain_size_ref = .125e-3;
00860    double h_star;
00861    double c_h = SQUALL_DEFAULT_C_H;
00862 
00863    if ( grain_size_in_m>grain_size_ref )
00864       h_star = c_h*(110. + 590.*pow( grain_size_ref/grain_size_in_m , 2.5 ) );
00865    else
00866       h_star = c_h*(500. + 200.*pow( grain_size_ref/grain_size_in_m , 0.6 ) );
00867 
00868    return h_star;
00869 }
00870 
00893 double *get_moveable_grains( double water_depth ,
00894                              double wave_height ,
00895                              double wave_period ,
00896                              Sed_sediment sed   ,
00897                              double *is_moveable )
00898 {
00899    double wave_length = 25.*wave_height;
00900    double u;
00901    double d;
00902    gssize n_grains;
00903 
00904    if ( sed )
00905       n_grains = sed_sediment_n_types( sed );
00906    else
00907       n_grains = sed_sediment_env_n_types( );
00908 
00909 //eh_make_note( water_depth*=5 );
00910 /*
00911    eh_require( wave_height>0 );
00912    eh_require( wave_period>0 );
00913    eh_require( water_depth>0 );
00914 */
00915 
00916    if ( !is_moveable )
00917       is_moveable = eh_new0( double , n_grains );
00918 
00919    if ( wave_height<= 0 || wave_period <= 0 || water_depth <= 0 )
00920    {
00921       gssize n;
00922       for ( n=0 ; n<n_grains ; n++ )
00923          is_moveable[n] = FALSE;
00924       return is_moveable;
00925    }
00926 
00927    u = M_PI*wave_height/(wave_period*sinh( 2*M_PI*water_depth/wave_length ));
00928 
00929    d = pow(   sed_rho_sea_water()
00930             / ((sed_rho_quartz()-sed_rho_sea_water())
00931             * sed_gravity()*.21) ,
00932             2. )
00933      * pow( u , 3 )
00934      * M_PI / wave_period;
00935 
00936    {
00937       double* grain_size = sed_sediment_property( NULL , &sed_type_grain_size_in_meters );
00938       gssize n;
00939 
00940       for ( n=0 ; n<n_grains ; n++ )
00941       {
00942          if ( grain_size[n] < d )
00943             is_moveable[n] = TRUE;
00944          else
00945             is_moveable[n] = FALSE;
00946 //if ( n!=0 )
00947 //   is_moveable[n] = 1;
00948       }
00949 
00950       eh_free( grain_size );
00951    }
00952    
00953    return is_moveable;
00954 }
00955 
00956 #include <math.h>
00957 
00958 double get_threshold_depth( double wave_length , double wave_height ,
00959                             double wave_period , double grain_size )
00960 {
00961    double rho   = sed_rho_sea_water();
00962    double rho_s = sed_rho_quartz();
00963    double g     = sed_gravity();
00964    double c_sq  = pow( rho/( .21*g*(rho_s-rho) ), 2.);
00965 
00966    return   wave_length/(2*M_PI)
00967           * asinh(   M_PI*wave_height/wave_period
00968                    * pow( M_PI*c_sq/(wave_period*grain_size) , 1./3 ) );
00969 }
00970 

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