📄 dynamics_vars.f90
字号:
#include <misc.h>#include <params.h>module dynamics_vars!BOP!! !MODULE: dynamics_vars --- Lin-Rood specific variables and methods!! !USES: use precision, only : r8 use pmgrid! !PUBLIC MEMBER FUNCTIONS: public dynamics_init, dynamics_clean! !PUBLIC DATA MEMBERS:!--------------------------------------------------------! Local variables specific to the Lin-Rood dynamical core!--------------------------------------------------------! Temporary variables to be removed at a later time real(r8), allocatable, save :: xyt(:,:,:) real(r8), allocatable, save :: yzt(:,:,:) real(r8), allocatable, save :: q3t(:,:,:,:)! Variables set by dynamics_init real(r8) dtime ! large time step integer iord ! order of LR scheme in X integer jord ! order of LR scheme in Y! General variables set by setrig real (r8) pi real(r8) dl ! Radians per (finite-volume) longitude real(r8) dp ! Radians per (finite-volume) latitude! Geometric arrays set by setrig real(r8) cosp(plat) ! Cos of latitude angle - volume mean real(r8) sinp(plat) ! Sin of latitude angle - volume mean real(r8) cose(plat) ! Cos of finite-volume edges real(r8) sine(plat) ! Sin of finite-volume edges real(r8) gw(plat) ! Delta sine! Variables set by set_eta real(r8) ptop ! pressure at top of atmosphere real(r8) pint ! pressure at sig-p interface real(r8) ak(plev+1) ! A's of the ETA-coordinate real(r8) bk(plev+1) ! B's of the ETA-coordinate integer ks ! Total number of pure-P layers! Variables set by rayf_init real(r8) rfac(plev) ! Rayleigh friction factor! Needed in Held-Suarez (hswf) set by hswf_init real (r8) sinp2(plat) real (r8) cosp2(plat) real (r8) cosp4(plat) real (r8) rf(plev)! Scalars set in dynpkg_init integer icd, jcd integer ng_c ! ghost zone needed by the c-gird dynamics integer ng_d ! ghost zone needed by the d-gird dynamics integer ng_s ! for certain arrays, max(ng_c+1,ng_d) real (r8) acap ! scaled polar cap area real (r8) rcap ! inverse of scaled polar cap area! Arrays initialized by dynpkg_init real(r8) coslon(plon) ! Cos of longitudes - volume center real(r8) sinlon(plon) ! Sin of longitudes - volume center real(r8) cosl5(plon) ! Cos of longitudes - volume center real(r8) sinl5(plon) ! Sin of longitudes - volume center real(r8) acosp(plat)! Scalars initialized by d_split integer ns ! total number of splits for Lagrangian dynamics! Constants -- these might be user-configurable at a later date logical rayf ! Rayleigh friction flag (off by default) parameter (rayf = .false.) ! off logical dcaf ! Dry convection flag (use only in ideal_phys case) parameter (dcaf = .false.) ! off!! !DESCRIPTION:!! This module provides variables which are specific to the Lin-Rood! dynamical core. Most of them were previously SAVE variables in ! different routines and were set with an "if (first)" statement.!! \begin{tabular}{|l|l|} \hline \hline! lr\_init & Initialize the Lin-Rood variables \\ \hline! lr\_clean & Deallocate all internal data structures \\ \hline ! \hline! \end{tabular}!! !REVISION HISTORY:! 01.06.06 Sawyer Consolidated from various code snippets! 01.07.12 Sawyer Removed CCM common blocks comtim.h and commap.h!! !BUGS:! o rayf_init and hswf_init should only be performed if they are used....! o Where does the value of ns0 come from??! o Should gw be here?? It's available in commap.h (physics variable)! and set in stepon. Needed in gmean.!!EOP!-----------------------------------------------------------------------contains!-----------------------------------------------------------------------!BOP! !IROUTINE: dynamics_init --- initialize the lin-rood dynamical core!! !INTERFACE: subroutine dynamics_init( dtime_in, iord_in, jord_in, nsplit_in )! !USES: use constituents, only : pcnst, pnats ! TEMPORARY!!!! implicit none! !INPUT PARAMETERS: real (r8), intent(in) :: dtime_in ! Large time step integer, intent(in) :: iord_in ! Order of LR scheme in X integer, intent(in) :: jord_in ! Order of LR scheme in Y integer, intent(in) :: nsplit_in ! Small time steps in dtime! !DESCRIPTION:!! Initialize Lin-Rood specific variables!! !REVISION HISTORY:!! 01.06.06 Sawyer Create!!EOP!-----------------------------------------------------------------------!BOC! dtime = dtime_in iord = iord_in jord = jord_in call setrig call set_eta if ( rayf ) call rayf_init call hswf_init call dynpkg_init call d_split(nsplit_in)#if defined( SPMD )! Temporary data structures allocate( xyt(beglonxy:endlonxy,beglatxy:endlatxy,plev) ) allocate( yzt(plon,beglat:endlat,beglev:endlev) ) allocate( q3t(plon,beglat:endlat,beglev:endlev,pcnst+pnats) )#endif return!EOC end subroutine dynamics_init!-----------------------------------------------------------------------!-----------------------------------------------------------------------!BOP! !IROUTINE: dynamics_clean -- clean up Lin-Rood-specific variables!! !INTERFACE: subroutine dynamics_clean implicit none! !DESCRIPTION:!! Clean up (deallocate) Lin-Rood-specific variables!! !REVISION HISTORY:!! 01.06.06 Sawyer Creation!!EOP!-----------------------------------------------------------------------!BOC!! Temporary data structures deallocate( q3t ) deallocate( yzt ) deallocate( xyt ) return!EOC end subroutine dynamics_clean!-----------------------------------------------------------------------!----------------------------------------------------------------------- !BOP! !ROUTINE: setrig --- Specify the grid attributes!! !INTERFACE: subroutine setrig! !USES: implicit none!! !DESCRIPTION:!! Specify the grid attributes, such as the spacing between! grid points in latitude and longitude, the sines and cosines of! latitude at cell midpoints and edges.!! !REVISION HISTORY: ! ??.??.?? Lin? Creation! 01.03.26 Sawyer Added ProTeX documentation!!EOP!-----------------------------------------------------------------------!BOC!! !LOCAL VARIABLES: integer j real*8 pi, ph5 ! This is to ensure 64-bit for any choice of r8 pi = 4.0 * atan(1.) dl = (pi+pi)/plon dp = pi/(plat-1) do j=2,plat ph5 = -0.5d0*pi + ((j-1)-0.5d0)*(pi/(plat-1)) sine(j) = sin(ph5) enddo do j=2,plat-1 gw(j) = sine(j+1) - sine(j) end do gw( 1) = 1. + sine(2) gw(plat) = 1. - sine(plat) cosp( 1) = 0. cosp(plat) = 0. do j=2,plat-1 cosp(j) = (sine(j+1)-sine(j)) / dp enddo! Define cosine at edges.. do j=2,plat cose(j) = 0.5 * (cosp(j-1) + cosp(j)) enddo cose(1) = cose(2) sinp( 1) = -1. sinp(plat) = 1. do j=2,plat-1 sinp(j) = 0.5 * (sine(j) + sine(j+1)) enddo return!EOC end subroutine setrig!-----------------------------------------------------------------------!----------------------------------------------------------------------- !BOP! !ROUTINE: set_eta --- Define vertical coordinate!! !INTERFACE: subroutine set_eta! !USES: implicit none!! !DESCRIPTION:!! Specify the vertical coordinate system. Currently this is a ! dual pressure - sigma coordinate system (???), which transitions at! level ks, but it could be just about anything reasonable.!! Choices for vertical resolutions are as follows:! \begin{tabular}{l}! NCAR: 18, 26, and 30 \\! NASA DAO: smoothed version of CAM's 30-level, 32, 55, 64, and 96 \\! New 32-layer setup with top at 0.4 mb for high horizontal! resolution runs. pint = 176.93 mb \\! Revised 55-level eta with pint at 176.93 mb SJL: 2000-03-20! \end{tabular}! ! !REVISION HISTORY: ! 98.01.15 Lin Creation! ongoing Lin Fine tuning! 01.03.26 Sawyer Added ProTeX documentation!!EOP!-----------------------------------------------------------------------!BOC!! !LOCAL VARIABLES:! NCAR specific real(r8) a18(19),b18(19) ! CCM3 real(r8) a26(27),b26(27) ! CAM real(r8) a30(31),b30(31) ! CCM3! NASA only real(r8) a30m(31),b30m(31) ! smoothed CAM 30-L real(r8) a32(33),b32(33) real(r8) a32old(33),b32old(33) real(r8) a55(56),b55(56) real(r8) a55old(56),b55old(56) real(r8) a64(65),b64(65) real(r8) a96(97),b96(97)
⌨️ 快捷键说明
复制代码
Ctrl + C
搜索代码
Ctrl + F
全屏模式
F11
切换主题
Ctrl + Shift + D
显示快捷键
?
增大字号
Ctrl + =
减小字号
Ctrl + -