module_kma2netcdf_interface.f90

References to this file elsewhere.
1 MODULE module_kma2netcdf_interface
2 
3    use module_kma_wave2grid
4 
5    USE module_domain
6    USE module_timing
7    USE module_driver_constants
8    USE module_configure
9    use module_kma_wave2grid
10    USE module_tiles
11 !  implicit none    !shc-wei
12 
13 CONTAINS
14 
15 SUBROUTINE kma2netcdf_interface ( grid, config_flags)
16 
17 !  IMPLICIT NONE    !shc-wei
18 
19 !--Input data.
20    TYPE (grid_config_rec_type)   :: config_flags
21    TYPE(domain), TARGET          :: grid
22 #ifdef DEREF_KLUDGE
23    INTEGER     :: sm31 , em31 , sm32 , em32 , sm33 , em33
24 
25    sm31             = grid%sm31
26    em31             = grid%em31
27    sm32             = grid%sm32
28    em32             = grid%em32
29    sm33             = grid%sm33
30    em33             = grid%em33
31 #endif
32 
33    call kma2netcdf_solver( grid, config_flags  &
34 
35 #include "em_actual_args.inc"
36 
37                )
38 
39 end SUBROUTINE kma2netcdf_interface
40 
41 
42 SUBROUTINE kma2netcdf_solver( grid, config_flags &
43 #include "em_dummy_args.inc"
44                  )
45 
46 !  IMPLICIT NONE   !shc-wei
47 
48 !--Input data.
49    TYPE (grid_config_rec_type)   :: config_flags
50    TYPE(domain), TARGET               :: grid
51 !  Definitions of dummy arguments to solve
52 
53 #include "em_dummy_decl.inc"
54 !    INCLUDE 'mpif.h'
55     real, allocatable      :: q(:,:,:) 
56     Integer                :: my_proc_id, ierr 
57     Integer                :: ii, jj,landmask_T213(428,215)  
58     real                   :: sfc_T213(428,215)
59 !---------------------------------------------------------------------------
60    INTEGER  :: ids,ide, jds,jde, kds,kde   ! domain dims.
61    INTEGER  :: ims,ime, jms,jme, kms,kme   ! memory dims.
62    INTEGER  :: ips,ipe, jps,jpe, kps,kpe   ! patch  dims.
63    INTEGER  :: its,ite, jts,jte, kts,kte   ! Tile   dims.
64    INTEGER  :: i
65    INTEGER  :: IGRD, JGRD, KGRD, JCAP, KMAX, INTVL                           
66 !rizvi add ------------------------------------------------------------------
67    NAMELIST /kma2netcdf_parm/ IGRD, JGRD, KGRD, JCAP, KMAX, INTVL                           
68 !
69       READ  (111, NML = kma2netcdf_parm, ERR = 8000)
70       close (111)
71       write(unit=*, fmt='(A,5(1x,/,5x,A,i6))')'kma2netcdf_parm namelist read: ',&
72       'IGRD= ',IGRD,'JGRD= ',JGRD,'KGRD= ',KGRD,'JCAP = ',JCAP,'KMAX= ',KMAX,'INTVL= ',INTVL
73       DPHI=180./(JGRD-1)
74       LMAX=KGRD-1
75       KLMAX=KMAX
76       MEND1 = JCAP + 1
77       NEND1 = JCAP + 1
78       JEND1 = JCAP + 1
79       IMAXG = IGRD 
80       JMAXG = JGRD - 1 
81 
82       IMAX  =360./DPHI+0.5
83       IOUT  =IMAX/INTVL  
84       JMAX  =180./DPHI+1.5
85       JOUT  =(JMAX-1)/INTVL+1
86       IMX   =IMAX+2
87       JOUTHF= (JOUT+1)/2
88       JMXGHF= (1+JMAXG)/2
89       KMX2  =KMAX*2  
90       LMX2 =LMAX*2
91       MNWAV =MEND1*(MEND1+1)/2
92 
93 
94 !rizvi add ------------------------------------------------------------------
95    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
96 
97   call copy_dims( grid, xp, &
98                    ids, ide, jds, jde, kds, kde, &
99                    ims, ime, jms, jme, kms, kme, &
100                    ips, ipe, jps, jpe, kps, kpe  )
101 !--Compute these starting and stopping locations for each tile and number of tiles.
102 
103    CALL set_tiles ( grid , ids , ide , jds , jde , ips , ipe , jps , jpe )
104 
105    call copy_tile_dims( grid, xp, its, ite, jts, jte, kts, kte )
106 
107    allocate (q(ims:ime,jms:jme,kms:kme))
108     
109 !   go to   100
110 #ifndef DEREF_KLUDGE
111      call W2GCONV(IGRD, JGRD, KGRD, JCAP, KMAX, INTVL , &
112              DPHI, LMAX, KLMAX, MEND1, NEND1, JEND1, IMAXG, JMAXG, &
113              IMAX, IOUT, JMAX, JOUT, IMX, JOUTHF, JMXGHF, KMX2, LMX2, MNWAV, &                            
114              ht, psfc, t_2, u_2, v_2, q, KMA_A, KMA_B, &
115              xp%ims, xp%jms, xp%kms, xp%ime, xp%jme, xp%kme,&
116              xp%ids, xp%jds, xp%kds, xp%ide, xp%jde, xp%kde,&
117              xp%its, xp%jts, xp%kts, xp%ite, xp%jte, xp%kte )
118 #else
119      call W2GCONV(IGRD, JGRD, KGRD, JCAP, KMAX, INTVL , &
120              DPHI, LMAX, KLMAX, MEND1, NEND1, JEND1, IMAXG, JMAXG, &
121              IMAX, IOUT, JMAX, JOUT, IMX, JOUTHF, JMXGHF, KMX2, LMX2, MNWAV, &                            
122      ht(ims,jms),psfc(ims,jms),t_2(ims,jms,kms),& 
123      u_2(ims,jms,kms), v_2(ims,jms,kms), q(ims,jms,kms),&
124      KMA_A(kms),KMA_B(kms), ims, jms, kms, ime, jme, kme, &
125                             ids, jds, kds, ide, jde, kde, &
126                             its, jts, kts, ite, jte, kte  )
127 #endif
128 
129 !  convert KMA pressure which is in hPa into Psacal in grid-array
130     psfc(its:ite,jts:jte) = 100. * psfc(its:ite,jts:jte)
131 100  continue
132  my_proc_id = 0
133 #ifdef DM_PARALLEL
134    call MPI_COMM_RANK( MPI_COMM_WORLD, my_proc_id, ierr )
135 #else
136    my_proc_id = 0
137 #endif
138 
139    moist(ims:ime,jms:jme,kms:kme,P_qv) = q(ims:ime,jms:jme,kms:kme)
140    deallocate (q)                 
141 ! Load landmask from KMA-original land mask for T213
142     if( JCAP == 213 ) then
143      open( UNIT = 151, file= 'KMA_landmask_428_215', status= 'old')
144        do jj=jds,jde
145         read(151,'(428i1)', err=9000) (landmask_T213(ii,jj),ii=1,428)
146        enddo
147        do jj=jts,jte
148         landmask(its:ite,jj) = landmask_T213(its:ite,jde-jj+1)
149        enddo
150        write(unit=*, fmt='(A)')'read successfully landmask'
151        close (151) 
152 ! Load U10 at T213 (428x215) resolution                          
153      open( UNIT = 151, file= 'nwpgr_UUUU.2007020100', status= 'old')
154         read(151,'(10e20.10)', err=9000) sfc_T213                     
155        do jj=jts,jte
156         u10(its:ite,jj) = sfc_T213(its:ite,jde-jj+1)
157        enddo
158        write(unit=*, fmt='(A)')'read successfully U10'
159        close (151) 
160 ! Load V10 at T213 (428x215) resolution                          
161      open( UNIT = 151, file= 'nwpgr_VVVV.2007020100', status= 'old')
162         read(151,'(10e20.10)', err=9000) sfc_T213                     
163        do jj=jts,jte
164         v10(its:ite,jj) = sfc_T213(its:ite,jde-jj+1)
165        enddo
166        write(unit=*, fmt='(A)')'read successfully V10'
167        close (151) 
168 ! Load T2  at T213 (428x215) resolution                          
169      open( UNIT = 151, file= 'nwpgr_TTTT.2007020100', status= 'old')
170         read(151,'(10e20.10)', err=9000) sfc_T213                     
171        do jj=jts,jte
172          t2(its:ite,jj) = sfc_T213(its:ite,jde-jj+1)
173        enddo
174        write(unit=*, fmt='(A)')'read successfully T2'
175        close (151) 
176 ! Load Q2  at T213 (428x215) resolution                          
177      open( UNIT = 151, file= 'nwpgr_QQQQ.2007020100', status= 'old')
178         read(151,'(10e20.10)', err=9000) sfc_T213                     
179        do jj=jts,jte
180          q2(its:ite,jj) = sfc_T213(its:ite,jde-jj+1)
181        enddo
182        write(unit=*, fmt='(A)')'read successfully Q2'
183        close (151) 
184 ! Load SST at T213 (428x215) resolution                          
185      open( UNIT = 151, file= 'nwpgr_SSTT.2007020100', status= 'old')
186         read(151,'(10e20.10)', err=9000) sfc_T213                     
187        do jj=jts,jte
188          sst(its:ite,jj) = sfc_T213(its:ite,jde-jj+1)
189        enddo
190        write(unit=*, fmt='(A)')'read successfully SST'
191        close (151) 
192 
193     else
194     write(unit=*, fmt='(A,i3)')'Surface data is not available for T',JCAP
195     endif
196 
197     write(unit=*, fmt='(A)')'Job done for kma2netcdf_solver' 
198     return
199 8000 write(unit=*, fmt='(A)')'read error on namelist unit 111'
200      stop
201 9000 write(unit=*, fmt='(A)')'read error on unit 151'
202      stop
203 END SUBROUTINE kma2netcdf_solver
204 
205 SUBROUTINE copy_dims(grid, xp, &
206                      ids, ide, jds, jde, kds, kde, &
207                      ims, ime, jms, jme, kms, kme, &
208                      ips, ipe, jps, jpe, kps, kpe )
209 !------------------------------------------------------------------------------
210 !  PURPOSE: Copy dimensioning information from grid structure.
211 !
212 !------------------------------------------------------------------------------
213 
214    USE module_domain
215    TYPE(domain), INTENT(IN)         :: grid
216    TYPE (xpose_type),INTENT(INOUT)  :: xp      ! Transpose variables.
217 
218    INTEGER,      INTENT(OUT)        :: ids,ide, jds,jde, kds,kde   ! domain dims.
219    INTEGER,      INTENT(OUT)        :: ims,ime, jms,jme, kms,kme   ! memory dims.
220    INTEGER,      INTENT(OUT)        :: ips,ipe, jps,jpe, kps,kpe   ! patch  dims.
221 
222 !--De-reference dimension information stored in the grid data structure.
223 
224    ids = grid%sd31 
225    ide = grid%ed31 - 1
226    jds = grid%sd32 
227    jde = grid%ed32 - 1
228    kds = grid%sd33 
229    kde = grid%ed33 - 1
230 
231    ims = grid%sm31 
232    ime = grid%em31 
233    jms = grid%sm32 
234    jme = grid%em32 
235    kms = grid%sm33 
236    kme = grid%em33 
237 
238    ips = grid%sp31 
239    ipe = grid%ep31 
240    jps = grid%sp32 
241    jpe = grid%ep32 
242    kps = grid%sp33 
243    kpe = grid%ep33 
244 
245 !--Indices for yz decomposition
246 
247    xp%idsx = grid%sd31
248    xp%idex = grid%ed31 - 1
249    xp%jdsx = grid%sd32
250    xp%jdex = grid%ed32 - 1
251    xp%kdsx = grid%sd33
252    xp%kdex = grid%ed33 - 1
253 
254    xp%imsx = grid%sm31x
255    xp%imex = grid%em31x
256    xp%jmsx = grid%sm32x
257    xp%jmex = grid%em32x
258    xp%kmsx = grid%sm33x
259    xp%kmex = grid%em33x
260 
261    xp%itsx = grid%sp31x
262    xp%itex = grid%ep31x
263    xp%jtsx = grid%sp32x
264    xp%jtex = grid%ep32x
265    xp%ktsx = grid%sp33x
266    xp%ktex = grid%ep33x
267 
268    xp%ipsx = grid%sp31x
269    xp%ipex = grid%ep31x
270    xp%jpsx = grid%sp32x
271    xp%jpex = grid%ep32x
272    xp%kpsx = grid%sp33x
273    xp%kpex = grid%ep33x
274 
275 !--Indices for xz decomposition
276 
277    xp%idsy = grid%sd31
278    xp%idey = grid%ed31 - 1
279    xp%jdsy = grid%sd32
280    xp%jdey = grid%ed32 - 1
281    xp%kdsy = grid%sd33
282    xp%kdey = grid%ed33 - 1
283 
284    xp%imsy = grid%sm31y
285    xp%imey = grid%em31y
286    xp%jmsy = grid%sm32y
287    xp%jmey = grid%em32y
288    xp%kmsy = grid%sm33y
289    xp%kmey = grid%em33y
290 
291    xp%itsy = grid%sp31y
292    xp%itey = grid%ep31y
293    xp%jtsy = grid%sp32y
294    xp%jtey = grid%ep32y
295    xp%ktsy = grid%sp33y
296    xp%ktey = grid%ep33y
297 
298    xp%ipsy = grid%sp31y
299    xp%ipey = grid%ep31y
300    xp%jpsy = grid%sp32y
301    xp%jpey = grid%ep32y
302    xp%kpsy = grid%sp33y
303    xp%kpey = grid%ep33y
304 
305    if(ipe > ide) ipe = ide
306    if(jpe > jde) jpe = jde
307    if(kpe > kde) kpe = kde
308 
309    ! Indices for yz decomposition
310 
311    if(xp%itex > ide) xp%itex = ide
312    if(xp%jtex > jde) xp%jtex = jde
313    if(xp%ktex > kde) xp%ktex = kde
314 
315    if(xp%ipex > ide) xp%ipex = ide
316    if(xp%jpex > jde) xp%jpex = jde
317    if(xp%kpex > kde) xp%kpex = kde
318 
319    ! Indices for xz decomposition
320 
321    if(xp%itey > ide) xp%itey = ide
322    if(xp%jtey > jde) xp%jtey = jde
323    if(xp%ktey > kde) xp%ktey = kde
324 
325    if(xp%ipey > ide) xp%ipey = ide
326    if(xp%jpey > jde) xp%jpey = jde
327    if(xp%kpey > kde) xp%kpey = kde
328 
329 !  Copy xpose dimensions from grid structure to xp structure.
330 
331 !--Indices for xy decomposition
332 
333    xp%ids = ids
334    xp%ide = ide
335    xp%jds = jds
336    xp%jde = jde
337    xp%kds = kds
338    xp%kde = kde
339 
340    xp%ims = ims
341    xp%ime = ime
342    xp%jms = jms
343    xp%jme = jme
344    xp%kms = kms
345    xp%kme = kme
346 
347    xp%ips = ips
348    xp%ipe = ipe
349    xp%jps = jps
350    xp%jpe = jpe
351    xp%kps = kps
352    xp%kpe = kpe
353 
354 END SUBROUTINE copy_dims
355 
356 SUBROUTINE copy_tile_dims( grid, xp, its, ite, jts, jte, kts, kte )
357 
358 !------------------------------------------------------------------------------
359 !  PURPOSE: Copy tile dimensions from grid structure.
360 !
361 !------------------------------------------------------------------------------
362 
363 !   USE module_domain
364    TYPE(domain), INTENT(IN)         :: grid
365    TYPE (xpose_type),INTENT(INOUT)  :: xp      ! Transpose variables.
366    INTEGER,      INTENT(OUT)        :: its,ite, jts,jte, kts,kte ! tile dims.
367 
368    INTEGER                  :: ij   ! Loop counter
369 
370 !  De-reference tile indices stored in the grid data structure.
371 
372    DO ij = 1 , grid%num_tiles
373      its = grid%i_start(ij)
374      ite = grid%i_end(ij)
375      jts = grid%j_start(ij)
376      jte = grid%j_end(ij)
377      kts = xp%kds
378      kte = xp%kde
379 
380      xp%its = its
381      xp%ite = ite
382      xp%jts = jts
383      xp%jte = jte
384      xp%kts = kts
385      xp%kte = kte
386 
387      if(xp%ite > xp%ide) xp%ite = xp%ide
388      if(xp%jte > xp%jde) xp%jte = xp%jde
389      if(xp%kte > xp%kde) xp%kte = xp%kde
390 
391      if(ite > xp%ide) ite = xp%ide
392      if(jte > xp%jde) jte = xp%jde
393      if(kte > xp%kde) kte = xp%kde
394 
395         write(unit=*, fmt='(/)')
396         write(unit=*, fmt='(A,i3,A,5x,3(i3,A,i3,5x))') 'Tile ',ij, &
397                 ' size:', its,':',ite, jts,':',jte, kts,':',kte
398    END DO
399 END SUBROUTINE copy_tile_dims
400 
401 END MODULE module_kma2netcdf_interface