!! SUBROUTINE ana_humid (ng, tile, model) ! !! git $Id$ !! svn $Id: ana_humid.h 1151 2023-02-09 03:08:53Z arango $ !!====================================================================== !! Copyright (c) 2002-2023 The ROMS/TOMS Group ! !! Licensed under a MIT/X style license ! !! See License_ROMS.md ! !======================================================================= ! ! ! This routine sets surface air humidity (moisture) using an ! ! analytical expression. There three types of humidity: ! ! ! ! 1) Absolute humidity: density of water vapor. ! ! 2) Specific humidity: ratio of the mass of water vapor to ! ! the mass of moist air cointaining the vapor (g/kg) ! ! 3) Relative humidity: ratio of the actual mixing ratio to ! ! saturation mixing ratio of the air at given temperature ! ! and pressure (percentage). ! ! ! !======================================================================= ! USE mod_param USE mod_forces USE mod_ncparam ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & __FILE__ ! #include "tile.h" ! CALL ana_humid_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & FORCES(ng) % Hair) ! ! Set analytical header file name used. ! #ifdef DISTRIBUTE IF (Lanafile) THEN #else IF (Lanafile.and.(tile.eq.0)) THEN #endif ANANAME( 9)=MyFile END IF ! RETURN END SUBROUTINE ana_humid ! !*********************************************************************** SUBROUTINE ana_humid_tile (ng, tile, model, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & Hair) !*********************************************************************** ! USE mod_param USE mod_scalars ! USE exchange_2d_mod, ONLY : exchange_r2d_tile #ifdef DISTRIBUTE USE mp_exchange_mod, ONLY : mp_exchange2d #endif ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile, model integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS ! #ifdef ASSUMED_SHAPE real(r8), intent(out) :: Hair(LBi:,LBj:) #else real(r8), intent(out) :: Hair(LBi:UBi,LBj:UBj) #endif ! ! Local variable declarations. ! integer :: i, j #include "set_bounds.h" ! !----------------------------------------------------------------------- ! Set analytical surface air humidity. !----------------------------------------------------------------------- ! #if defined BENCHMARK DO j=JstrT,JendT DO i=IstrT,IendT Hair(i,j)=0.8_r8 END DO END DO #elif defined BL_TEST DO j=JstrT,JendT DO i=IstrT,IendT Hair(i,j)=0.776_r8 END DO END DO #else ana_humidity.h: no values provided for Hair. #endif ! ! Exchange boundary data. ! IF (EWperiodic(ng).or.NSperiodic(ng)) THEN CALL exchange_r2d_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & Hair) END IF #ifdef DISTRIBUTE CALL mp_exchange2d (ng, tile, model, 1, & & LBi, UBi, LBj, UBj, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & Hair) #endif ! RETURN END SUBROUTINE ana_humid_tile