updated top-level README and version_decl for V4.4.2 (#1795)
[WRF.git] / external / io_pio / read_bdy_routines.F90
blobbd8c1729fdbc5a92d76f3fa93df800908737c050
1 !------------------------------------------------------------------
2 !$Id$
3 !------------------------------------------------------------------
5 subroutine transRg2l(ds1,de1,ds2,de2,ds3,de3, &
6                      ms1,me1,ms2,me2,ms3,me3, &
7                      ps1,pe1,ps2,pe2,ps3,pe3, &
8                      dlocal, dglobal)
9   integer           ,intent(in)  :: ds1,de1,ds2,de2,ds3,de3
10   integer           ,intent(in)  :: ms1,me1,ms2,me2,ms3,me3
11   integer           ,intent(in)  :: ps1,pe1,ps2,pe2,ps3,pe3
12   real              ,intent(out) ::  dlocal(ms1:me1,ms2:me2,ms3:me3)
13   real              ,intent(in)  :: dglobal(ds1:de1,ds2:de2,ds3:de3)
14   integer                        :: i,j,k
16 !$OMP PARALLEL DO SCHEDULE(RUNTIME) PRIVATE(i,j,k)
17   do k=ps3,pe3
18     do j=ps2,pe2
19       do i=ps1,pe1
20           dlocal(i,j,k) = dglobal(i,j,k)
21       enddo
22     enddo
23   enddo
24 !$OMP END PARALLEL DO
25   return
26 end subroutine transRg2l
28 subroutine transDg2l(ds1,de1,ds2,de2,ds3,de3, &
29                      ms1,me1,ms2,me2,ms3,me3, &
30                      ps1,pe1,ps2,pe2,ps3,pe3, &
31                      dlocal, dglobal)
32   integer           ,intent(in)  :: ds1,de1,ds2,de2,ds3,de3
33   integer           ,intent(in)  :: ms1,me1,ms2,me2,ms3,me3
34   integer           ,intent(in)  :: ps1,pe1,ps2,pe2,ps3,pe3
35   real*8            ,intent(out) ::  dlocal(ms1:me1,ms2:me2,ms3:me3)
36   real*8            ,intent(in)  :: dglobal(ds1:de1,ds2:de2,ds3:de3)
37   integer                        :: i,j,k
39 !$OMP PARALLEL DO SCHEDULE(RUNTIME) PRIVATE(i,j,k)
40   do k=ps3,pe3
41     do j=ps2,pe2
42       do i=ps1,pe1
43           dlocal(i,j,k) = dglobal(i,j,k)
44       enddo
45     enddo
46   enddo
47 !$OMP END PARALLEL DO
48   return
49 end subroutine transDg2l
51 subroutine read_bdy_RealFieldIO(DH,NDim,Dimens,MemoryStart,MemoryEnd, &
52                                 PatchStart,PatchEnd,Data,Status)
53   use pio
54   use pio_kinds
55   use wrf_data_pio
56   use pio_routines
57   implicit none
58   include 'wrf_status_codes.h'
59   type(wrf_data_handle)                    :: DH
60   integer                     ,intent(in)  :: NDim
61   integer, dimension(*)       ,intent(in)  :: Dimens
62   integer, dimension(*)       ,intent(in)  :: MemoryStart, MemoryEnd
63   integer, dimension(*)       ,intent(in)  :: PatchStart,  PatchEnd
64   real, dimension(*)          ,intent(out) :: Data
65   integer                     ,intent(out) :: Status
66   integer,dimension(4)                     :: Ones
67   integer                                  :: stat
68   integer                                  :: i,k,n,nloc,nglb
69   real,dimension(Dimens(1)*Dimens(2)*Dimens(3)) :: buffer
71  !write(unit=0, fmt='(//3a,i6)') 'file: ',__FILE__,', line', __LINE__
72  !write(unit=0, fmt='(4x,a,i3,3a,i2)') 'DH%VarNames(', DH%CurrentVariable, ') = ', &
73  !      trim(DH%VarNames(DH%CurrentVariable)), ', NDim = ', NDim
74  !write(unit=0, fmt='(4x,5(a,i6))') &
75  !     'Dimens(1)=', Dimens(1), &
76  !     ', MemoryStart(1)=', MemoryStart(1), &
77  !     ', MemoryEnd(1)=', MemoryEnd(1), &
78  !     ', PatchStart(1)=', PatchStart(1), &
79  !     ', PatchEnd(1)=', PatchEnd(1)
81  !call pio_setdebuglevel(10)
83   Ones = 1
85   if(3 == NDim) then
86      stat = pio_get_var(DH%file_handle,DH%descVar(DH%CurrentVariable),Ones,Dimens(1:4),buffer)
87   else if(2 == NDim) then
88      stat = pio_get_var(DH%file_handle,DH%descVar(DH%CurrentVariable),Ones(1:3),Dimens(1:3),buffer)
89   else
90      stat = -1
91   end if
92   call netcdf_err(stat,Status)
93   if(Status /= WRF_NO_ERR) then
94     write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
95     call wrf_debug ( WARN , msg)
96   endif
98   if(3 == Ndim) then
99     call transRg2l(1,Dimens(1),1,Dimens(2),1,Dimens(3), &
100                    MemoryStart(1),MemoryEnd(1),MemoryStart(2),MemoryEnd(2),MemoryStart(3),MemoryEnd(3), &
101                     PatchStart(1), PatchEnd(1), PatchStart(2), PatchEnd(2), PatchStart(3), PatchEnd(3), &
102                    Data, buffer)
103   else if(2 == Ndim) then
104     call transRg2l(1,Dimens(1),1,Dimens(2),1,Dimens(3), &
105                    MemoryStart(1),MemoryEnd(1),MemoryStart(2),MemoryEnd(2),1,1, &
106                     PatchStart(1), PatchEnd(1), PatchStart(2), PatchEnd(2),1,1, &
107                    Data, buffer)
108   else 
109     write(unit=0, fmt='(/3a,i6)') 'file: ',__FILE__,', line', __LINE__
110     write(unit=0, fmt='(a,i6)') 'Do not know how handle NDim = ', NDim
111     write(unit=0, fmt='(4x,a,i4,a,i3,4x,a,i3)') &
112          'DH%vartype(', DH%CurrentVariable, ') =', DH%vartype(DH%CurrentVariable), &
113          'BDY_VAR =', BDY_VAR
114     write(unit=0, fmt='(4x,a,i3,3a,i2)') 'DH%VarNames(', DH%CurrentVariable, ') = ', &
115           trim(DH%VarNames(DH%CurrentVariable)), ', NDim = ', NDim
116   end if
117 end subroutine read_bdy_RealFieldIO
119 subroutine read_bdy_DoubleFieldIO(DH,NDim,Dimens,MemoryStart,MemoryEnd, &
120                                   PatchStart,PatchEnd,Data,Status)
121   use pio
122   use pio_kinds
123   use wrf_data_pio
124   use pio_routines
125   implicit none
126   include 'wrf_status_codes.h'
127   type(wrf_data_handle)                    :: DH
128   integer                     ,intent(in)  :: NDim
129   integer, dimension(*)       ,intent(in)  :: Dimens
130   integer, dimension(*)       ,intent(in)  :: MemoryStart, MemoryEnd
131   integer, dimension(*)       ,intent(in)  :: PatchStart,  PatchEnd
132   real*8,  dimension(*)       ,intent(out) :: Data
133   integer                     ,intent(out) :: Status
134   integer,dimension(4)                     :: Ones
135   integer                                  :: stat
136   integer                                  :: i,k,n,nloc,nglb
137   real*8,dimension(Dimens(1)*Dimens(2)*Dimens(3)) :: buffer
139  !write(unit=0, fmt='(//3a,i6)') 'file: ',__FILE__,', line', __LINE__
140  !write(unit=0, fmt='(4x,a,i3,3a,i2)') 'DH%VarNames(', DH%CurrentVariable, ') = ', &
141  !      trim(DH%VarNames(DH%CurrentVariable)), ', NDim = ', NDim
142  !write(unit=0, fmt='(4x,5(a,i6))') &
143  !     'Dimens(1)=', Dimens(1), &
144  !     ', MemoryStart(1)=', MemoryStart(1), &
145  !     ', MemoryEnd(1)=', MemoryEnd(1), &
146  !     ', PatchStart(1)=', PatchStart(1), &
147  !     ', PatchEnd(1)=', PatchEnd(1)
149  !call pio_setdebuglevel(10)
151   Ones = 1
153   if(3 == NDim) then
154      stat = pio_get_var(DH%file_handle,DH%descVar(DH%CurrentVariable),Ones,Dimens(1:4),buffer)
155   else if(2 == NDim) then
156      stat = pio_get_var(DH%file_handle,DH%descVar(DH%CurrentVariable),Ones(1:3),Dimens(1:3),buffer)
157   else
158      stat = -1
159   end if
160   call netcdf_err(stat,Status)
161   if(Status /= WRF_NO_ERR) then
162     write(msg,*) 'NetCDF error in ',__FILE__,', line', __LINE__
163     call wrf_debug ( WARN , msg)
164   endif
166   if(3 == Ndim) then
167     call transDg2l(1,Dimens(1),1,Dimens(2),1,Dimens(3), &
168                    MemoryStart(1),MemoryEnd(1),MemoryStart(2),MemoryEnd(2),MemoryStart(3),MemoryEnd(3), &
169                     PatchStart(1), PatchEnd(1), PatchStart(2), PatchEnd(2), PatchStart(3), PatchEnd(3), &
170                    Data, buffer)
171   else if(2 == Ndim) then
172     call transDg2l(1,Dimens(1),1,Dimens(2),1,Dimens(3), &
173                    MemoryStart(1),MemoryEnd(1),MemoryStart(2),MemoryEnd(2),1,1, &
174                     PatchStart(1), PatchEnd(1), PatchStart(2), PatchEnd(2),1,1, &
175                    Data, buffer)
176   else 
177     write(unit=0, fmt='(/3a,i6)') 'file: ',__FILE__,', line', __LINE__
178     write(unit=0, fmt='(a,i6)') 'Do not know how handle NDim = ', NDim
179     write(unit=0, fmt='(4x,a,i4,a,i3,4x,a,i3)') &
180          'DH%vartype(', DH%CurrentVariable, ') =', DH%vartype(DH%CurrentVariable), &
181          'BDY_VAR =', BDY_VAR
182     write(unit=0, fmt='(4x,a,i3,3a,i2)') 'DH%VarNames(', DH%CurrentVariable, ') = ', &
183           trim(DH%VarNames(DH%CurrentVariable)), ', NDim = ', NDim
184   end if
185 end subroutine read_bdy_DoubleFieldIO