source: codes/icosagcm/trunk/src/caldyn_gcm.f90 @ 183

Last change on this file since 183 was 183, checked in by dubos, 11 years ago

moved Boussniesq calculation of pressure into compute_caldyn_horiz

File size: 37.9 KB
Line 
1MODULE caldyn_gcm_mod
2  USE icosa
3  USE transfert_mod
4  PRIVATE
5
6  INTEGER, PARAMETER :: energy=1, enstrophy=2
7  TYPE(t_field),POINTER :: f_out_u(:), f_qu(:), f_qv(:)
8  REAL(rstd),POINTER :: out_u(:,:), p(:,:), qu(:,:)
9
10  TYPE(t_field),POINTER :: f_buf_i(:), f_buf_ulon(:), f_buf_ulat(:), f_buf_u3d(:)
11  TYPE(t_field),POINTER :: f_buf_v(:), f_buf_s(:), f_buf_p(:)
12
13! temporary shared variable for caldyn
14  TYPE(t_field),POINTER :: f_theta(:)
15  TYPE(t_field),POINTER :: f_pk(:)
16  TYPE(t_field),POINTER :: f_geopot(:)
17  TYPE(t_field),POINTER :: f_wwuu(:)
18  TYPE(t_field),POINTER :: f_planetvel(:)
19   
20  INTEGER :: caldyn_conserv
21 
22  TYPE(t_message) :: req_ps, req_mass, req_theta_rhodz, req_u, req_qu
23
24  PUBLIC init_caldyn, caldyn_BC, caldyn, write_output_fields, &
25       req_ps, req_mass
26
27CONTAINS
28 
29  SUBROUTINE init_caldyn
30    USE icosa
31    USE exner_mod
32    USE mpipara
33    IMPLICIT NONE
34    CHARACTER(len=255) :: def
35    INTEGER            :: ind
36    REAL(rstd),POINTER :: planetvel(:)
37 
38    def='enstrophy'
39    CALL getin('caldyn_conserv',def)
40    SELECT CASE(TRIM(def))
41    CASE('energy')
42       caldyn_conserv=energy
43    CASE('enstrophy')
44       caldyn_conserv=enstrophy
45    CASE DEFAULT
46       IF (is_mpi_root) PRINT *,'Bad selector for variable caldyn_conserv : <', &
47            TRIM(def),'> options are <energy>, <enstrophy>'
48       STOP
49    END SELECT
50    IF (is_mpi_root) PRINT *, 'caldyn_conserv=',def
51
52    CALL allocate_caldyn
53
54    DO ind=1,ndomain
55       CALL swap_dimensions(ind)
56       CALL swap_geometry(ind)
57       planetvel=f_planetvel(ind)
58       CALL compute_planetvel(planetvel)
59    END DO
60
61  END SUBROUTINE init_caldyn
62
63  SUBROUTINE allocate_caldyn
64  USE icosa
65  IMPLICIT NONE
66
67    CALL allocate_field(f_out_u,field_u,type_real,llm) 
68    CALL allocate_field(f_qu,field_u,type_real,llm) 
69    CALL allocate_field(f_qv,field_z,type_real,llm) 
70 
71    CALL allocate_field(f_buf_i,   field_t,type_real,llm)
72    CALL allocate_field(f_buf_p,   field_t,type_real,llm+1) 
73    CALL allocate_field(f_buf_u3d, field_t,type_real,3,llm)  ! 3D vel at cell centers
74    CALL allocate_field(f_buf_ulon,field_t,type_real,llm)
75    CALL allocate_field(f_buf_ulat,field_t,type_real,llm)
76    CALL allocate_field(f_buf_v,   field_z,type_real,llm)
77    CALL allocate_field(f_buf_s,   field_t,type_real)
78
79    CALL allocate_field(f_theta, field_t,type_real,llm,  name='theta')   ! potential temperature
80    CALL allocate_field(f_pk,    field_t,type_real,llm,  name='pk')
81    CALL allocate_field(f_geopot,field_t,type_real,llm+1,name='geopot')  ! geopotential
82    CALL allocate_field(f_wwuu,  field_u,type_real,llm+1,name='wwuu')
83    CALL allocate_field(f_planetvel,  field_u,type_real, name='planetvel') ! planetary velocity at r=a
84
85  END SUBROUTINE allocate_caldyn
86
87  SUBROUTINE caldyn_BC(f_phis, f_wflux)
88    USE icosa
89    USE mpipara
90    USE omp_para
91    IMPLICIT NONE
92    TYPE(t_field),POINTER :: f_phis(:)
93    TYPE(t_field),POINTER :: f_wflux(:)
94    REAL(rstd),POINTER  :: phis(:)
95    REAL(rstd),POINTER  :: wflux(:,:)
96    REAL(rstd),POINTER  :: geopot(:,:)
97    REAL(rstd),POINTER  :: wwuu(:,:)
98
99    INTEGER :: ind,i,j,ij,l
100
101    IF (omp_first) THEN
102       DO ind=1,ndomain
103          CALL swap_dimensions(ind)
104          CALL swap_geometry(ind)
105          geopot=f_geopot(ind)
106          phis=f_phis(ind)
107          wflux=f_wflux(ind)
108          wwuu=f_wwuu(ind)
109         
110          DO ij=ij_begin_ext,ij_end_ext
111              ! lower BCs : geopot=phis, wflux=0, wwuu=0
112              geopot(ij,1) = phis(ij)
113              wflux(ij,1) = 0.
114              wwuu(ij+u_right,1)=0   
115              wwuu(ij+u_lup,1)=0   
116              wwuu(ij+u_ldown,1)=0
117              ! top BCs : wflux=0, wwuu=0
118              wflux(ij,llm+1)  = 0.
119              wwuu(ij+u_right,llm+1)=0   
120              wwuu(ij+u_lup,llm+1)=0   
121              wwuu(ij+u_ldown,llm+1)=0
122          ENDDO
123       END DO
124    ENDIF
125
126    !$OMP BARRIER
127  END SUBROUTINE caldyn_BC
128   
129  SUBROUTINE caldyn(write_out,f_phis, f_ps, f_mass, f_theta_rhodz, f_u, f_q, &
130       f_hflux, f_wflux, f_dps, f_dmass, f_dtheta_rhodz, f_du)
131    USE icosa
132    USE disvert_mod, ONLY : caldyn_eta, eta_mass
133    USE vorticity_mod
134    USE kinetic_mod
135    USE theta2theta_rhodz_mod
136    USE mpipara
137    USE trace
138    USE omp_para
139    USE output_field_mod
140    IMPLICIT NONE
141    LOGICAL,INTENT(IN)    :: write_out
142    TYPE(t_field),POINTER :: f_phis(:)
143    TYPE(t_field),POINTER :: f_ps(:)
144    TYPE(t_field),POINTER :: f_mass(:)
145    TYPE(t_field),POINTER :: f_theta_rhodz(:)
146    TYPE(t_field),POINTER :: f_u(:)
147    TYPE(t_field),POINTER :: f_q(:)
148    TYPE(t_field),POINTER :: f_hflux(:), f_wflux(:)
149    TYPE(t_field),POINTER :: f_dps(:)
150    TYPE(t_field),POINTER :: f_dmass(:)
151    TYPE(t_field),POINTER :: f_dtheta_rhodz(:)
152    TYPE(t_field),POINTER :: f_du(:)
153   
154    REAL(rstd),POINTER :: ps(:), dps(:)
155    REAL(rstd),POINTER :: mass(:,:), theta_rhodz(:,:), dtheta_rhodz(:,:)
156    REAL(rstd),POINTER :: u(:,:), du(:,:), hflux(:,:), wflux(:,:)
157    REAL(rstd),POINTER :: qu(:,:)
158    REAL(rstd),POINTER :: qv(:,:)
159
160! temporary shared variable
161    REAL(rstd),POINTER  :: theta(:,:) 
162    REAL(rstd),POINTER  :: pk(:,:)
163    REAL(rstd),POINTER  :: geopot(:,:)
164    REAL(rstd),POINTER  :: convm(:,:) 
165    REAL(rstd),POINTER  :: wwuu(:,:)
166       
167    INTEGER :: ind
168    LOGICAL,SAVE :: first=.TRUE.
169!$OMP THREADPRIVATE(first)
170   
171    ! MPI messages need to be sent at first call to caldyn
172    ! This is needed only once : the next ones will be sent by timeloop
173    IF (first) THEN
174      first=.FALSE.
175      IF(caldyn_eta==eta_mass) THEN
176         CALL init_message(f_ps,req_i1,req_ps)
177      ELSE
178         CALL init_message(f_mass,req_i1,req_mass)
179      END IF
180      CALL init_message(f_theta_rhodz,req_i1,req_theta_rhodz)
181      CALL init_message(f_u,req_e1_vect,req_u)
182      CALL init_message(f_qu,req_e1_scal,req_qu)
183      IF(caldyn_eta==eta_mass) THEN
184         CALL send_message(f_ps,req_ps) 
185      ELSE
186         CALL send_message(f_mass,req_mass) 
187      END IF
188    ENDIF
189   
190    CALL trace_start("caldyn")
191
192    CALL send_message(f_u,req_u)
193    CALL send_message(f_theta_rhodz,req_theta_rhodz) 
194
195    SELECT CASE(caldyn_conserv)
196    CASE(energy) ! energy-conserving
197       DO ind=1,ndomain
198          CALL swap_dimensions(ind)
199          CALL swap_geometry(ind)
200          ps=f_ps(ind)
201          u=f_u(ind)
202          theta_rhodz = f_theta_rhodz(ind)
203          mass=f_mass(ind)
204          theta = f_theta(ind)
205          qu=f_qu(ind)
206          qv=f_qv(ind)
207          CALL compute_pvort(ps,u,theta_rhodz, mass,theta,qu,qv)
208       ENDDO
209
210       CALL send_message(f_qu,req_qu)
211
212       DO ind=1,ndomain
213          CALL swap_dimensions(ind)
214          CALL swap_geometry(ind)
215          ps=f_ps(ind)
216          u=f_u(ind)
217          theta_rhodz=f_theta_rhodz(ind)
218          mass=f_mass(ind)
219          theta = f_theta(ind)
220          qu=f_qu(ind)
221          pk = f_pk(ind)
222          geopot = f_geopot(ind) 
223          CALL compute_geopot(ps,mass,theta, pk,geopot)
224          hflux=f_hflux(ind)
225          convm = f_dmass(ind)
226          dtheta_rhodz=f_dtheta_rhodz(ind)
227          du=f_du(ind)
228          CALL compute_caldyn_horiz(u,mass,qu,theta,pk,geopot, hflux,convm,dtheta_rhodz,du)
229          IF(caldyn_eta==eta_mass) THEN
230             wflux=f_wflux(ind)
231             wwuu=f_wwuu(ind)
232             dps=f_dps(ind)
233             CALL compute_caldyn_vert(u,theta,mass,convm, wflux,wwuu, dps, dtheta_rhodz, du)
234          END IF
235       ENDDO       
236       
237    CASE(enstrophy) ! enstrophy-conserving
238       DO ind=1,ndomain
239          CALL swap_dimensions(ind)
240          CALL swap_geometry(ind)
241          ps=f_ps(ind)
242          u=f_u(ind)
243          theta_rhodz=f_theta_rhodz(ind)
244          mass=f_mass(ind)
245          theta = f_theta(ind)
246          qu=f_qu(ind)
247          qv=f_qv(ind)
248          CALL compute_pvort(ps,u,theta_rhodz, mass,theta,qu,qv)
249          pk = f_pk(ind)
250          geopot = f_geopot(ind) 
251          CALL compute_geopot(ps,mass,theta, pk,geopot)
252          hflux=f_hflux(ind)
253          convm = f_dmass(ind)
254          dtheta_rhodz=f_dtheta_rhodz(ind)
255          du=f_du(ind)
256          CALL compute_caldyn_horiz(u,mass,qu,theta,pk,geopot, hflux,convm,dtheta_rhodz,du)
257          IF(caldyn_eta==eta_mass) THEN
258             wflux=f_wflux(ind)
259             wwuu=f_wwuu(ind)
260             dps=f_dps(ind)
261             CALL compute_caldyn_vert(u,theta,mass,convm, wflux,wwuu, dps, dtheta_rhodz, du)
262          END IF
263       ENDDO
264       
265    CASE DEFAULT
266       STOP
267    END SELECT
268
269!$OMP BARRIER
270    IF (write_out) THEN
271
272       IF (is_mpi_root) PRINT *,'CALL write_output_fields'
273
274! ---> for openMP test to fix later
275!       CALL write_output_fields(f_ps, f_phis, f_dps, f_u, f_theta_rhodz, f_q, &
276!            f_buf_i, f_buf_v, f_buf_u3d, f_buf_ulon, f_buf_ulat, f_buf_s, f_buf_p)
277
278       CALL output_field("ps",f_ps)
279       CALL output_field("dps",f_dps)
280       CALL output_field("mass",f_mass)
281       CALL output_field("dmass",f_dmass)
282       CALL output_field("vort",f_qv)
283       CALL output_field("theta",f_theta)
284       CALL output_field("exner",f_pk)
285       CALL output_field("pv",f_qv)
286 
287    END IF
288   
289    !    CALL check_mass_conservation(f_ps,f_dps)
290    CALL trace_end("caldyn")
291!$OMP BARRIER
292   
293END SUBROUTINE caldyn
294
295SUBROUTINE compute_planetvel(planetvel)
296  USE wind_mod
297  REAL(rstd),INTENT(OUT)  :: planetvel(iim*3*jjm)
298  REAL(rstd) :: ulon(iim*3*jjm)
299  REAL(rstd) :: ulat(iim*3*jjm)
300  REAL(rstd) :: lon,lat
301  INTEGER :: ij
302  DO ij=ij_begin_ext,ij_end_ext
303     CALL xyz2lonlat(xyz_e(ij+u_right,:),lon,lat)
304     ulon(ij+u_right)=a*omega*cos(lat)
305     ulat(ij+u_right)=0
306
307     CALL xyz2lonlat(xyz_e(ij+u_lup,:),lon,lat)
308     ulon(ij+u_lup)=a*omega*cos(lat)
309     ulat(ij+u_lup)=0
310
311     CALL xyz2lonlat(xyz_e(ij+u_ldown,:),lon,lat)
312     ulon(ij+u_ldown)=a*omega*cos(lat)
313     ulat(ij+u_ldown)=0
314  END DO
315  CALL compute_wind2D_perp_from_lonlat_compound(ulon, ulat, planetvel)
316END SUBROUTINE compute_planetvel
317   
318SUBROUTINE compute_pvort(ps,u,theta_rhodz, rhodz,theta,qu,qv)
319  USE icosa
320  USE disvert_mod, ONLY : mass_dak, mass_dbk, caldyn_eta, eta_mass
321  USE exner_mod
322  USE trace
323  USE omp_para
324  IMPLICIT NONE
325  REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)
326  REAL(rstd),INTENT(IN)  :: ps(iim*jjm)
327  REAL(rstd),INTENT(IN)  :: theta_rhodz(iim*jjm,llm)
328  REAL(rstd),INTENT(INOUT) :: rhodz(iim*jjm,llm)
329  REAL(rstd),INTENT(OUT) :: theta(iim*jjm,llm)
330  REAL(rstd),INTENT(OUT) :: qu(iim*3*jjm,llm)
331  REAL(rstd),INTENT(OUT) :: qv(iim*2*jjm,llm)
332 
333  INTEGER :: i,j,ij,l
334  REAL(rstd) :: etav,hv, m
335!  REAL(rstd) :: qv(2*iim*jjm,llm)     ! potential velocity 
336 
337  CALL trace_start("compute_pvort") 
338
339  IF(caldyn_eta==eta_mass) THEN
340     CALL wait_message(req_ps) 
341  ELSE
342     CALL wait_message(req_mass)
343  END IF
344  CALL wait_message(req_theta_rhodz) 
345
346  IF(caldyn_eta==eta_mass) THEN ! Compute mass & theta
347     DO l = ll_begin,ll_end
348        CALL test_message(req_u) 
349!DIR$ SIMD
350        DO ij=ij_begin_ext,ij_end_ext
351           m = ( mass_dak(l)+ps(ij)*mass_dbk(l) )/g
352           rhodz(ij,l) = m
353           theta(ij,l) = theta_rhodz(ij,l)/rhodz(ij,l)
354        ENDDO
355     ENDDO
356  ELSE ! Compute only theta
357     DO l = ll_begin,ll_end
358        CALL test_message(req_u) 
359!DIR$ SIMD
360       DO ij=ij_begin_ext,ij_end_ext
361         theta(ij,l) = theta_rhodz(ij,l)/rhodz(ij,l)
362       ENDDO
363     ENDDO
364  END IF
365
366  CALL wait_message(req_u)   
367 
368!!! Compute shallow-water potential vorticity
369  DO l = ll_begin,ll_end
370!DIR$ SIMD
371     DO ij=ij_begin_ext,ij_end_ext
372          etav= 1./Av(ij+z_up)*(  ne_rup        * u(ij+u_rup,l)        * de(ij+u_rup)         &
373                                + ne_left * u(ij+t_rup+u_left,l) * de(ij+t_rup+u_left)  &
374                                - ne_lup        * u(ij+u_lup,l)        * de(ij+u_lup) )                               
375
376          hv =  Riv2(ij,vup)          * rhodz(ij,l)            &
377              + Riv2(ij+t_rup,vldown) * rhodz(ij+t_rup,l)     &
378              + Riv2(ij+t_lup,vrdown) * rhodz(ij+t_lup,l)
379     
380          qv(ij+z_up,l) = ( etav+fv(ij+z_up) )/hv
381         
382          etav = 1./Av(ij+z_down)*(  ne_ldown         * u(ij+u_ldown,l)          * de(ij+u_ldown)          &
383                                   + ne_right * u(ij+t_ldown+u_right,l)  * de(ij+t_ldown+u_right)  &
384                                   - ne_rdown         * u(ij+u_rdown,l)          * de(ij+u_rdown) )
385       
386          hv =  Riv2(ij,vdown)        * rhodz(ij,l)          &
387              + Riv2(ij+t_ldown,vrup) * rhodz(ij+t_ldown,l)  &
388              + Riv2(ij+t_rdown,vlup) * rhodz(ij+t_rdown,l)
389                       
390          qv(ij+z_down,l) =( etav+fv(ij+z_down) )/hv
391
392      ENDDO
393
394!DIR$ SIMD
395      DO ij=ij_begin,ij_end
396         qu(ij+u_right,l) = 0.5*(qv(ij+z_rdown,l)+qv(ij+z_rup,l)) 
397         qu(ij+u_lup,l) = 0.5*(qv(ij+z_up,l)+qv(ij+z_lup,l)) 
398         qu(ij+u_ldown,l) = 0.5*(qv(ij+z_ldown,l)+qv(ij+z_down,l)) 
399      END DO
400     
401   ENDDO
402
403   CALL trace_end("compute_pvort")
404                                                       
405  END SUBROUTINE compute_pvort
406 
407  SUBROUTINE compute_geopot(ps,rhodz,theta, pk,geopot) 
408  USE icosa
409  USE disvert_mod
410  USE exner_mod
411  USE trace
412  USE omp_para
413  IMPLICIT NONE
414    REAL(rstd),INTENT(INOUT) :: ps(iim*jjm)
415    REAL(rstd),INTENT(IN)    :: rhodz(iim*jjm,llm)
416    REAL(rstd),INTENT(IN)    :: theta(iim*jjm,llm)    ! potential temperature
417    REAL(rstd),INTENT(INOUT) :: pk(iim*jjm,llm)       ! Exner function
418    REAL(rstd),INTENT(INOUT) :: geopot(iim*jjm,llm+1) ! geopotential
419   
420    INTEGER :: i,j,ij,l
421    REAL(rstd) :: p_ik, exner_ik
422
423    CALL trace_start("compute_geopot")
424
425    IF(caldyn_eta==eta_mass) THEN
426
427!!! Compute exner function and geopotential
428       DO l = 1,llm
429          !$OMP DO SCHEDULE(STATIC)
430          !DIR$ SIMD
431          DO ij=ij_begin_ext,ij_end_ext         
432             p_ik = ptop + mass_ak(l) + mass_bk(l)*ps(ij) ! FIXME : leave ps for the moment ; change ps to Ms later
433             !         p_ik = ptop + g*(mass_ak(l)+ mass_bk(l)*ps(i,j))
434             exner_ik = cpp * (p_ik/preff) ** kappa
435             pk(ij,l) = exner_ik
436             ! specific volume v = kappa*theta*pi/p = dphi/g/rhodz
437             geopot(ij,l+1) = geopot(ij,l) + (g*kappa)*rhodz(ij,l)*theta(ij,l)*exner_ik/p_ik
438          ENDDO
439       ENDDO
440
441    ELSE 
442       ! We are using a Lagrangian vertical coordinate
443       ! Pressure must be computed first top-down (temporarily stored in pk)
444       ! Then Exner pressure and geopotential are computed bottom-up
445       ! Notice that the computation below should work also when caldyn_eta=eta_mass         
446
447       IF(boussinesq) THEN ! compute only geopotential : pressure pk will be computed in compute_caldyn_horiz
448          ! specific volume 1 = dphi/g/rhodz
449          DO l = 1,llm
450             !$OMP DO SCHEDULE(STATIC)
451             !DIR$ SIMD
452             DO ij=ij_begin_ext,ij_end_ext         
453                geopot(ij,l+1) = geopot(ij,l) + g*rhodz(ij,l)
454             ENDDO
455          ENDDO
456       ELSE ! non-Boussinesq, compute geopotential and Exner pressure
457          ! uppermost layer
458          !DIR$ SIMD
459          DO ij=ij_begin_ext,ij_end_ext         
460             pk(ij,llm) = ptop + (.5*g)*rhodz(ij,llm)
461          END DO
462          ! other layers
463          DO l = llm-1, 1, -1
464             !$OMP DO SCHEDULE(STATIC)
465             !DIR$ SIMD
466             DO ij=ij_begin_ext,ij_end_ext         
467                pk(ij,l) = pk(ij,l+1) + (.5*g)*(rhodz(ij,l)+rhodz(ij,l+1))
468             END DO
469          END DO
470          ! surface pressure (for diagnostics)
471          DO ij=ij_begin_ext,ij_end_ext         
472             ps(ij) = pk(ij,1) + (.5*g)*rhodz(ij,1)
473          END DO
474
475          ! specific volume v = kappa*theta*pi/p = dphi/g/rhodz
476          DO l = 1,llm
477             !$OMP DO SCHEDULE(STATIC)
478             !DIR$ SIMD
479             DO ij=ij_begin_ext,ij_end_ext         
480                p_ik = pk(ij,l)
481                exner_ik = cpp * (p_ik/preff) ** kappa
482                geopot(ij,l+1) = geopot(ij,l) + (g*kappa)*rhodz(ij,l)*theta(ij,l)*exner_ik/p_ik 
483                pk(ij,l) = exner_ik
484             ENDDO
485          ENDDO
486       END IF
487
488    END IF
489
490  CALL trace_end("compute_geopot")
491
492  END SUBROUTINE compute_geopot
493
494  SUBROUTINE compute_caldyn_horiz(u,rhodz,qu,theta,pk,geopot, hflux,convm, dtheta_rhodz, du)
495  USE icosa
496  USE disvert_mod
497  USE exner_mod
498  USE trace
499  USE omp_para
500  IMPLICIT NONE
501    REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)    ! prognostic "velocity"
502    REAL(rstd),INTENT(IN)  :: rhodz(iim*jjm,llm)
503    REAL(rstd),INTENT(IN)  :: qu(iim*3*jjm,llm)
504    REAL(rstd),INTENT(IN)  :: theta(iim*jjm,llm)  ! potential temperature
505    REAL(rstd),INTENT(INOUT) :: pk(iim*jjm,llm) ! Exner function
506    REAL(rstd),INTENT(IN)  :: geopot(iim*jjm,llm+1)    ! geopotential
507
508    REAL(rstd),INTENT(OUT) :: hflux(iim*3*jjm,llm) ! hflux in kg/s
509    REAL(rstd),INTENT(OUT) :: convm(iim*jjm,llm)  ! mass flux convergence
510    REAL(rstd),INTENT(OUT) :: dtheta_rhodz(iim*jjm,llm)
511    REAL(rstd),INTENT(OUT) :: du(iim*3*jjm,llm)
512
513    REAL(rstd) :: cor_NT(iim*jjm,llm)  ! NT coriolis force u.(du/dPhi)
514    REAL(rstd) :: urel(3*iim*jjm,llm)  ! relative velocity
515    REAL(rstd) :: Ftheta(3*iim*jjm,llm) ! theta flux
516    REAL(rstd) :: berni(iim*jjm,llm)  ! Bernoulli function
517   
518    INTEGER :: i,j,ij,l
519    REAL(rstd) :: ww,uu
520
521    CALL trace_start("compute_caldyn_horiz")
522
523    CALL wait_message(req_theta_rhodz) 
524
525  DO l = ll_begin, ll_end
526!!!  Compute mass and theta fluxes
527    IF (caldyn_conserv==energy) CALL test_message(req_qu) 
528!DIR$ SIMD
529    DO ij=ij_begin_ext,ij_end_ext         
530        hflux(ij+u_right,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_right,l))*u(ij+u_right,l)*le(ij+u_right)
531        hflux(ij+u_lup,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_lup,l))*u(ij+u_lup,l)*le(ij+u_lup)
532        hflux(ij+u_ldown,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_ldown,l))*u(ij+u_ldown,l)*le(ij+u_ldown)
533
534        Ftheta(ij+u_right,l)=0.5*(theta(ij,l)+theta(ij+t_right,l))*hflux(ij+u_right,l)
535        Ftheta(ij+u_lup,l)=0.5*(theta(ij,l)+theta(ij+t_lup,l))*hflux(ij+u_lup,l)
536        Ftheta(ij+u_ldown,l)=0.5*(theta(ij,l)+theta(ij+t_ldown,l))*hflux(ij+u_ldown,l)
537    ENDDO
538
539!!! compute horizontal divergence of fluxes
540!DIR$ SIMD
541    DO ij=ij_begin,ij_end         
542        ! convm = -div(mass flux), sign convention as in Ringler et al. 2012, eq. 21
543        convm(ij,l)= -1./Ai(ij)*(ne_right*hflux(ij+u_right,l)   +  &
544                                ne_rup*hflux(ij+u_rup,l)       +  & 
545                                ne_lup*hflux(ij+u_lup,l)       +  & 
546                                ne_left*hflux(ij+u_left,l)     +  & 
547                                ne_ldown*hflux(ij+u_ldown,l)   +  & 
548                                ne_rdown*hflux(ij+u_rdown,l))   
549
550        ! signe ? attention d (rho theta dz)
551        ! dtheta_rhodz =  -div(flux.theta)
552        dtheta_rhodz(ij,l)=-1./Ai(ij)*(ne_right*Ftheta(ij+u_right,l)    +  &
553                                 ne_rup*Ftheta(ij+u_rup,l)        +  & 
554                                 ne_lup*Ftheta(ij+u_lup,l)        +  & 
555                                 ne_left*Ftheta(ij+u_left,l)      +  & 
556                                 ne_ldown*Ftheta(ij+u_ldown,l)    +  & 
557                                 ne_rdown*Ftheta(ij+u_rdown,l))
558    ENDDO
559
560 END DO
561
562!!! Compute potential vorticity (Coriolis) contribution to du
563 
564    SELECT CASE(caldyn_conserv)
565    CASE(energy) ! energy-conserving TRiSK
566
567       CALL wait_message(req_qu)
568
569        DO l=ll_begin,ll_end
570!DIR$ SIMD
571          DO ij=ij_begin,ij_end         
572   
573                uu = wee(ij+u_right,1,1)*hflux(ij+u_rup,l)*(qu(ij+u_right,l)+qu(ij+u_rup,l))+                             &
574                     wee(ij+u_right,2,1)*hflux(ij+u_lup,l)*(qu(ij+u_right,l)+qu(ij+u_lup,l))+                             &
575                     wee(ij+u_right,3,1)*hflux(ij+u_left,l)*(qu(ij+u_right,l)+qu(ij+u_left,l))+                           &
576                     wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)*(qu(ij+u_right,l)+qu(ij+u_ldown,l))+                         &
577                     wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)*(qu(ij+u_right,l)+qu(ij+u_rdown,l))+                         & 
578                     wee(ij+u_right,1,2)*hflux(ij+t_right+u_ldown,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_ldown,l))+         &
579                     wee(ij+u_right,2,2)*hflux(ij+t_right+u_rdown,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_rdown,l))+         &
580                     wee(ij+u_right,3,2)*hflux(ij+t_right+u_right,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_right,l))+         &
581                     wee(ij+u_right,4,2)*hflux(ij+t_right+u_rup,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_rup,l))+             &
582                     wee(ij+u_right,5,2)*hflux(ij+t_right+u_lup,l)*(qu(ij+u_right,l)+qu(ij+t_right+u_lup,l))
583                du(ij+u_right,l) = .5*uu/de(ij+u_right)
584               
585                uu = wee(ij+u_lup,1,1)*hflux(ij+u_left,l)*(qu(ij+u_lup,l)+qu(ij+u_left,l)) +                        &
586                     wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)*(qu(ij+u_lup,l)+qu(ij+u_ldown,l)) +                      &
587                     wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)*(qu(ij+u_lup,l)+qu(ij+u_rdown,l)) +                      &
588                     wee(ij+u_lup,4,1)*hflux(ij+u_right,l)*(qu(ij+u_lup,l)+qu(ij+u_right,l)) +                      &
589                     wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)*(qu(ij+u_lup,l)+qu(ij+u_rup,l)) +                          & 
590                     wee(ij+u_lup,1,2)*hflux(ij+t_lup+u_right,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_right,l)) +          &
591                     wee(ij+u_lup,2,2)*hflux(ij+t_lup+u_rup,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_rup,l)) +              &
592                     wee(ij+u_lup,3,2)*hflux(ij+t_lup+u_lup,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_lup,l)) +              &
593                     wee(ij+u_lup,4,2)*hflux(ij+t_lup+u_left,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_left,l)) +            &
594                     wee(ij+u_lup,5,2)*hflux(ij+t_lup+u_ldown,l)*(qu(ij+u_lup,l)+qu(ij+t_lup+u_ldown,l))
595                du(ij+u_lup,l) = .5*uu/de(ij+u_lup)
596
597               
598                uu = wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)*(qu(ij+u_ldown,l)+qu(ij+u_rdown,l)) +                      &
599                     wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)*(qu(ij+u_ldown,l)+qu(ij+u_right,l)) +                      &
600                     wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)*(qu(ij+u_ldown,l)+qu(ij+u_rup,l)) +                          &
601                     wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)*(qu(ij+u_ldown,l)+qu(ij+u_lup,l)) +                          &
602                     wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)*(qu(ij+u_ldown,l)+qu(ij+u_left,l)) +                        & 
603                     wee(ij+u_ldown,1,2)*hflux(ij+t_ldown+u_lup,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_lup,l)) +          &
604                     wee(ij+u_ldown,2,2)*hflux(ij+t_ldown+u_left,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_left,l)) +        &
605                     wee(ij+u_ldown,3,2)*hflux(ij+t_ldown+u_ldown,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_ldown,l)) +      &
606                     wee(ij+u_ldown,4,2)*hflux(ij+t_ldown+u_rdown,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_rdown,l)) +      &
607                     wee(ij+u_ldown,5,2)*hflux(ij+t_ldown+u_right,l)*(qu(ij+u_ldown,l)+qu(ij+t_ldown+u_right,l))
608                du(ij+u_ldown,l) = .5*uu/de(ij+u_ldown)
609               
610          ENDDO
611       ENDDO
612       
613    CASE(enstrophy) ! enstrophy-conserving TRiSK
614 
615        DO l=ll_begin,ll_end
616!DIR$ SIMD
617          DO ij=ij_begin,ij_end         
618
619                uu = wee(ij+u_right,1,1)*hflux(ij+u_rup,l)+                           &
620                     wee(ij+u_right,2,1)*hflux(ij+u_lup,l)+                           &
621                     wee(ij+u_right,3,1)*hflux(ij+u_left,l)+                          &
622                     wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)+                         &
623                     wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)+                         & 
624                     wee(ij+u_right,1,2)*hflux(ij+t_right+u_ldown,l)+                 &
625                     wee(ij+u_right,2,2)*hflux(ij+t_right+u_rdown,l)+                 &
626                     wee(ij+u_right,3,2)*hflux(ij+t_right+u_right,l)+                 &
627                     wee(ij+u_right,4,2)*hflux(ij+t_right+u_rup,l)+                   &
628                     wee(ij+u_right,5,2)*hflux(ij+t_right+u_lup,l)
629                du(ij+u_right,l) = qu(ij+u_right,l)*uu/de(ij+u_right)
630               
631               
632                uu = wee(ij+u_lup,1,1)*hflux(ij+u_left,l)+                        &
633                     wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)+                       &
634                     wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)+                       &
635                     wee(ij+u_lup,4,1)*hflux(ij+u_right,l)+                       &
636                     wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)+                         & 
637                     wee(ij+u_lup,1,2)*hflux(ij+t_lup+u_right,l)+                 &
638                     wee(ij+u_lup,2,2)*hflux(ij+t_lup+u_rup,l)+                   &
639                     wee(ij+u_lup,3,2)*hflux(ij+t_lup+u_lup,l)+                   &
640                     wee(ij+u_lup,4,2)*hflux(ij+t_lup+u_left,l)+                  &
641                     wee(ij+u_lup,5,2)*hflux(ij+t_lup+u_ldown,l)
642                du(ij+u_lup,l) = qu(ij+u_lup,l)*uu/de(ij+u_lup)
643
644                uu = wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)+                      &
645                     wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)+                      &
646                     wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)+                        &
647                     wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)+                        &
648                     wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)+                       & 
649                     wee(ij+u_ldown,1,2)*hflux(ij+t_ldown+u_lup,l)+                &
650                     wee(ij+u_ldown,2,2)*hflux(ij+t_ldown+u_left,l)+               &
651                     wee(ij+u_ldown,3,2)*hflux(ij+t_ldown+u_ldown,l)+              &
652                     wee(ij+u_ldown,4,2)*hflux(ij+t_ldown+u_rdown,l)+              &
653                     wee(ij+u_ldown,5,2)*hflux(ij+t_ldown+u_right,l)
654                du(ij+u_ldown,l) = qu(ij+u_ldown,l)*uu/de(ij+u_ldown)
655     
656          ENDDO
657       ENDDO
658       
659    CASE DEFAULT
660       STOP
661    END SELECT
662 
663!!! Compute bernouilli term = Kinetic Energy + geopotential
664    IF(boussinesq) THEN
665       ! first use hydrostatic balance with theta*rhodz to find pk (Lagrange multiplier=pressure)
666       ! uppermost layer
667       !DIR$ SIMD
668       DO ij=ij_begin_ext,ij_end_ext         
669          pk(ij,llm) = ptop + (.5*g)*theta(ij,llm)*rhodz(ij,llm)
670       END DO
671       ! other layers
672       DO l = llm-1, 1, -1
673          !$OMP DO SCHEDULE(STATIC)
674          !DIR$ SIMD
675          DO ij=ij_begin_ext,ij_end_ext         
676             pk(ij,l) = pk(ij,l+1) + (.5*g)*(theta(ij,l)*rhodz(ij,l)+theta(ij,l+1)*rhodz(ij,l+1))
677          END DO
678       END DO
679       ! surface pressure (for diagnostics) FIXME
680       ! DO ij=ij_begin_ext,ij_end_ext         
681       !   ps(ij) = pk(ij,1) + (.5*g)*theta(ij,1)*rhodz(ij,1)
682       ! END DO
683       ! now pk contains the Lagrange multiplier (pressure)
684
685       DO l=ll_begin,ll_end
686!DIR$ SIMD
687          DO ij=ij_begin,ij_end         
688                berni(ij,l) = pk(ij,l) &
689                     + 1/(4*Ai(ij))*(le(ij+u_right)*de(ij+u_right)*u(ij+u_right,l)**2 +    &
690                                     le(ij+u_rup)*de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
691                                     le(ij+u_lup)*de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
692                                     le(ij+u_left)*de(ij+u_left)*u(ij+u_left,l)**2 +       &
693                                     le(ij+u_ldown)*de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
694                                     le(ij+u_rdown)*de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) 
695                ! from now on pk contains the vertically-averaged geopotential
696                pk(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1))
697          ENDDO
698       ENDDO
699
700    ELSE ! compressible
701
702       DO l=ll_begin,ll_end
703!DIR$ SIMD
704            DO ij=ij_begin,ij_end         
705               
706                berni(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1)) &
707                     + 1/(4*Ai(ij))*(le(ij+u_right)*de(ij+u_right)*u(ij+u_right,l)**2 +    &
708                                     le(ij+u_rup)*de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
709                                     le(ij+u_lup)*de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
710                                     le(ij+u_left)*de(ij+u_left)*u(ij+u_left,l)**2 +       &
711                                     le(ij+u_ldown)*de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
712                                     le(ij+u_rdown)*de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) 
713          ENDDO
714       ENDDO
715
716    END IF ! Boussinesq/compressible
717
718!!! Add gradients of Bernoulli and Exner functions to du
719    DO l=ll_begin,ll_end
720!DIR$ SIMD
721       DO ij=ij_begin,ij_end         
722       
723             du(ij+u_right,l) = du(ij+u_right,l) + 1/de(ij+u_right) * (             &
724                               0.5*(theta(ij,l)+theta(ij+t_right,l))                &
725                                  *( ne_right*pk(ij,l)+ne_left*pk(ij+t_right,l))    &
726                                   + ne_right*berni(ij,l)+ne_left*berni(ij+t_right,l) )
727       
728       
729             du(ij+u_lup,l) = du(ij+u_lup,l) +  1/de(ij+u_lup) * (                  &
730                               0.5*(theta(ij,l)+theta(ij+t_lup,l))                  &
731                                  *( ne_lup*pk(ij,l)+ne_rdown*pk(ij+t_lup,l))       &
732                                   + ne_lup*berni(ij,l)+ne_rdown*berni(ij+t_lup,l) )
733               
734             du(ij+u_ldown,l) = du(ij+u_ldown,l) + 1/de(ij+u_ldown) * (             &
735                               0.5*(theta(ij,l)+theta(ij+t_ldown,l))                &
736                                  *( ne_ldown*pk(ij,l)+ne_rup*pk(ij+t_ldown,l))     &
737                                   + ne_ldown*berni(ij,l)+ne_rup*berni(ij+t_ldown,l) )
738
739       ENDDO
740    ENDDO
741 
742    CALL trace_end("compute_caldyn_horiz")
743
744END SUBROUTINE compute_caldyn_horiz
745
746SUBROUTINE compute_caldyn_vert(u,theta,rhodz,convm, wflux,wwuu, dps,dtheta_rhodz,du)
747  USE icosa
748  USE disvert_mod
749  USE exner_mod
750  USE trace
751  USE omp_para
752  IMPLICIT NONE
753    REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)
754    REAL(rstd),INTENT(IN)  :: theta(iim*jjm,llm)
755    REAL(rstd),INTENT(IN)  :: rhodz(iim*jjm,llm)
756    REAL(rstd),INTENT(INOUT)  :: convm(iim*jjm,llm)  ! mass flux convergence
757
758    REAL(rstd),INTENT(OUT) :: wflux(iim*jjm,llm+1) ! vertical mass flux (kg/m2/s)
759    REAL(rstd),INTENT(OUT) :: wwuu(iim*3*jjm,llm+1)
760    REAL(rstd),INTENT(OUT) :: du(iim*3*jjm,llm)
761    REAL(rstd),INTENT(OUT) :: dtheta_rhodz(iim*jjm,llm)
762    REAL(rstd),INTENT(OUT) :: dps(iim*jjm)
763
764! temporary variable   
765    INTEGER :: i,j,ij,l
766    REAL(rstd) :: p_ik, exner_ik
767
768
769  CALL trace_start("compute_caldyn_vert")
770
771!$OMP BARRIER   
772!!! cumulate mass flux convergence from top to bottom
773  DO  l = llm-1, 1, -1
774    IF (caldyn_conserv==energy) CALL test_message(req_qu) 
775!$OMP DO SCHEDULE(STATIC)
776!DIR$ SIMD
777    DO ij=ij_begin,ij_end         
778        convm(ij,l) = convm(ij,l) + convm(ij,l+1)
779    ENDDO
780  ENDDO
781
782! IMPLICIT FLUSH on convm
783!!!!!!!!!!!!!!!!!!!!!!!!!
784
785! compute dps
786  IF (omp_first) THEN
787!DIR$ SIMD
788     DO ij=ij_begin,ij_end         
789        ! dps/dt = -int(div flux)dz
790        dps(ij) = convm(ij,1) * g 
791    ENDDO
792  ENDIF
793 
794!!! Compute vertical mass flux (l=1,llm+1 done by caldyn_BC)
795  DO l=ll_beginp1,ll_end
796    IF (caldyn_conserv==energy) CALL test_message(req_qu) 
797!DIR$ SIMD
798      DO ij=ij_begin,ij_end         
799        ! w = int(z,ztop,div(flux)dz) + B(eta)dps/dt
800        ! => w>0 for upward transport
801        wflux( ij, l ) = bp(l) * convm( ij, 1 ) - convm( ij, l ) 
802    ENDDO
803  ENDDO
804
805  DO l=ll_begin,ll_endm1
806!DIR$ SIMD
807     DO ij=ij_begin,ij_end         
808        dtheta_rhodz(ij, l ) = dtheta_rhodz(ij, l )  -   0.5 * (  wflux(ij,l+1) * (theta(ij,l) + theta(ij,l+1)))
809    ENDDO
810  ENDDO
811 
812  DO l=ll_beginp1,ll_end
813!DIR$ SIMD
814      DO ij=ij_begin,ij_end         
815        dtheta_rhodz(ij, l ) = dtheta_rhodz(ij, l )  +   0.5 * ( wflux(ij,l  ) * (theta(ij,l-1) + theta(ij,l) ) )
816    ENDDO
817  ENDDO
818 
819! Compute vertical transport
820  DO l=ll_beginp1,ll_end
821!DIR$ SIMD
822      DO ij=ij_begin,ij_end         
823        wwuu(ij+u_right,l) = 0.5*( wflux(ij,l) + wflux(ij+t_right,l)) * (u(ij+u_right,l) - u(ij+u_right,l-1))
824        wwuu(ij+u_lup,l) = 0.5* ( wflux(ij,l) + wflux(ij+t_lup,l)) * (u(ij+u_lup,l) - u(ij+u_lup,l-1))
825        wwuu(ij+u_ldown,l) = 0.5*( wflux(ij,l) + wflux(ij+t_ldown,l)) * (u(ij+u_ldown,l) - u(ij+u_ldown,l-1))
826     ENDDO
827   ENDDO
828
829 !--> flush wwuu 
830 !$OMP BARRIER
831
832! Add vertical transport to du
833  DO l=ll_begin,ll_end
834!DIR$ SIMD
835     DO ij=ij_begin,ij_end         
836        du(ij+u_right, l )   = du(ij+u_right,l)  - (wwuu(ij+u_right,l+1)+ wwuu(ij+u_right,l)) / (rhodz(ij,l)+rhodz(ij+t_right,l))
837        du(ij+u_lup, l )     = du(ij+u_lup,l)    - (wwuu(ij+u_lup,l+1)  + wwuu(ij+u_lup,l))   / (rhodz(ij,l)+rhodz(ij+t_lup,l))
838        du(ij+u_ldown, l )   = du(ij+u_ldown,l)  - (wwuu(ij+u_ldown,l+1)+ wwuu(ij+u_ldown,l)) / (rhodz(ij,l)+rhodz(ij+t_ldown,l))
839    ENDDO
840  ENDDO     
841 
842  CALL trace_end("compute_caldyn_vert")
843
844  END SUBROUTINE compute_caldyn_vert
845
846!-------------------------------- Diagnostics ----------------------------
847
848  SUBROUTINE check_mass_conservation(f_ps,f_dps)
849  USE icosa
850  USE mpipara
851  IMPLICIT NONE
852    TYPE(t_field),POINTER :: f_ps(:)
853    TYPE(t_field),POINTER :: f_dps(:)
854    REAL(rstd),POINTER :: ps(:)
855    REAL(rstd),POINTER :: dps(:)
856    REAL(rstd) :: mass_tot,dmass_tot
857    INTEGER :: ind,i,j,ij
858   
859    mass_tot=0
860    dmass_tot=0
861   
862    CALL transfert_request(f_dps,req_i1)
863    CALL transfert_request(f_ps,req_i1)
864
865    DO ind=1,ndomain
866      CALL swap_dimensions(ind)
867      CALL swap_geometry(ind)
868
869      ps=f_ps(ind)
870      dps=f_dps(ind)
871
872      DO j=jj_begin,jj_end
873        DO i=ii_begin,ii_end
874          ij=(j-1)*iim+i
875          IF (domain(ind)%own(i,j)) THEN
876            mass_tot=mass_tot+ps(ij)*Ai(ij)/g
877            dmass_tot=dmass_tot+dps(ij)*Ai(ij)/g
878          ENDIF
879        ENDDO
880      ENDDO
881   
882    ENDDO
883    IF (is_mpi_root) PRINT*, "mass_tot ", mass_tot,"      dmass_tot ",dmass_tot       
884
885  END SUBROUTINE check_mass_conservation 
886 
887  SUBROUTINE write_output_fields(f_ps, f_phis, f_dps, f_u, f_theta_rhodz, f_q, &
888       f_buf_i, f_buf_v, f_buf_i3, f_buf1_i, f_buf2_i, f_buf_s, f_buf_p)
889    USE icosa
890    USE vorticity_mod
891    USE theta2theta_rhodz_mod
892    USE pression_mod
893    USE omega_mod
894    USE write_field
895    USE vertical_interp_mod
896    USE wind_mod
897    TYPE(t_field),POINTER :: f_ps(:), f_phis(:), f_u(:), f_theta_rhodz(:), f_q(:), f_dps(:), &
898         f_buf_i(:), f_buf_v(:), f_buf_i3(:), f_buf1_i(:), f_buf2_i(:), f_buf_s(:), f_buf_p(:)
899   
900    REAL(rstd) :: out_pression_lev
901    CHARACTER(LEN=255) :: str_pression
902    CHARACTER(LEN=255) :: physics_type
903   
904    out_pression_level=0
905    CALL getin("out_pression_level",out_pression_level) 
906    WRITE(str_pression,*) INT(out_pression_level/100)
907    str_pression=ADJUSTL(str_pression)
908   
909    CALL writefield("ps",f_ps)
910    CALL writefield("dps",f_dps)
911    CALL writefield("phis",f_phis)
912    CALL vorticity(f_u,f_buf_v)
913    CALL writefield("vort",f_buf_v)
914
915    CALL w_omega(f_ps, f_u, f_buf_i)
916    CALL writefield('omega', f_buf_i)
917    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
918      CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
919      CALL writefield("omega"//TRIM(str_pression),f_buf_s)
920    ENDIF
921   
922    ! Temperature
923!    CALL theta_rhodz2temperature(f_ps,f_theta_rhodz,f_buf_i) ; ! FIXME
924     
925    CALL getin('physics',physics_type)
926    IF (TRIM(physics_type)=='dcmip') THEN
927      CALL Tv2T(f_buf_i,f_q,f_buf1_i) 
928      CALL writefield("T",f_buf1_i)
929      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
930        CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
931        CALL writefield("T"//TRIM(str_pression),f_buf_s)
932      ENDIF
933    ELSE
934      CALL writefield("T",f_buf_i)
935      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
936        CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
937        CALL writefield("T"//TRIM(str_pression),f_buf_s)
938      ENDIF
939    ENDIF
940   
941    ! velocity components
942    CALL un2ulonlat(f_u, f_buf1_i, f_buf2_i)
943    CALL writefield("ulon",f_buf1_i)
944    CALL writefield("ulat",f_buf2_i)
945
946    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
947      CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
948      CALL writefield("ulon"//TRIM(str_pression),f_buf_s)
949      CALL vertical_interp(f_ps,f_buf2_i,f_buf_s,out_pression_level)
950      CALL writefield("ulat"//TRIM(str_pression),f_buf_s)
951    ENDIF
952   
953    ! geopotential ! FIXME
954    CALL thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_buf_s,f_buf_p,f_buf1_i,f_buf2_i,f_buf_i)
955    CALL writefield("p",f_buf_p)
956    CALL writefield("phi",f_geopot)   ! geopotential
957    CALL writefield("theta",f_buf1_i) ! potential temperature
958    CALL writefield("pk",f_buf2_i)    ! Exner pressure
959 
960  END SUBROUTINE write_output_fields
961 
962  SUBROUTINE thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_pks,f_p,f_theta,f_pk,f_phi) 
963    USE field_mod
964    USE pression_mod
965    USE exner_mod
966    USE geopotential_mod
967    USE theta2theta_rhodz_mod
968    TYPE(t_field), POINTER :: f_ps(:), f_phis(:), f_theta_rhodz(:), &  ! IN
969         f_pks(:), f_p(:), f_theta(:), f_pk(:), f_phi(:)               ! OUT
970    REAL(rstd),POINTER :: pk(:,:), p(:,:), theta(:,:), theta_rhodz(:,:), &
971         phi(:,:), phis(:), ps(:), pks(:)
972    INTEGER :: ind
973
974    DO ind=1,ndomain
975       CALL swap_dimensions(ind)
976       CALL swap_geometry(ind)
977       ps = f_ps(ind)
978       p  = f_p(ind)
979       CALL compute_pression(ps,p,0)
980       pk = f_pk(ind)
981       pks = f_pks(ind)
982       CALL compute_exner(ps,p,pks,pk,0)
983       theta_rhodz = f_theta_rhodz(ind)
984       theta = f_theta(ind)
985       CALL compute_theta_rhodz2theta(ps, theta_rhodz,theta,0)
986       phis = f_phis(ind)
987       phi = f_phi(ind)
988       CALL compute_geopotential(phis,pks,pk,theta,phi,0)
989    END DO
990
991  END SUBROUTINE thetarhodz2geopot
992 
993  SUBROUTINE Tv2T(f_Tv, f_q, f_T)
994  USE icosa
995  IMPLICIT NONE
996    TYPE(t_field), POINTER :: f_TV(:)
997    TYPE(t_field), POINTER :: f_q(:)
998    TYPE(t_field), POINTER :: f_T(:)
999   
1000    REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:)
1001    INTEGER :: ind
1002   
1003    DO ind=1,ndomain
1004       CALL swap_dimensions(ind)
1005       CALL swap_geometry(ind)
1006       Tv=f_Tv(ind)
1007       q=f_q(ind)
1008       T=f_T(ind)
1009       T=Tv/(1+0.608*q(:,:,1))
1010    END DO
1011   
1012  END SUBROUTINE Tv2T
1013 
1014END MODULE caldyn_gcm_mod
Note: See TracBrowser for help on using the repository browser.