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

Last change on this file since 295 was 295, checked in by ymipsl, 10 years ago

Merging OpenMP parallisme mode : by subdomain and on vertical level.
This feature is actually experimental but may be retro-compatible with the last method based only on subdomain

YM

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