source: LMDZ5/branches/LF-private/libf/phydev/iophy.F90 @ 5394

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

Updating/adapting iophy in phydev (part 2/2).
EM

File size: 13.5 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 
14
15#ifdef CPP_XIOS
16! interfaces for both IOIPSL and XIOS
17  INTERFACE histwrite_phy
18    MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy,histwrite2d_xios,histwrite3d_xios
19  END INTERFACE
20#else
21! interfaces for IOIPSL
22  INTERFACE histwrite_phy
23    MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy
24  END INTERFACE
25#endif
26
27#ifdef CPP_XIOS
28! interfaces for both IOIPSL and XIOS
29  INTERFACE histbeg_phy_all
30    MODULE PROCEDURE histbeg_phy, histbeg_phyxios
31  END INTERFACE
32#else
33! interfaces for IOIPSL
34  INTERFACE histbeg_phy_all
35    MODULE PROCEDURE histbeg_phy
36  END INTERFACE
37#endif
38
39contains
40
41  subroutine init_iophy_new(rlat,rlon)
42  USE dimphy
43  USE mod_phys_lmdz_para
44  USE mod_grid_phy_lmdz
45  USE ioipsl
46  implicit none
47  include 'dimensions.h'   
48    real,dimension(klon),intent(in) :: rlon
49    real,dimension(klon),intent(in) :: rlat
50
51    REAL,dimension(klon_glo)        :: rlat_glo
52    REAL,dimension(klon_glo)        :: rlon_glo
53   
54    INTEGER,DIMENSION(2) :: ddid
55    INTEGER,DIMENSION(2) :: dsg
56    INTEGER,DIMENSION(2) :: dsl
57    INTEGER,DIMENSION(2) :: dpf
58    INTEGER,DIMENSION(2) :: dpl
59    INTEGER,DIMENSION(2) :: dhs
60    INTEGER,DIMENSION(2) :: dhe
61    INTEGER :: i   
62
63    CALL gather(rlat,rlat_glo)
64    CALL bcast(rlat_glo)
65    CALL gather(rlon,rlon_glo)
66    CALL bcast(rlon_glo)
67   
68!$OMP MASTER 
69    ALLOCATE(io_lat(jjm+1-1/(iim*jjm)))
70    io_lat(1)=rlat_glo(1)
71    io_lat(jjm+1-1/(iim*jjm))=rlat_glo(klon_glo)
72    IF ((iim*jjm) > 1) then
73      DO i=2,jjm
74        io_lat(i)=rlat_glo(2+(i-2)*iim)
75      ENDDO
76    ENDIF
77
78    ALLOCATE(io_lon(iim))
79    io_lon(:)=rlon_glo(2-1/(iim*jjm):iim+1-1/(iim*jjm))
80!! (I) dtnb   : total number of domains
81!! (I) dnb    : domain number
82!! (I) did(:) : distributed dimensions identifiers
83!!              (up to 5 dimensions are supported)
84!! (I) dsg(:) : total number of points for each dimension
85!! (I) dsl(:) : local number of points for each dimension
86!! (I) dpf(:) : position of first local point for each dimension
87!! (I) dpl(:) : position of last local point for each dimension
88!! (I) dhs(:) : start halo size for each dimension
89!! (I) dhe(:) : end halo size for each dimension
90!! (C) cdnm   : Model domain definition name.
91!!              The names actually supported are :
92!!              "BOX", "APPLE", "ORANGE".
93!!              These names are case insensitive.
94    ddid=(/ 1,2 /)
95    dsg=(/ iim, jjm+1-1/(iim*jjm) /)
96    dsl=(/ iim, jj_nb /)
97    dpf=(/ 1,jj_begin /)
98    dpl=(/ iim, jj_end /)
99    dhs=(/ ii_begin-1,0 /)
100    if (mpi_rank==mpi_size-1) then
101      dhe=(/0,0/)
102    else
103      dhe=(/ iim-ii_end,0 /) 
104    endif
105   
106#ifndef CPP_NO_IOIPSL
107    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
108                      'APPLE',phys_domain_id)
109#endif
110#ifdef CPP_XIOS
111    !Pour els soucis en MPI, réglage du masque:
112    IF (mpi_rank == 0) THEN
113        data_ibegin = 0
114    ELSE
115        data_ibegin = ii_begin - 1
116    END IF
117
118    IF (mpi_rank == mpi_size-1) THEN
119        data_iend = nbp_lon
120    ELSE
121        data_iend = ii_end + 1
122    END IF
123
124    WRITE(*,*) "TOTO mpirank=",mpi_rank,"iibeg=",ii_begin , "jjbeg=",jj_begin,"jjnb=",jj_nb,"jjend=",jj_end
125
126    !On initialise le domaine xios, maintenant que tout est connu:
127    !SUBROUTINE wxios_domain_param(dom_id, is_sequential, ni, nj, ni_glo, nj_glo,        &
128    !                                ibegin, iend, jbegin, jend,                         &
129    !                                data_ni, data_ibegin,                               &
130    !                                io_lat, io_lon)
131    CALL wxios_domain_param("dom_glo", is_sequential, nbp_lon, jj_nb, nbp_lon, nbp_lat, &
132                            1, nbp_lon, ii_begin, ii_end, jj_begin, jj_end,             &
133                            klon_mpi+2*(nbp_lon-1), data_ibegin, data_iend,             &
134                            io_lat, io_lon)
135#endif
136!$OMP END MASTER
137     
138  END SUBROUTINE init_iophy_new
139 
140!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
141 
142  subroutine histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
143  USE dimphy
144  USE mod_phys_lmdz_para
145  use ioipsl
146  use write_field
147  implicit none
148  include 'dimensions.h'
149   
150    character*(*), intent(IN) :: name
151    integer, intent(in) :: itau0
152    real,intent(in) :: zjulian
153    real,intent(in) :: dtime
154    integer,intent(out) :: nhori
155    integer,intent(out) :: nid_day
156
157!$OMP MASTER   
158    if (is_sequential) then
159      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
160                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
161    else
162      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
163                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
164    endif
165!$OMP END MASTER
166 
167  end subroutine histbeg_phy
168
169!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
170
171#ifdef CPP_XIOS
172
173! SUBROUTINE histbeg_phyxios(name,itau0,zjulian,dtime,ffreq,lev,nhori,nid_day)
174 SUBROUTINE histbeg_phyxios(name,ffreq,lev)
175  USE dimphy
176  USE mod_phys_lmdz_para
177!  use ioipsl
178  use write_field
179  IMPLICIT NONE
180  include 'dimensions.h'
181   
182    character*(*), INTENT(IN) :: name
183!    integer, INTENT(IN) :: itau0
184!    REAL,INTENT(IN) :: zjulian
185!    REAL,INTENT(IN) :: dtime
186    character(LEN=*), INTENT(IN) :: ffreq
187    INTEGER,INTENT(IN) :: lev
188!    integer,intent(out) :: nhori
189!    integer,intent(out) :: nid_day
190
191!$OMP MASTER   
192
193    ! ug OMP en chantier...
194    IF((.NOT. is_using_mpi) .OR. is_mpi_root) THEN
195        ! ug Création du fichier
196        CALL wxios_add_file(name, ffreq, lev)
197    END IF
198
199!$OMP END MASTER
200 
201  END SUBROUTINE histbeg_phyxios
202
203#endif
204
205!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
206 
207  subroutine histwrite2d_phy(nid,lpoint,name,itau,field)
208  USE dimphy
209  USE mod_phys_lmdz_para
210  USE ioipsl
211  implicit none
212  include 'dimensions.h'
213   
214    integer,intent(in) :: nid
215    logical,intent(in) :: lpoint
216    character*(*), intent(IN) :: name
217    integer, intent(in) :: itau
218    real,dimension(:),intent(in) :: field
219    REAL,dimension(klon_mpi) :: buffer_omp
220    INTEGER, allocatable, dimension(:) :: index2d
221    REAL :: Field2d(iim,jj_nb)
222
223    integer :: ip
224    real,allocatable,dimension(:) :: fieldok
225
226    IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first dimension not equal to klon',1)
227   
228    CALL Gather_omp(field,buffer_omp)   
229!$OMP MASTER
230    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
231    if(.NOT.lpoint) THEN
232     ALLOCATE(index2d(iim*jj_nb))
233     ALLOCATE(fieldok(iim*jj_nb))
234     CALL histwrite(nid,name,itau,Field2d,iim*jj_nb,index2d)
235    else
236     ALLOCATE(fieldok(npstn))
237     ALLOCATE(index2d(npstn))
238
239     if(is_sequential) then
240!     klon_mpi_begin=1
241!     klon_mpi_end=klon
242      DO ip=1, npstn
243       fieldok(ip)=buffer_omp(nptabij(ip))
244      ENDDO
245     else
246      DO ip=1, npstn
247!     print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
248       IF(nptabij(ip).GE.klon_mpi_begin.AND. &
249          nptabij(ip).LE.klon_mpi_end) THEN
250         fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
251       ENDIF
252      ENDDO
253     endif
254     CALL histwrite(nid,name,itau,fieldok,npstn,index2d)
255!
256    endif
257    deallocate(index2d)
258    deallocate(fieldok)
259!$OMP END MASTER   
260  end subroutine histwrite2d_phy
261
262!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
263
264  subroutine histwrite3d_phy(nid,lpoint,name,itau,field)
265  USE dimphy
266  USE mod_phys_lmdz_para
267
268  use ioipsl
269  implicit none
270  include 'dimensions.h'
271   
272    integer,intent(in) :: nid
273    logical,intent(in) :: lpoint
274    character*(*), intent(IN) :: name
275    integer, intent(in) :: itau
276    real,dimension(:,:),intent(in) :: field  ! --> field(klon,:)
277    REAL,dimension(klon_mpi,size(field,2)) :: buffer_omp
278    REAL :: Field3d(iim,jj_nb,size(field,2))
279    INTEGER :: ip, n, nlev
280    INTEGER, ALLOCATABLE, dimension(:) :: index3d
281    real,allocatable, dimension(:,:) :: fieldok
282
283    IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first dimension not equal to klon',1)
284    nlev=size(field,2)
285
286!   print*,'hist3d_phy mpi_rank npstn=',mpi_rank,npstn
287
288!   DO ip=1, npstn
289!    print*,'hist3d_phy mpi_rank nptabij',mpi_rank,nptabij(ip)
290!   ENDDO
291
292    CALL Gather_omp(field,buffer_omp)
293!$OMP MASTER
294    CALL grid1Dto2D_mpi(buffer_omp,field3d)
295    if(.NOT.lpoint) THEN
296     ALLOCATE(index3d(iim*jj_nb*nlev))
297     ALLOCATE(fieldok(iim*jj_nb,nlev))
298     CALL histwrite(nid,name,itau,Field3d,iim*jj_nb*nlev,index3d)
299    else
300      nlev=size(field,2)
301      ALLOCATE(index3d(npstn*nlev))
302      ALLOCATE(fieldok(npstn,nlev))
303
304      if(is_sequential) then
305!      klon_mpi_begin=1
306!      klon_mpi_end=klon
307       DO n=1, nlev
308       DO ip=1, npstn
309        fieldok(ip,n)=buffer_omp(nptabij(ip),n)
310       ENDDO
311       ENDDO
312      else
313       DO n=1, nlev
314       DO ip=1, npstn
315        IF(nptabij(ip).GE.klon_mpi_begin.AND. &
316         nptabij(ip).LE.klon_mpi_end) THEN
317         fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
318        ENDIF
319       ENDDO
320       ENDDO
321      endif
322      CALL histwrite(nid,name,itau,fieldok,npstn*nlev,index3d)
323    endif
324  deallocate(index3d)
325  deallocate(fieldok)
326!$OMP END MASTER   
327  end subroutine histwrite3d_phy
328
329!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
330
331! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
332#ifdef CPP_XIOS
333  SUBROUTINE histwrite2d_xios(field_name,field)
334  USE dimphy
335  USE mod_phys_lmdz_para
336  USE wxios
337
338
339  IMPLICIT NONE
340  INCLUDE 'dimensions.h'
341  INCLUDE 'iniprint.h'
342
343    CHARACTER(LEN=*), INTENT(IN) :: field_name
344    REAL, DIMENSION(:), INTENT(IN) :: field
345     
346    REAL,DIMENSION(klon_mpi) :: buffer_omp
347    INTEGER, allocatable, DIMENSION(:) :: index2d
348    REAL :: Field2d(iim,jj_nb)
349
350    INTEGER :: ip
351    REAL, ALLOCATABLE, DIMENSION(:) :: fieldok
352
353    IF (prt_level >= 9) WRITE(lunout,*)'Begin histrwrite2d_xios ',field_name
354
355    !Et sinon on.... écrit
356    IF (SIZE(field)/=klon) CALL abort_gcm('iophy::histwrite2d_xios','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   
362!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
363!ATTENTION, STATIONS PAS GEREES !
364!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
365    !IF(.NOT.clef_stations(iff)) THEN
366    IF (.TRUE.) THEN
367        ALLOCATE(index2d(iim*jj_nb))
368        ALLOCATE(fieldok(iim*jj_nb))
369
370
371        CALL wxios_write_2D(field_name, Field2d)
372
373    ELSE
374        ALLOCATE(fieldok(npstn))
375        ALLOCATE(index2d(npstn))
376
377        IF (is_sequential) THEN
378            DO ip=1, npstn
379                fieldok(ip)=buffer_omp(nptabij(ip))
380            ENDDO
381        ELSE
382            DO ip=1, npstn
383                PRINT*,'histwrite2d_xios is_sequential npstn ip namenptabij',npstn,ip,field_name,nptabij(ip)
384                IF(nptabij(ip).GE.klon_mpi_begin.AND. &
385                nptabij(ip).LE.klon_mpi_end) THEN
386                    fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
387                ENDIF
388            ENDDO
389        ENDIF
390
391    ENDIF
392                 
393    deallocate(index2d)
394    deallocate(fieldok)
395!$OMP END MASTER   
396
397  IF (prt_level >= 9) WRITE(lunout,*)'End histrwrite2d_xios ',field_name
398  END SUBROUTINE histwrite2d_xios
399#endif
400
401!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
402
403! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
404#ifdef CPP_XIOS
405  SUBROUTINE histwrite3d_xios(field_name, field)
406  USE dimphy
407  USE mod_phys_lmdz_para
408  USE wxios
409
410
411  IMPLICIT NONE
412  INCLUDE 'dimensions.h'
413  INCLUDE 'iniprint.h'
414
415    CHARACTER(LEN=*), INTENT(IN) :: field_name
416    REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
417
418    REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp
419    REAL :: Field3d(iim,jj_nb,SIZE(field,2))
420    INTEGER :: ip, n, nlev
421    INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
422    REAL,ALLOCATABLE, DIMENSION(:,:) :: fieldok
423
424  IF (prt_level >= 9) write(lunout,*)'Begin histrwrite3d_xios ',field_name
425
426    !Et on.... écrit
427    IF (SIZE(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
428    nlev=SIZE(field,2)
429
430
431    CALL Gather_omp(field,buffer_omp)
432!$OMP MASTER
433    CALL grid1Dto2D_mpi(buffer_omp,field3d)
434
435    IF(.TRUE.)THEN
436        ALLOCATE(index3d(iim*jj_nb*nlev))
437        ALLOCATE(fieldok(iim*jj_nb,nlev))
438        CALL wxios_write_3D(field_name, Field3d(:,:,1:klev))
439                       
440    ELSE
441        nlev=size(field,2)
442        ALLOCATE(index3d(npstn*nlev))
443        ALLOCATE(fieldok(npstn,nlev))
444
445        IF (is_sequential) THEN
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    ENDIF
462    deallocate(index3d)
463    deallocate(fieldok)
464!$OMP END MASTER   
465
466  IF (prt_level >= 9) write(lunout,*)'End histrwrite3d_xios ',field_name
467  END SUBROUTINE histwrite3d_xios
468#endif
469
470end module iophy
Note: See TracBrowser for help on using the repository browser.