source: LMDZ5/branches/IPSLCM6.0.10/libf/phylmd/iophy.F90 @ 5445

Last change on this file since 5445 was 2669, checked in by Laurent Fairhead, 8 years ago

Merged trunk changes r2640:2664 into testing branch

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