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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 178 - (hide annotations)
Fri Mar 11 18:47:26 2016 UTC (8 years, 2 months ago) by guez
File size: 3903 byte(s)
Moved variables date0, deltat, datasz_max, ncvar_ids, point, buff_pos,
buffer, regular from module histcom_var to modules where they are
defined.

Removed procedure ioipslmpp, useless for a sequential program.

Added argument datasz_max to histwrite_real (to avoid circular
dependency with histwrite).

Removed useless variables and computations everywhere.

Changed real litteral constants from default kind to double precision
in lwb, lwu, lwvn, sw1s, swtt, swtt1, swu.

Removed unused arguments: paer of sw, sw1s, sw2s, swclr; pcldsw of
sw1s, sw2s; pdsig, prayl of swr; co2_ppm of clmain, clqh; tsol of
transp_lay; nsrf of screenp; kcrit and kknu of gwstress; pstd of
orosetup.

Added output of relative humidity.

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

  ViewVC Help
Powered by ViewVC 1.1.21