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

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

devel : more cleanup and reorganization in dynamics/

File size: 9.2 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    ngrid = physics_inout%ngrid
74    ! Input
75    ALLOCATE(physics_inout%Ai(ngrid))
76    ALLOCATE(physics_inout%lon(ngrid))
77    ALLOCATE(physics_inout%lat(ngrid))
78    ALLOCATE(physics_inout%phis(ngrid))
79    ALLOCATE(physics_inout%p(ngrid,llm+1))
80    ALLOCATE(physics_inout%pk(ngrid,llm))
81    ALLOCATE(physics_inout%Temp(ngrid,llm))
82    ALLOCATE(physics_inout%ulon(ngrid,llm))
83    ALLOCATE(physics_inout%ulat(ngrid,llm))
84    ALLOCATE(physics_inout%q(ngrid,llm,nqtot))
85    ! Output (tendencies)
86    ALLOCATE(physics_inout%dTemp(ngrid,llm))
87    ALLOCATE(physics_inout%dulon(ngrid,llm))
88    ALLOCATE(physics_inout%dulat(ngrid,llm))
89    ALLOCATE(physics_inout%dq(ngrid,llm,nqtot))
90    ! Physics-specific data
91    ALLOCATE(precl_packed(ngrid))
92    CALL allocate_field(f_precl, field_t,type_real)
93    CALL allocate_field(f_Q1, field_t,type_real,llm)
94    CALL allocate_field(f_Q2, field_t,type_real,llm)
95    CALL allocate_field(f_PS, field_t,type_real)
96    CALL allocate_field(f_rhodz, field_t,type_real,llm)
97    CALL allocate_field(f_Q1_col_int, field_t,type_real)
98    CALL allocate_field(f_Q2_col_int, field_t,type_real)
99
100    PRINT *, 'init_physics_new', SIZE(physics_inout%Ai)
101  END SUBROUTINE init_physics
102
103  SUBROUTINE full_physics
104    USE physics_interface_mod
105    CALL compute_physics(physics_inout%ngrid, physics_inout%dt_phys, &
106         physics_inout%lon, physics_inout%lat, physics_inout%p, physics_inout%pk, physics_inout%Temp, & 
107         physics_inout%ulon, physics_inout%ulat, physics_inout%q(:,:,1:5), &
108         physics_inout%dTemp, physics_inout%dulon, physics_inout%dulat, &
109         physics_inout%dq(:,:,1:5), precl_packed)
110  END SUBROUTINE full_physics
111
112  SUBROUTINE write_physics
113    USE output_field_mod
114    USE physics_interface_mod
115    use disvert_mod
116    REAL(rstd), POINTER :: Q1(:,:)
117    REAL(rstd), POINTER :: Q2(:,:)
118    REAL(rstd), POINTER :: PS(:)
119    REAL(rstd), POINTER :: rhodz(:,:)
120    REAL(rstd), POINTER :: Q1_col_int(:)
121    REAL(rstd), POINTER :: Q2_col_int(:)
122   
123   
124    CALL unpack_field(f_precl, precl_packed)
125    CALL output_field("precl",f_precl)
126   
127    CALL unpack_field(f_Q1,physics_inout%q(:,:,4))
128    CALL unpack_field(f_Q2,physics_inout%q(:,:,5))
129    CALL unpack_field(f_ps,physics_inout%p(:,1))
130   
131    DO ind=1,ndomain
132      IF (.NOT. assigned_domain(ind)) CYCLE
133        Q1=f_Q1(ind)
134        Q2=f_Q2(ind)
135        Q1_col_int=f_Q1_col_int(ind)
136        Q2_col_int=f_Q2_col_int(ind)
137        PS=f_PS(ind)
138        rhodz=f_rhodz(ind)
139        DO l=1,llm
140          rhodz(:,l)= ( ap(l)-ap(l+1) + (bp(l)-bp(l+1))*ps(:) )/g   
141        ENDDO
142        Q1_col_int=SUM(Q1*rhodz,2)/SUM(rhodz,2)
143        Q2_col_int=SUM(Q2*rhodz,2)/SUM(rhodz,2)
144    ENDDO
145     
146    CALL output_field("Q1_col_int",f_Q1_col_int)
147    CALL output_field("Q2_col_int",f_Q2_col_int)
148   
149  END SUBROUTINE write_physics
150
151  SUBROUTINE compute_physics(ngrid,dt_phys,lon, lat, p, pk, Temp,u,v,q, dTemp,du,dv,dq, precl)
152    USE icosa
153    USE dcmip2016_simple_physics_mod
154    USE dcmip2016_kessler_physic_mod
155    USE earth_const
156    USE terminator
157    IMPLICIT NONE
158    INTEGER    :: ngrid
159    REAL(rstd) :: lat(ngrid)
160    REAL(rstd) :: lon(ngrid)
161    REAL(rstd) :: ps(ngrid)
162    REAL(rstd) :: precl(ngrid)
163    ! arguments with bottom-up indexing (DYNAMICO)
164    REAL(rstd) :: p(ngrid,llm+1)
165    REAL(rstd) :: pk(ngrid,llm)
166    REAL(rstd) :: Temp(ngrid,llm)
167    REAL(rstd) :: u(ngrid,llm)
168    REAL(rstd) :: v(ngrid,llm)
169    REAL(rstd) :: q(ngrid,llm,5)
170    REAL(rstd) :: dTemp(ngrid,llm)
171    REAL(rstd) :: du(ngrid,llm)
172    REAL(rstd) :: dv(ngrid,llm)
173    REAL(rstd) :: dq(ngrid,llm,5)
174    ! local arrays with top-down vertical indexing (DCMIP)
175    REAL(rstd) :: pint(ngrid,llm+1)
176    REAL(rstd) :: pmid(ngrid,llm)
177    REAL(rstd) :: pdel(ngrid,llm)
178    REAL(rstd) :: Tfi(ngrid,llm)
179    REAL(rstd) :: ufi(ngrid,llm)
180    REAL(rstd) :: vfi(ngrid,llm)
181    REAL(rstd) :: qfi(ngrid,llm,5)
182
183    REAL(rstd)  :: rho(llm), z(llm), theta(llm), qv(llm),qc(llm),qr(llm)
184    REAL(rstd)  :: lastz
185    REAL(rstd)  :: dcl1,dcl2
186     INTEGER :: l,ll,ij
187    REAL(rstd) :: dt_phys, inv_dt
188    INTEGER :: simple_physic_testcase
189   
190    ! prepare input fields and mirror vertical index     
191    ps(:) = p(:,1) ! surface pressure
192
193    DO l=1,llm+1
194      DO ij=1,ngrid
195          pint(ij,l)=p(ij,llm+2-l)
196      ENDDO
197    ENDDO
198
199    DO l=1,llm
200       ll=llm+1-l
201       DO ij=1,ngrid
202          pmid(ij,l)=0.5*(pint(ij,l)+pint(ij,l+1)) ! Pressure inside layers
203          pdel(ij,l)=pint(ij,l+1)-pint(ij,l)       ! Pressure difference between two layers
204          ufi(ij,l)=u(ij,ll)
205          vfi(ij,l)=v(ij,ll)
206          qfi(ij,l,:)=q(ij,ll,:)
207          IF (physics_thermo==thermo_fake_moist) THEN
208            Tfi(ij,l)=Temp(ij,ll)/(1+0.608*qfi(ij,l,1)) 
209          ELSE
210            Tfi(ij,l)=Temp(ij,ll)
211          ENDIF
212       ENDDO
213    ENDDO
214   
215    precl=0.
216    IF (testcase==moist_baroclinic_full .OR. testcase==cyclone  ) THEN
217      IF (testcase==moist_baroclinic_full) simple_physic_testcase=1
218      IF (testcase==cyclone) simple_physic_testcase=0
219      CALL simple_physics(ngrid, llm, dt_phys, lat, tfi, qfi(:,:,1) , ufi, vfi, pmid, pint, pdel, 1./pdel, ps, precl, &
220                          simple_physic_testcase, .FALSE., PBL)
221    ENDIF
222
223 
224    IF (testcase==moist_baroclinic_full .OR. testcase==moist_baroclinic_kessler .OR. testcase==cyclone .OR. testcase==supercell ) THEN
225       DO ij=1,ngrid
226          lastz=0 
227          DO l=1,llm
228           ll=llm+1-l
229           rho(l) = pmid(ij,ll)/(287*Tfi(ij,ll)*(1+0.608*qfi(ij,ll,1)))
230           z(l)=lastz
231           lastz=lastz+ (p(ij,l)-p(ij,l+1)) /g / rho(l)
232           theta(l)= Tfi(ij,ll) / ( pk(ij,l) / cpp)
233          ENDDO
234         
235          DO l=1,llm-1
236           z(l)= 0.5*(z(l)+z(l+1))
237          ENDDO
238          z(llm)=z(llm)+(z(llm)-z(llm-1))
239         
240          qv(:)=max(qfi(ij,llm:1:-1,1),0.)
241          qc(:)=max(qfi(ij,llm:1:-1,2),0.)
242          qr(:)=max(qfi(ij,llm:1:-1,3),0.)
243         
244          CALL KESSLER(theta(:), qv, qc, qr, rho(:),  &
245                       pk(ij,:)/cpp, dt_phys, z(:), llm, precl(ij)) 
246         
247         
248          DO l=1,llm
249           ll=llm+1-l
250           Tfi(ij,ll) = theta(l)  * ( pk(ij,l) / cpp)
251          ENDDO
252
253          qfi(ij,:,1)=qv(llm:1:-1)
254          qfi(ij,:,2)=qc(llm:1:-1)
255          qfi(ij,:,3)=qr(llm:1:-1)
256
257       ENDDO
258    ENDIF
259   
260    DO l=1,llm
261      ll=llm+1-l
262      DO ij=1,ngrid
263        CALL  tendency_terminator( lat(ij), lon(ij), qfi(ij,ll,4), qfi(ij,ll,5), dt_phys, dcl1, dcl2)
264        qfi(ij,ll,4)=qfi(ij,ll,4)+ dt_phys*dcl1
265        qfi(ij,ll,5)=qfi(ij,ll,5)+ dt_phys*dcl2
266      ENDDO
267    ENDDO
268   
269   
270    ! Mirror vertical index and compute tendencies
271    inv_dt = 1./dt_phys
272    DO l=1,llm
273       ll=llm+1-l
274       DO ij=1,ngrid
275          IF (physics_thermo==thermo_fake_moist) THEN
276            dTemp(ij,l) = inv_dt * ( Tfi(ij,ll)*(1+0.608*qfi(ij,ll,1)) - Temp(ij,l) )
277          ELSE
278            dTemp(ij,l) = inv_dt * ( Tfi(ij,ll) - Temp(ij,l) )
279          ENDIF
280         
281          du(ij,l) = inv_dt * (ufi(ij,ll) - u(ij,l))
282          dv(ij,l) = inv_dt * (vfi(ij,ll) - v(ij,l))
283          dq(ij,l,:)  = inv_dt * (qfi(ij,ll,:) - q(ij,l,:))
284       ENDDO
285    ENDDO
286
287  END SUBROUTINE compute_physics
288   
289END MODULE physics_dcmip2016_mod
290
291
Note: See TracBrowser for help on using the repository browser.