1 | MODULE physics_dcmip2016_mod |
---|
2 | USE ICOSA |
---|
3 | PRIVATE |
---|
4 | |
---|
5 | INTEGER,SAVE :: testcase |
---|
6 | !$OMP THREADPRIVATE(testcase) |
---|
7 | |
---|
8 | TYPE(t_field),POINTER :: f_out_i(:) |
---|
9 | REAL(rstd),POINTER :: out_i(:,:) |
---|
10 | |
---|
11 | TYPE(t_field),POINTER :: f_precl(:) |
---|
12 | REAL(rstd),ALLOCATABLE :: precl_packed(:) |
---|
13 | |
---|
14 | PUBLIC :: init_physics, full_physics, write_physics |
---|
15 | |
---|
16 | INTEGER, PARAMETER :: dry_baroclinic=0 |
---|
17 | INTEGER, PARAMETER :: moist_baroclinic=1 |
---|
18 | INTEGER, PARAMETER :: cyclone=2 |
---|
19 | INTEGER, PARAMETER :: supercell=3 |
---|
20 | |
---|
21 | |
---|
22 | CONTAINS |
---|
23 | |
---|
24 | SUBROUTINE init_physics |
---|
25 | USE physics_interface_mod |
---|
26 | IMPLICIT NONE |
---|
27 | INTEGER :: ngrid |
---|
28 | CHARACTER(LEN=255) :: testcase_str |
---|
29 | |
---|
30 | CALL getin("physics_dcmip2016",testcase_str) |
---|
31 | |
---|
32 | SELECT CASE (TRIM(testcase_str)) |
---|
33 | CASE ('dry_baroclinic') |
---|
34 | testcase=dry_baroclinic |
---|
35 | CASE ('moist_baroclinic') |
---|
36 | testcase=moist_baroclinic |
---|
37 | CASE ('cyclone') |
---|
38 | testcase=cyclone |
---|
39 | CASE ('supercell') |
---|
40 | testcase=supercell |
---|
41 | CASE DEFAULT |
---|
42 | PRINT*, 'Bad selector for dcmip2016 test case <', testcase_str, & |
---|
43 | '> options are <dry_baroclinic>, <moist_baroclinic>, <cyclone>, <supercell>' |
---|
44 | STOP |
---|
45 | END SELECT |
---|
46 | |
---|
47 | |
---|
48 | ngrid = physics_inout%ngrid |
---|
49 | ! Input |
---|
50 | ALLOCATE(physics_inout%Ai(ngrid)) |
---|
51 | ALLOCATE(physics_inout%lon(ngrid)) |
---|
52 | ALLOCATE(physics_inout%lat(ngrid)) |
---|
53 | ALLOCATE(physics_inout%phis(ngrid)) |
---|
54 | ALLOCATE(physics_inout%p(ngrid,llm+1)) |
---|
55 | ALLOCATE(physics_inout%pk(ngrid,llm)) |
---|
56 | ALLOCATE(physics_inout%Temp(ngrid,llm)) |
---|
57 | ALLOCATE(physics_inout%ulon(ngrid,llm)) |
---|
58 | ALLOCATE(physics_inout%ulat(ngrid,llm)) |
---|
59 | ALLOCATE(physics_inout%q(ngrid,llm,nqtot)) |
---|
60 | ! Output (tendencies) |
---|
61 | ALLOCATE(physics_inout%dTemp(ngrid,llm)) |
---|
62 | ALLOCATE(physics_inout%dulon(ngrid,llm)) |
---|
63 | ALLOCATE(physics_inout%dulat(ngrid,llm)) |
---|
64 | ALLOCATE(physics_inout%dq(ngrid,llm,nqtot)) |
---|
65 | ! Physics-specific data |
---|
66 | ALLOCATE(precl_packed(ngrid)) |
---|
67 | CALL allocate_field(f_precl, field_t,type_real) |
---|
68 | |
---|
69 | PRINT *, 'init_physics_new', SIZE(physics_inout%Ai) |
---|
70 | END SUBROUTINE init_physics |
---|
71 | |
---|
72 | SUBROUTINE full_physics |
---|
73 | USE physics_interface_mod |
---|
74 | CALL compute_physics(physics_inout%ngrid, physics_inout%dt_phys, & |
---|
75 | physics_inout%lon, physics_inout%lat, physics_inout%p, physics_inout%pk, physics_inout%Temp, & |
---|
76 | physics_inout%ulon, physics_inout%ulat, physics_inout%q(:,:,1:5), & |
---|
77 | physics_inout%dTemp, physics_inout%dulon, physics_inout%dulat, & |
---|
78 | physics_inout%dq(:,:,1:5), precl_packed) |
---|
79 | END SUBROUTINE full_physics |
---|
80 | |
---|
81 | SUBROUTINE write_physics |
---|
82 | USE output_field_mod |
---|
83 | USE physics_interface_mod |
---|
84 | CALL unpack_field(f_precl, precl_packed) |
---|
85 | CALL output_field("precl",f_precl) |
---|
86 | |
---|
87 | END SUBROUTINE write_physics |
---|
88 | |
---|
89 | SUBROUTINE compute_physics(ngrid,dt_phys,lon, lat, p, pk, Temp,u,v,q, dTemp,du,dv,dq, precl) |
---|
90 | USE icosa |
---|
91 | USE dcmip2016_simple_physics_mod |
---|
92 | USE dcmip2016_kessler_physic_mod |
---|
93 | USE earth_const |
---|
94 | USE terminator |
---|
95 | IMPLICIT NONE |
---|
96 | INTEGER :: ngrid |
---|
97 | REAL(rstd) :: lat(ngrid) |
---|
98 | REAL(rstd) :: lon(ngrid) |
---|
99 | REAL(rstd) :: ps(ngrid) |
---|
100 | REAL(rstd) :: precl(ngrid) |
---|
101 | ! arguments with bottom-up indexing (DYNAMICO) |
---|
102 | REAL(rstd) :: p(ngrid,llm+1) |
---|
103 | REAL(rstd) :: pk(ngrid,llm) |
---|
104 | REAL(rstd) :: Temp(ngrid,llm) |
---|
105 | REAL(rstd) :: u(ngrid,llm) |
---|
106 | REAL(rstd) :: v(ngrid,llm) |
---|
107 | REAL(rstd) :: q(ngrid,llm,5) |
---|
108 | REAL(rstd) :: dTemp(ngrid,llm) |
---|
109 | REAL(rstd) :: du(ngrid,llm) |
---|
110 | REAL(rstd) :: dv(ngrid,llm) |
---|
111 | REAL(rstd) :: dq(ngrid,llm,5) |
---|
112 | ! local arrays with top-down vertical indexing (DCMIP) |
---|
113 | REAL(rstd) :: pint(ngrid,llm+1) |
---|
114 | REAL(rstd) :: pmid(ngrid,llm) |
---|
115 | REAL(rstd) :: pdel(ngrid,llm) |
---|
116 | REAL(rstd) :: Tfi(ngrid,llm) |
---|
117 | REAL(rstd) :: ufi(ngrid,llm) |
---|
118 | REAL(rstd) :: vfi(ngrid,llm) |
---|
119 | REAL(rstd) :: qfi(ngrid,llm,5) |
---|
120 | |
---|
121 | REAL(rstd) :: rho(llm), z(llm), theta(llm), qv(llm),qc(llm),qr(llm) |
---|
122 | REAL(rstd) :: lastz |
---|
123 | REAL(rstd) :: dcl1,dcl2 |
---|
124 | INTEGER :: l,ll,ij |
---|
125 | REAL(rstd) :: dt_phys, inv_dt |
---|
126 | INTEGER :: simple_physic_testcase |
---|
127 | |
---|
128 | ! prepare input fields and mirror vertical index |
---|
129 | ps(:) = p(:,1) ! surface pressure |
---|
130 | |
---|
131 | DO l=1,llm+1 |
---|
132 | DO ij=1,ngrid |
---|
133 | pint(ij,l)=p(ij,llm+2-l) |
---|
134 | ENDDO |
---|
135 | ENDDO |
---|
136 | |
---|
137 | DO l=1,llm |
---|
138 | ll=llm+1-l |
---|
139 | DO ij=1,ngrid |
---|
140 | pmid(ij,l)=0.5*(pint(ij,l)+pint(ij,l+1)) ! Pressure inside layers |
---|
141 | pdel(ij,l)=pint(ij,l+1)-pint(ij,l) ! Pressure difference between two layers |
---|
142 | ufi(ij,l)=u(ij,ll) |
---|
143 | vfi(ij,l)=v(ij,ll) |
---|
144 | qfi(ij,l,:)=q(ij,ll,:) |
---|
145 | IF (physics_thermo==thermo_fake_moist) THEN |
---|
146 | Tfi(ij,l)=Temp(ij,ll)/(1+0.608*qfi(ij,l,1)) |
---|
147 | ELSE |
---|
148 | Tfi(ij,l)=Temp(ij,ll) |
---|
149 | ENDIF |
---|
150 | ENDDO |
---|
151 | ENDDO |
---|
152 | |
---|
153 | |
---|
154 | IF (testcase==moist_baroclinic .OR. testcase==cyclone ) THEN |
---|
155 | IF (testcase==moist_baroclinic) simple_physic_testcase=1 |
---|
156 | IF (testcase==cyclone) simple_physic_testcase=0 |
---|
157 | CALL simple_physics(ngrid, llm, dt_phys, lat, tfi, qfi(:,:,1) , ufi, vfi, pmid, pint, pdel, 1./pdel, ps, precl, & |
---|
158 | simple_physic_testcase, .FALSE., .FALSE.) |
---|
159 | ENDIF |
---|
160 | |
---|
161 | |
---|
162 | IF (testcase==moist_baroclinic .OR. testcase==cyclone .OR. testcase==supercell ) THEN |
---|
163 | DO ij=1,ngrid |
---|
164 | lastz=0 |
---|
165 | DO l=1,llm |
---|
166 | ll=llm+1-l |
---|
167 | rho(l) = pmid(ij,ll)/(287*Temp(ij,l)) |
---|
168 | z(l)=lastz |
---|
169 | lastz=lastz+ (p(ij,l)-p(ij,l+1)) /g / rho(l) |
---|
170 | theta(l)= Tfi(ij,ll) / ( pk(ij,l) / cpp) |
---|
171 | ENDDO |
---|
172 | |
---|
173 | qv(:)=qfi(ij,llm:1:-1,1) |
---|
174 | qc(:)=qfi(ij,llm:1:-1,2) |
---|
175 | qr(:)=qfi(ij,llm:1:-1,3) |
---|
176 | |
---|
177 | CALL KESSLER(theta(:), qv, qc, qr, rho(:), & |
---|
178 | pk(ij,:)/cpp, dt_phys, z(:), llm, precl(ij)) |
---|
179 | |
---|
180 | |
---|
181 | DO l=1,llm |
---|
182 | ll=llm+1-l |
---|
183 | Tfi(ij,ll) = theta(l) * ( pk(ij,l) / cpp) |
---|
184 | ENDDO |
---|
185 | |
---|
186 | qfi(ij,:,1)=qv(llm:1:-1) |
---|
187 | qfi(ij,:,2)=qc(llm:1:-1) |
---|
188 | qfi(ij,:,3)=qr(llm:1:-1) |
---|
189 | |
---|
190 | ENDDO |
---|
191 | ENDIF |
---|
192 | |
---|
193 | DO l=1,llm |
---|
194 | ll=llm+1-l |
---|
195 | DO ij=1,ngrid |
---|
196 | CALL tendency_terminator( lat(ij), lon(ij), qfi(ij,ll,4), qfi(ij,ll,5), dt_phys, dcl1, dcl2) |
---|
197 | qfi(ij,ll,4)=qfi(ij,ll,4)+ dt_phys*dcl1 |
---|
198 | qfi(ij,ll,5)=qfi(ij,ll,5)+ dt_phys*dcl2 |
---|
199 | ENDDO |
---|
200 | ENDDO |
---|
201 | |
---|
202 | |
---|
203 | ! Mirror vertical index and compute tendencies |
---|
204 | inv_dt = 1./dt_phys |
---|
205 | DO l=1,llm |
---|
206 | ll=llm+1-l |
---|
207 | DO ij=1,ngrid |
---|
208 | IF (physics_thermo==thermo_fake_moist) THEN |
---|
209 | dTemp(ij,l) = inv_dt * ( Tfi(ij,ll)*(1+0.608*qfi(ij,ll,1)) - Temp(ij,l) ) |
---|
210 | ELSE |
---|
211 | dTemp(ij,l) = inv_dt * ( Tfi(ij,ll) - Temp(ij,l) ) |
---|
212 | ENDIF |
---|
213 | |
---|
214 | du(ij,l) = inv_dt * (ufi(ij,ll) - u(ij,l)) |
---|
215 | dv(ij,l) = inv_dt * (vfi(ij,ll) - v(ij,l)) |
---|
216 | dq(ij,l,:) = inv_dt * (qfi(ij,ll,:) - q(ij,l,:)) |
---|
217 | ENDDO |
---|
218 | ENDDO |
---|
219 | |
---|
220 | END SUBROUTINE compute_physics |
---|
221 | |
---|
222 | END MODULE physics_dcmip2016_mod |
---|
223 | |
---|
224 | |
---|