source: LMDZ5/trunk/libf/phylmd/iophy.F90 @ 2854

Last change on this file since 2854 was 2854, checked in by oboucher, 7 years ago

Introducing dry AOD diagnostics for the total aerosols and specieswise
The calculations are only performed if the diagnostics are requested

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