source: LMDZ5/branches/LMDZ5_SPLA/libf/phylmd/read_pstoke.F90 @ 4934

Last change on this file since 4934 was 1992, checked in by lguez, 11 years ago

Converted to free source form files in libf/phylmd which were still in
fixed source form. The conversion was done using the polish mode of
the NAG Fortran Compiler.

In addition to converting to free source form, the processing of the
files also:

-- indented the code (including comments);

-- set Fortran keywords to uppercase, and set all other identifiers
to lower case;

-- added qualifiers to end statements (for example "end subroutine
conflx", instead of "end");

-- changed the terminating statements of all DO loops so that each
loop ends with an ENDDO statement (instead of a labeled continue).

-- replaced #include by include.

  • 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
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 14.3 KB
Line 
1
2! $Id: read_pstoke.F90 1992 2014-03-05 13:19:12Z abarral $
3
4
5
6SUBROUTINE read_pstoke(irec, zrec, zklono, zklevo, airefi, phisfi, t, mfu, &
7    mfd, en_u, de_u, en_d, de_d, coefh, fm_therm, en_therm, frac_impa, &
8    frac_nucl, pyu1, pyv1, ftsol, psrf)
9
10  ! ******************************************************************************
11  ! Frederic HOURDIN, Abderrahmane IDELKADI
12  ! Lecture des parametres physique stockes online necessaires pour
13  ! recalculer offline le transport de traceurs sur une grille 2x plus fine
14  ! que
15  ! celle online
16  ! A FAIRE : une seule routine au lieu de 2 (lectflux, redecoupe)!
17  ! ******************************************************************************
18
19  USE netcdf
20  USE dimphy
21  USE control_mod
22  USE indice_sol_mod
23
24  IMPLICIT NONE
25
26  include "netcdf.inc"
27  include "dimensions.h"
28  include "paramet.h"
29  include "comconst.h"
30  include "comgeom.h"
31  include "temps.h"
32  include "ener.h"
33  include "logic.h"
34  include "description.h"
35  include "serre.h"
36  ! ccc#include "dimphy.h"
37
38  INTEGER klono, klevo, imo, jmo
39  PARAMETER (imo=iim/2, jmo=(jjm+1)/2)
40  PARAMETER (klono=(jmo-1)*imo+2, klevo=llm)
41  REAL phisfi(klono)
42  REAL phisfi2(imo, jmo+1), airefi2(imo, jmo+1)
43
44  REAL mfu(klono, klevo), mfd(klono, klevo)
45  REAL en_u(klono, klevo), de_u(klono, klevo)
46  REAL en_d(klono, klevo), de_d(klono, klevo)
47  REAL coefh(klono, klevo)
48  REAL fm_therm(klono, klevo), en_therm(klono, klevo)
49
50  REAL mfu2(imo, jmo+1, klevo), mfd2(imo, jmo+1, klevo)
51  REAL en_u2(imo, jmo+1, klevo), de_u2(imo, jmo+1, klevo)
52  REAL en_d2(imo, jmo+1, klevo), de_d2(imo, jmo+1, klevo)
53  REAL coefh2(imo, jmo+1, klevo)
54  REAL fm_therm2(imo, jmo+1, klevo)
55  REAL en_therm2(imo, jmo+1, klevo)
56
57  REAL pl(klevo)
58  INTEGER irec
59  INTEGER xid, yid, zid, tid
60  REAL zrec, zklono, zklevo, zim, zjm
61  INTEGER ncrec, ncklono, ncklevo, ncim, ncjm
62
63  REAL airefi(klono)
64  CHARACTER *20 namedim
65
66  ! !! attention !!
67  ! attention il y a aussi le pb de def klono
68  ! dim de phis??
69
70
71  REAL frac_impa(klono, klevo), frac_nucl(klono, klevo)
72  REAL frac_impa2(imo, jmo+1, klevo), frac_nucl2(imo, jmo+1, klevo)
73  REAL pyu1(klono), pyv1(klono)
74  REAL pyu12(imo, jmo+1), pyv12(imo, jmo+1)
75  REAL ftsol(klono, nbsrf)
76  REAL psrf(klono, nbsrf)
77  REAL ftsol1(klono), ftsol2(klono), ftsol3(klono), ftsol4(klono)
78  REAL psrf1(klono), psrf2(klono), psrf3(klono), psrf4(klono)
79  REAL ftsol12(imo, jmo+1), ftsol22(imo, jmo+1), ftsol32(imo, jmo+1), &
80    ftsol42(imo, jmo+1)
81  REAL psrf12(imo, jmo+1), psrf22(imo, jmo+1), psrf32(imo, jmo+1), &
82    psrf42(imo, jmo+1)
83  REAL t(klono, klevo)
84  REAL t2(imo, jmo+1, klevo)
85  INTEGER ncidp
86  SAVE ncidp
87  INTEGER varidt
88  INTEGER varidmfu, varidmfd, varidps, varidenu, variddeu
89  INTEGER varidend, varidded, varidch, varidfi, varidfn
90  INTEGER varidfmth, varidenth
91  INTEGER varidyu1, varidyv1, varidpl, varidai, varididvt
92  INTEGER varidfts1, varidfts2, varidfts3, varidfts4
93  INTEGER varidpsr1, varidpsr2, varidpsr3, varidpsr4
94  SAVE varidmfu, varidmfd, varidps, varidenu, variddeu
95  SAVE varidend, varidded, varidch, varidfi, varidfn
96  SAVE varidfmth, varidenth
97  SAVE varidyu1, varidyv1, varidpl, varidai, varididvt
98  SAVE varidfts1, varidfts2, varidfts3, varidfts4
99  SAVE varidpsr1, varidpsr2, varidpsr3, varidpsr4
100  SAVE varidt
101
102  INTEGER l, i
103  INTEGER start(4), count(4), status
104  REAL rcode
105  LOGICAL first
106  SAVE first
107  DATA first/.TRUE./
108
109
110
111  ! ---------------------------------------------
112  ! Initialisation de la lecture des fichiers
113  ! ---------------------------------------------
114
115  IF (irec==0) THEN
116
117    rcode = nf90_open('phystoke.nc', nf90_nowrite, ncidp)
118
119    rcode = nf90_inq_varid(ncidp, 'phis', varidps)
120    PRINT *, 'ncidp,varidps', ncidp, varidps
121
122    rcode = nf90_inq_varid(ncidp, 'sig_s', varidpl)
123    PRINT *, 'ncidp,varidpl', ncidp, varidpl
124
125    rcode = nf90_inq_varid(ncidp, 'aire', varidai)
126    PRINT *, 'ncidp,varidai', ncidp, varidai
127
128    ! A FAIRE: Es-il necessaire de stocke t?
129    rcode = nf90_inq_varid(ncidp, 't', varidt)
130    PRINT *, 'ncidp,varidt', ncidp, varidt
131
132    rcode = nf90_inq_varid(ncidp, 'mfu', varidmfu)
133    PRINT *, 'ncidp,varidmfu', ncidp, varidmfu
134
135    rcode = nf90_inq_varid(ncidp, 'mfd', varidmfd)
136    PRINT *, 'ncidp,varidmfd', ncidp, varidmfd
137
138    rcode = nf90_inq_varid(ncidp, 'en_u', varidenu)
139    PRINT *, 'ncidp,varidenu', ncidp, varidenu
140
141    rcode = nf90_inq_varid(ncidp, 'de_u', variddeu)
142    PRINT *, 'ncidp,variddeu', ncidp, variddeu
143
144    rcode = nf90_inq_varid(ncidp, 'en_d', varidend)
145    PRINT *, 'ncidp,varidend', ncidp, varidend
146
147    rcode = nf90_inq_varid(ncidp, 'de_d', varidded)
148    PRINT *, 'ncidp,varidded', ncidp, varidded
149
150    rcode = nf90_inq_varid(ncidp, 'coefh', varidch)
151    PRINT *, 'ncidp,varidch', ncidp, varidch
152
153    ! abder (pour thermiques)
154    rcode = nf90_inq_varid(ncidp, 'fm_th', varidfmth)
155    PRINT *, 'ncidp,varidfmth', ncidp, varidfmth
156
157    rcode = nf90_inq_varid(ncidp, 'en_th', varidenth)
158    PRINT *, 'ncidp,varidenth', ncidp, varidenth
159
160    rcode = nf90_inq_varid(ncidp, 'frac_impa', varidfi)
161    PRINT *, 'ncidp,varidfi', ncidp, varidfi
162
163    rcode = nf90_inq_varid(ncidp, 'frac_nucl', varidfn)
164    PRINT *, 'ncidp,varidfn', ncidp, varidfn
165
166    rcode = nf90_inq_varid(ncidp, 'pyu1', varidyu1)
167    PRINT *, 'ncidp,varidyu1', ncidp, varidyu1
168
169    rcode = nf90_inq_varid(ncidp, 'pyv1', varidyv1)
170    PRINT *, 'ncidp,varidyv1', ncidp, varidyv1
171
172    rcode = nf90_inq_varid(ncidp, 'ftsol1', varidfts1)
173    PRINT *, 'ncidp,varidfts1', ncidp, varidfts1
174
175    rcode = nf90_inq_varid(ncidp, 'ftsol2', varidfts2)
176    PRINT *, 'ncidp,varidfts2', ncidp, varidfts2
177
178    rcode = nf90_inq_varid(ncidp, 'ftsol3', varidfts3)
179    PRINT *, 'ncidp,varidfts3', ncidp, varidfts3
180
181    rcode = nf90_inq_varid(ncidp, 'ftsol4', varidfts4)
182    PRINT *, 'ncidp,varidfts4', ncidp, varidfts4
183
184    rcode = nf90_inq_varid(ncidp, 'psrf1', varidpsr1)
185    PRINT *, 'ncidp,varidpsr1', ncidp, varidpsr1
186
187    rcode = nf90_inq_varid(ncidp, 'psrf2', varidpsr2)
188    PRINT *, 'ncidp,varidpsr2', ncidp, varidpsr2
189
190    rcode = nf90_inq_varid(ncidp, 'psrf3', varidpsr3)
191    PRINT *, 'ncidp,varidpsr3', ncidp, varidpsr3
192
193    rcode = nf90_inq_varid(ncidp, 'psrf4', varidpsr4)
194    PRINT *, 'ncidp,varidpsr4', ncidp, varidpsr4
195
196    ! ID pour les dimensions
197
198    status = nf_inq_dimid(ncidp, 'y', yid)
199    status = nf_inq_dimid(ncidp, 'x', xid)
200    status = nf_inq_dimid(ncidp, 'sig_s', zid)
201    status = nf_inq_dimid(ncidp, 'time_counter', tid)
202
203    ! lecture des dimensions
204
205    status = nf_inq_dim(ncidp, yid, namedim, ncjm)
206    status = nf_inq_dim(ncidp, xid, namedim, ncim)
207    status = nf_inq_dim(ncidp, zid, namedim, ncklevo)
208    status = nf_inq_dim(ncidp, tid, namedim, ncrec)
209
210    zrec = ncrec
211    zklevo = ncklevo
212    zim = ncim
213    zjm = ncjm
214
215    zklono = zim*(zjm-2) + 2
216
217    WRITE (*, *) 'read_pstoke : zrec = ', zrec
218    WRITE (*, *) 'read_pstoke : zklevo = ', zklevo
219    WRITE (*, *) 'read_pstoke : zim = ', zim
220    WRITE (*, *) 'read_pstoke : zjm = ', zjm
221    WRITE (*, *) 'read_pstoke : zklono = ', zklono
222
223    ! niveaux de pression
224#ifdef NC_DOUBLE
225    status = nf_get_vara_double(ncidp, varidpl, 1, zklevo, pl)
226#else
227    status = nf_get_vara_real(ncidp, varidpl, 1, zklevo, pl)
228#endif
229
230    ! lecture de aire et phis
231
232    start(1) = 1
233    start(2) = 1
234    start(3) = 1
235    start(4) = 0
236
237    count(1) = zim
238    count(2) = zjm
239    count(3) = 1
240    count(4) = 0
241
242    ! phis
243#ifdef NC_DOUBLE
244    status = nf_get_vara_double(ncidp, varidps, start, count, phisfi2)
245#else
246    status = nf_get_vara_real(ncidp, varidps, start, count, phisfi2)
247#endif
248    CALL gr_ecrit_fi(1, klono, imo, jmo+1, phisfi2, phisfi)
249
250    ! aire
251#ifdef NC_DOUBLE
252    status = nf_get_vara_double(ncidp, varidai, start, count, airefi2)
253#else
254    status = nf_get_vara_real(ncidp, varidai, start, count, airefi2)
255#endif
256    CALL gr_ecrit_fi(1, klono, imo, jmo+1, airefi2, airefi)
257  ELSE
258
259    PRINT *, 'ok1'
260
261    ! ---------------------
262    ! lecture des champs
263    ! ---------------------
264
265    PRINT *, 'WARNING!!! Il n y a pas de test de coherence'
266    PRINT *, 'sur le nombre de niveaux verticaux dans le fichier nc'
267
268    start(1) = 1
269    start(2) = 1
270    start(3) = 1
271    start(4) = irec
272
273    count(1) = zim
274    count(2) = zjm
275    count(3) = zklevo
276    count(4) = 1
277
278
279    ! *** Lessivage******************************************************
280    ! frac_impa
281#ifdef NC_DOUBLE
282    status = nf_get_vara_double(ncidp, varidfi, start, count, frac_impa2)
283#else
284    status = nf_get_vara_real(ncidp, varidfi, start, count, frac_impa2)
285#endif
286    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, frac_impa2, frac_impa)
287
288    ! frac_nucl
289#ifdef NC_DOUBLE
290    status = nf_get_vara_double(ncidp, varidfn, start, count, frac_nucl2)
291#else
292    status = nf_get_vara_real(ncidp, varidfn, start, count, frac_nucl2)
293#endif
294    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, frac_nucl2, frac_nucl)
295
296    ! *** Temperature ******************************************************
297    ! abder t
298#ifdef NC_DOUBLE
299    status = nf_get_vara_double(ncidp, varidt, start, count, t2)
300#else
301    status = nf_get_vara_real(ncidp, varidt, start, count, t2)
302#endif
303    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, t2, t)
304
305    ! *** Flux pour le calcul de la convection TIEDTK ***********************
306    ! mfu
307#ifdef NC_DOUBLE
308    status = nf_get_vara_double(ncidp, varidmfu, start, count, mfu2)
309#else
310    status = nf_get_vara_real(ncidp, varidmfu, start, count, mfu2)
311#endif
312    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, mfu2, mfu)
313
314    ! mfd
315#ifdef NC_DOUBLE
316    status = nf_get_vara_double(ncidp, varidmfd, start, count, mfd2)
317#else
318    status = nf_get_vara_real(ncidp, varidmfd, start, count, mfd2)
319#endif
320    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, mfd2, mfd)
321
322    ! en_u
323#ifdef NC_DOUBLE
324    status = nf_get_vara_double(ncidp, varidenu, start, count, en_u2)
325#else
326    status = nf_get_vara_real(ncidp, varidenu, start, count, en_u2)
327#endif
328    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_u2, en_u)
329
330    ! de_u
331#ifdef NC_DOUBLE
332    status = nf_get_vara_double(ncidp, variddeu, start, count, de_u2)
333#else
334    status = nf_get_vara_real(ncidp, variddeu, start, count, de_u2)
335#endif
336    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, de_u2, de_u)
337
338    ! en_d
339#ifdef NC_DOUBLE
340    status = nf_get_vara_double(ncidp, varidend, start, count, en_d2)
341#else
342    status = nf_get_vara_real(ncidp, varidend, start, count, en_d2)
343#endif
344    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_d2, en_d)
345
346    ! de_d
347#ifdef NC_DOUBLE
348    status = nf_get_vara_double(ncidp, varidded, start, count, de_d2)
349#else
350    status = nf_get_vara_real(ncidp, varidded, start, count, de_d2)
351#endif
352    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, de_d2, de_d)
353
354    ! **** Coeffecient du mellange
355    ! turbulent**********************************
356    ! coefh
357#ifdef NC_DOUBLE
358    status = nf_get_vara_double(ncidp, varidch, start, count, coefh2)
359#else
360    status = nf_get_vara_real(ncidp, varidch, start, count, coefh2)
361#endif
362    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, coefh2, coefh)
363
364    ! *** Flux ascendant et entrant pour les
365    ! Thermiques************************
366    ! abder thermiques
367#ifdef NC_DOUBLE
368    status = nf_get_vara_double(ncidp, varidfmth, start, count, fm_therm2)
369#else
370    status = nf_get_vara_real(ncidp, varidfmth, start, count, fm_therm2)
371#endif
372    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, fm_therm2, fm_therm)
373
374#ifdef NC_DOUBLE
375    status = nf_get_vara_double(ncidp, varidenth, start, count, en_therm2)
376#else
377    status = nf_get_vara_real(ncidp, varidenth, start, count, en_therm2)
378#endif
379    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_therm2, en_therm)
380
381    ! *** Vitesses aux sol
382    ! ******************************************************
383    start(3) = irec
384    start(4) = 0
385    count(3) = 1
386    count(4) = 0
387    ! pyu1
388#ifdef NC_DOUBLE
389    status = nf_get_vara_double(ncidp, varidyu1, start, count, pyu12)
390#else
391    status = nf_get_vara_real(ncidp, varidyu1, start, count, pyu12)
392#endif
393    CALL gr_ecrit_fi(1, klono, imo, jmo+1, pyu12, pyu1)
394
395    ! pyv1
396#ifdef NC_DOUBLE
397    status = nf_get_vara_double(ncidp, varidyv1, start, count, pyv12)
398#else
399    status = nf_get_vara_real(ncidp, varidyv1, start, count, pyv12)
400#endif
401    CALL gr_ecrit_fi(1, klono, imo, jmo+1, pyv12, pyv1)
402
403    ! *** Temperature au sol ********************************************
404    ! ftsol1
405#ifdef NC_DOUBLE
406    status = nf_get_vara_double(ncidp, varidfts1, start, count, ftsol12)
407#else
408    status = nf_get_vara_real(ncidp, varidfts1, start, count, ftsol12)
409#endif
410    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol12, ftsol1)
411
412    ! ftsol2
413#ifdef NC_DOUBLE
414    status = nf_get_vara_double(ncidp, varidfts2, start, count, ftsol22)
415#else
416    status = nf_get_vara_real(ncidp, varidfts2, start, count, ftsol22)
417#endif
418    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol22, ftsol2)
419
420    ! ftsol3
421#ifdef NC_DOUBLE
422    status = nf_get_vara_double(ncidp, varidfts3, start, count, ftsol32)
423#else
424    status = nf_get_vara_real(ncidp, varidfts3, start, count, ftsol32)
425#endif
426    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol32, ftsol3)
427
428    ! ftsol4
429#ifdef NC_DOUBLE
430    status = nf_get_vara_double(ncidp, varidfts4, start, count, ftsol42)
431#else
432    status = nf_get_vara_real(ncidp, varidfts4, start, count, ftsol42)
433#endif
434    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol42, ftsol4)
435
436    ! *** Nature du sol **************************************************
437    ! psrf1
438#ifdef NC_DOUBLE
439    status = nf_get_vara_double(ncidp, varidpsr1, start, count, psrf12)
440#else
441    status = nf_get_vara_real(ncidp, varidpsr1, start, count, psrf12)
442#endif
443    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf12, psrf1)
444
445    ! psrf2
446#ifdef NC_DOUBLE
447    status = nf_get_vara_double(ncidp, varidpsr2, start, count, psrf22)
448#else
449    status = nf_get_vara_real(ncidp, varidpsr2, start, count, psrf22)
450#endif
451    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf22, psrf2)
452
453    ! psrf3
454#ifdef NC_DOUBLE
455    status = nf_get_vara_double(ncidp, varidpsr3, start, count, psrf32)
456#else
457    status = nf_get_vara_real(ncidp, varidpsr3, start, count, psrf32)
458#endif
459    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf32, psrf3)
460
461    ! psrf4
462#ifdef NC_DOUBLE
463    status = nf_get_vara_double(ncidp, varidpsr4, start, count, psrf42)
464#else
465    status = nf_get_vara_real(ncidp, varidpsr4, start, count, psrf42)
466#endif
467    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf42, psrf4)
468
469    DO i = 1, klono
470
471      psrf(i, 1) = psrf1(i)
472      psrf(i, 2) = psrf2(i)
473      psrf(i, 3) = psrf3(i)
474      psrf(i, 4) = psrf4(i)
475
476      ftsol(i, 1) = ftsol1(i)
477      ftsol(i, 2) = ftsol2(i)
478      ftsol(i, 3) = ftsol3(i)
479      ftsol(i, 4) = ftsol4(i)
480
481    END DO
482
483  END IF
484
485  RETURN
486
487END SUBROUTINE read_pstoke
488
Note: See TracBrowser for help on using the repository browser.