source: codes/icosagcm/trunk/src/dcmip/terminator.f90

Last change on this file was 548, checked in by dubos, 7 years ago

trunk : reorganize source tree

File size: 4.0 KB
Line 
1!================================================================================================
2! This is the "toy" chemistry module.
3!================================================================================================
4
5MODULE Terminator
6
7implicit none
8private
9save
10
11public :: initial_value_Terminator       ! initialize cl and cl2
12public :: tendency_Terminator             ! interface to tendency computation
13
14!integer, parameter :: 8 = selected_real_kind (12)
15
16real(8), parameter :: cly_constant = 4.e-6_8
17real(8), parameter :: pi = 3.1415926535897932384626433832795028841971693993751058209749445923078164_8
18real(8), parameter :: half_pi = pi*0.5_8
19real(8), parameter :: degrees_to_radians = pi/180.0_8
20real(8), parameter :: k1_lat_center =   20.d0*degrees_to_radians
21real(8), parameter :: k1_lon_center =  300.d0*degrees_to_radians
22
23contains
24
25!===============================================================================
26!  Solar photolysis rate and recombination rate
27!===============================================================================
28
29subroutine k_vals( lat, lon, k1, k2 )
30
31!-----------------------------------------------------------------------
32! Arguments:
33!-----------------------------------------------------------------------
34real(8), intent(in)    :: lat, lon  ! latitude and longitude, radians
35real(8), intent(out)   :: k1, k2    ! reaction rates
36
37k1 = 1.0_8*max(0.d0,sin(lat)*sin(k1_lat_center) + cos(lat)*cos(k1_lat_center)*cos(lon-k1_lon_center))
38k2 = 1._8
39
40return
41
42end subroutine k_vals
43
44!===============================================================================
45!  Tendencies of cl and cl2
46!===============================================================================
47
48subroutine tendency_Terminator( lat, lon, cl, cl2, dt, cl_f, cl2_f )
49
50!-----------------------------------------------------------------------
51! Arguments:
52!-----------------------------------------------------------------------
53
54real(8), intent(in)    :: lat, lon  ! latitude and longitude, degrees
55real(8), intent(in)    :: cl, cl2   ! molar mixing ratio of cl and cl2
56real(8), intent(in)    :: dt        ! size of physics time step
57
58real(8), intent(out)   :: cl_f, cl2_f  ! time rate of change of cl and cl2
59
60!-----------------------------------------------------------------------
61! Local variables
62!-----------------------------------------------------------------------
63
64real(8) :: r, det, expdt, el ! useful algebaic quantities used in the computation
65real(8) :: k1, k2            ! reaction rates
66real(8) :: cly               ! quantity that should be conseved
67
68call k_vals( lat, lon, k1, k2 )
69
70r = k1 / (4._8*k2)
71cly = cl + 2._8* cl2
72
73det = sqrt( r*r + 2._8*r*cly )
74expdt = exp( -4._8*k2*det*dt )
75
76if ( abs(det * k2 * dt) .gt. 1e-16 ) then
77el = (1._8 - expdt) /det /dt
78else
79el = 4._8*k2
80endif
81
82cl_f  = -el * (cl - det + r)*(cl + det + r) / (1._8 + expdt + dt*el*(cl + r))
83cl2_f = -cl_f / 2._8
84
85return
86
87end subroutine tendency_Terminator
88
89!===============================================================================
90!  Compute initial values
91!===============================================================================
92
93subroutine initial_value_Terminator( lat, lon, cl, cl2 )
94
95!-----------------------------------------------------------------------
96! Arguments:
97!-----------------------------------------------------------------------
98
99real(8), intent(in)  :: lat, lon  ! latitude and longitude, degrees
100real(8), intent(out) :: cl, cl2   ! molar mixing ratio of cl and cl2
101
102!-----------------------------------------------------------------------
103! Local variables
104!-----------------------------------------------------------------------
105
106real(8) :: r, det  ! useful algebraic forms
107real(8) :: k1, k2  ! reaction rates
108
109call k_vals( lat, lon, k1, k2 )
110
111r = k1 / (4._8*k2)
112det = sqrt(r*r + 2._8*cly_constant*r)
113
114cl  = (det-r)
115cl2 = cly_constant/2._8 - (det-r)/2._8
116
117return
118
119end subroutine initial_value_Terminator
120
121end module Terminator
122
123
124
Note: See TracBrowser for help on using the repository browser.