/[lmdze]/trunk/Sources/dyn3d/test_disvert.f
ViewVC logotype

Annotation of /trunk/Sources/dyn3d/test_disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 107 - (hide annotations)
Thu Sep 11 15:09:15 2014 UTC (9 years, 8 months ago) by guez
Original Path: trunk/dyn3d/test_disvert.f
File size: 3474 byte(s)
Imported procedure grilles_gcm_sub from LMDZ. Had then to transform
local variable phis of etat to argument.

Replaced calls to lnblnk by calls to trim.

Removed arguments nlat, klevel and griscal of filtreg. Replaced
integer arguments ifiltre and iaire by logical arguments direct and
intensive.

Changed default values of guide_t and guide_q to false.

1 guez 99 module test_disvert_m
2    
3     implicit none
4    
5     contains
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, and writes the vertical coordinates to files
13     ! "half_level.csv" and "full_level.csv" for a reference surface
14     ! pressure. We arbitrarily choose to test ngrid values of the
15     ! surface pressure, which sample possible values on Earth.
16    
17     use abort_gcm_m, only: abort_gcm
18     use comconst, only: kappa, cpp
19 guez 107 use disvert_m, only: ap, bp, preff, presnivs, s
20 guez 99 use dimens_m, only: llm
21     use exner_hyb_m, only: exner_hyb
22     use jumble, only: new_unit
23    
24     integer i, unit, l
25    
26     integer, parameter:: ngrid = 8
27    
28     real p(ngrid, 1, llm + 1) ! pressure at half-level, in Pa
29     real z(llm) ! pressure-altitude at half-level (km)
30     real pks(ngrid, 1) ! exner function at the surface, in J K-1 kg-1
31     real pk(ngrid, 1, llm) ! exner function at full level, in J K-1 kg-1
32     real ps(ngrid, 1) ! surface pressure, in Pa
33     real p_lay(ngrid, llm) ! pressure at full level, in Pa
34     real, parameter:: delta_ps = 6e4 / (ngrid - 2) ! in Pa
35    
36     !---------------------
37    
38     print *, "Call sequence information: test_disvert"
39    
40     ps(:, 1) = (/(5e4 + delta_ps * i, i = 0, ngrid - 2), preff/)
41     forall (l = 1: llm + 1) p(:, 1, l) = ap(l) + bp(l) * ps(:, 1)
42     call exner_hyb(ps, p, pks, pk)
43     p_lay = preff * (pk(:, 1, :) / cpp)**(1. / kappa)
44    
45     ! Are pressure values in the right order?
46     if (any(p(:, 1, :llm) <= p_lay .or. p_lay <= p(:, 1, 2:))) then
47     ! List details and stop:
48     do l = 1, llm
49     do i = 1, ngrid
50     if (p(i, 1, l) <= p_lay(i, l)) then
51     print 1000, "ps = ", ps(i, 1) / 100., "hPa, p(level ", l, &
52     ") = ", p(i, 1, l) / 100., " hPa <= p(layer ", l, ") = ", &
53     p_lay(i, l) / 100., " hPa"
54     end if
55     if (p_lay(i, l) <= p(i, 1, l+1)) then
56     print 1000, "ps = ", ps(i, 1) / 100., &
57     "hPa, p(layer ", l, ") = ", p_lay(i, l) / 100., &
58     " hPa <= p(level ", l + 1, ") = ", &
59     p(i, 1, 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    
66     ! Write distribution for the reference surface pressure (index ngrid):
67    
68     z = 7. * log(preff / p(ngrid, 1, :llm))
69     call new_unit(unit)
70    
71     open(unit, file="half_level.csv", status="replace", action="write")
72     ! Title line:
73     write(unit, fmt=*) '"ap (hPa)" "bp" "s" "pressure (hPa)" "z (km)"'
74     do l = 1, llm
75     write(unit, fmt=*) ap(l) / 100., bp(l), s(l), p(ngrid, 1, l) / 100., z(l)
76     end do
77     close(unit)
78     print *, 'The file "half_level.csv" has been created.'
79    
80     open(unit, file="full_level.csv", status="replace", action="write")
81     ! Title line:
82     write(unit, fmt=*) &
83     '"pressure (hPa)" "z (km)" "presnivs (hPa)" "delta z (km)"'
84     do l = 1, llm - 1
85     write(unit, fmt=*) p_lay(ngrid, l) / 100., &
86     7. * log(preff / p_lay(ngrid, l)), presnivs(l) / 100., z(l+1) - z(l)
87     end do
88     write(unit, fmt=*) p_lay(ngrid, llm) / 100., &
89     7. * log(preff / p_lay(ngrid, llm)), presnivs(llm) / 100.
90     close(unit)
91     print *, 'The file "full_level.csv" has been created.'
92    
93     1000 format (3(a, g10.4, a, i0))
94    
95     end subroutine test_disvert
96    
97     end module test_disvert_m

  ViewVC Help
Powered by ViewVC 1.1.21