source: LMDZ6/branches/DYNAMICO-conv/libf/phylmd/iophy.F90

Last change on this file was 3411, checked in by Laurent Fairhead, 6 years ago

Undoing merge with trunk (r3356) to properly register Yann's latest modifications

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