source: LMDZ6/branches/Optimisation_LMDZ/libf/phylmd/iophy.F90 @ 5423

Last change on this file since 5423 was 3755, checked in by adurocher, 4 years ago

Fixed compilation error without RRTM or without XIOS

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