source: LMDZ6/branches/contrails/libf/phylmd/iophy.F90 @ 5450

Last change on this file since 5450 was 5310, checked in by abarral, 7 weeks ago

unify abort_gcm
rename wxios -> wxios_mod

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