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