source: trunk/LMDZ.COMMON/libf/evolution/phys_constants.F90 @ 3985

Last change on this file since 3985 was 3985, checked in by jbclement, 6 days ago

PEM:
Addition of a module "phys_constants" to read and store physical parameter of the planet properly, i.e. without going through the module "comcstfi_h" and/or "comconst_mod".
JBC

File size: 2.3 KB
Line 
1MODULE phys_constants
2
3implicit none
4
5real, parameter :: pi = 4.*atan(1.) ! PI = 3.14159...
6
7real :: g     ! Gravity (m/s2)
8real :: r     ! Reduced gas constant,r = 8.314511/(mugaz/1000.0)
9real :: mugaz ! Molar mass of the atmosphere (g/mol)
10real :: rad   ! Radius of the planet (m)
11real :: cpp   ! Cp of the atmosphere
12real :: rcp   ! r/cpp
13
14!=======================================================================
15contains
16!=======================================================================
17
18SUBROUTINE read_constants(fname)
19
20use netcdf
21
22implicit none
23
24! Arguments
25!----------
26character(*), intent(in) :: fname
27
28! Local variables
29!----------------
30real, dimension(:), allocatable :: controle
31integer :: ncid           ! File ID
32integer :: varid_controle ! Variable ID for 'controle'
33integer :: dimid_index    ! Dimension ID for 'index'
34integer :: nindex         ! Size of dimension 'index'
35integer :: ierr           ! Return codes
36
37! Code
38!-----
39! Open the NetCDF file
40ierr = nf90_open(trim(fname),NF90_NOWRITE,ncid)
41if (ierr /= nf90_noerr) then
42    write(*,*) "Error opening file:", trim(nf90_strerror(ierr))
43    error stop
44endif
45
46! Get the dimension size of 'index'
47ierr = nf90_inq_dimid(ncid,"index",dimid_index)
48if (ierr /= nf90_noerr) then
49    write(*,*) "Error getting dimid 'index':", trim(nf90_strerror(ierr))
50    error stop
51endif
52ierr = nf90_inquire_dimension(ncid,dimid_index,len = nindex)
53if (ierr /= nf90_noerr) then
54    write(*,*) "Error getting dimension length:", trim(nf90_strerror(ierr))
55    error stop
56endif
57
58! Get the variable ID for 'controle'
59allocate(controle(nindex))
60ierr = nf90_inq_varid(ncid,"controle",varid_controle)
61if (ierr /= nf90_noerr) then
62    write(*,*) "Error getting variable ID 'controle':", trim(nf90_strerror(ierr))
63    error stop
64endif
65ierr = nf90_get_var(ncid,varid_controle,controle)
66if (ierr /= nf90_noerr) then
67    write(*,*) "Error reading 'controle':", trim(nf90_strerror(ierr))
68    error stop
69endif
70
71! Close the file
72ierr = nf90_close(ncid)
73if (ierr /= nf90_noerr) then
74    write(*,*) "Error closing file:", trim(nf90_strerror(ierr))
75    error stop
76endif
77
78! Initialize the constants with 'controle' data
79rad   = controle(5)
80g     = controle(7)
81mugaz = controle(8)
82rcp   = controle(9)
83r     = 8.314511*1000./mugaz
84cpp   = r/rcp
85deallocate(controle)
86
87END SUBROUTINE read_constants
88
89END MODULE phys_constants
Note: See TracBrowser for help on using the repository browser.