source: LMDZ6/branches/Amaury_dev/libf/phylmd/press_coefoz_m.F90 @ 5409

Last change on this file since 5409 was 5160, checked in by abarral, 5 months ago

Put .h into modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 2.1 KB
Line 
1! $Id$
2module press_coefoz_m
3
4  IMPLICIT NONE
5
6  REAL, ALLOCATABLE, save:: plev(:)
7  ! (pressure level of Mobidic input data, converted to Pa, in strictly
8  ! ascending order)
9
10  REAL, ALLOCATABLE, save:: press_in_edg(:)
11  ! (edges of pressure intervals for Mobidic input data, in Pa, in strictly
12  ! ascending order)
13
14CONTAINS
15
16  SUBROUTINE press_coefoz
17
18    ! This procedure is called once per "gcm" run.
19    ! A single thread of the root process reads the pressure levels
20    ! from "coefoz_LMDZ.nc" and broadcasts them to the other processes.
21
22    ! We assume that, in "coefoz_LMDZ.nc", the pressure levels are in hPa
23    ! and in strictly ascending order.
24
25    USE netcdf95, ONLY: nf95_open, nf95_close, nf95_gw_var, nf95_inq_varid
26    USE netcdf, ONLY: nf90_nowrite
27
28    USE lmdz_phys_mpi_data, ONLY: is_mpi_root
29    USE lmdz_phys_mpi_transfert, ONLY: bcast_mpi ! broadcast
30
31    ! Variables local to the procedure:
32    INTEGER ncid, varid ! for NetCDF
33    INTEGER n_plev ! number of pressure levels in the input data
34    INTEGER k
35
36    !---------------------------------------
37
38    !$omp single
39    PRINT *, "Call sequence information: press_coefoz"
40
41    IF (is_mpi_root) THEN
42       CALL nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
43
44       CALL nf95_inq_varid(ncid, "plev", varid)
45       CALL nf95_gw_var(ncid, varid, plev)
46       ! Convert from hPa to Pa because "paprs" and "pplay" are in Pa:
47       plev = plev * 100.
48       n_plev = size(plev)
49
50       CALL nf95_close(ncid)
51    end if
52
53    CALL bcast_mpi(n_plev)
54    IF (.NOT. is_mpi_root) allocate(plev(n_plev))
55    CALL bcast_mpi(plev)
56   
57    ! Compute edges of pressure intervals:
58    allocate(press_in_edg(n_plev + 1))
59    IF (is_mpi_root) THEN
60       press_in_edg(1) = 0.
61       ! We choose edges halfway in logarithm:
62       DO k = 2,n_plev
63          press_in_edg(k) = SQRT(plev(k - 1) * plev(k))
64       ENDDO
65       press_in_edg(n_plev + 1) = huge(0.)
66       ! (infinity, but any value guaranteed to be greater than the
67       ! surface pressure would do)
68    end if
69    CALL bcast_mpi(press_in_edg)
70    !$omp end single
71
72  END SUBROUTINE  press_coefoz
73
74END MODULE press_coefoz_m
Note: See TracBrowser for help on using the repository browser.