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

Last change on this file since 179 was 179, checked in by mtort, 11 years ago

Prepared structure for non-traditional shallow-water equations

File size: 36.8 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,u, 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          CALL compute_pvort(ps,u,theta_rhodz, mass,theta,qu,qv)
248          pk = f_pk(ind)
249          geopot = f_geopot(ind) 
250          CALL compute_geopot(ps,mass,theta,u, pk,geopot)
251          hflux=f_hflux(ind)
252          convm = f_dmass(ind)
253          dtheta_rhodz=f_dtheta_rhodz(ind)
254          du=f_du(ind)
255          CALL compute_caldyn_horiz(u,mass,qu,theta,pk,geopot, hflux,convm,dtheta_rhodz,du)
256          IF(caldyn_eta==eta_mass) THEN
257             wflux=f_wflux(ind)
258             wwuu=f_wwuu(ind)
259             dps=f_dps(ind)
260             CALL compute_caldyn_vert(u,theta,mass,convm, wflux,wwuu, dps, dtheta_rhodz, du)
261          END IF
262       ENDDO
263       
264    CASE DEFAULT
265       STOP
266    END SELECT
267
268!$OMP BARRIER
269    IF (write_out) THEN
270
271       IF (is_mpi_root) PRINT *,'CALL write_output_fields'
272
273! ---> for openMP test to fix later
274!       CALL write_output_fields(f_ps, f_phis, f_dps, f_u, f_theta_rhodz, f_q, &
275!            f_buf_i, f_buf_v, f_buf_u3d, f_buf_ulon, f_buf_ulat, f_buf_s, f_buf_p)
276
277       CALL output_field("ps",f_ps)
278       CALL output_field("dps",f_dps)
279       CALL output_field("mass",f_mass)
280       CALL output_field("dmass",f_dmass)
281       CALL output_field("vort",f_qv)
282       CALL output_field("theta",f_theta)
283       CALL output_field("exner",f_pk)
284       CALL output_field("pv",f_qv)
285 
286    END IF
287   
288    !    CALL check_mass_conservation(f_ps,f_dps)
289    CALL trace_end("caldyn")
290!$OMP BARRIER
291   
292END SUBROUTINE caldyn
293
294SUBROUTINE compute_planetvel(planetvel)
295  USE wind_mod
296  REAL(rstd),INTENT(OUT)  :: planetvel(iim*3*jjm)
297  REAL(rstd) :: ulon(iim*3*jjm)
298  REAL(rstd) :: ulat(iim*3*jjm)
299  REAL(rstd) :: lon,lat
300  INTEGER :: ij
301  DO ij=ij_begin_ext,ij_end_ext
302     CALL xyz2lonlat(xyz_e(ij+u_right,:),lon,lat)
303     ulon(ij+u_right)=a*omega*cos(lat)
304     ulat(ij+u_right)=0
305
306     CALL xyz2lonlat(xyz_e(ij+u_lup,:),lon,lat)
307     ulon(ij+u_lup)=a*omega*cos(lat)
308     ulat(ij+u_lup)=0
309
310     CALL xyz2lonlat(xyz_e(ij+u_ldown,:),lon,lat)
311     ulon(ij+u_ldown)=a*omega*cos(lat)
312     ulat(ij+u_ldown)=0
313  END DO
314  CALL compute_wind2D_perp_from_lonlat_compound(ulon, ulat, planetvel)
315END SUBROUTINE compute_planetvel
316   
317SUBROUTINE compute_pvort(ps,u,theta_rhodz, rhodz,theta,qu,qv)
318  USE icosa
319  USE disvert_mod, ONLY : mass_dak, mass_dbk, caldyn_eta, eta_mass
320  USE exner_mod
321  USE trace
322  USE omp_para
323  IMPLICIT NONE
324  REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)
325  REAL(rstd),INTENT(IN)  :: ps(iim*jjm)
326  REAL(rstd),INTENT(IN)  :: theta_rhodz(iim*jjm,llm)
327  REAL(rstd),INTENT(INOUT) :: rhodz(iim*jjm,llm)
328  REAL(rstd),INTENT(OUT) :: theta(iim*jjm,llm)
329  REAL(rstd),INTENT(OUT) :: qu(iim*3*jjm,llm)
330  REAL(rstd),INTENT(OUT) :: qv(iim*2*jjm,llm)
331 
332  INTEGER :: i,j,ij,l
333  REAL(rstd) :: etav,hv, m
334!  REAL(rstd) :: qv(2*iim*jjm,llm)     ! potential velocity 
335 
336  CALL trace_start("compute_pvort") 
337
338  IF(caldyn_eta==eta_mass) THEN
339     CALL wait_message(req_ps) 
340  ELSE
341     CALL wait_message(req_mass)
342  END IF
343  CALL wait_message(req_theta_rhodz) 
344
345  IF(caldyn_eta==eta_mass) THEN ! Compute mass & theta
346     DO l = ll_begin,ll_end
347        CALL test_message(req_u) 
348!DIR$ SIMD
349        DO ij=ij_begin_ext,ij_end_ext
350           m = ( mass_dak(l)+ps(ij)*mass_dbk(l) )/g
351           rhodz(ij,l) = m
352           theta(ij,l) = theta_rhodz(ij,l)/rhodz(ij,l)
353        ENDDO
354     ENDDO
355  ELSE ! Compute only theta
356     DO l = ll_begin,ll_end
357        CALL test_message(req_u) 
358!DIR$ SIMD
359       DO ij=ij_begin_ext,ij_end_ext
360         theta(ij,l) = theta_rhodz(ij,l)/rhodz(ij,l)
361       ENDDO
362     ENDDO
363  END IF
364
365  CALL wait_message(req_u)   
366 
367!!! Compute shallow-water potential vorticity
368  DO l = ll_begin,ll_end
369!DIR$ SIMD
370     DO ij=ij_begin_ext,ij_end_ext
371          etav= 1./Av(ij+z_up)*(  ne_rup        * u(ij+u_rup,l)        * de(ij+u_rup)         &
372                                + ne_left * u(ij+t_rup+u_left,l) * de(ij+t_rup+u_left)  &
373                                - ne_lup        * u(ij+u_lup,l)        * de(ij+u_lup) )                               
374
375          hv =  Riv2(ij,vup)          * rhodz(ij,l)            &
376              + Riv2(ij+t_rup,vldown) * rhodz(ij+t_rup,l)     &
377              + Riv2(ij+t_lup,vrdown) * rhodz(ij+t_lup,l)
378     
379          qv(ij+z_up,l) = ( etav+fv(ij+z_up) )/hv
380         
381          etav = 1./Av(ij+z_down)*(  ne_ldown         * u(ij+u_ldown,l)          * de(ij+u_ldown)          &
382                                   + ne_right * u(ij+t_ldown+u_right,l)  * de(ij+t_ldown+u_right)  &
383                                   - ne_rdown         * u(ij+u_rdown,l)          * de(ij+u_rdown) )
384       
385          hv =  Riv2(ij,vdown)        * rhodz(ij,l)          &
386              + Riv2(ij+t_ldown,vrup) * rhodz(ij+t_ldown,l)  &
387              + Riv2(ij+t_rdown,vlup) * rhodz(ij+t_rdown,l)
388                       
389          qv(ij+z_down,l) =( etav+fv(ij+z_down) )/hv
390
391      ENDDO
392
393!DIR$ SIMD
394      DO ij=ij_begin,ij_end
395         qu(ij+u_right,l) = 0.5*(qv(ij+z_rdown,l)+qv(ij+z_rup,l)) 
396         qu(ij+u_lup,l) = 0.5*(qv(ij+z_up,l)+qv(ij+z_lup,l)) 
397         qu(ij+u_ldown,l) = 0.5*(qv(ij+z_ldown,l)+qv(ij+z_down,l)) 
398      END DO
399     
400   ENDDO
401
402   CALL trace_end("compute_pvort")
403                                                       
404  END SUBROUTINE compute_pvort
405 
406  SUBROUTINE compute_geopot(ps,rhodz,theta,u, pk,geopot) 
407  USE icosa
408  USE disvert_mod
409  USE exner_mod
410  USE trace
411  USE omp_para
412  IMPLICIT NONE
413    REAL(rstd),INTENT(INOUT) :: ps(iim*jjm)
414    REAL(rstd),INTENT(IN)    :: rhodz(iim*jjm,llm)
415    REAL(rstd),INTENT(IN)    :: theta(iim*jjm,llm)    ! potential temperature
416    REAL(rstd),INTENT(IN)    :: u(3*iim*jjm,llm)      ! prognostic velocity
417    REAL(rstd),INTENT(INOUT) :: pk(iim*jjm,llm)       ! Exner function
418    REAL(rstd),INTENT(INOUT) :: geopot(iim*jjm,llm+1) ! geopotential
419   
420    REAL(rstd) :: urel(3*iim*jjm,llm)  ! relative velocity
421    REAL(rstd) :: cor_NT(iim*jjm,llm)  ! NT coriolis force u.(du/dPhi)
422    INTEGER :: i,j,ij,l
423    REAL(rstd) :: p_ik, exner_ik
424
425    CALL trace_start("compute_geopot")
426
427    IF(caldyn_eta==eta_mass) THEN
428   
429!!! Compute exner function and geopotential
430       DO l = 1,llm
431          !$OMP DO SCHEDULE(STATIC)
432!DIR$ SIMD
433         DO ij=ij_begin_ext,ij_end_ext         
434            p_ik = ptop + mass_ak(l) + mass_bk(l)*ps(ij) ! FIXME : leave ps for the moment ; change ps to Ms later
435            !         p_ik = ptop + g*(mass_ak(l)+ mass_bk(l)*ps(i,j))
436            exner_ik = cpp * (p_ik/preff) ** kappa
437            pk(ij,l) = exner_ik
438            ! specific volume v = kappa*theta*pi/p = dphi/g/rhodz
439            geopot(ij,l+1) = geopot(ij,l) + (g*kappa)*rhodz(ij,l)*theta(ij,l)*exner_ik/p_ik
440          ENDDO
441       ENDDO
442
443    ELSE 
444       ! We are using a Lagrangian vertical coordinate
445       ! Pressure must be computed first top-down (temporarily stored in pk)
446       ! Then Exner pressure and geopotential are computed bottom-up
447       ! Notice that the computation below should work also when caldyn_eta=eta_mass         
448
449       ! uppermost layer
450!DIR$ SIMD
451       DO ij=ij_begin_ext,ij_end_ext         
452         pk(ij,llm) = ptop + (.5*g)*rhodz(ij,llm)
453       END DO
454       ! other layers
455       DO l = llm-1, 1, -1
456          !$OMP DO SCHEDULE(STATIC)
457!DIR$ SIMD
458          DO ij=ij_begin_ext,ij_end_ext         
459             pk(ij,l) = pk(ij,l+1) + (.5*g)*(rhodz(ij,l)+rhodz(ij,l+1))
460          END DO
461       END DO
462       ! surface pressure (for diagnostics)
463       DO ij=ij_begin_ext,ij_end_ext         
464         ps(ij) = pk(ij,1) + (.5*g)*rhodz(ij,1)
465       END DO
466
467       IF(boussinesq) THEN ! specific volume 1 = dphi/g/rhodz
468          DO l = 1,llm
469             !$OMP DO SCHEDULE(STATIC)
470!DIR$ SIMD
471             DO ij=ij_begin_ext,ij_end_ext         
472               geopot(ij,l+1) = geopot(ij,l) + g*rhodz(ij,l)
473             ENDDO
474          ENDDO
475       ELSE ! 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)
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) :: urel(3*iim*jjm,llm)  ! relative velocity
514    REAL(rstd) :: Ftheta(3*iim*jjm,llm) ! theta flux
515    REAL(rstd) :: berni(iim*jjm,llm)  ! Bernoulli function
516   
517    INTEGER :: i,j,ij,l
518    REAL(rstd) :: ww,uu
519
520    CALL trace_start("compute_caldyn_horiz")
521
522    CALL wait_message(req_theta_rhodz) 
523
524  DO l = ll_begin, ll_end
525!!!  Compute mass and theta fluxes
526    IF (caldyn_conserv==energy) CALL test_message(req_qu) 
527!DIR$ SIMD
528    DO ij=ij_begin_ext,ij_end_ext         
529        hflux(ij+u_right,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_right,l))*u(ij+u_right,l)*le(ij+u_right)
530        hflux(ij+u_lup,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_lup,l))*u(ij+u_lup,l)*le(ij+u_lup)
531        hflux(ij+u_ldown,l)=0.5*(rhodz(ij,l)+rhodz(ij+t_ldown,l))*u(ij+u_ldown,l)*le(ij+u_ldown)
532
533        Ftheta(ij+u_right,l)=0.5*(theta(ij,l)+theta(ij+t_right,l))*hflux(ij+u_right,l)
534        Ftheta(ij+u_lup,l)=0.5*(theta(ij,l)+theta(ij+t_lup,l))*hflux(ij+u_lup,l)
535        Ftheta(ij+u_ldown,l)=0.5*(theta(ij,l)+theta(ij+t_ldown,l))*hflux(ij+u_ldown,l)
536    ENDDO
537
538!!! compute horizontal divergence of fluxes
539!DIR$ SIMD
540    DO ij=ij_begin,ij_end         
541        ! convm = -div(mass flux), sign convention as in Ringler et al. 2012, eq. 21
542        convm(ij,l)= -1./Ai(ij)*(ne_right*hflux(ij+u_right,l)   +  &
543                                ne_rup*hflux(ij+u_rup,l)       +  & 
544                                ne_lup*hflux(ij+u_lup,l)       +  & 
545                                ne_left*hflux(ij+u_left,l)     +  & 
546                                ne_ldown*hflux(ij+u_ldown,l)   +  & 
547                                ne_rdown*hflux(ij+u_rdown,l))   
548
549        ! signe ? attention d (rho theta dz)
550        ! dtheta_rhodz =  -div(flux.theta)
551        dtheta_rhodz(ij,l)=-1./Ai(ij)*(ne_right*Ftheta(ij+u_right,l)    +  &
552                                 ne_rup*Ftheta(ij+u_rup,l)        +  & 
553                                 ne_lup*Ftheta(ij+u_lup,l)        +  & 
554                                 ne_left*Ftheta(ij+u_left,l)      +  & 
555                                 ne_ldown*Ftheta(ij+u_ldown,l)    +  & 
556                                 ne_rdown*Ftheta(ij+u_rdown,l))
557    ENDDO
558
559 END DO
560
561!!! Compute potential vorticity (Coriolis) contribution to du
562 
563    SELECT CASE(caldyn_conserv)
564    CASE(energy) ! energy-conserving TRiSK
565
566       CALL wait_message(req_qu)
567
568        DO l=ll_begin,ll_end
569!DIR$ SIMD
570          DO ij=ij_begin,ij_end         
571   
572                uu = wee(ij+u_right,1,1)*hflux(ij+u_rup,l)*(qu(ij+u_right,l)+qu(ij+u_rup,l))+                             &
573                     wee(ij+u_right,2,1)*hflux(ij+u_lup,l)*(qu(ij+u_right,l)+qu(ij+u_lup,l))+                             &
574                     wee(ij+u_right,3,1)*hflux(ij+u_left,l)*(qu(ij+u_right,l)+qu(ij+u_left,l))+                           &
575                     wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)*(qu(ij+u_right,l)+qu(ij+u_ldown,l))+                         &
576                     wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)*(qu(ij+u_right,l)+qu(ij+u_rdown,l))+                         & 
577                     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))+         &
578                     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))+         &
579                     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))+         &
580                     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))+             &
581                     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))
582                du(ij+u_right,l) = .5*uu/de(ij+u_right)
583               
584                uu = wee(ij+u_lup,1,1)*hflux(ij+u_left,l)*(qu(ij+u_lup,l)+qu(ij+u_left,l)) +                        &
585                     wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)*(qu(ij+u_lup,l)+qu(ij+u_ldown,l)) +                      &
586                     wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)*(qu(ij+u_lup,l)+qu(ij+u_rdown,l)) +                      &
587                     wee(ij+u_lup,4,1)*hflux(ij+u_right,l)*(qu(ij+u_lup,l)+qu(ij+u_right,l)) +                      &
588                     wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)*(qu(ij+u_lup,l)+qu(ij+u_rup,l)) +                          & 
589                     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)) +          &
590                     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)) +              &
591                     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)) +              &
592                     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)) +            &
593                     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))
594                du(ij+u_lup,l) = .5*uu/de(ij+u_lup)
595
596               
597                uu = wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)*(qu(ij+u_ldown,l)+qu(ij+u_rdown,l)) +                      &
598                     wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)*(qu(ij+u_ldown,l)+qu(ij+u_right,l)) +                      &
599                     wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)*(qu(ij+u_ldown,l)+qu(ij+u_rup,l)) +                          &
600                     wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)*(qu(ij+u_ldown,l)+qu(ij+u_lup,l)) +                          &
601                     wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)*(qu(ij+u_ldown,l)+qu(ij+u_left,l)) +                        & 
602                     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)) +          &
603                     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)) +        &
604                     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)) +      &
605                     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)) +      &
606                     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))
607                du(ij+u_ldown,l) = .5*uu/de(ij+u_ldown)
608               
609          ENDDO
610       ENDDO
611       
612    CASE(enstrophy) ! enstrophy-conserving TRiSK
613 
614        DO l=ll_begin,ll_end
615!DIR$ SIMD
616          DO ij=ij_begin,ij_end         
617
618                uu = wee(ij+u_right,1,1)*hflux(ij+u_rup,l)+                           &
619                     wee(ij+u_right,2,1)*hflux(ij+u_lup,l)+                           &
620                     wee(ij+u_right,3,1)*hflux(ij+u_left,l)+                          &
621                     wee(ij+u_right,4,1)*hflux(ij+u_ldown,l)+                         &
622                     wee(ij+u_right,5,1)*hflux(ij+u_rdown,l)+                         & 
623                     wee(ij+u_right,1,2)*hflux(ij+t_right+u_ldown,l)+                 &
624                     wee(ij+u_right,2,2)*hflux(ij+t_right+u_rdown,l)+                 &
625                     wee(ij+u_right,3,2)*hflux(ij+t_right+u_right,l)+                 &
626                     wee(ij+u_right,4,2)*hflux(ij+t_right+u_rup,l)+                   &
627                     wee(ij+u_right,5,2)*hflux(ij+t_right+u_lup,l)
628                du(ij+u_right,l) = qu(ij+u_right,l)*uu/de(ij+u_right)
629               
630               
631                uu = wee(ij+u_lup,1,1)*hflux(ij+u_left,l)+                        &
632                     wee(ij+u_lup,2,1)*hflux(ij+u_ldown,l)+                       &
633                     wee(ij+u_lup,3,1)*hflux(ij+u_rdown,l)+                       &
634                     wee(ij+u_lup,4,1)*hflux(ij+u_right,l)+                       &
635                     wee(ij+u_lup,5,1)*hflux(ij+u_rup,l)+                         & 
636                     wee(ij+u_lup,1,2)*hflux(ij+t_lup+u_right,l)+                 &
637                     wee(ij+u_lup,2,2)*hflux(ij+t_lup+u_rup,l)+                   &
638                     wee(ij+u_lup,3,2)*hflux(ij+t_lup+u_lup,l)+                   &
639                     wee(ij+u_lup,4,2)*hflux(ij+t_lup+u_left,l)+                  &
640                     wee(ij+u_lup,5,2)*hflux(ij+t_lup+u_ldown,l)
641                du(ij+u_lup,l) = qu(ij+u_lup,l)*uu/de(ij+u_lup)
642
643                uu = wee(ij+u_ldown,1,1)*hflux(ij+u_rdown,l)+                      &
644                     wee(ij+u_ldown,2,1)*hflux(ij+u_right,l)+                      &
645                     wee(ij+u_ldown,3,1)*hflux(ij+u_rup,l)+                        &
646                     wee(ij+u_ldown,4,1)*hflux(ij+u_lup,l)+                        &
647                     wee(ij+u_ldown,5,1)*hflux(ij+u_left,l)+                       & 
648                     wee(ij+u_ldown,1,2)*hflux(ij+t_ldown+u_lup,l)+                &
649                     wee(ij+u_ldown,2,2)*hflux(ij+t_ldown+u_left,l)+               &
650                     wee(ij+u_ldown,3,2)*hflux(ij+t_ldown+u_ldown,l)+              &
651                     wee(ij+u_ldown,4,2)*hflux(ij+t_ldown+u_rdown,l)+              &
652                     wee(ij+u_ldown,5,2)*hflux(ij+t_ldown+u_right,l)
653                du(ij+u_ldown,l) = qu(ij+u_ldown,l)*uu/de(ij+u_ldown)
654     
655          ENDDO
656       ENDDO
657       
658    CASE DEFAULT
659       STOP
660    END SELECT
661 
662!!! Compute bernouilli term = Kinetic Energy + geopotential
663    IF(boussinesq) THEN
664
665       DO l=ll_begin,ll_end
666!DIR$ SIMD
667          DO ij=ij_begin,ij_end         
668               
669                berni(ij,l) = pk(ij,l) &
670                     + 1/(4*Ai(ij))*(le(ij+u_right)*de(ij+u_right)*u(ij+u_right,l)**2 +    &
671                                     le(ij+u_rup)*de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
672                                     le(ij+u_lup)*de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
673                                     le(ij+u_left)*de(ij+u_left)*u(ij+u_left,l)**2 +       &
674                                     le(ij+u_ldown)*de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
675                                     le(ij+u_rdown)*de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) 
676                pk(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1))
677          ENDDO
678       ENDDO
679
680    ELSE
681
682       DO l=ll_begin,ll_end
683!DIR$ SIMD
684            DO ij=ij_begin,ij_end         
685               
686                berni(ij,l) = .5*(geopot(ij,l)+geopot(ij,l+1)) &
687                     + 1/(4*Ai(ij))*(le(ij+u_right)*de(ij+u_right)*u(ij+u_right,l)**2 +    &
688                                     le(ij+u_rup)*de(ij+u_rup)*u(ij+u_rup,l)**2 +          &
689                                     le(ij+u_lup)*de(ij+u_lup)*u(ij+u_lup,l)**2 +          &
690                                     le(ij+u_left)*de(ij+u_left)*u(ij+u_left,l)**2 +       &
691                                     le(ij+u_ldown)*de(ij+u_ldown)*u(ij+u_ldown,l)**2 +    &
692                                     le(ij+u_rdown)*de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) 
693          ENDDO
694       ENDDO
695
696    END IF ! Boussinesq/compressible
697
698!!! Add gradients of Bernoulli and Exner functions to du
699    DO l=ll_begin,ll_end
700!DIR$ SIMD
701       DO ij=ij_begin,ij_end         
702       
703             du(ij+u_right,l) = du(ij+u_right,l) + 1/de(ij+u_right) * (             &
704                               0.5*(theta(ij,l)+theta(ij+t_right,l))                &
705                                  *( ne_right*pk(ij,l)+ne_left*pk(ij+t_right,l))    &
706                                   + ne_right*berni(ij,l)+ne_left*berni(ij+t_right,l) )
707       
708       
709             du(ij+u_lup,l) = du(ij+u_lup,l) +  1/de(ij+u_lup) * (                  &
710                               0.5*(theta(ij,l)+theta(ij+t_lup,l))                  &
711                                  *( ne_lup*pk(ij,l)+ne_rdown*pk(ij+t_lup,l))       &
712                                   + ne_lup*berni(ij,l)+ne_rdown*berni(ij+t_lup,l) )
713               
714             du(ij+u_ldown,l) = du(ij+u_ldown,l) + 1/de(ij+u_ldown) * (             &
715                               0.5*(theta(ij,l)+theta(ij+t_ldown,l))                &
716                                  *( ne_ldown*pk(ij,l)+ne_rup*pk(ij+t_ldown,l))     &
717                                   + ne_ldown*berni(ij,l)+ne_rup*berni(ij+t_ldown,l) )
718
719       ENDDO
720    ENDDO
721 
722    CALL trace_end("compute_caldyn_horiz")
723
724END SUBROUTINE compute_caldyn_horiz
725
726SUBROUTINE compute_caldyn_vert(u,theta,rhodz,convm, wflux,wwuu, dps,dtheta_rhodz,du)
727  USE icosa
728  USE disvert_mod
729  USE exner_mod
730  USE trace
731  USE omp_para
732  IMPLICIT NONE
733    REAL(rstd),INTENT(IN)  :: u(iim*3*jjm,llm)
734    REAL(rstd),INTENT(IN)  :: theta(iim*jjm,llm)
735    REAL(rstd),INTENT(IN)  :: rhodz(iim*jjm,llm)
736    REAL(rstd),INTENT(INOUT)  :: convm(iim*jjm,llm)  ! mass flux convergence
737
738    REAL(rstd),INTENT(OUT) :: wflux(iim*jjm,llm+1) ! vertical mass flux (kg/m2/s)
739    REAL(rstd),INTENT(OUT) :: wwuu(iim*3*jjm,llm+1)
740    REAL(rstd),INTENT(OUT) :: du(iim*3*jjm,llm)
741    REAL(rstd),INTENT(OUT) :: dtheta_rhodz(iim*jjm,llm)
742    REAL(rstd),INTENT(OUT) :: dps(iim*jjm)
743
744! temporary variable   
745    INTEGER :: i,j,ij,l
746    REAL(rstd) :: p_ik, exner_ik
747
748
749  CALL trace_start("compute_caldyn_vert")
750
751!$OMP BARRIER   
752!!! cumulate mass flux convergence from top to bottom
753  DO  l = llm-1, 1, -1
754    IF (caldyn_conserv==energy) CALL test_message(req_qu) 
755!$OMP DO SCHEDULE(STATIC)
756!DIR$ SIMD
757    DO ij=ij_begin,ij_end         
758        convm(ij,l) = convm(ij,l) + convm(ij,l+1)
759    ENDDO
760  ENDDO
761
762! IMPLICIT FLUSH on convm
763!!!!!!!!!!!!!!!!!!!!!!!!!
764
765! compute dps
766  IF (omp_first) THEN
767!DIR$ SIMD
768     DO ij=ij_begin,ij_end         
769        ! dps/dt = -int(div flux)dz
770        dps(ij) = convm(ij,1) * g 
771    ENDDO
772  ENDIF
773 
774!!! Compute vertical mass flux (l=1,llm+1 done by caldyn_BC)
775  DO l=ll_beginp1,ll_end
776    IF (caldyn_conserv==energy) CALL test_message(req_qu) 
777!DIR$ SIMD
778      DO ij=ij_begin,ij_end         
779        ! w = int(z,ztop,div(flux)dz) + B(eta)dps/dt
780        ! => w>0 for upward transport
781        wflux( ij, l ) = bp(l) * convm( ij, 1 ) - convm( ij, l ) 
782    ENDDO
783  ENDDO
784
785  DO l=ll_begin,ll_endm1
786!DIR$ SIMD
787     DO ij=ij_begin,ij_end         
788        dtheta_rhodz(ij, l ) = dtheta_rhodz(ij, l )  -   0.5 * (  wflux(ij,l+1) * (theta(ij,l) + theta(ij,l+1)))
789    ENDDO
790  ENDDO
791 
792  DO l=ll_beginp1,ll_end
793!DIR$ SIMD
794      DO ij=ij_begin,ij_end         
795        dtheta_rhodz(ij, l ) = dtheta_rhodz(ij, l )  +   0.5 * ( wflux(ij,l  ) * (theta(ij,l-1) + theta(ij,l) ) )
796    ENDDO
797  ENDDO
798 
799! Compute vertical transport
800  DO l=ll_beginp1,ll_end
801!DIR$ SIMD
802      DO ij=ij_begin,ij_end         
803        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))
804        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))
805        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))
806     ENDDO
807   ENDDO
808
809 !--> flush wwuu 
810 !$OMP BARRIER
811
812! Add vertical transport to du
813  DO l=ll_begin,ll_end
814!DIR$ SIMD
815     DO ij=ij_begin,ij_end         
816        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))
817        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))
818        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))
819    ENDDO
820  ENDDO     
821 
822  CALL trace_end("compute_caldyn_vert")
823
824  END SUBROUTINE compute_caldyn_vert
825
826!-------------------------------- Diagnostics ----------------------------
827
828  SUBROUTINE check_mass_conservation(f_ps,f_dps)
829  USE icosa
830  USE mpipara
831  IMPLICIT NONE
832    TYPE(t_field),POINTER :: f_ps(:)
833    TYPE(t_field),POINTER :: f_dps(:)
834    REAL(rstd),POINTER :: ps(:)
835    REAL(rstd),POINTER :: dps(:)
836    REAL(rstd) :: mass_tot,dmass_tot
837    INTEGER :: ind,i,j,ij
838   
839    mass_tot=0
840    dmass_tot=0
841   
842    CALL transfert_request(f_dps,req_i1)
843    CALL transfert_request(f_ps,req_i1)
844
845    DO ind=1,ndomain
846      CALL swap_dimensions(ind)
847      CALL swap_geometry(ind)
848
849      ps=f_ps(ind)
850      dps=f_dps(ind)
851
852      DO j=jj_begin,jj_end
853        DO i=ii_begin,ii_end
854          ij=(j-1)*iim+i
855          IF (domain(ind)%own(i,j)) THEN
856            mass_tot=mass_tot+ps(ij)*Ai(ij)/g
857            dmass_tot=dmass_tot+dps(ij)*Ai(ij)/g
858          ENDIF
859        ENDDO
860      ENDDO
861   
862    ENDDO
863    IF (is_mpi_root) PRINT*, "mass_tot ", mass_tot,"      dmass_tot ",dmass_tot       
864
865  END SUBROUTINE check_mass_conservation 
866 
867  SUBROUTINE write_output_fields(f_ps, f_phis, f_dps, f_u, f_theta_rhodz, f_q, &
868       f_buf_i, f_buf_v, f_buf_i3, f_buf1_i, f_buf2_i, f_buf_s, f_buf_p)
869    USE icosa
870    USE vorticity_mod
871    USE theta2theta_rhodz_mod
872    USE pression_mod
873    USE omega_mod
874    USE write_field
875    USE vertical_interp_mod
876    USE wind_mod
877    TYPE(t_field),POINTER :: f_ps(:), f_phis(:), f_u(:), f_theta_rhodz(:), f_q(:), f_dps(:), &
878         f_buf_i(:), f_buf_v(:), f_buf_i3(:), f_buf1_i(:), f_buf2_i(:), f_buf_s(:), f_buf_p(:)
879   
880    REAL(rstd) :: out_pression_lev
881    CHARACTER(LEN=255) :: str_pression
882    CHARACTER(LEN=255) :: physics_type
883   
884    out_pression_level=0
885    CALL getin("out_pression_level",out_pression_level) 
886    WRITE(str_pression,*) INT(out_pression_level/100)
887    str_pression=ADJUSTL(str_pression)
888   
889    CALL writefield("ps",f_ps)
890    CALL writefield("dps",f_dps)
891    CALL writefield("phis",f_phis)
892    CALL vorticity(f_u,f_buf_v)
893    CALL writefield("vort",f_buf_v)
894
895    CALL w_omega(f_ps, f_u, f_buf_i)
896    CALL writefield('omega', f_buf_i)
897    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
898      CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
899      CALL writefield("omega"//TRIM(str_pression),f_buf_s)
900    ENDIF
901   
902    ! Temperature
903!    CALL theta_rhodz2temperature(f_ps,f_theta_rhodz,f_buf_i) ; ! FIXME
904     
905    CALL getin('physics',physics_type)
906    IF (TRIM(physics_type)=='dcmip') THEN
907      CALL Tv2T(f_buf_i,f_q,f_buf1_i) 
908      CALL writefield("T",f_buf1_i)
909      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
910        CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
911        CALL writefield("T"//TRIM(str_pression),f_buf_s)
912      ENDIF
913    ELSE
914      CALL writefield("T",f_buf_i)
915      IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
916        CALL vertical_interp(f_ps,f_buf_i,f_buf_s,out_pression_level)
917        CALL writefield("T"//TRIM(str_pression),f_buf_s)
918      ENDIF
919    ENDIF
920   
921    ! velocity components
922    CALL un2ulonlat(f_u, f_buf1_i, f_buf2_i)
923    CALL writefield("ulon",f_buf1_i)
924    CALL writefield("ulat",f_buf2_i)
925
926    IF (out_pression_level<=preff .AND. out_pression_level > 0) THEN
927      CALL vertical_interp(f_ps,f_buf1_i,f_buf_s,out_pression_level)
928      CALL writefield("ulon"//TRIM(str_pression),f_buf_s)
929      CALL vertical_interp(f_ps,f_buf2_i,f_buf_s,out_pression_level)
930      CALL writefield("ulat"//TRIM(str_pression),f_buf_s)
931    ENDIF
932   
933    ! geopotential ! FIXME
934    CALL thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_buf_s,f_buf_p,f_buf1_i,f_buf2_i,f_buf_i)
935    CALL writefield("p",f_buf_p)
936    CALL writefield("phi",f_geopot)   ! geopotential
937    CALL writefield("theta",f_buf1_i) ! potential temperature
938    CALL writefield("pk",f_buf2_i)    ! Exner pressure
939 
940  END SUBROUTINE write_output_fields
941 
942  SUBROUTINE thetarhodz2geopot(f_ps,f_phis,f_theta_rhodz, f_pks,f_p,f_theta,f_pk,f_phi) 
943    USE field_mod
944    USE pression_mod
945    USE exner_mod
946    USE geopotential_mod
947    USE theta2theta_rhodz_mod
948    TYPE(t_field), POINTER :: f_ps(:), f_phis(:), f_theta_rhodz(:), &  ! IN
949         f_pks(:), f_p(:), f_theta(:), f_pk(:), f_phi(:)               ! OUT
950    REAL(rstd),POINTER :: pk(:,:), p(:,:), theta(:,:), theta_rhodz(:,:), &
951         phi(:,:), phis(:), ps(:), pks(:)
952    INTEGER :: ind
953
954    DO ind=1,ndomain
955       CALL swap_dimensions(ind)
956       CALL swap_geometry(ind)
957       ps = f_ps(ind)
958       p  = f_p(ind)
959       CALL compute_pression(ps,p,0)
960       pk = f_pk(ind)
961       pks = f_pks(ind)
962       CALL compute_exner(ps,p,pks,pk,0)
963       theta_rhodz = f_theta_rhodz(ind)
964       theta = f_theta(ind)
965       CALL compute_theta_rhodz2theta(ps, theta_rhodz,theta,0)
966       phis = f_phis(ind)
967       phi = f_phi(ind)
968       CALL compute_geopotential(phis,pks,pk,theta,phi,0)
969    END DO
970
971  END SUBROUTINE thetarhodz2geopot
972 
973  SUBROUTINE Tv2T(f_Tv, f_q, f_T)
974  USE icosa
975  IMPLICIT NONE
976    TYPE(t_field), POINTER :: f_TV(:)
977    TYPE(t_field), POINTER :: f_q(:)
978    TYPE(t_field), POINTER :: f_T(:)
979   
980    REAL(rstd),POINTER :: Tv(:,:), q(:,:,:), T(:,:)
981    INTEGER :: ind
982   
983    DO ind=1,ndomain
984       CALL swap_dimensions(ind)
985       CALL swap_geometry(ind)
986       Tv=f_Tv(ind)
987       q=f_q(ind)
988       T=f_T(ind)
989       T=Tv/(1+0.608*q(:,:,1))
990    END DO
991   
992  END SUBROUTINE Tv2T
993 
994END MODULE caldyn_gcm_mod
Note: See TracBrowser for help on using the repository browser.