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

Contents of /trunk/dyn3d/test_disvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 99 - (show annotations)
Wed Jul 2 18:39:15 2014 UTC (9 years, 10 months ago) by guez
File size: 3517 byte(s)
Created procedure test_disvert (following LMDZ). Added procedures
hybrid and funcd in module disvert_m. Upgraded compute_ab from
internal procedure of disvert to module procedure. Added variables y,
ya in module disvert_m. Upgraded s from local variable of procedure
disvert to module variable.

Renamed allowed value of variable vert_sampling in procedure disvert
from "read" to "read_hybrid". Added possibility to read pressure
values, value "read_pressure". Replaced vertical distribution for
value "param" by the distribution "strato_correct" from LMDZ (but kept
the value "param"). In case "tropo", replaced 1 by dsigmin (following
LMDZ). In case "strato", replaced 0.3 by dsigmin (following LMDZ).

Changed computation of bp in procedure compute_ab.

Removed debugindex case in clmain. Removed useless argument rlon of
procedure clmain. Removed useless variables ytaux, ytauy of procedure
clmain.

Removed intermediary variables tsol, qsol, tsolsrf, tslab in procedure
etat0.

Removed variable ok_veget:. coupling with the model Orchid is not
possible. Removed variable ocean: modeling an ocean slab is not
possible.

Removed useless variables tmp_rriv and tmp_rcoa from module
interface_surf.

Moved initialization of variables da, mp, phi in procedure physiq to
to inside the test iflag_con >= 3.

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

  ViewVC Help
Powered by ViewVC 1.1.21