!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !! !! !! GNU General Public License !! !! !! !! This file is part of the Flexible Modeling System (FMS). !! !! !! !! FMS is free software; you can redistribute it and/or modify !! !! it and are expected to follow the terms of the GNU General Public !! !! License as published by the Free Software Foundation. !! !! !! !! FMS is distributed in the hope that it will be useful, !! !! but WITHOUT ANY WARRANTY; without even the implied warranty of !! !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the !! !! GNU General Public License for more details. !! !! !! !! You should have received a copy of the GNU General Public License !! !! along with FMS; if not, write to: !! !! Free Software Foundation, Inc. !! !! 59 Temple Place, Suite 330 !! !! Boston, MA 02111-1307 USA !! !! or see: !! !! http://www.gnu.org/licenses/gpl.txt !! !! !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! ! ! Diagnostics for use in the Energy Balance Model ! diagnostic tag is 'ebm' ! Add geopotential height to diagnostic fields ! ========================================================================================== module ebm_diagnostics_mod ! ! B.L. Samuels ! ! ! Zhi Liang ! ! ! Diagnostics for use in the Energy Balance Model ! ! ! Diagnostics for use in the Energy Balance Model ! ! ! diagnostic tag is 'ebm' ! use fms_mod, only: write_version_number use constants_mod, only: RADIAN use transforms_mod, only: get_grid_boundaries, get_deg_lon use transforms_mod, only: get_deg_lat, grid_domain use diag_manager_mod, only: diag_axis_init, register_diag_field use diag_manager_mod, only: register_static_field, send_data use time_manager_mod, only: time_type implicit none private public :: ebm_diagnostics_init, ebm_diagnostics_up, ebm_diagnostics_down character(len=84), parameter :: version = '$Id: ebm_diagnostics.F90,v 11.0 2004/09/28 19:10:27 fms Exp $' character(len=84), parameter :: tagname='$Name: mom4p1_pubrel_dec2009_nnz $' !----------------------------------------------------------------------- !------------------------- axis names ---------------------------------- character(len=8) :: axiset = 'ebm' character(len=8) :: mod_name = 'ebm' !----------------------------------------------------------------------- ! axis and field identifiers for the diag manager integer :: id_t, id_q, id_prec, id_lprec, id_fprec, id_dt_qg integer :: id_olr_toa, id_olr_boa, id_ilr_surf, id_nsw_toa, id_nsw_surf contains !----------------------------------------------------------------------------------------------------------------- ! ! ! ! setup/write netcdf metadata and static fields ! ! ! setup/write netcdf metadata and static fields ! ! ! ! Current/Init time. ! ! ! longitude of spectral atmopsheric grid ! ! ! latitude of spectral atmopsheric grid ! ! ! axes identifiers ! subroutine ebm_diagnostics_init(Time, lon_max, lat_max, axis_id) type(time_type), intent(in) :: Time integer, intent(in) :: lon_max, lat_max integer, intent(out),dimension(:) :: axis_id real, dimension(lon_max ) :: lon real, dimension(lon_max+1) :: lonb real, dimension(lat_max ) :: lat real, dimension(lat_max+1) :: latb integer :: axes_2d(2) integer :: log_unit, id_lonb, id_lon, id_latb, id_lat integer :: namelist_unit, ierr, io logical :: used !--- write out version information ------------------------------------- call write_version_number (version, tagname) call get_grid_boundaries(lonb,latb,global=.true.) call get_deg_lon(lon) call get_deg_lat(lat) id_lonb=diag_axis_init('lonb', RADIAN*lonb, 'degrees_E', 'x', 'longitude edges', set_name=axiset, Domain2=grid_domain) id_latb=diag_axis_init('latb', RADIAN*latb, 'degrees_N', 'y', 'latitude edges', set_name=axiset, Domain2=grid_domain) id_lon =diag_axis_init('lon', lon, 'degrees_E', 'x', 'longitude', set_name=axiset, Domain2=grid_domain, edges=id_lonb) id_lat =diag_axis_init('lat', lat, 'degrees_N', 'y', 'latitude', set_name=axiset, Domain2=grid_domain, edges=id_latb) axis_id(1) = id_lon axis_id(2) = id_lat axes_2d = axis_id(1:2) id_t = register_diag_field(mod_name, 'temp' , axes_2d, Time, 'temperature' , 'deg_k' ) id_q = register_diag_field(mod_name, 'sphum' , axes_2d, Time, 'specific humidity', 'none' ) id_prec = register_diag_field(mod_name, 'prec' , axes_2d, Time, 'precipitation' , 'Kg/(m2-s)') id_lprec = register_diag_field(mod_name, 'lprec',axes_2d, Time, 'rain' , 'Kg/(m2-s)') id_fprec = register_diag_field(mod_name, 'fprec' , axes_2d, Time, 'snow' ,' Kg/(m2-s)') id_dt_qg = register_diag_field(mod_name, 'dt_qg' , axes_2d, Time, 'advection' ,' 1/s') id_olr_toa = register_diag_field(mod_name, 'olr_toa', axes_2d, Time, 'outgoing longwave (to space)', 'W/m2') id_olr_boa = register_diag_field(mod_name, 'olr_boa', axes_2d, Time, 'outgoing longwave (to ocean)', 'W/m2') id_ilr_surf = register_diag_field(mod_name, 'ilr_surf', axes_2d, Time, 'incoming longwave (from ocean)', 'W/m2') id_nsw_toa = register_diag_field(mod_name, 'nsw_toa', axes_2d, Time, 'net shortwave (toa)', 'W/m2') id_nsw_surf = register_diag_field(mod_name, 'nsw_surf', axes_2d, Time, 'net shortwave (to ocean)', 'W/m2') return end subroutine ebm_diagnostics_init ! ! ! ! atmospheric temperature ! ! ! specific humidity ! ! ! total precipitation per unit ocean area ! ! ! rain per unit ocean area ! ! ! snow per unit ocean area ! ! ! advection ! ! ! outgoing longwave (to space) ! ! ! outgoing longwave (to ocean) ! ! ! incoming longwave (from ocean) ! ! ! net shortwave (toa) ! ! ! net shortwave (to ocean) ! ! !-------------------------------------------------------------------------------------------- ! ! ! ! Diagnostics passed from the "ocean" surface up thru the atmosphere ! ! ! Diagnostics passed from the "ocean" surface up thru the atmosphere ! ! ! ! ! time at the current time level (tau) ! ! ! atmospheric temperature ( deg_k ) ! ! ! specific humidity ( no units ) ! ! ! liquid precipitation rate (rain) in kg/m2/s ! ! ! frozen precipitation rate (snow) in kg/m2/s ! ! ! rate of change in specific humidity (1/s) ! subroutine ebm_diagnostics_up(Time,tg, qg, lprec, fprec, dt_qg ) type(time_type), intent(in) :: Time real, intent(in), dimension(:,:) :: tg, qg, lprec, fprec, dt_qg logical :: used if(id_t > 0) used = send_data(id_t , tg , time) if(id_q > 0) used = send_data(id_q , qg , time) if(id_prec > 0) used = send_data(id_prec , lprec+fprec , time) if(id_lprec > 0) used = send_data(id_lprec , lprec , time) if(id_fprec > 0) used = send_data(id_fprec , fprec , time) if(id_dt_qg > 0) used = send_data(id_dt_qg , dt_qg , time) return end subroutine ebm_diagnostics_up ! !-------------------------------------------------------------------------------------------- ! ! ! ! Diagnostics passed down from the atmosphere to the "ocean" surface ! ! ! Diagnostics passed down from the atmosphere to the "ocean" surface ! ! ! ! ! time at the current time level (tau) ! subroutine ebm_diagnostics_down(Time,olr_toa,olr_boa,ilr_surf,nsw_toa,nsw_surf) type(time_type), intent(in) :: Time real, intent(in), dimension(:,:) :: olr_toa,olr_boa,ilr_surf,nsw_toa,nsw_surf logical :: used if(id_olr_toa > 0) used = send_data(id_olr_toa ,olr_toa ,time) if(id_olr_boa > 0) used = send_data(id_olr_boa ,olr_boa ,time) if(id_ilr_surf > 0) used = send_data(id_ilr_surf,ilr_surf ,time) if(id_nsw_toa > 0) used = send_data(id_nsw_toa ,nsw_toa ,time) if(id_nsw_surf > 0) used = send_data(id_nsw_surf,nsw_surf ,time) return end subroutine ebm_diagnostics_down ! !-------------------------------------------------------------------------------------------- end module ebm_diagnostics_mod