source: codes/icosagcm/trunk/src/dcmip/physics_dcmip2016.f90

Last change on this file was 899, checked in by adurocher, 5 years ago

trunk : Fixed GCC warnings

Fixed iso c bindings
fixed warnings with -Wall -Wno-aliasing -Wno-unused -Wno-unused-dummy-argument -Wno-maybe-uninitialized -Wno-tabs warnings
Removed all unused variables (-Wunused-variable)
vector%dot_product is now dot_product_3d to avoid compilation warning "dot_product shadows intrinsic" with GCC

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