Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Addition of an example unit test for users to learn from #581

Merged
Show file tree
Hide file tree
Changes from all 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
7 changes: 4 additions & 3 deletions sorc/chgres_cube.fd/grib2_util.F90
Original file line number Diff line number Diff line change
Expand Up @@ -17,18 +17,19 @@ module grib2_util

implicit none

public :: rh2spfh
public :: convert_omega
public :: rh2spfh_gfs
public :: fpvsnew


contains

!> Convert relative humidity to specific humidity.
!> Calculation of saturation water vapor pressure is based on
!> Brock and Richardson 2001 (Meterological Measurement
!> Systems, p. 86, equation 5.1)
!!
!! @param[inout] rh_sphum rel humidity on input. spec hum on output.
!! @param[inout] rh_sphum rel humidity (%) on input. spec hum (kg/kg) on output.
!! @param[in] p pressure in Pa
!! @param[in] t temperature
!! @author Larissa Reames
Expand Down Expand Up @@ -59,7 +60,7 @@ subroutine rh2spfh(rh_sphum,p,t)
!print *, 'q = ', sphum

!if (P .eq. 100000.0) THEN
! print *, 'T = ', T, ' RH = ', RH, ' P = ', P, ' es = ', es, ' e = ', e, ' q = ', sphum
!print *, 'T = ', t, ' RH = ', rh, ' P = ', p, ' es = ', es, ' e = ', e, ' q = ', rh_sphum
!end if

end subroutine RH2SPFH
Expand Down
13 changes: 12 additions & 1 deletion tests/README.md
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
## HOW TO CREATE A UNIT TEST.

Unit tests should test only small parts of a program. For example,
a single routine or function.
a single routine or function, or multiple closely-linked routines.

Source code for a test shall be placed in a program
specific directory under ./tests.
Expand Down Expand Up @@ -86,6 +86,17 @@ To run the tests locally, invoke one of the following after the `make install`:
For the parallel tests to run locally, update the machine-specific
"mpi_exec" script under [cmake](../cmake) with your run account and queue.

### EXAMPLE UNIT TEST:

An example unit test, ftst_example.F90, exists in the tests/chgres_cube directory
and is currently compiled and run with existing tests

This simple test checks whether the routine rh2spfh in chgres_cube.fd/grib2_util.F90
is working correctly, and contains detailed comments explaining what each section
of the test does. It also prompts the user to create additional lines of code to
test one more subroutine from grib2_util.F90. Use this to get started understanding
the unit test framework.

### QUESTIONS

Please contact the repository managers: https://github.com/NOAA-EMC/UFS_UTILS/wiki
4 changes: 4 additions & 0 deletions tests/chgres_cube/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -199,3 +199,7 @@ add_mpi_test(chgres_cube-ftst_surface_search_many
add_executable(ftst_read_vcoord ftst_read_vcoord.F90)
target_link_libraries(ftst_read_vcoord chgres_cube_lib)
add_test(NAME chgres_cube-ftst_read_vcoord COMMAND ftst_read_vcoord)

add_executable(ftst_example ftst_example.F90)
target_link_libraries(ftst_example chgres_cube_lib)
add_test(NAME chgres_cube-ftst_example COMMAND ftst_example)
99 changes: 99 additions & 0 deletions tests/chgres_cube/ftst_example.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,99 @@
! Unit test for rh2spfh in grib2_utils to be used as an example for users
! learning how to write unit tests. Users are prompted to add an additional
! test for convert_omega
! Larissa Reames OU/CIMMS/NOAA/NSSL/FRDD

program ftst_example

use esmf
use model_grid, only : i_input, j_input
use grib2_util, only : rh2spfh

implicit none

real(esmf_kind_r4), allocatable :: rh_spfh(:,:), &
rh_orig(:,:), &
spfh_returned(:,:), &
spfh_correct(:,:)
real(esmf_kind_r8),allocatable :: t(:,:)
real(esmf_kind_r8) :: p
real,parameter :: EPS = 1.0E-6

i_input = 2
j_input = 2
allocate(rh_spfh(i_input,j_input))
allocate(rh_orig(i_input,j_input))
allocate(spfh_returned(i_input,j_input))
allocate(spfh_correct(i_input,j_input))
allocate(t(i_input,j_input))

! -------------------------------------------------------------------------
! Set constants/values to be passed to the unit test. In this case it's a
! set of single values, but it could be more complicated, like an
! n-dimensional array or ESMF objects.
! -------------------------------------------------------------------------

rh_spfh(:,:) = 50.0 ! Relative humidity (%)
p = 100000.0 ! Pressure (Pa)
t(:,:) = 300.0 ! Temperature (K)
spfh_correct(:,:) = 10.978297E-3 ! Correct specific humidity value (kg/kg)

print*, "Starting Unit Testing rh2spfh."

!-------------------------------------------------------------------------
! Execute testing below by running target function rh2spfh and providing
! values set above
!-------------------------------------------------------------------------

rh_orig = rh_spfh ! Save the original RH value for posterity
call rh2spfh(rh_spfh,p,t)
spfh_returned = rh_spfh ! Rename the returned value for clarity

!-------------------------------------------------------------------------
! Check the returned value against what we know to be the correct answer.
! If the correct result is returned (within a certain small tolerance),
! then the test passes and the called routine is working as expected. If the
! incorrect value is passed back, the test fails and an error is returned.
!-------------------------------------------------------------------------

if ( any(abs(spfh_returned - spfh_correct) .gt. EPS)) then
stop 1
print*, "RH2SPFH TEST FAILED."
endif
edwardhartnett marked this conversation as resolved.
Show resolved Hide resolved

!-------------------------------------------------------------------------
! If you are trying to debug a test failure, code like the commented
! section below might prove useful.
!-------------------------------------------------------------------------

! if ( any(abs(spfh_returned - spfh_correct) .lt. EPS)) then
! print*, "RH2SPFH TEST PASSED. SUCCESS!"
! else
! print*, "RH2SPFH TEST FAILED."
! print*, "TEST RETURNED VALUE OF ", spfh_returned
! print*, "RETURNED VALUE EXPECT TO BE ", spfh_correct
! stop 1
! endif

!-------------------------------------------------------------------------
! Make sure to deallocate any and all allocatable arrays that you use
!-------------------------------------------------------------------------

deallocate(rh_spfh,spfh_correct,rh_orig,spfh_returned,t)

! -------------------------------------------------------------------------
GeorgeGayno-NOAA marked this conversation as resolved.
Show resolved Hide resolved
! You can test multiple subroutines (units) in any test file. This would
! be a good place to test the other subroutine in grib2_util,
! convert_omega. Make note of the difference in variable size for this
! routine. You don't have to pass an array of size you'd normally
! encounter for these variables (like 200x200x64), just choose a small
! size with the proper dimensionality, say 3x3x2, and fill it with values
! typical of the various arrays. You can check the returned array element-
! by-element, or use the any() command to check all elements at once. Make
! sure to provide a helpful failure message indicating where the failure
! occured and what the returned/expected values were at that location. Also,
! don't forget to deallocate any allocatable arrays as this will cause a
! failure when compiling the test with gnu compilers and address sanitizers.
! -------------------------------------------------------------------------

end program ftst_example