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

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

trunk : backported commits r582-r598 (transport diagnostics)

File size: 13.7 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_Fel(:), f_buf_uh(:), & ! horizontal velocity, different from prognostic velocity if NH
8       f_buf_ulon(:), f_buf_ulat(:)
9  TYPE(t_field),POINTER, SAVE :: f_buf1_i(:), f_buf2_i(:)
10  TYPE(t_field),POINTER, SAVE :: f_buf_v(:), f_buf_s(:), f_buf_p(:)
11  TYPE(t_field),POINTER, SAVE :: f_pmid(:)
12
13! temporary shared variable for caldyn
14  TYPE(t_field),POINTER, SAVE :: f_theta(:)
15
16  PUBLIC init_observable, write_output_fields_basic, f_theta
17
18CONTAINS
19 
20  SUBROUTINE init_observable
21    CALL allocate_field(f_buf_i,   field_t,type_real,llm,name="buffer_i")
22    CALL allocate_field(f_buf1_i,   field_t,type_real,llm,name="buffer1_i")
23    CALL allocate_field(f_buf2_i,   field_t,type_real,llm,name="buffer2_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') ! 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 output_field_mod
41    USE omp_para
42    USE time_mod
43    USE xios_mod
44    USE earth_const
45    USE pression_mod
46    USE vertical_interp_mod
47    USE theta2theta_rhodz_mod
48    USE omega_mod
49    USE diagflux_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 transfert_request(f_epotfluxt,req_e1_vect) 
169          CALL flux_centered_lonlat(1./(itau_out*dt) , f_epotfluxt, f_buf_ulon, f_buf_ulat)
170          CALL output_field("epot_t", f_epot)
171          CALL output_field("epotflux_lon",f_buf_ulon)
172          CALL output_field("epotflux_lat",f_buf_ulat)
173
174          CALL transfert_request(f_ekinfluxt,req_e1_vect) 
175          CALL flux_centered_lonlat(1./(itau_out*dt) , f_ekinfluxt, f_buf_ulon, f_buf_ulat)
176          CALL output_field("ekin_t", f_ekin)
177          CALL output_field("ekinflux_lon",f_buf_ulon)
178          CALL output_field("ekinflux_lat",f_buf_ulat)
179
180          CALL transfert_request(f_enthalpyfluxt,req_e1_vect) 
181          CALL flux_centered_lonlat(1./(itau_out*dt) , f_enthalpyfluxt, f_buf_ulon, f_buf_ulat)
182          CALL output_field("enthalpy_t", f_enthalpy)
183          CALL output_field("enthalpyflux_lon",f_buf_ulon)
184          CALL output_field("enthalpyflux_lat",f_buf_ulat)
185
186          CALL qflux_centered_lonlat(1./(itau_out*dt) , f_qfluxt, f_qfluxt_lon, f_qfluxt_lat)
187          CALL output_field("qmass_t", f_qmasst)
188          CALL output_field("qflux_lon",f_qfluxt_lon)
189          CALL output_field("qflux_lat",f_qfluxt_lat)
190          CALL zero_qfluxt  ! restart accumulating fluxes
191       END IF
192    END IF
193  END SUBROUTINE write_output_fields_basic
194 
195 !------------------- Conversion from prognostic to observable variables ------------------
196
197  SUBROUTINE progonostic_vel_to_horiz(f_geopot, f_ps, f_rhodz, f_u, f_W, f_uh, f_uz)
198    USE disvert_mod
199    TYPE(t_field), POINTER :: f_geopot(:), f_ps(:), f_rhodz(:), &
200         f_u(:), f_W(:), f_uz(:), &  ! IN
201         f_uh(:)                         ! OUT
202    REAL(rstd),POINTER :: geopot(:,:), ps(:), rhodz(:,:), u(:,:), W(:,:), uh(:,:), uz(:,:), F_el(:,:)
203    INTEGER :: ind
204   
205    DO ind=1,ndomain
206       IF (.NOT. assigned_domain(ind)) CYCLE
207       CALL swap_dimensions(ind)
208       CALL swap_geometry(ind)
209       geopot = f_geopot(ind)
210       rhodz = f_rhodz(ind)
211       u = f_u(ind)
212       W = f_W(ind)
213       uh  = f_uh(ind)
214       F_el  = f_buf_Fel(ind)
215       IF(caldyn_eta==eta_mass) THEN
216          ps=f_ps(ind)
217          CALL compute_rhodz(.TRUE., ps, rhodz)
218       END IF
219       uz = f_uz(ind)
220       !$OMP BARRIER
221       CALL compute_prognostic_vel_to_horiz(geopot,rhodz,u,W, F_el, uh,uz)
222       !$OMP BARRIER
223    END DO
224  END SUBROUTINE progonostic_vel_to_horiz
225 
226  SUBROUTINE compute_prognostic_vel_to_horiz(Phi, rhodz, u, W, F_el, uh, uz)
227    USE omp_para
228    REAL(rstd), INTENT(IN) :: Phi(iim*jjm,llm+1)
229    REAL(rstd), INTENT(IN) :: rhodz(iim*jjm,llm)
230    REAL(rstd), INTENT(IN) :: u(3*iim*jjm,llm)
231    REAL(rstd), INTENT(IN) :: W(iim*jjm,llm+1)
232    REAL(rstd), INTENT(OUT) :: uh(3*iim*jjm,llm)
233    REAL(rstd), INTENT(OUT) :: uz(iim*jjm,llm)
234    INTEGER :: ij,l
235    REAL(rstd) :: F_el(3*iim*jjm,llm+1)
236    REAL(rstd) :: uu_right, uu_lup, uu_ldown, W_el, DePhil
237    ! NB : u and uh are not in DEC form, they are normal components   
238    ! => we must divide by de
239    IF(hydrostatic) THEN
240       uh(:,:)=u(:,:)
241       uz(:,:)=0.
242    ELSE
243       DO l=ll_begin, ll_endp1 ! compute on l levels (interfaces)
244          DO ij=ij_begin_ext, ij_end_ext
245             ! Compute on edge 'right'
246             W_el   = .5*( W(ij,l)+W(ij+t_right,l) )
247             DePhil = ne_right*(Phi(ij+t_right,l)-Phi(ij,l))
248             F_el(ij+u_right,l)   = DePhil*W_el/de(ij+u_right)
249             ! Compute on edge 'lup'
250             W_el   = .5*( W(ij,l)+W(ij+t_lup,l) )
251             DePhil = ne_lup*(Phi(ij+t_lup,l)-Phi(ij,l))
252             F_el(ij+u_lup,l)   = DePhil*W_el/de(ij+u_lup)
253             ! Compute on edge 'ldown'
254             W_el   = .5*( W(ij,l)+W(ij+t_ldown,l) )
255             DePhil = ne_ldown*(Phi(ij+t_ldown,l)-Phi(ij,l))
256             F_el(ij+u_ldown,l) = DePhil*W_el/de(ij+u_ldown)
257          END DO
258       END DO
259
260       ! We need a barrier here because we compute F_el above and do a vertical average below
261       !$OMP BARRIER
262
263       DO l=ll_begin, ll_end ! compute on k levels (full levels)
264          DO ij=ij_begin_ext, ij_end_ext
265             ! w = vertical momentum = g^-2*dPhi/dt = uz/g
266             uz(ij,l) = (.5*g)*(W(ij,l)+W(ij,l+1))/rhodz(ij,l)
267             ! uh = u-w.grad(Phi) = u - uz.grad(z)
268             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))
269             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))
270             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))
271          END DO
272       END DO
273
274    END IF
275  END SUBROUTINE compute_prognostic_vel_to_horiz
276
277  SUBROUTINE diagnose_temperature(f_pmid,f_q,f_temp)
278    USE icosa
279    USE pression_mod
280    IMPLICIT NONE
281    TYPE(t_field), POINTER :: f_pmid(:)           ! IN
282    TYPE(t_field), POINTER :: f_q(:)            ! IN
283    TYPE(t_field), POINTER :: f_temp(:)         ! INOUT
284   
285    REAL(rstd), POINTER :: pmid(:,:)
286    REAL(rstd), POINTER :: q(:,:,:)
287    REAL(rstd), POINTER :: temp(:,:)
288    INTEGER :: ind
289   
290    DO ind=1,ndomain
291       IF (.NOT. assigned_domain(ind)) CYCLE
292       CALL swap_dimensions(ind)
293       CALL swap_geometry(ind)
294       pmid=f_pmid(ind)
295       q=f_q(ind)
296       temp=f_temp(ind)
297       CALL compute_diagnose_temp(pmid,q,temp)
298    END DO
299  END SUBROUTINE diagnose_temperature
300 
301  SUBROUTINE compute_diagnose_temp(pmid,q,temp)
302    USE omp_para
303    USE pression_mod
304    REAL(rstd),INTENT(IN)    :: pmid(iim*jjm,llm)
305    REAL(rstd),INTENT(IN)    :: q(iim*jjm,llm,nqtot)
306    REAL(rstd),INTENT(INOUT) :: temp(iim*jjm,llm)
307
308    REAL(rstd) :: Rd, p_ik, theta_ik, temp_ik, qv, chi, Rmix
309    INTEGER :: ij,l
310
311    Rd = kappa*cpp
312    DO l=ll_begin,ll_end
313       DO ij=ij_begin,ij_end
314          p_ik = pmid(ij,l)
315          theta_ik = temp(ij,l)
316          qv = q(ij,l,1) ! water vaper mixing ratio = mv/md
317          SELECT CASE(caldyn_thermo)
318          CASE(thermo_theta)
319             temp_ik = theta_ik*((p_ik/preff)**kappa)
320          CASE(thermo_entropy)
321             temp_ik = Treff*exp((theta_ik + Rd*log(p_ik/preff))/cpp)
322          CASE(thermo_moist)
323             Rmix = Rd+qv*Rv
324             chi = ( theta_ik + Rmix*log(p_ik/preff) ) / (cpp + qv*cppv) ! log(T/Treff)
325             temp_ik = Treff*exp(chi)
326          END SELECT
327          IF(physics_thermo==thermo_fake_moist) temp_ik=temp_ik/(1+0.608*qv)
328          temp(ij,l)=temp_ik
329       END DO
330    END DO
331  END SUBROUTINE compute_diagnose_temp
332
333  SUBROUTINE Tv2T(f_Tv, f_q, f_T)
334    TYPE(t_field), POINTER :: f_TV(:)
335    TYPE(t_field), POINTER :: f_q(:)
336    TYPE(t_field), POINTER :: f_T(:)
337   
338    REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:)
339    INTEGER :: ind
340   
341    DO ind=1,ndomain
342       IF (.NOT. assigned_domain(ind)) CYCLE
343       CALL swap_dimensions(ind)
344       CALL swap_geometry(ind)
345       Tv=f_Tv(ind)
346       T=f_T(ind)
347       SELECT CASE(physics_thermo)
348       CASE(thermo_dry)
349          T=Tv
350       CASE(thermo_fake_moist)
351          q=f_q(ind)
352          T=Tv/(1+0.608*q(:,:,1))
353       END SELECT
354    END DO
355  END SUBROUTINE Tv2T
356
357  SUBROUTINE divide_by_mass(iq, f_mass, f_theta_rhodz, f_theta)
358    INTEGER, INTENT(IN) :: iq
359    TYPE(t_field), POINTER :: f_mass(:), f_theta_rhodz(:), f_theta(:)
360    REAL(rstd), POINTER :: mass(:,:), theta_rhodz(:,:,:), theta(:,:)
361    INTEGER :: ind
362    DO ind=1,ndomain
363       IF (.NOT. assigned_domain(ind)) CYCLE
364       CALL swap_dimensions(ind)
365       CALL swap_geometry(ind)
366       mass=f_mass(ind)
367       theta_rhodz=f_theta_rhodz(ind)
368       theta=f_theta(ind)
369       theta(:,:) = theta_rhodz(:,:,iq) / mass(:,:)
370    END DO
371  END SUBROUTINE divide_by_mass
372   
373END MODULE observable_mod
Note: See TracBrowser for help on using the repository browser.