source: flexpart.git/src/getrb.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 2.0 KB
Line 
1subroutine getrb(nc,ustar,nyl,diffh2o,reldiff,rb)
2  !                 i    i    i     i       i    o
3  !*****************************************************************************
4  !                                                                            *
5  !  Calculation of the quasilaminar sublayer resistance to dry deposition.    *
6  !                                                                            *
7  !      AUTHOR: Andreas Stohl, 20 May 1995                                    *
8  !                                                                            *
9  !*****************************************************************************
10  !                                                                            *
11  ! Variables:                                                                 *
12  ! rb(ncmax)       sublayer resistance                                        *
13  ! schmidt         Schmidt number                                             *
14  ! ustar [m/s]     friction velocity                                          *
15  ! diffh20 [m2/s]  diffusivity of water vapor in air                          *
16  ! reldiff         diffusivity relative to H2O                                *
17  !                                                                            *
18  ! Constants:                                                                 *
19  ! karman          von Karman constant                                        *
20  ! pr              Prandtl number                                             *
21  !                                                                            *
22  !*****************************************************************************
23
24  use par_mod
25
26  implicit none
27
28  real :: ustar,diffh2o,rb(maxspec),schmidt,nyl
29  real :: reldiff(maxspec)
30  integer :: ic,nc
31  real,parameter :: pr=0.72
32
33  do ic=1,nc
34    if (reldiff(ic).gt.0.) then
35      schmidt=nyl/diffh2o*reldiff(ic)
36      rb(ic)=2.0*(schmidt/pr)**0.67/(karman*ustar)
37    endif
38  end do
39
40end subroutine getrb
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG