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