Skip to content

AeroBulk is a modern-FORTRAN-based package/library that gathers state-of-the-art aerodynamic bulk formulae algorithms used to compute turbulent air-sea fluxes of momentum, heat and freshwater.

License

Notifications You must be signed in to change notification settings

brodeau/aerobulk

Folders and files

NameName
Last commit message
Last commit date

Latest commit

 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 

Repository files navigation

DOI

AeroBulk is a FORTRAN90-based library and suite of tools (including a C++ interface) that feature state of the art parameterizations to estimate turbulent air-sea fluxes by means of the traditional aerodynamic bulk formulae.

These turbulent fluxes, namely, wind stress, evaporation (latent heat flux) and sensible heat flux, are estimated using the sea surface temperature (bulk or skin), and the near-surface atmospheric surface state: wind speed, air temperature and humidity. If the cool-skin/warm-layer schemes need to be called to estimate the skin temperature, surface downwelling shortwave and longwave radiative fluxes are required.

> Bulk formula and their parameterizations

AeroBulk relies on the following traditional set of bulk formula to estimate turbulent fluxes at the air-sea interface:

$$ \vec{\tau} = \rho\ C_D \ \vec{U}_z \ U_B $$

$$ E = \rho \ C_E \ \big[ q_s - q_z \big] \ U_B $$

$$ Q_L = -L_v \ E $$

$$ Q_H = \rho \ C_H \ C_P \ \big[ \theta_z - T_s \big] \ U_B $$

where $\rho$ is the density of air. The $z$ subscript relates to the reference height above the air-sea interface (generally z=10m). $\vec{U}_z$ is the wind speed vector at the reference height. $U_B$ is the bulk scalar wind speed at the reference height (very close to $|\vec{U}_z|$ in most cases). $\theta_z$ and $q_z$ are the potential temperature and specific humidity of air at the reference height, respectively. $T_s$ and $q_s$ are the absolute (= potential) temperature and specific humidity of air immediately at the air-sea interface (z=0), respectively; if the cool-skin/warm-layer scheme is used, these two are deduced from the skin temperature, otherwise they are deduced from the bulk SST (default).

Any decent level of accuracy from this set of formula can only be achieved through a consistent estimate of the 3 bulk transfer coefficients: $C_D$, $C_E$, and $C_H$ (drag, evaporation, and sensible heat coefficients). In AeroBulk, these bulk coefficients can be estimated thanks to a collection of bulk parameterizations a.k.a bulk algorithms, which relate the value of these coefficients to the near-surface atmospheric stability, the wind speed, and (ideally) the roughness of the sea surface.

The following figure provides a schematic overview on the way turbulent fluxes are computed in AeroBulk:

 

Currently, in AeroBulk, 5 bulk parameterizations are available to compute $C_D$, $C_E$, and $C_H$ used in the bulk formula:

In the COARE and ECMWF algorithms, a cool-skin/warm layer scheme is included and can be activated if the input sea-surface temperature to be used is the bulk SST (usually measured a few tenths of meters below the surface). Activation of these cool-skin/warm layer schemes requires the surface downwelling shortwave and longwave radiative flux components to be provided. Other parameterizations, such as NCAR, are meant to be used with the bulk SST, and do not feature a cool-skin/warm layer scheme.

Beside bulk algorithms, AeroBulk also provides a collection of functions (module mod_phymbl.f90) to accurately estimate relevant atmospheric parameters such as: density of air, different expressions of the humidity of air, viscosity of air, specific humidity at saturation, Obukhov length, bulk Richardson number, wind gustiness, etc...

The focus in AeroBulk is readability, efficiency, and portability towards modern ocean & atmosphere GCMs (Fortran 90, set of modules and a library).

Example of a set of figures generated with one of AeroBulk diagnostic Python scripts:

AeroBulk Approach

Fig 1/ Comparison of the stability correction profiles $\Psi(\zeta)$ as used in 4 different bulk algorithms.

 

Fig 2/ Comparison of the neutral drag (thick lines) and evaporation coefficients (thinner lines) as a function of the neutral wind speed at 10m.

 

> Compilation of AeroBulk

Normally, AeroBulk should compile fine with GNU gfortran (version 7 and above) and Intel ifort. Prior to typing make, make sure that a file or symbolic link named make.macro, which contains architecture specific compilation paths and options, is present in the root directory of AeroBulk. Template examples for make.macro can be found in the arch/ directory.

Compile the library and basic executables:

make

Compile the suite of test / diagnostic executables:

make test

Compile the C++ interface:

make cpp

Using aerobulk routines and functions in your GCM

GCMs written in Fortran

Include compilation flag for the compiler to locate the *.mod module files:

