source: flexpart.git/src/boundcond_domainfill_mpi.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 27.9 KB
Line 
1subroutine boundcond_domainfill(itime,loutend)
2!                                  i      i
3!*****************************************************************************
4!                                                                            *
5! Particles are created by this subroutine continuously throughout the       *
6! simulation at the boundaries of the domain-filling box.                    *
7! All particles carry the same amount of mass which alltogether comprises the*
8! mass of air within the box, which remains (more or less) constant.         *
9!                                                                            *
10!     Author: A. Stohl                                                       *
11!                                                                            *
12!     16 October 2002                                                        *
13!                                                                            *
14!*****************************************************************************
15!                                                                            *
16! Variables:                                                                 *
17!                                                                            *
18! nx_we(2)       grid indices for western and eastern boundary of domain-    *
19!                filling trajectory calculations                             *
20! ny_sn(2)       grid indices for southern and northern boundary of domain-  *
21!                filling trajectory calculations                             *
22!                                                                            *
23!*****************************************************************************
24!  CHANGES                                                                   
25!    08/2016 eso: MPI version:                                               
26!
27!   -Root process release particles and distributes to other processes.
28!    Temporary arrays are used, also for the non-root (receiving) processes.
29!   -The scheme can be improved by having all processes report numpart
30!    (keeping track of how many particles have left the domain), so that
31!    a proportional amount of new particles can be distributed (however
32!    we have a separate function called from timemanager that will
33!    redistribute particles among processes if there are imbalances)     
34!*****************************************************************************
35
36  use point_mod
37  use par_mod
38  use com_mod
39  use random_mod, only: ran1
40  use mpi_mod
41
42  implicit none
43
44  real :: dz,dz1,dz2,dt1,dt2,dtt,ylat,xm,cosfact,accmasst
45  integer :: itime,in,indz,indzp,i,loutend
46  integer :: j,k,ix,jy,m,indzh,indexh,minpart,ipart,mmass
47  integer :: numactiveparticles, numpart_total, rel_counter
48  integer,allocatable,dimension(:) ::  numrel_mpi !, numactiveparticles_mpi
49
50  real :: windl(2),rhol(2)
51  real :: windhl(2),rhohl(2)
52  real :: windx,rhox
53  real :: deltaz,boundarea,fluxofmass
54
55  integer :: ixm,ixp,jym,jyp,indzm,mm
56  real :: pvpart,ddx,ddy,rddx,rddy,p1,p2,p3,p4,y1(2),yh1(2)
57
58  integer :: idummy = -11
59  integer :: mtag
60  logical :: first_call=.true.
61! Sizes of temporary arrays are maxpartfract*maxpart. Increase maxpartfract if
62! needed. 
63  real,parameter :: maxpartfract=0.1
64  integer :: tmp_size = int(maxpartfract*maxpart)
65
66! Use different seed for each process
67  if (first_call) then
68    idummy=idummy+mp_seed
69    first_call=.false.
70  end if
71
72
73! If domain-filling is global, no boundary conditions are needed
74!***************************************************************
75
76  if (gdomainfill) return
77
78  accmasst=0.
79  numactiveparticles=0
80! Keep track of active particles on each process
81  allocate(numrel_mpi(0:mp_partgroup_np-1))
82! numactiveparticles_mpi(0:mp_partgroup_np-1)
83
84! New particles to be released on each process
85  numrel_mpi(:)=0
86
87! Terminate trajectories that have left the domain, if domain-filling
88! trajectory calculation domain is not global. Done for all processes
89!********************************************************************
90
91  do i=1,numpart
92    if (itra1(i).eq.itime) then
93      if ((ytra1(i).gt.real(ny_sn(2))).or. &
94           (ytra1(i).lt.real(ny_sn(1)))) itra1(i)=-999999999
95      if (((.not.xglobal).or.(nx_we(2).ne.(nx-2))).and. &
96           ((xtra1(i).lt.real(nx_we(1))).or. &
97           (xtra1(i).gt.real(nx_we(2))))) itra1(i)=-999999999
98    endif
99    if (itra1(i).ne.-999999999) numactiveparticles= &
100         numactiveparticles+1
101  end do
102!  numactiveparticles_mpi(mp_partid) = numactiveparticles
103
104
105! Collect number of active particles from all processes
106  ! call MPI_Allgather(numactiveparticles, 1, MPI_INTEGER, &
107  !      &numactiveparticles_mpi, 1, MPI_INTEGER, mp_comm_used, mp_ierr)
108
109
110! Total number of new releases
111  numpart_total = 0
112
113
114! This section only done by root process
115!***************************************
116
117  if (lroot) then
118
119! Use separate arrays for newly released particles
120!*************************************************
121
122    allocate(itra1_tmp(tmp_size),npoint_tmp(tmp_size),nclass_tmp(tmp_size),&
123         & idt_tmp(tmp_size),itramem_tmp(tmp_size),itrasplit_tmp(tmp_size),&
124         & xtra1_tmp(tmp_size),ytra1_tmp(tmp_size),ztra1_tmp(tmp_size),&
125         & xmass1_tmp(tmp_size, maxspec))
126
127! Initialize all particles as non-existent   
128    itra1_tmp(:)=-999999999
129
130! Determine auxiliary variables for time interpolation
131!*****************************************************
132
133    dt1=real(itime-memtime(1))
134    dt2=real(memtime(2)-itime)
135    dtt=1./(dt1+dt2)
136
137! Initialize auxiliary variable used to search for vacant storage space
138!**********************************************************************
139
140    minpart=1
141
142!***************************************
143! Western and eastern boundary condition
144!***************************************
145
146! Loop from south to north
147!*************************
148
149    do jy=ny_sn(1),ny_sn(2)
150
151! Loop over western (index 1) and eastern (index 2) boundary
152!***********************************************************
153
154      do k=1,2
155
156! Loop over all release locations in a column
157!********************************************
158
159        do j=1,numcolumn_we(k,jy)
160
161! Determine, for each release location, the area of the corresponding boundary
162!*****************************************************************************
163
164          if (j.eq.1) then
165            deltaz=(zcolumn_we(k,jy,2)+zcolumn_we(k,jy,1))/2.
166          else if (j.eq.numcolumn_we(k,jy)) then
167!        deltaz=height(nz)-(zcolumn_we(k,jy,j-1)+
168!    +        zcolumn_we(k,jy,j))/2.
169! In order to avoid taking a very high column for very many particles,
170! use the deltaz from one particle below instead
171            deltaz=(zcolumn_we(k,jy,j)-zcolumn_we(k,jy,j-2))/2.
172          else
173            deltaz=(zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))/2.
174          endif
175          if ((jy.eq.ny_sn(1)).or.(jy.eq.ny_sn(2))) then
176            boundarea=deltaz*111198.5/2.*dy
177          else
178            boundarea=deltaz*111198.5*dy
179          endif
180
181
182! Interpolate the wind velocity and density to the release location
183!******************************************************************
184
185! Determine the model level below the release position
186!*****************************************************
187
188          do i=2,nz
189            if (height(i).gt.zcolumn_we(k,jy,j)) then
190              indz=i-1
191              indzp=i
192              goto 6
193            endif
194          end do
1956         continue
196
197! Vertical distance to the level below and above current position
198!****************************************************************
199
200          dz1=zcolumn_we(k,jy,j)-height(indz)
201          dz2=height(indzp)-zcolumn_we(k,jy,j)
202          dz=1./(dz1+dz2)
203
204! Vertical and temporal interpolation
205!************************************
206
207          do m=1,2
208            indexh=memind(m)
209            do in=1,2
210              indzh=indz+in-1
211              windl(in)=uu(nx_we(k),jy,indzh,indexh)
212              rhol(in)=rho(nx_we(k),jy,indzh,indexh)
213            end do
214
215            windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz
216            rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz
217          end do
218
219          windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt
220          rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt
221
222! Calculate mass flux
223!********************
224
225          fluxofmass=windx*rhox*boundarea*real(lsynctime)
226
227
228! If the mass flux is directed into the domain, add it to previous mass fluxes;
229! if it is out of the domain, set accumulated mass flux to zero
230!******************************************************************************
231
232          if (k.eq.1) then
233            if (fluxofmass.ge.0.) then
234              acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+fluxofmass
235            else
236              acc_mass_we(k,jy,j)=0.
237            endif
238          else
239            if (fluxofmass.le.0.) then
240              acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)+abs(fluxofmass)
241            else
242              acc_mass_we(k,jy,j)=0.
243            endif
244          endif
245          accmasst=accmasst+acc_mass_we(k,jy,j)
246
247! If the accumulated mass exceeds half the mass that each particle shall carry,
248! one (or more) particle(s) is (are) released and the accumulated mass is
249! reduced by the mass of this (these) particle(s)
250!******************************************************************************
251
252          if (acc_mass_we(k,jy,j).ge.xmassperparticle/2.) then
253            mmass=int((acc_mass_we(k,jy,j)+xmassperparticle/2.)/ &
254                 xmassperparticle)
255            acc_mass_we(k,jy,j)=acc_mass_we(k,jy,j)- &
256                 real(mmass)*xmassperparticle
257          else
258            mmass=0
259          endif
260
261          do m=1,mmass
262            do ipart=minpart,maxpart
263
264! If a vacant storage space is found, attribute everything to this array element
265! TODO: for the MPI version this test can be removed, as all
266!       elements in _tmp arrays are initialized to zero
267!*****************************************************************************
268
269              if (itra1_tmp(ipart).ne.itime) then
270
271! Assign particle positions
272!**************************
273
274                xtra1_tmp(ipart)=real(nx_we(k))
275                if (jy.eq.ny_sn(1)) then
276                  ytra1_tmp(ipart)=real(jy)+0.5*ran1(idummy)
277                else if (jy.eq.ny_sn(2)) then
278                  ytra1_tmp(ipart)=real(jy)-0.5*ran1(idummy)
279                else
280                  ytra1_tmp(ipart)=real(jy)+(ran1(idummy)-.5)
281                endif
282                if (j.eq.1) then
283                  ztra1_tmp(ipart)=zcolumn_we(k,jy,1)+(zcolumn_we(k,jy,2)- &
284                       zcolumn_we(k,jy,1))/4.
285                else if (j.eq.numcolumn_we(k,jy)) then
286                  ztra1_tmp(ipart)=(2.*zcolumn_we(k,jy,j)+ &
287                       zcolumn_we(k,jy,j-1)+height(nz))/4.
288                else
289                  ztra1_tmp(ipart)=zcolumn_we(k,jy,j-1)+ran1(idummy)* &
290                       (zcolumn_we(k,jy,j+1)-zcolumn_we(k,jy,j-1))
291                endif
292
293! Interpolate PV to the particle position
294!****************************************
295                ixm=int(xtra1_tmp(ipart))
296                jym=int(ytra1_tmp(ipart))
297                ixp=ixm+1
298                jyp=jym+1
299                ddx=xtra1_tmp(ipart)-real(ixm)
300                ddy=ytra1_tmp(ipart)-real(jym)
301                rddx=1.-ddx
302                rddy=1.-ddy
303                p1=rddx*rddy
304                p2=ddx*rddy
305                p3=rddx*ddy
306                p4=ddx*ddy
307                do i=2,nz
308                  if (height(i).gt.ztra1_tmp(ipart)) then
309                    indzm=i-1
310                    indzp=i
311                    goto 26
312                  endif
313                end do
31426              continue
315                dz1=ztra1_tmp(ipart)-height(indzm)
316                dz2=height(indzp)-ztra1_tmp(ipart)
317                dz=1./(dz1+dz2)
318                do mm=1,2
319                  indexh=memind(mm)
320                  do in=1,2
321                    indzh=indzm+in-1
322                    y1(in)=p1*pv(ixm,jym,indzh,indexh) &
323                         +p2*pv(ixp,jym,indzh,indexh) &
324                         +p3*pv(ixm,jyp,indzh,indexh) &
325                         +p4*pv(ixp,jyp,indzh,indexh)
326                  end do
327                  yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz
328                end do
329                pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt
330                ylat=ylat0+ytra1_tmp(ipart)*dy
331                if (ylat.lt.0.) pvpart=-1.*pvpart
332
333
334! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere
335!*****************************************************************************
336
337                if (((ztra1_tmp(ipart).gt.3000.).and. &
338                     (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then
339                  nclass_tmp(ipart)=min(int(ran1(idummy)* &
340                       real(nclassunc))+1,nclassunc)
341                  numactiveparticles=numactiveparticles+1
342                  numparticlecount=numparticlecount+1
343                  npoint_tmp(ipart)=numparticlecount
344                  idt_tmp(ipart)=mintime
345                  itra1_tmp(ipart)=itime
346                  itramem_tmp(ipart)=itra1_tmp(ipart)
347                  itrasplit_tmp(ipart)=itra1_tmp(ipart)+ldirect*itsplit
348                  xmass1_tmp(ipart,1)=xmassperparticle
349                  if (mdomainfill.eq.2) xmass1_tmp(ipart,1)= &
350                       xmass1_tmp(ipart,1)*pvpart*48./29.*ozonescale/10.**9
351                else
352                  goto 71
353                endif
354
355
356! Increase numpart, if necessary
357!*******************************
358
359                numpart_total=max(numpart_total,ipart)
360                goto 73      ! Storage space has been found, stop searching
361              endif
362            end do
363            if (ipart.gt.tmp_size) &
364                 stop 'boundcond_domainfill_mpi.f90: too many particles required'
36573          minpart=ipart+1
36671          continue
367          end do
368
369
370        end do
371      end do
372    end do
373
374
375!*****************************************
376! Southern and northern boundary condition
377!*****************************************
378
379! Loop from west to east
380!***********************
381
382    do ix=nx_we(1),nx_we(2)
383
384! Loop over southern (index 1) and northern (index 2) boundary
385!*************************************************************
386
387      do k=1,2
388        ylat=ylat0+real(ny_sn(k))*dy
389        cosfact=cos(ylat*pi180)
390
391! Loop over all release locations in a column
392!********************************************
393
394        do j=1,numcolumn_sn(k,ix)
395
396! Determine, for each release location, the area of the corresponding boundary
397!*****************************************************************************
398
399          if (j.eq.1) then
400            deltaz=(zcolumn_sn(k,ix,2)+zcolumn_sn(k,ix,1))/2.
401          else if (j.eq.numcolumn_sn(k,ix)) then
402!        deltaz=height(nz)-(zcolumn_sn(k,ix,j-1)+
403!    +        zcolumn_sn(k,ix,j))/2.
404! In order to avoid taking a very high column for very many particles,
405! use the deltaz from one particle below instead
406            deltaz=(zcolumn_sn(k,ix,j)-zcolumn_sn(k,ix,j-2))/2.
407          else
408            deltaz=(zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))/2.
409          endif
410          if ((ix.eq.nx_we(1)).or.(ix.eq.nx_we(2))) then
411            boundarea=deltaz*111198.5/2.*cosfact*dx
412          else
413            boundarea=deltaz*111198.5*cosfact*dx
414          endif
415
416
417! Interpolate the wind velocity and density to the release location
418!******************************************************************
419
420! Determine the model level below the release position
421!*****************************************************
422
423          do i=2,nz
424            if (height(i).gt.zcolumn_sn(k,ix,j)) then
425              indz=i-1
426              indzp=i
427              goto 16
428            endif
429          end do
43016        continue
431
432! Vertical distance to the level below and above current position
433!****************************************************************
434
435          dz1=zcolumn_sn(k,ix,j)-height(indz)
436          dz2=height(indzp)-zcolumn_sn(k,ix,j)
437          dz=1./(dz1+dz2)
438
439! Vertical and temporal interpolation
440!************************************
441
442          do m=1,2
443            indexh=memind(m)
444            do in=1,2
445              indzh=indz+in-1
446              windl(in)=vv(ix,ny_sn(k),indzh,indexh)
447              rhol(in)=rho(ix,ny_sn(k),indzh,indexh)
448            end do
449
450            windhl(m)=(dz2*windl(1)+dz1*windl(2))*dz
451            rhohl(m)=(dz2*rhol(1)+dz1*rhol(2))*dz
452          end do
453
454          windx=(windhl(1)*dt2+windhl(2)*dt1)*dtt
455          rhox=(rhohl(1)*dt2+rhohl(2)*dt1)*dtt
456
457! Calculate mass flux
458!********************
459
460          fluxofmass=windx*rhox*boundarea*real(lsynctime)
461
462! If the mass flux is directed into the domain, add it to previous mass fluxes;
463! if it is out of the domain, set accumulated mass flux to zero
464!******************************************************************************
465
466          if (k.eq.1) then
467            if (fluxofmass.ge.0.) then
468              acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+fluxofmass
469            else
470              acc_mass_sn(k,ix,j)=0.
471            endif
472          else
473            if (fluxofmass.le.0.) then
474              acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)+abs(fluxofmass)
475            else
476              acc_mass_sn(k,ix,j)=0.
477            endif
478          endif
479          accmasst=accmasst+acc_mass_sn(k,ix,j)
480
481! If the accumulated mass exceeds half the mass that each particle shall carry,
482! one (or more) particle(s) is (are) released and the accumulated mass is
483! reduced by the mass of this (these) particle(s)
484!******************************************************************************
485
486          if (acc_mass_sn(k,ix,j).ge.xmassperparticle/2.) then
487            mmass=int((acc_mass_sn(k,ix,j)+xmassperparticle/2.)/ &
488                 xmassperparticle)
489            acc_mass_sn(k,ix,j)=acc_mass_sn(k,ix,j)- &
490                 real(mmass)*xmassperparticle
491          else
492            mmass=0
493          endif
494
495          do m=1,mmass
496            do ipart=minpart,maxpart
497
498! If a vacant storage space is found, attribute everything to this array element
499!*****************************************************************************
500
501              if (itra1_tmp(ipart).ne.itime) then
502
503! Assign particle positions
504!**************************
505
506                ytra1_tmp(ipart)=real(ny_sn(k))
507                if (ix.eq.nx_we(1)) then
508                  xtra1_tmp(ipart)=real(ix)+0.5*ran1(idummy)
509                else if (ix.eq.nx_we(2)) then
510                  xtra1_tmp(ipart)=real(ix)-0.5*ran1(idummy)
511                else
512                  xtra1_tmp(ipart)=real(ix)+(ran1(idummy)-.5)
513                endif
514                if (j.eq.1) then
515                  ztra1_tmp(ipart)=zcolumn_sn(k,ix,1)+(zcolumn_sn(k,ix,2)- &
516                       zcolumn_sn(k,ix,1))/4.
517                else if (j.eq.numcolumn_sn(k,ix)) then
518                  ztra1_tmp(ipart)=(2.*zcolumn_sn(k,ix,j)+ &
519                       zcolumn_sn(k,ix,j-1)+height(nz))/4.
520                else
521                  ztra1_tmp(ipart)=zcolumn_sn(k,ix,j-1)+ran1(idummy)* &
522                       (zcolumn_sn(k,ix,j+1)-zcolumn_sn(k,ix,j-1))
523                endif
524
525
526! Interpolate PV to the particle position
527!****************************************
528                ixm=int(xtra1_tmp(ipart))
529                jym=int(ytra1_tmp(ipart))
530                ixp=ixm+1
531                jyp=jym+1
532                ddx=xtra1_tmp(ipart)-real(ixm)
533                ddy=ytra1_tmp(ipart)-real(jym)
534                rddx=1.-ddx
535                rddy=1.-ddy
536                p1=rddx*rddy
537                p2=ddx*rddy
538                p3=rddx*ddy
539                p4=ddx*ddy
540                do i=2,nz
541                  if (height(i).gt.ztra1_tmp(ipart)) then
542                    indzm=i-1
543                    indzp=i
544                    goto 126
545                  endif
546                end do
547126             continue
548                dz1=ztra1_tmp(ipart)-height(indzm)
549                dz2=height(indzp)-ztra1_tmp(ipart)
550                dz=1./(dz1+dz2)
551                do mm=1,2
552                  indexh=memind(mm)
553                  do in=1,2
554                    indzh=indzm+in-1
555                    y1(in)=p1*pv(ixm,jym,indzh,indexh) &
556                         +p2*pv(ixp,jym,indzh,indexh) &
557                         +p3*pv(ixm,jyp,indzh,indexh) &
558                         +p4*pv(ixp,jyp,indzh,indexh)
559                  end do
560                  yh1(mm)=(dz2*y1(1)+dz1*y1(2))*dz
561                end do
562                pvpart=(yh1(1)*dt2+yh1(2)*dt1)*dtt
563                if (ylat.lt.0.) pvpart=-1.*pvpart
564
565
566! For domain-filling option 2 (stratospheric O3), do the rest only in the stratosphere
567!*****************************************************************************
568
569                if (((ztra1_tmp(ipart).gt.3000.).and. &
570                     (pvpart.gt.pvcrit)).or.(mdomainfill.eq.1)) then
571                  nclass_tmp(ipart)=min(int(ran1(idummy)* &
572                       real(nclassunc))+1,nclassunc)
573                  numactiveparticles=numactiveparticles+1
574                  numparticlecount=numparticlecount+1
575                  npoint_tmp(ipart)=numparticlecount
576                  idt_tmp(ipart)=mintime
577                  itra1_tmp(ipart)=itime
578                  itramem_tmp(ipart)=itra1_tmp(ipart)
579                  itrasplit_tmp(ipart)=itra1_tmp(ipart)+ldirect*itsplit
580                  xmass1_tmp(ipart,1)=xmassperparticle
581                  if (mdomainfill.eq.2) xmass1_tmp(ipart,1)= &
582                       xmass1_tmp(ipart,1)*pvpart*48./29.*ozonescale/10.**9
583                else
584                  goto 171
585                endif
586
587
588! Increase numpart, if necessary
589!*******************************
590                numpart_total=max(numpart_total,ipart)
591                goto 173      ! Storage space has been found, stop searching
592              endif
593            end do
594            if (ipart.gt.tmp_size) &
595                 stop 'boundcond_domainfill.f: too many particles required'
596173         minpart=ipart+1
597171         continue
598          end do
599
600
601        end do
602      end do
603    end do
604
605
606! xm=0.
607! do i=1,numpart_total
608!   if (itra1_tmp(i).eq.itime) xm=xm+xmass1(i,1)
609! end do
610
611!write(*,*) itime,numactiveparticles,numparticlecount,numpart,
612!    +xm,accmasst,xm+accmasst
613
614  end if ! if lroot
615
616! Distribute the number of particles to be released
617! *************************************************
618  call MPI_Bcast(numpart_total, 1, MPI_INTEGER, id_root, mp_comm_used, mp_ierr)
619
620  do i=0, mp_partgroup_np-1
621    numrel_mpi(i) = numpart_total/mp_partgroup_np
622    if (i.lt.mod(numpart_total,mp_partgroup_np)) numrel_mpi(i) = numrel_mpi(i) + 1
623  end do
624
625! Allocate temporary arrays for receiving processes
626  if (.not.lroot) then
627    allocate(itra1_tmp(numrel_mpi(mp_partid)),&
628         & npoint_tmp(numrel_mpi(mp_partid)),&
629         & nclass_tmp(numrel_mpi(mp_partid)),&
630         & idt_tmp(numrel_mpi(mp_partid)),&
631         & itramem_tmp(numrel_mpi(mp_partid)),&
632         & itrasplit_tmp(numrel_mpi(mp_partid)),&
633         & xtra1_tmp(numrel_mpi(mp_partid)),&
634         & ytra1_tmp(numrel_mpi(mp_partid)),&
635         & ztra1_tmp(numrel_mpi(mp_partid)),&
636         & xmass1_tmp(numrel_mpi(mp_partid),maxspec))
637     
638! Initialize all particles as non-existent   
639    itra1_tmp(:)=-999999999
640  end if
641
642! Distribute particles
643! Keep track of released particles so far
644  rel_counter = 0
645  mtag = 1000
646
647  do i=0, mp_partgroup_np-1
648
649! For root process, nothing to do except update release count
650    if (i.eq.0) then
651      rel_counter = rel_counter + numrel_mpi(i)
652      cycle
653    end if
654     
655! Send particles from root to non-root processes
656    if (lroot.and.numrel_mpi(i).gt.0) then
657
658      call MPI_SEND(nclass_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
659           &numrel_mpi(i),MPI_INTEGER,i,mtag+1*i,mp_comm_used,mp_ierr)
660
661      call MPI_SEND(npoint_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
662           &numrel_mpi(i),MPI_INTEGER,i,mtag+2*i,mp_comm_used,mp_ierr)
663
664      call MPI_SEND(itra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
665           &numrel_mpi(i),MPI_INTEGER,i,mtag+3*i,mp_comm_used,mp_ierr)
666
667      call MPI_SEND(idt_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
668           &numrel_mpi(i),MPI_INTEGER,i,mtag+4*i,mp_comm_used,mp_ierr)
669
670      call MPI_SEND(itramem_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
671           &numrel_mpi(i),MPI_INTEGER,i,mtag+5*i,mp_comm_used,mp_ierr)
672
673      call MPI_SEND(itrasplit_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
674           &numrel_mpi(i),MPI_INTEGER,i,mtag+6*i,mp_comm_used,mp_ierr)
675
676      call MPI_SEND(xtra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
677           &numrel_mpi(i),mp_dp,i,mtag+7*i,mp_comm_used,mp_ierr)
678
679      call MPI_SEND(ytra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
680           &numrel_mpi(i),mp_dp,i,mtag+8*i,mp_comm_used,mp_ierr)
681
682      call MPI_SEND(ztra1_tmp(rel_counter+1:rel_counter+numrel_mpi(i)),&
683           &numrel_mpi(i),mp_sp,i,mtag+9*i,mp_comm_used,mp_ierr)
684
685      do j=1,nspec
686        call MPI_SEND(xmass1_tmp(rel_counter+1:rel_counter+numrel_mpi(i),j),&
687             &numrel_mpi(i),mp_sp,i,mtag+(9+j)*i,mp_comm_used,mp_ierr)
688      end do
689
690! Non-root processes issue receive requests
691    else if (i.eq.mp_partid.and.numrel_mpi(i).gt.0) then
692      call MPI_RECV(nclass_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
693           &MPI_INTEGER,id_root,mtag+1*i,mp_comm_used,mp_status,mp_ierr)
694
695      call MPI_RECV(npoint_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
696           &MPI_INTEGER,id_root,mtag+2*i,mp_comm_used,mp_status,mp_ierr)
697
698      call MPI_RECV(itra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
699           &MPI_INTEGER,id_root,mtag+3*i,mp_comm_used,mp_status,mp_ierr)
700
701      call MPI_RECV(idt_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
702           &MPI_INTEGER,id_root,mtag+4*i,mp_comm_used,mp_status,mp_ierr)
703
704      call MPI_RECV(itramem_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
705           &MPI_INTEGER,id_root,mtag+5*i,mp_comm_used,mp_status,mp_ierr)
706
707      call MPI_RECV(itrasplit_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
708           &MPI_INTEGER,id_root,mtag+6*i,mp_comm_used,mp_status,mp_ierr)
709
710      call MPI_RECV(xtra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
711           &mp_dp,id_root,mtag+7*i,mp_comm_used,mp_status,mp_ierr)
712
713      call MPI_RECV(ytra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
714           &mp_dp,id_root,mtag+8*i,mp_comm_used,mp_status,mp_ierr)
715
716      call MPI_RECV(ztra1_tmp(1:numrel_mpi(i)),numrel_mpi(i),&
717           &mp_sp,id_root,mtag+9*i,mp_comm_used,mp_status,mp_ierr)
718
719      do j=1,nspec
720        call MPI_RECV(xmass1_tmp(1:numrel_mpi(i),j),numrel_mpi(i),&
721             &mp_sp,id_root,mtag+(9+j)*i,mp_comm_used,mp_status,mp_ierr)
722
723      end do
724    end if
725    rel_counter = rel_counter + numrel_mpi(i)
726  end do
727
728! Find free storage space for the new particles.
729! This section is independent of the redistribution scheme used
730! ********************************************************************
731
732! Keep track of released particles so far
733  minpart=1
734   
735! The algorithm should be correct also for root process
736  do i=1, numrel_mpi(mp_partid)
737    do ipart=minpart, maxpart
738      if (itra1(ipart).ne.itime) then
739        itra1(ipart) = itra1_tmp(i)
740        npoint(ipart) = npoint_tmp(i)
741        nclass(ipart) = nclass_tmp(i)
742        idt(ipart) = idt_tmp(i)
743        itramem(ipart) = itramem_tmp(i)
744        itrasplit(ipart) = itrasplit_tmp(i)
745        xtra1(ipart) = xtra1_tmp(i)
746        ytra1(ipart) = ytra1_tmp(i)
747        ztra1(ipart) = ztra1_tmp(i)
748        xmass1(ipart,:) = xmass1_tmp(i,:)
749! Increase numpart, if necessary
750        numpart=max(numpart,ipart)
751        goto 200 ! Storage space has been found, stop searching
752      end if
753    end do
754200 minpart=ipart+1
755  end do
756
757! If particles shall be dumped, then accumulated masses at the domain boundaries
758! must be dumped, too, to be used for later runs
759!*****************************************************************************
760
761  if ((ipout.gt.0).and.(itime.eq.loutend)) then
762    if (lroot) then
763      call mpif_mtime('iotime',0)
764      open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', &
765           form='unformatted')
766      write(unitboundcond) numcolumn_we,numcolumn_sn, &
767           zcolumn_we,zcolumn_sn,acc_mass_we,acc_mass_sn
768      close(unitboundcond)
769      call mpif_mtime('iotime',1)
770    end if
771  endif
772
773! Deallocate temporary arrays
774  deallocate(itra1_tmp,npoint_tmp,nclass_tmp,idt_tmp,itramem_tmp,itrasplit_tmp,&
775       & xtra1_tmp,ytra1_tmp,ztra1_tmp,xmass1_tmp,numrel_mpi)
776! numactiveparticles_mpi
777
778
779end subroutine boundcond_domainfill
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG