Skip to content

Commit

Permalink
chgres_cube: Add check of canopy water content for GEFS v12 GRIB2 data (
Browse files Browse the repository at this point in the history
ufs-community#482)

New routine 'check_cnwat' to correct erroneously high canopy moisture content 
values. Create unit test for this new routine.

Fixes ufs-community#481
  • Loading branch information
LarissaReames-NOAA authored May 17, 2021
1 parent 8995dc5 commit 42e39ae
Show file tree
Hide file tree
Showing 2 changed files with 80 additions and 10 deletions.
26 changes: 26 additions & 0 deletions sorc/chgres_cube.fd/input_data.F90
Original file line number Diff line number Diff line change
Expand Up @@ -56,6 +56,7 @@ module input_data
type(esmf_field), public :: ps_input_grid !< surface pressure
type(esmf_field), public :: terrain_input_grid !< terrain height
type(esmf_field), public :: temp_input_grid !< temperature

type(esmf_field), public :: u_input_grid !< u/v wind at grid
type(esmf_field), public :: v_input_grid !< box center
type(esmf_field), public :: wind_input_grid !< 3-component wind
Expand Down Expand Up @@ -133,6 +134,7 @@ module input_data
public :: read_input_nst_data
public :: cleanup_input_nst_data
public :: check_soilt
public :: check_cnwat
public :: quicksort
public :: convert_winds

Expand Down Expand Up @@ -5239,6 +5241,7 @@ subroutine read_input_sfc_grib2_file(localpet)
dummy2d(:,:) = 0.0_esmf_kind_r4
endif
endif
call check_cnwat(dummy2d)
dummy2d_8= real(dummy2d,esmf_kind_r8)
print*,'cnwat ',maxval(dummy2d),minval(dummy2d)
endif
Expand Down Expand Up @@ -6639,4 +6642,27 @@ subroutine check_soilt(soilt, landmask, skint)
enddo
enddo
end subroutine check_soilt

!> When using GEFS data, some points on the target grid have
!> unreasonable canpy moisture content, so zero out any
!> locations with unrealistic canopy moisture values (>0.5).
!!
!! @param cnwat [input] 2-dimensional canopy moisture content
!! @author Larissa Reames CIMMS/NSSL

subroutine check_cnwat(cnwat)
implicit none
real(esmf_kind_r4), intent(inout) :: cnwat(i_input,j_input)

real(esmf_kind_r4) :: max_cnwat = 0.5

integer :: i, j

do i = 1,i_input
do j = 1,j_input
if (cnwat(i,j) .gt. max_cnwat) cnwat(i,j) = 0.0_esmf_kind_r4
enddo
enddo
end subroutine check_cnwat

end module input_data
64 changes: 54 additions & 10 deletions tests/chgres_cube/ftst_sfc_input_data.F90
Original file line number Diff line number Diff line change
@@ -1,7 +1,10 @@
program test_sfc_input_data

! Test sfc_input_data. For now only one test is performed
! on the subroutine soilt_check using a sample 3x3 matrix.
! Test sfc_input_data. Tests include:
!
! subroutine soilt_check using a sample 3x3 matrix
!
! subroutine cnwat_check using a sample 3x3 matrix
!
! author: Larissa Reames (larissa.reames@noaa.gov)

Expand All @@ -19,9 +22,16 @@ program test_sfc_input_data
soilt_updated(:,:,:), &
soilt_correct(:,:,:)
real(esmf_kind_r8), allocatable :: skint(:,:)
real(esmf_kind_r4), allocatable :: cnwat_bad(:,:), &
cnwat_updated(:,:), &
cnwat_correct(:,:)

call mpi_init(ierr)

!--------------------------------------------------------!
!------------------- CHECK_SOILT TEST -------------------!
!--------------------------------------------------------!

lsoil_input = 2
i_input = 3
j_input = 3
Expand Down Expand Up @@ -80,20 +90,54 @@ program test_sfc_input_data
call check_soilt(soilt_updated,mask,skint)

if (any(soilt_updated /= soilt_correct)) then
print*,'TEST FAILED '
print*,'ARRAY SHOULD BE:', soilt_correct
print*,'ARRAY FROM TEST:', soilt_updated
print*,'SOILT TEST FAILED '
print*,'SOILT ARRAY SHOULD BE:', soilt_correct
print*,'SOILT ARRAY FROM TEST:', soilt_updated
stop 2
else
print*, 'SOILT TEST SUCCEEDED. CONTINUING'
endif
deallocate(mask,skint,soilt_bad,soilt_updated,soilt_correct)

deallocate(mask)
deallocate(skint)
deallocate(soilt_bad)
deallocate(soilt_updated)
deallocate(soilt_correct)
!--------------------------------------------------------!
!------------------- CHECK_CNWAT TEST -------------------!
!--------------------------------------------------------!

i_input = 3
j_input = 3

allocate(cnwat_bad(i_input,j_input))
allocate(cnwat_updated(i_input,j_input))
allocate(cnwat_correct(i_input,j_input))

! Canopy water content array with some values that are too high.
! This will be the input array for check_cnwat
cnwat_bad = reshape((/999.0, 0.3, 0.4, 0.1, 0.3, 999.0, 999.0, 0.3, 0.15/), &
(/3,3/))

! Corrected canopy water content array that should be passed
! back from check_cnwat.
cnwat_correct = reshape((/0.0, 0.3, 0.4, 0.1, 0.3, 0.0, 0.0, 0.3, 0.15/), &
(/3,3/))

print*,"Starting test of check_cnwat subroutine."

cnwat_updated = cnwat_bad
call check_cnwat(cnwat_updated)

if (any(cnwat_updated /= cnwat_correct)) then
print*,'CNWAT TEST FAILED '
print*,'CNWAT ARRAY SHOULD BE:', cnwat_correct
print*,'CNWAT ARRAY FROM TEST:', cnwat_updated
stop 2
else
print*, 'CNWAT TEST SUCCEEDED.'
endif
deallocate(cnwat_bad,cnwat_updated,cnwat_correct)

call mpi_finalize(ierr)

print*,"SUCCESS!"


end program test_sfc_input_data

0 comments on commit 42e39ae

Please sign in to comment.