source: LMDZ4/trunk/libf/phylmd/iophy.F90 @ 931

Last change on this file since 931 was 931, checked in by lmdzadmin, 17 years ago

Oubli mise a jour iophy.F90 AI
IM

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 4.8 KB
Line 
1!
2! $Header$
3!
4module iophy
5 
6  REAL,private,allocatable,dimension(:,:),save :: tmp_tab2d
7  REAL,private,allocatable,dimension(:,:,:),save :: tmp_tab3d
8  INTEGER,private,allocatable,dimension(:),save :: ndex2d
9  INTEGER,private,allocatable,dimension(:),save :: ndex3d
10! abd  REAL,private,allocatable,dimension(:),save :: io_lat
11! abd  REAL,private,allocatable,dimension(:),save :: io_lon
12  REAL,allocatable,dimension(:),save :: io_lat
13  REAL,allocatable,dimension(:),save :: io_lon
14  INTEGER, save :: phys_domain_id
15 
16  INTERFACE histwrite_phy
17    MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy
18  END INTERFACE
19
20contains
21
22  subroutine init_iophy(lat,lon)
23  USE dimphy
24  USE mod_phys_lmdz_para
25  use ioipsl
26  implicit none
27  include 'dimensions.h'   
28    real,dimension(iim),intent(in) :: lon
29    real,dimension(jjm+1-1/iim),intent(in) :: lat
30
31    INTEGER,DIMENSION(2) :: ddid
32    INTEGER,DIMENSION(2) :: dsg
33    INTEGER,DIMENSION(2) :: dsl
34    INTEGER,DIMENSION(2) :: dpf
35    INTEGER,DIMENSION(2) :: dpl
36    INTEGER,DIMENSION(2) :: dhs
37    INTEGER,DIMENSION(2) :: dhe
38
39!$OMP MASTER 
40    allocate(io_lat(jjm+1-1/iim))
41    io_lat(:)=lat(:)
42    allocate(io_lon(iim))
43    io_lon(:)=lon(:)
44    allocate(tmp_tab2d(iim,jj_nb))
45    allocate(tmp_tab3d(iim,jj_nb,klev))
46    allocate(ndex2d(iim*jj_nb))
47    allocate(ndex3d(iim*jj_nb*klev))
48    ndex2d(:)=0
49    ndex3d(:)=0
50   
51    ddid=(/ 1,2 /)
52    dsg=(/ iim, jjm+1-1/iim /)
53    dsl=(/ iim, jj_nb /)
54    dpf=(/ 1,jj_begin /)
55    dpl=(/ iim, jj_end /)
56    dhs=(/ ii_begin-1,0 /)
57    if (mpi_rank==mpi_size-1) then
58      dhe=(/0,0/)
59    else
60      dhe=(/ iim-ii_end,0 /) 
61    endif
62   
63    call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
64                      'APPLE',phys_domain_id)
65
66!$OMP END MASTER
67     
68  end subroutine init_iophy
69 
70  subroutine histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day)
71  USE dimphy
72  USE mod_phys_lmdz_para
73  use ioipsl
74  use write_field
75  implicit none
76  include 'dimensions.h'
77   
78    character*(*), intent(IN) :: name
79    integer, intent(in) :: itau0
80    real,intent(in) :: zjulian
81    real,intent(in) :: dtime
82    integer,intent(out) :: nhori
83    integer,intent(out) :: nid_day
84
85!$OMP MASTER   
86    if (is_sequential) then
87      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
88                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day)
89    else
90      call histbeg(name,iim,io_lon, jj_nb,io_lat(jj_begin:jj_end), &
91                   1,iim,1,jj_nb,itau0, zjulian, dtime, nhori, nid_day,phys_domain_id)
92    endif
93!$OMP END MASTER
94 
95  end subroutine histbeg_phy
96 
97  subroutine histwrite2d_phy(nid,name,itau,field)
98  USE dimphy
99  USE mod_phys_lmdz_para
100  USE ioipsl
101  implicit none
102  include 'dimensions.h'
103   
104    integer,intent(in) :: nid
105    character*(*), intent(IN) :: name
106    integer, intent(in) :: itau
107    real,dimension(klon),intent(in) :: field
108   
109    REAL,dimension(klon_mpi) :: buffer_omp
110   
111    CALL Gather_omp(field,buffer_omp)   
112!$OMP MASTER
113    CALL grid1Dto2D_mpi(buffer_omp,tmp_tab2d)
114    CALL histwrite(nid,name,itau,tmp_tab2d,iim*jj_nb,ndex2d)
115!$OMP END MASTER   
116  end subroutine histwrite2d_phy
117 
118  subroutine histwrite3d_phy(nid,name,itau,field)
119  USE dimphy
120  USE mod_phys_lmdz_para
121
122  use ioipsl
123  implicit none
124  include 'dimensions.h'
125   
126    integer,intent(in) :: nid
127    character*(*), intent(IN) :: name
128    integer, intent(in) :: itau
129    real,dimension(klon,klev),intent(in) :: field
130
131    REAL,dimension(klon_mpi,klev) :: buffer_omp
132   
133    CALL Gather_omp(field,buffer_omp)
134!$OMP MASTER
135    CALL grid1Dto2D_mpi(buffer_omp,tmp_tab3d)
136    CALL histwrite(nid,name,itau,tmp_tab3d,iim*jj_nb*klev,ndex3d)
137!$OMP END MASTER   
138  end subroutine histwrite3d_phy
139 
140 
141!  subroutine phy2dyn(field_phy,field_dyn,nlev)
142!  USE dimphy_old
143!  implicit none
144!  include 'dimensions.h'
145
146!    real,dimension(klon_mpi,nlev),intent(in) :: field_phy
147!    real,dimension(iim,jjphy_nb,nlev),intent(out) :: field_dyn
148!    integer,intent(in) :: nlev
149!   
150!    integer :: next
151!    integer :: j,l
152!   
153!      do l=1,nlev
154!               
155!       if (jjphy_begin==jjphy_end) then
156!         field_dyn(:,1,l)=0.
157!         field_dyn(iiphy_begin:iiphy_end,1,l)=field_phy(1:klon_mpi,l)
158!       else
159!       
160!        if (jjphy_begin==1) then
161!           field_dyn(:,1,l)=field_phy(1,l)
162!           next=2
163!        else
164!          field_dyn(:,1,l)=0.
165!          next=iim-iiphy_begin+2
166!          field_dyn(iiphy_begin:iim,1,l)=field_phy(1:next-1,l)   
167!        endif
168!       
169!         do j=2,jjphy_nb-1
170!           field_dyn(:,j,l)=field_phy(next:next+iim-1,l)
171!           next=next+iim
172!         enddo
173!         
174!         if (jjphy_end==jjm+1-1/iim) then
175!             field_dyn(:,jjphy_nb,l)=field_phy(klon_mpi,l)
176!         else
177!          field_dyn(:,jjphy_nb,l)=0.
178!          field_dyn(1:iiphy_end,jjphy_nb,l)=field_phy(next:next+iiphy_end-1,l)   
179!         endif
180!         
181!       endif
182!     
183!     enddo
184!       
185!    end subroutine phy2dyn         
186 
187         
188end module iophy
Note: See TracBrowser for help on using the repository browser.