Changeset 6b22af9 in flexpart.git for src


Ignore:
Timestamp:
Nov 20, 2015, 5:41:29 PM (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:
a1f4dd6
Parents:
df967a99
Message:

Local branch for espen, working with OpenMP version of readwind

Location:
src
Files:
8 edited

Legend:

Unmodified
Added
Removed
  • src/com_mod.f90

    r43225d1 r6b22af9  
    732732  integer :: count_clock, count_clock0,  count_rate, count_max
    733733  real    :: tins
    734   logical, parameter :: nmlout=.false.
     734  logical, parameter :: nmlout=.true.
    735735
    736736  ! These variables are used to avoid having separate versions of
  • src/concoutput.f90

    rb069789 r6b22af9  
    118118! This fixes a bug where the dates file kept growing across multiple runs
    119119
    120 ! If 'dates' file exists, make a backup
     120! If 'dates' file exists in output directory, make a backup
    121121  inquire(file=path(2)(1:length(2))//'dates', exist=ldates_file)
    122122  if (ldates_file.and.init) then
     
    129129      if (ierr.ne.0) exit
    130130      write(unit=unittmp, fmt='(a)', iostat=ierr, advance='yes') trim(dates_char)
    131 !      if (ierr.ne.0) write(*,*) "Write error, ", ierr
    132131    end do
    133132    close(unit=unitdates)
  • src/makefile

    radf46ae r6b22af9  
    3939
    4040ifeq ($(gcc), 4.9)
    41 # Compiled libraries under user user ~flexpart, gfortran v4.9
     41# Compiled libraries under users ~flexpart (or ~espen), gfortran v4.9
    4242        ROOT_DIR = /homevip/flexpart/
    4343#       ROOT_DIR = /homevip/espen/
     
    7070LIBS = -lgrib_api_f90 -lgrib_api -lm -ljasper -lnetcdff -llapack  # -fopenmp # -llapack
    7171
    72 FFLAGS   = -I$(INCPATH1) -I$(INCPATH2) -O$(O_LEV) -g -m64 -mcmodel=medium -fconvert=little-endian -frecord-marker=4 -fmessage-length=0 -flto=jobserver -O$(O_LEV) -mtune=native $(FUSER) # -march=native
     72FFLAGS   = -I$(INCPATH1) -I$(INCPATH2) -O$(O_LEV) -g -m64 -mcmodel=medium -fconvert=little-endian -frecord-marker=4 -fmessage-length=0 -flto=jobserver -O$(O_LEV) $(FUSER) # -march=native
    7373
    7474DBGFLAGS = -I$(INCPATH1) -I$(INCPATH2) -O$(O_LEV_DBG) -g3 -ggdb3 -m64 -mcmodel=medium -fconvert=little-endian -frecord-marker=4 -fmessage-length=0 -flto=jobserver -O$(O_LEV_DBG) -fbacktrace -Warray-bounds  -Wall -fcheck=all $(FUSER)  #  -ffpe-trap=invalid,overflow,denormal,underflow,zero -fdump-core
     
    420420timemanager_mpi.o: com_mod.o flux_mod.o mpi_mod.o oh_mod.o outg_mod.o \
    421421        par_mod.o point_mod.o unc_mod.o xmass_mod.o
    422 verttransform.o: cmapf_mod.o com_mod.o par_mod.o
     422verttransform.o: cmapf_mod.o com_mod.o par_mod.o 
    423423verttransform_gfs.o: cmapf_mod.o com_mod.o par_mod.o
    424424verttransform_nests.o: com_mod.o par_mod.o
  • src/mpi_mod.f90

    rf55fdce r6b22af9  
    105105! If setting this to .false., numwfmem must be set to 3
    106106!===============================================================================
    107   logical :: lmp_sync=.true.
     107  logical :: lmp_sync=.false.
    108108!===============================================================================
    109109
     
    119119  logical, parameter :: mp_dbg_out = .false.
    120120  logical, parameter :: mp_time_barrier=.true.
    121   logical, parameter :: mp_measure_time=.false.
     121  logical, parameter :: mp_measure_time=.true.
    122122  logical, parameter :: mp_exact_numpart=.true.
    123123
     
    214214      lmp_sync=.true. ! :DBG: eso fix this...
    215215    end if
    216 ! TODO: Add warnings for unimplemented flexpart features
     216
     217! TODO: Add more warnings for unimplemented flexpart features
    217218
    218219! Set ID of process that calls getfield/readwind.
     
    434435! Distribute particle variables from pid0 to all processes.
    435436! Called from timemanager
    436 !
     437! *NOT IN USE* at the moment, but can be useful for debugging
    437438!
    438439!***********************************************************************
     
    588589! concentrations/writing output grids.                                         *
    589590!                                                                              *
    590 ! Currently not in use, as each process calculates concentrations              *
    591 ! separately.                                                                  *
     591! Currently *not in use*, as each process calculates concentrations            *
     592! separately, but can be useful for debugging                                  *
    592593!                                                                              *
    593594! To use this routine (together with mpif_tm_send_vars) to transfer data       *
     
    595596! (in timemanager_mpi.f90)                                                     *
    596597!                                                                              *
    597 !      if (lroot) tot_numpart=numpart ! root stores total numpart            *
     598!      if (lroot) tot_numpart=numpart ! root stores total numpart              *
    598599!      call mpif_tm_send_dims                                                  *
    599600!      if (numpart>1) then                                                     *
     
    605606!                                                                              *
    606607!      if (numpart>0) then                                                     *
    607 !          if (lroot) numpart = tot_numpart                                  *
     608!          if (lroot) numpart = tot_numpart                                    *
    608609!       call mpif_tm_collect_vars                                              *
    609610!      end if                                                                  *
    610611!                                                                              *
    611 !    Then a section that begins with "if (lroot) ..." will behave like       *
    612 !    serial flexpart                                                           *
     612! Then a section that begins with "if (lroot) ..." will behave like            *
     613! serial flexpart                                                              *
    613614!                                                                              *
    614615!*******************************************************************************
     
    20662067
    20672068
     2069  subroutine set_fields_synthetic()
     2070!*******************************************************************************
     2071! DESCRIPTION
     2072!   Set all meteorological fields to synthetic (usually constant/homogeneous)
     2073!   values.
     2074!   Used for validation and error-checking
     2075!
     2076! NOTE
     2077!   This version uses asynchronious communications.
     2078!
     2079! VARIABLES
     2080!   
     2081!
     2082! TODO
     2083!
     2084!*******************************************************************************
     2085    use com_mod
     2086
     2087    implicit none
     2088
     2089    integer, parameter :: li=1, ui=3 ! wfmem indices (i.e, operate on entire field)
     2090
     2091
     2092! Variables transferred at initialization only
     2093!*********************************************
     2094!    readclouds=readclouds_
     2095    oro(:,:)=0.0
     2096    excessoro(:,:)=0.0
     2097    lsm(:,:)=0.0
     2098    xlanduse(:,:,:)=0.0
     2099!    wftime
     2100!    numbwf
     2101!    nmixz
     2102!    height
     2103
     2104! Time-varying fields:
     2105    uu(:,:,:,li:ui) = 10.0
     2106    vv(:,:,:,li:ui) = 0.0
     2107    uupol(:,:,:,li:ui) = 10.0 ! TODO check if ok
     2108    vvpol(:,:,:,li:ui)=0.0
     2109    ww(:,:,:,li:ui)=0.
     2110    tt(:,:,:,li:ui)=300.
     2111    rho(:,:,:,li:ui)=1.3
     2112    drhodz(:,:,:,li:ui)=.0
     2113    tth(:,:,:,li:ui)=0.0
     2114    qvh(:,:,:,li:ui)=1.0
     2115    qv(:,:,:,li:ui)=1.0
     2116
     2117    pv(:,:,:,li:ui)=1.0
     2118    clouds(:,:,:,li:ui)=0.0
     2119
     2120    clwc(:,:,:,li:ui)=0.0
     2121    ciwc(:,:,:,li:ui)=0.0
     2122 
     2123! 2D fields
     2124
     2125    cloudsh(:,:,li:ui)=0.0
     2126    vdep(:,:,:,li:ui)=0.0
     2127    ps(:,:,:,li:ui)=1.0e5
     2128    sd(:,:,:,li:ui)=0.0
     2129    tcc(:,:,:,li:ui)=0.0
     2130    tt2(:,:,:,li:ui)=300.
     2131    td2(:,:,:,li:ui)=250.
     2132    lsprec(:,:,:,li:ui)=0.0
     2133    convprec(:,:,:,li:ui)=0.0
     2134    ustar(:,:,:,li:ui)=1.0
     2135    wstar(:,:,:,li:ui)=1.0
     2136    hmix(:,:,:,li:ui)=10000.
     2137    tropopause(:,:,:,li:ui)=10000.
     2138    oli(:,:,:,li:ui)=0.01
     2139    z0=1.0
     2140 
     2141  end subroutine set_fields_synthetic
     2142
    20682143end module mpi_mod
  • src/outgrid_init.f90

    r8a65cb0 r6b22af9  
    205205    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
    206206  if (ldirect.gt.0) then
    207   allocate(wetgridunc(0:numxgrid-1,0:numygrid-1,maxspec, &
     207    allocate(wetgridunc(0:numxgrid-1,0:numygrid-1,maxspec, &
    208208       maxpointspec_act,nclassunc,maxageclass),stat=stat)
    209     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
     209    if (stat.ne.0) write(*,*)'ERROR: could not allocate wetgridunc'
    210210  allocate(drygridunc(0:numxgrid-1,0:numygrid-1,maxspec, &
    211211       maxpointspec_act,nclassunc,maxageclass),stat=stat)
    212     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
     212    if (stat.ne.0) write(*,*)'ERROR: could not allocate drygridunc'
    213213  endif
    214214
  • src/par_mod.f90

    radf46ae r6b22af9  
    210210  !**************************************************
    211211
    212   integer,parameter :: maxpart=60000000
     212  integer,parameter :: maxpart=40000000
    213213!  integer,parameter :: maxpart=60000000
    214214!  integer,parameter :: maxpart=120000000
    215   integer,parameter :: maxspec=6
     215  integer,parameter :: maxspec=1
    216216
    217217  ! maxpart                 Maximum number of particles
     
    226226  ! ---------
    227227  integer,parameter :: maxwf=50000, maxtable=1000, numclass=13, ni=11
    228   integer,parameter :: numwfmem=2 ! Serial version/MPI with 2 fields
    229   !integer,parameter :: numwfmem=3 ! MPI with 3 fields
     228  !integer,parameter :: numwfmem=2 ! Serial version/MPI with 2 fields
     229  integer,parameter :: numwfmem=3 ! MPI with 3 fields
    230230
    231231  ! maxwf                   maximum number of wind fields to be used for simulation
     
    254254  !*********************************
    255255
    256   integer,parameter :: maxrand=120000000
    257 !  integer,parameter :: maxrand=2000000
    258 !  integer,parameter :: maxrand=20
    259 
     256!  integer,parameter :: maxrand=120000000
     257  integer,parameter :: maxrand=200000000
     258!
    260259  ! maxrand                 number of random numbers used
    261260 
  • src/pathnames

    r5f9d14a r6b22af9  
    11../options/
    2 ../output/
     2../out_clwc_leo/
    33/
    4 /xnilu_wrk/flex_wrk/WIND_FIELDS/AVAILABLE_ECMWF_OPER_fields_global
     4/xnilu_wrk/flex_wrk/zhg/ECMWF_CLWC/Availables
    55============================================
    66
     7../options_nikos_backwards/
     8/xnilu_wrk/flex_wrk/WIND_FIELDS/AVAILABLE_ECMWF_OPER_fields_global
     9/xnilu_wrk/flex_wrk/WIND_FIELDS/NCEP/2010_data_05/AVAILABLE
     10
     11
     12 
  • src/timemanager_mpi.f90

    rf55fdce r6b22af9  
    218218      call mpif_gf_send_vars(memstat)
    219219      call mpif_gf_send_vars_nest(memstat)
    220 ! Version 2  (lmp_sync=.false.) is also used whenever 2 new fields are read (in which
    221 ! case async send/recv is impossible.
     220! Version 2  (lmp_sync=.false., see below) is also used whenever 2 new fields are
     221! read (as at first time step), in which case async send/recv is impossible.
    222222    else if (.not.lmp_sync.and.lmp_use_reader.and.memstat.ge.32) then
    223223      call mpif_gf_send_vars(memstat)
     
    233233! READER PROCESS:
    234234      if (memstat.gt.0..and.memstat.lt.32.and.lmp_use_reader.and.lmpreader) then
     235        if (mp_dev_mode) write(*,*) 'Reader process: calling mpif_gf_send_vars_async'
    235236        call mpif_gf_send_vars_async(memstat)
    236237      end if
     
    249250      ! at next time step. Issue receive request anyway, cancel at mpif_gf_request
    250251      if (memstat.gt.0.and.lmp_use_reader.and..not.lmpreader) then
     252        if (mp_dev_mode) write(*,*) 'Receiving process: calling mpif_gf_send_vars_async. PID: ', mp_pid
    251253        call mpif_gf_recv_vars_async(memstat)
    252254      end if
     
    256258    if (mp_measure_time.and..not.(lmpreader.and.lmp_use_reader)) call mpif_mtime('getfields',1)
    257259
     260! For validation and tests: call the function below to set all fields to simple
     261! homogeneous values
     262    if (mp_dev_mode) call set_fields_synthetic
    258263
    259264!*******************************************************************************
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG