source: LMDZ5/branches/LMDZ5_AR5/libf/phylmd/iophy.F90 @ 5447

Last change on this file since 5447 was 1539, checked in by musat, 14 years ago

Ajouts CFMIP2/CMIP5

  • 6eme fichier de sortie "stations" histstn.nc qui necessite 2 fichiers (voir DefLists?): npCFMIP_param.data(_*) contenant le nombre de points (120 pour simulations AMIP, 73 pour aqua) pointlocations.txt(_*) contenant le numero, les coordonnees (lon,lat) et le nom de chaque station
  • flag LOGICAL dans tous les appels histwrite_phy pour pouvoir sortir le fichier histstn.nc

NB: 1) les flags de type phys_ que l'on met dans le physiq.def_L* pour ajouter plus de sorties

necessitent dorenavant 6 valeurs, la 6eme correspondant au fichier histstn.nc

2) par defaut le fichier histstn.nc ne sort pas; pour le sortir ajouter les lignes suivantes

dans physiq.def_L*

### Type de fichier : global (n) ou stations (y)
phys_out_filestations = n n n n n y

  • introduction de 2 jeux de flags pour les taux des GES; taux actuels avec suffixes _act, taux futurs avec "_per" avec 2 appels au rayonnement si taux "_per" different des taux "_act" (utiles pour diags. CFMIP 4CO2)
  • flags "betaCRF" pour calculs CRF pour experiences sensibilite proprietes optiques eau liquide nuageuse avec initialisations par defaut; sinon besoin de fichier beta_crf.data

