From e6dd89a370879740b8e5a85fbe7bbdd754236be6 Mon Sep 17 00:00:00 2001 From: Fredrik Jansson Date: Thu, 20 Aug 2020 14:51:16 +0200 Subject: [PATCH] fix #52, making modtimedepsv allocations consistent with modtimedep and removing redundant broadcast --- src/modtimedepsv.f90 | 18 +++++++++++++----- 1 file changed, 13 insertions(+), 5 deletions(-) diff --git a/src/modtimedepsv.f90 b/src/modtimedepsv.f90 index 243ed829..4d16b278 100644 --- a/src/modtimedepsv.f90 +++ b/src/modtimedepsv.f90 @@ -42,8 +42,8 @@ module modtimedepsv logical :: ltimedepsvz = .false. !< Switch for large scale forcings logical :: ltimedepsvsurf = .true. !< Switch for surface fluxes - integer, parameter :: kflux = 100 - integer, parameter :: kls = 100 + integer :: kflux + integer :: kls real, allocatable :: timesvsurf (:) real, allocatable :: svst (:,:) !< Time dependent surface scalar concentration @@ -56,8 +56,9 @@ module modtimedepsv !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine inittimedepsv use mpi - use modmpi, only :myid,my_real,mpi_logical,mpierr,comm3d - use modglobal,only :cexpnr,kmax,k1,ifinput,runtime,nsv + use modmpi, only :myid,my_real,mpi_logical,mpierr,comm3d + use modglobal, only :cexpnr,kmax,k1,ifinput,runtime,nsv + use modtestbed, only :ltestbed,ntnudge,tb_time implicit none character (80):: chmess @@ -69,6 +70,14 @@ subroutine inittimedepsv if (nsv==0 .or. .not.ltimedepsv ) return + if (ltestbed) then + kflux = ntnudge + kls = ntnudge + else + kflux = 10000 + kls = 10000 + end if + allocate(height(k1)) allocate(timesvsurf (0:kflux)) allocate(svst (kflux,nsv)) @@ -150,7 +159,6 @@ subroutine inittimedepsv call MPI_BCAST(timesvsurf(1:kflux),kflux,MY_REAL,0,comm3d,mpierr) - call MPI_BCAST(timesvz(1:kflux),kflux,MY_REAL,0,comm3d,mpierr) call MPI_BCAST(svst ,kflux*nsv,MY_REAL,0,comm3d,mpierr) call MPI_BCAST(timesvz(1:kls) ,kls,MY_REAL ,0,comm3d,mpierr) call MPI_BCAST(ltimedepsvsurf ,1,MPI_LOGICAL,0,comm3d,mpierr)