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

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

trunk : backported r600-603 from devel

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