source: LMDZ6/trunk/libf/phylmd/iophy.F90 @ 3457

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

Wrong order in call sequence of XIOS context initialisation meant that the coupled model hung.
LF

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