$(FC) $(FFLAGS) -I<path_to_aerobulk>/mod ...

Linking to AeroBulk's library:

... -L<path_to_aerobulk>/lib -laerobulk

GCMs written in C++

Include compilation flag for the compiler to locate the *.hpp header file:

$(CXX) $(CXXLAGS) -I<path_to_aerobulk>/include ...

Linking to AeroBulk's library:

... -L<path_to_aerobulk>/lib -laerobulk_cxx -laerobulk

 

> Giving AeroBulk a first try in interactive "toy mode"

Check that bin/aerobulk_toy.x has been compiled and execute it:

 ./bin/aerobulk_toy.x

You will be interactively prompted for different sea-surface and ABL related parameters, such as SST, air temperature and humidity, wind speed, etc. Then aerobulk_toy.x will compute all the turbulent air-sea fluxes with all the algorithms available (as well as third-party diagnostics of the ABL), and will print it in the form of a summary table.

List of command line options for aerobulk_toy.x:

-p   => Ask for sea-level pressure, otherwise assume 1010 hPa

-r   => Ask for relative humidity rather than specific humidity

-S   => Use the Cool Skin Warm Layer parameterization to compute
        and use the skin temperature instead of the bulk SST
        only in COARE and ECMWF

-N   => Force neutral stability in surface atmospheric layer
        -> will not ask for air temp. at zt, instead will only
        ask for relative humidity at zt and will force the air
        temperature at zt that yields a neutral-stability!

-h   => Show this message

Example of an output obtained with the following setup:

  • zu = 10 m
  • zt = 2 m
  • SST = 22°C
  • tzt = 20°C
  • qzt = 12 g/kg
  • Uzu = 5 m/s
 ==============================================================================================
    Algorithm:      coare3p0  |  coare3p6  |    ncar    |    ecmwf   | andreas
 ==============================================================================================
       C_D     =    1.1954       1.0775       1.2038       1.2862       1.0167      [10^-3]
       C_E     =    1.3345       1.3729       1.3618       1.3143       1.1565      [10^-3]
       C_H     =    1.3345       1.3729       1.2776       1.2635       1.1103      [10^-3]

       z_0     =   4.40936E-05  2.19285E-05  4.49880E-05  6.98835E-05  1.56119E-05  [m]
       u*      =   0.17578      0.16672      0.17348      0.18192       0.1594      [m/s]
       L       =   -20.383      -16.919      -20.494      -24.029      -18.558      [m]
       Ri_bulk =  -3.78706E-02 -3.79537E-02 -3.90686E-02 -3.79799E-02 -3.87826E-02  [-]

                  *** Neutral stability: **
       UN10    =    5.4192       5.4311       5.3396       5.3992       5.3289      [m/s]
       C_D_N   =    1.0521      0.94234       1.0555       1.1353       0.8950      [10^-3]
       C_E_N   =    1.1077       1.1119       1.1241       1.1064       0.9600      [10^-3]
       C_H_N   =    1.1077       1.1119       1.0624       1.0680       0.9260      [10^-3]

 Equ. Charn p. =   1.10033E-02  4.23674E-03  1.15475E-02  1.80029E-02  2.02298E-03

  Wind stress  =    36.198       32.596       35.849       38.860       30.275      [mN/m^2]
  Evaporation  =    3.1059       3.1925       3.1215       3.0538       2.6270      [mm/day]
     QL        =   -88.031      -90.486      -88.473      -86.555      -74.459      [W/m^2]
     QH        =   -17.833      -18.330      -16.730      -16.805      -14.441      [W/m^2]

 

> Computing turbulent fluxes with AeroBulk

AeroBulk can also directly compute the 3 turbulent fluxes by means of the aerobulk_model() routine of module mod_aerobulk (mod_aerobulk.f90):

   PROGRAM TEST_FLUX
       USE mod_aerobulk
       ...
       CALL AEROBULK_MODEL( jt, Nt, calgo, zt, zu, sst, t_zt, hum_zt, U_zu, V_zu, SLP, &
       &                    Qe, Qh, Tau_x, Tau_y, Evap                                 &
       &                   [, Niter=N, rad_sw=Rsw, rad_lw=Rlw] )
       ...
   END PROGRAM TEST_FLUX

