/[lmdze]/trunk/Sources/dyn3d/Guide/Read_reanalyse/read_reanalyse.f
ViewVC logotype

Contents of /trunk/Sources/dyn3d/Guide/Read_reanalyse/read_reanalyse.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 172 - (show annotations)
Wed Sep 30 15:59:14 2015 UTC (8 years, 8 months ago) by guez
File size: 3861 byte(s)
Just indented correctbid and nat2gcm.

The procedure read_reanalyse just reads the next time slab every time
it is called. No use keeping track of the time index in the calling
procedure, guide. It is simpler to do it in read_reanalyse. Also
simpler to read the number of vertical levels in read_reanalyse than
in guide, since we have already in read_reanalyse the input of
pressure levels. We then have to make the arrays containing reanalyses
static allocatable instead of automatic. Also only read pressure
levels at the first call of read_reanalyse instead of at every call.

masserea2 not used in guide. Remove it down the chain in
read_reanalyse and reanalyse2nat.

1 module read_reanalyse_m
2
3 IMPLICIT NONE
4
5 contains
6
7 subroutine read_reanalyse(psi, u, v, t, q)
8
9 ! From LMDZ4/libf/dyn3d/read_reanalyse.F, version 1.3, 2005/04/15 12:31:21
10
11 USE conf_guide_m, ONLY: guide_q, guide_t, guide_u, guide_v, ncep
12 USE dimens_m, ONLY: iim, jjm, llm
13 USE netcdf, ONLY: nf90_nowrite
14 USE netcdf95, ONLY: nf95_get_var, nf95_inq_dimid, nf95_inq_varid, &
15 nf95_inquire_dimension, nf95_open
16 USE paramet_m, ONLY: iip1, jjp1
17 use reanalyse2nat_m, only: reanalyse2nat
18
19 real, intent(in):: psi(:, :) ! (iip1, jjp1)
20 real, intent(out):: u(:, :, :) ! (iip1, jjp1, llm)
21 real, intent(out):: v(:, :, :) ! (iip1, jjm, llm)
22 real, intent(out):: t(:, :, :), q(:, :, :) ! (iip1, jjp1, llm)
23
24 ! Local:
25 integer, save:: nlevnc
26 integer:: timestep = 0
27 real pk(iip1, jjp1, llm)
28 integer, save:: ncidu, varidu, ncidv, varidv, ncidt, varidt, ncidQ, varidQ
29 integer ncid, varid, dimid
30 real, allocatable, save:: unc(:, :, :) ! (iip1, jjp1, nlevnc)
31 real, allocatable, save:: vnc(:, :, :) ! (iip1, jjm, nlevnc)
32 real, allocatable, save:: tnc(:, :, :), Qnc(:, :, :) ! (iip1, jjp1, nlevnc)
33 real, allocatable, save:: pl(:) ! (nlevnc)
34 logical:: first = .true.
35 character(len = 8) name
36
37 ! -----------------------------------------------------------------
38
39 ! Initialisation de la lecture des fichiers
40
41 if (first) then
42 print *, 'Intitialisation de read reanalsye'
43
44 ! Vent zonal
45 if (guide_u) then
46 call nf95_open('u.nc', nf90_nowrite, ncidu)
47 call nf95_inq_varid(ncidu, 'UWND', varidu)
48 endif
49
50 ! Vent meridien
51 if (guide_v) then
52 call nf95_open('v.nc', nf90_nowrite, ncidv)
53 call nf95_inq_varid(ncidv, 'VWND', varidv)
54 endif
55
56 ! Temperature
57 if (guide_T) then
58 call nf95_open('T.nc', nf90_nowrite, ncidt)
59 call nf95_inq_varid(ncidt, 'AIR', varidt)
60 endif
61
62 ! Humidite
63 if (guide_Q) then
64 call nf95_open('hur.nc', nf90_nowrite, ncidQ)
65 call nf95_inq_varid(ncidQ, 'RH', varidQ)
66 endif
67
68 ! Coordonn\'ee verticale :
69
70 if (guide_u) then
71 ncid = ncidu
72 else if (guide_v) then
73 ncid = ncidv
74 else if (guide_T) then
75 ncid = ncidt
76 else
77 ncid = ncidq
78 end if
79
80 name = merge('LEVEL ', 'PRESSURE', ncep)
81 call nf95_inq_dimid(ncid, name, dimid)
82 call nf95_inquire_dimension(ncid, dimid, nclen = nlevnc)
83 call nf95_inq_varid(ncid, name, varid)
84 PRINT *, 'nlevnc = ', nlevnc
85 allocate(unc(iip1, jjp1, nlevnc), vnc(iip1, jjm, nlevnc))
86 allocate(tnc(iip1, jjp1, nlevnc), Qnc(iip1, jjp1, nlevnc), pl(nlevnc))
87 call NF95_GET_VAR(ncid, varid, pl)
88 pl = 100. * pl ! passage en pascal
89 first = .false.
90 endif
91
92 ! lecture des champs u, v, T
93
94 timestep = timestep + 1
95 unc = 0.
96 vnc = 0.
97 tnc = 0.
98 Qnc = 0.
99
100 ! Vent zonal
101 if (guide_u) then
102 call NF95_GET_VAR(ncidu, varidu, unc, start = (/1, 1, 1, timestep/))
103 call correctbid(iim, jjp1 * nlevnc, unc)
104 endif
105
106 ! Temperature
107 if (guide_T) then
108 call NF95_GET_VAR(ncidt, varidt, tnc, start = (/1, 1, 1, timestep/))
109 call correctbid(iim, jjp1 * nlevnc, tnc)
110 endif
111
112 ! Humidite
113 if (guide_Q) then
114 call NF95_GET_VAR(ncidQ, varidQ, Qnc, start = (/1, 1, 1, timestep/))
115 call correctbid(iim, jjp1 * nlevnc, Qnc)
116 endif
117
118 ! Vent meridien
119 if (guide_v) then
120 call NF95_GET_VAR(ncidv, varidv, vnc, start = (/1, 1, 1, timestep/))
121 call correctbid(iim, jjm * nlevnc, vnc)
122 endif
123
124 call reanalyse2nat(nlevnc, psi, unc, vnc, tnc, Qnc, pl, u, v, t, Q, pk)
125 call nat2gcm(u, v, t, Q, pk, u, v, t, Q)
126
127 end subroutine read_reanalyse
128
129 end module read_reanalyse_m

  ViewVC Help
Powered by ViewVC 1.1.21