source: LMDZ6/trunk/libf/phylmdiso/phyredem.F90 @ 3927

Last change on this file since 3927 was 3927, checked in by Laurent Fairhead, 3 years ago

Initial import of the physics wih isotopes from Camille Risi
CR

File size: 19.4 KB
Line 
1!
2! $Id: phyredem.F90 3506 2019-05-16 14:38:11Z ymeurdesoif $
3!
4SUBROUTINE phyredem (fichnom)
5!
6!-------------------------------------------------------------------------------
7! Author: Z.X. Li (LMD/CNRS), 1993/08/18
8!-------------------------------------------------------------------------------
9! Purpose: Write restart state for physics.
10!-------------------------------------------------------------------------------
11  USE dimphy, ONLY: klon, klev
12  USE fonte_neige_mod,  ONLY : fonte_neige_final
13  USE pbl_surface_mod,  ONLY : pbl_surface_final
14  USE phys_state_var_mod, ONLY: radpas, zmasq, pctsrf, ftsol, falb_dir,      &
15                                falb_dif, qsol, fevap, radsol, solsw, sollw, &
16                                sollwdown, rain_fall, snow_fall, z0m, z0h,   &
17                                agesno, zmea, zstd, zsig, zgam, zthe, zpic,  &
18                                zval, rugoro, t_ancien, q_ancien,            &
19                                prw_ancien, prlw_ancien, prsw_ancien,        &
20                                ql_ancien, qs_ancien,  u_ancien,             &
21                                v_ancien, clwcon, rnebcon, ratqs, pbl_tke,   &
22                                wake_delta_pbl_tke, zmax0, f0, sig1, w01,    &
23                                wake_deltat, wake_deltaq, wake_s, wake_dens, &
24                                wake_cstar,                                  &
25                                wake_pe, wake_fip, fm_therm, entr_therm,     &
26                                detr_therm, ale_bl, ale_bl_trig, alp_bl,     &
27                                ale_wake, ale_bl_stat,                       &
28                                du_gwd_rando, du_gwd_front, u10m, v10m,      &
29                                treedrg
30#ifdef ISO
31  USE phys_state_var_mod, ONLY: xtsol, fxtevap,xtrain_fall, xtsnow_fall,     &
32                                xt_ancien, xtl_ancien, xts_ancien,           &
33                                wake_deltaxt                             
34#endif
35  USE geometry_mod, ONLY : longitude_deg, latitude_deg
36  USE iostart, ONLY: open_restartphy, close_restartphy, enddef_restartphy, put_field, put_var
37  USE traclmdz_mod, ONLY : traclmdz_to_restart
38  USE infotrac_phy, ONLY: type_trac, niadv, tname, nbtr, nqo,itr_indice
39#ifdef ISO
40  USE infotrac_phy, ONLY: itr_indice,niso,ntraciso
41#ifdef ISOVERIF
42  USE isotopes_verif_mod
43#endif
44#endif
45  USE carbon_cycle_mod, ONLY : carbon_cycle_cpl, co2_send
46  USE indice_sol_mod, ONLY: nbsrf, is_oce, is_sic, is_ter, is_lic, epsfra
47  USE surface_data, ONLY: type_ocean, version_ocean
48  USE ocean_slab_mod, ONLY : nslay, tslab, seaice, tice, fsic
49  USE time_phylmdz_mod, ONLY: annee_ref, day_end, itau_phy, pdtphys
50
51  IMPLICIT none
52
53  include "dimsoil.h"
54  include "clesphys.h"
55  include "thermcell.h"
56  include "compbl.h"
57  !======================================================================
58  CHARACTER*(*) fichnom
59
60  ! les variables globales ecrites dans le fichier restart
61
62  REAL tsoil(klon, nsoilmx, nbsrf)
63  REAL qsurf(klon, nbsrf)
64  REAL snow(klon, nbsrf)
65  real fder(klon)
66  REAL run_off_lic_0(klon)
67  REAL trs(klon, nbtr)
68#ifdef ISO
69  REAL xtsnow(niso,klon, nbsrf)
70  REAL xtrun_off_lic_0(niso,klon)
71  REAL Rland_ice(niso,klon)
72#endif
73  INTEGER iq ! C Risi
74
75  INTEGER nid, nvarid, idim1, idim2, idim3
76  INTEGER ierr
77  INTEGER length
78  PARAMETER (length=100)
79  REAL tab_cntrl(length)
80
81  INTEGER isoil, nsrf,isw
82  CHARACTER (len=2) :: str2
83  CHARACTER (len=256) :: nam, lnam
84  INTEGER           :: it, iiq, pass
85
86  !======================================================================
87
88  ! Get variables which will be written to restart file from module
89  ! pbl_surface_mod
90  CALL pbl_surface_final(fder, snow, qsurf,  tsoil &
91#ifdef ISO
92       ,xtsnow,Rland_ice &
93#endif       
94       )
95
96  ! Get a variable calculated in module fonte_neige_mod
97  CALL fonte_neige_final(run_off_lic_0 &
98#ifdef ISO
99       ,xtrun_off_lic_0 &
100#endif       
101       )
102
103  !======================================================================
104
105  CALL open_restartphy(fichnom)
106
107 
108  DO ierr = 1, length
109     tab_cntrl(ierr) = 0.0
110  ENDDO
111  tab_cntrl(1) = pdtphys
112  tab_cntrl(2) = radpas
113  ! co2_ppm : current value of atmospheric CO2
114  tab_cntrl(3) = co2_ppm
115  tab_cntrl(4) = solaire
116  tab_cntrl(5) = iflag_con
117  tab_cntrl(6) = nbapp_rad
118
119  IF( iflag_cycle_diurne.GE.1 ) tab_cntrl( 7 ) = iflag_cycle_diurne
120  IF(   soil_model ) tab_cntrl( 8 ) = 1.
121  IF(     new_oliq ) tab_cntrl( 9 ) = 1.
122  IF(     ok_orodr ) tab_cntrl(10 ) = 1.
123  IF(     ok_orolf ) tab_cntrl(11 ) = 1.
124
125  tab_cntrl(13) = day_end
126  tab_cntrl(14) = annee_ref
127  tab_cntrl(15) = itau_phy
128
129  ! co2_ppm0 : initial value of atmospheric CO2
130  tab_cntrl(16) = co2_ppm0
131
132  DO pass=1,2   ! pass=1 netcdf definition ; pass=2 netcdf write
133 
134    CALL put_var(pass, "controle", "Parametres de controle", tab_cntrl)
135
136    CALL put_field(pass,"longitude", &
137         "Longitudes de la grille physique", longitude_deg)
138
139    CALL put_field(pass,"latitude", "Latitudes de la grille physique", latitude_deg)
140
141    ! PB ajout du masque terre/mer
142
143    CALL put_field(pass,"masque", "masque terre mer", zmasq)
144
145    ! BP ajout des fraction de chaque sous-surface
146
147    ! Get last fractions from slab ocean
148    IF (type_ocean == 'slab' .AND. version_ocean == "sicINT") THEN
149        WHERE (1.-zmasq(:).GT.EPSFRA)
150            pctsrf(:,is_oce)=(1.-fsic(:))*(1.-zmasq(:))
151            pctsrf(:,is_sic)=fsic(:)*(1.-zmasq(:))
152        END WHERE
153    END IF
154
155    ! 1. fraction de terre
156
157    CALL put_field(pass,"FTER", "fraction de continent", pctsrf(:, is_ter))
158
159    ! 2. Fraction de glace de terre
160
161    CALL put_field(pass,"FLIC", "fraction glace de terre", pctsrf(:, is_lic))
162
163    ! 3. fraction ocean
164
165    CALL put_field(pass,"FOCE", "fraction ocean", pctsrf(:, is_oce))
166
167    ! 4. Fraction glace de mer
168
169    CALL put_field(pass,"FSIC", "fraction glace mer", pctsrf(:, is_sic))
170
171    IF(nbsrf>99) THEN
172      PRINT*, "Trop de sous-mailles";  CALL abort_physic("phyredem", "", 1)
173    END IF
174    IF(nsoilmx>99) THEN
175      PRINT*, "Trop de sous-surfaces"; CALL abort_physic("phyredem", "", 1)
176    END IF
177    IF(nsw>99) THEN
178      PRINT*, "Trop de bandes"; CALL abort_physic("phyredem", "", 1)
179    END IF
180
181    CALL put_field_srf1(pass,"TS","Temperature",ftsol(:,:))
182
183! ================== Albedo =======================================
184    print*,'PHYREDEM NOUVEAU'
185    CALL put_field_srf2(pass,"A_dir_SW","Albedo direct",falb_dir(:,:,:))
186    CALL put_field_srf2(pass,"A_dif_SW","Albedo diffus",falb_dif(:,:,:))
187
188    CALL put_field_srf1(pass,"U10M", "u a 10m", u10m)
189
190    CALL put_field_srf1(pass,"V10M", "v a 10m", v10m)
191
192
193! ================== Tsoil =========================================
194    CALL put_field_srf2(pass,"Tsoil","Temperature",tsoil(:,:,:))
195!FC
196!  CALL put_field_srf2("treedrg","freinage arbres",treedrg(:,:,:))
197    CALL put_field(pass,"treedrg_ter","freinage arbres",treedrg(:,:,is_ter))
198
199
200    CALL put_field_srf1(pass,"QS"  , "Humidite",qsurf(:,:))
201
202    CALL put_field     (pass,"QSOL", "Eau dans le sol (mm)", qsol)
203
204    CALL put_field_srf1(pass,"EVAP", "Evaporation", fevap(:,:))
205
206    CALL put_field_srf1(pass,"SNOW", "Neige", snow(:,:))
207
208    CALL put_field(pass,"RADS", "Rayonnement net a la surface", radsol)
209
210    CALL put_field(pass,"solsw", "Rayonnement solaire a la surface", solsw)
211
212    CALL put_field(pass,"sollw", "Rayonnement IF a la surface", sollw)
213
214    CALL put_field(pass,"sollwdown", "Rayonnement down IF a la surface", sollwdown)
215
216    CALL put_field(pass,"fder", "Derive de flux", fder)
217
218    CALL put_field(pass,"rain_f", "precipitation liquide", rain_fall)
219
220    CALL put_field(pass,"snow_f", "precipitation solide", snow_fall)
221
222    CALL put_field_srf1(pass,"Z0m", "rugosite", z0m(:,:))
223
224    CALL put_field_srf1(pass,"Z0h", "rugosite", z0h(:,:))
225
226    CALL put_field_srf1(pass,"AGESNO", "Age de la neige", agesno(:,:))
227
228    CALL put_field(pass,"ZMEA", "ZMEA", zmea)
229
230    CALL put_field(pass,"ZSTD", "ZSTD", zstd)
231
232    CALL put_field(pass,"ZSIG", "ZSIG", zsig)
233
234    CALL put_field(pass,"ZGAM", "ZGAM", zgam)
235
236    CALL put_field(pass,"ZTHE", "ZTHE", zthe)
237
238    CALL put_field(pass,"ZPIC", "ZPIC", zpic)
239
240    CALL put_field(pass,"ZVAL", "ZVAL", zval)
241
242    CALL put_field(pass,"RUGSREL", "RUGSREL", rugoro)
243
244    CALL put_field(pass,"TANCIEN", "TANCIEN", t_ancien)
245
246    CALL put_field(pass,"QANCIEN", "QANCIEN", q_ancien)
247
248    CALL put_field(pass,"QLANCIEN", "QLANCIEN", ql_ancien)
249
250    CALL put_field(pass,"QSANCIEN", "QSANCIEN", qs_ancien)
251
252    CALL put_field(pass,"PRWANCIEN", "PRWANCIEN", prw_ancien)
253
254    CALL put_field(pass,"PRLWANCIEN", "PRLWANCIEN", prlw_ancien)
255
256    CALL put_field(pass,"PRSWANCIEN", "PRSWANCIEN", prsw_ancien)
257
258    CALL put_field(pass,"UANCIEN", "UANCIEN", u_ancien)
259
260    CALL put_field(pass,"VANCIEN", "VANCIEN", v_ancien)
261
262    CALL put_field(pass,"CLWCON", "Eau liquide convective", clwcon)
263
264    CALL put_field(pass,"RNEBCON", "Nebulosite convective", rnebcon)
265
266    CALL put_field(pass,"RATQS", "Ratqs", ratqs)
267
268    ! run_off_lic_0
269
270    CALL put_field(pass,"RUNOFFLIC0", "Runofflic0", run_off_lic_0)
271
272    ! DEB TKE PBL !
273
274    IF (iflag_pbl>1) then
275      CALL put_field_srf3(pass,"TKE", "Energ. Cineti. Turb.", &
276           pbl_tke(:,:,:))
277      CALL put_field_srf3(pass,"DELTATKE", "Del TKE wk/env.", &
278           wake_delta_pbl_tke(:,:,:))
279    END IF
280
281    ! FIN TKE PBL !
282    !IM ajout zmax0, f0, sig1, w01
283    !IM wake_deltat, wake_deltaq, wake_s, wake_cstar, wake_pe, wake_fip
284
285    CALL put_field(pass,"ZMAX0", "ZMAX0", zmax0)
286
287    CALL put_field(pass,"F0", "F0", f0)
288
289    CALL put_field(pass,"sig1", "sig1 Emanuel", sig1)
290
291    CALL put_field(pass,"w01", "w01 Emanuel", w01)
292
293    ! wake_deltat
294    CALL put_field(pass,"WAKE_DELTAT", "WAKE_DELTAT", wake_deltat)
295
296    CALL put_field(pass,"WAKE_DELTAQ", "WAKE_DELTAQ", wake_deltaq)
297
298    CALL put_field(pass,"WAKE_S", "Wake frac. area", wake_s)
299
300    CALL put_field(pass,"WAKE_DENS", "Wake num. /unit area", wake_dens)
301
302    CALL put_field(pass,"WAKE_CSTAR", "WAKE_CSTAR", wake_cstar)
303
304    CALL put_field(pass,"WAKE_PE", "WAKE_PE", wake_pe)
305
306    CALL put_field(pass,"WAKE_FIP", "WAKE_FIP", wake_fip)
307
308    ! thermiques
309
310    CALL put_field(pass,"FM_THERM", "FM_THERM", fm_therm)
311
312    CALL put_field(pass,"ENTR_THERM", "ENTR_THERM", entr_therm)
313
314    CALL put_field(pass,"DETR_THERM", "DETR_THERM", detr_therm)
315
316    CALL put_field(pass,"ALE_BL", "ALE_BL", ale_bl)
317
318    CALL put_field(pass,"ALE_BL_TRIG", "ALE_BL_TRIG", ale_bl_trig)
319
320    CALL put_field(pass,"ALP_BL", "ALP_BL", alp_bl)
321
322    CALL put_field(pass,"ALE_WAKE", "ALE_WAKE", ale_wake)
323
324    CALL put_field(pass,"ALE_BL_STAT", "ALE_BL_STAT", ale_bl_stat)
325
326
327    ! trs from traclmdz_mod
328    IF (type_trac == 'lmdz') THEN
329       CALL traclmdz_to_restart(trs)
330       DO it=1, nbtr
331!!        iiq=niadv(it+2)                                                           ! jyg
332          !iiq=niadv(it+nqo) ! C Risi: on efface pour remplacer:
333          iq=itr_indice(it)                                                           ! jyg
334          iiq=niadv(iq)                                                           ! jyg
335          CALL put_field(pass,"trs_"//tname(iiq), "", trs(:, it))
336       END DO
337       IF (carbon_cycle_cpl) THEN
338          IF (.NOT. ALLOCATED(co2_send)) THEN
339             ! This is the case of create_etat0_limit, ce0l
340             ALLOCATE(co2_send(klon))
341             co2_send(:) = co2_ppm0
342          END IF
343          CALL put_field(pass,"co2_send", "co2_ppm for coupling", co2_send)
344       END IF
345    END IF
346
347    ! Restart variables for Slab ocean
348    IF (type_ocean == 'slab') THEN
349        IF (nslay.EQ.1) THEN
350          CALL put_field(pass,"tslab", "Slab ocean temperature", tslab)
351        ELSE
352          DO it=1,nslay
353            WRITE(str2,'(i2.2)') it
354            CALL put_field(pass,"tslab"//str2, "Slab ocean temperature", tslab(:,it))
355          END DO
356        END IF
357        IF (version_ocean == 'sicINT') THEN
358            CALL put_field(pass,"seaice", "Slab seaice (kg/m2)", seaice)
359            CALL put_field(pass,"slab_tice", "Slab sea ice temperature", tice)
360        END IF
361    END IF
362
363    if (ok_gwd_rando) call put_field(pass,"du_gwd_rando", &
364         "tendency on zonal wind due to flott gravity waves", du_gwd_rando)
365
366    IF (.not. ok_hines .and. ok_gwd_rando) call put_field(pass,"du_gwd_front", &
367         "tendency on zonal wind due to acama gravity waves", du_gwd_front)
368
369#ifdef ISO
370      write(*,*) 'phyredem 342'
371      call phyisoredem (pass, &
372     &           xtsnow, &
373     &           xtrun_off_lic_0,Rland_ice, &
374     &           run_off_lic_0)
375#endif
376
377    IF (pass==1) CALL enddef_restartphy
378    IF (pass==2) CALL close_restartphy
379  ENDDO ! DO pass=1,2   ! pass=1 netcdf definition ; pass=2 netcdf write
380 
381
382 
383  !$OMP BARRIER
384
385
386  CONTAINS
387
388
389SUBROUTINE put_field_srf1(pass,nam,lnam,field)
390
391  IMPLICIT NONE
392  INTEGER, INTENT(IN)            :: pass
393  CHARACTER(LEN=*),  INTENT(IN) :: nam, lnam
394  REAL,              INTENT(IN) :: field(:,:)
395  CHARACTER(LEN=256) :: nm, lm, str
396  DO nsrf = 1, SIZE(field,2)
397    WRITE(str, '(i2.2)') nsrf
398    nm=TRIM(nam)//TRIM(str)
399    lm=TRIM(lnam)//" de surface No. "//TRIM(str)
400    CALL put_field(pass,nm,lm,field(:,nsrf))
401  END DO
402
403END SUBROUTINE put_field_srf1
404
405
406SUBROUTINE put_field_srf2(pass,nam,lnam,field)
407
408  IMPLICIT NONE
409  INTEGER, INTENT(IN)            :: pass
410  CHARACTER(LEN=*),  INTENT(IN) :: nam, lnam
411  REAL,              INTENT(IN) :: field(:,:,:)
412  CHARACTER(LEN=256) :: nm, lm, str
413  DO nsrf = 1, SIZE(field,3)
414    DO isoil=1, SIZE(field,2)
415      WRITE(str, '(i2.2,"srf",i2.2)')isoil,nsrf
416!      WRITE(lunout,*)"PHYREDEM ",TRIM(nam)//TRIM(str)
417      nm=TRIM(nam)//TRIM(str)
418      lm=TRIM(lnam)//" du sol No. "//TRIM(str)
419      CALL put_field(pass,nm,lm,field(:,isoil,nsrf))
420    END DO
421  END DO
422
423END SUBROUTINE put_field_srf2
424
425
426SUBROUTINE put_field_srf3(pass,nam,lnam,field)
427
428  IMPLICIT NONE
429  INTEGER, INTENT(IN)            :: pass
430  CHARACTER(LEN=*),  INTENT(IN) :: nam, lnam
431  REAL,              INTENT(IN) :: field(:,:,:)
432  CHARACTER(LEN=256) :: nm, lm, str
433  DO nsrf = 1, SIZE(field,3)
434    WRITE(str, '(i2.2)') nsrf
435    nm=TRIM(nam)//TRIM(str)
436    lm=TRIM(lnam)//TRIM(str)
437    CALL put_field(pass,nm,lm,field(:,1:klev+1,nsrf))
438  END DO
439
440END SUBROUTINE put_field_srf3
441
442#ifdef ISO
443! je voulais mettre cette subroutine dans isotopes_mod, mais elle a besoin de put_field_srf1 qui est contenue dans la subroutine phyredem. Si on veut mettre cette routine dans isotopes_mod, il faudrait convertir ce fichier en module pour pouvoir en appeler des routines
444
445      SUBROUTINE phyisoredem (pass, &
446     &          xtsnow, &
447     &          xtrun_off_lic_0,Rland_ice, &
448     &          run_off_lic_0)
449      USE dimphy
450      !USE mod_grid_phy_lmdz
451      !USE mod_phys_lmdz_para
452      USE phys_state_var_mod, ONLY: q_ancien,xt_ancien,wake_deltaq,wake_deltaxt, &
453        xtrain_fall,xtsnow_fall, ql_ancien,xtl_ancien,qs_ancien,xts_ancien, &
454        xtsol,fxtevap
455      USE infotrac_phy,ONLY: niso, ntraciso
456      !USE control_mod
457      USE indice_sol_mod, ONLY: nbsrf
458      USE iostart, ONLY: put_field
459      USE isotopes_mod, ONLY: striso,iso_eau
460#ifdef ISOVERIF
461      USE isotopes_verif_mod
462#endif
463#ifdef ISOTRAC
464    use isotrac_mod, only: index_zone,index_iso,strtrac
465#endif
466      !USE phyredem, ONLY: put_field_srf1
467
468        implicit none
469
470        ! equivalent isotopique de phyredem
471
472#include "dimsoil.h"
473#include "clesphys.h"
474#include "thermcell.h"
475#include "compbl.h"     
476      ! inputs
477      !REAL xtsol(niso,klon)
478      REAL xtsnow(niso,klon,nbsrf)
479      !REAL xtevap(ntraciso,klon,nbsrf)     
480      REAL xtrun_off_lic_0(niso,klon)
481      REAL Rland_ice(niso,klon)
482      real run_off_lic_0(klon)
483      integer, intent(in) :: pass
484
485      ! locals
486      real iso_tmp(klon)
487      real iso_tmp_lonlev(klon,klev)
488      real iso_tmp_lonsrf(klon,nbsrf)
489      integer i,ixt,k,nsrf
490      INTEGER nid, nvarid
491      INTEGER ierr
492      CHARACTER*7 str7
493      CHARACTER*2 str2
494      CHARACTER*50 striso_sortie
495      integer lnblnk
496#ifdef ISOTRAC
497      integer iiso,izone
498#endif     
499
500      write(*,*) 'phyisoredem 41: entrée'
501#ifdef ISOVERIF
502     if (iso_eau.gt.0) then
503      do k=1,klev
504        do i=1,klon
505           call iso_verif_egalite(xt_ancien(iso_eau,i,k),q_ancien(i,k), &
506     &           'phyisoredem 50a')
507           call iso_verif_egalite(xtl_ancien(iso_eau,i,k),ql_ancien(i,k), &
508     &           'phyisoredem 50b')
509           call iso_verif_egalite(xts_ancien(iso_eau,i,k),qs_ancien(i,k), &
510     &           'phyisoredem 50c')
511         
512        enddo !do i=1,klon
513      enddo !do k=1,klev
514      do i=1,klon
515        DO nsrf = 1, nbsrf
516           call iso_verif_egalite(fxtevap(iso_eau,i,nsrf),fevap(i,nsrf), &
517     &           'phyisoredem 50d')
518        enddo !DO nsrf = 1, nbsrf
519       enddo
520      endif !if (iso_eau.gt.0) then
521      do i=1,klon
522       do ixt=1,niso
523        call iso_verif_noNaN(xtsol(ixt,i),'phyisoredem 72')
524       enddo !do ixt=1,niso
525      enddo !do i=1,klon
526#ifdef ISOTRAC       
527      do k=1,klev
528        do i=1,klon 
529          call iso_verif_traceur(xt_ancien(1,i,k), &
530     &                   'phyisoredem 60') 
531        enddo !do i=1,klon
532      enddo !do k=1,kle
533#endif
534#endif
535
536   do ixt=1,ntraciso
537
538     if (ixt.le.niso) then
539        striso_sortie=striso(ixt)
540     else
541#ifdef ISOTRAC
542        iiso=index_iso(ixt)
543        izone=index_zone(ixt)       
544        striso_sortie=striso(iiso)//strtrac(izone)
545#else
546        write(*,*) 'phyredem 546: ixt,ntraciso=', ixt,ntraciso
547        stop
548#endif
549     endif !if (ixt.le.niso) then
550      write(*,*) 'phyredem 550: ixt,striso_sortie=',ixt,striso_sortie(1:lnblnk(striso_sortie))
551     
552      iso_tmp_lonsrf(:,:)=fxtevap(ixt,:,:)
553      CALL put_field_srf1(pass,"XTEVAP"//striso_sortie(1:lnblnk(striso_sortie)), &
554     &           "Evaporation de surface",iso_tmp_lonsrf)
555
556      iso_tmp_lonsrf(:,:)=xtsnow(ixt,:,:)
557      CALL put_field_srf1(pass,"XTSNOW"//striso_sortie(1:lnblnk(striso_sortie)), &
558     &           "NEIGE",iso_tmp_lonsrf)       
559
560      iso_tmp(:)=xtrain_fall(ixt,:)
561      CALL put_field(pass,"xtrain_f"//striso_sortie(1:lnblnk(striso_sortie)), &
562     &   "precipitation liquide",iso_tmp)
563
564      iso_tmp(:)=xtsnow_fall(ixt,:)
565      CALL put_field(pass,"xtsnow_f"//striso_sortie(1:lnblnk(striso_sortie)), &
566     &    "precipitation solide",iso_tmp)
567
568      iso_tmp_lonlev(:,:)=xt_ancien(ixt,:,:)
569      CALL put_field(pass,"XTANCIEN"//striso_sortie(1:lnblnk(striso_sortie)), &
570     &   "QANCIEN",iso_tmp_lonlev)
571
572      iso_tmp_lonlev(:,:)=xtl_ancien(ixt,:,:)
573      CALL put_field(pass,"XTLANCIEN"//striso_sortie(1:lnblnk(striso_sortie)), &
574     &   "QLANCIEN",iso_tmp_lonlev)
575
576      iso_tmp_lonlev(:,:)=xts_ancien(ixt,:,:)
577      CALL put_field(pass,"XTSANCIEN"//striso_sortie(1:lnblnk(striso_sortie)), &
578     &   "QSANCIEN",iso_tmp_lonlev)
579
580      iso_tmp_lonlev(:,:)=wake_deltaxt(ixt,:,:)
581      CALL put_field(pass,"WAKE_DELTAXT"//striso_sortie(1:lnblnk(striso_sortie)), &
582     &   "WAKE_DELTAQ",iso_tmp_lonlev)
583
584      iso_tmp(:)=xtrun_off_lic_0(ixt,:)
585      CALL put_field(pass,"XTRUNOFFLIC0"//striso_sortie(1:lnblnk(striso_sortie)), &
586     &   "Runofflic0",iso_tmp)
587
588      iso_tmp_lonlev(:,:)=wake_deltaxt(ixt,:,:)
589      CALL put_field(pass,"WAKE_DELTAXT"//striso_sortie(1:lnblnk(striso_sortie)), &
590     &   "WAKE_DELTAXT",iso_tmp_lonlev)
591
592      ! variables seulement pour niso:
593      if (ixt.le.niso) then
594
595      iso_tmp(:)=xtsol(ixt,:)
596      CALL put_field(pass,"XTSOL"//striso_sortie(1:lnblnk(striso_sortie)), &
597     &   "Eau dans le sol (mm)",iso_tmp)
598
599      iso_tmp(:)=Rland_ice(ixt,:)
600      CALL put_field(pass,"Rland_ice"//striso_sortie(1:lnblnk(striso_sortie)), &
601     &   "ratio land ice",iso_tmp)
602
603      endif ! if (ixt.le.niso) then
604
605      enddo !do ixt=1,niso
606
607      write(*,*) 'phyisoredem 261: sortie'
608      END SUBROUTINE phyisoredem
609#endif
610
611END SUBROUTINE phyredem
Note: See TracBrowser for help on using the repository browser.