Changeset 32b49c3 in flexpart.git


Ignore:
Timestamp:
Apr 13, 2016, 11:57:08 AM (8 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:
341f4b7
Parents:
e31b3b5
Message:

Parallel version can now save/restart simulations with IPOUT/IPIN

Location:
src
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • src/FLEXPART_MPI.f90

    rdb712a8 r32b49c3  
    317317      print*,'call readpartpositions'
    318318    endif
    319     call readpartpositions
     319    ! readwind process skips this step
     320    if (lmp_use_reader.and..not.lmpreader) call readpartpositions
    320321  else
    321322    if (verbosity.gt.0 .and. lroot) then
  • src/makefile

    r1c3c778 r32b49c3  
    1414#
    1515#  USAGE
    16 #    Compile serial FLEXPART (ECMWF)a
     16#    Compile serial FLEXPART (ECMWF)
    1717#      make [-j] ecmwf
    1818#
  • src/mpi_mod.f90

    re31b3b5 r32b49c3  
    165165!***********************************************************************
    166166    use par_mod, only: maxpart, numwfmem, dep_prec
    167     use com_mod, only: mpi_mode
     167    use com_mod, only: mpi_mode, verbosity
    168168
    169169    implicit none
     
    386386! Dynamic allocation of arrays (in serial code these have size maxpart)
    387387!
     388! Not in use anymore, moved to com_mod for interoperability with serial version
     389!
     390! TODO: remove
    388391!*******************************************************************************
    389392    use com_mod
     
    10431046!***********************************************************************
    10441047    use com_mod
    1045     use par_mod, only: numwfmem
     1048    use par_mod,only: numwfmem
    10461049
    10471050    implicit none
     
    12261229    integer :: d2s1 = nxmax*nymax
    12271230    integer :: d2s2 = nxmax*nymax*maxspec
    1228     integer :: d1s1 = maxwf
     1231!    integer :: d1s1 = maxwf
    12291232
    12301233!*******************************************************************************
     
    14301433! last read was asynchronous to index 2, 3 is free
    14311434      mind=3
     1435    else
     1436! illegal state
     1437      write(*,*) 'mpi_mod> FLEXPART ERROR: Illegal memstat value. Exiting.'
     1438      stop
    14321439    end if
    14331440
     
    16021609!
    16031610!*******************************************************************************
    1604     use com_mod, only: readclouds
     1611!    use com_mod,only: readclouds
    16051612
    16061613    implicit none
    16071614
    16081615
    1609     integer :: n_req
    1610     integer :: i
     1616    integer :: n_req !,i
    16111617
    16121618!***********************************************************************
     
    19581964    implicit none
    19591965
    1960     integer :: ip,j,r
     1966    integer :: ip !,j,r
    19611967
    19621968!***********************************************************************
  • src/par_mod.f90

    r8ed5f11 r32b49c3  
    186186  !**************************************************
    187187
    188   integer,parameter :: maxpart=10000
     188  integer,parameter :: maxpart=4000000
    189189  integer,parameter :: maxspec=1
    190190  real,parameter :: minmass=0.0001
     
    229229  !*********************************
    230230
    231 !  integer,parameter :: maxrand=120000000
    232231  integer,parameter :: maxrand=200000000
    233 !
     232
    234233  ! maxrand                 number of random numbers used
    235234 
  • src/readpartpositions_mpi.f90

    r5f9d14a r32b49c3  
    4545  use com_mod
    4646  use random_mod, only: ran1
    47   use mpi_mod, only: mp_seed
     47  use mpi_mod !, only: mp_seed
    4848
    4949  implicit none
    5050
    51   integer :: ibdatein,ibtimein,nspecin,itimein,numpointin,i,j,ix
     51  integer :: ibdatein,ibtimein,nspecin,itimein,numpointin,i,j,ix,ip
    5252  integer :: id1,id2,it1,it2
     53  integer :: addone,numparticlecount_all,numpart_all,lbnd,ubnd
    5354  real :: xlonin,ylatin,topo,hmixi,pvi,qvi,rhoi,tri,tti
    5455  character :: specin*7
     
    5758  integer :: idummy = -8
    5859
     60  ! These variables are allocated at the root process for all particles in file.
     61  real,dimension(maxpart) :: xtra1_all,ytra1_all,ztra1_all
     62  real,dimension(maxpart,maxspec) :: xmass1_all
     63  integer,dimension(maxpart) :: npoint_all,itramem_all
     64
    5965  ! Different seed for each process
    6066  idummy=idummy+mp_seed
    6167
    6268  numparticlecount=0
     69  numparticlecount_all=0
     70  numpart_all=0
    6371
    6472  ! Open header file of dumped particle data
    65   !*****************************************
    66 
    67   open(unitpartin,file=path(2)(1:length(2))//'header', &
    68        form='unformatted',err=998)
    69 
    70   read(unitpartin) ibdatein,ibtimein
    71   read(unitpartin)
    72   read(unitpartin)
    73 
    74   read(unitpartin)
    75   read(unitpartin)
    76   read(unitpartin) nspecin
    77   nspecin=nspecin/3
    78   if ((ldirect.eq.1).and.(nspec.ne.nspecin)) goto 997
    79 
    80   do i=1,nspecin
    81     read(unitpartin)
    82     read(unitpartin)
    83     read(unitpartin) j,specin
    84     if ((ldirect.eq.1).and.(species(i)(1:7).ne.specin)) goto 996
     73  ! Each MPI process sequentially access file (just in case)
     74  !*********************************************************
     75
     76  do ip=0, mp_partgroup_np-1
     77    call mpif_mpi_barrier
     78    open(unitpartin,file=path(2)(1:length(2))//'header', &
     79         form='unformatted',err=998)
     80
     81    read(unitpartin) ibdatein,ibtimein
     82    read(unitpartin)
     83    read(unitpartin)
     84
     85    read(unitpartin)
     86    read(unitpartin)
     87    read(unitpartin) nspecin
     88    nspecin=nspecin/3
     89    if ((ldirect.eq.1).and.(nspec.ne.nspecin)) goto 997
     90
     91    do i=1,nspecin
     92      read(unitpartin)
     93      read(unitpartin)
     94      read(unitpartin) j,specin
     95      if ((ldirect.eq.1).and.(species(i)(1:7).ne.specin)) goto 996
     96    end do
     97   
     98    read(unitpartin) numpointin
     99    if (numpointin.ne.numpoint) then ! goto 995
     100
     101! eso 2016: moved this warning here to avoid out-of-block goto
     102!995   write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
     103      write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
     104      write(*,*) ' #### NUMBER OF RELEASE LOCATIONS DOES NOT     #### '
     105      write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
     106!  stop
     107      goto 999
     108    end if
     109
     110999 continue
     111    do i=1,numpointin
     112      read(unitpartin)
     113      read(unitpartin)
     114      read(unitpartin)
     115      read(unitpartin)
     116      do j=1,nspec
     117        read(unitpartin)
     118        read(unitpartin)
     119        read(unitpartin)
     120      end do
     121    end do
     122    read(unitpartin)
     123    read(unitpartin)
     124
     125    do ix=0,numxgrid-1
     126      read(unitpartin)
     127    end do
     128    close(unitpartin)
     129
     130
     131  ! Open data file of dumped particle data
     132  ! All processes read the whole file
     133  !***************************************
     134
     135  open(unitpartin,file=path(2)(1:length(2))//'partposit_end', &
     136         form='unformatted',err=998)
     137   
     138100 read(unitpartin,end=99) itimein
     139    i=0
     140200 i=i+1
     141    read(unitpartin) npoint_all(i),xlonin,ylatin,ztra1_all(i),itramem_all(i), &
     142         topo,pvi,qvi,rhoi,hmixi,tri,tti,(xmass1_all(i,j),j=1,nspec)
     143   
     144    if (xlonin.eq.-9999.9) goto 100
     145    xtra1_all(i)=(xlonin-xlon0)/dx
     146    ytra1_all(i)=(ylatin-ylat0)/dy
     147    numparticlecount_all=max(numparticlecount_all,npoint(i))
     148    goto 200
     149
     15099  numpart_all=i-1
     151   
     152    close(unitpartin)
     153
    85154  end do
    86155
    87   read(unitpartin) numpointin
    88   if (numpointin.ne.numpoint) goto 995
    89 999 continue
    90   do i=1,numpointin
    91     read(unitpartin)
    92     read(unitpartin)
    93     read(unitpartin)
    94     read(unitpartin)
    95     do j=1,nspec
    96       read(unitpartin)
    97       read(unitpartin)
    98       read(unitpartin)
    99     end do
     156
     157  ! Distribute particles among processes
     158  !**************************************
     159  lbnd=1
     160  ubnd=0
     161
     162  do ip=0, mp_partgroup_np-1
     163   
     164! Extra particle distributed in case remainder in the division
     165    if (ip.lt.mod(numpart_all,mp_partgroup_np)) then
     166      addone=1
     167    else
     168      addone=0
     169    end if
     170
     171    if (ip==0) then
     172      ubnd=ubnd + numpart_all/mp_partgroup_np + addone
     173    else
     174      ubnd=lbnd + numpart_all/mp_partgroup_np + addone - 1
     175    end if
     176   
     177    if (ip==mp_pid) then
     178     
     179      numparticlecount=numparticlecount_all/mp_partgroup_np+addone
     180      numpart=numpart_all/mp_partgroup_np+addone
     181
     182      xtra1(1:numpart) = xtra1_all(lbnd:ubnd)
     183      ytra1(1:numpart) = ytra1_all(lbnd:ubnd)
     184      ztra1(1:numpart) = ztra1_all(lbnd:ubnd)
     185      xmass1(1:numpart,:) = xmass1_all(lbnd:ubnd,:)
     186      npoint(1:numpart) = npoint_all(lbnd:ubnd)
     187      itramem(1:numpart) = itramem_all(lbnd:ubnd)
     188
     189    end if
     190
     191    lbnd=ubnd+1
     192
    100193  end do
    101   read(unitpartin)
    102   read(unitpartin)
    103 
    104   do ix=0,numxgrid-1
    105     read(unitpartin)
    106   end do
    107 
    108 
    109   ! Open data file of dumped particle data
    110   !***************************************
    111 
    112   close(unitpartin)
    113   open(unitpartin,file=path(2)(1:length(2))//'partposit_end', &
    114        form='unformatted',err=998)
    115 
    116 
    117 100   read(unitpartin,end=99) itimein
    118     i=0
    119 200   i=i+1
    120     read(unitpartin) npoint(i),xlonin,ylatin,ztra1(i),itramem(i), &
    121          topo,pvi,qvi,rhoi,hmixi,tri,tti,(xmass1(i,j),j=1,nspec)
    122 
    123     if (xlonin.eq.-9999.9) goto 100
    124     xtra1(i)=(xlonin-xlon0)/dx
    125     ytra1(i)=(ylatin-ylat0)/dy
    126     numparticlecount=max(numparticlecount,npoint(i))
    127     goto 200
    128 
    129 99   numpart=i-1
    130 
    131   close(unitpartin)
    132194
    133195  julin=juldate(ibdatein,ibtimein)+real(itimein,kind=dp)/86400._dp
     
    156218  stop
    157219
    158 !995   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
    159 995   write(*,*) ' #### FLEXPART MODEL WARNING IN READPARTPOSITIONS#### '
    160   write(*,*) ' #### NUMBER OF RELEASE LOCATIONS DOES NOT     #### '
    161   write(*,*) ' #### AGREE WITH CURRENT SETTINGS!             #### '
    162 !  stop
    163   goto 999
    164 
    165220996   write(*,*) ' #### FLEXPART MODEL ERROR IN READPARTPOSITIONS#### '
    166221  write(*,*) ' #### SPECIES NAMES TO BE READ IN DO NOT       #### '
  • src/readspecies.f90

    rec7fc72 r32b49c3  
    151151!  write(*,*) density(pos_spec)
    152152    read(unitspecies,'(e18.1)',end=22) dquer(pos_spec)
    153     write(*,*) 'dquer(pos_spec):', dquer(pos_spec)
     153!    write(*,*) 'dquer(pos_spec):', dquer(pos_spec)
    154154    read(unitspecies,'(e18.1)',end=22) dsigma(pos_spec)
    155155!  write(*,*) dsigma(pos_spec)
  • src/timemanager_mpi.f90

    re31b3b5 r32b49c3  
    104104
    105105  logical :: reqv_state=.false. ! .true. if waiting for a MPI_Irecv to complete
    106   integer :: j,ks,kp,l,n,itime=0,nstop,nstop1,memstat=0,mind
     106  integer :: j,ks,kp,l,n,itime=0,nstop,nstop1,memstat=0 !,mind
    107107! integer :: ksp
    108108  integer :: ip
     
    115115
    116116  real(sp) :: gridtotalunc
    117   real(dep_prec) :: drygridtotalunc,wetgridtotalunc,drydeposit(maxspec)
     117  real(dep_prec) :: drygridtotalunc=0_dep_prec,wetgridtotalunc=0_dep_prec,&
     118       & drydeposit(maxspec)=0_dep_prec
    118119  real :: xold,yold,zold,xmassfract
    119120  real, parameter :: e_inv = 1.0/exp(1.0)
     
    137138
    138139!  itime=0
    139   if (lroot) then
    140   !  write(*,45) itime,numpart*mp_partgroup_np,gridtotalunc,wetgridtotalunc,drygridtotalunc
    141     write(*,46) float(itime)/3600,itime,numpart*mp_partgroup_np
     140  if (lroot.or.mp_dev_mode) then
     141    write(*,45) itime,numpart*mp_partgroup_np,gridtotalunc,wetgridtotalunc,drygridtotalunc
     142  !  write(*,46) float(itime)/3600,itime,numpart*mp_partgroup_np
    142143   
    143144    if (verbosity.gt.0) then
     
    275276! For validation and tests: call the function below to set all fields to simple
    276277! homogeneous values
    277     if (mp_dev_mode) call set_fields_synthetic
     278!    if (mp_dev_mode) call set_fields_synthetic
    278279
    279280!*******************************************************************************
     
    554555       
    555556        !CGZ-lifetime: output species lifetime
    556         if (lroot) then
     557        if (lroot.or.mp_dev_mode) then
    557558        !   write(*,*) 'Overview species lifetime in days', &
    558559        !        real((species_lifetime(:,1)/species_lifetime(:,2))/real(3600.0*24.0))
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG