1 | MODULE drag_noro_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | |
---|
5 | CONTAINS |
---|
6 | |
---|
7 | SUBROUTINE drag_noro (ngrid,nlayer,ptimestep,pplay,pplev,pvar, psig, pgam, pthe, & |
---|
8 | kgwd,kgwdim,kdx,ktest,t, u, v, & |
---|
9 | ! output |
---|
10 | pulow, pvlow, pustr, pvstr, d_t, d_u, d_v) |
---|
11 | |
---|
12 | !-------------------------------------------------------------------------------- |
---|
13 | ! MODULE contains DRAG_NORO SUBROUNTINE for sub-grid scale orographic sheme. |
---|
14 | ! 1) Initialization the variables |
---|
15 | ! 2) Overturn the levels and computes geopotential height |
---|
16 | ! 3) Call the and ORODRAG subroutine to updates the tendencies |
---|
17 | ! 4) Overturn back to the normal level of the tendencies and transfer it into |
---|
18 | ! increments and sum the oro-gw stress. |
---|
19 | ! the scheme has been called. |
---|
20 | ! Z.X. Li + F.Lott (LMD/CNRS/ENS) 1995-02-01 |
---|
21 | ! |
---|
22 | ! UPDATE: J.Liu 03/03/2022 Translate the code into .F90. The name of the |
---|
23 | ! variables are uniformed. Comments are made. |
---|
24 | ! Unused variables are deleted. |
---|
25 | !------------------------------------------------------------------------------- |
---|
26 | |
---|
27 | use dimradmars_mod, only: ndomainsz |
---|
28 | USE orodrag_mod, ONLY: orodrag |
---|
29 | USE comcstfi_h, ONLY: g, r |
---|
30 | |
---|
31 | IMPLICIT none |
---|
32 | |
---|
33 | ! 0. DECLARATIONS: |
---|
34 | |
---|
35 | ! 0.1 Input: |
---|
36 | INTEGER, intent(in):: ngrid ! Number of atmospheric columns [only nd] |
---|
37 | INTEGER, intent(in):: nlayer ! Number of atmospheric layers |
---|
38 | REAL, intent(in):: ptimestep ! Time step of the Physics(s) |
---|
39 | real, intent(in):: pplay(ndomainsz,nlayer) ! Pressure at full levels(Pa) |
---|
40 | real, intent(in):: pplev(ndomainsz,nlayer+1) ! Pressure at 1/2 levels(Pa) |
---|
41 | REAL, intent(in):: pvar(ndomainsz) ! sub-grid scale standard deviation |
---|
42 | REAL, intent(in):: psig(ndomainsz) ! sub-grid scale standard slope |
---|
43 | REAL, intent(in):: pgam(ndomainsz) ! sub-grid scale standard anisotropy |
---|
44 | REAL, intent(in):: pthe(ndomainsz) ! sub-grid scale principal axes angle |
---|
45 | REAL, intent(in):: u(ndomainsz,nlayer) ! Zonal wind at full levels(m/s) |
---|
46 | REAL, intent(in):: v(ndomainsz,nlayer) ! Meridional wind at full levels(m/s) |
---|
47 | REAL, intent(in):: t(ndomainsz,nlayer) ! Temperature at full levels(m/s) |
---|
48 | INTEGER, intent(in):: kgwd ! Number of points at which the scheme is called |
---|
49 | INTEGER, intent(in):: kgwdim ! kgwdim=MAX(1,kgwd) |
---|
50 | INTEGER, intent(in):: kdx(ndomainsz) ! Points at which to call the scheme |
---|
51 | INTEGER, intent(in):: ktest(ndomainsz) ! Map of calling points |
---|
52 | |
---|
53 | ! 0.2 Output: |
---|
54 | REAL, intent(out):: pulow(ndomainsz) ! Low level zonal wind |
---|
55 | REAL, intent(out):: pvlow(ndomainsz) ! Low level meridional wind |
---|
56 | REAL, intent(out):: pustr(ndomainsz) ! Low level zonal stress |
---|
57 | REAL, intent(out):: pvstr(ndomainsz) ! Low level meridional stress |
---|
58 | REAL, intent(out):: d_t(ndomainsz,nlayer) ! Temperature increment(K) due to orographic gravity waves |
---|
59 | REAL, intent(out):: d_u(ndomainsz,nlayer) ! Zonal increment(m/s) due to orographic gravity waves |
---|
60 | REAL, intent(out):: d_v(ndomainsz,nlayer) ! Meridional increment(m/s) due to orographic gravity waves |
---|
61 | |
---|
62 | ! 0.3 Variables locales: |
---|
63 | INTEGER i, k |
---|
64 | REAL inv_pplev(ndomainsz,nlayer+1) ! Inversed (by inverse the pplev pressure) pressure at 1/2 levels |
---|
65 | REAL inv_pplay(ndomainsz,nlayer) ! Inversed (by inverse the pplay pressure) pressure at full levels |
---|
66 | REAL zgeom(ndomainsz,nlayer) ! Geopotetial height (Inversed ??) |
---|
67 | REAL inv_pt(ndomainsz,nlayer) ! Inversed (by inverse the t) temperature (K) at full levels |
---|
68 | REAL inv_pu(ndomainsz,nlayer) ! Inversed (by inverse the u) zonal wind (m/s) at full levels |
---|
69 | REAL inv_pv(ndomainsz,nlayer) ! Inversed (by inverse the v) meridional wind (m/s) at full levels |
---|
70 | REAL pdtdt(ndomainsz,nlayer) ! Temperature tendency outputs from main routine ORODRAG |
---|
71 | REAL pdudt(ndomainsz,nlayer) ! Zonal winds tendency outputs from main routine ORODRAG |
---|
72 | REAL pdvdt(ndomainsz,nlayer) ! Meridional winds tendency outputs from main routine ORODRAG |
---|
73 | |
---|
74 | !----------------------------------------------------------------------------------------------------------------------- |
---|
75 | ! 1. INITIALISATIONS |
---|
76 | !----------------------------------------------------------------------------------------------------------------------- |
---|
77 | ! 1.1 Initialize the input/output variables (for security purpose) |
---|
78 | DO i = 1,ngrid |
---|
79 | pulow(i) = 0.0 |
---|
80 | pvlow(i) = 0.0 |
---|
81 | pustr(i) = 0.0 |
---|
82 | pvstr(i) = 0.0 |
---|
83 | ENDDO |
---|
84 | |
---|
85 | DO k = 1, nlayer |
---|
86 | DO i = 1, ngrid |
---|
87 | d_t(i,k) = 0.0 |
---|
88 | d_u(i,k) = 0.0 |
---|
89 | d_v(i,k) = 0.0 |
---|
90 | pdudt(i,k)=0.0 |
---|
91 | pdvdt(i,k)=0.0 |
---|
92 | pdtdt(i,k)=0.0 |
---|
93 | ENDDO |
---|
94 | ENDDO |
---|
95 | |
---|
96 | ! 1.2 Prepare the ininv_put varibales (Attention: The order of vertical levels increases from top to bottom ) |
---|
97 | ! Here the levels are overturned, that is, the surface becomes to the top and the top becomes to the surface |
---|
98 | ! and so forth |
---|
99 | DO k = 1, nlayer |
---|
100 | DO i = 1, ngrid |
---|
101 | inv_pt(i,k) = t(i,nlayer-k+1) |
---|
102 | inv_pu(i,k) = u(i,nlayer-k+1) |
---|
103 | inv_pv(i,k) = v(i,nlayer-k+1) |
---|
104 | inv_pplay(i,k) = pplay(i,nlayer-k+1) |
---|
105 | inv_pplev(i,k) = pplev(i,nlayer+1-k+1) |
---|
106 | ENDDO |
---|
107 | endDO |
---|
108 | |
---|
109 | DO i = 1, ngrid |
---|
110 | inv_pplev(i,nlayer+1) = pplev(i,1) |
---|
111 | ENDDO |
---|
112 | |
---|
113 | !calculate g*dz by R*T*[LN(P(z)/P(z+dz))] |
---|
114 | DO i = 1, ngrid |
---|
115 | zgeom(i,nlayer) = r * inv_pt(i,nlayer)* LOG(inv_pplev(i,nlayer+1)/inv_pplay(i,nlayer)) |
---|
116 | ENDDO |
---|
117 | |
---|
118 | ! sum g*dz from surface to top to get geopotential height |
---|
119 | DO k = nlayer-1, 1, -1 |
---|
120 | DO i = 1, ngrid |
---|
121 | zgeom(i,k) = zgeom(i,k+1) + r * (inv_pt(i,k)+inv_pt(i,k+1))/2.0 * LOG(inv_pplay(i,k+1)/inv_pplay(i,k)) |
---|
122 | ENDDO |
---|
123 | endDO |
---|
124 | |
---|
125 | !----------------------------------------------------------------------------------------------------------------------- |
---|
126 | ! 2. CALL the Main Rountine OORDRAG |
---|
127 | !----------------------------------------------------------------------------------------------------------------------- |
---|
128 | ! 2.1 call the main rountine to get the tendencies |
---|
129 | CALL ORODRAG(ngrid,nlayer,kgwd,kgwdim,kdx,ktest,ptimestep,inv_pplev,inv_pplay,zgeom,& |
---|
130 | inv_pt, inv_pu, inv_pv, pvar, psig, pgam, pthe, & |
---|
131 | ! output: |
---|
132 | pulow,pvlow,pdudt,pdvdt,pdtdt) |
---|
133 | |
---|
134 | ! 2.2 Transfer tendency into increment by multiply the physical time steps |
---|
135 | ! maybe in the future here cancel multiply the ptimestep thus it will not divide ptimestep again in the main rountine |
---|
136 | ! calldrag_noro_mod.F90 |
---|
137 | DO k = 1, nlayer |
---|
138 | DO i = 1, ngrid |
---|
139 | d_u(i,nlayer+1-k) = ptimestep*pdudt(i,k) |
---|
140 | d_v(i,nlayer+1-k) = ptimestep*pdvdt(i,k) |
---|
141 | d_t(i,nlayer+1-k) = ptimestep*pdtdt(i,k) |
---|
142 | pustr(i) = pustr(i)+g*pdudt(i,k)*(inv_pplev(i,k+1)-inv_pplev(i,k)) |
---|
143 | pvstr(i) = pvstr(i)+g*pdvdt(i,k)*(inv_pplev(i,k+1)-inv_pplev(i,k)) |
---|
144 | ENDDO |
---|
145 | ENDDO |
---|
146 | |
---|
147 | END SUBROUTINE drag_noro |
---|
148 | |
---|
149 | END MODULE drag_noro_mod |
---|