source: trunk/libf/phylmd/iophy.F90 @ 1

Last change on this file since 1 was 1, checked in by emillour, 14 years ago

Import initial LMDZ5

File size: 5.1 KB
RevLine 
[1]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 
12  INTERFACE histwrite_phy
13    MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy
14  END INTERFACE
15
16
17contains
18
19  subroutine init_iophy_new(rlat,rlon)
20  USE dimphy
21  USE mod_phys_lmdz_para
22  USE mod_grid_phy_lmdz
23  USE ioipsl
24  implicit none
25  include 'dimensions.h'   
26    real,dimension(klon),intent(in) :: rlon
27    real,dimension(klon),intent(in) :: rlat
28
29    REAL,dimension(klon_glo)        :: rlat_glo
30    REAL,dimension(klon_glo)        :: rlon_glo
31   
32    INTEGER,DIMENSION(2) :: ddid
33    INTEGER,DIMENSION(2) :: dsg
34    INTEGER,DIMENSION(2) :: dsl
35    INTEGER,DIMENSION(2) :: dpf
36    INTEGER,DIMENSION(2) :: dpl
37    INTEGER,DIMENSION(2) :: dhs
38    INTEGER,DIMENSION(2) :: dhe
39    INTEGER :: i   
40
41    CALL gather(rlat,rlat_glo)
42    CALL bcast(rlat_glo)
43    CALL gather(rlon,rlon_glo)
44    CALL bcast(rlon_glo)
45   
46!$OMP MASTER 
47    ALLOCATE(io_lat(jjm+1-1/iim))
48    io_lat(1)=rlat_glo(1)
49    io_lat(jjm+1-1/iim)=rlat_glo(klon_glo)
50    IF (iim > 1) then
51      DO i=2,jjm
52        io_lat(i)=rlat_glo(2+(i-2)*iim)
53      ENDDO
54    ENDIF
55
56    ALLOCATE(io_lon(iim))
57    io_lon(:)=rlon_glo(2-1/iim:iim+1-1/iim)
58
59    ddid=(/ 1,2 /)
60    dsg=(/ iim, jjm+1-1/iim /)
61    dsl=(/ iim, jj_nb /)
62    dpf=(/ 1,jj_begin /)
63    dpl=(/ iim, jj_end /)
64    dhs=(/ ii_begin-1,0 /)
65    if (mpi_rank==mpi_size-1) then
66      dhe=(/0,0/)
67    else
68      dhe=(/ iim-ii_end,0 /) 
69    endif
70   
71    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
72                      'APPLE',phys_domain_id)
73
74!$OMP END MASTER
75     
76  end subroutine init_iophy_new
77
78  subroutine init_iophy(lat,lon)
79  USE dimphy
80  USE mod_phys_lmdz_para
81  use ioipsl
82  implicit none
83  include 'dimensions.h'   
84    real,dimension(iim),intent(in) :: lon
85    real,dimension(jjm+1-1/iim),intent(in) :: lat
86
87    INTEGER,DIMENSION(2) :: ddid
88    INTEGER,DIMENSION(2) :: dsg
89    INTEGER,DIMENSION(2) :: dsl
90    INTEGER,DIMENSION(2) :: dpf
91    INTEGER,DIMENSION(2) :: dpl
92    INTEGER,DIMENSION(2) :: dhs
93    INTEGER,DIMENSION(2) :: dhe
94
95!$OMP MASTER 
96    allocate(io_lat(jjm+1-1/iim))
97    io_lat(:)=lat(:)
98    allocate(io_lon(iim))
99    io_lon(:)=lon(:)
100   
101    ddid=(/ 1,2 /)
102    dsg=(/ iim, jjm+1-1/iim /)
103    dsl=(/ iim, jj_nb /)
104    dpf=(/ 1,jj_begin /)
105    dpl=(/ iim, jj_end /)
106    dhs=(/ ii_begin-1,0 /)
107    if (mpi_rank==mpi_size-1) then
108      dhe=(/0,0/)
109    else
110      dhe=(/ iim-ii_end,0 /) 
111    endif
112   
113    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
114                      'APPLE',phys_domain_id)
115
116!$OMP END MASTER
117     
118  end subroutine init_iophy
119 
120  subroutine histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
121  USE dimphy
122  USE mod_phys_lmdz_para
123  use ioipsl
124  use write_field
125  implicit none
126  include 'dimensions.h'
127   
128    character*(*), intent(IN) :: name
129    integer, intent(in) :: itau0
130    real,intent(in) :: zjulian
131    real,intent(in) :: dtime
132    integer,intent(out) :: nhori
133    integer,intent(out) :: nid_day
134
135!$OMP MASTER   
136    if (is_sequential) then
137      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
138                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
139    else
140      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
141                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
142    endif
143!$OMP END MASTER
144 
145  end subroutine histbeg_phy
146 
147  subroutine histwrite2d_phy(nid,name,itau,field)
148  USE dimphy
149  USE mod_phys_lmdz_para
150  USE ioipsl
151  implicit none
152  include 'dimensions.h'
153   
154    integer,intent(in) :: nid
155    character*(*), intent(IN) :: name
156    integer, intent(in) :: itau
157    real,dimension(:),intent(in) :: field
158    REAL,dimension(klon_mpi) :: buffer_omp
159    INTEGER :: index2d(iim*jj_nb)
160    REAL :: Field2d(iim,jj_nb)
161
162    IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first dimension not equal to klon',1)
163   
164    CALL Gather_omp(field,buffer_omp)   
165!$OMP MASTER
166    CALL grid1Dto2D_mpi(buffer_omp,Field2d)
167    CALL histwrite(nid,name,itau,Field2d,iim*jj_nb,index2d)
168!$OMP END MASTER   
169  end subroutine histwrite2d_phy
170
171
172 
173  subroutine histwrite3d_phy(nid,name,itau,field)
174  USE dimphy
175  USE mod_phys_lmdz_para
176
177  use ioipsl
178  implicit none
179  include 'dimensions.h'
180   
181    integer,intent(in) :: nid
182    character*(*), intent(IN) :: name
183    integer, intent(in) :: itau
184    real,dimension(:,:),intent(in) :: field  ! --> field(klon,:)
185    REAL,dimension(klon_mpi,size(field,2)) :: buffer_omp
186    INTEGER :: nlev
187    INTEGER :: index3d(iim*jj_nb*size(field,2))
188    REAL :: Field3d(iim,jj_nb,size(field,2))
189   
190    IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first dimension not equal to klon',1)
191    nlev=size(field,2)
192   
193    CALL Gather_omp(field,buffer_omp)
194!$OMP MASTER
195    CALL grid1Dto2D_mpi(buffer_omp,field3d)
196    CALL histwrite(nid,name,itau,Field3d,iim*jj_nb*nlev,index3d)
197!$OMP END MASTER   
198  end subroutine histwrite3d_phy
199 
200 
201
202end module iophy
Note: See TracBrowser for help on using the repository browser.