source: LMDZ5/branches/LMDZ5_SPLA/libf/phylmd/test_disvert_m.F90 @ 2265

Last change on this file since 2265 was 2022, checked in by lguez, 11 years ago

Forgot copyright property on new file.

  • 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 
1module test_disvert_m
2
3  implicit none
4
5contains
6
7  subroutine test_disvert
8
9    ! Author: Lionel GUEZ
10
11    ! This procedure tests the order of pressure values at half-levels
12    ! and full levels. We arbitrarily choose to test ngrid values of
13    ! the surface pressure, which sample possible values on Earth.
14
15    use exner_hyb_m, only: exner_hyb
16
17    ! For llm:
18    include "dimensions.h"
19
20    ! For ap, bp, preff:
21    include "comvert.h"
22
23    ! For kappa, cpp:
24    include "comconst.h"
25
26    ! Local:
27    integer l, i
28    integer, parameter:: ngrid = 7
29    real p(ngrid, llm + 1) ! pressure at half-level, in Pa
30    real pks(ngrid) ! exner function at the surface, in J K-1 kg-1
31    real pk(ngrid, llm) ! exner function at full level, in J K-1 kg-1
32    real ps(ngrid) ! surface pressure, in Pa
33    real p_lay(ngrid, llm) ! pressure at full level, in Pa
34    real delta_ps ! in Pa
35
36    !---------------------
37
38    print *, "Call sequence information: test_disvert"
39
40    delta_ps = 6e4 / (ngrid - 1)
41    ps = (/(5e4 + delta_ps * i, i = 0, ngrid - 1)/)
42    forall (l = 1: llm + 1) p(:, l) = ap(l) + bp(l) * ps
43    call exner_hyb(ngrid, ps, p, pks, pk)
44    p_lay = preff * (pk / cpp)**(1. / kappa)
45
46    ! Are pressure values in the right order?
47    if (any(p(:, :llm) <= p_lay .or. p_lay <= p(:, 2:))) then
48       ! List details and stop:
49       do l = 1, llm
50          do i = 1, ngrid
51             if (p(i, l) <= p_lay(i, l)) then
52                print 1000, "p(", i, ", ", l, ") = ", p(i, l) / 100., &
53                     " hPa <= p_lay(", i, ", ", l, ") = ", &
54                     p_lay(i, l) / 100., " hPa"
55             end if
56             if (p_lay(i, l) <= p(i, l + 1)) then
57                print 1000, "p_lay(", i, ", ", l, ") = ", p_lay(i, l) / 100., &
58                     " hPa <= p(", i, ", ", l + 1, ") = ", &
59                     p(i, l + 1) / 100., " hPa"
60             end if
61          end do
62       end do
63       call abort_gcm("test_disvert", "bad order of pressure values", 1)
64    end if
65
661000 format (2(a, 2(i0, a), g13.6), a)
67
68  end subroutine test_disvert
69
70end module test_disvert_m
Note: See TracBrowser for help on using the repository browser.