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

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

Déplacement de nombreuses variables de physiq.F vers phys_local_var_mod.
UG
................................
Moving of numerous vars from physiq.F to phys_local_var_mod.
UG

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 19.8 KB
Line 
1!
2! $Header$
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  INTERFACE histwrite_phy
18    MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy,histwrite2d_phy_old,histwrite3d_phy_old
19  END INTERFACE
20
21  INTERFACE histbeg_phy_all
22    MODULE PROCEDURE histbeg_phy,histbeg_phy_points
23  END INTERFACE
24
25
26CONTAINS
27
28! ug Routine pour définir itau_iophy depuis phys_output_write_mod:
29  SUBROUTINE set_itau_iophy(ito)
30      IMPLICIT NONE
31      INTEGER, INTENT(IN) :: ito
32      itau_iophy = ito
33  END SUBROUTINE
34
35  SUBROUTINE init_iophy_new(rlat,rlon)
36  USE dimphy
37  USE mod_phys_lmdz_para
38  USE mod_grid_phy_lmdz
39  USE ioipsl
40  IMPLICIT NONE
41  INCLUDE 'dimensions.h'   
42    REAL,DIMENSION(klon),INTENT(IN) :: rlon
43    REAL,DIMENSION(klon),INTENT(IN) :: rlat
44
45    REAL,DIMENSION(klon_glo)        :: rlat_glo
46    REAL,DIMENSION(klon_glo)        :: rlon_glo
47   
48    INTEGER,DIMENSION(2) :: ddid
49    INTEGER,DIMENSION(2) :: dsg
50    INTEGER,DIMENSION(2) :: dsl
51    INTEGER,DIMENSION(2) :: dpf
52    INTEGER,DIMENSION(2) :: dpl
53    INTEGER,DIMENSION(2) :: dhs
54    INTEGER,DIMENSION(2) :: dhe
55    INTEGER :: i   
56
57    CALL gather(rlat,rlat_glo)
58    CALL bcast(rlat_glo)
59    CALL gather(rlon,rlon_glo)
60    CALL bcast(rlon_glo)
61   
62!$OMP MASTER 
63    ALLOCATE(io_lat(jjm+1-1/(iim*jjm)))
64    io_lat(1)=rlat_glo(1)
65    io_lat(jjm+1-1/(iim*jjm))=rlat_glo(klon_glo)
66    IF ((iim*jjm) > 1) then
67      DO i=2,jjm
68        io_lat(i)=rlat_glo(2+(i-2)*iim)
69      ENDDO
70    ENDIF
71
72    ALLOCATE(io_lon(iim))
73    io_lon(:)=rlon_glo(2-1/(iim*jjm):iim+1-1/(iim*jjm))
74
75    ddid=(/ 1,2 /)
76    dsg=(/ iim, jjm+1-1/(iim*jjm) /)
77    dsl=(/ iim, jj_nb /)
78    dpf=(/ 1,jj_begin /)
79    dpl=(/ iim, jj_end /)
80    dhs=(/ ii_begin-1,0 /)
81    IF (mpi_rank==mpi_size-1) THEN
82      dhe=(/0,0/)
83    ELSE
84      dhe=(/ iim-ii_end,0 /) 
85    ENDIF
86   
87    CALL flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
88                      'APPLE',phys_domain_id)
89
90!$OMP END MASTER
91     
92  END SUBROUTINE init_iophy_new
93
94  SUBROUTINE init_iophy(lat,lon)
95  USE dimphy
96  USE mod_phys_lmdz_para
97  USE ioipsl
98  IMPLICIT NONE
99  INCLUDE 'dimensions.h'   
100    REAL,DIMENSION(iim),INTENT(IN) :: lon
101    REAL,DIMENSION(jjm+1-1/(iim*jjm)),INTENT(IN) :: lat
102
103    INTEGER,DIMENSION(2) :: ddid
104    INTEGER,DIMENSION(2) :: dsg
105    INTEGER,DIMENSION(2) :: dsl
106    INTEGER,DIMENSION(2) :: dpf
107    INTEGER,DIMENSION(2) :: dpl
108    INTEGER,DIMENSION(2) :: dhs
109    INTEGER,DIMENSION(2) :: dhe
110
111!$OMP MASTER 
112    allocate(io_lat(jjm+1-1/(iim*jjm)))
113    io_lat(:)=lat(:)
114    allocate(io_lon(iim))
115    io_lon(:)=lon(:)
116   
117    ddid=(/ 1,2 /)
118    dsg=(/ iim, jjm+1-1/(iim*jjm) /)
119    dsl=(/ iim, jj_nb /)
120    dpf=(/ 1,jj_begin /)
121    dpl=(/ iim, jj_end /)
122    dhs=(/ ii_begin-1,0 /)
123    if (mpi_rank==mpi_size-1) then
124      dhe=(/0,0/)
125    else
126      dhe=(/ iim-ii_end,0 /) 
127    endif
128   
129    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
130                      'APPLE',phys_domain_id)
131
132!$OMP END MASTER
133     
134  end SUBROUTINE init_iophy
135 
136  SUBROUTINE histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
137  USE dimphy
138  USE mod_phys_lmdz_para
139  use ioipsl
140  use write_field
141  IMPLICIT NONE
142  include 'dimensions.h'
143   
144    character*(*), INTENT(IN) :: name
145    integer, INTENT(IN) :: itau0
146    REAL,INTENT(IN) :: zjulian
147    REAL,INTENT(IN) :: dtime
148    integer,intent(out) :: nhori
149    integer,intent(out) :: nid_day
150
151!$OMP MASTER   
152    if (is_sequential) then
153      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
154                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
155    else
156      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
157                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
158    endif
159!$OMP END MASTER
160 
161  END SUBROUTINE histbeg_phy
162
163  SUBROUTINE histbeg_phy_points(rlon,rlat,pim,tabij,ipt,jpt, &
164             plon,plat,plon_bounds,plat_bounds, &
165             nname,itau0,zjulian,dtime,nnhori,nnid_day)
166  USE dimphy
167  USE mod_phys_lmdz_para
168  USE mod_grid_phy_lmdz
169  use ioipsl
170  use write_field
171  IMPLICIT NONE
172  include 'dimensions.h'
173
174    REAL,DIMENSION(klon),INTENT(IN) :: rlon
175    REAL,DIMENSION(klon),INTENT(IN) :: rlat
176    integer, INTENT(IN) :: itau0
177    REAL,INTENT(IN) :: zjulian
178    REAL,INTENT(IN) :: dtime
179    integer, INTENT(IN) :: pim
180    integer, intent(out) :: nnhori
181    character(len=20), INTENT(IN) :: nname
182    INTEGER, intent(out) :: nnid_day
183    integer :: i
184    REAL,DIMENSION(klon_glo)        :: rlat_glo
185    REAL,DIMENSION(klon_glo)        :: rlon_glo
186    INTEGER, DIMENSION(pim), INTENT(IN)  :: tabij
187    REAL,DIMENSION(pim), INTENT(IN) :: plat, plon
188    INTEGER,DIMENSION(pim), INTENT(IN) :: ipt, jpt
189    REAL,DIMENSION(pim,2), intent(out) :: plat_bounds, plon_bounds
190
191    INTEGER, SAVE :: tabprocbeg, tabprocend
192!$OMP THREADPRIVATE(tabprocbeg, tabprocend)
193    INTEGER :: ip
194    INTEGER, PARAMETER :: nip=1
195    INTEGER :: npproc
196    REAL, allocatable, DIMENSION(:) :: npplat, npplon
197    REAL, allocatable, DIMENSION(:,:) :: npplat_bounds, npplon_bounds
198    INTEGER, PARAMETER :: jjmp1=jjm+1-1/jjm
199    REAL, DIMENSION(iim,jjmp1) :: zx_lon, zx_lat
200
201    CALL gather(rlat,rlat_glo)
202    CALL bcast(rlat_glo)
203    CALL gather(rlon,rlon_glo)
204    CALL bcast(rlon_glo)
205
206!$OMP MASTER
207    DO i=1,pim
208
209!    print*,'CFMIP_iophy i tabij lon lat',i,tabij(i),plon(i),plat(i)
210
211     plon_bounds(i,1)=rlon_glo(tabij(i)-1)
212     plon_bounds(i,2)=rlon_glo(tabij(i)+1)
213     if(plon_bounds(i,2).LE.0..AND.plon_bounds(i,1).GE.0.) THEN
214      if(rlon_glo(tabij(i)).GE.0.) THEN
215       plon_bounds(i,2)=-1*plon_bounds(i,2)
216      endif
217     endif
218     if(plon_bounds(i,2).GE.0..AND.plon_bounds(i,1).LE.0.) THEN
219      if(rlon_glo(tabij(i)).LE.0.) THEN
220       plon_bounds(i,2)=-1*plon_bounds(i,2)
221      endif
222     endif
223!
224     IF ( tabij(i).LE.iim) THEN
225      plat_bounds(i,1)=rlat_glo(tabij(i))
226     ELSE
227      plat_bounds(i,1)=rlat_glo(tabij(i)-iim)
228     ENDIF
229     plat_bounds(i,2)=rlat_glo(tabij(i)+iim)
230!
231!    print*,'CFMIP_iophy point i lon lon_bds',i,plon_bounds(i,1),rlon_glo(tabij(i)),plon_bounds(i,2)
232!    print*,'CFMIP_iophy point i lat lat_bds',i,plat_bounds(i,1),rlat_glo(tabij(i)),plat_bounds(i,2)
233!
234    ENDDO
235    if (is_sequential) then
236
237     npstn=pim
238     IF(.NOT. ALLOCATED(nptabij)) THEN
239      ALLOCATE(nptabij(pim))
240     ENDIF
241     DO i=1,pim
242      nptabij(i)=tabij(i)
243     ENDDO
244
245       CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlon_glo,zx_lon)
246       if ((iim*jjm).gt.1) then
247       DO i = 1, iim
248         zx_lon(i,1) = rlon_glo(i+1)
249         zx_lon(i,jjmp1) = rlon_glo(i+1)
250       ENDDO
251       endif
252       CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlat_glo,zx_lat)
253
254    DO i=1,pim
255!    print*,'CFMIP_iophy i tabij lon lat',i,tabij(i),plon(i),plat(i)
256
257     plon_bounds(i,1)=zx_lon(ipt(i)-1,jpt(i))
258     plon_bounds(i,2)=zx_lon(ipt(i)+1,jpt(i))
259
260     if (ipt(i).EQ.1) then
261      plon_bounds(i,1)=zx_lon(iim,jpt(i))
262      plon_bounds(i,2)=360.+zx_lon(ipt(i)+1,jpt(i))
263     endif
264 
265     if (ipt(i).EQ.iim) then
266      plon_bounds(i,2)=360.+zx_lon(1,jpt(i))
267     endif
268
269     plat_bounds(i,1)=zx_lat(ipt(i),jpt(i)-1)
270     plat_bounds(i,2)=zx_lat(ipt(i),jpt(i)+1)
271
272     if (jpt(i).EQ.1) then
273      plat_bounds(i,1)=zx_lat(ipt(i),1)+0.001
274      plat_bounds(i,2)=zx_lat(ipt(i),1)-0.001
275     endif
276 
277     if (jpt(i).EQ.jjmp1) then
278      plat_bounds(i,1)=zx_lat(ipt(i),jjmp1)+0.001
279      plat_bounds(i,2)=zx_lat(ipt(i),jjmp1)-0.001
280     endif
281!
282!    print*,'CFMIP_iophy point i lon lon_bds',i,plon_bounds(i,1),rlon(tabij(i)),plon_bounds(i,2)
283!    print*,'CFMIP_iophy point i lat lat_bds',i,plat_bounds(i,1),rlat(tabij(i)),plat_bounds(i,2)
284!
285    ENDDO
286!    print*,'iophy is_sequential nname, nnhori, nnid_day=',trim(nname),nnhori,nnid_day
287     call histbeg(nname,pim,plon,plon_bounds, &
288                           plat,plat_bounds, &
289                           itau0, zjulian, dtime, nnhori, nnid_day)
290    else
291     npproc=0
292     DO ip=1, pim
293      tabprocbeg=klon_mpi_begin
294      tabprocend=klon_mpi_end
295      IF(tabij(ip).GE.tabprocbeg.AND.tabij(ip).LE.tabprocend) THEN
296       npproc=npproc+1
297       npstn=npproc
298      ENDIF
299     ENDDO
300!    print*,'CFMIP_iophy mpi_rank npstn',mpi_rank,npstn
301     IF(.NOT. ALLOCATED(nptabij)) THEN
302      ALLOCATE(nptabij(npstn))
303      ALLOCATE(npplon(npstn), npplat(npstn))
304      ALLOCATE(npplon_bounds(npstn,2), npplat_bounds(npstn,2))
305     ENDIF
306     npproc=0
307     DO ip=1, pim
308      IF(tabij(ip).GE.tabprocbeg.AND.tabij(ip).LE.tabprocend) THEN
309       npproc=npproc+1
310       nptabij(npproc)=tabij(ip)
311!      print*,'mpi_rank npproc ip plon plat tabij=',mpi_rank,npproc,ip, &
312!      plon(ip),plat(ip),tabij(ip)
313       npplon(npproc)=plon(ip)
314       npplat(npproc)=plat(ip)
315       npplon_bounds(npproc,1)=plon_bounds(ip,1)
316       npplon_bounds(npproc,2)=plon_bounds(ip,2)
317       npplat_bounds(npproc,1)=plat_bounds(ip,1)
318       npplat_bounds(npproc,2)=plat_bounds(ip,2)
319!!!
320!!! print qui sert a reordonner les points stations selon l'ordre CFMIP
321!!! ne pas enlever
322        print*,'iophy_mpi rank ip lon lat',mpi_rank,ip,plon(ip),plat(ip)
323!!!
324      ENDIF
325     ENDDO
326     call histbeg(nname,npstn,npplon,npplon_bounds, &
327                            npplat,npplat_bounds, &
328                            itau0,zjulian,dtime,nnhori,nnid_day,phys_domain_id)
329    endif
330!$OMP END MASTER
331
332  end SUBROUTINE histbeg_phy_points
333 
334  SUBROUTINE histwrite2d_phy_old(nid,lpoint,name,itau,field)
335  USE dimphy
336  USE mod_phys_lmdz_para
337  USE phys_output_var_mod
338  USE ioipsl
339  IMPLICIT NONE
340  include 'dimensions.h'
341  include 'iniprint.h'
342   
343    integer,INTENT(IN) :: nid
344    logical,INTENT(IN) :: lpoint
345    character*(*), INTENT(IN) :: name
346    integer, INTENT(IN) :: itau
347    REAL,DIMENSION(:),INTENT(IN) :: field
348    REAL,DIMENSION(klon_mpi) :: buffer_omp
349    INTEGER, allocatable, DIMENSION(:) :: index2d
350    REAL :: Field2d(iim,jj_nb)
351
352    integer :: ip
353    REAL,allocatable,DIMENSION(:) :: fieldok
354
355
356    IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first DIMENSION not equal to klon',1)
357   
358    CALL Gather_omp(field,buffer_omp)   
359!$OMP MASTER
360    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
361    if(.NOT.lpoint) THEN
362     ALLOCATE(index2d(iim*jj_nb))
363     ALLOCATE(fieldok(iim*jj_nb))
364     IF (prt_level >= 9) write(lunout,*)'Sending ',name,' to IOIPSL'
365     CALL histwrite(nid,name,itau,Field2d,iim*jj_nb,index2d)
366     IF (prt_level >= 9) write(lunout,*)'Finished sending ',name,' to IOIPSL'
367    else
368     ALLOCATE(fieldok(npstn))
369     ALLOCATE(index2d(npstn))
370
371     if(is_sequential) then
372!     klon_mpi_begin=1
373!     klon_mpi_end=klon
374      DO ip=1, npstn
375       fieldok(ip)=buffer_omp(nptabij(ip))
376      ENDDO
377     else
378      DO ip=1, npstn
379!     print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
380       IF(nptabij(ip).GE.klon_mpi_begin.AND. &
381          nptabij(ip).LE.klon_mpi_end) THEN
382         fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
383       ENDIF
384      ENDDO
385     endif
386     IF (prt_level >= 9) write(lunout,*)'Sending ',name,' to IOIPSL'
387     CALL histwrite(nid,name,itau,fieldok,npstn,index2d)
388     IF (prt_level >= 9) write(lunout,*)'Finished sending ',name,' to IOIPSL'
389!
390    endif
391    deallocate(index2d)
392    deallocate(fieldok)
393!$OMP END MASTER   
394
395 
396  end SUBROUTINE histwrite2d_phy_old
397
398  SUBROUTINE histwrite3d_phy_old(nid,lpoint,name,itau,field)
399  USE dimphy
400  USE mod_phys_lmdz_para
401  USE phys_output_var_mod
402
403  use ioipsl
404  IMPLICIT NONE
405  include 'dimensions.h'
406  include 'iniprint.h'
407   
408    integer,INTENT(IN) :: nid
409    logical,INTENT(IN) :: lpoint
410    character*(*), INTENT(IN) :: name
411    integer, INTENT(IN) :: itau
412    REAL,DIMENSION(:,:),INTENT(IN) :: field  ! --> field(klon,:)
413    REAL,DIMENSION(klon_mpi,size(field,2)) :: buffer_omp
414    REAL :: Field3d(iim,jj_nb,size(field,2))
415    INTEGER :: ip, n, nlev
416    INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
417    REAL,allocatable, DIMENSION(:,:) :: fieldok
418
419
420    IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
421    nlev=size(field,2)
422
423!   print*,'hist3d_phy mpi_rank npstn=',mpi_rank,npstn
424
425!   DO ip=1, npstn
426!    print*,'hist3d_phy mpi_rank nptabij',mpi_rank,nptabij(ip)
427!   ENDDO
428
429    CALL Gather_omp(field,buffer_omp)
430!$OMP MASTER
431    CALL grid1Dto2D_mpi(buffer_omp,field3d)
432    if(.NOT.lpoint) THEN
433     ALLOCATE(index3d(iim*jj_nb*nlev))
434     ALLOCATE(fieldok(iim*jj_nb,nlev))
435     IF (prt_level >= 9) write(lunout,*)'Sending ',name,' to IOIPSL'
436     CALL histwrite(nid,name,itau,Field3d,iim*jj_nb*nlev,index3d)
437     IF (prt_level >= 9) write(lunout,*)'Finished sending ',name,' to IOIPSL'
438   else
439      nlev=size(field,2)
440      ALLOCATE(index3d(npstn*nlev))
441      ALLOCATE(fieldok(npstn,nlev))
442
443      if(is_sequential) then
444!      klon_mpi_begin=1
445!      klon_mpi_end=klon
446       DO n=1, nlev
447       DO ip=1, npstn
448        fieldok(ip,n)=buffer_omp(nptabij(ip),n)
449       ENDDO
450       ENDDO
451      else
452       DO n=1, nlev
453       DO ip=1, npstn
454        IF(nptabij(ip).GE.klon_mpi_begin.AND. &
455         nptabij(ip).LE.klon_mpi_end) THEN
456         fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
457        ENDIF
458       ENDDO
459       ENDDO
460      endif
461      IF (prt_level >= 9) write(lunout,*)'Sending ',name,' to IOIPSL'
462      CALL histwrite(nid,name,itau,fieldok,npstn*nlev,index3d)
463      IF (prt_level >= 9) write(lunout,*)'Finished sending ',name,' to IOIPSL'
464    endif
465  deallocate(index3d)
466  deallocate(fieldok)
467!$OMP END MASTER   
468
469  end SUBROUTINE histwrite3d_phy_old
470
471
472! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE
473  SUBROUTINE histwrite2d_phy(var,field, STD_iff)
474  USE dimphy
475  USE mod_phys_lmdz_para
476  USE ioipsl
477!Pour avoir nfiles, nidfiles tout ça tout ça...
478  USE phys_output_var_mod
479 
480 
481
482#ifdef CPP_XIOS
483!  USE WXIOS
484#endif
485
486  IMPLICIT NONE
487  include 'dimensions.h'
488   
489!    integer,INTENT(IN) :: nid
490!    logical,INTENT(IN) :: lpoint
491!    character*(*), INTENT(IN) :: name
492!    integer, INTENT(IN) :: itau
493!    REAL,DIMENSION(:),INTENT(IN) :: field
494
495      TYPE(ctrl_out), INTENT(IN) :: var
496      REAL, DIMENSION(:), INTENT(IN) :: field
497      INTEGER, INTENT(IN), OPTIONAL :: STD_iff ! ug RUSTINE POUR LES STD LEVS.....
498     
499      INTEGER :: iff, iff_beg, iff_end
500     
501    REAL,DIMENSION(klon_mpi) :: buffer_omp
502    INTEGER, allocatable, DIMENSION(:) :: index2d
503    REAL :: Field2d(iim,jj_nb)
504
505    INTEGER :: ip
506    REAL, ALLOCATABLE, DIMENSION(:) :: fieldok
507
508! ug RUSTINE POUR LES STD LEVS.....
509      IF (PRESENT(STD_iff)) THEN
510            iff_beg = STD_iff
511            iff_end = STD_iff
512      ELSE
513            iff_beg = 1
514            iff_end = nfiles
515      END IF
516
517    IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first DIMENSION not equal to klon',1)
518   
519    CALL Gather_omp(field,buffer_omp)   
520!$OMP MASTER
521    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
522   
523! La boucle sur les fichiers:
524      DO iff=iff_beg, iff_end
525            IF (var%flag(iff) <= lev_files(iff) .AND. clef_files(iff)) THEN
526   
527                  IF(.NOT.clef_stations(iff)) THEN
528                        ALLOCATE(index2d(iim*jj_nb))
529                        ALLOCATE(fieldok(iim*jj_nb))
530     
531                        CALL histwrite(nid_files(iff),var%name,itau_iophy,Field2d,iim*jj_nb,index2d)
532#ifdef CPP_XIOS
533!                        IF (iff .EQ. 1) THEN
534!                              CALL wxios_write_2D(var%name, Field2d)
535!                        ENDIF
536#endif
537                  ELSE
538                        ALLOCATE(fieldok(npstn))
539                        ALLOCATE(index2d(npstn))
540
541                        IF (is_sequential) THEN
542!                            klon_mpi_begin=1
543!                             klon_mpi_end=klon
544                              DO ip=1, npstn
545                                    fieldok(ip)=buffer_omp(nptabij(ip))
546                              ENDDO
547                             ELSE
548                              DO ip=1, npstn
549!                                   print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
550                                     IF(nptabij(ip).GE.klon_mpi_begin.AND. &
551                                        nptabij(ip).LE.klon_mpi_end) THEN
552                                       fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
553                                     ENDIF
554                              ENDDO
555                       ENDIF
556     
557                       CALL histwrite(nid_files(iff),var%name,itau_iophy,fieldok,npstn,index2d)
558                  ENDIF
559                 
560                deallocate(index2d)
561                deallocate(fieldok)
562            ENDIF !levfiles
563      ENDDO
564!$OMP END MASTER   
565
566  END SUBROUTINE histwrite2d_phy
567
568
569! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE
570  SUBROUTINE histwrite3d_phy(var, field)
571  USE dimphy
572  USE mod_phys_lmdz_para
573
574  use ioipsl
575!Pour avoir nfiles, nidfiles tout ça tout ça...
576  USE phys_output_var_mod 
577 
578
579#ifdef CPP_XIOS
580! USE WXIOS
581#endif
582
583
584  IMPLICIT NONE
585  include 'dimensions.h'
586   
587!    integer,INTENT(IN) :: nid
588!    logical,INTENT(IN) :: lpoint
589!    character*(*), INTENT(IN) :: name
590!    integer, INTENT(IN) :: itau
591!    REAL,DIMENSION(:,:),INTENT(IN) :: field  ! --> field(klon,:)
592
593      TYPE(ctrl_out), INTENT(IN) :: var
594      REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
595
596
597    REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp
598    REAL :: Field3d(iim,jj_nb,SIZE(field,2))
599    INTEGER :: ip, n, nlev, iff
600    INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
601    REAL,ALLOCATABLE, DIMENSION(:,:) :: fieldok
602
603    IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
604    nlev=size(field,2)
605
606!   print*,'hist3d_phy mpi_rank npstn=',mpi_rank,npstn
607
608!   DO ip=1, npstn
609!    print*,'hist3d_phy mpi_rank nptabij',mpi_rank,nptabij(ip)
610!   ENDDO
611
612    CALL Gather_omp(field,buffer_omp)
613!$OMP MASTER
614    CALL grid1Dto2D_mpi(buffer_omp,field3d)
615
616
617! BOUCLE SUR LES FICHIERS
618      DO iff=1, nfiles
619            IF (var%flag(iff) <= lev_files(iff) .AND. clef_files(iff)) THEN
620                IF (.NOT.clef_stations(iff)) THEN
621                        ALLOCATE(index3d(iim*jj_nb*nlev))
622                        ALLOCATE(fieldok(iim*jj_nb,nlev))
623                        CALL histwrite(nid_files(iff),var%name,itau_iophy,Field3d,iim*jj_nb*nlev,index3d)
624#ifdef CPP_XIOS
625!                        IF (iff .EQ. 1) THEN
626!                              CALL wxios_write_3D(var%name, Field3d(:,:,1:klev))
627!                        ENDIF
628#endif
629                       
630                ELSE
631                        nlev=size(field,2)
632                        ALLOCATE(index3d(npstn*nlev))
633                        ALLOCATE(fieldok(npstn,nlev))
634
635                        IF (is_sequential) THEN
636            !                  klon_mpi_begin=1
637            !                  klon_mpi_end=klon
638                              DO n=1, nlev
639                                    DO ip=1, npstn
640                                          fieldok(ip,n)=buffer_omp(nptabij(ip),n)
641                                    ENDDO
642                              ENDDO
643                        ELSE
644                              DO n=1, nlev
645                                    DO ip=1, npstn
646                                                IF(nptabij(ip).GE.klon_mpi_begin.AND. &
647                                                      nptabij(ip).LE.klon_mpi_end) THEN
648                                                fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
649                                          ENDIF
650                                    ENDDO
651                              ENDDO
652                        ENDIF
653                        CALL histwrite(nid_files(iff),var%name,itau_iophy,fieldok,npstn*nlev,index3d)
654                  ENDIF
655                  deallocate(index3d)
656                  deallocate(fieldok)
657            ENDIF
658      ENDDO
659!$OMP END MASTER   
660  END SUBROUTINE histwrite3d_phy
661 
662end module iophy
Note: See TracBrowser for help on using the repository browser.