source: LMDZ5/branches/LF-private/libf/phylmd/iophy.F90 @ 5442

Last change on this file since 5442 was 1852, checked in by Ehouarn Millour, 11 years ago

Implémentation des sorties XIOS dans LMDZ. Activation via -cpp CPP_XIOS.
ATTENTION: un problème de raccord subsiste en mode MPI !
UG
................................
Adding XIOS output to LMDZ. Activated by the CPP_XIOS key.
WARNING: buggy for now in MPI mode.
UG

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