source: LMDZ5/trunk/libf/phylmd/read_pstoke0.F90 @ 2216

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