Changeset c8fc724 in flexpart.git


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.

Location:
src
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • src/FLEXPART.f90

    r4c64400 rc8fc724  
    414414
    415415! NIK 16.02.2005
    416   write(*,*) '**********************************************'
    417   write(*,*) 'Total number of occurences of below-cloud scavenging', tot_blc_count
    418   write(*,*) 'Total number of occurences of in-cloud    scavenging', tot_inc_count
    419   write(*,*) '**********************************************'
    420 
     416  do i=1,nspec
     417    write(*,*) '**********************************************'
     418    write(*,*) 'Scavenging statistics for species ', species(i), ':'
     419    write(*,*) 'Total number of occurences of below-cloud scavenging', &
     420         & tot_blc_count(i)
     421    write(*,*) 'Total number of occurences of in-cloud    scavenging', &
     422         & tot_inc_count(i)
     423    write(*,*) '**********************************************'
     424  end do
     425 
    421426  write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE&
    422427       &XPART MODEL RUN!'
  • 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&
  • src/com_mod.f90

    r4c64400 rc8fc724  
    140140
    141141!NIK 16.02.2015
    142   integer(selected_int_kind(16)) :: tot_blc_count=0, tot_inc_count=0
     142  integer(selected_int_kind(16)), dimension(maxspec) :: tot_blc_count=0, &
     143       &tot_inc_count=0
    143144
    144145
     
    576577  real :: dxoutn,dyoutn,outlon0n,outlat0n,xoutshiftn,youtshiftn
    577578  !real outheight(maxzgrid),outheighthalf(maxzgrid)
    578   logical :: DEP,DRYDEP,DRYDEPSPEC(maxspec),WETDEP,OHREA,ASSSPEC
     579  logical :: DEP,DRYDEP,DRYDEPSPEC(maxspec),WETDEP,WETDEPSPEC(maxspec),&
     580       & OHREA,ASSSPEC
    579581
    580582  ! numxgrid,numygrid       number of grid points in x,y-direction
     
    593595  ! DRYDEPSPEC              .true., if dry deposition is switched on for that species
    594596  ! WETDEP                  .true., if wet deposition is switched on
     597  ! WETDEPSPEC              .true., if wet deposition is switched on for that species
    595598  ! OHREA                   .true., if OH reaction is switched on
    596599  ! ASSSPEC                 .true., if there are two species asscoiated
  • src/readreleases.f90

    r4c64400 rc8fc724  
    276276  do i=1,maxspec
    277277    DRYDEPSPEC(i)=.false.
     278    WETDEPSPEC(i)=.false.
    278279  end do
    279280
     
    369370         &(dquer(i).gt.0. .and. (crain_aero(i) .gt. 0. .or. csnow_aero(i).gt.0.)))  then
    370371      WETDEP=.true.
     372      WETDEPSPEC(i)=.true.
    371373      if (lroot) then
    372374        write (*,*) '  Below-cloud scavenging: ON'
     
    380382    if (dquer(i).gt.0..and.(ccn_aero(i).gt.0. .or. in_aero(i).gt.0.))  then
    381383      WETDEP=.true.
     384      WETDEPSPEC(i)=.true.
    382385      if (lroot) then
    383386        write (*,*) '  In-cloud scavenging: ON'
     
    399402    endif
    400403
    401   end do
     404  end do ! end loop over species
    402405
    403406  if (WETDEP.or.DRYDEP) DEP=.true.
  • src/readspecies.f90

    r62e65c7 rc8fc724  
    197197    weta_gas(pos_spec)=pweta_gas
    198198    wetb_gas(pos_spec)=pwetb_gas
    199     crain_aero=pcrain_aero
    200     csnow_aero=pcsnow_aero
     199    crain_aero(pos_spec)=pcrain_aero
     200    csnow_aero(pos_spec)=pcsnow_aero
    201201    ccn_aero(pos_spec)=pccn_aero
    202202    in_aero(pos_spec)=pin_aero
  • src/wetdepo.f90

    r05cf28d rc8fc724  
    9090  real :: frac_act, liq_frac, dquer_m
    9191
    92   integer :: blc_count, inc_count
     92  integer(selected_int_kind(16)), dimension(nspec) :: blc_count, inc_count
    9393  real    :: Si_dummy, wetscav_dummy
    9494  logical :: readclouds_this_nest
     
    107107!************************
    108108
    109   blc_count=0
    110   inc_count=0
     109  blc_count(:)=0
     110  inc_count(:)=0
    111111
    112112  do jpart=1,numpart
     
    256256    do ks=1,nspec      ! loop over species
    257257      wetdeposit(ks)=0.
    258       wetscav=0.   
     258      wetscav=0.
     259
     260! Cycle loop if wet deposition for the species is off
     261      if (WETDEPSPEC(ks).eqv..false.) cycle
    259262
    260263      if (ngrid.gt.0) then
     
    274277        if ((dquer(ks).le.0.).and.(weta_gas(ks).gt.0..or.wetb_gas(ks).gt.0.)) then
    275278          !        if (weta(ks).gt.0. .or. wetb(ks).gt.0.) then
    276           blc_count=blc_count+1
     279          blc_count(ks)=blc_count(ks)+1
    277280          wetscav=weta_gas(ks)*prec(1)**wetb_gas(ks)
    278281
     
    280283!*********************************************************************************
    281284        else if ((dquer(ks).gt.0.).and.(crain_aero(ks).gt.0..or.csnow_aero(ks).gt.0.)) then
    282           blc_count=blc_count+1
     285          blc_count(ks)=blc_count(ks)+1
    283286
    284287!NIK 17.02.2015
     
    320323! given in species file, or if gas and positive Henry's constant
    321324        if ((ccn_aero(ks).gt.0. .or. in_aero(ks).gt.0.).or.(henry(ks).gt.0.and.dquer(ks).le.0)) then
    322           inc_count=inc_count+1
     325          inc_count(ks)=inc_count(ks)+1
    323326! if negative coefficients (turned off) set to zero for use in equation
    324327          if (ccn_aero(ks).lt.0.) ccn_aero(ks)=0.
     
    432435
    433436! count the total number of below-cloud and in-cloud occurences:
    434   tot_blc_count=tot_blc_count+blc_count
    435   tot_inc_count=tot_inc_count+inc_count
     437  tot_blc_count(1:nspec)=tot_blc_count(1:nspec)+blc_count(1:nspec)
     438  tot_inc_count(1:nspec)=tot_inc_count(1:nspec)+inc_count(1:nspec)
    436439
    437440end subroutine wetdepo
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG