source: LMDZ6/branches/Amaury_dev/libf/phydev/iophy.F90 @ 5224

Last change on this file since 5224 was 5135, checked in by abarral, 4 months ago

Put iotd* into lmdz_iotd.f90

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:keywords set to Id
File size: 13.2 KB
Line 
1! $Id: iophy.F90 5135 2024-07-26 17:20:23Z abarral $
2
3module iophy
4
5  ! abd  REAL,PRIVATE,ALLOCATABLE,DIMENSION(:),save :: io_lat
6  ! abd  REAL,PRIVATE,ALLOCATABLE,DIMENSION(:),save :: io_lon
7  REAL, ALLOCATABLE, DIMENSION(:), save :: io_lat
8  REAL, ALLOCATABLE, DIMENSION(:), save :: io_lon
9  INTEGER, save :: phys_domain_id
10  INTEGER, save :: npstn
11  INTEGER, ALLOCATABLE, DIMENSION(:), save :: nptabij
12
13
14  ! interfaces for both IOIPSL and XIOS
15  INTERFACE histwrite_phy
16    MODULE PROCEDURE histwrite2d_phy, histwrite3d_phy, histwrite2d_xios, histwrite3d_xios
17  END INTERFACE
18
19  ! interfaces for both IOIPSL and XIOS
20  INTERFACE histbeg_phy_all
21    MODULE PROCEDURE histbeg_phy, histbeg_phyxios
22  END INTERFACE
23
24CONTAINS
25
26  SUBROUTINE init_iophy_new(rlat, rlon)
27    USE dimphy, ONLY: klon
28    USE lmdz_phys_para, ONLY: gather, bcast, &
29            jj_nb, jj_begin, jj_end, ii_begin, ii_end, &
30            mpi_size, mpi_rank, klon_mpi, &
31            is_sequential, is_south_pole_dyn
32    USE lmdz_grid_phy, ONLY: nbp_lon, nbp_lat, klon_glo
33    USE lmdz_print_control, ONLY: lunout, prt_level
34    USE lmdz_grid_phy, ONLY: nbp_lon, nbp_lat
35    USE ioipsl, ONLY: flio_dom_set
36    USE lmdz_wxios, ONLY: wxios_domain_param, using_xios
37    IMPLICIT NONE
38    REAL, DIMENSION(klon), INTENT(IN) :: rlon
39    REAL, DIMENSION(klon), INTENT(IN) :: rlat
40
41    REAL, DIMENSION(klon_glo) :: rlat_glo
42    REAL, DIMENSION(klon_glo) :: rlon_glo
43
44    INTEGER, DIMENSION(2) :: ddid
45    INTEGER, DIMENSION(2) :: dsg
46    INTEGER, DIMENSION(2) :: dsl
47    INTEGER, DIMENSION(2) :: dpf
48    INTEGER, DIMENSION(2) :: dpl
49    INTEGER, DIMENSION(2) :: dhs
50    INTEGER, DIMENSION(2) :: dhe
51    INTEGER :: i
52    INTEGER :: data_ibegin, data_iend
53
54    CALL gather(rlat, rlat_glo)
55    CALL bcast(rlat_glo)
56    CALL gather(rlon, rlon_glo)
57    CALL bcast(rlon_glo)
58
59    !$OMP MASTER
60    ALLOCATE(io_lat(nbp_lat))
61    io_lat(1) = rlat_glo(1)
62    io_lat(nbp_lat) = rlat_glo(klon_glo)
63    IF ((nbp_lon * nbp_lat) > 1) THEN
64      DO i = 2, nbp_lat - 1
65        io_lat(i) = rlat_glo(2 + (i - 2) * nbp_lon)
66      ENDDO
67    ENDIF
68
69    ALLOCATE(io_lon(nbp_lon))
70    IF ((nbp_lon * nbp_lat) > 1) THEN
71      io_lon(:) = rlon_glo(2:nbp_lon + 1)
72    ELSE
73      io_lon(1) = rlon_glo(1)
74    ENDIF
75    !! (I) dtnb   : total number of domains
76    !! (I) dnb    : domain number
77    !! (I) did(:) : distributed dimensions identifiers
78    !!              (up to 5 dimensions are supported)
79    !! (I) dsg(:) : total number of points for each dimension
80    !! (I) dsl(:) : local number of points for each dimension
81    !! (I) dpf(:) : position of first local point for each dimension
82    !! (I) dpl(:) : position of last local point for each dimension
83    !! (I) dhs(:) : start halo size for each dimension
84    !! (I) dhe(:) : end halo size for each dimension
85    !! (C) cdnm   : Model domain definition name.
86    !!              The names actually supported are :
87    !!              "BOX", "APPLE", "ORANGE".
88    !!              These names are case insensitive.
89    ddid = (/ 1, 2 /)
90    dsg = (/ nbp_lon, nbp_lat /)
91    dsl = (/ nbp_lon, jj_nb /)
92    dpf = (/ 1, jj_begin /)
93    dpl = (/ nbp_lon, jj_end /)
94    dhs = (/ ii_begin - 1, 0 /)
95    IF (mpi_rank==mpi_size - 1) THEN
96      dhe = (/0, 0/)
97    else
98      dhe = (/ nbp_lon - ii_end, 0 /)
99    endif
100
101#ifndef CPP_IOIPSL_NO_OUTPUT
102    CALL flio_dom_set(mpi_size, mpi_rank, ddid, dsg, dsl, dpf, dpl, dhs, dhe, &
103            'APPLE', phys_domain_id)
104#endif
105    IF (using_xios) THEN
106      ! Set values for the mask:
107      IF (mpi_rank == 0) THEN
108        data_ibegin = 0
109      ELSE
110        data_ibegin = ii_begin - 1
111      END IF
112
113      IF (mpi_rank == mpi_size - 1) THEN
114        data_iend = nbp_lon
115      ELSE
116        data_iend = ii_end + 1
117      END IF
118
119      IF (prt_level>=10) THEN
120        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
121        WRITE(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " nbp_lon=", nbp_lon, " nbp_lat=", nbp_lat
122        WRITE(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " data_ibegin=", data_ibegin, " data_iend=", data_iend
123        WRITE(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " data_ibegin=", data_ibegin, " data_iend=", data_iend
124        WRITE(lunout, *) "init_iophy_new: mpirank=", mpi_rank, " is_south_pole=", is_south_pole_dyn
125      endif
126
127      ! Initialize the XIOS domain coreesponding to this process:
128      CALL wxios_domain_param("dom_glo")
129      !      CALL wxios_domain_param("dom_glo", is_sequential, nbp_lon, jj_nb, nbp_lon, nbp_lat, &
130      !                              1, nbp_lon, ii_begin, ii_end, jj_begin, jj_end,             &
131      !                              klon_mpi+2*(nbp_lon-1), data_ibegin, data_iend,             &
132      !                              io_lat, io_lon,is_south_pole_dyn,mpi_rank)
133    ENDIF
134    !$OMP END MASTER
135
136  END SUBROUTINE init_iophy_new
137
138  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
139
140  SUBROUTINE histbeg_phy(name, itau0, zjulian, dtime, nhori, nid_day)
141    USE lmdz_phys_para, ONLY: is_sequential, jj_begin, jj_end, jj_nb
142    USE ioipsl, ONLY: histbeg
143    USE lmdz_print_control, ONLY: prt_level, lunout
144    USE lmdz_grid_phy, ONLY: nbp_lon
145    IMPLICIT NONE
146
147    CHARACTER*(*), INTENT(IN) :: name
148    INTEGER, INTENT(IN) :: itau0
149    REAL, INTENT(IN) :: zjulian
150    REAL, INTENT(IN) :: dtime
151    INTEGER, INTENT(OUT) :: nhori
152    INTEGER, INTENT(OUT) :: nid_day
153
154    !$OMP MASTER
155    IF (is_sequential) THEN
156      CALL histbeg(name, nbp_lon, io_lon, jj_nb, io_lat(jj_begin:jj_end), &
157              1, nbp_lon, 1, jj_nb, itau0, zjulian, dtime, nhori, nid_day)
158    else
159      CALL histbeg(name, nbp_lon, io_lon, jj_nb, io_lat(jj_begin:jj_end), &
160              1, nbp_lon, 1, jj_nb, itau0, zjulian, dtime, nhori, nid_day, phys_domain_id)
161    endif
162    !$OMP END MASTER
163
164  END SUBROUTINE  histbeg_phy
165
166  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
167
168
169  ! SUBROUTINE histbeg_phyxios(name,itau0,zjulian,dtime,ffreq,lev,nhori,nid_day)
170  SUBROUTINE histbeg_phyxios(name, ffreq, lev)
171    USE lmdz_phys_para, ONLY: is_using_mpi, is_mpi_root
172    USE lmdz_wxios, ONLY: wxios_add_file
173    IMPLICIT NONE
174
175    CHARACTER*(*), INTENT(IN) :: name
176    !    INTEGER, INTENT(IN) :: itau0
177    !    REAL,INTENT(IN) :: zjulian
178    !    REAL,INTENT(IN) :: dtime
179    character(LEN = *), INTENT(IN) :: ffreq
180    INTEGER, INTENT(IN) :: lev
181    !    INTEGER,INTENT(OUT) :: nhori
182    !    INTEGER,INTENT(OUT) :: nid_day
183
184    !$OMP MASTER
185
186    ! ug OMP en chantier...
187    IF((.NOT. is_using_mpi) .OR. is_mpi_root) THEN
188      ! ug Création du fichier
189      CALL wxios_add_file(name, ffreq, lev)
190    END IF
191
192    !$OMP END MASTER
193
194  END SUBROUTINE histbeg_phyxios
195
196
197  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
198
199  SUBROUTINE histwrite2d_phy(nid, lpoint, name, itau, field)
200    USE dimphy, ONLY: klon
201    USE lmdz_phys_para, ONLY: Gather_omp, grid1Dto2D_mpi, &
202            is_sequential, klon_mpi_begin, klon_mpi_end, &
203            jj_nb, klon_mpi
204    USE ioipsl, ONLY: histwrite
205    USE lmdz_grid_phy, ONLY: nbp_lon
206    USE lmdz_abort_physic, ONLY: abort_physic
207    IMPLICIT NONE
208
209    INTEGER, INTENT(IN) :: nid
210    logical, INTENT(IN) :: lpoint
211    CHARACTER*(*), INTENT(IN) :: name
212    INTEGER, INTENT(IN) :: itau
213    REAL, DIMENSION(:), INTENT(IN) :: field
214    REAL, DIMENSION(klon_mpi) :: buffer_omp
215    INTEGER, ALLOCATABLE, DIMENSION(:) :: index2d
216    REAL :: Field2d(nbp_lon, jj_nb)
217
218    INTEGER :: ip
219    REAL, ALLOCATABLE, DIMENSION(:) :: fieldok
220
221    IF (size(field)/=klon) CALL abort_physic('iophy::histwrite2d', 'Field first dimension not equal to klon', 1)
222
223    CALL Gather_omp(field, buffer_omp)
224    !$OMP MASTER
225    CALL grid1Dto2D_mpi(buffer_omp, Field2d)
226    IF(.NOT.lpoint) THEN
227      ALLOCATE(index2d(nbp_lon * jj_nb))
228      ALLOCATE(fieldok(nbp_lon * jj_nb))
229      CALL histwrite(nid, name, itau, Field2d, nbp_lon * jj_nb, index2d)
230    else
231      ALLOCATE(fieldok(npstn))
232      ALLOCATE(index2d(npstn))
233
234      IF(is_sequential) THEN
235        !     klon_mpi_begin=1
236        !     klon_mpi_end=klon
237        DO ip = 1, npstn
238          fieldok(ip) = buffer_omp(nptabij(ip))
239        ENDDO
240      else
241        DO ip = 1, npstn
242          !     PRINT*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
243          IF(nptabij(ip)>=klon_mpi_begin.AND. &
244                  nptabij(ip)<=klon_mpi_end) THEN
245            fieldok(ip) = buffer_omp(nptabij(ip) - klon_mpi_begin + 1)
246          ENDIF
247        ENDDO
248      endif
249      CALL histwrite(nid, name, itau, fieldok, npstn, index2d)
250
251    endif
252    deallocate(index2d)
253    deallocate(fieldok)
254    !$OMP END MASTER
255  END SUBROUTINE  histwrite2d_phy
256
257  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
258
259  SUBROUTINE histwrite3d_phy(nid, lpoint, name, itau, field)
260    USE dimphy, ONLY: klon
261    USE lmdz_phys_para, ONLY: Gather_omp, grid1Dto2D_mpi, &
262            is_sequential, klon_mpi_begin, klon_mpi_end, &
263            jj_nb, klon_mpi
264    USE ioipsl, ONLY: histwrite
265    USE lmdz_grid_phy, ONLY: nbp_lon
266    USE lmdz_abort_physic, ONLY: abort_physic
267    IMPLICIT NONE
268
269    INTEGER, INTENT(IN) :: nid
270    logical, INTENT(IN) :: lpoint
271    CHARACTER*(*), INTENT(IN) :: name
272    INTEGER, INTENT(IN) :: itau
273    REAL, DIMENSION(:, :), INTENT(IN) :: field  ! --> field(klon,:)
274    REAL, DIMENSION(klon_mpi, size(field, 2)) :: buffer_omp
275    REAL :: Field3d(nbp_lon, jj_nb, size(field, 2))
276    INTEGER :: ip, n, nlev
277    INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d
278    REAL, ALLOCATABLE, DIMENSION(:, :) :: fieldok
279
280    IF (size(field, 1)/=klon) CALL abort_physic('iophy::histwrite3d', 'Field first dimension not equal to klon', 1)
281    nlev = size(field, 2)
282
283    CALL Gather_omp(field, buffer_omp)
284    !$OMP MASTER
285    CALL grid1Dto2D_mpi(buffer_omp, field3d)
286    IF(.NOT.lpoint) THEN
287      ALLOCATE(index3d(nbp_lon * jj_nb * nlev))
288      ALLOCATE(fieldok(nbp_lon * jj_nb, nlev))
289      CALL histwrite(nid, name, itau, Field3d, nbp_lon * jj_nb * nlev, index3d)
290    else
291      nlev = size(field, 2)
292      ALLOCATE(index3d(npstn * nlev))
293      ALLOCATE(fieldok(npstn, nlev))
294
295      IF(is_sequential) THEN
296        !      klon_mpi_begin=1
297        !      klon_mpi_end=klon
298        DO n = 1, nlev
299          DO ip = 1, npstn
300            fieldok(ip, n) = buffer_omp(nptabij(ip), n)
301          ENDDO
302        ENDDO
303      else
304        DO n = 1, nlev
305          DO ip = 1, npstn
306            IF(nptabij(ip)>=klon_mpi_begin.AND. &
307                    nptabij(ip)<=klon_mpi_end) THEN
308              fieldok(ip, n) = buffer_omp(nptabij(ip) - klon_mpi_begin + 1, n)
309            ENDIF
310          ENDDO
311        ENDDO
312      endif
313      CALL histwrite(nid, name, itau, fieldok, npstn * nlev, index3d)
314    endif
315    deallocate(index3d)
316    deallocate(fieldok)
317    !$OMP END MASTER
318  END SUBROUTINE  histwrite3d_phy
319
320  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
321
322  ! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
323
324  SUBROUTINE histwrite2d_xios(field_name, field)
325    USE dimphy, ONLY: klon
326    USE lmdz_phys_para, ONLY: gather_omp, grid1Dto2D_mpi, &
327            jj_nb, klon_mpi
328    USE lmdz_xios, ONLY: xios_send_field
329    USE lmdz_print_control, ONLY: prt_level, lunout
330    USE lmdz_grid_phy, ONLY: nbp_lon
331    USE lmdz_abort_physic, ONLY: abort_physic
332    IMPLICIT NONE
333
334    CHARACTER(LEN = *), INTENT(IN) :: field_name
335    REAL, DIMENSION(:), INTENT(IN) :: field
336
337    REAL, DIMENSION(klon_mpi) :: buffer_omp
338    REAL :: Field2d(nbp_lon, jj_nb)
339
340    IF (prt_level >= 10) WRITE(lunout, *)'Begin histrwrite2d_xios ', trim(field_name)
341
342    IF (SIZE(field)/=klon) CALL abort_physic('iophy::histwrite2d_xios', 'Field first DIMENSION not equal to klon', 1)
343
344    CALL Gather_omp(field, buffer_omp)
345    !$OMP MASTER
346    CALL grid1Dto2D_mpi(buffer_omp, Field2d)
347
348    CALL xios_send_field(field_name, Field2d)
349    !$OMP END MASTER
350
351    IF (prt_level >= 10) WRITE(lunout, *)'End histrwrite2d_xios ', trim(field_name)
352  END SUBROUTINE histwrite2d_xios
353
354  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
355
356  ! VERSION DES HISTWRITE DEDIEES AU TOUT-XIOS-XML DEJA UTILISEE DANS PHYDEV
357
358  SUBROUTINE histwrite3d_xios(field_name, field)
359    USE dimphy, ONLY: klon, klev
360    USE lmdz_phys_para, ONLY: gather_omp, grid1Dto2D_mpi, &
361            jj_nb, klon_mpi
362    USE lmdz_xios, ONLY: xios_send_field
363    USE lmdz_print_control, ONLY: prt_level, lunout
364    USE lmdz_grid_phy, ONLY: nbp_lon
365    USE lmdz_abort_physic, ONLY: abort_physic
366    IMPLICIT NONE
367
368    CHARACTER(LEN = *), INTENT(IN) :: field_name
369    REAL, DIMENSION(:, :), INTENT(IN) :: field ! --> field(klon,:)
370
371    REAL, DIMENSION(klon_mpi, SIZE(field, 2)) :: buffer_omp
372    REAL :: Field3d(nbp_lon, jj_nb, SIZE(field, 2))
373    INTEGER :: ip, n, nlev
374
375    IF (prt_level >= 10) WRITE(lunout, *)'Begin histrwrite3d_xios ', trim(field_name)
376
377    !Et on.... écrit
378    IF (SIZE(field, 1)/=klon) CALL abort_physic('iophy::histwrite3d', 'Field first DIMENSION not equal to klon', 1)
379    nlev = SIZE(field, 2)
380
381    CALL Gather_omp(field, buffer_omp)
382    !$OMP MASTER
383    CALL grid1Dto2D_mpi(buffer_omp, field3d)
384
385    CALL xios_send_field(field_name, Field3d(:, :, 1:nlev))
386    !$OMP END MASTER
387
388    IF (prt_level >= 10) WRITE(lunout, *)'End histrwrite3d_xios ', trim(field_name)
389  END SUBROUTINE histwrite3d_xios
390
391END MODULE iophy
Note: See TracBrowser for help on using the repository browser.