source: LMDZ6/branches/LMDZ_ECRad/libf/phylmd/iophy.F90 @ 5403

Last change on this file since 5403 was 4727, checked in by idelkadi, 14 months ago

Merged trunk changes -r4488:4726 LMDZ_ECRad branch

  • 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 Id
File size: 52.7 KB
Line 
1!
2! $Id: iophy.F90 4727 2023-10-19 14:02:57Z abarral $
3!
4MODULE iophy
5
6! abd  REAL,private,allocatable,DIMENSION(:),save :: io_lat
7! abd  REAL,private,allocatable,DIMENSION(:),save :: io_lon
8  REAL,ALLOCATABLE,DIMENSION(:),SAVE :: io_lat
9  REAL,ALLOCATABLE,DIMENSION(:),SAVE :: io_lon
10  INTEGER, SAVE :: phys_domain_id
11  INTEGER, SAVE :: npstn
12  INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: nptabij
13  INTEGER, SAVE :: itau_iophy
14  LOGICAL :: check_dim = .false.
15
16!$OMP THREADPRIVATE(itau_iophy)
17
18  INTERFACE histwrite_phy
19    MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy,histwrite2d_phy_old,histwrite3d_phy_old,histwrite2d_xios,histwrite3d_xios,histwrite0d_xios
20  END INTERFACE
21
22  INTERFACE histbeg_phy_all
23    MODULE PROCEDURE histbeg_phy,histbeg_phyxios,histbeg_phy_points
24  END INTERFACE
25
26
27CONTAINS
28
29! ug Routine pour définir itau_iophy depuis phys_output_write_mod:
30  SUBROUTINE set_itau_iophy(ito)
31    IMPLICIT NONE
32    INTEGER, INTENT(IN) :: ito
33    itau_iophy = ito
34  END SUBROUTINE
35
36  SUBROUTINE init_iophy_new(rlat,rlon)
37
38    USE dimphy, ONLY: klon
39    USE mod_phys_lmdz_para, ONLY: gather, bcast, &
40                                  jj_nb, jj_begin, jj_end, ii_begin, ii_end, &
41                                  mpi_size, mpi_rank, klon_mpi, &
42                                is_sequential, is_south_pole_dyn
43  USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, klon_glo, grid_type, unstructured
44  USE print_control_mod, ONLY: prt_level,lunout
45#ifdef CPP_IOIPSL
46    USE ioipsl, ONLY: flio_dom_set
47#endif
48  use wxios, ONLY: wxios_domain_param, wxios_domain_param_unstructured, wxios_context_init, using_xios
49    IMPLICIT NONE
50    REAL,DIMENSION(klon),INTENT(IN) :: rlon
51    REAL,DIMENSION(klon),INTENT(IN) :: rlat
52
53    REAL,DIMENSION(klon_glo)        :: rlat_glo
54    REAL,DIMENSION(klon_glo)        :: rlon_glo
55   
56    INTEGER,DIMENSION(2) :: ddid
57    INTEGER,DIMENSION(2) :: dsg
58    INTEGER,DIMENSION(2) :: dsl
59    INTEGER,DIMENSION(2) :: dpf
60    INTEGER,DIMENSION(2) :: dpl
61    INTEGER,DIMENSION(2) :: dhs
62    INTEGER,DIMENSION(2) :: dhe
63    INTEGER :: i   
64    INTEGER :: data_ibegin, data_iend
65
66    if (using_xios)  CALL wxios_context_init
67   
68
69    IF (grid_type==unstructured) THEN
70   
71      CALL wxios_domain_param_unstructured("dom_glo")
72
73    ELSE
74
75      CALL gather(rlat,rlat_glo)
76      CALL bcast(rlat_glo)
77      CALL gather(rlon,rlon_glo)
78      CALL bcast(rlon_glo)
79   
80!$OMP MASTER 
81    ALLOCATE(io_lat(nbp_lat))
82    IF (klon_glo == 1) THEN
83      io_lat(1)=rlat_glo(1)
84    ELSE
85      io_lat(1)=rlat_glo(1)
86      io_lat(nbp_lat)=rlat_glo(klon_glo)
87      DO i=2,nbp_lat-1
88        io_lat(i)=rlat_glo(2+(i-2)*nbp_lon)
89      ENDDO
90    ENDIF
91
92    ALLOCATE(io_lon(nbp_lon))
93    IF (klon_glo == 1) THEN
94      io_lon(1)=rlon_glo(1)
95    ELSE
96      io_lon(1:nbp_lon)=rlon_glo(2:nbp_lon+1)
97    ENDIF
98
99!! (I) dtnb   : total number of domains
100!! (I) dnb    : domain number
101!! (I) did(:) : distributed dimensions identifiers
102!!              (up to 5 dimensions are supported)
103!! (I) dsg(:) : total number of points for each dimension
104!! (I) dsl(:) : local number of points for each dimension
105!! (I) dpf(:) : position of first local point for each dimension
106!! (I) dpl(:) : position of last local point for each dimension
107!! (I) dhs(:) : start halo size for each dimension
108!! (I) dhe(:) : end halo size for each dimension
109!! (C) cdnm   : Model domain definition name.
110!!              The names actually supported are :
111!!              "BOX", "APPLE", "ORANGE".
112!!              These names are case insensitive.
113
114    ddid=(/ 1,2 /)
115    dsg=(/ nbp_lon, nbp_lat /)
116    dsl=(/ nbp_lon, jj_nb /)
117    dpf=(/ 1,jj_begin /)
118    dpl=(/ nbp_lon, jj_end /)
119    dhs=(/ ii_begin-1,0 /)
120    IF (mpi_rank==mpi_size-1) THEN
121      dhe=(/0,0/)
122    ELSE
123      dhe=(/ nbp_lon-ii_end,0 /) 
124    ENDIF
125
126#ifndef CPP_IOIPSL_NO_OUTPUT   
127    CALL flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
128                      'APPLE',phys_domain_id)
129#endif
130    IF (using_xios) THEN 
131      ! Set values for the mask:
132      IF (mpi_rank == 0) THEN
133          data_ibegin = 0
134      ELSE
135          data_ibegin = ii_begin - 1
136      END IF
137
138      IF (mpi_rank == mpi_size-1) THEN
139          data_iend = nbp_lon
140      ELSE
141          data_iend = ii_end + 1
142      END IF
143
144      IF (prt_level>=10) THEN
145        write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," iibegin=",ii_begin , " ii_end=",ii_end," jjbegin=",jj_begin," jj_nb=",jj_nb," jj_end=",jj_end
146        write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," nbp_lon=",nbp_lon," nbp_lat=",nbp_lat
147        write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," data_ibegin=",data_ibegin," data_iend=",data_iend
148        write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," data_ibegin=",data_ibegin," data_iend=",data_iend
149        write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," is_south_pole_dyn=",is_south_pole_dyn
150      ENDIF
151
152      ! Initialize the XIOS domain coreesponding to this process:
153    ENDIF
154!$OMP END MASTER
155
156      IF (using_xios) CALL wxios_domain_param("dom_glo")
157     
158    ENDIF
159     
160  END SUBROUTINE init_iophy_new
161
162
163  SUBROUTINE init_iophy(lat,lon)
164
165    USE mod_phys_lmdz_para, ONLY: jj_begin, jj_end, ii_begin, ii_end, jj_nb, &
166                                  mpi_size, mpi_rank
167    USE ioipsl, ONLY: flio_dom_set
168    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
169
170    IMPLICIT NONE
171
172    REAL,DIMENSION(nbp_lon),INTENT(IN) :: lon
173    REAL,DIMENSION(nbp_lat),INTENT(IN) :: lat
174
175    INTEGER,DIMENSION(2) :: ddid
176    INTEGER,DIMENSION(2) :: dsg
177    INTEGER,DIMENSION(2) :: dsl
178    INTEGER,DIMENSION(2) :: dpf
179    INTEGER,DIMENSION(2) :: dpl
180    INTEGER,DIMENSION(2) :: dhs
181    INTEGER,DIMENSION(2) :: dhe
182
183!$OMP MASTER 
184    ALLOCATE(io_lat(nbp_lat))
185    io_lat(:)=lat(:)
186    ALLOCATE(io_lon(nbp_lon))
187    io_lon(:)=lon(:)
188   
189    ddid=(/ 1,2 /)
190    dsg=(/ nbp_lon, nbp_lat /)
191    dsl=(/ nbp_lon, jj_nb /)
192    dpf=(/ 1,jj_begin /)
193    dpl=(/ nbp_lon, jj_end /)
194    dhs=(/ ii_begin-1,0 /)
195    IF (mpi_rank==mpi_size-1) THEN
196      dhe=(/0,0/)
197    ELSE
198      dhe=(/ nbp_lon-ii_end,0 /) 
199    ENDIF
200   
201#ifndef CPP_IOIPSL_NO_OUTPUT
202    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
203                      'APPLE',phys_domain_id)
204#endif
205!$OMP END MASTER
206     
207  END SUBROUTINE init_iophy
208
209 SUBROUTINE histbeg_phyxios(name,itau0,zjulian,dtime,ffreq,lev,nhori,nid_day)
210!  USE dimphy
211  USE mod_phys_lmdz_para, ONLY: is_sequential, is_using_mpi, is_mpi_root, &
212                                jj_begin, jj_end, jj_nb
213  USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
214  USE ioipsl, ONLY: histbeg
215  USE wxios, ONLY: wxios_add_file, using_xios
216  IMPLICIT NONE
217  INCLUDE 'clesphys.h'
218   
219  CHARACTER*(*), INTENT(IN) :: name
220  INTEGER, INTENT(IN) :: itau0
221  REAL,INTENT(IN) :: zjulian
222  REAL,INTENT(IN) :: dtime
223  CHARACTER(LEN=*), INTENT(IN) :: ffreq
224  INTEGER,INTENT(IN) :: lev
225  INTEGER,INTENT(OUT) :: nhori
226  INTEGER,INTENT(OUT) :: nid_day
227
228!$OMP MASTER   
229  IF (is_sequential) THEN
230    CALL histbeg(name,nbp_lon,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
231                 1,nbp_lon,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
232  ELSE
233    CALL histbeg(name,nbp_lon,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
234                 1,nbp_lon,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
235  ENDIF
236
237  IF (using_xios) THEN 
238    ! ug OMP en chantier...
239    IF((.NOT. is_using_mpi) .OR. is_mpi_root) THEN
240        ! ug Création du fichier
241      IF (.not. ok_all_xml) THEN
242        CALL wxios_add_file(name, ffreq, lev)
243      ENDIF
244    ENDIF
245  ENDIF
246!$OMP END MASTER
247 
248  END SUBROUTINE histbeg_phyxios
249 
250  SUBROUTINE histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
251
252  USE mod_phys_lmdz_para, ONLY: jj_begin, jj_end, jj_nb, is_sequential
253  USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
254  USE ioipsl, ONLY: histbeg
255
256  IMPLICIT NONE
257   
258    CHARACTER*(*), INTENT(IN) :: name
259    INTEGER, INTENT(IN) :: itau0
260    REAL,INTENT(IN) :: zjulian
261    REAL,INTENT(IN) :: dtime
262    INTEGER,INTENT(OUT) :: nhori
263    INTEGER,INTENT(OUT) :: nid_day
264
265!$OMP MASTER   
266#ifndef CPP_IOIPSL_NO_OUTPUT
267    IF (is_sequential) THEN
268      CALL histbeg(name,nbp_lon,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
269                   1,nbp_lon,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
270    ELSE
271      CALL histbeg(name,nbp_lon,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
272                   1,nbp_lon,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
273    ENDIF
274#endif
275!$OMP END MASTER
276 
277  END SUBROUTINE histbeg_phy
278
279
280  SUBROUTINE histbeg_phy_points(rlon,rlat,pim,tabij,ipt,jpt, &
281             plon,plat,plon_bounds,plat_bounds, &
282             nname,itau0,zjulian,dtime,nnhori,nnid_day)
283  USE dimphy, ONLY: klon
284  USE mod_phys_lmdz_para, ONLY: gather, bcast, &
285                                is_sequential, klon_mpi_begin, klon_mpi_end, &
286                                mpi_rank
287  USE mod_grid_phy_lmdz, ONLY: klon_glo, nbp_lon, nbp_lat, grid1dTo2d_glo
288  USE ioipsl, ONLY: histbeg
289
290  IMPLICIT NONE
291
292    REAL,DIMENSION(klon),INTENT(IN) :: rlon
293    REAL,DIMENSION(klon),INTENT(IN) :: rlat
294    INTEGER, INTENT(IN) :: itau0
295    REAL,INTENT(IN) :: zjulian
296    REAL,INTENT(IN) :: dtime
297    INTEGER, INTENT(IN) :: pim
298    INTEGER, intent(out) :: nnhori
299    CHARACTER(len=20), INTENT(IN) :: nname
300    INTEGER, INTENT(OUT) :: nnid_day
301    INTEGER :: i
302    REAL,DIMENSION(klon_glo)        :: rlat_glo
303    REAL,DIMENSION(klon_glo)        :: rlon_glo
304    INTEGER, DIMENSION(pim), INTENT(IN)  :: tabij
305    REAL,DIMENSION(pim), INTENT(IN) :: plat, plon
306    INTEGER,DIMENSION(pim), INTENT(IN) :: ipt, jpt
307    REAL,DIMENSION(pim,2), intent(out) :: plat_bounds, plon_bounds
308
309    INTEGER, SAVE :: tabprocbeg, tabprocend
310!$OMP THREADPRIVATE(tabprocbeg, tabprocend)
311    INTEGER :: ip
312    INTEGER, PARAMETER :: nip=1
313    INTEGER :: npproc
314    REAL, allocatable, DIMENSION(:) :: npplat, npplon
315    REAL, allocatable, DIMENSION(:,:) :: npplat_bounds, npplon_bounds
316    REAL, DIMENSION(nbp_lon,nbp_lat) :: zx_lon, zx_lat
317
318    CALL gather(rlat,rlat_glo)
319    CALL bcast(rlat_glo)
320    CALL gather(rlon,rlon_glo)
321    CALL bcast(rlon_glo)
322
323!$OMP MASTER
324    DO i=1,pim
325
326!    print*,'CFMIP_iophy i tabij lon lat',i,tabij(i),plon(i),plat(i)
327
328     plon_bounds(i,1)=rlon_glo(tabij(i)-1)
329     plon_bounds(i,2)=rlon_glo(tabij(i)+1)
330     IF (plon_bounds(i,2).LE.0..AND.plon_bounds(i,1).GE.0.) THEN
331      IF (rlon_glo(tabij(i)).GE.0.) THEN
332       plon_bounds(i,2)=-1*plon_bounds(i,2)
333      ENDIF
334     ENDIF
335     IF (plon_bounds(i,2).GE.0..AND.plon_bounds(i,1).LE.0.) THEN
336      IF (rlon_glo(tabij(i)).LE.0.) THEN
337       plon_bounds(i,2)=-1*plon_bounds(i,2)
338      ENDIF
339     ENDIF
340!
341     IF ( tabij(i).LE.nbp_lon) THEN
342      plat_bounds(i,1)=rlat_glo(tabij(i))
343     ELSE
344      plat_bounds(i,1)=rlat_glo(tabij(i)-nbp_lon)
345     ENDIF
346     plat_bounds(i,2)=rlat_glo(tabij(i)+nbp_lon)
347!
348!    print*,'CFMIP_iophy point i lon lon_bds',i,plon_bounds(i,1),rlon_glo(tabij(i)),plon_bounds(i,2)
349!    print*,'CFMIP_iophy point i lat lat_bds',i,plat_bounds(i,1),rlat_glo(tabij(i)),plat_bounds(i,2)
350!
351    ENDDO
352    if (is_sequential) then
353
354     npstn=pim
355     IF(.NOT. ALLOCATED(nptabij)) THEN
356      ALLOCATE(nptabij(pim))
357     ENDIF
358     DO i=1,pim
359      nptabij(i)=tabij(i)
360     ENDDO
361
362       CALL grid1dTo2d_glo(rlon_glo,zx_lon)
363       IF ((nbp_lon*nbp_lat).GT.1) THEN
364       DO i = 1, nbp_lon
365         zx_lon(i,1) = rlon_glo(i+1)
366         zx_lon(i,nbp_lat) = rlon_glo(i+1)
367       ENDDO
368       ENDIF
369       CALL grid1dTo2d_glo(rlat_glo,zx_lat)
370
371    DO i=1,pim
372!    print*,'CFMIP_iophy i tabij lon lat',i,tabij(i),plon(i),plat(i)
373
374     plon_bounds(i,1)=zx_lon(ipt(i)-1,jpt(i))
375     plon_bounds(i,2)=zx_lon(ipt(i)+1,jpt(i))
376
377     IF (ipt(i).EQ.1) THEN
378      plon_bounds(i,1)=zx_lon(nbp_lon,jpt(i))
379      plon_bounds(i,2)=360.+zx_lon(ipt(i)+1,jpt(i))
380     ENDIF
381 
382     IF (ipt(i).EQ.nbp_lon) THEN
383      plon_bounds(i,2)=360.+zx_lon(1,jpt(i))
384     ENDIF
385
386     plat_bounds(i,1)=zx_lat(ipt(i),jpt(i)-1)
387     plat_bounds(i,2)=zx_lat(ipt(i),jpt(i)+1)
388
389     IF (jpt(i).EQ.1) THEN
390      plat_bounds(i,1)=zx_lat(ipt(i),1)+0.001
391      plat_bounds(i,2)=zx_lat(ipt(i),1)-0.001
392     ENDIF
393 
394     IF (jpt(i).EQ.nbp_lat) THEN
395      plat_bounds(i,1)=zx_lat(ipt(i),nbp_lat)+0.001
396      plat_bounds(i,2)=zx_lat(ipt(i),nbp_lat)-0.001
397     ENDIF
398!
399!    print*,'CFMIP_iophy point i lon lon_bds',i,plon_bounds(i,1),rlon(tabij(i)),plon_bounds(i,2)
400!    print*,'CFMIP_iophy point i lat lat_bds',i,plat_bounds(i,1),rlat(tabij(i)),plat_bounds(i,2)
401!
402    ENDDO
403
404#ifndef CPP_IOIPSL_NO_OUTPUT
405     CALL histbeg(nname,pim,plon,plon_bounds, &
406                           plat,plat_bounds, &
407                           itau0, zjulian, dtime, nnhori, nnid_day)
408#endif
409    ELSE
410     npproc=0
411     DO ip=1, pim
412      tabprocbeg=klon_mpi_begin
413      tabprocend=klon_mpi_end
414      IF(tabij(ip).GE.tabprocbeg.AND.tabij(ip).LE.tabprocend) THEN
415       npproc=npproc+1
416       npstn=npproc
417      ENDIF
418     ENDDO
419!    print*,'CFMIP_iophy mpi_rank npstn',mpi_rank,npstn
420     IF(.NOT. ALLOCATED(nptabij)) THEN
421      ALLOCATE(nptabij(npstn))
422      ALLOCATE(npplon(npstn), npplat(npstn))
423      ALLOCATE(npplon_bounds(npstn,2), npplat_bounds(npstn,2))
424     ENDIF
425     npproc=0
426     DO ip=1, pim
427      IF(tabij(ip).GE.tabprocbeg.AND.tabij(ip).LE.tabprocend) THEN
428       npproc=npproc+1
429       nptabij(npproc)=tabij(ip)
430!      print*,'mpi_rank npproc ip plon plat tabij=',mpi_rank,npproc,ip, &
431!      plon(ip),plat(ip),tabij(ip)
432       npplon(npproc)=plon(ip)
433       npplat(npproc)=plat(ip)
434       npplon_bounds(npproc,1)=plon_bounds(ip,1)
435       npplon_bounds(npproc,2)=plon_bounds(ip,2)
436       npplat_bounds(npproc,1)=plat_bounds(ip,1)
437       npplat_bounds(npproc,2)=plat_bounds(ip,2)
438!!!
439!!! print qui sert a reordonner les points stations selon l'ordre CFMIP
440!!! ne pas enlever
441        print*,'iophy_mpi rank ip lon lat',mpi_rank,ip,plon(ip),plat(ip)
442!!!
443      ENDIF
444     ENDDO
445#ifndef CPP_IOIPSL_NO_OUTPUT
446     CALL histbeg(nname,npstn,npplon,npplon_bounds, &
447                            npplat,npplat_bounds, &
448                            itau0,zjulian,dtime,nnhori,nnid_day,phys_domain_id)
449#endif
450    ENDIF
451!$OMP END MASTER
452
453  END SUBROUTINE histbeg_phy_points
454
455
456  SUBROUTINE histdef2d_old (iff,lpoint,flag_var,nomvar,titrevar,unitvar)
457
458    USE ioipsl, ONLY: histdef
459    USE mod_phys_lmdz_para, ONLY: jj_nb, is_master
460    USE phys_output_var_mod, ONLY: type_ecri, zoutm, zdtime_moy, lev_files, &
461                                   nid_files, nhorim, swaero_diag, dryaod_diag, nfiles, &
462                                   ok_4xCO2atm
463    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
464    USE aero_mod, ONLY : naero_tot, name_aero_tau
465    USE print_control_mod, ONLY: prt_level,lunout
466
467    IMPLICIT NONE
468
469    INCLUDE "clesphys.h"
470
471    INTEGER                          :: iff
472    INTEGER                          :: naero
473    LOGICAL                          :: lpoint
474    INTEGER, DIMENSION(nfiles)       :: flag_var
475    CHARACTER(LEN=20)                :: nomvar
476    CHARACTER(LEN=*)                 :: titrevar
477    CHARACTER(LEN=*)                 :: unitvar
478
479    REAL zstophym
480
481    IF (type_ecri(iff)=='inst(X)'.OR.type_ecri(iff)=='once') THEN
482       zstophym=zoutm(iff)
483    ELSE
484       zstophym=zdtime_moy
485    ENDIF
486    IF (check_dim .AND. is_master) WRITE(lunout,*)'histdef2d_old for ', nomvar
487    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
488    CALL conf_physoutputs(nomvar,flag_var)
489
490    IF(.NOT.lpoint) THEN 
491       IF ( flag_var(iff)<=lev_files(iff) ) THEN
492          CALL histdef (nid_files(iff),nomvar,titrevar,unitvar, &
493               nbp_lon,jj_nb,nhorim(iff), 1,1,1, -99, 32, &
494               type_ecri(iff), zstophym,zoutm(iff))               
495       ENDIF
496    ELSE
497       IF ( flag_var(iff)<=lev_files(iff) ) THEN
498          CALL histdef (nid_files(iff),nomvar,titrevar,unitvar, &
499               npstn,1,nhorim(iff), 1,1,1, -99, 32, &
500               type_ecri(iff), zstophym,zoutm(iff))               
501       ENDIF
502    ENDIF
503
504    ! Set swaero_diag=true if at least one of the concerned variables are defined
505    IF (nomvar=='topswad' .OR. nomvar=='topswad0' .OR. nomvar=='solswad' .OR. nomvar=='solswad0' .OR. &
506        nomvar=='toplwad' .OR. nomvar=='toplwad0' .OR. nomvar=='sollwad' .OR. nomvar=='sollwad0' .OR. &
507        nomvar=='topswai' .OR. nomvar=='solswai' ) THEN
508       IF  ( flag_var(iff)<=lev_files(iff) ) swaero_diag=.TRUE.
509    ENDIF
510
511    ! Set dryaod_diag=true if at least one of the concerned variables are defined
512    IF (nomvar=='dryod550aer') THEN
513      IF  ( flag_var(iff)<=lev_files(iff) ) dryaod_diag=.TRUE.
514    ENDIF
515    DO naero = 1, naero_tot-1
516      IF (nomvar=='dryod550_'//name_aero_tau(naero)) THEN
517        IF  ( flag_var(iff)<=lev_files(iff) ) dryaod_diag=.TRUE.
518      ENDIF
519    ENDDO
520
521    ! Set ok_4xCO2atm=true if at least one of the concerned variables are
522    ! defined
523    IF (nomvar=='rsut4co2'.OR.nomvar=='rlut4co2'.OR.nomvar=='rsutcs4co2' &
524        .OR. nomvar=='rlutcs4co2'.OR.nomvar=='rsu4co2'.OR.nomvar=='rsucs4co2' &
525        .OR.nomvar=='rsu4co2'.OR.nomvar=='rsucs4co2'.OR.nomvar=='rsd4co2'.OR. &
526        nomvar=='rsdcs4co2'.OR.nomvar=='rlu4co2'.OR.nomvar=='rlucs4co2'.OR.&
527        nomvar=='rld4co2'.OR.nomvar=='rldcs4co2') THEN
528        IF ( flag_var(iff)<=lev_files(iff) ) ok_4xCO2atm=.TRUE.
529    ENDIF
530  END SUBROUTINE histdef2d_old
531
532  SUBROUTINE histdef3d_old (iff,lpoint,flag_var,nomvar,titrevar,unitvar)
533
534    USE ioipsl, ONLY: histdef
535    USE dimphy, ONLY: klev
536    USE mod_phys_lmdz_para, ONLY: jj_nb, is_master
537    USE phys_output_var_mod, ONLY: type_ecri, zoutm, lev_files, nid_files, &
538                                   nhorim, zdtime_moy, levmin, levmax, &
539                                   nvertm, nfiles
540    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
541    USE print_control_mod, ONLY: prt_level,lunout
542    IMPLICIT NONE
543
544    INCLUDE "clesphys.h"
545
546    INTEGER                          :: iff
547    LOGICAL                          :: lpoint
548    INTEGER, DIMENSION(nfiles)       :: flag_var
549    CHARACTER(LEN=20)                 :: nomvar
550    CHARACTER(LEN=*)                 :: titrevar
551    CHARACTER(LEN=*)                 :: unitvar
552
553    REAL zstophym
554
555    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
556    CALL conf_physoutputs(nomvar,flag_var)
557
558    IF (check_dim .AND. is_master) WRITE(lunout,*)'histdef3d_old for ', nomvar
559
560    IF (type_ecri(iff)=='inst(X)'.OR.type_ecri(iff)=='once') THEN
561       zstophym=zoutm(iff)
562    ELSE
563       zstophym=zdtime_moy
564    ENDIF
565
566    IF(.NOT.lpoint) THEN
567       IF ( flag_var(iff)<=lev_files(iff) ) THEN
568          CALL histdef (nid_files(iff), nomvar, titrevar, unitvar, &
569               nbp_lon, jj_nb, nhorim(iff), klev, levmin(iff), &
570               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, type_ecri(iff), &
571               zstophym, zoutm(iff))
572       ENDIF
573    ELSE
574       IF ( flag_var(iff)<=lev_files(iff) ) THEN
575          CALL histdef (nid_files(iff), nomvar, titrevar, unitvar, &
576               npstn,1,nhorim(iff), klev, levmin(iff), &
577               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, &
578               type_ecri(iff), zstophym,zoutm(iff))
579       ENDIF
580    ENDIF
581  END SUBROUTINE histdef3d_old
582
583  SUBROUTINE histdef2d (iff,var)
584
585    USE ioipsl, ONLY: histdef
586    USE mod_phys_lmdz_para, ONLY: jj_nb, is_master
587    USE phys_output_var_mod, ONLY: ctrl_out, type_ecri_files, zoutm, zdtime_moy, &
588                                   clef_stations, phys_out_filenames, lev_files, &
589                                   nid_files, nhorim, swaerofree_diag, swaero_diag, dryaod_diag,&
590                                   ok_4xCO2atm
591    USE print_control_mod, ONLY: prt_level,lunout
592    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
593    USE aero_mod, ONLY : naero_tot, name_aero_tau
594    USE wxios, ONLY: wxios_add_field_to_file, using_xios
595    USE print_control_mod, ONLY: prt_level,lunout
596    IMPLICIT NONE
597
598    INCLUDE "clesphys.h"
599
600    INTEGER                          :: iff
601    INTEGER                          :: naero
602    TYPE(ctrl_out)                   :: var
603
604    REAL zstophym
605    CHARACTER(LEN=20) :: typeecrit
606
607    IF (check_dim .AND. is_master) WRITE(lunout,*)'histdef2d for ', var%name
608
609    ! ug On récupère le type écrit de la structure:
610    !       Assez moche, à refaire si meilleure méthode...
611    IF (INDEX(var%type_ecrit(iff), "once") > 0) THEN
612       typeecrit = 'once'
613    ELSE IF(INDEX(var%type_ecrit(iff), "t_min") > 0) THEN
614       typeecrit = 't_min(X)'
615    ELSE IF(INDEX(var%type_ecrit(iff), "t_max") > 0) THEN
616       typeecrit = 't_max(X)'
617    ELSE IF(INDEX(var%type_ecrit(iff), "inst") > 0) THEN
618       typeecrit = 'inst(X)'
619    ELSE
620       typeecrit = type_ecri_files(iff)
621    ENDIF
622
623    IF (typeecrit=='inst(X)'.OR.typeecrit=='once') THEN
624       zstophym=zoutm(iff)
625    ELSE
626       zstophym=zdtime_moy
627    ENDIF
628
629    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
630    CALL conf_physoutputs(var%name, var%flag)
631
632    IF(.NOT.clef_stations(iff)) THEN 
633
634      IF (using_xios) THEN 
635        IF (.not. ok_all_xml) THEN
636          IF ( var%flag(iff)<=lev_files(iff) ) THEN
637            CALL wxios_add_field_to_file(var%name, 2, iff, phys_out_filenames(iff), &
638            var%description, var%unit, var%flag(iff), typeecrit)
639            IF (prt_level >= 10) THEN
640              WRITE(lunout,*) 'histdef2d: call wxios_add_field_to_file var%name iff: ', &
641                              trim(var%name),iff
642            ENDIF
643          ENDIF
644        ENDIF
645      ENDIF
646#ifndef CPP_IOIPSL_NO_OUTPUT
647
648       IF ( var%flag(iff)<=lev_files(iff) ) THEN
649          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
650               nbp_lon,jj_nb,nhorim(iff), 1,1,1, -99, 32, &
651               typeecrit, zstophym,zoutm(iff))               
652       ENDIF
653    ELSE
654       IF ( var%flag(iff)<=lev_files(iff)) THEN
655          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
656               npstn,1,nhorim(iff), 1,1,1, -99, 32, &
657               typeecrit, zstophym,zoutm(iff))               
658       ENDIF
659#endif
660    ENDIF
661
662    ! Set swaero_diag=true if at least one of the concerned variables are defined
663    !--OB 30/05/2016 use wider set of variables
664    IF ( var%name=='topswad' .OR. var%name=='topswad0' .OR. var%name=='solswad' .OR. var%name=='solswad0' .OR. &
665         var%name=='topswai' .OR. var%name=='solswai'  .OR. ( iflag_rrtm==1 .AND. (                            &
666         var%name=='toplwad' .OR. var%name=='toplwad0' .OR. var%name=='sollwad' .OR. var%name=='sollwad0' .OR. &
667         var%name=='toplwai' .OR. var%name=='sollwai'  ) ) ) THEN
668       IF  ( var%flag(iff)<=lev_files(iff) ) swaero_diag=.TRUE.
669    ENDIF
670
671    ! Set swaerofree_diag=true if at least one of the concerned variables are defined
672    IF (var%name=='SWupTOAcleanclr' .OR. var%name=='SWupSFCcleanclr' .OR. var%name=='SWdnSFCcleanclr' .OR. &
673        var%name=='LWupTOAcleanclr' .OR. var%name=='LWdnSFCcleanclr' ) THEN
674       IF  ( var%flag(iff)<=lev_files(iff) ) swaerofree_diag=.TRUE.
675    ENDIF
676
677    ! set dryaod_dry=true if at least one of the concerned variables are defined
678    IF (var%name=='dryod550aer') THEN
679      IF  ( var%flag(iff)<=lev_files(iff) ) dryaod_diag=.TRUE.
680    ENDIF
681    !
682    DO naero = 1, naero_tot-1
683      IF (var%name=='dryod550_'//name_aero_tau(naero)) THEN
684        IF  ( var%flag(iff)<=lev_files(iff) ) dryaod_diag=.TRUE.
685      ENDIF
686    ENDDO
687    ! Set ok_4xCO2atm=true if at least one of the concerned variables are
688    ! defined
689    IF (var%name=='rsut4co2'.OR.var%name=='rlut4co2'.OR.var%name=='rsutcs4co2' &
690        .OR. var%name=='rlutcs4co2'.OR.var%name=='rsu4co2'.OR.var%name=='rsucs4co2' &
691        .OR.var%name=='rsu4co2'.OR.var%name=='rsucs4co2'.OR.var%name=='rsd4co2'.OR. &
692        var%name=='rsdcs4co2'.OR.var%name=='rlu4co2'.OR.var%name=='rlucs4co2'.OR.&
693        var%name=='rld4co2'.OR.var%name=='rldcs4co2') THEN
694        IF ( var%flag(iff)<=lev_files(iff) ) ok_4xCO2atm=.TRUE.
695    ENDIF
696  END SUBROUTINE histdef2d
697
698  SUBROUTINE histdef3d (iff,var)
699
700    USE ioipsl, ONLY: histdef
701    USE dimphy, ONLY: klev
702    USE mod_phys_lmdz_para, ONLY: jj_nb, is_master
703    USE phys_output_var_mod, ONLY: ctrl_out, type_ecri_files, zoutm, zdtime_moy, &
704                                   clef_stations, phys_out_filenames, lev_files, &
705                                   nid_files, nhorim, swaerofree_diag, levmin, &
706                                   levmax, nvertm
707    USE print_control_mod, ONLY: prt_level,lunout
708    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
709    USE wxios, ONLY: wxios_add_field_to_file, using_xios
710    USE print_control_mod, ONLY: prt_level,lunout
711    IMPLICIT NONE
712
713    INCLUDE "clesphys.h"
714
715    INTEGER                          :: iff
716    TYPE(ctrl_out)                   :: var
717
718    REAL zstophym
719    CHARACTER(LEN=20) :: typeecrit
720
721    IF (check_dim .AND. is_master) WRITE(lunout,*)'histdef3d for ', var%name
722
723    ! ug On récupère le type écrit de la structure:
724    !       Assez moche, à refaire si meilleure méthode...
725    IF (INDEX(var%type_ecrit(iff), "once") > 0) THEN
726       typeecrit = 'once'
727    ELSE IF(INDEX(var%type_ecrit(iff), "t_min") > 0) THEN
728       typeecrit = 't_min(X)'
729    ELSE IF(INDEX(var%type_ecrit(iff), "t_max") > 0) THEN
730       typeecrit = 't_max(X)'
731    ELSE IF(INDEX(var%type_ecrit(iff), "inst") > 0) THEN
732       typeecrit = 'inst(X)'
733    ELSE
734       typeecrit = type_ecri_files(iff)
735    ENDIF
736
737    ! Appel a la lecture des noms et niveau d'ecriture des variables dans output.def
738    CALL conf_physoutputs(var%name,var%flag)
739
740    IF (typeecrit=='inst(X)'.OR.typeecrit=='once') THEN
741       zstophym=zoutm(iff)
742    ELSE
743       zstophym=zdtime_moy
744    ENDIF
745
746    IF(.NOT.clef_stations(iff)) THEN
747
748      IF (using_xios) THEN 
749
750        IF (.not. ok_all_xml) THEN
751          IF ( var%flag(iff)<=lev_files(iff) ) THEN
752          CALL wxios_add_field_to_file(var%name, 3, iff, phys_out_filenames(iff), &
753          var%description, var%unit, var%flag(iff), typeecrit)
754            IF (prt_level >= 10) THEN
755              WRITE(lunout,*) 'histdef3d: call wxios_add_field_to_file var%name iff: ', &
756                              trim(var%name),iff
757            ENDIF
758          ENDIF
759        ENDIF
760      ENDIF
761#ifndef CPP_IOIPSL_NO_OUTPUT
762
763       IF ( var%flag(iff)<=lev_files(iff) ) THEN
764          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
765               nbp_lon, jj_nb, nhorim(iff), klev, levmin(iff), &
766               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, typeecrit, &
767               zstophym, zoutm(iff))
768       ENDIF
769    ELSE
770       IF ( var%flag(iff)<=lev_files(iff)) THEN
771          CALL histdef (nid_files(iff), var%name, var%description, var%unit, &
772               npstn,1,nhorim(iff), klev, levmin(iff), &
773               levmax(iff)-levmin(iff)+1, nvertm(iff), 32, &
774               typeecrit, zstophym,zoutm(iff))
775       ENDIF
776#endif
777    ENDIF
778
779    ! Set swaerofree_diag=true if at least one of the concerned variables are defined
780    IF (var%name=='rsucsaf' .OR. var%name=='rsdcsaf') THEN
781       IF  ( var%flag(iff)<=lev_files(iff) ) swaerofree_diag=.TRUE.
782    ENDIF
783
784  END SUBROUTINE histdef3d
785
786  SUBROUTINE conf_physoutputs(nam_var,flag_var)
787!!! Lecture des noms et niveau de sortie des variables dans output.def
788    !   en utilisant les routines getin de IOIPSL 
789    USE ioipsl, ONLY: getin
790    USE phys_output_var_mod, ONLY: nfiles
791    USE print_control_mod, ONLY: prt_level,lunout
792    IMPLICIT NONE
793
794    CHARACTER(LEN=*), INTENT(INOUT) :: nam_var
795    INTEGER,          INTENT(INOUT) :: flag_var(nfiles)
796
797    IF(prt_level>10) WRITE(lunout,*)'Avant getin: nam_var flag_var ',nam_var,flag_var(:)
798    CALL getin('flag_'//nam_var,flag_var)
799    CALL getin('name_'//nam_var,nam_var)
800    IF(prt_level>10) WRITE(lunout,*)'Apres getin: nam_var flag_var ',nam_var,flag_var(:)
801
802  END SUBROUTINE conf_physoutputs
803
804 
805  SUBROUTINE histwrite2d_phy_old(nid,lpoint,name,itau,field)
806
807    USE dimphy, ONLY: klon
808    USE mod_phys_lmdz_para, ONLY: Gather_omp, grid1Dto2D_mpi, &
809                                  is_sequential, klon_mpi_begin, klon_mpi_end, &
810                                  jj_nb, klon_mpi, is_master
811    USE ioipsl, ONLY: histwrite
812    USE print_control_mod, ONLY: prt_level,lunout
813    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
814
815    IMPLICIT NONE
816   
817    INTEGER,INTENT(IN) :: nid
818    LOGICAL,INTENT(IN) :: lpoint
819    CHARACTER*(*), INTENT(IN) :: name
820    INTEGER, INTENT(IN) :: itau
821    REAL,DIMENSION(:),INTENT(IN) :: field
822    REAL,DIMENSION(klon_mpi) :: buffer_omp
823    INTEGER, allocatable, DIMENSION(:) :: index2d
824    REAL :: Field2d(nbp_lon,jj_nb)
825
826    INTEGER :: ip
827    REAL,ALLOCATABLE,DIMENSION(:) :: fieldok
828
829    IF (size(field)/=klon) CALL abort_physic('iophy::histwrite2d','Field first DIMENSION not equal to klon',1)
830    IF (check_dim .AND. is_master) WRITE(lunout,*)'histwrite2d_phy_old for ', name
831
832    CALL Gather_omp(field,buffer_omp)   
833!$OMP MASTER
834    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
835    IF (.NOT.lpoint) THEN
836     ALLOCATE(index2d(nbp_lon*jj_nb))
837     ALLOCATE(fieldok(nbp_lon*jj_nb))
838     IF (prt_level >= 10) write(lunout,*)'Sending ',name,' to IOIPSL'
839     CALL histwrite(nid,name,itau,Field2d,nbp_lon*jj_nb,index2d)
840     IF (prt_level >= 10) write(lunout,*)'Finished sending ',name,' to IOIPSL'
841    ELSE
842     ALLOCATE(fieldok(npstn))
843     ALLOCATE(index2d(npstn))
844
845     IF (is_sequential) THEN
846!     klon_mpi_begin=1
847!     klon_mpi_end=klon
848      DO ip=1, npstn
849       fieldok(ip)=buffer_omp(nptabij(ip))
850      ENDDO
851     ELSE
852      DO ip=1, npstn
853!     print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
854       IF(nptabij(ip).GE.klon_mpi_begin.AND. &
855          nptabij(ip).LE.klon_mpi_end) THEN
856         fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
857       ENDIF
858      ENDDO
859     ENDIF
860     IF (prt_level >= 10) write(lunout,*)'Sending ',name,' to IOIPSL'
861     CALL histwrite(nid,name,itau,fieldok,npstn,index2d)
862     IF (prt_level >= 10) write(lunout,*)'Finished sending ',name,' to IOIPSL'
863!
864    ENDIF
865    DEALLOCATE(index2d)
866    DEALLOCATE(fieldok)
867!$OMP END MASTER   
868
869 
870  END SUBROUTINE histwrite2d_phy_old
871
872  SUBROUTINE histwrite3d_phy_old(nid,lpoint,name,itau,field)
873
874    USE dimphy, ONLY: klon
875    USE mod_phys_lmdz_para, ONLY: Gather_omp, grid1Dto2D_mpi, &
876                                  is_sequential, klon_mpi_begin, klon_mpi_end, &
877                                  jj_nb, klon_mpi, is_master
878    USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat
879    USE ioipsl, ONLY: histwrite
880    USE print_control_mod, ONLY: prt_level,lunout
881
882    IMPLICIT NONE
883   
884    INTEGER,INTENT(IN) :: nid
885    LOGICAL,INTENT(IN) :: lpoint
886    CHARACTER*(*), INTENT(IN) :: name
887    INTEGER, INTENT(IN) :: itau
888    REAL,DIMENSION(:,:),INTENT(IN) :: field  ! --> field(klon,:)
889    REAL,DIMENSION(klon_mpi,size(field,2)) :: buffer_omp
890    REAL :: Field3d(nbp_lon,jj_nb,size(field,2))
891    INTEGER :: ip, n, nlev
892    INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
893    REAL,allocatable, DIMENSION(:,:) :: fieldok
894
895    IF (check_dim .AND. is_master) WRITE(lunout,*)'histwrite3d_phy_old for ', name
896
897    IF (size(field,1)/=klon) CALL abort_physic('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
898    nlev=size(field,2)
899
900    CALL Gather_omp(field,buffer_omp)
901!$OMP MASTER
902    CALL grid1Dto2D_mpi(buffer_omp,field3d)
903    IF (.NOT.lpoint) THEN
904      ALLOCATE(index3d(nbp_lon*jj_nb*nlev))
905      ALLOCATE(fieldok(nbp_lon*jj_nb,nlev))
906      IF (prt_level >= 10) write(lunout,*)'Sending ',name,' to IOIPSL'
907      CALL histwrite(nid,name,itau,Field3d,nbp_lon*jj_nb*nlev,index3d)
908      IF (prt_level >= 10) write(lunout,*)'Finished sending ',name,' to IOIPSL'
909    ELSE
910      nlev=size(field,2)
911      ALLOCATE(index3d(npstn*nlev))
912      ALLOCATE(fieldok(npstn,nlev))
913
914      IF (is_sequential) THEN
915!       klon_mpi_begin=1
916!       klon_mpi_end=klon
917        DO n=1, nlev
918        DO ip=1, npstn
919          fieldok(ip,n)=buffer_omp(nptabij(ip),n)
920        ENDDO
921        ENDDO
922      ELSE
923        DO n=1, nlev
924        DO ip=1, npstn
925          IF(nptabij(ip).GE.klon_mpi_begin.AND. &
926           nptabij(ip).LE.klon_mpi_end) THEN
927           fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
928          ENDIF
929        ENDDO
930        ENDDO
931      ENDIF
932      IF (prt_level >= 10) write(lunout,*)'Sending ',name,' to IOIPSL'
933      CALL histwrite(nid,name,itau,fieldok,npstn*nlev,index3d)
934      IF (prt_level >= 10) write(lunout,*)'Finished sending ',name,' to IOIPSL'
935    ENDIF
936    DEALLOCATE(index3d)
937    DEALLOCATE(fieldok)
938!$OMP END MASTER   
939
940  END SUBROUTINE histwrite3d_phy_old
941
942
943! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE
944  SUBROUTINE histwrite2d_phy(var,field, STD_iff)
945
946  USE mod_phys_lmdz_omp_transfert, ONLY: bcast_omp
947  USE dimphy, ONLY: klon, klev
948  USE mod_phys_lmdz_para, ONLY: gather_omp, grid1dto2d_mpi, &
949                                jj_nb, klon_mpi, klon_mpi_begin, &
950                                klon_mpi_end, is_sequential, is_master
951  USE ioipsl, ONLY: histwrite
952  USE phys_output_var_mod, ONLY: ctrl_out, clef_files, lev_files, &
953                                 nfiles, vars_defined, clef_stations, &
954                                 nid_files, swaerofree_diag, swaero_diag, dryaod_diag, ok_4xCO2atm
955  USE print_control_mod, ONLY: prt_level,lunout
956  USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat, grid_type, unstructured, regular_lonlat
957  USE lmdz_xios, ONLY: xios_send_field, xios_field_is_active, using_xios
958  USE print_control_mod, ONLY: lunout, prt_level
959
960  IMPLICIT NONE
961  INCLUDE 'clesphys.h'
962
963  TYPE(ctrl_out), INTENT(IN) :: var
964  REAL, DIMENSION(:), INTENT(IN) :: field
965  INTEGER, INTENT(IN), OPTIONAL :: STD_iff ! ug RUSTINE POUR LES STD LEVS.....
966     
967  INTEGER :: iff, iff_beg, iff_end
968  LOGICAL, SAVE  :: firstx
969!$OMP THREADPRIVATE(firstx)
970
971  REAL,DIMENSION(klon_mpi) :: buffer_omp
972  INTEGER, allocatable, DIMENSION(:) :: index2d
973  REAL :: Field2d(nbp_lon,jj_nb)
974
975  INTEGER :: ip
976  REAL, ALLOCATABLE, DIMENSION(:) :: fieldok
977  logical, save :: is_active = .true.
978
979  IF (check_dim .AND. is_master) WRITE(lunout,*)'histwrite2d_phy for ',trim(var%name)
980
981  IF (prt_level >= 10) THEN
982    WRITE(lunout,*)'Begin histwrite2d_phy for ',trim(var%name)
983  ENDIF
984
985! ug RUSTINE POUR LES STD LEVS.....
986  IF (PRESENT(STD_iff)) THEN
987        iff_beg = STD_iff
988        iff_end = STD_iff
989  ELSE
990        iff_beg = 1
991        iff_end = nfiles
992  ENDIF
993
994  ! On regarde si on est dans la phase de définition ou d'écriture:
995  IF (.NOT.vars_defined) THEN
996!$OMP MASTER
997      !Si phase de définition.... on définit
998      IF (.not. ok_all_xml) THEN
999      IF (prt_level >= 10) THEN
1000      write(lunout,*)"histwrite2d_phy: .not.vars_defined ; time to define ", &
1001                     trim(var%name)
1002      ENDIF
1003      DO iff=iff_beg, iff_end
1004         IF (clef_files(iff)) THEN
1005            CALL histdef2d(iff, var)
1006         ENDIF
1007      ENDDO
1008      ENDIF
1009!$OMP END MASTER
1010!--broadcasting the flags that have been changed in histdef2d on OMP masters
1011      CALL bcast_omp(swaero_diag)
1012      CALL bcast_omp(swaerofree_diag)
1013      CALL bcast_omp(dryaod_diag)
1014      CALL bcast_omp(ok_4xCO2atm)
1015
1016  ELSE
1017    IF (using_xios) THEN
1018      IF (ok_all_xml) THEN
1019        !$omp barrier
1020        !$omp master
1021        is_active = xios_field_is_active(var%name, at_current_timestep_arg=.false.)
1022        !$omp end master
1023        !$omp barrier
1024        IF(.not. is_active) RETURN
1025      ENDIF
1026    ENDIF
1027
1028    !Et sinon on.... écrit
1029    IF (SIZE(field)/=klon .AND. SIZE(field)/=klev .AND. SIZE(field)/=klev+1) &
1030         CALL abort_physic('iophy::histwrite2d_phy',&
1031         'Field first DIMENSION not equal to klon/klev',1)
1032    IF (prt_level >= 10) THEn
1033      WRITE (lunout,*)"histwrite2d_phy: .not.vars_defined ; time to gather and write ", trim(var%name)
1034    ENDIF
1035   
1036   
1037    IF (SIZE(field) == klon) then
1038        CALL Gather_omp(field,buffer_omp)
1039    ELSE
1040        buffer_omp(:)=0.
1041    ENDIF
1042!$OMP MASTER
1043    IF (grid_type==regular_lonlat) CALL grid1Dto2D_mpi(buffer_omp,Field2d)
1044
1045! La boucle sur les fichiers:
1046      firstx=.true.
1047
1048      IF (ok_all_xml) THEN
1049        IF (using_xios) THEN 
1050          IF (prt_level >= 10) THEN
1051             write(lunout,*)'Dans iophy histwrite2D,var%name ', trim(var%name)                       
1052          ENDIF
1053         
1054          IF (grid_type==regular_lonlat) THEN
1055            IF (SIZE(field) == klon) then
1056              CALL xios_send_field(var%name, Field2d)
1057            ELSE
1058               CALL xios_send_field(var%name, field)
1059            ENDIF
1060          ELSE IF (grid_type==unstructured) THEN
1061            IF (SIZE(field) == klon) then
1062              CALL xios_send_field(var%name, buffer_omp)
1063            ELSE
1064               CALL xios_send_field(var%name, field)
1065            ENDIF
1066
1067          ENDIF
1068          IF (prt_level >= 10) THEN
1069             write(lunout,*)'Dans iophy histwrite2D,var%name apres xios_send ',&
1070                             trim(var%name)                       
1071          ENDIF
1072        ELSE
1073          CALL abort_physic ('iophy','cannot have ok_all_xml = .T. without CPP_XIOS defined' ,1)
1074        ENDIF
1075      ELSE 
1076        DO iff=iff_beg, iff_end
1077            IF (var%flag(iff) <= lev_files(iff) .AND. clef_files(iff)) THEN
1078
1079              IF (using_xios) THEN
1080                IF (firstx) THEN
1081                   IF (prt_level >= 10) THEN
1082                      write(lunout,*)'Dans iophy histwrite2D,iff,var%name ',&
1083                                     iff,trim(var%name)                       
1084                      write(lunout,*)"histwrite2d_phy:.NOT.clef_stations(iff)and iff==iff_beg, call xios_send_field"
1085                   ENDIF
1086                   IF (grid_type==regular_lonlat) THEN
1087                     IF (SIZE(field) == klon) then
1088                        CALL xios_send_field(var%name, Field2d)
1089                     ELSE
1090                        CALL xios_send_field(var%name, field)
1091                     ENDIF
1092                   ELSE IF (grid_type==unstructured) THEN
1093                     IF (SIZE(field) == klon) then
1094                       CALL xios_send_field(var%name, buffer_omp)
1095                     ELSE
1096                       CALL xios_send_field(var%name, field)
1097                     ENDIF
1098                   ENDIF
1099
1100                   firstx=.false.
1101                ENDIF
1102              ENDIF
1103
1104                  IF (.NOT.clef_stations(iff)) THEN
1105                        ALLOCATE(index2d(nbp_lon*jj_nb))
1106                        ALLOCATE(fieldok(nbp_lon*jj_nb))
1107#ifndef CPP_IOIPSL_NO_OUTPUT
1108                        CALL histwrite(nid_files(iff),var%name,itau_iophy,Field2d,nbp_lon*jj_nb,index2d)
1109#endif
1110!    IF (using_xios) THEN 
1111!                        IF (iff == iff_beg) THEN
1112!                          IF (prt_level >= 10) THEN
1113!                            write(lunout,*)"histwrite2d_phy: .NOT.clef_stations(iff) and iff==iff_beg, call xios_send_field"
1114!                          ENDIF
1115!                          CALL xios_send_field(var%name, Field2d)
1116!                        ENDIF
1117!    ENDIF
1118                  ELSE
1119                        ALLOCATE(fieldok(npstn))
1120                        ALLOCATE(index2d(npstn))
1121
1122                        IF (is_sequential) THEN
1123                          DO ip=1, npstn
1124                            fieldok(ip)=buffer_omp(nptabij(ip))
1125                          ENDDO
1126                        ELSE
1127                              DO ip=1, npstn
1128                                write(lunout,*)'histwrite2d_phy is_sequential npstn ip namenptabij',npstn,ip,var%name,nptabij(ip)
1129                                     IF(nptabij(ip).GE.klon_mpi_begin.AND. &
1130                                        nptabij(ip).LE.klon_mpi_end) THEN
1131                                       fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
1132                                     ENDIF
1133                              ENDDO
1134                       ENDIF ! of IF (is_sequential)
1135#ifndef CPP_IOIPSL_NO_OUTPUT
1136                       IF (prt_level >= 10) THEN
1137                         write(lunout,*)"histwrite2d_phy: clef_stations(iff) and iff==iff_beg, call wxios_write_2D"
1138                       ENDIF
1139                       CALL histwrite(nid_files(iff),var%name,itau_iophy,fieldok,npstn,index2d)
1140#endif
1141                  ENDIF ! of IF(.NOT.clef_stations(iff))
1142                 
1143                DEALLOCATE(index2d)
1144                DEALLOCATE(fieldok)
1145            ENDIF !levfiles
1146        ENDDO ! of DO iff=iff_beg, iff_end
1147      ENDIF
1148!$OMP END MASTER   
1149  ENDIF ! vars_defined
1150
1151  IF (prt_level >= 10) WRITE(lunout,*)'End histwrite2d_phy ',trim(var%name)
1152
1153  END SUBROUTINE histwrite2d_phy
1154
1155
1156! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE
1157  SUBROUTINE histwrite3d_phy(var, field, STD_iff)
1158
1159  USE mod_phys_lmdz_omp_transfert, ONLY: bcast_omp
1160  USE dimphy, ONLY: klon, klev
1161  USE mod_phys_lmdz_para, ONLY: gather_omp, grid1dto2d_mpi, &
1162                                jj_nb, klon_mpi, klon_mpi_begin, &
1163                                klon_mpi_end, is_sequential, is_master
1164  USE ioipsl, ONLY: histwrite
1165  USE phys_output_var_mod, ONLY: ctrl_out, clef_files, lev_files, &
1166                                 nfiles, vars_defined, clef_stations, &
1167                                 nid_files, swaerofree_diag
1168  USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat, grid_type, regular_lonlat, unstructured
1169  USE lmdz_xios, ONLY: xios_send_field, xios_field_is_active, using_xios
1170  USE print_control_mod, ONLY: prt_level,lunout
1171
1172  IMPLICIT NONE
1173  INCLUDE 'clesphys.h'
1174
1175  TYPE(ctrl_out), INTENT(IN) :: var
1176  REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
1177  INTEGER, INTENT(IN), OPTIONAL :: STD_iff ! ug RUSTINE POUR LES STD LEVS.....
1178     
1179  INTEGER :: iff, iff_beg, iff_end
1180  LOGICAL, SAVE  :: firstx
1181!$OMP THREADPRIVATE(firstx)
1182  REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp
1183  REAL :: Field3d(nbp_lon,jj_nb,SIZE(field,2))
1184  INTEGER :: ip, n, nlev, nlevx
1185  INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
1186  REAL,ALLOCATABLE, DIMENSION(:,:) :: fieldok
1187  logical, save :: is_active = .true.
1188
1189  IF (check_dim .AND. is_master) WRITE(lunout,*)'histwrite3d_phy for ', trim(var%name)
1190
1191  IF (prt_level >= 10) write(lunout,*)'Begin histrwrite3d ',var%name
1192
1193! ug RUSTINE POUR LES STD LEVS.....
1194      IF (PRESENT(STD_iff)) THEN
1195            iff_beg = STD_iff
1196            iff_end = STD_iff
1197      ELSE
1198            iff_beg = 1
1199            iff_end = nfiles
1200      ENDIF
1201
1202  ! On regarde si on est dans la phase de définition ou d'écriture:
1203  IF (.NOT.vars_defined) THEN
1204      !Si phase de définition.... on définit
1205!$OMP MASTER
1206      DO iff=iff_beg, iff_end
1207        IF (clef_files(iff)) THEN
1208          CALL histdef3d(iff, var)
1209        ENDIF
1210      ENDDO
1211!$OMP END MASTER
1212!--broadcasting the flag that have been changed in histdef3d on OMP masters
1213      CALL bcast_omp(swaerofree_diag)
1214  ELSE
1215    IF (using_xios) THEN 
1216      IF (ok_all_xml) THEN
1217        !$omp barrier
1218        !$omp master
1219        is_active = xios_field_is_active(var%name, at_current_timestep_arg=.false.)
1220        !$omp end master
1221        !$omp barrier
1222        IF(.not. is_active) RETURN
1223      ENDIF
1224    ENDIF
1225
1226    !Et sinon on.... écrit
1227    IF (SIZE(field,1)/=klon .AND. SIZE(field,1)/=klev &
1228         .AND. SIZE(field,1)/=klev+1) &
1229         CALL abort_physic('iophy::histwrite3d_phy', &
1230         'Field first DIMENSION not equal to klon/klev',1)
1231
1232    nlev=SIZE(field,2)
1233    nlevx=nlev
1234!    IF (nlev.EQ.klev+1) THEN
1235!        nlevx=klev
1236!    ELSE
1237!        nlevx=nlev
1238!    ENDIF
1239
1240    IF (SIZE(field,1) == klon) then
1241        CALL Gather_omp(field,buffer_omp)
1242    ELSE
1243        buffer_omp(:,:)=0.
1244    ENDIF
1245!$OMP MASTER
1246    IF (grid_type==regular_lonlat) CALL grid1Dto2D_mpi(buffer_omp,field3d)
1247
1248! BOUCLE SUR LES FICHIERS
1249    firstx=.true.
1250
1251    IF (ok_all_xml) THEN
1252      IF (using_xios) THEN 
1253          IF (prt_level >= 10) THEN
1254             write(lunout,*)'Dans iophy histwrite3D,var%name ',&
1255                             trim(var%name)                       
1256          ENDIF
1257          IF (grid_type==regular_lonlat) THEN
1258            IF (SIZE(field,1) == klon) then
1259               CALL xios_send_field(var%name, Field3d(:,:,1:nlevx))
1260            ELSE
1261               CALL xios_send_field(var%name, field)
1262            ENDIF
1263          ELSE IF (grid_type==unstructured) THEN
1264            IF (SIZE(field,1) == klon) then
1265               CALL xios_send_field(var%name, buffer_omp(:,1:nlevx))
1266            ELSE
1267               CALL xios_send_field(var%name, field)
1268            ENDIF
1269        ENDIF
1270
1271      ELSE
1272        CALL abort_physic ('iophy','cannot have ok_all_xml = .T. without CPP_XIOS defined' ,1)
1273      ENDIF
1274    ELSE 
1275
1276      DO iff=iff_beg, iff_end
1277          IF (var%flag(iff) <= lev_files(iff) .AND. clef_files(iff)) THEN
1278            IF (using_xios) THEN 
1279              IF (firstx) THEN
1280                IF (prt_level >= 10) THEN
1281                  write(lunout,*)'Dans iophy, histwrite3D iff nlev klev firstx', &
1282                                  iff,nlev,klev, firstx                       
1283                  write(lunout,*)'histwrite3d_phy: call xios_send_field for ', &
1284                                  trim(var%name), ' with iim jjm nlevx = ', &
1285                                  nbp_lon,jj_nb,nlevx
1286                ENDIF
1287                IF (grid_type==regular_lonlat) THEN
1288                  IF (SIZE(field,1) == klon) then
1289                      CALL xios_send_field(var%name, Field3d(:,:,1:nlevx))
1290                  ELSE
1291                       CALL xios_send_field(var%name, field)
1292                  ENDIF
1293                ELSE IF (grid_type==unstructured) THEN
1294                  IF (SIZE(field,1) == klon) then
1295                     CALL xios_send_field(var%name, buffer_omp(:,1:nlevx))
1296                  ELSE
1297                     CALL xios_send_field(var%name, field)
1298                  ENDIF
1299                ENDIF
1300
1301                firstx=.false.
1302              ENDIF
1303            ENDIF
1304           
1305            IF (.NOT.clef_stations(iff)) THEN
1306                      ALLOCATE(index3d(nbp_lon*jj_nb*nlev))
1307                      ALLOCATE(fieldok(nbp_lon*jj_nb,nlev))
1308
1309#ifndef CPP_IOIPSL_NO_OUTPUT
1310                        CALL histwrite(nid_files(iff),var%name,itau_iophy,Field3d,nbp_lon*jj_nb*nlev,index3d)
1311#endif
1312
1313!    IF (using_xios) THEN 
1314!                        IF (iff == 1) THEN
1315!                              CALL xios_send_field(var%name, Field3d(:,:,1:klev))
1316!                        ENDIF
1317!    ENDIF
1318!                       
1319              ELSE
1320                        nlev=size(field,2)
1321                        ALLOCATE(index3d(npstn*nlev))
1322                        ALLOCATE(fieldok(npstn,nlev))
1323
1324                        IF (is_sequential) THEN
1325                              DO n=1, nlev
1326                                    DO ip=1, npstn
1327                                          fieldok(ip,n)=buffer_omp(nptabij(ip),n)
1328                                    ENDDO
1329                              ENDDO
1330                        ELSE
1331                              DO n=1, nlev
1332                                    DO ip=1, npstn
1333                                                IF(nptabij(ip).GE.klon_mpi_begin.AND. &
1334                                                      nptabij(ip).LE.klon_mpi_end) THEN
1335                                                fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
1336                                          ENDIF
1337                                    ENDDO
1338                              ENDDO
1339                        ENDIF
1340#ifndef CPP_IOIPSL_NO_OUTPUT
1341                        CALL histwrite(nid_files(iff),var%name,itau_iophy,fieldok,npstn*nlev,index3d)
1342#endif
1343              ENDIF
1344              DEALLOCATE(index3d)
1345              DEALLOCATE(fieldok)
1346          ENDIF
1347      ENDDO
1348    ENDIF
1349!$OMP END MASTER   
1350  ENDIF ! vars_defined
1351
1352  IF (prt_level >= 10) write(lunout,*)'End histrwrite3d ',var%name
1353
1354  END SUBROUTINE histwrite3d_phy
1355 
1356
1357! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
1358  SUBROUTINE histwrite2d_xios(field_name,field)
1359
1360  USE dimphy, ONLY: klon, klev
1361  USE mod_phys_lmdz_para, ONLY: gather_omp, grid1Dto2D_mpi, &
1362                                is_sequential, klon_mpi_begin, klon_mpi_end, &
1363                                jj_nb, klon_mpi, is_master
1364  USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat, grid_type, unstructured
1365  USE lmdz_xios, ONLY: xios_send_field
1366  USE print_control_mod, ONLY: prt_level,lunout
1367
1368  IMPLICIT NONE
1369
1370  CHARACTER(LEN=*), INTENT(IN) :: field_name
1371  REAL, DIMENSION(:), INTENT(IN) :: field
1372     
1373  REAL,DIMENSION(klon_mpi) :: buffer_omp
1374  INTEGER, allocatable, DIMENSION(:) :: index2d
1375  REAL :: Field2d(nbp_lon,jj_nb)
1376
1377  INTEGER :: ip
1378  REAL, ALLOCATABLE, DIMENSION(:) :: fieldok
1379
1380  IF (check_dim .AND. is_master) WRITE(lunout,*)'histwrite2d_xios for ', field_name
1381
1382  IF (prt_level >= 10) WRITE(lunout,*)'Begin histrwrite2d_xios ',field_name
1383
1384    !Et sinon on.... écrit
1385    IF (SIZE(field)/=klon .AND. SIZE(field)/=klev .AND. SIZE(field)/=klev+1) CALL abort_physic('iophy::histwrite2d_xios','Field first DIMENSION not equal to klon/klev',1)
1386    IF (SIZE(field) == klev .OR. SIZE(field) == klev+1) then
1387!$OMP MASTER
1388
1389        CALL xios_send_field(field_name,field)
1390!$OMP END MASTER   
1391  ELSE
1392        CALL Gather_omp(field,buffer_omp)   
1393!$OMP MASTER
1394
1395      IF (grid_type==unstructured) THEN
1396 
1397        CALL xios_send_field(field_name, buffer_omp)
1398
1399      ELSE
1400
1401        CALL grid1Dto2D_mpi(buffer_omp,Field2d)
1402   
1403!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1404!ATTENTION, STATIONS PAS GEREES !
1405!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1406    !IF(.NOT.clef_stations(iff)) THEN
1407        IF (.TRUE.) THEN
1408   
1409            CALL xios_send_field(field_name, Field2d)
1410   
1411        ELSE
1412            ALLOCATE(fieldok(npstn))
1413            ALLOCATE(index2d(npstn))
1414   
1415            IF (is_sequential) THEN
1416                DO ip=1, npstn
1417                    fieldok(ip)=buffer_omp(nptabij(ip))
1418                ENDDO
1419            ELSE
1420                DO ip=1, npstn
1421                    PRINT*,'histwrite2d_xios is_sequential npstn ip namenptabij',npstn,ip,field_name,nptabij(ip)
1422                    IF(nptabij(ip).GE.klon_mpi_begin.AND. &
1423                    nptabij(ip).LE.klon_mpi_end) THEN
1424                        fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
1425                    ENDIF
1426                ENDDO
1427            ENDIF
1428            DEALLOCATE(index2d)
1429            DEALLOCATE(fieldok)
1430   
1431        ENDIF                 
1432      ENDIF
1433!$OMP END MASTER   
1434  ENDIF
1435
1436  IF (prt_level >= 10) WRITE(lunout,*)'End histrwrite2d_xios ',field_name
1437  END SUBROUTINE histwrite2d_xios
1438
1439
1440! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE
1441  SUBROUTINE histwrite3d_xios(field_name, field)
1442
1443  USE dimphy, ONLY: klon, klev
1444  USE mod_phys_lmdz_para, ONLY: gather_omp, grid1Dto2D_mpi, &
1445                                is_sequential, klon_mpi_begin, klon_mpi_end, &
1446                                jj_nb, klon_mpi, is_master
1447  USE lmdz_xios, ONLY: xios_send_field
1448  USE mod_grid_phy_lmdz, ONLY : nbp_lon, nbp_lat, grid_type, unstructured
1449  USE print_control_mod, ONLY: prt_level,lunout
1450
1451  IMPLICIT NONE
1452
1453  CHARACTER(LEN=*), INTENT(IN) :: field_name
1454  REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
1455
1456  REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp
1457  REAL :: Field3d(nbp_lon,jj_nb,SIZE(field,2))
1458  INTEGER :: ip, n, nlev
1459  INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
1460  REAL,ALLOCATABLE, DIMENSION(:,:) :: fieldok
1461
1462  IF (check_dim .AND. is_master) WRITE(lunout,*)'histwrite3d_xios for ', field_name
1463
1464  IF (prt_level >= 10) write(lunout,*)'Begin histrwrite3d_xios ',field_name
1465
1466    !Et on.... écrit
1467    IF (SIZE(field,1)/=klon .AND. SIZE(field,1)/=klev .AND. SIZE(field,1)/=klev+1) then
1468      write(lunout,*)' histrwrite3d_xios ', field_name, SIZE(field)
1469      CALL abort_physic('iophy::histwrite3d_xios','Field first DIMENSION not equal to klon/klev',1)
1470    ENDIF
1471   
1472    IF (SIZE(field,1) == klev .OR. SIZE(field,1) == klev+1) then
1473!$OMP MASTER
1474        CALL xios_send_field(field_name,field)
1475!$OMP END MASTER   
1476  ELSE
1477        nlev=SIZE(field,2)
1478
1479
1480        CALL Gather_omp(field,buffer_omp)
1481!$OMP MASTER
1482
1483    IF (grid_type==unstructured) THEN
1484
1485      CALL xios_send_field(field_name, buffer_omp(:,1:nlev))
1486
1487    ELSE
1488        CALL grid1Dto2D_mpi(buffer_omp,field3d)
1489
1490!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1491!ATTENTION, STATIONS PAS GEREES !
1492!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1493    !IF (.NOT.clef_stations(iff)) THEN
1494        IF(.TRUE.)THEN
1495
1496           CALL xios_send_field(field_name, Field3d(:,:,1:nlev))
1497                           
1498        ELSE
1499            nlev=size(field,2)
1500            ALLOCATE(index3d(npstn*nlev))
1501            ALLOCATE(fieldok(npstn,nlev))
1502   
1503            IF (is_sequential) THEN
1504                DO n=1, nlev
1505                    DO ip=1, npstn
1506                        fieldok(ip,n)=buffer_omp(nptabij(ip),n)
1507                    ENDDO
1508                ENDDO
1509            ELSE
1510                DO n=1, nlev
1511                    DO ip=1, npstn
1512                        IF(nptabij(ip).GE.klon_mpi_begin.AND. &
1513                        nptabij(ip).LE.klon_mpi_end) THEN
1514                            fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
1515                        ENDIF
1516                    ENDDO
1517                ENDDO
1518            ENDIF
1519            DEALLOCATE(index3d)
1520            DEALLOCATE(fieldok)
1521        ENDIF
1522      ENDIF
1523!$OMP END MASTER   
1524  ENDIF
1525
1526  IF (prt_level >= 10) write(lunout,*)'End histrwrite3d_xios ',field_name
1527
1528  END SUBROUTINE histwrite3d_xios
1529
1530  SUBROUTINE histwrite0d_xios(field_name, field)
1531  USE lmdz_xios, ONLY: xios_send_field
1532  USE mod_phys_lmdz_para, ONLY: is_master
1533  USE print_control_mod, ONLY: prt_level,lunout
1534  IMPLICIT NONE
1535
1536  CHARACTER(LEN=*), INTENT(IN) :: field_name
1537  REAL, INTENT(IN) :: field ! --> scalar
1538
1539  IF (check_dim .AND. is_master) WRITE(lunout,*)'histwrite0d_xios for ', field_name
1540
1541!$OMP MASTER
1542  CALL xios_send_field(field_name, field)
1543!$OMP END MASTER
1544
1545  END SUBROUTINE histwrite0d_xios
1546
1547END MODULE iophy
Note: See TracBrowser for help on using the repository browser.