Skip to content
Merged
Show file tree
Hide file tree
Changes from 8 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
30 changes: 21 additions & 9 deletions src/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -13,18 +13,20 @@ project(clm_tests Fortran C)

include(CIME_utils)


# find needed external packages
# NetCDF is required -- because PIO and NetCDF are required by the standard default ESMF libraries
find_package(NetCDF 4.7.4 REQUIRED Fortran)

# The following - for finding ESMF - is copied from the share CMakeLists.txt
if (DEFINED ENV{ESMF_ROOT})
if(DEFINED ENV{ESMF_ROOT})
list(APPEND CMAKE_MODULE_PATH $ENV{ESMF_ROOT}/cmake)
endif()

find_package(ESMF REQUIRED)

# This adds include directories needed for ESMF
set(CMAKE_Fortran_FLAGS "${CMAKE_Fortran_FLAGS} ${ESMF_F90COMPILEPATHS} ")

# This (which is *not* done in the share CMakeLists.txt) adds all directories and
# libraries needed when linking ESMF, including any dependencies of ESMF. (But note that
# this does *not* include the "-lesmf" itself). In particular, note that this includes any
Expand All @@ -43,12 +45,12 @@ add_subdirectory(${SRCROOT}/share/src csm_share)
add_subdirectory(${SRCROOT}/share/unit_test_stubs/util csm_share_stubs)

# Add files needed from CMEPS
list ( APPEND drv_sources_needed
list(APPEND drv_sources_needed
${SRCROOT}/components/cmeps/cesm/nuopc_cap_share/glc_elevclass_mod.F90
${SRCROOT}/components/cmeps/cesm/nuopc_cap_share/shr_dust_emis_mod.F90
${SRCROOT}/components/cmeps/cesm/nuopc_cap_share/shr_expr_parser_mod.F90
${SRCROOT}/components/cmeps/cesm/nuopc_cap_share/shr_fire_emis_mod.F90
)
)

# Add CLM source directories
add_subdirectory(${CLM_ROOT}/src/utils clm_utils)
Expand All @@ -60,6 +62,14 @@ add_subdirectory(${CLM_ROOT}/src/main clm_main)
add_subdirectory(${CLM_ROOT}/src/init_interp clm_init_interp)
add_subdirectory(${CLM_ROOT}/src/self_tests clm_self_tests)

# Add FATES source directories
add_subdirectory(${CLM_ROOT}/src/fates/main fates_main)
add_subdirectory(${CLM_ROOT}/src/fates/biogeochem fates_biogeochem)
add_subdirectory(${CLM_ROOT}/src/fates/biogeophys fates_biogeophys)
add_subdirectory(${CLM_ROOT}/src/fates/parteh fates_parteh)
add_subdirectory(${CLM_ROOT}/src/fates/fire fates_fire)
add_subdirectory(${CLM_ROOT}/src/fates/radiation fates_radiation)

# Add general unit test directories (stubbed out files, etc.)
add_subdirectory(unit_test_stubs)
add_subdirectory(unit_test_shr)
Expand All @@ -69,10 +79,11 @@ add_subdirectory(unit_test_shr)
# TODO: this should be moved into a general-purpose function in Sourcelist_utils.
# Then each removal could be replaced with a single call, like:
# remove_source_file(${share_sources} "shr_mpi_mod.F90")
foreach (sourcefile ${share_sources})
foreach(sourcefile ${share_sources})
# Remove shr_mpi_mod from share_sources.
# This is needed because we want to use the mock shr_mpi_mod in place of the real one
string(REGEX MATCH "shr_mpi_mod.F90" match_found ${sourcefile})

if(match_found)
list(REMOVE_ITEM share_sources ${sourcefile})
endif()
Expand All @@ -83,6 +94,7 @@ foreach (sourcefile ${share_sources})
# error message, "Cannot open module file 'pio.mod'") on a Mac without a pre-built PIO
# (where ESMF was built with its internal PIO).
string(REGEX MATCH "shr_pio_mod.F90" match_found ${sourcefile})

if(match_found)
list(REMOVE_ITEM share_sources ${sourcefile})
endif()
Expand All @@ -93,19 +105,20 @@ endforeach()
add_library(csm_share ${share_sources} ${drv_sources_needed})
declare_generated_dependencies(csm_share "${share_genf90_sources}")
add_library(clm ${clm_sources})
add_library(fates ${fates_sources})
declare_generated_dependencies(clm "${clm_genf90_sources}")
add_dependencies(clm csm_share ESMF)
add_dependencies(clm csm_share ESMF fates)

# We need to look for header files here, in order to pick up shr_assert.h
include_directories(${SRCROOT}/share/include)

# Tell cmake to look for libraries & mod files here, because this is where we built libraries
include_directories(${CMAKE_CURRENT_BINARY_DIR})
include_directories (${NETCDF}/include)
include_directories(${NETCDF}/include)

