source: codes/icosagcm/trunk/src/diagnostics/observable.f90

Last change on this file was 1024, checked in by millour, 4 years ago

Bug fix: qv computation is only possible if using at least 1 tracer.
EM

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