source: LMDZ5/trunk/libf/phydev/iophy.F90 @ 1985

Last change on this file since 1985 was 1907, checked in by lguez, 11 years ago

Added a copyright property to every file of the distribution, except
for the fcm files (which have their own copyright). Use svn propget on
a file to see the copyright. For instance:

$ svn propget copyright libf/phylmd/physiq.F90
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

Also added the files defining the CeCILL version 2 license, in French
and English, at the top of the LMDZ tree.

  • 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
File size: 12.7 KB
Line 
1!
2! $Id: $
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, only: klon
43  USE mod_phys_lmdz_para, only: gather, bcast, &
44                                jj_nb, jj_begin, jj_end, ii_begin, ii_end, &
45                                mpi_size, mpi_rank, klon_mpi, &
46                                is_sequential, is_south_pole
47  USE mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, klon_glo
48#ifdef CPP_IOIPSL
49  USE ioipsl, only: flio_dom_set
50#endif
51#ifdef CPP_XIOS
52  use wxios, only: wxios_domain_param
53#endif
54  implicit none
55  include 'dimensions.h'
56  include 'iniprint.h'
57    real,dimension(klon),intent(in) :: rlon
58    real,dimension(klon),intent(in) :: rlat
59
60    REAL,dimension(klon_glo)        :: rlat_glo
61    REAL,dimension(klon_glo)        :: rlon_glo
62   
63    INTEGER,DIMENSION(2) :: ddid
64    INTEGER,DIMENSION(2) :: dsg
65    INTEGER,DIMENSION(2) :: dsl
66    INTEGER,DIMENSION(2) :: dpf
67    INTEGER,DIMENSION(2) :: dpl
68    INTEGER,DIMENSION(2) :: dhs
69    INTEGER,DIMENSION(2) :: dhe
70    INTEGER :: i   
71    integer :: data_ibegin,data_iend
72
73    CALL gather(rlat,rlat_glo)
74    CALL bcast(rlat_glo)
75    CALL gather(rlon,rlon_glo)
76    CALL bcast(rlon_glo)
77   
78!$OMP MASTER 
79    ALLOCATE(io_lat(jjm+1-1/(iim*jjm)))
80    io_lat(1)=rlat_glo(1)
81    io_lat(jjm+1-1/(iim*jjm))=rlat_glo(klon_glo)
82    IF ((iim*jjm) > 1) then
83      DO i=2,jjm
84        io_lat(i)=rlat_glo(2+(i-2)*iim)
85      ENDDO
86    ENDIF
87
88    ALLOCATE(io_lon(iim))
89    io_lon(:)=rlon_glo(2-1/(iim*jjm):iim+1-1/(iim*jjm))
90!! (I) dtnb   : total number of domains
91!! (I) dnb    : domain number
92!! (I) did(:) : distributed dimensions identifiers
93!!              (up to 5 dimensions are supported)
94!! (I) dsg(:) : total number of points for each dimension
95!! (I) dsl(:) : local number of points for each dimension
96!! (I) dpf(:) : position of first local point for each dimension
97!! (I) dpl(:) : position of last local point for each dimension
98!! (I) dhs(:) : start halo size for each dimension
99!! (I) dhe(:) : end halo size for each dimension
100!! (C) cdnm   : Model domain definition name.
101!!              The names actually supported are :
102!!              "BOX", "APPLE", "ORANGE".
103!!              These names are case insensitive.
104    ddid=(/ 1,2 /)
105    dsg=(/ iim, jjm+1-1/(iim*jjm) /)
106    dsl=(/ iim, jj_nb /)
107    dpf=(/ 1,jj_begin /)
108    dpl=(/ iim, jj_end /)
109    dhs=(/ ii_begin-1,0 /)
110    if (mpi_rank==mpi_size-1) then
111      dhe=(/0,0/)
112    else
113      dhe=(/ iim-ii_end,0 /) 
114    endif
115   
116#ifndef CPP_NO_IOIPSL
117    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
118                      'APPLE',phys_domain_id)
119#endif
120#ifdef CPP_XIOS
121    ! Set values for the mask:
122    IF (mpi_rank == 0) THEN
123        data_ibegin = 0
124    ELSE
125        data_ibegin = ii_begin - 1
126    END IF
127
128    IF (mpi_rank == mpi_size-1) THEN
129        data_iend = nbp_lon
130    ELSE
131        data_iend = ii_end + 1
132    END IF
133
134    if (prt_level>=10) then
135      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
136      write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," nbp_lon=",nbp_lon," nbp_lat=",nbp_lat
137      write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," data_ibegin=",data_ibegin," data_iend=",data_iend
138      write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," data_ibegin=",data_ibegin," data_iend=",data_iend
139      write(lunout,*) "init_iophy_new: mpirank=",mpi_rank," is_south_pole=",is_south_pole
140    endif
141
142    ! Initialize the XIOS domain coreesponding to this process:
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,is_south_pole,mpi_rank)
147#endif
148!$OMP END MASTER
149     
150  END SUBROUTINE init_iophy_new
151 
152!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
153 
154  subroutine histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
155  USE mod_phys_lmdz_para, only: is_sequential, jj_begin, jj_end, jj_nb
156  use ioipsl, only: histbeg
157  implicit none
158  include 'dimensions.h'
159   
160    character*(*), intent(IN) :: name
161    integer, intent(in) :: itau0
162    real,intent(in) :: zjulian
163    real,intent(in) :: dtime
164    integer,intent(out) :: nhori
165    integer,intent(out) :: nid_day
166
167!$OMP MASTER   
168    if (is_sequential) then
169      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
170                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
171    else
172      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
173                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
174    endif
175!$OMP END MASTER
176 
177  end subroutine histbeg_phy
178
179!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
180
181#ifdef CPP_XIOS
182
183! SUBROUTINE histbeg_phyxios(name,itau0,zjulian,dtime,ffreq,lev,nhori,nid_day)
184 SUBROUTINE histbeg_phyxios(name,ffreq,lev)
185  USE mod_phys_lmdz_para, only: is_using_mpi, is_mpi_root
186  use wxios, only: wxios_add_file
187  IMPLICIT NONE
188  include 'dimensions.h'
189   
190    character*(*), INTENT(IN) :: name
191!    integer, INTENT(IN) :: itau0
192!    REAL,INTENT(IN) :: zjulian
193!    REAL,INTENT(IN) :: dtime
194    character(LEN=*), INTENT(IN) :: ffreq
195    INTEGER,INTENT(IN) :: lev
196!    integer,intent(out) :: nhori
197!    integer,intent(out) :: nid_day
198
199!$OMP MASTER   
200
201    ! ug OMP en chantier...
202    IF((.NOT. is_using_mpi) .OR. is_mpi_root) THEN
203        ! ug Création du fichier
204        CALL wxios_add_file(name, ffreq, lev)
205    END IF
206
207!$OMP END MASTER
208 
209  END SUBROUTINE histbeg_phyxios
210
211#endif
212
213!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
214 
215  subroutine histwrite2d_phy(nid,lpoint,name,itau,field)
216  USE dimphy, only: klon
217  USE mod_phys_lmdz_para, only: Gather_omp, grid1Dto2D_mpi, &
218                                is_sequential, klon_mpi_begin, klon_mpi_end, &
219                                jj_nb, klon_mpi
220  USE ioipsl, only: histwrite
221  implicit none
222  include 'dimensions.h'
223   
224    integer,intent(in) :: nid
225    logical,intent(in) :: lpoint
226    character*(*), intent(IN) :: name
227    integer, intent(in) :: itau
228    real,dimension(:),intent(in) :: field
229    REAL,dimension(klon_mpi) :: buffer_omp
230    INTEGER, allocatable, dimension(:) :: index2d
231    REAL :: Field2d(iim,jj_nb)
232
233    integer :: ip
234    real,allocatable,dimension(:) :: fieldok
235
236    IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first dimension not equal to klon',1)
237   
238    CALL Gather_omp(field,buffer_omp)   
239!$OMP MASTER
240    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
241    if(.NOT.lpoint) THEN
242     ALLOCATE(index2d(iim*jj_nb))
243     ALLOCATE(fieldok(iim*jj_nb))
244     CALL histwrite(nid,name,itau,Field2d,iim*jj_nb,index2d)
245    else
246     ALLOCATE(fieldok(npstn))
247     ALLOCATE(index2d(npstn))
248
249     if(is_sequential) then
250!     klon_mpi_begin=1
251!     klon_mpi_end=klon
252      DO ip=1, npstn
253       fieldok(ip)=buffer_omp(nptabij(ip))
254      ENDDO
255     else
256      DO ip=1, npstn
257!     print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
258       IF(nptabij(ip).GE.klon_mpi_begin.AND. &
259          nptabij(ip).LE.klon_mpi_end) THEN
260         fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
261       ENDIF
262      ENDDO
263     endif
264     CALL histwrite(nid,name,itau,fieldok,npstn,index2d)
265!
266    endif
267    deallocate(index2d)
268    deallocate(fieldok)
269!$OMP END MASTER   
270  end subroutine histwrite2d_phy
271
272!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
273
274  subroutine histwrite3d_phy(nid,lpoint,name,itau,field)
275  USE dimphy, only: klon
276  USE mod_phys_lmdz_para, only: Gather_omp, grid1Dto2D_mpi, &
277                                is_sequential, klon_mpi_begin, klon_mpi_end, &
278                                jj_nb, klon_mpi
279  USE ioipsl, only: histwrite
280  implicit none
281  include 'dimensions.h'
282   
283    integer,intent(in) :: nid
284    logical,intent(in) :: lpoint
285    character*(*), intent(IN) :: name
286    integer, intent(in) :: itau
287    real,dimension(:,:),intent(in) :: field  ! --> field(klon,:)
288    REAL,dimension(klon_mpi,size(field,2)) :: buffer_omp
289    REAL :: Field3d(iim,jj_nb,size(field,2))
290    INTEGER :: ip, n, nlev
291    INTEGER, ALLOCATABLE, dimension(:) :: index3d
292    real,allocatable, dimension(:,:) :: fieldok
293
294    IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first dimension not equal to klon',1)
295    nlev=size(field,2)
296
297    CALL Gather_omp(field,buffer_omp)
298!$OMP MASTER
299    CALL grid1Dto2D_mpi(buffer_omp,field3d)
300    if(.NOT.lpoint) THEN
301     ALLOCATE(index3d(iim*jj_nb*nlev))
302     ALLOCATE(fieldok(iim*jj_nb,nlev))
303     CALL histwrite(nid,name,itau,Field3d,iim*jj_nb*nlev,index3d)
304    else
305      nlev=size(field,2)
306      ALLOCATE(index3d(npstn*nlev))
307      ALLOCATE(fieldok(npstn,nlev))
308
309      if(is_sequential) then
310!      klon_mpi_begin=1
311!      klon_mpi_end=klon
312       DO n=1, nlev
313       DO ip=1, npstn
314        fieldok(ip,n)=buffer_omp(nptabij(ip),n)
315       ENDDO
316       ENDDO
317      else
318       DO n=1, nlev
319       DO ip=1, npstn
320        IF(nptabij(ip).GE.klon_mpi_begin.AND. &
321         nptabij(ip).LE.klon_mpi_end) THEN
322         fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
323        ENDIF
324       ENDDO
325       ENDDO
326      endif
327      CALL histwrite(nid,name,itau,fieldok,npstn*nlev,index3d)
328    endif
329  deallocate(index3d)
330  deallocate(fieldok)
331!$OMP END MASTER   
332  end subroutine histwrite3d_phy
333
334!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
335
336! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
337#ifdef CPP_XIOS
338  SUBROUTINE histwrite2d_xios(field_name,field)
339  USE dimphy, only: klon
340  USE mod_phys_lmdz_para, only: gather_omp, grid1Dto2D_mpi, &
341                                jj_nb, klon_mpi
342  USE wxios, only: wxios_write_2D
343
344
345  IMPLICIT NONE
346  INCLUDE 'dimensions.h'
347  INCLUDE 'iniprint.h'
348
349    CHARACTER(LEN=*), INTENT(IN) :: field_name
350    REAL, DIMENSION(:), INTENT(IN) :: field
351     
352    REAL,DIMENSION(klon_mpi) :: buffer_omp
353    REAL :: Field2d(iim,jj_nb)
354
355    IF (prt_level >= 10) WRITE(lunout,*)'Begin histrwrite2d_xios ',trim(field_name)
356
357    IF (SIZE(field)/=klon) CALL abort_gcm('iophy::histwrite2d_xios','Field first DIMENSION not equal to klon',1)
358   
359    CALL Gather_omp(field,buffer_omp)   
360!$OMP MASTER
361    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
362   
363    CALL wxios_write_2D(field_name, Field2d)
364!$OMP END MASTER   
365
366    IF (prt_level >= 10) WRITE(lunout,*)'End histrwrite2d_xios ',trim(field_name)
367  END SUBROUTINE histwrite2d_xios
368#endif
369
370!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
371
372! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
373#ifdef CPP_XIOS
374  SUBROUTINE histwrite3d_xios(field_name, field)
375  USE dimphy, only: klon, klev
376  USE mod_phys_lmdz_para, only: gather_omp, grid1Dto2D_mpi, &
377                                jj_nb, klon_mpi
378  USE wxios, only: wxios_write_3D
379
380
381  IMPLICIT NONE
382  INCLUDE 'dimensions.h'
383  INCLUDE 'iniprint.h'
384
385    CHARACTER(LEN=*), INTENT(IN) :: field_name
386    REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:)
387
388    REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp
389    REAL :: Field3d(iim,jj_nb,SIZE(field,2))
390    INTEGER :: ip, n, nlev
391
392  IF (prt_level >= 10) write(lunout,*)'Begin histrwrite3d_xios ',trim(field_name)
393
394    !Et on.... écrit
395    IF (SIZE(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first DIMENSION not equal to klon',1)
396    nlev=SIZE(field,2)
397
398
399    CALL Gather_omp(field,buffer_omp)
400!$OMP MASTER
401    CALL grid1Dto2D_mpi(buffer_omp,field3d)
402
403    CALL wxios_write_3D(field_name, Field3d(:,:,1:klev))
404!$OMP END MASTER   
405
406    IF (prt_level >= 10) write(lunout,*)'End histrwrite3d_xios ',trim(field_name)
407  END SUBROUTINE histwrite3d_xios
408#endif
409
410end module iophy
Note: See TracBrowser for help on using the repository browser.