UPP (upp-srw-2.2.0)
Loading...
Searching...
No Matches
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.
137! The X decomposition is specified by the third
138! argument and the Y decoposition is specified by
139! the fourth argument. The product of the third and fourth arguments
140! must be num_procs and the third and fourth arguments must be integral
141! factors of num_procs.
142!
143! for the special case of 1D decomposition, numx is set to 1 and the
144! fourth argument becomes the number of MPI ranks for the job. numx=1
145! makes the code fully compatible with the old 1D decomposition.
146
147
148 call para_range2(im,jm,numx,num_procs/numx,me,ista,iend,jsta,jend)
149
150 jsta_m = jsta
151 jsta_m2 = jsta
152 jend_m = jend
153 jend_m2 = jend
154 ista_m = ista
155 ista_m2 = ista
156 iend_m = iend
157 iend_m2 = iend
158
159 if (me<numx)then
160 jsta_m=2
161 jsta_m2=3
162 end if
163
164 if(mod(me,numx)==0)then
165 ista_m=2
166 ista_m2=3
167 end if
168
169 if (me>=(num_procs-numx))then
170 jend_m=jm-1
171 jend_m2=jm-2
172 end if
173
174 if(mod(me+1,numx)==0)then
175 iend_m=im-1
176 iend_m2=im-2
177 end if
178
179 102 format(6i10,a20)
180
181!
182 if ( me == 0 ) then
183 idn = mpi_proc_null
184 end if
185 if ( me == num_procs - 1 ) then
186 iup = mpi_proc_null
187 end if
188!
189! GWV. Array of i/j coordinates for bookkeeping tests. Not used in
190! calculations but to check if scatter,gather, and exchanges are doing as
191! expected. Both real and integer arrays are sent. Integer will be needed
192! for very large domains because real mantissas overflow and both coordinates'
193! information can't be packed into a real mantisa. Real is easier to use
194! because the datatype is the same as for actual data
195
196 allocate(icoords(im,jm))
197 allocate(rcoords(im,jm))
198 allocate(ibuff(im*jm))
199 allocate(rbuff(im*jm))
200 do j=1,jm
201 do i=1,im
202 icoords(i,j)=10000*i+j ! both I and J information is in each element
203 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
204 end do
205 end do
206
207! end COORDS test
208
209! counts, disps for gatherv and scatterv
210
211 isum=1
212 allocate(isxa(0:num_procs-1) )
213 allocate(jsxa(0:num_procs-1) )
214 allocate(iexa(0:num_procs-1) )
215 allocate(jexa(0:num_procs-1) )
216 do i = 0, num_procs - 1
217 call para_range2(im,jm,numx,num_procs/numx,i,isx,iex,jsx,jex)
218 icnt(i) = ((jex-jsx)+1)*((iex-isx)+1)
219 isxa(i)=isx
220 iexa(i)=iex
221 jsxa(i)=jsx
222 jexa(i)=jex
223
224 idsp(i)=isumm
225 isumm=isumm+icnt(i)
226 if(jsx .eq. 1 .or. jex .eq. jm) then
227 icnt2(i) = (iex-isx+1)
228 else
229 icnt2(i)=0
230 endif
231 idsp2(i)=isumm2
232 if(jsx .eq. 1 .or. jex .eq. jm) isumm2=isumm2+(iex-isx+1)
233
234! GWV Create send buffer for scatter. This is now needed because we are no
235! longer sending contiguous slices of the im,jm full state arrays to the
236! processors with scatter. Instead we are sending a slice of I and a slice of J
237! and so have to reshape the send buffer below to make it contiguous groups of
238! isx:iex,jsx:jex arrays
239
240 do jj=jsx,jex
241 do ii=isx,iex
242 ibuff(isum)=icoords(ii,jj)
243 rbuff(isum)=rcoords(ii,jj)
244 isum=isum+1
245 end do
246 end do
247
248 end do ! enddo of num_procs
249!
250! extraction limits -- set to two rows
251!
252 jsta_2l = max(jsta - 2, 1 )
253 jend_2u = min(jend + 2, jm )
254 if(modelname=='GFS') then
255 ista_2l=max(ista-2,0)
256 iend_2u=min(iend+2,im+1)
257 else
258 ista_2l=max(ista-2,1)
259 iend_2u=min(iend+2,im)
260 endif
261
262! special for c-grid v
263 jvend_2u = min(jend + 2, jm+1 )
264!
265! NEW neighbors
266
267 ileft = me - 1
268 iright = me + 1
269 iup=mpi_proc_null
270 idn=mpi_proc_null
271
272 !if(mod(me,numx) .eq. 0) print *,' LEFT POINT',me
273 !if(mod(me+1,numx) .eq. 0) print *,' RIGHT POINT',me
274 if(mod(me,numx) .eq. 0) ileft=mpi_proc_null
275 if(mod(me,numx) .eq. 0) ileftb=me+numx-1
276 if(mod(me+1,numx) .eq. 0 .or. me .eq. num_procs-1) iright=mpi_proc_null
277 if(mod(me+1,numx) .eq. 0 .or. me .eq. num_procs-1) irightb=me-numx+1
278 if(me .ge. numx) idn=me-numx
279 if(me+1 .le. num_procs-numx) iup=me+numx
280
281 !print 102,me,ileft,iright,iup,idn,num_procs,'GWVX BOUNDS'
282
283! allocate arrays
284
285 ibsize = ( (iend-ista) +1) * ( (jend-jsta)+1)
286 allocate(ibcoords(ista_2l:iend_2u,jsta_2l:jend_2u))
287 allocate(rbcoords(ista_2l:iend_2u,jsta_2l:jend_2u))
288 allocate(ibufs(ibsize))
289 allocate(rbufs(ibsize))
290 call mpi_scatterv(ibuff,icnt,idsp,mpi_integer &
291 ,ibufs,icnt(me),mpi_integer ,0,mpi_comm_world,j)
292 call mpi_scatterv(rbuff,icnt,idsp,mpi_real &
293 ,rbufs,icnt(me),mpi_real ,0,mpi_comm_world,j)
294
295!
296!GWV reshape the receive subdomain
297
298 isum=1
299 do j=jsta,jend
300 do i=ista,iend
301 ibcoords(i,j)=ibufs(isum)
302 rbcoords(i,j)=rbufs(isum)
303 isum=isum+1
304 end do
305 end do
306
307!GWV end reshape
308 do j=jsta,jend
309 do i=ista,iend
310 ii=ibcoords(i,j)/10000
311 jj=ibcoords( i,j)-(ii*10000)
312 if(ii .ne. i .or. jj .ne. j) then
313 print *,i,j,ii,jj,ibcoords(i,j),' GWVX FAIL '
314 else
315 continue
316 endif
317 end do
318 end do
319
320 allocate(ipoles(im,2),ipole(ista:iend))
321 allocate(rpoles(im,2),rpole(ista:iend))
322 ipole=9900000
323 ipoles=-999999999
324
325 do i=ista,iend
326 if(me .lt. num_procs/2. .and. jsta_2l .le. 1 .and. icnt2(me) .gt. 0) ipole(i)=ibcoords(i,1)
327 if(me .lt. num_procs/2. .and. jsta_2l .le. 1 .and. icnt2(me) .gt. 0) rpole(i)=rbcoords(i,1)
328 if(me .gt. num_procs/2. .and. jend_2u .ge. jm .and. icnt2(me) .gt. 0) ipole(i)=ibcoords(i,jm)
329 if(me .gt. num_procs/2. .and. jend_2u .ge. jm .and. icnt2(me) .gt. 0) rpole(i)=rbcoords(i,jm)
330
331! check code to be removed upon debugging
332 if(me .lt. num_procs/2. .and. jsx .eq. 1) then
333 continue
334 endif
335 if(me .gt. num_procs/2. .and. jend_2u .ge. jm) then
336 continue
337 endif
338 end do ! end check code
339
340! test pole gather
341 !print 105,' GWVX GATHER DISP ',icnt2(me),idsp2(me),me
342 105 format(a30,3i12)
343
344 call mpi_gatherv(ipole(ista),icnt2(me),mpi_integer, ipoles,icnt2,idsp2,mpi_integer,0,mpi_comm_world, ierr )
345 call mpi_gatherv(rpole(ista),icnt2(me),mpi_real , rpoles,icnt2,idsp2,mpi_real ,0,mpi_comm_world, ierr )
346
347 if(me .eq. 0) then
348 do j=1,2
349 do i=1,im
350 ii=rpoles(i,j)/4000
351 jj=rpoles(i,j) -ii*4000
352 if(ii .ne. i .or. jj .ne. 1 .and. jj .ne. jm ) then
353 write(*,169)' IPOLES BAD POINT',rpoles(i,j),ii,i,jj,' jm= ',jm
354 else
355 continue
356! write(*,169)' IPOLES GOOD POINT',rpoles(i,j),ii,i,jj,' jm= ',jm
357 endif
358 end do
359 end do
360 endif
361
362 107 format(a20,10i10)
363 169 format(a25,f20.1,3i10,a10,4i10)
364!
365! print *, ' me, jsta_2l, jend_2u = ',me,jsta_2l, jend_2u, &
366! 'jvend_2u=',jvend_2u,'im=',im,'jm=',jm,'lm=',lm, &
367! 'lp1=',lp1
368! write(*,'(A,5I10)') 'MPI_FIRST me,jsta,jend,ista,iend,=',me,jsta,jend,ista,iend
369
370 end
371
372! subroutine sub(a)
373! return
374! end
375
376
377
378 subroutine fullpole(a,rpoles)
379
380 use ctlblk_mod, only: num_procs, jend, iup, jsta, idn, mpi_comm_comp, im,modelname,numx,&
381 icoords,ibcoords,rbcoords,bufs,ibufs,me, &
382 jsta_2l,jend_2u,ileft,iright,ista_2l,iend_2u,ista,iend,jm,icnt2,idsp2
383!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
384 implicit none
385!
386 include 'mpif.h'
387!
388 real,intent(inout) :: a ( ista_2l:iend_2u,jsta_2l:jend_2u ),rpoles(im,2)
389 real, allocatable :: rpole(:)
390
391 integer status(MPI_STATUS_SIZE)
392 integer ierr
393 integer size,ubound,lbound
394 integer i,ii,jj, ibl,ibu,jbl,jbu,icc,jcc
395 integer ifirst
396 data ifirst/0/
397 integer iwest,ieast
398 data iwest,ieast/0,0/
399 allocate(rpole(ista:iend))
400
401 do i=ista,iend
402 if(me .lt. num_procs/2. .and. jsta_2l .le. 1 .and. icnt2(me) .gt. 0) rpole(i)=a(i,1)
403 if(me .ge. num_procs/2. .and. jend_2u .ge. jm .and. icnt2(me) .gt. 0) rpole(i)=a(i,jm)
404 end do
405
406 call mpi_allgatherv(rpole(ista),icnt2(me),mpi_real,rpoles,icnt2,idsp2,mpi_real, mpi_comm_comp,ierr)
407
408 call mpi_barrier(mpi_comm_comp,ierr)
409 ifirst=1
410
411 end
412
subroutine mpi_first()
SUBPROGRAM: MPI_FIRST SET UP MESSGAE PASSING INFO PRGRMMR: TUCCILLO ORG: IBM.
Definition MPI_FIRST.f:41