source: codes/icosagcm/trunk/src/observable.f90 @ 401

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

Prepare DCMIP2016 output by XIOS2

YM

File size: 12.3 KB
Line 
1MODULE observable_mod
2  USE icosa
3  IMPLICIT NONE
4  PRIVATE
5
6  TYPE(t_field),POINTER, SAVE :: f_buf_i(:), &
7       f_buf_uh(:), & ! horizontal velocity, different from prognostic velocity if NH
8       f_buf_ulon(:), f_buf_ulat(:), &
9       f_buf_u3d(:) ! unused, remove ?
10  TYPE(t_field),POINTER, SAVE :: f_buf1_i(:), f_buf2_i(:)
11  TYPE(t_field),POINTER, SAVE :: f_buf_v(:), f_buf_s(:), f_buf_p(:)
12  TYPE(t_field),POINTER, SAVE :: f_pmid(:)
13
14! temporary shared variable for caldyn
15  TYPE(t_field),POINTER, SAVE :: f_theta(:)
16
17  PUBLIC init_observable, write_output_fields_basic, f_theta
18  LOGICAL,SAVE :: first_output=.TRUE.
19!$OMP THREADPRIVATE(first_output)
20 
21CONTAINS
22 
23  SUBROUTINE init_observable
24    CALL allocate_field(f_buf_i,   field_t,type_real,llm,name="buffer_i")
25    CALL allocate_field(f_buf1_i,   field_t,type_real,llm,name="buffer1_i")
26    CALL allocate_field(f_buf2_i,   field_t,type_real,llm,name="buffer2_i")
27    CALL allocate_field(f_buf_p,   field_t,type_real,llm+1) 
28    CALL allocate_field(f_buf_u3d, field_t,type_real,3,llm)  ! 3D vel at cell centers
29    CALL allocate_field(f_buf_ulon,field_t,type_real,llm, name="buf_ulon")
30    CALL allocate_field(f_buf_ulat,field_t,type_real,llm, name="buf_ulat")
31    CALL allocate_field(f_buf_uh,  field_u,type_real,llm, name="buf_uh")
32    CALL allocate_field(f_buf_v,   field_z,type_real,llm, name="buf_v")
33    CALL allocate_field(f_buf_s,   field_t,type_real, name="buf_s")
34
35    CALL allocate_field(f_theta, field_t,type_real,llm,  name='theta')   ! potential temperature
36    CALL allocate_field(f_pmid,  field_t,type_real,llm,  name='pmid')   ! mid layer pressure
37  END SUBROUTINE init_observable
38 
39  SUBROUTINE write_output_fields_basic(f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
40    USE wind_mod
41    USE output_field_mod
42    USE omp_para
43    USE time_mod
44    USE xios
45    USE disvert_mod
46    USE earth_const
47    USE pression_mod
48    USE vertical_interp_mod
49    USE theta2theta_rhodz_mod
50    USE wind_mod
51    USE omega_mod
52   
53    TYPE(t_field),POINTER :: f_ps(:), f_mass(:), f_geopot(:), f_theta_rhodz(:), f_u(:), f_W(:), f_q(:)
54!    IF (is_master) PRINT *,'CALL write_output_fields_basic'
55    REAL :: scalar(1)
56    REAL :: mid_ap(llm)
57    REAL :: mid_bp(llm)
58    INTEGER :: l
59
60    IF (first_output) THEN
61      scalar(1)=dt
62      CALL xios_send_field("timestep", scalar)
63      scalar(1)=preff
64      CALL xios_send_field("preff", scalar)
65      CALL xios_send_field("ap",ap)
66      CALL xios_send_field("bp",bp)
67      DO l=1,llm
68        mid_ap(l)=(ap(l)+ap(l+1))/2
69        mid_bp(l)=(bp(l)+bp(l+1))/2
70      ENDDO
71      CALL xios_send_field("mid_ap",mid_ap)
72      CALL xios_send_field("mid_bp",mid_bp)
73       
74      first_output=.FALSE.
75    ENDIF
76   
77    CALL progonostic_vel_to_horiz(f_geopot, f_ps, f_mass, f_u, f_W, f_buf_uh, f_buf_i)
78    CALL transfert_request(f_buf_uh,req_e1_vect) 
79    CALL output_field("uz",f_buf_i)
80    CALL vertical_interp(f_ps,f_buf_i,f_buf_s,85000.)
81    CALL output_field("w850",f_buf_s)
82    CALL vertical_interp(f_ps,f_buf_i,f_buf_s,50000.)
83    CALL output_field("w500",f_buf_s)
84   
85   
86    CALL un2ulonlat(f_buf_uh, f_buf_ulon, f_buf_ulat)
87    CALL output_field("ulon",f_buf_ulon)
88    CALL output_field("ulat",f_buf_ulat)
89    CALL output_field("ps",f_ps)
90    CALL output_field("Ai",geom%Ai)
91
92    !       CALL output_field("dps",f_dps)
93    CALL output_field("mass",f_mass)
94    CALL output_field("geopot",f_geopot)
95    !       CALL output_field("dmass",f_dmass)
96    !       CALL output_field("vort",f_qv)
97   
98   
99    CALL theta_rhodz2temperature(f_ps,f_theta_rhodz,f_buf_i) ; 
100    CALL Tv2T(f_buf_i,f_q,f_buf1_i) 
101    CALL output_field("temp",f_buf_i)
102    CALL vertical_interp(f_ps,f_buf_i,f_buf_s,85000.)
103    CALL output_field("t850",f_buf_s)
104    CALL vertical_interp(f_ps,f_buf_i,f_buf_s,50000.)
105    CALL output_field("t500",f_buf_s)
106    CALL vertical_interp(f_ps,f_buf_i,f_buf_s,preff)
107    CALL output_field("SST",f_buf_s)
108
109
110    CALL extract_slice(f_theta_rhodz, f_buf_i,1)
111    CALL output_field("theta",f_buf_i)
112           
113    !       CALL output_field("exner",f_pk)
114    !       CALL output_field("pv",f_qv)
115    CALL output_field("q",f_q)
116    CALL pression_mid(f_ps, f_pmid)
117    CALL output_field("p",f_pmid)
118
119    CALL vertical_interp(f_ps,f_buf_ulon,f_buf_s,85000.)
120    CALL output_field("u850",f_buf_s)
121    CALL vertical_interp(f_ps,f_buf_ulon,f_buf_s,50000.)
122    CALL output_field("u500",f_buf_s)
123
124    CALL vertical_interp(f_ps,f_buf_ulat,f_buf_s,85000.)
125    CALL output_field("v850",f_buf_s)
126    CALL vertical_interp(f_ps,f_buf_ulat,f_buf_s,50000.)
127    CALL output_field("v500",f_buf_s)
128   
129    CALL w_omega(f_ps, f_u, f_buf_i)
130    CALL output_field("omega",f_buf_i)
131    CALL vertical_interp(f_ps,f_buf_i,f_buf_s,85000.)
132    CALL output_field("omega850",f_buf_s)
133    CALL vertical_interp(f_ps,f_buf_i,f_buf_s,50000.)
134    CALL output_field("omega500",f_buf_s)
135
136
137   
138  END SUBROUTINE write_output_fields_basic
139 
140  SUBROUTINE write_output_fields(f_ps, f_phis, f_dps, f_u, f_theta_rhodz, f_q, &
141       f_buf_i, f_buf_v, f_buf_i3, f_buf1_i, f_buf2_i, f_buf_s, f_buf_p)
142    USE vorticity_mod
143    USE theta2theta_rhodz_mod
144    USE pression_mod
145    USE omega_mod
146    USE write_field_mod
147    USE vertical_interp_mod
148    USE wind_mod
149    TYPE(t_field),POINTER :: f_ps(:), f_phis(:), f_u(:), f_theta_rhodz(:), f_q(:), f_dps(:), &
150         f_buf_i(:), f_buf_v(:), f_buf_i3(:), f_buf1_i(:), f_buf2_i(:), f_buf_s(:), f_buf_p(:)
151   
152    REAL(rstd) :: out_pression_level
153    CHARACTER(LEN=255) :: str_pression
154    CHARACTER(LEN=255) :: physics_type
155   
156    out_pression_level=0.
157    CALL getin("out_pression_level",out_pression_level) 
158    WRITE(str_pression,*) INT(out_pression_level/100)
159    str_pression=ADJUSTL(str_pression)
160   
161    CALL writefield("ps",f_ps)
162    CALL writefield("dps",f_dps)
163    CALL writefield("phis",f_phis)
164    CALL vorticity(f_u,f_buf_v)
165    CALL writefield("vort",f_buf_v)
166
167    CALL w_omega(f_ps, f_u, f_buf_i)
168    CALL writefield('omega', f_buf_i)
169    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
170      CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
171      CALL writefield("omega"//TRIM(str_pression),f_buf_s)
172    ENDIF
173   
174    ! Temperature
175!    CALL theta_rhodz2temperature(f_ps,f_theta_rhodz,f_buf_i) ; ! FIXME
176     
177    CALL getin('physics',physics_type)
178    IF (TRIM(physics_type)=='dcmip') THEN
179      CALL Tv2T(f_buf_i,f_q,f_buf1_i) 
180      CALL writefield("T",f_buf1_i)
181      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
182        CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
183        CALL writefield("T"//TRIM(str_pression),f_buf_s)
184      ENDIF
185    ELSE
186      CALL writefield("T",f_buf_i)
187      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
188        CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
189        CALL writefield("T"//TRIM(str_pression),f_buf_s)
190      ENDIF
191    ENDIF
192   
193    ! velocity components
194    CALL un2ulonlat(f_u, f_buf1_i, f_buf2_i)
195    CALL writefield("ulon",f_buf1_i)
196    CALL writefield("ulat",f_buf2_i)
197
198    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
199      CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
200      CALL writefield("ulon"//TRIM(str_pression),f_buf_s)
201      CALL vertical_interp(f_ps,f_buf2_i,f_buf_s,out_pression_level)
202      CALL writefield("ulat"//TRIM(str_pression),f_buf_s)
203    ENDIF
204   
205    ! geopotential ! FIXME
206    CALL thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_buf_s,f_buf_p,f_buf1_i,f_buf2_i,f_buf_i)
207    CALL writefield("p",f_buf_p)
208!    CALL writefield("phi",f_geopot)   ! geopotential
209    CALL writefield("theta",f_buf1_i) ! potential temperature
210    CALL writefield("pk",f_buf2_i)    ! Exner pressure
211 
212  END SUBROUTINE write_output_fields
213
214!------------------- Conversion from prognostic to observable variables ------------------
215
216  SUBROUTINE progonostic_vel_to_horiz(f_geopot, f_ps, f_rhodz, f_u, f_W, f_uh, f_uz)
217    USE disvert_mod
218    TYPE(t_field), POINTER :: f_geopot(:), f_ps(:), f_rhodz(:), &
219         f_u(:), f_W(:), f_uz(:), &  ! IN
220         f_uh(:)                         ! OUT
221    REAL(rstd),POINTER :: geopot(:,:), ps(:), rhodz(:,:), u(:,:), W(:,:), uh(:,:), uz(:,:)
222    INTEGER :: ind
223   
224    DO ind=1,ndomain
225       IF (.NOT. assigned_domain(ind)) CYCLE
226       CALL swap_dimensions(ind)
227       CALL swap_geometry(ind)
228       geopot = f_geopot(ind)
229       rhodz = f_rhodz(ind)
230       u = f_u(ind)
231       W = f_W(ind)
232       uh  = f_uh(ind)
233       IF(caldyn_eta==eta_mass) THEN
234          ps=f_ps(ind)
235          CALL compute_rhodz(.TRUE., ps, rhodz)
236       END IF
237       uz = f_uz(ind)
238       CALL compute_prognostic_vel_to_horiz(geopot,rhodz,u,W,uh,uz)
239    END DO
240  END SUBROUTINE progonostic_vel_to_horiz
241 
242  SUBROUTINE compute_prognostic_vel_to_horiz(Phi, rhodz, u, W, uh, uz)
243    USE omp_para
244    REAL(rstd), INTENT(IN) :: Phi(iim*jjm,llm+1)
245    REAL(rstd), INTENT(IN) :: rhodz(iim*jjm,llm)
246    REAL(rstd), INTENT(IN) :: u(3*iim*jjm,llm)
247    REAL(rstd), INTENT(IN) :: W(iim*jjm,llm+1)
248    REAL(rstd), INTENT(OUT) :: uh(3*iim*jjm,llm)
249    REAL(rstd), INTENT(OUT) :: uz(iim*jjm,llm)
250    INTEGER :: ij,l
251    REAL(rstd) :: F_el(3*iim*jjm,llm+1)
252    REAL(rstd) :: uu_right, uu_lup, uu_ldown, W_el, DePhil
253    ! NB : u and uh are not in DEC form, they are normal components   
254    ! => we must divide by de
255    IF(hydrostatic) THEN
256       uh(:,:)=u(:,:)
257       uz(:,:)=0.
258    ELSE
259       DO l=ll_begin, ll_endp1 ! compute on l levels (interfaces)
260          DO ij=ij_begin_ext, ij_end_ext
261             ! Compute on edge 'right'
262             W_el   = .5*( W(ij,l)+W(ij+t_right,l) )
263             DePhil = ne_right*(Phi(ij+t_right,l)-Phi(ij,l))
264             F_el(ij+u_right,l)   = DePhil*W_el/de(ij+u_right)
265             ! Compute on edge 'lup'
266             W_el   = .5*( W(ij,l)+W(ij+t_lup,l) )
267             DePhil = ne_lup*(Phi(ij+t_lup,l)-Phi(ij,l))
268             F_el(ij+u_lup,l)   = DePhil*W_el/de(ij+u_lup)
269             ! Compute on edge 'ldown'
270             W_el   = .5*( W(ij,l)+W(ij+t_ldown,l) )
271             DePhil = ne_ldown*(Phi(ij+t_ldown,l)-Phi(ij,l))
272             F_el(ij+u_ldown,l) = DePhil*W_el/de(ij+u_ldown)
273          END DO
274       END DO
275
276       DO l=ll_begin, ll_end ! compute on k levels (full levels)
277          DO ij=ij_begin_ext, ij_end_ext
278             ! w = vertical momentum = g^-2*dPhi/dt = uz/g
279             uz(ij,l) = (.5*g)*(W(ij,l)+W(ij,l+1))/rhodz(ij,l)
280             ! uh = u-w.grad(Phi) = u - uz.grad(z)
281             uh(ij+u_right,l) = u(ij+u_right,l) - (F_el(ij+u_right,l)+F_el(ij+u_right,l+1)) / (rhodz(ij,l)+rhodz(ij+t_right,l))
282             uh(ij+u_lup,l)   = u(ij+u_lup,l)   - (F_el(ij+u_lup,l)+F_el(ij+u_lup,l+1))     / (rhodz(ij,l)+rhodz(ij+t_lup,l))
283             uh(ij+u_ldown,l) = u(ij+u_ldown,l) - (F_el(ij+u_ldown,l)+F_el(ij+u_ldown,l+1)) / (rhodz(ij,l)+rhodz(ij+t_ldown,l))
284          END DO
285       END DO
286
287    END IF
288  END SUBROUTINE compute_prognostic_vel_to_horiz
289
290  SUBROUTINE thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_pks,f_p,f_theta,f_pk,f_phi) 
291    USE field_mod
292    USE pression_mod
293    USE exner_mod
294    USE geopotential_mod
295    USE theta2theta_rhodz_mod
296    TYPE(t_field), POINTER :: f_ps(:), f_phis(:), f_theta_rhodz(:), &  ! IN
297         f_pks(:), f_p(:), f_theta(:), f_pk(:), f_phi(:)               ! OUT
298    REAL(rstd),POINTER :: pk(:,:), p(:,:), theta(:,:), theta_rhodz(:,:,:), &
299         phi(:,:), phis(:), ps(:), pks(:)
300    INTEGER :: ind
301   
302    DO ind=1,ndomain
303       IF (.NOT. assigned_domain(ind)) CYCLE
304       CALL swap_dimensions(ind)
305       CALL swap_geometry(ind)
306       ps = f_ps(ind)
307       p  = f_p(ind)
308!$OMP BARRIER
309       CALL compute_pression(ps,p,0)
310       pk = f_pk(ind)
311       pks = f_pks(ind)
312!$OMP BARRIER
313       CALL compute_exner(ps,p,pks,pk,0)
314!$OMP BARRIER
315       theta_rhodz = f_theta_rhodz(ind)
316       theta = f_theta(ind)
317       CALL compute_theta_rhodz2theta(ps, theta_rhodz(:,:,1),theta,0)
318       phis = f_phis(ind)
319       phi = f_phi(ind)
320       CALL compute_geopotential(phis,pks,pk,theta,phi,0)
321    END DO
322
323  END SUBROUTINE thetarhodz2geopot
324 
325  SUBROUTINE Tv2T(f_Tv, f_q, f_T)
326    TYPE(t_field), POINTER :: f_TV(:)
327    TYPE(t_field), POINTER :: f_q(:)
328    TYPE(t_field), POINTER :: f_T(:)
329   
330    REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:)
331    INTEGER :: ind
332   
333    DO ind=1,ndomain
334       IF (.NOT. assigned_domain(ind)) CYCLE
335       CALL swap_dimensions(ind)
336       CALL swap_geometry(ind)
337       Tv=f_Tv(ind)
338       q=f_q(ind)
339       T=f_T(ind)
340       T=Tv/(1+0.608*q(:,:,1))
341    END DO
342   
343  END SUBROUTINE Tv2T
344 
345END MODULE observable_mod
Note: See TracBrowser for help on using the repository browser.