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

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

trunk : upgrading to devel

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