source: LMDZ6/trunk/libf/dynphy_lonlat/phylmd/test_disvert_m.f90 @ 5297

Last change on this file since 5297 was 5271, checked in by abarral, 7 days ago

Move dimensions.h into a module
Nb: doesn't compile yet

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