# Directories and libraries to include in the link step
link_directories(${CMAKE_CURRENT_BINARY_DIR})
link_libraries( netcdf esmf )
link_libraries(netcdf esmf)

# Add the test directories
# Note: it's possible that these could be added by each source directory that
Expand All @@ -125,4 +138,3 @@ add_subdirectory(${CLM_ROOT}/src/self_tests/test clm_self_tests_test)
# Add driver unit test directories
# (these should be moved to the appropriate submodule)
add_subdirectory(${CLM_ROOT}/src/drv_test drv_test)

1 change: 1 addition & 0 deletions src/main/test/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -8,4 +8,5 @@ add_subdirectory(filter_test)
add_subdirectory(initVertical_test)
add_subdirectory(ncdio_utils_test)
add_subdirectory(topo_test)
add_subdirectory(decomp_test)
add_subdirectory(abortutils_test)
8 changes: 8 additions & 0 deletions src/main/test/decomp_test/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -0,0 +1,8 @@
set(pfunit_sources
test_decompMod.pf)

add_pfunit_ctest(decomp
TEST_SOURCES "${pfunit_sources}"
LINK_LIBRARIES clm csm_share esmf
EXTRA_FINALIZE unittest_finalize_esmf
EXTRA_USE unittestInitializeAndFinalize)
106 changes: 106 additions & 0 deletions src/main/test/decomp_test/test_decompMod.pf
Original file line number Diff line number Diff line change
@@ -0,0 +1,106 @@
module test_decompMod

! Tests of decompMod

use funit
use decompMod
use shr_kind_mod , only : r8 => shr_kind_r8

implicit none

@TestCase
type, extends(TestCase) :: TestDecompMod
contains
procedure :: setUp
procedure :: tearDown
procedure :: create_simpleSingleDecomp
end type TestDecompMod

integer, parameter :: ni = 2
integer, parameter :: nj = 2

contains

! ========================================================================
! Helper routines
! ========================================================================

subroutine setUp(this)
class(TestDecompMod), intent(inout) :: this

call this%create_simpleSingleDecomp()
end subroutine setUp

subroutine tearDown(this)
class(TestDecompMod), intent(inout) :: this

call decompmod_clean()

end subroutine tearDown

subroutine create_simpleSingleDecomp(this)
use spmdMod, only : iam
class(TestDecompMod), intent(inout) :: this

integer :: clump_pproc
! TOTO: When decompMod has it's own allocate method that could be used here
nclumps = 1
clump_pproc = nclumps
allocate(procinfo%cid(clump_pproc))
allocate(clumps(nclumps))
! Set the procinfo and clumps values
! TOD: Use initialization method when available (currently in decompInitMod)
procinfo%cid = 1
procinfo%ncells = ni*nj
procinfo%begg = 1
procinfo%endg = procinfo%ncells
procinfo%nclumps = nclumps
clumps(:)%owner = iam
clumps(:)%begg = 1
clumps(:)%endg = procinfo%ncells

end subroutine create_simpleSingleDecomp
! ========================================================================
! Begin tests
! ========================================================================

@Test
subroutine test_get_clump_bounds(this)
class(TestDecompMod), intent(inout) :: this

type(bounds_type) :: bounds
integer :: n

do n = 1, procinfo%nclumps
call get_clump_bounds(n, bounds)
@assertEqual(bounds%level, bounds_level_clump)
@assertEqual(bounds%clump_index, n)
end do
end subroutine test_get_clump_bounds

@Test
subroutine test_get_proc_bounds(this)
class(TestDecompMod), intent(inout) :: this

type(bounds_type) :: bounds

! Add optional argument, just to test that it can handle it
call get_proc_bounds(bounds, allow_call_from_threaded_region=.true.)
@assertEqual(bounds%level, bounds_level_proc)
@assertEqual(bounds%clump_index, -1)
end subroutine test_get_proc_bounds

@Test
subroutine test_proc_clump_bounds_equal(this)
class(TestDecompMod), intent(inout) :: this

type(bounds_type) :: bounds_clump, bounds_proc

@assertTrue(procinfo%nclumps == 1)
call get_clump_bounds(1, bounds_clump)
call get_proc_bounds(bounds_proc)
@assertEqual(bounds_proc%begg, bounds_clump%begg)
@assertEqual(bounds_proc%endg, bounds_clump%endg)
end subroutine test_proc_clump_bounds_equal

end module test_decompMod
1 change: 1 addition & 0 deletions src/utils/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -22,6 +22,7 @@ list(APPEND clm_sources
SparseMatrixMultiplyMod.F90
IssueFixedMetadataHandler.F90
NumericsMod.F90
spmdMod.F90
)

sourcelist_to_parent(clm_sources)