Ajout flag LOGICAL lCOSP necessaire pour sortir un fichier stations
IM

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 12.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  INTERFACE histwrite_phy
15    MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy
16  END INTERFACE
17
18  INTERFACE histbeg_phy_all
19    MODULE PROCEDURE histbeg_phy,histbeg_phy_points
20  END INTERFACE
21
22
23contains
24
25  subroutine init_iophy_new(rlat,rlon)
26  USE dimphy
27  USE mod_phys_lmdz_para
28  USE mod_grid_phy_lmdz
29  USE ioipsl
30  implicit none
31  include 'dimensions.h'   
32    real,dimension(klon),intent(in) :: rlon
33    real,dimension(klon),intent(in) :: rlat
34
35    REAL,dimension(klon_glo)        :: rlat_glo
36    REAL,dimension(klon_glo)        :: rlon_glo
37   
38    INTEGER,DIMENSION(2) :: ddid
39    INTEGER,DIMENSION(2) :: dsg
40    INTEGER,DIMENSION(2) :: dsl
41    INTEGER,DIMENSION(2) :: dpf
42    INTEGER,DIMENSION(2) :: dpl
43    INTEGER,DIMENSION(2) :: dhs
44    INTEGER,DIMENSION(2) :: dhe
45    INTEGER :: i   
46
47    CALL gather(rlat,rlat_glo)
48    CALL bcast(rlat_glo)
49    CALL gather(rlon,rlon_glo)
50    CALL bcast(rlon_glo)
51   
52!$OMP MASTER 
53    ALLOCATE(io_lat(jjm+1-1/iim))
54    io_lat(1)=rlat_glo(1)
55    io_lat(jjm+1-1/iim)=rlat_glo(klon_glo)
56    IF (iim > 1) then
57      DO i=2,jjm
58        io_lat(i)=rlat_glo(2+(i-2)*iim)
59      ENDDO
60    ENDIF
61
62    ALLOCATE(io_lon(iim))
63    io_lon(:)=rlon_glo(2-1/iim:iim+1-1/iim)
64
65    ddid=(/ 1,2 /)
66    dsg=(/ iim, jjm+1-1/iim /)
67    dsl=(/ iim, jj_nb /)
68    dpf=(/ 1,jj_begin /)
69    dpl=(/ iim, jj_end /)
70    dhs=(/ ii_begin-1,0 /)
71    if (mpi_rank==mpi_size-1) then
72      dhe=(/0,0/)
73    else
74      dhe=(/ iim-ii_end,0 /) 
75    endif
76   
77    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
78                      'APPLE',phys_domain_id)
79
80!$OMP END MASTER
81     
82  end subroutine init_iophy_new
83
84  subroutine init_iophy(lat,lon)
85  USE dimphy
86  USE mod_phys_lmdz_para
87  use ioipsl
88  implicit none
89  include 'dimensions.h'   
90    real,dimension(iim),intent(in) :: lon
91    real,dimension(jjm+1-1/iim),intent(in) :: lat
92
93    INTEGER,DIMENSION(2) :: ddid
94    INTEGER,DIMENSION(2) :: dsg
95    INTEGER,DIMENSION(2) :: dsl
96    INTEGER,DIMENSION(2) :: dpf
97    INTEGER,DIMENSION(2) :: dpl
98    INTEGER,DIMENSION(2) :: dhs
99    INTEGER,DIMENSION(2) :: dhe
100
101!$OMP MASTER 
102    allocate(io_lat(jjm+1-1/iim))
103    io_lat(:)=lat(:)
104    allocate(io_lon(iim))
105    io_lon(:)=lon(:)
106   
107    ddid=(/ 1,2 /)
108    dsg=(/ iim, jjm+1-1/iim /)
109    dsl=(/ iim, jj_nb /)
110    dpf=(/ 1,jj_begin /)
111    dpl=(/ iim, jj_end /)
112    dhs=(/ ii_begin-1,0 /)
113    if (mpi_rank==mpi_size-1) then
114      dhe=(/0,0/)
115    else
116      dhe=(/ iim-ii_end,0 /) 
117    endif
118   
119    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
120                      'APPLE',phys_domain_id)
121
122!$OMP END MASTER
123     
124  end subroutine init_iophy
125 
126  subroutine histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
127  USE dimphy
128  USE mod_phys_lmdz_para
129  use ioipsl
130  use write_field
131  implicit none
132  include 'dimensions.h'
133   
134    character*(*), intent(IN) :: name
135    integer, intent(in) :: itau0
136    real,intent(in) :: zjulian
137    real,intent(in) :: dtime
138    integer,intent(out) :: nhori
139    integer,intent(out) :: nid_day
140
141!$OMP MASTER   
142    if (is_sequential) then
143      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
144                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
145    else
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,phys_domain_id)
148    endif
149!$OMP END MASTER
150 
151  end subroutine histbeg_phy
152
153  subroutine histbeg_phy_points(rlon,rlat,pim,tabij,ipt,jpt, &
154             plon,plat,plon_bounds,plat_bounds, &
155             nname,itau0,zjulian,dtime,nnhori,nnid_day)
156  USE dimphy
157  USE mod_phys_lmdz_para
158  USE mod_grid_phy_lmdz
159  use ioipsl
160  use write_field
161  implicit none
162  include 'dimensions.h'
163
164    real,dimension(klon),intent(in) :: rlon
165    real,dimension(klon),intent(in) :: rlat
166    integer, intent(in) :: itau0
167    real,intent(in) :: zjulian
168    real,intent(in) :: dtime
169    integer, intent(in) :: pim
170    integer, intent(out) :: nnhori
171    character(len=20), intent(in) :: nname
172    INTEGER, intent(out) :: nnid_day
173    integer :: i
174    REAL,dimension(klon_glo)        :: rlat_glo
175    REAL,dimension(klon_glo)        :: rlon_glo
176    INTEGER, DIMENSION(pim), intent(in)  :: tabij
177    REAL,dimension(pim), intent(in) :: plat, plon
178    INTEGER,dimension(pim), intent(in) :: ipt, jpt
179    REAL,dimension(pim,2), intent(out) :: plat_bounds, plon_bounds
180
181    INTEGER, SAVE :: tabprocbeg, tabprocend
182!$OMP THREADPRIVATE(tabprocbeg, tabprocend)
183    INTEGER :: ip
184    INTEGER, PARAMETER :: nip=1
185    INTEGER :: npproc
186    REAL, allocatable, dimension(:) :: npplat, npplon
187    REAL, allocatable, dimension(:,:) :: npplat_bounds, npplon_bounds
188    INTEGER, PARAMETER :: jjmp1=jjm+1-1/jjm
189    REAL, dimension(iim,jjmp1) :: zx_lon, zx_lat
190
191    CALL gather(rlat,rlat_glo)
192    CALL bcast(rlat_glo)
193    CALL gather(rlon,rlon_glo)
194    CALL bcast(rlon_glo)
195
196!$OMP MASTER
197    DO i=1,pim
198
199!    print*,'CFMIP_iophy i tabij lon lat',i,tabij(i),plon(i),plat(i)
200
201     plon_bounds(i,1)=rlon_glo(tabij(i)-1)
202     plon_bounds(i,2)=rlon_glo(tabij(i)+1)
203     if(plon_bounds(i,2).LE.0..AND.plon_bounds(i,1).GE.0.) THEN
204      if(rlon_glo(tabij(i)).GE.0.) THEN
205       plon_bounds(i,2)=-1*plon_bounds(i,2)
206      endif
207     endif
208     if(plon_bounds(i,2).GE.0..AND.plon_bounds(i,1).LE.0.) THEN
209      if(rlon_glo(tabij(i)).LE.0.) THEN
210       plon_bounds(i,2)=-1*plon_bounds(i,2)
211      endif
212     endif
213!
214     IF ( tabij(i).LE.iim) THEN
215      plat_bounds(i,1)=rlat_glo(tabij(i))
216     ELSE
217      plat_bounds(i,1)=rlat_glo(tabij(i)-iim)
218     ENDIF
219     plat_bounds(i,2)=rlat_glo(tabij(i)+iim)
220!
221!    print*,'CFMIP_iophy point i lon lon_bds',i,plon_bounds(i,1),rlon_glo(tabij(i)),plon_bounds(i,2)
222!    print*,'CFMIP_iophy point i lat lat_bds',i,plat_bounds(i,1),rlat_glo(tabij(i)),plat_bounds(i,2)
223!
224    ENDDO
225    if (is_sequential) then
226
227     npstn=pim
228     IF(.NOT. ALLOCATED(nptabij)) THEN
229      ALLOCATE(nptabij(pim))
230     ENDIF
231     DO i=1,pim
232      nptabij(i)=tabij(i)
233     ENDDO
234
235       CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlon_glo,zx_lon)
236       if (iim.gt.1) then
237       DO i = 1, iim
238         zx_lon(i,1) = rlon_glo(i+1)
239         zx_lon(i,jjmp1) = rlon_glo(i+1)
240       ENDDO
241       endif
242       CALL gr_fi_ecrit(1,klon,iim,jjmp1,rlat_glo,zx_lat)
243
244    DO i=1,pim
245!    print*,'CFMIP_iophy i tabij lon lat',i,tabij(i),plon(i),plat(i)
246
247     plon_bounds(i,1)=zx_lon(ipt(i)-1,jpt(i))
248     plon_bounds(i,2)=zx_lon(ipt(i)+1,jpt(i))
249
250     if (ipt(i).EQ.1) then
251      plon_bounds(i,1)=zx_lon(iim,jpt(i))
252      plon_bounds(i,2)=360.+zx_lon(ipt(i)+1,jpt(i))
253     endif
254 
255     if (ipt(i).EQ.iim) then
256      plon_bounds(i,2)=360.+zx_lon(1,jpt(i))
257     endif
258
259     plat_bounds(i,1)=zx_lat(ipt(i),jpt(i)-1)
260     plat_bounds(i,2)=zx_lat(ipt(i),jpt(i)+1)
261
262     if (jpt(i).EQ.1) then
263      plat_bounds(i,1)=zx_lat(ipt(i),1)+0.001
264      plat_bounds(i,2)=zx_lat(ipt(i),1)-0.001
265     endif
266 
267     if (jpt(i).EQ.jjmp1) then
268      plat_bounds(i,1)=zx_lat(ipt(i),jjmp1)+0.001
269      plat_bounds(i,2)=zx_lat(ipt(i),jjmp1)-0.001
270     endif
271!
272!    print*,'CFMIP_iophy point i lon lon_bds',i,plon_bounds(i,1),rlon(tabij(i)),plon_bounds(i,2)
273!    print*,'CFMIP_iophy point i lat lat_bds',i,plat_bounds(i,1),rlat(tabij(i)),plat_bounds(i,2)
274!
275    ENDDO
276!    print*,'iophy is_sequential nname, nnhori, nnid_day=',trim(nname),nnhori,nnid_day
277     call histbeg(nname,pim,plon,plon_bounds, &
278                           plat,plat_bounds, &
279                           itau0, zjulian, dtime, nnhori, nnid_day)
280    else
281     npproc=0
282     DO ip=1, pim
283      tabprocbeg=klon_mpi_begin
284      tabprocend=klon_mpi_end
285      IF(tabij(ip).GE.tabprocbeg.AND.tabij(ip).LE.tabprocend) THEN
286       npproc=npproc+1
287       npstn=npproc
288      ENDIF
289     ENDDO
290!    print*,'CFMIP_iophy mpi_rank npstn',mpi_rank,npstn
291     IF(.NOT. ALLOCATED(nptabij)) THEN
292      ALLOCATE(nptabij(npstn))
293      ALLOCATE(npplon(npstn), npplat(npstn))
294      ALLOCATE(npplon_bounds(npstn,2), npplat_bounds(npstn,2))
295     ENDIF
296     npproc=0
297     DO ip=1, pim
298      IF(tabij(ip).GE.tabprocbeg.AND.tabij(ip).LE.tabprocend) THEN
299       npproc=npproc+1
300       nptabij(npproc)=tabij(ip)
301!      print*,'mpi_rank npproc ip plon plat tabij=',mpi_rank,npproc,ip, &
302!      plon(ip),plat(ip),tabij(ip)
303       npplon(npproc)=plon(ip)
304       npplat(npproc)=plat(ip)
305       npplon_bounds(npproc,1)=plon_bounds(ip,1)
306       npplon_bounds(npproc,2)=plon_bounds(ip,2)
307       npplat_bounds(npproc,1)=plat_bounds(ip,1)
308       npplat_bounds(npproc,2)=plat_bounds(ip,2)
309!!!
310!!! print qui sert a reordonner les points stations selon l'ordre CFMIP
311!!! ne pas enlever
312        print*,'iophy_mpi rank ip lon lat',mpi_rank,ip,plon(ip),plat(ip)
313!!!
314      ENDIF
315     ENDDO
316     call histbeg(nname,npstn,npplon,npplon_bounds, &
317                            npplat,npplat_bounds, &
318                            itau0,zjulian,dtime,nnhori,nnid_day,phys_domain_id)
319    endif
320!$OMP END MASTER
321
322  end subroutine histbeg_phy_points
323 
324  subroutine histwrite2d_phy(nid,lpoint,name,itau,field)
325  USE dimphy
326  USE mod_phys_lmdz_para
327  USE ioipsl
328  implicit none
329  include 'dimensions.h'
330   
331    integer,intent(in) :: nid
332    logical,intent(in) :: lpoint
333    character*(*), intent(IN) :: name
334    integer, intent(in) :: itau
335    real,dimension(:),intent(in) :: field
336    REAL,dimension(klon_mpi) :: buffer_omp
337    INTEGER, allocatable, dimension(:) :: index2d
338    REAL :: Field2d(iim,jj_nb)
339
340    integer :: ip
341    real,allocatable,dimension(:) :: fieldok
342
343    IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first dimension not equal to klon',1)
344   
345    CALL Gather_omp(field,buffer_omp)   
346!$OMP MASTER
347    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
348    if(.NOT.lpoint) THEN
349     ALLOCATE(index2d(iim*jj_nb))
350     ALLOCATE(fieldok(iim*jj_nb))
351     CALL histwrite(nid,name,itau,Field2d,iim*jj_nb,index2d)
352    else
353     ALLOCATE(fieldok(npstn))
354     ALLOCATE(index2d(npstn))
355
356     if(is_sequential) then
357!     klon_mpi_begin=1
358!     klon_mpi_end=klon
359      DO ip=1, npstn
360       fieldok(ip)=buffer_omp(nptabij(ip))
361      ENDDO
362     else
363      DO ip=1, npstn
364!     print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip)
365       IF(nptabij(ip).GE.klon_mpi_begin.AND. &
366          nptabij(ip).LE.klon_mpi_end) THEN
367         fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1)
368       ENDIF
369      ENDDO
370     endif
371     CALL histwrite(nid,name,itau,fieldok,npstn,index2d)
372!
373    endif
374    deallocate(index2d)
375    deallocate(fieldok)
376!$OMP END MASTER   
377  end subroutine histwrite2d_phy
378
379  subroutine histwrite3d_phy(nid,lpoint,name,itau,field)
380  USE dimphy
381  USE mod_phys_lmdz_para
382
383  use ioipsl
384  implicit none
385  include 'dimensions.h'
386   
387    integer,intent(in) :: nid
388    logical,intent(in) :: lpoint
389    character*(*), intent(IN) :: name
390    integer, intent(in) :: itau
391    real,dimension(:,:),intent(in) :: field  ! --> field(klon,:)
392    REAL,dimension(klon_mpi,size(field,2)) :: buffer_omp
393    REAL :: Field3d(iim,jj_nb,size(field,2))
394    INTEGER :: ip, n, nlev
395    INTEGER, ALLOCATABLE, dimension(:) :: index3d
396    real,allocatable, dimension(:,:) :: fieldok
397
398    IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first dimension not equal to klon',1)
399    nlev=size(field,2)
400
401!   print*,'hist3d_phy mpi_rank npstn=',mpi_rank,npstn
402
403!   DO ip=1, npstn
404!    print*,'hist3d_phy mpi_rank nptabij',mpi_rank,nptabij(ip)
405!   ENDDO
406
407    CALL Gather_omp(field,buffer_omp)
408!$OMP MASTER
409    CALL grid1Dto2D_mpi(buffer_omp,field3d)
410    if(.NOT.lpoint) THEN
411     ALLOCATE(index3d(iim*jj_nb*nlev))
412     ALLOCATE(fieldok(iim*jj_nb,nlev))
413     CALL histwrite(nid,name,itau,Field3d,iim*jj_nb*nlev,index3d)
414    else
415      nlev=size(field,2)
416      ALLOCATE(index3d(npstn*nlev))
417      ALLOCATE(fieldok(npstn,nlev))
418
419      if(is_sequential) then
420!      klon_mpi_begin=1
421!      klon_mpi_end=klon
422       DO n=1, nlev
423       DO ip=1, npstn
424        fieldok(ip,n)=buffer_omp(nptabij(ip),n)
425       ENDDO
426       ENDDO
427      else
428       DO n=1, nlev
429       DO ip=1, npstn
430        IF(nptabij(ip).GE.klon_mpi_begin.AND. &
431         nptabij(ip).LE.klon_mpi_end) THEN
432         fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n)
433        ENDIF
434       ENDDO
435       ENDDO
436      endif
437      CALL histwrite(nid,name,itau,fieldok,npstn*nlev,index3d)
438    endif
439  deallocate(index3d)
440  deallocate(fieldok)
441!$OMP END MASTER   
442  end subroutine histwrite3d_phy
443 
444end module iophy
Note: See TracBrowser for help on using the repository browser.