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

Last change on this file since 4231 was 4170, checked in by dcugnet, 2 years ago

The variable "types_trac" is the equivalent of "type_trac" in case multiple sections must be read
and used in "tracer.def" file.
Tests on the "type_trac" were replaced with tests on the vector "types_trac".
Most of the time, there are two components: 'lmdz' and a second one. The later has priority on 'lmdz'
and must be used for the tests. For more components, care must be taken to execute specific parts
of the code on the right tracers ; the tracers(:)%component has been created in that respect.

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