/[lmdze]/trunk/Sources/dyn3d/comconst.f
ViewVC logotype

Diff of /trunk/Sources/dyn3d/comconst.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

trunk/libf/dyn3d/comconst.f90 revision 3 by guez, Wed Feb 27 13:16:39 2008 UTC trunk/dyn3d/comconst.f revision 82 by guez, Wed Mar 5 14:57:53 2014 UTC
# Line 1  Line 1 
1  module comconst  module comconst
2    
3    use dimens_m, only: jjm    use nr_util, only: pi
4    
5    implicit none    implicit none
6    
   INTEGER im, jm, lllm, imp1  
   integer, parameter:: jmp1 = jjm + 1  
   integer lllmm1, lllmp1, lcl  
   REAL dtvr ! time step for dynamics (in s)  
7    real, parameter:: daysec = 86400. ! number of seconds per day    real, parameter:: daysec = 86400. ! number of seconds per day
8    REAL pi, dtphys, dtdiss  
9      REAL dtvr ! time step for dynamics, in s
10      REAL dtphys ! time step for physics, in s
11    
12    real, parameter:: rad = 6371229. ! radius of the Earth (in m)    real, parameter:: rad = 6371229. ! radius of the Earth (in m)
13    real r  
14    real, parameter:: cpp = 1004.70885, kappa = 0.2857143    real, parameter:: cpp = 1004.70885
15    REAL cotot, unsim    ! specific heat capacity at constant pressure of dry air, in J K-1 kg-1
16    
17      real, parameter:: kappa = 0.2857143 ! r / cpp
18    
19      real, parameter:: r = cpp * kappa
20      ! specific ideal gas constant for dry air, in J K-1 kg-1
21    
22    real, parameter:: g = 9.8 ! acceleration of gravity (in m s-2)    real, parameter:: g = 9.8 ! acceleration of gravity (in m s-2)
   real omeg ! angular speed of rotation of the Earth (in rad s-1)  
23    
24    private jjm    real, parameter:: omeg = 2 * pi / daysec
25      ! angular speed of rotation of the Earth (in rad s-1)
26    
27      private pi
28    
29  contains  contains
30    
31    subroutine initialize    SUBROUTINE iniconst
32    
33        ! From dyn3d/iniconst.F,v 1.1.1.1 2004/05/19 12:53:05
34        ! P. Le Van
35    
36        USE conf_gcm_m, ONLY: day_step, iphysiq
37    
38        IMPLICIT NONE
39    
40        !-----------------------------------------------------------------------
41    
42        dtvr = daysec / real(day_step)
43        dtphys  = iphysiq * dtvr
44    
45      print *, "Call sequence information: initialize"      print *, 'dtvr = ', dtvr
46      pi     =  Acos(-1.)      print *, 'dtphys = ', dtphys
47      omeg   = 2 * pi / daysec      PRINT *, 'cpp = ', cpp
48        PRINT *, 'R = ', r
49        PRINT *, 'kappa = ', kappa
50    
51    end subroutine initialize    END SUBROUTINE iniconst
52    
53  end module comconst  end module comconst

Legend:
Removed from v.3  
changed lines
  Added in v.82

  ViewVC Help
Powered by ViewVC 1.1.21