source: codes/icosagcm/trunk/src/physics_dcmip2016.f90 @ 397

Last change on this file since 397 was 397, checked in by ymipsl, 8 years ago

Prepare DCMIP2016 output by XIOS2

YM

File size: 6.6 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),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 
22CONTAINS
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    precl_packed(:)=0.
68    CALL allocate_field(f_precl, field_t,type_real)
69
70    PRINT *, 'init_physics_new', SIZE(physics_inout%Ai)
71  END SUBROUTINE init_physics
72
73  SUBROUTINE full_physics
74    USE physics_interface_mod
75    CALL compute_physics(physics_inout%ngrid, physics_inout%dt_phys, &
76         physics_inout%lon, physics_inout%lat, physics_inout%p, physics_inout%pk, physics_inout%Temp, & 
77         physics_inout%ulon, physics_inout%ulat, physics_inout%q(:,:,1:5), &
78         physics_inout%dTemp, physics_inout%dulon, physics_inout%dulat, &
79         physics_inout%dq(:,:,1:5), precl_packed)
80  END SUBROUTINE full_physics
81
82  SUBROUTINE write_physics
83    USE output_field_mod
84    USE physics_interface_mod
85    CALL unpack_field(f_precl, precl_packed)
86    CALL output_field("precl",f_precl)
87    precl_packed(:)=0.
88   
89  END SUBROUTINE write_physics
90
91  SUBROUTINE compute_physics(ngrid,dt_phys,lon, lat, p, pk, Temp,u,v,q, dTemp,du,dv,dq, precl)
92    USE icosa
93    USE dcmip2016_simple_physics_mod
94    USE dcmip2016_kessler_physic_mod
95    USE terminator
96    IMPLICIT NONE
97    INTEGER    :: ngrid
98    REAL(rstd) :: lat(ngrid)
99    REAL(rstd) :: lon(ngrid)
100    REAL(rstd) :: ps(ngrid)
101    REAL(rstd) :: precl(ngrid)
102    ! arguments with bottom-up indexing (DYNAMICO)
103    REAL(rstd) :: p(ngrid,llm+1)
104    REAL(rstd) :: pk(ngrid,llm)
105    REAL(rstd) :: Temp(ngrid,llm)
106    REAL(rstd) :: u(ngrid,llm)
107    REAL(rstd) :: v(ngrid,llm)
108    REAL(rstd) :: q(ngrid,llm,5)
109    REAL(rstd) :: dTemp(ngrid,llm)
110    REAL(rstd) :: du(ngrid,llm)
111    REAL(rstd) :: dv(ngrid,llm)
112    REAL(rstd) :: dq(ngrid,llm,5)
113    ! local arrays with top-down vertical indexing (DCMIP)
114    REAL(rstd) :: pint(ngrid,llm+1)
115    REAL(rstd) :: pmid(ngrid,llm)
116    REAL(rstd) :: pdel(ngrid,llm)
117    REAL(rstd) :: Tfi(ngrid,llm)
118    REAL(rstd) :: ufi(ngrid,llm)
119    REAL(rstd) :: vfi(ngrid,llm)
120    REAL(rstd) :: qfi(ngrid,llm,5)
121
122    REAL(rstd)  :: rho(llm), z(llm), theta(llm), qv(llm),qc(llm),qr(llm)
123    REAL(rstd)  :: lastz
124    REAL(rstd)  :: dcl1,dcl2
125     INTEGER :: l,ll,ij
126    REAL(rstd) :: dt_phys, inv_dt
127    INTEGER :: simple_physic_testcase
128   
129    ! prepare input fields and mirror vertical index     
130    ps(:) = p(:,1) ! surface pressure
131
132    DO l=1,llm+1
133      DO ij=1,ngrid
134          pint(ij,l)=p(ij,llm+2-l)
135      ENDDO
136    ENDDO
137
138    DO l=1,llm
139       ll=llm+1-l
140       DO ij=1,ngrid
141          pmid(ij,l)=0.5*(pint(ij,l)+pint(ij,l+1)) ! Pressure inside layers
142          pdel(ij,l)=pint(ij,l+1)-pint(ij,l)       ! Pressure difference between two layers
143          ufi(ij,l)=u(ij,ll)
144          vfi(ij,l)=v(ij,ll)
145          qfi(ij,l,:)=q(ij,ll,:)
146          Tfi(ij,l)=Temp(ij,ll)/(1+0.608*qfi(ij,l,1)) 
147       ENDDO
148    ENDDO
149   
150   
151    IF (testcase==moist_baroclinic .OR. testcase==cyclone  ) THEN
152      IF (testcase==moist_baroclinic) simple_physic_testcase=1
153      IF (testcase==cyclone) simple_physic_testcase=0
154      CALL simple_physics(ngrid, llm, dt_phys, lat, tfi, qfi(:,:,1) , ufi, vfi, pmid, pint, pdel, 1./pdel, ps, precl, &
155                          simple_physic_testcase, .TRUE., .FALSE.)
156    ENDIF
157
158 
159    IF (testcase==moist_baroclinic .OR. testcase==cyclone .OR. testcase==supercell ) THEN
160       DO ij=1,ngrid
161          lastz=0 
162          DO l=1,llm
163           ll=llm+1-l
164           rho(l) = pmid(ij,ll)/(287*Temp(ij,l))
165           z(l)=lastz+ (p(ij,l)-p(ij,l+1)) /g / rho(l)
166           lastz=z(l)
167!           theta(l)= Tfi(ij,ll)*(1+0.608*qfi(ij,ll,1)) / ( pk(ij,l) / cpp)
168           theta(l)= Tfi(ij,ll) / ( pk(ij,l) / cpp)
169          ENDDO
170
171          qv(:)=qfi(ij,llm:1:-1,1)
172          qc(:)=qfi(ij,llm:1:-1,2)
173          qr(:)=qfi(ij,llm:1:-1,3)
174         
175!          CALL KESSLER(theta(:), qv, qc, qr, rho(:),  &
176!                       pk(ij,:), dt_phys, z(:), llm, precl(ij))
177         
178         
179          DO l=1,llm
180           ll=llm+1-l
181!           Tfi(ij,ll) = theta(l) /(1+0.608*qfi(ij,ll,1)) * ( pk(ij,l) / cpp)
182           Tfi(ij,ll) = theta(l)  * ( pk(ij,l) / cpp)
183          ENDDO
184
185          qfi(ij,:,1)=qv(llm:1:-1)
186          qfi(ij,:,2)=qc(llm:1:-1)
187          qfi(ij,:,3)=qr(llm:1:-1)
188
189       ENDDO
190    ENDIF
191   
192    DO l=1,llm
193      ll=llm+1-l
194      DO ij=1,ngrid
195        CALL  tendency_terminator( lat(ij), lon(ij), qfi(ij,ll,4), qfi(ij,ll,5), dt_phys, dcl1, dcl2)
196        qfi(ij,ll,4)=qfi(ij,ll,4)+ dt_phys*dcl1
197        qfi(ij,ll,5)=qfi(ij,ll,5)+ dt_phys*dcl2
198      ENDDO
199    ENDDO
200   
201   
202    ! Mirror vertical index and compute tendencies
203    inv_dt = 1./dt_phys
204    DO l=1,llm
205       ll=llm+1-l
206       DO ij=1,ngrid
207          dTemp(ij,l) = inv_dt * ( Tfi(ij,ll)*(1+0.608*qfi(ij,ll,1)) - Temp(ij,l) )
208          du(ij,l) = inv_dt * (ufi(ij,ll) - u(ij,l))
209          dv(ij,l) = inv_dt * (vfi(ij,ll) - v(ij,l))
210          dq(ij,l,:)  = inv_dt * (qfi(ij,ll,:) - q(ij,l,:))
211       ENDDO
212    ENDDO
213
214  END SUBROUTINE compute_physics
215   
216END MODULE physics_dcmip2016_mod
217
218
Note: See TracBrowser for help on using the repository browser.