source: LMDZ5/trunk/libf/phylmd/read_pstoke.F90 @ 2172

Last change on this file since 2172 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
RevLine 
[1992]1
[1403]2! $Id: read_pstoke.F90 1992 2014-03-05 13:19:12Z jescribano $
[524]3
4
[1403]5
[1992]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)
[524]9
[1992]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  ! ******************************************************************************
[524]18
[1992]19  USE netcdf
20  USE dimphy
21  USE control_mod
22  USE indice_sol_mod
[524]23
[1992]24  IMPLICIT NONE
[524]25
[1992]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"
[524]37
[1992]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)
[524]43
[1992]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)
[524]49
[1992]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)
[524]56
[1992]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
[524]62
[1992]63  REAL airefi(klono)
64  CHARACTER *20 namedim
[524]65
[1992]66  ! !! attention !!
67  ! attention il y a aussi le pb de def klono
68  ! dim de phis??
[524]69
70
[1992]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
[524]101
[1992]102  INTEGER l, i
103  INTEGER start(4), count(4), status
104  REAL rcode
105  LOGICAL first
106  SAVE first
107  DATA first/.TRUE./
[524]108
109
110
[1992]111  ! ---------------------------------------------
112  ! Initialisation de la lecture des fichiers
113  ! ---------------------------------------------
[541]114
[1992]115  IF (irec==0) THEN
[524]116
[1992]117    rcode = nf90_open('phystoke.nc', nf90_nowrite, ncidp)
[524]118
[1992]119    rcode = nf90_inq_varid(ncidp, 'phis', varidps)
120    PRINT *, 'ncidp,varidps', ncidp, varidps
[524]121
[1992]122    rcode = nf90_inq_varid(ncidp, 'sig_s', varidpl)
123    PRINT *, 'ncidp,varidpl', ncidp, varidpl
[524]124
[1992]125    rcode = nf90_inq_varid(ncidp, 'aire', varidai)
126    PRINT *, 'ncidp,varidai', ncidp, varidai
[541]127
[1992]128    ! A FAIRE: Es-il necessaire de stocke t?
129    rcode = nf90_inq_varid(ncidp, 't', varidt)
130    PRINT *, 'ncidp,varidt', ncidp, varidt
[541]131
[1992]132    rcode = nf90_inq_varid(ncidp, 'mfu', varidmfu)
133    PRINT *, 'ncidp,varidmfu', ncidp, varidmfu
[524]134
[1992]135    rcode = nf90_inq_varid(ncidp, 'mfd', varidmfd)
136    PRINT *, 'ncidp,varidmfd', ncidp, varidmfd
[524]137
[1992]138    rcode = nf90_inq_varid(ncidp, 'en_u', varidenu)
139    PRINT *, 'ncidp,varidenu', ncidp, varidenu
[524]140
[1992]141    rcode = nf90_inq_varid(ncidp, 'de_u', variddeu)
142    PRINT *, 'ncidp,variddeu', ncidp, variddeu
[524]143
[1992]144    rcode = nf90_inq_varid(ncidp, 'en_d', varidend)
145    PRINT *, 'ncidp,varidend', ncidp, varidend
[524]146
[1992]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
[541]224#ifdef NC_DOUBLE
[1992]225    status = nf_get_vara_double(ncidp, varidpl, 1, zklevo, pl)
[541]226#else
[1992]227    status = nf_get_vara_real(ncidp, varidpl, 1, zklevo, pl)
[541]228#endif
[524]229
[1992]230    ! lecture de aire et phis
[524]231
[1992]232    start(1) = 1
233    start(2) = 1
234    start(3) = 1
235    start(4) = 0
[524]236
[1992]237    count(1) = zim
238    count(2) = zjm
239    count(3) = 1
240    count(4) = 0
241
242    ! phis
[541]243#ifdef NC_DOUBLE
[1992]244    status = nf_get_vara_double(ncidp, varidps, start, count, phisfi2)
[541]245#else
[1992]246    status = nf_get_vara_real(ncidp, varidps, start, count, phisfi2)
[541]247#endif
[1992]248    CALL gr_ecrit_fi(1, klono, imo, jmo+1, phisfi2, phisfi)
[524]249
[1992]250    ! aire
[541]251#ifdef NC_DOUBLE
[1992]252    status = nf_get_vara_double(ncidp, varidai, start, count, airefi2)
[541]253#else
[1992]254    status = nf_get_vara_real(ncidp, varidai, start, count, airefi2)
[541]255#endif
[1992]256    CALL gr_ecrit_fi(1, klono, imo, jmo+1, airefi2, airefi)
257  ELSE
[524]258
[1992]259    PRINT *, 'ok1'
[524]260
[1992]261    ! ---------------------
262    ! lecture des champs
263    ! ---------------------
[524]264
[1992]265    PRINT *, 'WARNING!!! Il n y a pas de test de coherence'
266    PRINT *, 'sur le nombre de niveaux verticaux dans le fichier nc'
[524]267
[1992]268    start(1) = 1
269    start(2) = 1
270    start(3) = 1
271    start(4) = irec
[524]272
[1992]273    count(1) = zim
274    count(2) = zjm
275    count(3) = zklevo
276    count(4) = 1
[541]277
[1992]278
279    ! *** Lessivage******************************************************
280    ! frac_impa
[541]281#ifdef NC_DOUBLE
[1992]282    status = nf_get_vara_double(ncidp, varidfi, start, count, frac_impa2)
[541]283#else
[1992]284    status = nf_get_vara_real(ncidp, varidfi, start, count, frac_impa2)
[541]285#endif
[1992]286    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, frac_impa2, frac_impa)
[524]287
[1992]288    ! frac_nucl
[541]289#ifdef NC_DOUBLE
[1992]290    status = nf_get_vara_double(ncidp, varidfn, start, count, frac_nucl2)
[541]291#else
[1992]292    status = nf_get_vara_real(ncidp, varidfn, start, count, frac_nucl2)
[541]293#endif
[1992]294    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, frac_nucl2, frac_nucl)
[524]295
[1992]296    ! *** Temperature ******************************************************
297    ! abder t
[541]298#ifdef NC_DOUBLE
[1992]299    status = nf_get_vara_double(ncidp, varidt, start, count, t2)
[541]300#else
[1992]301    status = nf_get_vara_real(ncidp, varidt, start, count, t2)
[541]302#endif
[1992]303    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, t2, t)
[524]304
[1992]305    ! *** Flux pour le calcul de la convection TIEDTK ***********************
306    ! mfu
[541]307#ifdef NC_DOUBLE
[1992]308    status = nf_get_vara_double(ncidp, varidmfu, start, count, mfu2)
[541]309#else
[1992]310    status = nf_get_vara_real(ncidp, varidmfu, start, count, mfu2)
[541]311#endif
[1992]312    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, mfu2, mfu)
[524]313
[1992]314    ! mfd
[541]315#ifdef NC_DOUBLE
[1992]316    status = nf_get_vara_double(ncidp, varidmfd, start, count, mfd2)
[541]317#else
[1992]318    status = nf_get_vara_real(ncidp, varidmfd, start, count, mfd2)
[541]319#endif
[1992]320    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, mfd2, mfd)
[524]321
[1992]322    ! en_u
[541]323#ifdef NC_DOUBLE
[1992]324    status = nf_get_vara_double(ncidp, varidenu, start, count, en_u2)
[541]325#else
[1992]326    status = nf_get_vara_real(ncidp, varidenu, start, count, en_u2)
[541]327#endif
[1992]328    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_u2, en_u)
[524]329
[1992]330    ! de_u
[541]331#ifdef NC_DOUBLE
[1992]332    status = nf_get_vara_double(ncidp, variddeu, start, count, de_u2)
[541]333#else
[1992]334    status = nf_get_vara_real(ncidp, variddeu, start, count, de_u2)
[541]335#endif
[1992]336    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, de_u2, de_u)
[524]337
[1992]338    ! en_d
[541]339#ifdef NC_DOUBLE
[1992]340    status = nf_get_vara_double(ncidp, varidend, start, count, en_d2)
[541]341#else
[1992]342    status = nf_get_vara_real(ncidp, varidend, start, count, en_d2)
[541]343#endif
[1992]344    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_d2, en_d)
[524]345
[1992]346    ! de_d
[541]347#ifdef NC_DOUBLE
[1992]348    status = nf_get_vara_double(ncidp, varidded, start, count, de_d2)
[541]349#else
[1992]350    status = nf_get_vara_real(ncidp, varidded, start, count, de_d2)
[541]351#endif
[1992]352    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, de_d2, de_d)
[524]353
[1992]354    ! **** Coeffecient du mellange
355    ! turbulent**********************************
356    ! coefh
[541]357#ifdef NC_DOUBLE
[1992]358    status = nf_get_vara_double(ncidp, varidch, start, count, coefh2)
[541]359#else
[1992]360    status = nf_get_vara_real(ncidp, varidch, start, count, coefh2)
[541]361#endif
[1992]362    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, coefh2, coefh)
[524]363
[1992]364    ! *** Flux ascendant et entrant pour les
365    ! Thermiques************************
366    ! abder thermiques
[541]367#ifdef NC_DOUBLE
[1992]368    status = nf_get_vara_double(ncidp, varidfmth, start, count, fm_therm2)
[541]369#else
[1992]370    status = nf_get_vara_real(ncidp, varidfmth, start, count, fm_therm2)
[541]371#endif
[1992]372    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, fm_therm2, fm_therm)
[541]373
374#ifdef NC_DOUBLE
[1992]375    status = nf_get_vara_double(ncidp, varidenth, start, count, en_therm2)
[541]376#else
[1992]377    status = nf_get_vara_real(ncidp, varidenth, start, count, en_therm2)
[541]378#endif
[1992]379    CALL gr_ecrit_fi(klevo, klono, imo, jmo+1, en_therm2, en_therm)
[541]380
[1992]381    ! *** Vitesses aux sol
382    ! ******************************************************
383    start(3) = irec
384    start(4) = 0
385    count(3) = 1
386    count(4) = 0
387    ! pyu1
[541]388#ifdef NC_DOUBLE
[1992]389    status = nf_get_vara_double(ncidp, varidyu1, start, count, pyu12)
[541]390#else
[1992]391    status = nf_get_vara_real(ncidp, varidyu1, start, count, pyu12)
[541]392#endif
[1992]393    CALL gr_ecrit_fi(1, klono, imo, jmo+1, pyu12, pyu1)
[524]394
[1992]395    ! pyv1
[541]396#ifdef NC_DOUBLE
[1992]397    status = nf_get_vara_double(ncidp, varidyv1, start, count, pyv12)
[541]398#else
[1992]399    status = nf_get_vara_real(ncidp, varidyv1, start, count, pyv12)
[541]400#endif
[1992]401    CALL gr_ecrit_fi(1, klono, imo, jmo+1, pyv12, pyv1)
[524]402
[1992]403    ! *** Temperature au sol ********************************************
404    ! ftsol1
[541]405#ifdef NC_DOUBLE
[1992]406    status = nf_get_vara_double(ncidp, varidfts1, start, count, ftsol12)
[541]407#else
[1992]408    status = nf_get_vara_real(ncidp, varidfts1, start, count, ftsol12)
[541]409#endif
[1992]410    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol12, ftsol1)
[524]411
[1992]412    ! ftsol2
[541]413#ifdef NC_DOUBLE
[1992]414    status = nf_get_vara_double(ncidp, varidfts2, start, count, ftsol22)
[541]415#else
[1992]416    status = nf_get_vara_real(ncidp, varidfts2, start, count, ftsol22)
[541]417#endif
[1992]418    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol22, ftsol2)
[524]419
[1992]420    ! ftsol3
[541]421#ifdef NC_DOUBLE
[1992]422    status = nf_get_vara_double(ncidp, varidfts3, start, count, ftsol32)
[541]423#else
[1992]424    status = nf_get_vara_real(ncidp, varidfts3, start, count, ftsol32)
[541]425#endif
[1992]426    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol32, ftsol3)
[524]427
[1992]428    ! ftsol4
[541]429#ifdef NC_DOUBLE
[1992]430    status = nf_get_vara_double(ncidp, varidfts4, start, count, ftsol42)
[541]431#else
[1992]432    status = nf_get_vara_real(ncidp, varidfts4, start, count, ftsol42)
[541]433#endif
[1992]434    CALL gr_ecrit_fi(1, klono, imo, jmo+1, ftsol42, ftsol4)
[524]435
[1992]436    ! *** Nature du sol **************************************************
437    ! psrf1
[541]438#ifdef NC_DOUBLE
[1992]439    status = nf_get_vara_double(ncidp, varidpsr1, start, count, psrf12)
[541]440#else
[1992]441    status = nf_get_vara_real(ncidp, varidpsr1, start, count, psrf12)
[541]442#endif
[1992]443    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf12, psrf1)
[524]444
[1992]445    ! psrf2
[541]446#ifdef NC_DOUBLE
[1992]447    status = nf_get_vara_double(ncidp, varidpsr2, start, count, psrf22)
[541]448#else
[1992]449    status = nf_get_vara_real(ncidp, varidpsr2, start, count, psrf22)
[541]450#endif
[1992]451    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf22, psrf2)
[524]452
[1992]453    ! psrf3
[541]454#ifdef NC_DOUBLE
[1992]455    status = nf_get_vara_double(ncidp, varidpsr3, start, count, psrf32)
[541]456#else
[1992]457    status = nf_get_vara_real(ncidp, varidpsr3, start, count, psrf32)
[541]458#endif
[1992]459    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf32, psrf3)
[524]460
[1992]461    ! psrf4
[541]462#ifdef NC_DOUBLE
[1992]463    status = nf_get_vara_double(ncidp, varidpsr4, start, count, psrf42)
[541]464#else
[1992]465    status = nf_get_vara_real(ncidp, varidpsr4, start, count, psrf42)
[541]466#endif
[1992]467    CALL gr_ecrit_fi(1, klono, imo, jmo+1, psrf42, psrf4)
[524]468
[1992]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.