source: codes/icosagcm/trunk/src/observable.f90 @ 418

Last change on this file since 418 was 418, checked in by ymipsl, 8 years ago

Supress spurious openmp directive

YM

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_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_v,   field_z,type_real,llm, name="buf_v")
31    CALL allocate_field(f_buf_s,   field_t,type_real, name="buf_s")
32
33    CALL allocate_field(f_theta, field_t,type_real,llm,nqdyn,  name='theta') ! potential temperature
34    CALL allocate_field(f_pmid,  field_t,type_real,llm,  name='pmid')       ! mid layer pressure
35  END SUBROUTINE init_observable
36
37  SUBROUTINE write_output_fields_basic(init, f_phis, f_ps, f_mass, f_geopot, f_theta_rhodz, f_u, f_W, f_q)
38    USE xios
39    USE disvert_mod
40    USE wind_mod
41    USE output_field_mod
42    USE omp_para
43    USE time_mod
44    USE xios
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       scalar(1)=dt
63       CALL xios_send_field("timestep", scalar)
64       scalar(1)=preff
65       CALL xios_send_field("preff", scalar)
66       CALL xios_send_field("ap",ap)
67       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       CALL xios_send_field("mid_ap",mid_ap)
73       CALL xios_send_field("mid_bp",mid_bp)
74
75       CALL output_field("phis",f_phis)
76       CALL output_field("Ai",geom%Ai)       
77    END IF
78
79    CALL divide_by_mass(1, f_mass, f_theta_rhodz, f_buf_i)
80    IF(init) THEN
81       CALL output_field("theta_init",f_buf_i)
82    ELSE
83       CALL output_field("theta",f_buf_i)
84    END IF
85
86    IF(nqdyn>1) THEN
87       CALL divide_by_mass(2, f_mass, f_theta_rhodz, f_buf_i)
88       IF(init) THEN
89          CALL output_field("dyn_q_init",f_buf_i)
90       ELSE
91          CALL output_field("dyn_q",f_buf_i)
92       END IF
93    END IF
94
95    CALL theta_rhodz2temperature(f_ps,f_theta_rhodz,f_buf_i) ; 
96    CALL Tv2T(f_buf_i,f_q,f_buf1_i)
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_ps,f_buf_i,f_buf_s,85000.)
102       CALL output_field("t850",f_buf_s)
103       CALL vertical_interp(f_ps,f_buf_i,f_buf_s,50000.)
104       CALL output_field("t500",f_buf_s)
105       CALL vertical_interp(f_ps,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    CALL pression_mid(f_ps, f_pmid)
113    IF(init) THEN
114       CALL output_field("uz_init",f_buf_i)
115       CALL output_field("ulon_init",f_buf_ulon)
116       CALL output_field("ulat_init",f_buf_ulat)
117       CALL output_field("p_init",f_pmid)
118       CALL output_field("ps_init",f_ps)
119       CALL output_field("mass_init",f_mass)
120       CALL output_field("geopot_init",f_geopot)
121       CALL output_field("q_init",f_q)
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_ps,f_buf_ulon,f_buf_s,85000.)
136       CALL output_field("u850",f_buf_s)
137       CALL vertical_interp(f_ps,f_buf_ulon,f_buf_s,50000.)
138       CALL output_field("u500",f_buf_s)
139       
140       CALL vertical_interp(f_ps,f_buf_ulat,f_buf_s,85000.)
141       CALL output_field("v850",f_buf_s)
142       CALL vertical_interp(f_ps,f_buf_ulat,f_buf_s,50000.)
143       CALL output_field("v500",f_buf_s)
144
145       CALL vertical_interp(f_ps,f_buf_i,f_buf_s,85000.)
146       CALL output_field("w850",f_buf_s)
147       CALL vertical_interp(f_ps,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_ps,f_buf_i,f_buf_s,85000.)
153       CALL output_field("omega850",f_buf_s)
154       CALL vertical_interp(f_ps,f_buf_i,f_buf_s,50000.)
155       CALL output_field("omega500",f_buf_s)
156    END IF
157
158  END SUBROUTINE write_output_fields_basic
159 
160  SUBROUTINE write_output_fields(f_ps, f_phis, f_dps, f_u, f_theta_rhodz, f_q, &
161       f_buf_i, f_buf_v, f_buf_i3, f_buf1_i, f_buf2_i, f_buf_s, f_buf_p)
162    USE vorticity_mod
163    USE theta2theta_rhodz_mod
164    USE pression_mod
165    USE omega_mod
166    USE write_field_mod
167    USE vertical_interp_mod
168    USE wind_mod
169    TYPE(t_field),POINTER :: f_ps(:), f_phis(:), f_u(:), f_theta_rhodz(:), f_q(:), f_dps(:), &
170         f_buf_i(:), f_buf_v(:), f_buf_i3(:), f_buf1_i(:), f_buf2_i(:), f_buf_s(:), f_buf_p(:)
171   
172    REAL(rstd) :: out_pression_level
173    CHARACTER(LEN=255) :: str_pression
174    CHARACTER(LEN=255) :: physics_type
175   
176    out_pression_level=0.
177    CALL getin("out_pression_level",out_pression_level) 
178    WRITE(str_pression,*) INT(out_pression_level/100)
179    str_pression=ADJUSTL(str_pression)
180   
181    CALL writefield("ps",f_ps)
182    CALL writefield("dps",f_dps)
183    CALL writefield("phis",f_phis)
184    CALL vorticity(f_u,f_buf_v)
185    CALL writefield("vort",f_buf_v)
186
187    CALL w_omega(f_ps, f_u, f_buf_i)
188    CALL writefield('omega', f_buf_i)
189    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
190      CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
191      CALL writefield("omega"//TRIM(str_pression),f_buf_s)
192    ENDIF
193   
194    ! Temperature
195!    CALL theta_rhodz2temperature(f_ps,f_theta_rhodz,f_buf_i) ; ! FIXME
196     
197    CALL getin('physics',physics_type)
198    IF (TRIM(physics_type)=='dcmip') THEN
199      CALL Tv2T(f_buf_i,f_q,f_buf1_i) 
200      CALL writefield("T",f_buf1_i)
201      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
202        CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
203        CALL writefield("T"//TRIM(str_pression),f_buf_s)
204      ENDIF
205    ELSE
206      CALL writefield("T",f_buf_i)
207      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
208        CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
209        CALL writefield("T"//TRIM(str_pression),f_buf_s)
210      ENDIF
211    ENDIF
212   
213    ! velocity components
214    CALL un2ulonlat(f_u, f_buf1_i, f_buf2_i)
215    CALL writefield("ulon",f_buf1_i)
216    CALL writefield("ulat",f_buf2_i)
217
218    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
219      CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
220      CALL writefield("ulon"//TRIM(str_pression),f_buf_s)
221      CALL vertical_interp(f_ps,f_buf2_i,f_buf_s,out_pression_level)
222      CALL writefield("ulat"//TRIM(str_pression),f_buf_s)
223    ENDIF
224   
225    ! geopotential ! FIXME
226    CALL thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_buf_s,f_buf_p,f_buf1_i,f_buf2_i,f_buf_i)
227    CALL writefield("p",f_buf_p)
228!    CALL writefield("phi",f_geopot)   ! geopotential
229    CALL writefield("theta",f_buf1_i) ! potential temperature
230    CALL writefield("pk",f_buf2_i)    ! Exner pressure
231 
232  END SUBROUTINE write_output_fields
233
234!------------------- Conversion from prognostic to observable variables ------------------
235
236  SUBROUTINE progonostic_vel_to_horiz(f_geopot, f_ps, f_rhodz, f_u, f_W, f_uh, f_uz)
237    USE disvert_mod
238    TYPE(t_field), POINTER :: f_geopot(:), f_ps(:), f_rhodz(:), &
239         f_u(:), f_W(:), f_uz(:), &  ! IN
240         f_uh(:)                         ! OUT
241    REAL(rstd),POINTER :: geopot(:,:), ps(:), rhodz(:,:), u(:,:), W(:,:), uh(:,:), uz(:,:)
242    INTEGER :: ind
243   
244    DO ind=1,ndomain
245       IF (.NOT. assigned_domain(ind)) CYCLE
246       CALL swap_dimensions(ind)
247       CALL swap_geometry(ind)
248       geopot = f_geopot(ind)
249       rhodz = f_rhodz(ind)
250       u = f_u(ind)
251       W = f_W(ind)
252       uh  = f_uh(ind)
253       IF(caldyn_eta==eta_mass) THEN
254          ps=f_ps(ind)
255          CALL compute_rhodz(.TRUE., ps, rhodz)
256       END IF
257       uz = f_uz(ind)
258       CALL compute_prognostic_vel_to_horiz(geopot,rhodz,u,W,uh,uz)
259    END DO
260  END SUBROUTINE progonostic_vel_to_horiz
261 
262  SUBROUTINE compute_prognostic_vel_to_horiz(Phi, rhodz, u, W, uh, uz)
263    USE omp_para
264    REAL(rstd), INTENT(IN) :: Phi(iim*jjm,llm+1)
265    REAL(rstd), INTENT(IN) :: rhodz(iim*jjm,llm)
266    REAL(rstd), INTENT(IN) :: u(3*iim*jjm,llm)
267    REAL(rstd), INTENT(IN) :: W(iim*jjm,llm+1)
268    REAL(rstd), INTENT(OUT) :: uh(3*iim*jjm,llm)
269    REAL(rstd), INTENT(OUT) :: uz(iim*jjm,llm)
270    INTEGER :: ij,l
271    REAL(rstd) :: F_el(3*iim*jjm,llm+1)
272    REAL(rstd) :: uu_right, uu_lup, uu_ldown, W_el, DePhil
273    ! NB : u and uh are not in DEC form, they are normal components   
274    ! => we must divide by de
275    IF(hydrostatic) THEN
276       uh(:,:)=u(:,:)
277       uz(:,:)=0.
278    ELSE
279       DO l=ll_begin, ll_endp1 ! compute on l levels (interfaces)
280          DO ij=ij_begin_ext, ij_end_ext
281             ! Compute on edge 'right'
282             W_el   = .5*( W(ij,l)+W(ij+t_right,l) )
283             DePhil = ne_right*(Phi(ij+t_right,l)-Phi(ij,l))
284             F_el(ij+u_right,l)   = DePhil*W_el/de(ij+u_right)
285             ! Compute on edge 'lup'
286             W_el   = .5*( W(ij,l)+W(ij+t_lup,l) )
287             DePhil = ne_lup*(Phi(ij+t_lup,l)-Phi(ij,l))
288             F_el(ij+u_lup,l)   = DePhil*W_el/de(ij+u_lup)
289             ! Compute on edge 'ldown'
290             W_el   = .5*( W(ij,l)+W(ij+t_ldown,l) )
291             DePhil = ne_ldown*(Phi(ij+t_ldown,l)-Phi(ij,l))
292             F_el(ij+u_ldown,l) = DePhil*W_el/de(ij+u_ldown)
293          END DO
294       END DO
295
296       DO l=ll_begin, ll_end ! compute on k levels (full levels)
297          DO ij=ij_begin_ext, ij_end_ext
298             ! w = vertical momentum = g^-2*dPhi/dt = uz/g
299             uz(ij,l) = (.5*g)*(W(ij,l)+W(ij,l+1))/rhodz(ij,l)
300             ! uh = u-w.grad(Phi) = u - uz.grad(z)
301             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))
302             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))
303             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))
304          END DO
305       END DO
306
307    END IF
308  END SUBROUTINE compute_prognostic_vel_to_horiz
309
310  SUBROUTINE thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_pks,f_p,f_theta,f_pk,f_phi) 
311    USE field_mod
312    USE pression_mod
313    USE exner_mod
314    USE geopotential_mod
315    USE theta2theta_rhodz_mod
316    TYPE(t_field), POINTER :: f_ps(:), f_phis(:), f_theta_rhodz(:), &  ! IN
317         f_pks(:), f_p(:), f_theta(:), f_pk(:), f_phi(:)               ! OUT
318    REAL(rstd),POINTER :: pk(:,:), p(:,:), theta(:,:), theta_rhodz(:,:,:), &
319         phi(:,:), phis(:), ps(:), pks(:)
320    INTEGER :: ind
321   
322    DO ind=1,ndomain
323       IF (.NOT. assigned_domain(ind)) CYCLE
324       CALL swap_dimensions(ind)
325       CALL swap_geometry(ind)
326       ps = f_ps(ind)
327       p  = f_p(ind)
328!$OMP BARRIER
329       CALL compute_pression(ps,p,0)
330       pk = f_pk(ind)
331       pks = f_pks(ind)
332!$OMP BARRIER
333       CALL compute_exner(ps,p,pks,pk,0)
334!$OMP BARRIER
335       theta_rhodz = f_theta_rhodz(ind)
336       theta = f_theta(ind)
337       CALL compute_theta_rhodz2theta(ps, theta_rhodz(:,:,1),theta,0)
338       phis = f_phis(ind)
339       phi = f_phi(ind)
340       CALL compute_geopotential(phis,pks,pk,theta,phi,0)
341    END DO
342
343  END SUBROUTINE thetarhodz2geopot
344 
345  SUBROUTINE Tv2T(f_Tv, f_q, f_T)
346    TYPE(t_field), POINTER :: f_TV(:)
347    TYPE(t_field), POINTER :: f_q(:)
348    TYPE(t_field), POINTER :: f_T(:)
349   
350    REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:)
351    INTEGER :: ind
352   
353    DO ind=1,ndomain
354       IF (.NOT. assigned_domain(ind)) CYCLE
355       CALL swap_dimensions(ind)
356       CALL swap_geometry(ind)
357       Tv=f_Tv(ind)
358       q=f_q(ind)
359       T=f_T(ind)
360       T=Tv/(1+0.608*q(:,:,1))
361    END DO
362   
363  END SUBROUTINE Tv2T
364
365  SUBROUTINE divide_by_mass(iq, f_mass, f_theta_rhodz, f_theta)
366    INTEGER, INTENT(IN) :: iq
367    TYPE(t_field), POINTER :: f_mass(:), f_theta_rhodz(:), f_theta(:)
368    REAL(rstd), POINTER :: mass(:,:), theta_rhodz(:,:,:), theta(:,:)
369    INTEGER :: ind
370    DO ind=1,ndomain
371       IF (.NOT. assigned_domain(ind)) CYCLE
372       CALL swap_dimensions(ind)
373       CALL swap_geometry(ind)
374       mass=f_mass(ind)
375       theta_rhodz=f_theta_rhodz(ind)
376       theta=f_theta(ind)
377       theta(:,:) = theta_rhodz(:,:,iq) / mass(:,:)
378    END DO
379  END SUBROUTINE divide_by_mass
380   
381END MODULE observable_mod
Note: See TracBrowser for help on using the repository browser.