source: codes/icosagcm/devel/src/diagnostics/observable.f90 @ 603

Last change on this file since 603 was 603, checked in by dubos, 7 years ago

devel : output theta and momentum fluxes

File size: 13.6 KB
Line 
1MODULE observable_mod
2  USE icosa
3  USE diagflux_mod
4  USE output_field_mod
5  IMPLICIT NONE
6  PRIVATE
7
8  TYPE(t_field),POINTER, SAVE :: f_buf_i(:), &
9       f_buf_Fel(:), f_buf_uh(:), & ! horizontal velocity, different from prognostic velocity if NH
10       f_buf_ulon(:), f_buf_ulat(:)
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, f_buf_i
18
19CONTAINS
20 
21  SUBROUTINE init_observable
22    CALL allocate_field(f_buf_i,   field_t,type_real,llm,name="buffer_i")
23    CALL allocate_field(f_buf_p,   field_t,type_real,llm+1) 
24    CALL allocate_field(f_buf_ulon,field_t,type_real,llm, name="buf_ulon")
25    CALL allocate_field(f_buf_ulat,field_t,type_real,llm, name="buf_ulat")
26    CALL allocate_field(f_buf_uh,  field_u,type_real,llm, name="buf_uh")
27    CALL allocate_field(f_buf_Fel, field_u,type_real,llm+1, name="buf_F_el")
28    CALL allocate_field(f_buf_v,   field_z,type_real,llm, name="buf_v")
29    CALL allocate_field(f_buf_s,   field_t,type_real, name="buf_s")
30
31    CALL allocate_field(f_theta, field_t,type_real,llm,nqdyn,  name='theta') ! potential temperature
32    CALL allocate_field(f_pmid,  field_t,type_real,llm,  name='pmid')       ! mid layer pressure
33  END SUBROUTINE init_observable
34
35  SUBROUTINE write_output_fields_basic(init, f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
36    USE xios_mod
37    USE disvert_mod
38    USE wind_mod
39    USE omp_para
40    USE time_mod
41    USE xios_mod
42    USE earth_const
43    USE pression_mod
44    USE vertical_interp_mod
45    USE theta2theta_rhodz_mod
46    USE omega_mod
47    LOGICAL, INTENT(IN) :: init
48    INTEGER :: l
49    REAL :: scalar(1)
50    REAL :: mid_ap(llm)
51    REAL :: mid_bp(llm)
52
53    TYPE(t_field),POINTER :: f_phis(:), 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
56    CALL transfert_request(f_ps,req_i1)
57   
58    IF(init) THEN
59       IF(is_master) PRINT *, 'Creating output files ...'
60       scalar(1)=dt
61       IF (is_omp_master) CALL xios_send_field("timestep", scalar)
62       scalar(1)=preff
63       IF (is_omp_master) CALL xios_send_field("preff", scalar)
64       IF (is_omp_master) CALL xios_send_field("ap",ap)
65       IF (is_omp_master) CALL xios_send_field("bp",bp)
66       DO l=1,llm
67          mid_ap(l)=(ap(l)+ap(l+1))/2
68          mid_bp(l)=(bp(l)+bp(l+1))/2
69       ENDDO
70       IF (is_omp_master) CALL xios_send_field("mid_ap",mid_ap)
71       IF (is_omp_master) CALL xios_send_field("mid_bp",mid_bp)
72
73       CALL output_field("phis",f_phis)
74       CALL output_field("Ai",geom%Ai) 
75       IF(is_master) PRINT *, '... done creating output files. Writing initial condition ...'
76    END IF
77
78    IF(nqdyn>1) THEN
79       CALL divide_by_mass(2, f_mass, f_theta_rhodz, f_buf_i)
80       IF(init) THEN
81          CALL output_field("dyn_q_init",f_buf_i)
82       ELSE
83          CALL output_field("dyn_q",f_buf_i)
84       END IF
85    END IF
86
87    CALL divide_by_mass(1, f_mass, f_theta_rhodz, f_buf_i)
88    IF(init) THEN
89       CALL output_field("theta_init",f_buf_i)
90    ELSE
91       CALL output_field("theta",f_buf_i)
92    END IF
93
94    CALL pression_mid(f_ps, f_pmid)
95    CALL diagnose_temperature(f_pmid, f_q, f_buf_i) ! f_buf_i : IN = theta, out = T
96
97    IF(init) THEN
98       CALL output_field("temp_init",f_buf_i)
99    ELSE
100       CALL output_field("temp",f_buf_i)
101       CALL vertical_interp(f_pmid,f_buf_i,f_buf_s,85000.)
102       CALL output_field("t850",f_buf_s)
103       CALL vertical_interp(f_pmid,f_buf_i,f_buf_s,50000.)
104       CALL output_field("t500",f_buf_s)
105       CALL vertical_interp(f_pmid,f_buf_i,f_buf_s,preff)
106       CALL output_field("SST",f_buf_s)       
107    END IF
108   
109    CALL progonostic_vel_to_horiz(f_geopot, f_ps, f_mass, f_u, f_W, f_buf_uh, f_buf_i)
110    CALL transfert_request(f_buf_uh,req_e1_vect) 
111    CALL un2ulonlat(f_buf_uh, f_buf_ulon, f_buf_ulat)
112    IF(init) THEN
113       CALL output_field("uz_init",f_buf_i)
114       CALL output_field("ulon_init",f_buf_ulon)
115       CALL output_field("ulat_init",f_buf_ulat)
116       CALL output_field("p_init",f_pmid)
117       CALL output_field("ps_init",f_ps)
118       CALL output_field("mass_init",f_mass)
119       CALL output_field("geopot_init",f_geopot)
120       CALL output_field("q_init",f_q)
121       IF(is_master) PRINT *, 'Done writing initial condition ...'
122    ELSE
123       CALL output_field("uz",f_buf_i)
124       CALL output_field("ulon",f_buf_ulon)
125       CALL output_field("ulat",f_buf_ulat)
126       CALL output_field("p",f_pmid)
127       CALL output_field("ps",f_ps)
128       CALL output_field("mass",f_mass)
129       CALL output_field("geopot",f_geopot)
130       CALL output_field("q",f_q)
131
132       !       CALL output_field("exner",f_pk)
133       !       CALL output_field("pv",f_qv)
134       
135       CALL vertical_interp(f_pmid,f_buf_ulon,f_buf_s,85000.)
136       CALL output_field("u850",f_buf_s)
137       CALL vertical_interp(f_pmid,f_buf_ulon,f_buf_s,50000.)
138       CALL output_field("u500",f_buf_s)
139       
140       CALL vertical_interp(f_pmid,f_buf_ulat,f_buf_s,85000.)
141       CALL output_field("v850",f_buf_s)
142       CALL vertical_interp(f_pmid,f_buf_ulat,f_buf_s,50000.)
143       CALL output_field("v500",f_buf_s)
144
145       CALL vertical_interp(f_pmid,f_buf_i,f_buf_s,85000.)
146       CALL output_field("w850",f_buf_s)
147       CALL vertical_interp(f_pmid,f_buf_i,f_buf_s,50000.)
148       CALL output_field("w500",f_buf_s)   
149
150       CALL w_omega(f_ps, f_u, f_buf_i)
151       CALL output_field("omega",f_buf_i)
152       CALL vertical_interp(f_pmid,f_buf_i,f_buf_s,85000.)
153       CALL output_field("omega850",f_buf_s)
154       CALL vertical_interp(f_pmid,f_buf_i,f_buf_s,50000.)
155       CALL output_field("omega500",f_buf_s)
156    END IF
157
158    IF(.NOT. init) THEN
159       IF(diagflux_on) THEN
160          CALL flux_centered_lonlat(1./(itau_out*dt) , f_massfluxt, f_buf_ulon, f_buf_ulat)
161          CALL output_field("mass_t", f_masst)
162          CALL output_field("massflux_lon",f_buf_ulon)
163          CALL output_field("massflux_lat",f_buf_ulat)
164
165          CALL output_energyflux(f_ulont, f_ulonfluxt, "ulon_t", "ulonflux_lon", "ulonflux_lat")
166          CALL output_energyflux(f_thetat, f_thetafluxt, "theta_t", "thetaflux_lon", "thetaflux_lat")
167          CALL output_energyflux(f_epot, f_epotfluxt, "epot_t", "epotflux_lon", "epotflux_lat")
168          CALL output_energyflux(f_ekin, f_ekinfluxt, "ekin_t", "ekinflux_lon", "ekinflux_lat")
169          CALL output_energyflux(f_enthalpy, f_enthalpyfluxt, "enthalpy_t", "enthalpyflux_lon", "enthalpyflux_lat")
170
171          CALL qflux_centered_lonlat(1./(itau_out*dt) , f_qfluxt, f_qfluxt_lon, f_qfluxt_lat)
172          CALL output_field("qmass_t", f_qmasst)
173          CALL output_field("qflux_lon",f_qfluxt_lon)
174          CALL output_field("qflux_lat",f_qfluxt_lat)
175          CALL zero_qfluxt  ! restart accumulating fluxes
176       END IF
177    END IF
178  END SUBROUTINE write_output_fields_basic
179
180  SUBROUTINE output_energyflux(f_energy, f_flux, name_energy, name_fluxlon, name_fluxlat)
181    TYPE(t_field), POINTER :: f_energy(:), f_flux(:)
182    CHARACTER(*), INTENT(IN) :: name_energy, name_fluxlon, name_fluxlat
183    CALL transfert_request(f_flux,req_e1_vect)
184    CALL flux_centered_lonlat(1./(itau_out*dt) , f_flux, f_buf_ulon, f_buf_ulat)
185    CALL output_field(name_energy,  f_energy)
186    CALL output_field(name_fluxlon, f_buf_ulon)
187    CALL output_field(name_fluxlat, f_buf_ulat)
188  END SUBROUTINE output_energyflux
189 
190 !------------------- Conversion from prognostic to observable variables ------------------
191
192  SUBROUTINE progonostic_vel_to_horiz(f_geopot, f_ps, f_rhodz, f_u, f_W, f_uh, f_uz)
193    USE disvert_mod
194    TYPE(t_field), POINTER :: f_geopot(:), f_ps(:), f_rhodz(:), &
195         f_u(:), f_W(:), f_uz(:), &  ! IN
196         f_uh(:)                         ! OUT
197    REAL(rstd),POINTER :: geopot(:,:), ps(:), rhodz(:,:), u(:,:), W(:,:), uh(:,:), uz(:,:), F_el(:,:)
198    INTEGER :: ind
199   
200    DO ind=1,ndomain
201       IF (.NOT. assigned_domain(ind)) CYCLE
202       CALL swap_dimensions(ind)
203       CALL swap_geometry(ind)
204       geopot = f_geopot(ind)
205       rhodz = f_rhodz(ind)
206       u = f_u(ind)
207       W = f_W(ind)
208       uh  = f_uh(ind)
209       F_el  = f_buf_Fel(ind)
210       IF(caldyn_eta==eta_mass) THEN
211          ps=f_ps(ind)
212          CALL compute_rhodz(.TRUE., ps, rhodz)
213       END IF
214       uz = f_uz(ind)
215       !$OMP BARRIER
216       CALL compute_prognostic_vel_to_horiz(geopot,rhodz,u,W, F_el, uh,uz)
217       !$OMP BARRIER
218    END DO
219  END SUBROUTINE progonostic_vel_to_horiz
220 
221  SUBROUTINE compute_prognostic_vel_to_horiz(Phi, rhodz, u, W, F_el, uh, uz)
222    USE omp_para
223    REAL(rstd), INTENT(IN) :: Phi(iim*jjm,llm+1)
224    REAL(rstd), INTENT(IN) :: rhodz(iim*jjm,llm)
225    REAL(rstd), INTENT(IN) :: u(3*iim*jjm,llm)
226    REAL(rstd), INTENT(IN) :: W(iim*jjm,llm+1)
227    REAL(rstd), INTENT(OUT) :: uh(3*iim*jjm,llm)
228    REAL(rstd), INTENT(OUT) :: uz(iim*jjm,llm)
229    INTEGER :: ij,l
230    REAL(rstd) :: F_el(3*iim*jjm,llm+1)
231    REAL(rstd) :: uu_right, uu_lup, uu_ldown, W_el, DePhil
232    ! NB : u and uh are not in DEC form, they are normal components   
233    ! => we must divide by de
234    IF(hydrostatic) THEN
235       uh(:,:)=u(:,:)
236       uz(:,:)=0.
237    ELSE
238       DO l=ll_begin, ll_endp1 ! compute on l levels (interfaces)
239          DO ij=ij_begin_ext, ij_end_ext
240             ! Compute on edge 'right'
241             W_el   = .5*( W(ij,l)+W(ij+t_right,l) )
242             DePhil = ne_right*(Phi(ij+t_right,l)-Phi(ij,l))
243             F_el(ij+u_right,l)   = DePhil*W_el/de(ij+u_right)
244             ! Compute on edge 'lup'
245             W_el   = .5*( W(ij,l)+W(ij+t_lup,l) )
246             DePhil = ne_lup*(Phi(ij+t_lup,l)-Phi(ij,l))
247             F_el(ij+u_lup,l)   = DePhil*W_el/de(ij+u_lup)
248             ! Compute on edge 'ldown'
249             W_el   = .5*( W(ij,l)+W(ij+t_ldown,l) )
250             DePhil = ne_ldown*(Phi(ij+t_ldown,l)-Phi(ij,l))
251             F_el(ij+u_ldown,l) = DePhil*W_el/de(ij+u_ldown)
252          END DO
253       END DO
254
255       ! We need a barrier here because we compute F_el above and do a vertical average below
256       !$OMP BARRIER
257
258       DO l=ll_begin, ll_end ! compute on k levels (full levels)
259          DO ij=ij_begin_ext, ij_end_ext
260             ! w = vertical momentum = g^-2*dPhi/dt = uz/g
261             uz(ij,l) = (.5*g)*(W(ij,l)+W(ij,l+1))/rhodz(ij,l)
262             ! uh = u-w.grad(Phi) = u - uz.grad(z)
263             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))
264             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))
265             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))
266          END DO
267       END DO
268
269    END IF
270  END SUBROUTINE compute_prognostic_vel_to_horiz
271
272  SUBROUTINE diagnose_temperature(f_pmid,f_q,f_temp)
273    USE icosa
274    USE pression_mod
275    IMPLICIT NONE
276    TYPE(t_field), POINTER :: f_pmid(:)           ! IN
277    TYPE(t_field), POINTER :: f_q(:)            ! IN
278    TYPE(t_field), POINTER :: f_temp(:)         ! INOUT
279   
280    REAL(rstd), POINTER :: pmid(:,:)
281    REAL(rstd), POINTER :: q(:,:,:)
282    REAL(rstd), POINTER :: temp(:,:)
283    INTEGER :: ind
284   
285    DO ind=1,ndomain
286       IF (.NOT. assigned_domain(ind)) CYCLE
287       CALL swap_dimensions(ind)
288       CALL swap_geometry(ind)
289       pmid=f_pmid(ind)
290       q=f_q(ind)
291       temp=f_temp(ind)
292       CALL compute_diagnose_temp(pmid,q,temp)
293    END DO
294  END SUBROUTINE diagnose_temperature
295 
296  SUBROUTINE compute_diagnose_temp(pmid,q,temp)
297    USE omp_para
298    USE pression_mod
299    REAL(rstd),INTENT(IN)    :: pmid(iim*jjm,llm)
300    REAL(rstd),INTENT(IN)    :: q(iim*jjm,llm,nqtot)
301    REAL(rstd),INTENT(INOUT) :: temp(iim*jjm,llm)
302
303    REAL(rstd) :: Rd, p_ik, theta_ik, temp_ik, qv, chi, Rmix
304    INTEGER :: ij,l
305
306    Rd = kappa*cpp
307    DO l=ll_begin,ll_end
308       DO ij=ij_begin,ij_end
309          p_ik = pmid(ij,l)
310          theta_ik = temp(ij,l)
311          qv = q(ij,l,1) ! water vaper mixing ratio = mv/md
312          SELECT CASE(caldyn_thermo)
313          CASE(thermo_theta)
314             temp_ik = theta_ik*((p_ik/preff)**kappa)
315          CASE(thermo_entropy)
316             temp_ik = Treff*exp((theta_ik + Rd*log(p_ik/preff))/cpp)
317          CASE(thermo_moist)
318             Rmix = Rd+qv*Rv
319             chi = ( theta_ik + Rmix*log(p_ik/preff) ) / (cpp + qv*cppv) ! log(T/Treff)
320             temp_ik = Treff*exp(chi)
321          END SELECT
322          IF(physics_thermo==thermo_fake_moist) temp_ik=temp_ik/(1+0.608*qv)
323          temp(ij,l)=temp_ik
324       END DO
325    END DO
326  END SUBROUTINE compute_diagnose_temp
327
328  SUBROUTINE Tv2T(f_Tv, f_q, f_T)
329    TYPE(t_field), POINTER :: f_TV(:)
330    TYPE(t_field), POINTER :: f_q(:)
331    TYPE(t_field), POINTER :: f_T(:)
332   
333    REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:)
334    INTEGER :: ind
335   
336    DO ind=1,ndomain
337       IF (.NOT. assigned_domain(ind)) CYCLE
338       CALL swap_dimensions(ind)
339       CALL swap_geometry(ind)
340       Tv=f_Tv(ind)
341       T=f_T(ind)
342       SELECT CASE(physics_thermo)
343       CASE(thermo_dry)
344          T=Tv
345       CASE(thermo_fake_moist)
346          q=f_q(ind)
347          T=Tv/(1+0.608*q(:,:,1))
348       END SELECT
349    END DO
350  END SUBROUTINE Tv2T
351
352  SUBROUTINE divide_by_mass(iq, f_mass, f_theta_rhodz, f_theta)
353    INTEGER, INTENT(IN) :: iq
354    TYPE(t_field), POINTER :: f_mass(:), f_theta_rhodz(:), f_theta(:)
355    REAL(rstd), POINTER :: mass(:,:), theta_rhodz(:,:,:), theta(:,:)
356    INTEGER :: ind
357    DO ind=1,ndomain
358       IF (.NOT. assigned_domain(ind)) CYCLE
359       CALL swap_dimensions(ind)
360       CALL swap_geometry(ind)
361       mass=f_mass(ind)
362       theta_rhodz=f_theta_rhodz(ind)
363       theta=f_theta(ind)
364       theta(:,:) = theta_rhodz(:,:,iq) / mass(:,:)
365    END DO
366  END SUBROUTINE divide_by_mass
367   
368END MODULE observable_mod
Note: See TracBrowser for help on using the repository browser.