source: LMDZ4/branches/LMDZ4_V3_patches/libf/phylmd/iophy.F90 @ 5434

Last change on this file since 5434 was 793, checked in by Laurent Fairhead, 18 years ago

Modifications suite a la transformation des fichiers include pour
qu'ils soient compatibles a la fois au format fixe et au format libre
Un bon nombre de fichiers *.inc du coup disparaissent
LF

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