source: LMDZ6/branches/contrails/libf/phylmdiso/phyredem.F90 @ 5450

Last change on this file since 5450 was 5296, checked in by abarral, 8 weeks ago

Turn compbl.h into a module
Move calcul_REGDYN.h to obsolete
Create phys_constants_mod.f90

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