source: codes/icosagcm/devel/src/dcmip/physics_dcmip2016.f90 @ 739

Last change on this file since 739 was 739, checked in by dubos, 6 years ago

devel : small cleanup in idealized physics

File size: 8.6 KB
Line 
1MODULE physics_dcmip2016_mod
2  USE ICOSA
3  USE caldyn_vars_mod
4  PRIVATE
5 
6  INTEGER,SAVE :: testcase
7!$OMP THREADPRIVATE(testcase)
8
9  TYPE(t_field),SAVE,POINTER :: f_out_i(:)
10  REAL(rstd),SAVE,POINTER :: out_i(:,:)
11
12  TYPE(t_field),SAVE,POINTER  :: f_precl(:)
13  REAL(rstd),SAVE,ALLOCATABLE :: precl_packed(:)
14!$OMP THREADPRIVATE(precl_packed) 
15
16  TYPE(t_field),SAVE,POINTER  :: f_Q1(:)
17  TYPE(t_field),SAVE,POINTER  :: f_Q2(:)
18  TYPE(t_field),SAVE,POINTER  :: f_PS(:)
19  TYPE(t_field),SAVE,POINTER  :: f_rhodz(:)
20  TYPE(t_field),SAVE,POINTER  :: f_Q1_col_int(:)
21  TYPE(t_field),SAVE,POINTER  :: f_Q2_col_int(:)
22  PUBLIC :: init_physics, full_physics, write_physics
23
24  INTEGER, PARAMETER :: dry_baroclinic=0
25  INTEGER, PARAMETER :: moist_baroclinic_full=1
26  INTEGER, PARAMETER :: moist_baroclinic_kessler=2
27  INTEGER, PARAMETER :: cyclone=3
28  INTEGER, PARAMETER :: supercell=4
29
30  LOGICAL,SAVE :: PBL   !  boundary layer
31                        !  True : George Bryan
32                        !  False : Reed & Jablonowsi
33  !$OMP THREADPRIVATE(PBL)
34CONTAINS
35
36  SUBROUTINE init_physics
37    USE physics_interface_mod
38    USE omp_para
39    IMPLICIT NONE
40    INTEGER :: ngrid
41    CHARACTER(LEN=255) :: testcase_str 
42   
43    testcase_str='undefined'
44    CALL getin("physics_dcmip2016",testcase_str)
45   
46    SELECT CASE (TRIM(testcase_str))
47      CASE ('dry_baroclinic')
48        testcase=dry_baroclinic
49      CASE ('moist_baroclinic_full') 
50        testcase=moist_baroclinic_full
51      CASE ('moist_baroclinic_kessler') 
52        testcase=moist_baroclinic_kessler
53      CASE ('cyclone') 
54        testcase=cyclone
55      CASE ('supercell') 
56        testcase=supercell
57      CASE DEFAULT
58         PRINT*, 'Bad selector for dcmip2016 test case <', testcase_str, &
59              '> options are <dry_baroclinic>, <moist_baroclinic_full>, <moist_baroclinic_kessler>, <cyclone>, <supercell>'
60         STOP
61    END SELECT
62
63    PBL=.FALSE.
64    CALL getin("physics_dcmip2016_PBL",PBL)
65    IF(is_master) THEN
66       IF(PBL) THEN
67          PRINT *, 'PBL=.TRUE., Bryan PBL activated'
68       ELSE
69          PRINT *, 'PBL=.FALSE., Reed & Jablonowski PBL activated'
70       END IF
71    END IF
72
73    ! Physics-specific data
74    ngrid = physics_inout%ngrid
75    ALLOCATE(precl_packed(ngrid))
76    CALL allocate_field(f_precl, field_t,type_real)
77    CALL allocate_field(f_Q1, field_t,type_real,llm)
78    CALL allocate_field(f_Q2, field_t,type_real,llm)
79    CALL allocate_field(f_PS, field_t,type_real)
80    CALL allocate_field(f_rhodz, field_t,type_real,llm)
81    CALL allocate_field(f_Q1_col_int, field_t,type_real)
82    CALL allocate_field(f_Q2_col_int, field_t,type_real)
83
84    PRINT *, 'init_physics_new', SIZE(physics_inout%Ai)
85  END SUBROUTINE init_physics
86
87  SUBROUTINE full_physics
88    USE physics_interface_mod
89    CALL compute_physics(physics_inout%ngrid, physics_inout%dt_phys, &
90         physics_inout%lon, physics_inout%lat, physics_inout%p, physics_inout%pk, physics_inout%Temp, & 
91         physics_inout%ulon, physics_inout%ulat, physics_inout%q(:,:,1:5), &
92         physics_inout%dTemp, physics_inout%dulon, physics_inout%dulat, &
93         physics_inout%dq(:,:,1:5), precl_packed)
94  END SUBROUTINE full_physics
95
96  SUBROUTINE write_physics
97    USE output_field_mod
98    USE physics_interface_mod
99    use disvert_mod
100    REAL(rstd), POINTER :: Q1(:,:)
101    REAL(rstd), POINTER :: Q2(:,:)
102    REAL(rstd), POINTER :: PS(:)
103    REAL(rstd), POINTER :: rhodz(:,:)
104    REAL(rstd), POINTER :: Q1_col_int(:)
105    REAL(rstd), POINTER :: Q2_col_int(:)
106   
107   
108    CALL unpack_field(f_precl, precl_packed)
109    CALL output_field("precl",f_precl)
110   
111    CALL unpack_field(f_Q1,physics_inout%q(:,:,4))
112    CALL unpack_field(f_Q2,physics_inout%q(:,:,5))
113    CALL unpack_field(f_ps,physics_inout%p(:,1))
114   
115    DO ind=1,ndomain
116      IF (.NOT. assigned_domain(ind)) CYCLE
117        Q1=f_Q1(ind)
118        Q2=f_Q2(ind)
119        Q1_col_int=f_Q1_col_int(ind)
120        Q2_col_int=f_Q2_col_int(ind)
121        PS=f_PS(ind)
122        rhodz=f_rhodz(ind)
123        DO l=1,llm
124          rhodz(:,l)= ( ap(l)-ap(l+1) + (bp(l)-bp(l+1))*ps(:) )/g   
125        ENDDO
126        Q1_col_int=SUM(Q1*rhodz,2)/SUM(rhodz,2)
127        Q2_col_int=SUM(Q2*rhodz,2)/SUM(rhodz,2)
128    ENDDO
129     
130    CALL output_field("Q1_col_int",f_Q1_col_int)
131    CALL output_field("Q2_col_int",f_Q2_col_int)
132   
133  END SUBROUTINE write_physics
134
135  SUBROUTINE compute_physics(ngrid,dt_phys,lon, lat, p, pk, Temp,u,v,q, dTemp,du,dv,dq, precl)
136    USE icosa
137    USE dcmip2016_simple_physics_mod
138    USE dcmip2016_kessler_physic_mod
139    USE earth_const
140    USE terminator
141    IMPLICIT NONE
142    INTEGER    :: ngrid
143    REAL(rstd) :: lat(ngrid)
144    REAL(rstd) :: lon(ngrid)
145    REAL(rstd) :: ps(ngrid)
146    REAL(rstd) :: precl(ngrid)
147    ! arguments with bottom-up indexing (DYNAMICO)
148    REAL(rstd) :: p(ngrid,llm+1)
149    REAL(rstd) :: pk(ngrid,llm)
150    REAL(rstd) :: Temp(ngrid,llm)
151    REAL(rstd) :: u(ngrid,llm)
152    REAL(rstd) :: v(ngrid,llm)
153    REAL(rstd) :: q(ngrid,llm,5)
154    REAL(rstd) :: dTemp(ngrid,llm)
155    REAL(rstd) :: du(ngrid,llm)
156    REAL(rstd) :: dv(ngrid,llm)
157    REAL(rstd) :: dq(ngrid,llm,5)
158    ! local arrays with top-down vertical indexing (DCMIP)
159    REAL(rstd) :: pint(ngrid,llm+1)
160    REAL(rstd) :: pmid(ngrid,llm)
161    REAL(rstd) :: pdel(ngrid,llm)
162    REAL(rstd) :: Tfi(ngrid,llm)
163    REAL(rstd) :: ufi(ngrid,llm)
164    REAL(rstd) :: vfi(ngrid,llm)
165    REAL(rstd) :: qfi(ngrid,llm,5)
166
167    REAL(rstd)  :: rho(llm), z(llm), theta(llm), qv(llm),qc(llm),qr(llm)
168    REAL(rstd)  :: lastz
169    REAL(rstd)  :: dcl1,dcl2
170     INTEGER :: l,ll,ij
171    REAL(rstd) :: dt_phys, inv_dt
172    INTEGER :: simple_physic_testcase
173   
174    ! prepare input fields and mirror vertical index     
175    ps(:) = p(:,1) ! surface pressure
176
177    DO l=1,llm+1
178      DO ij=1,ngrid
179          pint(ij,l)=p(ij,llm+2-l)
180      ENDDO
181    ENDDO
182
183    DO l=1,llm
184       ll=llm+1-l
185       DO ij=1,ngrid
186          pmid(ij,l)=0.5*(pint(ij,l)+pint(ij,l+1)) ! Pressure inside layers
187          pdel(ij,l)=pint(ij,l+1)-pint(ij,l)       ! Pressure difference between two layers
188          ufi(ij,l)=u(ij,ll)
189          vfi(ij,l)=v(ij,ll)
190          qfi(ij,l,:)=q(ij,ll,:)
191          IF (physics_thermo==thermo_fake_moist) THEN
192            Tfi(ij,l)=Temp(ij,ll)/(1+0.608*qfi(ij,l,1)) 
193          ELSE
194            Tfi(ij,l)=Temp(ij,ll)
195          ENDIF
196       ENDDO
197    ENDDO
198   
199    precl=0.
200    IF (testcase==moist_baroclinic_full .OR. testcase==cyclone  ) THEN
201      IF (testcase==moist_baroclinic_full) simple_physic_testcase=1
202      IF (testcase==cyclone) simple_physic_testcase=0
203      CALL simple_physics(ngrid, llm, dt_phys, lat, tfi, qfi(:,:,1) , ufi, vfi, pmid, pint, pdel, 1./pdel, ps, precl, &
204                          simple_physic_testcase, .FALSE., PBL)
205    ENDIF
206
207 
208    IF (testcase==moist_baroclinic_full .OR. testcase==moist_baroclinic_kessler .OR. testcase==cyclone .OR. testcase==supercell ) THEN
209       DO ij=1,ngrid
210          lastz=0 
211          DO l=1,llm
212           ll=llm+1-l
213           rho(l) = pmid(ij,ll)/(287*Tfi(ij,ll)*(1+0.608*qfi(ij,ll,1)))
214           z(l)=lastz
215           lastz=lastz+ (p(ij,l)-p(ij,l+1)) /g / rho(l)
216           theta(l)= Tfi(ij,ll) / ( pk(ij,l) / cpp)
217          ENDDO
218         
219          DO l=1,llm-1
220           z(l)= 0.5*(z(l)+z(l+1))
221          ENDDO
222          z(llm)=z(llm)+(z(llm)-z(llm-1))
223         
224          qv(:)=max(qfi(ij,llm:1:-1,1),0.)
225          qc(:)=max(qfi(ij,llm:1:-1,2),0.)
226          qr(:)=max(qfi(ij,llm:1:-1,3),0.)
227         
228          CALL KESSLER(theta(:), qv, qc, qr, rho(:),  &
229                       pk(ij,:)/cpp, dt_phys, z(:), llm, precl(ij)) 
230         
231         
232          DO l=1,llm
233           ll=llm+1-l
234           Tfi(ij,ll) = theta(l)  * ( pk(ij,l) / cpp)
235          ENDDO
236
237          qfi(ij,:,1)=qv(llm:1:-1)
238          qfi(ij,:,2)=qc(llm:1:-1)
239          qfi(ij,:,3)=qr(llm:1:-1)
240
241       ENDDO
242    ENDIF
243   
244    DO l=1,llm
245      ll=llm+1-l
246      DO ij=1,ngrid
247        CALL  tendency_terminator( lat(ij), lon(ij), qfi(ij,ll,4), qfi(ij,ll,5), dt_phys, dcl1, dcl2)
248        qfi(ij,ll,4)=qfi(ij,ll,4)+ dt_phys*dcl1
249        qfi(ij,ll,5)=qfi(ij,ll,5)+ dt_phys*dcl2
250      ENDDO
251    ENDDO
252   
253   
254    ! Mirror vertical index and compute tendencies
255    inv_dt = 1./dt_phys
256    DO l=1,llm
257       ll=llm+1-l
258       DO ij=1,ngrid
259          IF (physics_thermo==thermo_fake_moist) THEN
260            dTemp(ij,l) = inv_dt * ( Tfi(ij,ll)*(1+0.608*qfi(ij,ll,1)) - Temp(ij,l) )
261          ELSE
262            dTemp(ij,l) = inv_dt * ( Tfi(ij,ll) - Temp(ij,l) )
263          ENDIF
264         
265          du(ij,l) = inv_dt * (ufi(ij,ll) - u(ij,l))
266          dv(ij,l) = inv_dt * (vfi(ij,ll) - v(ij,l))
267          dq(ij,l,:)  = inv_dt * (qfi(ij,ll,:) - q(ij,l,:))
268       ENDDO
269    ENDDO
270
271  END SUBROUTINE compute_physics
272   
273END MODULE physics_dcmip2016_mod
274
275
Note: See TracBrowser for help on using the repository browser.