Changeset c8fc724 in flexpart.git for src/FLEXPART_MPI.f90


Ignore:
Timestamp:
Feb 2, 2017, 2:23:04 PM (7 years ago)
Author:
Espen Sollum ATMOS <eso@…>
Branches:
master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug, univie
Children:
93786a1
Parents:
4c64400
Message:

Bugfix: using namelist species files would sometimes set incorrect below-cloud scavenging parameters when using more than 1 species. Update: per-species scavenging statistics are written out at end of simulation.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/FLEXPART_MPI.f90

    r4c64400 rc8fc724  
    462462! NIK 16.02.2005
    463463  if (lroot) then
    464     call MPI_Reduce(MPI_IN_PLACE, tot_blc_count, 1, MPI_INTEGER8, MPI_SUM, id_root, &
     464    call MPI_Reduce(MPI_IN_PLACE, tot_blc_count, nspec, MPI_INTEGER8, MPI_SUM, id_root, &
    465465         & mp_comm_used, mp_ierr)
    466     call MPI_Reduce(MPI_IN_PLACE, tot_inc_count, 1, MPI_INTEGER8, MPI_SUM, id_root, &
     466    call MPI_Reduce(MPI_IN_PLACE, tot_inc_count, nspec, MPI_INTEGER8, MPI_SUM, id_root, &
    467467         & mp_comm_used, mp_ierr)
    468468  else
    469469    if (mp_partgroup_pid.ge.0) then ! Skip for readwind process
    470       call MPI_Reduce(tot_blc_count, 0, 1, MPI_INTEGER8, MPI_SUM, id_root, &
     470      call MPI_Reduce(tot_blc_count, 0, nspec, MPI_INTEGER8, MPI_SUM, id_root, &
    471471           & mp_comm_used, mp_ierr)
    472       call MPI_Reduce(tot_inc_count, 0, 1, MPI_INTEGER8, MPI_SUM, id_root, &
     472      call MPI_Reduce(tot_inc_count, 0, nspec, MPI_INTEGER8, MPI_SUM, id_root, &
    473473           & mp_comm_used, mp_ierr)
    474474    end if
     
    476476
    477477  if (lroot) then
    478     write(*,*) '**********************************************'
    479     write(*,*) 'Total number of occurences of below-cloud scavenging', &
    480          & tot_blc_count
    481     write(*,*) 'Total number of occurences of in-cloud    scavenging', &
    482          & tot_inc_count
    483     write(*,*) '**********************************************'
     478    do i=1,nspec
     479      write(*,*) '**********************************************'
     480      write(*,*) 'Scavenging statistics for species ', species(i), ':'
     481      write(*,*) 'Total number of occurences of below-cloud scavenging', &
     482           & tot_blc_count(i)
     483      write(*,*) 'Total number of occurences of in-cloud    scavenging', &
     484           & tot_inc_count(i)
     485      write(*,*) '**********************************************'
     486    end do
    484487
    485488    write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE&
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG