[1154] | 1 | ! $Id$ |
---|
| 2 | module open_climoz_m |
---|
| 3 | |
---|
| 4 | implicit none |
---|
| 5 | |
---|
| 6 | contains |
---|
| 7 | |
---|
| 8 | subroutine open_climoz(ncid, press_in_edg) |
---|
| 9 | |
---|
| 10 | ! This procedure should be called once per "gcm" run, by a single |
---|
| 11 | ! thread of each MPI process. |
---|
| 12 | ! The root process opens "climoz_LMDZ.nc", reads the pressure |
---|
| 13 | ! levels and broadcasts them to the other processes. |
---|
| 14 | |
---|
| 15 | ! We assume that, in "climoz_LMDZ.nc", the pressure levels are in hPa |
---|
| 16 | ! and in strictly ascending order. |
---|
| 17 | |
---|
| 18 | use netcdf95, only: nf95_open, nf95_close, nf95_gw_var, nf95_inq_varid |
---|
| 19 | use netcdf, only: nf90_nowrite |
---|
| 20 | |
---|
| 21 | use mod_phys_lmdz_mpi_data, only: is_mpi_root |
---|
| 22 | use mod_phys_lmdz_mpi_transfert, only: bcast_mpi ! broadcast |
---|
| 23 | |
---|
| 24 | integer, intent(out):: ncid ! of "climoz_LMDZ.nc" |
---|
| 25 | |
---|
| 26 | real, pointer:: press_in_edg(:) |
---|
| 27 | ! (edges of pressure intervals for ozone climatology, in Pa, in strictly |
---|
| 28 | ! ascending order) |
---|
| 29 | |
---|
| 30 | ! Variables local to the procedure: |
---|
| 31 | |
---|
| 32 | real, pointer:: plev(:) |
---|
| 33 | ! (pressure levels for ozone climatology, converted to Pa, in strictly |
---|
| 34 | ! ascending order) |
---|
| 35 | |
---|
| 36 | integer varid ! for NetCDF |
---|
| 37 | integer n_plev ! number of pressure levels in the input data |
---|
| 38 | integer k |
---|
| 39 | |
---|
| 40 | !--------------------------------------- |
---|
| 41 | |
---|
| 42 | print *, "Call sequence information: open_climoz" |
---|
| 43 | |
---|
| 44 | if (is_mpi_root) then |
---|
| 45 | call nf95_open("climoz_LMDZ.nc", nf90_nowrite, ncid) |
---|
| 46 | |
---|
| 47 | call nf95_inq_varid(ncid, "plev", varid) |
---|
| 48 | call nf95_gw_var(ncid, varid, plev) |
---|
| 49 | ! Convert from hPa to Pa because "paprs" and "pplay" are in Pa: |
---|
| 50 | plev = plev * 100. |
---|
| 51 | n_plev = size(plev) |
---|
| 52 | end if |
---|
| 53 | |
---|
| 54 | call bcast_mpi(n_plev) |
---|
| 55 | if (.not. is_mpi_root) allocate(plev(n_plev)) |
---|
| 56 | call bcast_mpi(plev) |
---|
| 57 | |
---|
| 58 | ! Compute edges of pressure intervals: |
---|
| 59 | allocate(press_in_edg(n_plev + 1)) |
---|
| 60 | if (is_mpi_root) then |
---|
| 61 | press_in_edg(1) = 0. |
---|
| 62 | ! We choose edges halfway in logarithm: |
---|
| 63 | forall (k = 2:n_plev) press_in_edg(k) = sqrt(plev(k - 1) * plev(k)) |
---|
| 64 | press_in_edg(n_plev + 1) = huge(0.) |
---|
| 65 | ! (infinity, but any value guaranteed to be greater than the |
---|
| 66 | ! surface pressure would do) |
---|
| 67 | end if |
---|
| 68 | call bcast_mpi(press_in_edg) |
---|
| 69 | deallocate(plev) ! pointer |
---|
| 70 | |
---|
| 71 | end subroutine open_climoz |
---|
| 72 | |
---|
| 73 | end module open_climoz_m |
---|