Changeset 328fdf9 in flexpart.git for src/init_domainfill_mpi.f90
- Timestamp:
- Apr 28, 2019, 7:53:42 PM (5 years ago)
- Branches:
- master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug
- Children:
- 0a94e13
- Parents:
- 77783e3
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
src/init_domainfill_mpi.f90
rb5127f9 r328fdf9 110 110 endif 111 111 112 ! Exit here if resuming a run from particle dump 113 !*********************************************** 114 if (gdomainfill.and.ipin.ne.0) return 115 112 116 ! Do not release particles twice (i.e., not at both in the leftmost and rightmost 113 117 ! grid cell) for a global domain … … 213 217 colmass(ix,jy)=(pp(1)-pp(nz))/ga*gridarea(jy) 214 218 colmasstotal=colmasstotal+colmass(ix,jy) 215 216 219 end do 217 220 end do … … 466 469 467 470 ! eso TODO: only needed for root process 468 if ( ipin.eq.1) then471 if ((ipin.eq.1).and.(.not.gdomainfill)) then 469 472 open(unitboundcond,file=path(2)(1:length(2))//'boundcond.bin', & 470 473 form='unformatted') … … 474 477 endif 475 478 476 numpart = numpart/mp_partgroup_np 477 if (mod(numpart,mp_partgroup_np).ne.0) numpart=numpart+1 478 479 else ! Allocate dummy arrays for receiving processes 480 allocate(itra1_tmp(nullsize),npoint_tmp(nullsize),nclass_tmp(nullsize),& 481 & idt_tmp(nullsize),itramem_tmp(nullsize),itrasplit_tmp(nullsize),& 482 & xtra1_tmp(nullsize),ytra1_tmp(nullsize),ztra1_tmp(nullsize),& 483 & xmass1_tmp(nullsize, nullsize)) 479 if (ipin.eq.0) then 480 numpart = numpart/mp_partgroup_np 481 if (mod(numpart,mp_partgroup_np).ne.0) numpart=numpart+1 482 end if 483 484 else ! Allocate dummy arrays for receiving processes 485 if (ipin.eq.0) then 486 allocate(itra1_tmp(nullsize),npoint_tmp(nullsize),nclass_tmp(nullsize),& 487 & idt_tmp(nullsize),itramem_tmp(nullsize),itrasplit_tmp(nullsize),& 488 & xtra1_tmp(nullsize),ytra1_tmp(nullsize),ztra1_tmp(nullsize),& 489 & xmass1_tmp(nullsize, nullsize)) 490 end if 484 491 485 end if ! end if(lroot) 492 end if ! end if(lroot) 493 486 494 487 495 488 496 ! Distribute particles to other processes (numpart is 'per-process', not total) 489 call MPI_Bcast(numpart, 1, MPI_INTEGER, id_root, mp_comm_used, mp_ierr) 490 ! eso TODO: xmassperparticle: not necessary to send 491 call MPI_Bcast(xmassperparticle, 1, mp_sp, id_root, mp_comm_used, mp_ierr)492 call mpif_send_part_properties(numpart)497 ! Only if not restarting from previous run 498 if (ipin.eq.0) then 499 call MPI_Bcast(numpart, 1, MPI_INTEGER, id_root, mp_comm_used, mp_ierr) 500 call mpif_send_part_properties(npart(1)/mp_partgroup_np) 493 501 494 502 ! Deallocate the temporary arrays used for all particles 495 deallocate(itra1_tmp,npoint_tmp,nclass_tmp,idt_tmp,itramem_tmp,&503 deallocate(itra1_tmp,npoint_tmp,nclass_tmp,idt_tmp,itramem_tmp,& 496 504 & itrasplit_tmp,xtra1_tmp,ytra1_tmp,ztra1_tmp,xmass1_tmp) 505 end if 497 506 498 507
Note: See TracChangeset
for help on using the changeset viewer.