INPUT ARGUMENTS:

  • jt : (Sc,int) current time step (to go from 1 to Nt)
  • Nt : (Sc,int) number of time records to go for
  • calgo : (String) algorithm to use (coare3p0/coare3p6/ncar/ecmwf)
  • zt : (Sc,real) height for temperature and spec. hum. of air [m]
  • zu : (Sc,real) height for wind speed (generally 10m) [m]
  • sst : (2D,real) SST [K]
  • t_zt : (2D,real) ABSOLUTE air temperature at zt [K]
  • hum_zt: (2D,real) humidity of air at zt given as specific [kg/kg], or relative [%], or dew-point [K]
  • U_zu : (2D,real) zonal scalar wind speed at 10m [m/s]
  • V_zu : (2D,real) meridional scalar wind speed at 10m [m/s]
  • SLP : (2D,real) sea-level pressure [Pa]

[ OPTIONAL INPUT ARGUMENT: ]

  • Niter : (Sc,int) number of iterations for the bulk algorithms (default is 5)
  • rad_sw : (2D,real) downw. shortwave rad. at surface (>0) [W/m^2]
  • rad_lw : (2D,real) downw. longwave rad. at surface (>0) [W/m^2]

(The presence of rad_sw and rad_sw triggers the use of the Cool-Skin Warm-Layer parameterization with COARE and ECMWF algorithms)

OUTPUT ARGUMENTS:

  • Qe : (2D,real) latent heat flux [W/m^2]
  • Qh : (2D,real) sensible heat flux [W/m^2]
  • Tau_x : (2D,real) zonal wind stress [N/m^2]
  • Tau_y : (2D,real) meridional wind stress [N/m^2]
  • Evap : (2D,real) evaporation [mm/s]

Example of a call, using COARE 3.0 algorithm with cool-skin warm-layer parameterization and 10 iterations:

CALL AEROBULK_MODEL( 1, 24, 'coare3p0', 2., 10., sst, t_zt, hum_zt, U_zu, V_zu, SLP, &
     &               Qe, Qh, Tau_x, Tau_y, Evap,                                     &
     &               Niter=10, rad_sw=Rsw, rad_lw=Rlw T_s=SSST )

 

> Computing transfer coefficients with AeroBulk

In AeroBulk, 5 different routines are available to compute the bulk transfer (a.k.a exchange) coefficients CD, CH and CE. Beside computing the transfer coefficients, these routines adjust air potential temperature and specific humidity from height zt to the reference height (wind) zu. They also return the bulk wind speed, which is the scalar wind speed at height zu with the potential inclusion of a gustiness contribution (in calm and unstable conditions).

> TURB_COARE3p0, transfer coefficients with COARE 3.0 (replace "3p0" by "3p6" to use COARE 3.6)

Use turb_coare3p0() of module mod_blk_coare3p0 (mod_blk_coare3p0.f90). Example of a call:

PROGRAM TEST_COEFF
    USE mod_const
    USE mod_aerobulk, ONLY: AEROBULK_INIT, AEROBULK_BYE
    USE mod_blk_coare3p0
    ...
    !! Initialization and sanity check:
    CALL AEROBULK_INIT( Nt, 'coare3p0', T_s, t_zt, q_zt, Ux_zu, Uy_zy, P )
    
    ...
    
    q_s = 0.98 * q_sat( T_s, P ) ! spec. hum. of saturation at the air-sea interface (z=0) [kg/kg]
    
    ...
    
    CALL TURB_COARE3P0( jt, zt, zu, T_s, t_zt, q_s, q_zt, U_zu, l_use_cool_skin, l_use_warm_layer, &
    &                   Cd, Ch, Ce, t_zu, q_zu, U_blk             &
    &                   [ , Qsw=Rsw, rad_lw=Rlw, slp=P ]          &
    &                   [ , xz0=z0, xu_star=u_s, xL=L ] )
    
    ...
    
    CALL AEROBULK_BYE()
    
END PROGRAM TEST_COEFF

INPUT ARGUMENTS:

  • jt : (Sc,int) current time step (to go from 1 to nitend==Nt)
  • zt : (Sc,real) height for air temperature and humidity [m]
  • zu : (Sc,real) height for wind speed (generally 10m) [m]
  • t_zt: (2D,real) POTENTIAL air temperature at zt [K]
  • q_zt: (2D,real) air spec. humidity of at zt [kg/kg]
  • U_zu: (2D,real) scalar wind speed at zu [m/s]

INPUT and OUTPUT ARGUMENTS:

  • T_s : (2D,real) surface temperature [K]
    • input: bulk SST
    • output: skin temperature or SST (unchanged)
  • q_s : (2D,real) surface satur. spec. humidity at T_s [kg/kg]
    • input: saturation at bulk SST (not needed if skin p. used)
    • output: saturation at skin temp. or at SST (unchanged)

[ OPTIONAL INPUT ARGUMENTS: ]

  • Qsw : (2D,real) net shortw. rad. at surface (>0, after albedo!) [W/m^2]
  • rad_lw : (2D,real) downw. longw. rad. at surface (>0) [W/m^2]
  • slp : (2D,real) sea-level pressure [Pa]

(The presence of these 3 optional input parameters triggers the use of the Cool-Skin Warm-Layer parameterization)

OUTPUT ARGUMENTS:

  • Cd : (2D,real) drag coefficient
  • Ch : (2D,real) sensible heat transfer coefficient
  • Ce : (2D,real) moisture transfer (evaporation) coefficient
  • t_zu : (2D,real) air pot. temperature adjusted at zu [K]
  • q_zu : (2D,real) air spec. humidity adjusted at zu [kg/kg]
  • Ublk : (2D,real) bulk wind speed at 10m [m/s]

[ OPTIONAL OUTPUT ARGUMENTS: ]

  • z0 : (2D,real) roughness length of the sea surface [m]
  • u_s : (2D,real) friction velocity [m/s]
  • L : (2D,real) Obukhov length [m]

> Some Examples

Using COARE 3.6 without the use of the cool-skin & warm-layer schemes, with air temperature and humidity provided at 2m and wind at 10m:

PROGRAM TEST_COEFF
    USE mod_const
    USE mod_aerobulk, ONLY: AEROBULK_INIT, AEROBULK_BYE
    USE mod_blk_coare3p0
    ...
    !! Initialization and sanity check:
    CALL AEROBULK_INIT( Nt, 'coare3p6', Ts, t2, q2, U10, V10, P0,  l_cswl=.FALSE. )
    ...
    DO jt = 1, Nt
        ...
        qs = 0.98 * q_sat( Ts, P0 ) ! spec. hum. of saturation at the air-sea interface (z=0) [kg/kg]
        ...
        CALL TURB_COARE3P6( jt, 2., 10., Ts, t2, qs, q2, SQRT(U10*U10+V10*V10), .false., .false., &
        &                   Cd, Ch, Ce, t10, q10, U_blk )
        ...
    END DO
    ...
END PROGRAM TEST_COEFF

In this case, Ts and qs, the surface temperature and saturation specific humidity won't be modified. The actual value of qs must be provided as an input.

Now the same but using the cool-skin & warm-layer schemes:

PROGRAM TEST_COEFF
    USE mod_const
    USE mod_aerobulk, ONLY: AEROBULK_INIT, AEROBULK_BYE
    USE mod_blk_coare3p0
    ...
    !! Initialization and sanity check:
    CALL AEROBULK_INIT( Nt, 'coare3p6', Ts, t2, q2, U10, V10, P0,  l_cswl=.TRUE. )
    ...
    DO jt = 1, Nt
        ...
        CALL TURB_COARE3P6( jt, 2., 10., Ts, t2, qs, q2, SQRT(U10*U10+V10*V10), .true., .true., &
        &                   Cd, Ch, Ce, t10, q10, U_blk,                       &
        &                   Qsw=Rsw, rad_lw=Rlw, slp=P0,                       &
        &                   isecday_utc=50400, plong=xlongitudes )
        ...
    END DO
    ...
END PROGRAM TEST_COEFF

We provide arrays of the NET solar flux to the ocean (Rsw), downwelling longwave (infrared) flux (Rlw), sea-level atmospheric pressure (P0), and longitudes (xlongitudes); as well as the current UTC time as seconds since 00:00 (here 2PM = 50400). Note: Ts is the bulk SST as an input and is updated to the skin temperature as an output! Value passed to qs as an input won't be used, however, as an output, qs is the saturation specific humidity at temperature Ts!

 

> Computing atmospheric state variables with AeroBulk

A selection of useful functions to estimate some atmospheric state variables in the marine boundary layer are available in the module mod_phymbl (mod_phymbl.f90).

Example for computing SSQ of Eq.(1) out of the SST and the SLP:

PROGRAM TEST_PHYMBL
    USE mod_const
    USE mod_phymbl
    ...

    SSQ(:,:) = 0.98*q_sat(SST, SLP)
    ...
END PROGRAM TEST_PHYMBL

 

> Acknowledging AeroBulk

To acknowledge/reference AeroBulk in your scientific work, please cite: Brodeau, L., B. Barnier, S. Gulev, and C. Woods, 2016: Climatologically significant effects of some approximations in the bulk parameterizations of turbulent air-sea fluxes. J. Phys. Oceanogr., 47 (1), 5–28, 10.1175/JPO-D-16-0169.1. [ DOI ]

About

AeroBulk is a modern-FORTRAN-based package/library that gathers state-of-the-art aerodynamic bulk formulae algorithms used to compute turbulent air-sea fluxes of momentum, heat and freshwater.

Topics

Resources

License

Stars

Watchers

Forks

Contributors 3

  •  
  •  
  •