updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / frame / module_intermediate_nmm.F
blob93394deb59f0ffb6c68516d7e451b7381ab24e9c
1 ! MODULE module_intermediate_nmm
3 ! This module contains routines that feed parent grid variables to the
4 ! intermediate grid when doing up-interpolation.  This is needed by
5 ! the new NMM interpolation routines, which require certain variables
6 ! on the target domain in order to do log(P)-space vertical
7 ! interpolation.
9 ! This module is also used during forcing (parent->nest boundary) to
10 ! copy variables to the intermediate domain that may not otherwise be
11 ! copied by the forcing routines.
13 ! Author: Samuel Trahan
15 ! History:
16 !   Aug 2012 - written by Sam Trahan for up-interpolation
17 !   Sep 2012 - updated to also work with forcing (parent->nest bdy)
19 module module_intermediate_nmm
20 #if (NMM_CORE == 1 && NMM_NEST==1)
21 contains
22   SUBROUTINE parent_to_inter_part1 ( grid, intermediate_grid, ngrid, config_flags )
23     USE module_state_description
24     USE module_domain, ONLY : domain, get_ijk_from_grid
25     USE module_configure, ONLY : grid_config_rec_type
26     USE module_dm, ONLY : ntasks, ntasks_x, ntasks_y, itrace, local_communicator, mytask, &
27          ipe_save, jpe_save, ips_save, jps_save, get_dm_max_halo_width,                   &
28          nest_pes_x, nest_pes_y,                                                          &
29          intercomm_active, nest_task_offsets,                                             &
30          mpi_comm_to_mom, mpi_comm_to_kid, which_kid!,                                     &
31          !push_communicators_for_domain,pop_communicators_for_domain
34     USE module_timing
35     IMPLICIT NONE
37     TYPE(domain), POINTER :: grid          
38     TYPE(domain), POINTER :: intermediate_grid
39     TYPE(domain), POINTER :: ngrid
40     INTEGER nlev, msize
41     INTEGER i,j,pig,pjg,cm,cn,nig,njg,retval,k,ioffset,ierr
42     INTEGER iparstrt,jparstrt,sw
43     TYPE (grid_config_rec_type)            :: config_flags
44     REAL xv(500)
45     INTEGER       ::          cids, cide, cjds, cjde, ckds, ckde,    &
46          cims, cime, cjms, cjme, ckms, ckme,    &
47          cips, cipe, cjps, cjpe, ckps, ckpe
48     INTEGER       ::          iids, iide, ijds, ijde, ikds, ikde,    &
49          iims, iime, ijms, ijme, ikms, ikme,    &
50          iips, iipe, ijps, ijpe, ikps, ikpe
51     INTEGER       ::          nids, nide, njds, njde, nkds, nkde,    &
52          nims, nime, njms, njme, nkms, nkme,    &
53          nips, nipe, njps, njpe, nkps, nkpe
55     INTEGER idim1,idim2,idim3,idim4,idim5,idim6,idim7
57     INTEGER icoord, jcoord, idim_cd, jdim_cd, pgr
58     INTEGER local_comm, myproc, nproc
59     INTEGER thisdomain_max_halo_width
61 !    CALL wrf_get_dm_communicator ( local_comm )
62     CALL wrf_get_myproc( myproc )
63     CALL wrf_get_nproc( nproc )
65     CALL get_ijk_from_grid (  grid ,                   &
66          cids, cide, cjds, cjde, ckds, ckde,    &
67          cims, cime, cjms, cjme, ckms, ckme,    &
68          cips, cipe, cjps, cjpe, ckps, ckpe    )
69     CALL get_ijk_from_grid (  intermediate_grid ,              &
70          iids, iide, ijds, ijde, ikds, ikde,    &
71          iims, iime, ijms, ijme, ikms, ikme,    &
72          iips, iipe, ijps, ijpe, ikps, ikpe    )
73     CALL get_ijk_from_grid (  ngrid ,              &
74          nids, nide, njds, njde, nkds, nkde,    &
75          nims, nime, njms, njme, nkms, nkme,    &
76          nips, nipe, njps, njpe, nkps, nkpe    )
78     CALL nl_get_parent_grid_ratio ( ngrid%id, pgr )
79     CALL nl_get_i_parent_start ( intermediate_grid%id, iparstrt )
80     CALL nl_get_j_parent_start ( intermediate_grid%id, jparstrt )
81     CALL nl_get_shw            ( intermediate_grid%id, sw )
82     icoord =    iparstrt - sw
83     jcoord =    jparstrt - sw
84     idim_cd = iide - iids + 1
85     jdim_cd = ijde - ijds + 1
87     nlev  = ckde - ckds + 1
89     CALL get_dm_max_halo_width ( ngrid%id , thisdomain_max_halo_width )
90     CALL wrf_dm_nestexchange_init
92     IF ( grid%active_this_task ) THEN
93       msize = 5
94       CALL rsl_lite_to_child_info( msize*4                               &
95          ,cips,cipe,cjps,cjpe                               &
96          ,iids,iide,ijds,ijde                               &
97          ,nids,nide,njds,njde                               &
98          ,pgr , sw                                          &
99          ,nest_task_offsets(ngrid%id)                       &
100          ,nest_pes_x(grid%id)                               &
101          ,nest_pes_y(grid%id)                               &
102          ,nest_pes_x(intermediate_grid%id)                  &
103          ,nest_pes_y(intermediate_grid%id)                  &
104          ,thisdomain_max_halo_width                         &
105          ,icoord,jcoord                                     &
106          ,idim_cd,jdim_cd                                   &
107          ,pig,pjg,retval )
108       DO while ( retval .eq. 1 )
109        IF ( SIZE(grid%hres_fis) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
110           xv(1)=grid%hres_fis(pig,pjg)
111           CALL rsl_lite_to_child_msg(4,xv)
112        ENDIF
113        IF ( SIZE(grid%sm) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
114           xv(1)=grid%sm(pig,pjg)
115           CALL rsl_lite_to_child_msg(4,xv)
116        ENDIF
117        IF ( SIZE(grid%pd) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
118           xv(1)=grid%pd(pig,pjg)
119           CALL rsl_lite_to_child_msg(4,xv)
120        ENDIF
121        IF ( SIZE(grid%fis) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
122           xv(1)=grid%fis(pig,pjg)
123           CALL rsl_lite_to_child_msg(4,xv)
124        ENDIF
125        CALL rsl_lite_to_child_info( msize*4                               &
126             ,cips,cipe,cjps,cjpe                               &
127             ,iids,iide,ijds,ijde                               &
128             ,nids,nide,njds,njde                               &
129             ,pgr , sw                                          &
130             ,nest_task_offsets(ngrid%id)                       &
131             ,nest_pes_x(grid%id)                               &
132             ,nest_pes_y(grid%id)                               &
133             ,nest_pes_x(intermediate_grid%id)                  &
134             ,nest_pes_y(intermediate_grid%id)                  &
135             ,thisdomain_max_halo_width                         &
136             ,icoord,jcoord                                     &
137             ,idim_cd,jdim_cd                                   &
138             ,pig,pjg,retval )
139       ENDDO
140     ENDIF ! grid%active_this_task
142     IF ( intercomm_active( grid%id ) ) THEN        ! I am parent
143       local_comm = mpi_comm_to_kid( which_kid(ngrid%id), grid%id )
144       ioffset = nest_task_offsets(ngrid%id)
145     ELSE IF ( intercomm_active( ngrid%id ) ) THEN  ! I am nest
146       local_comm = mpi_comm_to_mom( ngrid%id )
147       ioffset = nest_task_offsets(ngrid%id)
148     ENDIF
150     IF ( grid%active_this_task .OR. ngrid%active_this_task ) THEN
151 #if defined(DM_PARALLEL) && !defined(STUBMPI)
152       CALL mpi_comm_rank(local_comm,myproc,ierr)
153       CALL mpi_comm_size(local_comm,nproc,ierr)
154 #endif
155       CALL rsl_lite_bcast_msgs( myproc, nest_pes_x(grid%id)*nest_pes_y(grid%id),         &
156                                           nest_pes_x(ngrid%id)*nest_pes_y(ngrid%id),       &
157                                           ioffset, local_comm )
158     ENDIF
160     RETURN
161   END SUBROUTINE parent_to_inter_part1
163   SUBROUTINE parent_to_inter_part2 ( grid, ngrid, config_flags )
164     USE module_state_description
165     USE module_domain, ONLY : domain, get_ijk_from_grid
166     USE module_configure, ONLY : grid_config_rec_type
167     USE module_dm, ONLY : ntasks, ntasks_x, ntasks_y, itrace, local_communicator, mytask, &
168          ipe_save, jpe_save, ips_save, jps_save, get_dm_max_halo_width!,                   &
169          !push_communicators_for_domain,pop_communicators_for_domain
171     USE module_comm_dm, ONLY : HALO_NMM_INT_UP_sub
172     IMPLICIT NONE
174     TYPE(domain), POINTER :: grid          
175     TYPE(domain), POINTER :: cgrid
176     TYPE(domain), POINTER :: ngrid
178     INTEGER nlev, msize
179     INTEGER i,j,pig,pjg,cm,cn,nig,njg,retval,k
180     TYPE (grid_config_rec_type)            :: config_flags
181     REAL xv(500)
182     INTEGER       ::          cids, cide, cjds, cjde, ckds, ckde,    &
183          cims, cime, cjms, cjme, ckms, ckme,    &
184          cips, cipe, cjps, cjpe, ckps, ckpe
185     INTEGER       ::          nids, nide, njds, njde, nkds, nkde,    &
186          nims, nime, njms, njme, nkms, nkme,    &
187          nips, nipe, njps, njpe, nkps, nkpe
188     INTEGER       ::          ids, ide, jds, jde, kds, kde,    &
189          ims, ime, jms, jme, kms, kme,    &
190          ips, ipe, jps, jpe, kps, kpe
192     INTEGER idim1,idim2,idim3,idim4,idim5,idim6,idim7
193     REAL  dummy_xs, dummy_xe, dummy_ys, dummy_ye
195     integer myproc
197     CALL get_ijk_from_grid (  grid ,                   &
198          cids, cide, cjds, cjde, ckds, ckde,    &
199          cims, cime, cjms, cjme, ckms, ckme,    &
200          cips, cipe, cjps, cjpe, ckps, ckpe    )
202 IF ( ngrid%active_this_task ) THEN
203     nlev  = ckde - ckds + 1 
204     !write(0,*) 'IN parent_to_inter_part2'
205     CALL rsl_lite_from_parent_info(pig,pjg,retval)
206     DO while ( retval .eq. 1 )
207     !write(0,*) 'top of loop'
208        IF ( SIZE(grid%hres_fis) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
209           CALL rsl_lite_from_parent_msg(4,xv)
210           grid%hres_fis(pig,pjg) = xv(1)
211        ENDIF
212        !write(0,*)'do sm'
213        IF ( SIZE(grid%sm) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
214           CALL rsl_lite_from_parent_msg(4,xv)
215           grid%sm(pig,pjg) = xv(1)
216        ENDIF
217        !write(0,*)'do pd'
218        IF ( SIZE(grid%pd) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
219           CALL rsl_lite_from_parent_msg(4,xv)
220           grid%pd(pig,pjg) = xv(1)
221        ENDIF
222        !write(0,*)'do fis'
223        IF ( SIZE(grid%fis) .GT. 1 ) THEN ! okay for intermediate_grid too. see comment in gen_comms.c
224           CALL rsl_lite_from_parent_msg(4,xv)
225           grid%fis(pig,pjg) = xv(1)
226        ENDIF
227        !write(0,*) 'call rsl_lite_from_parent_info'
228        CALL rsl_lite_from_parent_info(pig,pjg,retval)
229        !write(0,*) 'back with retval=',retval
230     ENDDO
233     !write(0,*) 'out of loop'
235     CALL get_ijk_from_grid (  grid ,              &
236          ids, ide, jds, jde, kds, kde,    &
237          ims, ime, jms, jme, kms, kme,    &
238          ips, ipe, jps, jpe, kps, kpe    )
240     CALL push_communicators_for_domain( grid%id )
241 #include "HALO_NMM_INT_UP.inc"
242     CALL pop_communicators_for_domain
244 ENDIF
245     CALL wrf_dm_nestexchange_init
246     RETURN
247   END SUBROUTINE parent_to_inter_part2
248 #endif
249 end module module_intermediate_nmm