Skip to content

Commit

Permalink
Include precipitation in the YOG conservation check.
Browse files Browse the repository at this point in the history
  • Loading branch information
jatkinson1000 committed Jan 9, 2025
1 parent 292fbfd commit 34c0ab3
Showing 1 changed file with 11 additions and 4 deletions.
15 changes: 11 additions & 4 deletions src/physics/cam/nn_interface_cam.F90
Original file line number Diff line number Diff line change
Expand Up @@ -74,14 +74,15 @@ subroutine nn_convection_flux_CAM_init(nn_filename, sounding_filename)

end subroutine nn_convection_flux_CAM_init

subroutine yog_conservation_check(p_int, t, qv, qc, qi, &
subroutine yog_conservation_check(p_int, t, qv, qc, qi, prec, &
ncol, nver)
! Use this function to check the energy and water conversation by outputting the
! respective sums.

real(8), intent(in) :: p_int(:,:) ! pressure on grid
real(8), intent(in) :: t(:,:) ! temperature
real(8), intent(in) :: qv(:,:), qc(:,:), qi(:, :) ! moisture components: water vapor, cloud water, cloud ice
real(8), intent(in) :: prec(:) !! precipitation that has come out of the column
real(8) :: pdel(ncol, nver-1) ! pressure in grid cell
real(8) :: se(ncol) ! sum of energy
real(8) :: sw(ncol) ! sum of water
Expand Down Expand Up @@ -112,6 +113,7 @@ subroutine yog_conservation_check(p_int, t, qv, qc, qi, &

write(iulog, *), "Energy sum: ", se
write(iulog, *), "Water sum: ", sw
write(iulog, *), "Water + prec: ", sw + prec

end subroutine yog_conservation_check

Expand Down Expand Up @@ -171,6 +173,8 @@ subroutine nn_convection_flux_CAM(pres_cam, pres_int_cam, pres_sfc_cam, &
! Variables for the energy checker calculations at the end of the code
real(8), dimension(:,:), allocatable :: tabs_cam_out
!! absolute temperature [K] to the CAM model
real(8), dimension(:), allocatable :: prec_cam_out
!! precipitation to the CAM model
real(8), dimension(:,:), allocatable :: qv_cam_out, qc_cam_out, qi_cam_out
!! moisture content [kg/kg] to the CAM model
integer :: ncol_chnk, nver_chnk
Expand All @@ -181,7 +185,7 @@ subroutine nn_convection_flux_CAM(pres_cam, pres_int_cam, pres_sfc_cam, &
precsfc(:)=0.
!-----------------------------------------------------
write(iulog, *), "Before conversion:"
call yog_conservation_check(pres_int_cam, tabs_cam, qv_cam, qc_cam, qi_cam, ncol, pver)
call yog_conservation_check(pres_int_cam, tabs_cam, qv_cam, qc_cam, qi_cam, precsfc, ncol, pver)
!-----------------------------------------------------

! Interpolate CAM variables to the SAM pressure levels
Expand Down Expand Up @@ -260,7 +264,7 @@ subroutine nn_convection_flux_CAM(pres_cam, pres_int_cam, pres_sfc_cam, &
presi_col(i, :) = presi(:)*100.0
end do

call yog_conservation_check(presi_col, tabs_sam, qv_sam, qc_sam, qi_sam, ncol, nrf)
call yog_conservation_check(presi_col, tabs_sam, qv_sam, qc_sam, qi_sam, precsfc, ncol, nrf)
!-----------------------------------------------------

! Convert back into CAM variable tendencies (diff div by dtn) on SAM grid
Expand Down Expand Up @@ -292,19 +296,22 @@ subroutine nn_convection_flux_CAM(pres_cam, pres_int_cam, pres_sfc_cam, &
allocate(qv_cam_out(ncol_chnk,nver_chnk))
allocate(qc_cam_out(ncol_chnk,nver_chnk))
allocate(qi_cam_out(ncol_chnk,nver_chnk))
allocate(prec_cam_out(ncol_chnk))

tabs_cam_out = tabs_cam + dtn * ds / cp_cam
qv_cam_out = qv_cam + dtn * dqv
qc_cam_out = qc_cam + dtn * dqc
qi_cam_out = qi_cam + dtn * dqi
prec_cam_out = dtn * precsfc

write(iulog, *), "After conversion back:"
call yog_conservation_check(pres_int_cam, tabs_cam_out, qv_cam_out, qc_cam_out, qi_cam_out, ncol, pver)
call yog_conservation_check(pres_int_cam, tabs_cam_out, qv_cam_out, qc_cam_out, qi_cam_out, prec_cam_out, ncol, pver)

deallocate(tabs_cam_out)
deallocate(qv_cam_out)
deallocate(qc_cam_out)
deallocate(qi_cam_out)
deallocate(prec_cam_out)



Expand Down

0 comments on commit 34c0ab3

Please sign in to comment.