UPP  11.0.0
 All Data Structures Files Functions Variables Pages
MPI_FIRST.f
Go to the documentation of this file.
1 
2 !
39 !@PROCESS NOEXTCHK
40  SUBROUTINE mpi_first()
41 
42 !
43  use vrbls4d, only: dust, salt, soot, waso, suso, no3, nh4, pp25, pp10
44  use vrbls3d, only: u, v, t, q, uh, vh, wh, pmid, pmidv, pint, alpint, zmid, &
45  zint, q2, omga, t_adj, ttnd, rswtt, rlwtt, exch_h, train, tcucn, &
46  el_pbl, cwm, f_ice, f_rain, f_rimef, qqw, qqi, qqr, qqs,qqg, qqni, qqnr, &
47  extcof55, cfr, dbz, dbzr, dbzi, dbzc, mcvg, nlice, nrain, o3, vdifftt, &
48  tcucns, vdiffmois, dconvmois, sconvmois, nradtt, o3vdiff, o3prod, &
49  o3tndy, mwpv, unknown, vdiffzacce, zgdrag, cnvctummixing, vdiffmacce, &
50  mgdrag, cnvctvmmixing, ncnvctcfrac, cnvctumflx, cnvctdmflx, cnvctdetmflx,&
51  cnvctzgdrag, cnvctmgdrag, icing_gfip, asy, ssa, duem, dusd, dudp, &
52  duwt, suem, susd, sudp, suwt, ocem, ocsd, ocdp, ocwt, bcem, bcsd, &
53  bcdp, bcwt, ssem, sssd, ssdp, sswt, ext, dpres, rhomid, effri, effrl, &
54  effrs
55  use vrbls2d, only: wspd10max, w_up_max, w_dn_max, w_mean, refd_max, up_heli_max, &
56  prate_max, fprate_max, swupt, &
57  up_heli_max16, grpl_max, up_heli, up_heli16, ltg1_max, ltg2_max, &
58  up_heli_min, up_heli_min16, up_heli_max02, up_heli_min02, up_heli_max03, &
59  up_heli_min03, rel_vort_max, rel_vort_max01, wspd10umax, wspd10vmax, &
60  refdm10c_max, hail_max2d, hail_maxk1, ltg3_max,rel_vort_maxhy1, &
61  nci_ltg, nca_ltg, nci_wq, nca_wq, nci_refd, &
62  u10, v10, tshltr, qshltr, mrshltr, smstav, ssroff, bgroff, &
63  nca_refd, vegfrc, acsnow, acsnom, cmc, sst, qz0, thz0, uz0, vz0, qs, ths,&
64  sno, snonc, snoavg, psfcavg, t10m, t10avg, akmsavg, akhsavg, u10max, &
65  v10max, u10h, v10h, akms, akhs, cuprec, acprec, ancprc, cuppt, &
66  rainc_bucket, rainnc_bucket, pcp_bucket, snow_bucket, qrmax, tmax, &
67  snownc, graupelnc, tsnow, qvg, qv2m, rswin, rlwin, rlwtoa, tg, sfcshx, &
68  fis, t500, cfracl, cfracm, cfrach, acfrst, acfrcv, hbot, potevp, &
69  sfclhx, htop, aswin, alwin, aswout, alwout, aswtoa, alwtoa, czen, czmean,&
70  sigt4, rswout, radot, ncfrst, ncfrcv, smstot, pctsno, pshltr, th10, &
71  q10, sr, prec, subshx, snopcx, sfcuvx, sfcevp, z0, ustar, pblh, mixht, &
72  twbs, qwbs, sfcexc, grnflx, soiltb, z1000, slp, pslp, f, albedo, albase, &
73  cldfra, cprate, cnvcfr, ivgtyp, hbotd, htopd, hbots, isltyp, htops, &
74  cldefi, islope, si, lspa, rswinc, vis, pd, mxsnal, epsr, sfcux, &
75  sfcvx, sfcuxi, sfcvxi, avgalbedo, avgcprate, avgprec, ptop, pbot, avgcfrach, avgcfracm, &
76  avgcfracl, avgtcdc, auvbin, auvbinc, ptopl, pbotl, ttopl, ptopm, &
77  pbotm, ttopm, ptoph, pboth, ttoph, sfcugs, sfcvgs, pblcfr, cldwork, &
78  gtaux, gtauy, mdltaux, mdltauy, runoff, maxtshltr, mintshltr, &
79  maxrhshltr, minrhshltr, dzice, alwinc, alwoutc, alwtoac, aswinc, &
80  aswoutc,aswtoac, aswintoa, smcwlt, suntime, fieldcapa, avisbeamswin, &
81  avisdiffswin, airbeamswin, airdiffswin, snowfall, dusmass, ducmass, &
82  dusmass25, susmass, sucmass, susmass25, sucmass25, ocsmass, occmass, &
83  ocsmass25, occmass25, bcsmass, bccmass, bcsmass25, bccmass25, &
84  sssmass, sscmass, sssmass25, sscmass25, ducmass25, &
85  dustcb, sscb, bccb, occb, sulfcb, dustallcb, ssallcb,dustpm,sspm, pp25cb,&
86  no3cb, nh4cb, dustpm10, pp10cb, maod, ti
87  use soil, only: smc, stc, sh2o, sldpth, rtdpth, sllevel
88  use masks, only: htm, vtm, hbm2, sm, sice, lmh, gdlat, gdlon, dx, dy, lmv
89  use ctlblk_mod, only: me, num_procs, jm, jsta, jend, jsta_m, jsta_m2,ista,iend , &
90  jend_m, jend_m2, iup, idn, icnt, im, idsp, jsta_2l, jend_2u,idsp2,icnt2, &
91  jvend_2u, lm, lp1, jsta_2l, jend_2u, nsoil, nbin_du, nbin_ss, &
92  nbin_bc, nbin_oc, nbin_su, nbin_no3, nbin_nh4, &
93  ista_m,iend_m,ista_m2,iend_m2, ista_m,iend_m,ista_m2,iend_m2, &
94  ileft,iright,ileftb,irightb,ibsize,ibsum, isxa,iexa,jsxa,jexa, &
95  icoords,ibcoords,bufs,ibufs, rbufs, rcoords,rbcoords, &
96  ista_2l, iend_2u,ivend_2u,numx,modelname
97 
98 !
99 ! use params_mod
100 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101  implicit none
102 
103  include 'mpif.h'
104 !
105  integer ierr,i,jsx,jex,isx,iex,j
106  integer size,ubound,lbound
107  integer isumm,isum ,ii,jj,isumm2
108  integer , allocatable :: ibuff(:)
109  real , allocatable :: rbuff(:)
110  integer, allocatable :: ipole(:),ipoles(:,:)
111  real , allocatable :: rpole(:),rpoles(:,:)
112 
113  isumm=0
114  isumm2=0
115 
116  if ( me == 0 ) then
117  write(*,*) ' NUM_PROCS,NUMX,NUMY = ',num_procs,numx,num_procs/numx
118  end if
119 
120  if ( num_procs > 1024 ) then
121  print *, ' too many MPI tasks, max is 1024, stopping'
122  call mpi_abort(mpi_comm_world,1,ierr)
123  stop
124  end if
125 
126 ! error check
127 
128  if ( num_procs > jm/2 ) then
129  print *, ' too many MPI tasks, max is ',jm/2,' stopping'
130  call mpi_abort(mpi_comm_world,1,ierr)
131  stop
132  end if
133 
134 ! global loop ranges
135 !
136 ! para_range2 supports a 2D decomposition. The rest of the post
137 ! supports 1D still and the call here is the special case where each
138 ! processor gets all of the longitudes in the latitude 1D subdomain
139 ! jsta:jend. The X decomposition will be specified by the third
140 ! argument (currently 1) and the Y decoposition will be specified by
141 ! the fourth argument (currently all of the ranks) When X is
142 ! subdivided the third and fourth arguments will have to be integral
143 ! factors of num_procs
144 
145  call para_range2(im,jm,numx,num_procs/numx,me,ista,iend,jsta,jend)
146 
147  jsta_m = jsta
148  jsta_m2 = jsta
149  jend_m = jend
150  jend_m2 = jend
151  ista_m = ista
152  ista_m2 = ista
153  iend_m = iend
154  iend_m2 = iend
155 
156  if (me<numx)then
157  jsta_m=2
158  jsta_m2=3
159  end if
160 
161  if(mod(me,numx)==0)then
162  ista_m=2
163  ista_m2=3
164  end if
165 
166  if (me>=(num_procs-numx))then
167  jend_m=jm-1
168  jend_m2=jm-2
169  end if
170 
171  if(mod(me+1,numx)==0)then
172  iend_m=im-1
173  iend_m2=im-2
174  end if
175 
176  102 format(6i10,a20)
177 
178 !
179  if ( me == 0 ) then
180  idn = mpi_proc_null
181  end if
182  if ( me == num_procs - 1 ) then
183  iup = mpi_proc_null
184  end if
185 !
186 ! GWV. Array of i/j coordinates for bookkeeping tests. Not used in
187 ! calculations but to check if scatter,gather, and exchanges are doing as
188 ! expected. Both real and integer arrays are sent. Integer will be needed
189 ! for very large domains because real mantissas overflow and both coordinates'
190 ! information can't be packed into a real mantisa. Real is easier to use
191 ! because the datatype is the same as for actual data
192 
193  allocate(icoords(im,jm))
194  allocate(rcoords(im,jm))
195  allocate(ibuff(im*jm))
196  allocate(rbuff(im*jm))
197  do j=1,jm
198  do i=1,im
199  icoords(i,j)=10000*i+j ! both I and J information is in each element
200  rcoords(i,j)=4000*i+j ! both I and J information is in each element but it overflows for large I I to 3600 is safe
201  end do
202  end do
203 
204 ! end COORDS test
205 
206 ! counts, disps for gatherv and scatterv
207 
208  isum=1
209  allocate(isxa(0:num_procs-1) )
210  allocate(jsxa(0:num_procs-1) )
211  allocate(iexa(0:num_procs-1) )
212  allocate(jexa(0:num_procs-1) )
213  do i = 0, num_procs - 1
214  call para_range2(im,jm,numx,num_procs/numx,i,isx,iex,jsx,jex)
215  icnt(i) = ((jex-jsx)+1)*((iex-isx)+1)
216  isxa(i)=isx
217  iexa(i)=iex
218  jsxa(i)=jsx
219  jexa(i)=jex
220 
221  idsp(i)=isumm
222  isumm=isumm+icnt(i)
223  if(jsx .eq. 1 .or. jex .eq. jm) then
224  icnt2(i) = (iex-isx+1)
225  else
226  icnt2(i)=0
227  endif
228  idsp2(i)=isumm2
229  if(jsx .eq. 1 .or. jex .eq. jm) isumm2=isumm2+(iex-isx+1)
230 
231 ! GWV Create send buffer for scatter. This is now needed because we are no
232 ! longer sending contiguous slices of the im,jm full state arrays to the
233 ! processors with scatter. Instead we are sending a slice of I and a slice of J
234 ! and so have to reshape the send buffer below to make it contiguous groups of
235 ! isx:iex,jsx:jex arrays
236 
237  do jj=jsx,jex
238  do ii=isx,iex
239  ibuff(isum)=icoords(ii,jj)
240  rbuff(isum)=rcoords(ii,jj)
241  isum=isum+1
242  end do
243  end do
244 
245  end do ! enddo of num_procs
246 !
247 ! extraction limits -- set to two rows
248 !
249  jsta_2l = max(jsta - 2, 1 )
250  jend_2u = min(jend + 2, jm )
251  if(modelname=='GFS') then
252  ista_2l=max(ista-2,0)
253  iend_2u=min(iend+2,im+1)
254  else
255  ista_2l=max(ista-2,1)
256  iend_2u=min(iend+2,im)
257  endif
258 
259 ! special for c-grid v
260  jvend_2u = min(jend + 2, jm+1 )
261 !
262 ! NEW neighbors
263 
264  ileft = me - 1
265  iright = me + 1
266  iup=mpi_proc_null
267  idn=mpi_proc_null
268 
269  if(mod(me,numx) .eq. 0) print *,' LEFT POINT',me
270  if(mod(me+1,numx) .eq. 0) print *,' RIGHT POINT',me
271  if(mod(me,numx) .eq. 0) ileft=mpi_proc_null
272  if(mod(me,numx) .eq. 0) ileftb=me+numx-1
273  if(mod(me+1,numx) .eq. 0 .or. me .eq. num_procs-1) iright=mpi_proc_null
274  if(mod(me+1,numx) .eq. 0 .or. me .eq. num_procs-1) irightb=me-numx+1
275  if(me .ge. numx) idn=me-numx
276  if(me+1 .le. num_procs-numx) iup=me+numx
277 
278  print 102,me,ileft,iright,iup,idn,num_procs,'GWVX BOUNDS'
279 
280 ! allocate arrays
281 
282  ibsize = ( (iend-ista) +1) * ( (jend-jsta)+1)
283  allocate(ibcoords(ista_2l:iend_2u,jsta_2l:jend_2u))
284  allocate(rbcoords(ista_2l:iend_2u,jsta_2l:jend_2u))
285  allocate(ibufs(ibsize))
286  allocate(rbufs(ibsize))
287  call mpi_scatterv(ibuff,icnt,idsp,mpi_integer &
288  ,ibufs,icnt(me),mpi_integer ,0,mpi_comm_world,j)
289  call mpi_scatterv(rbuff,icnt,idsp,mpi_real &
290  ,rbufs,icnt(me),mpi_real ,0,mpi_comm_world,j)
291 
292 !
293 !GWV reshape the receive subdomain
294 
295  isum=1
296  do j=jsta,jend
297  do i=ista,iend
298  ibcoords(i,j)=ibufs(isum)
299  rbcoords(i,j)=rbufs(isum)
300  isum=isum+1
301  end do
302  end do
303 
304 !GWV end reshape
305  do j=jsta,jend
306  do i=ista,iend
307  ii=ibcoords(i,j)/10000
308  jj=ibcoords( i,j)-(ii*10000)
309  if(ii .ne. i .or. jj .ne. j) then
310  print *,i,j,ii,jj,ibcoords(i,j),' GWVX FAIL '
311  else
312  continue
313  endif
314  end do
315  end do
316 
317  allocate(ipoles(im,2),ipole(ista:iend))
318  allocate(rpoles(im,2),rpole(ista:iend))
319  ipole=9900000
320  ipoles=-999999999
321 
322  do i=ista,iend
323  if(me .lt. num_procs/2. .and. jsta_2l .le. 1 .and. icnt2(me) .gt. 0) ipole(i)=ibcoords(i,1)
324  if(me .lt. num_procs/2. .and. jsta_2l .le. 1 .and. icnt2(me) .gt. 0) rpole(i)=rbcoords(i,1)
325  if(me .gt. num_procs/2. .and. jend_2u .ge. jm .and. icnt2(me) .gt. 0) ipole(i)=ibcoords(i,jm)
326  if(me .gt. num_procs/2. .and. jend_2u .ge. jm .and. icnt2(me) .gt. 0) rpole(i)=rbcoords(i,jm)
327 
328 ! check code to be removed upon debugging
329  if(me .lt. num_procs/2. .and. jsx .eq. 1) then
330  continue
331  endif
332  if(me .gt. num_procs/2. .and. jend_2u .ge. jm) then
333  continue
334  endif
335  end do ! end check code
336 
337 ! test pole gather
338  print 105,' GWVX GATHER DISP ',icnt2(me),idsp2(me),me
339  105 format(a30,3i12)
340 
341  call mpi_gatherv(ipole(ista),icnt2(me),mpi_integer, ipoles,icnt2,idsp2,mpi_integer,0,mpi_comm_world, ierr )
342  call mpi_gatherv(rpole(ista),icnt2(me),mpi_real , rpoles,icnt2,idsp2,mpi_real ,0,mpi_comm_world, ierr )
343 
344  if(me .eq. 0) then
345  do j=1,2
346  do i=1,im
347  ii=rpoles(i,j)/4000
348  jj=rpoles(i,j) -ii*4000
349  if(ii .ne. i .or. jj .ne. 1 .and. jj .ne. jm ) then
350  write(*,169)' IPOLES BAD POINT',rpoles(i,j),ii,i,jj,' jm= ',jm
351  else
352  continue
353 ! write(*,169)' IPOLES GOOD POINT',rpoles(i,j),ii,i,jj,' jm= ',jm
354  endif
355  end do
356  end do
357  endif
358 
359  107 format(a20,10i10)
360  169 format(a25,f20.1,3i10,a10,4i10)
361 !
362  print *, ' me, jsta_2l, jend_2u = ',me,jsta_2l, jend_2u, &
363  'jvend_2u=',jvend_2u,'im=',im,'jm=',jm,'lm=',lm, &
364  'lp1=',lp1
365  write(*,'(A,5I10)') 'MPI_FIRST me,jsta,jend,ista,iend,=',me,jsta,jend,ista,iend
366 
367  end
368 
369 ! subroutine sub(a)
370 ! return
371 ! end
372 
373 
374 
375  subroutine fullpole(a,rpoles)
376 
377  use ctlblk_mod, only: num_procs, jend, iup, jsta, idn, mpi_comm_comp, im,modelname,numx,&
378  icoords,ibcoords,rbcoords,bufs,ibufs,me, &
379  jsta_2l,jend_2u,ileft,iright,ista_2l,iend_2u,ista,iend,jm,icnt2,idsp2
380 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
381  implicit none
382 !
383  include 'mpif.h'
384 !
385  real,intent(inout) :: a ( ista_2l:iend_2u,jsta_2l:jend_2u ),rpoles(im,2)
386  real, allocatable :: rpole(:)
387 
388  integer status(mpi_status_size)
389  integer ierr
390  integer size,ubound,lbound
391  integer i,ii,jj, ibl,ibu,jbl,jbu,icc,jcc
392  integer ifirst
393  data ifirst/0/
394  integer iwest,ieast
395  data iwest,ieast/0,0/
396  allocate(rpole(ista:iend))
397 
398  do i=ista,iend
399  if(me .lt. num_procs/2. .and. jsta_2l .le. 1 .and. icnt2(me) .gt. 0) rpole(i)=a(i,1)
400  if(me .ge. num_procs/2. .and. jend_2u .ge. jm .and. icnt2(me) .gt. 0) rpole(i)=a(i,jm)
401  end do
402 
403  call mpi_allgatherv(rpole(ista),icnt2(me),mpi_real,rpoles,icnt2,idsp2,mpi_real, mpi_comm_comp,ierr)
404 
405  call mpi_barrier(mpi_comm_comp,ierr)
406  ifirst=1
407 
408  end
409 
Definition: MASKS_mod.f:1
Definition: SOIL_mod.f:1
subroutine mpi_first()
SUBPROGRAM: MPI_FIRST SET UP MESSGAE PASSING INFO PRGRMMR: TUCCILLO ORG: IBM.
Definition: MPI_FIRST.f:40