source: codes/icosagcm/trunk/src/timeloop_gcm.f90 @ 339

Last change on this file since 339 was 339, checked in by dubos, 9 years ago

Kinnmark (1984) RK2.5 scheme

File size: 19.1 KB
Line 
1MODULE timeloop_gcm_mod
2  USE transfert_mod
3  USE icosa
4  PRIVATE
5 
6  PUBLIC :: init_timeloop, timeloop
7
8  INTEGER, PARAMETER :: euler=1, rk4=2, mlf=3, rk25=4
9  INTEGER, PARAMETER :: itau_sync=10
10  REAL(rstd), DIMENSION(4), PARAMETER :: coef_rk4 = (/ .25, 1./3., .5, 1. /)
11  REAL(rstd), DIMENSION(5), PARAMETER :: coef_rk25 = (/ .25, 1./6., 3./8., .5, 1. /)
12
13  TYPE(t_message),SAVE :: req_ps0, req_mass0, req_theta_rhodz0, req_u0, req_q0
14
15  TYPE(t_field),POINTER,SAVE :: f_q(:)
16  TYPE(t_field),POINTER,SAVE :: f_rhodz(:), f_mass(:), f_massm1(:), f_massm2(:), f_dmass(:)
17  TYPE(t_field),POINTER,SAVE :: f_phis(:), f_ps(:),f_psm1(:), f_psm2(:), f_dps(:)
18  TYPE(t_field),POINTER,SAVE :: f_u(:),f_um1(:),f_um2(:), f_du(:)
19  TYPE(t_field),POINTER,SAVE :: f_theta_rhodz(:),f_theta_rhodzm1(:),f_theta_rhodzm2(:), f_dtheta_rhodz(:)
20  TYPE(t_field),POINTER,SAVE :: f_hflux(:), f_wflux(:), f_hfluxt(:), f_wfluxt(:) 
21
22  INTEGER,SAVE :: nb_stage, matsuno_period, scheme   
23!$OMP THREADPRIVATE(nb_stage, matsuno_period, scheme)
24
25CONTAINS
26 
27  SUBROUTINE init_timeloop
28  USE icosa
29  USE dissip_gcm_mod
30  USE caldyn_mod
31  USE etat0_mod
32  USE disvert_mod
33  USE guided_mod
34  USE advect_tracer_mod
35  USE physics_mod
36  USE mpipara
37  USE omp_para
38  USE trace
39  USE transfert_mod
40  USE check_conserve_mod
41  USE output_field_mod
42  USE write_field
43  USE theta2theta_rhodz_mod
44  USE sponge_mod
45  IMPLICIT NONE
46
47    CHARACTER(len=255) :: def
48
49
50   IF (xios_output) itau_out=1
51   IF (.NOT. enable_io) itau_out=HUGE(itau_out)
52
53! Time-independant orography
54    CALL allocate_field(f_phis,field_t,type_real,name='phis')
55! Trends
56    CALL allocate_field(f_du,field_u,type_real,llm,name='du')
57    CALL allocate_field(f_dtheta_rhodz,field_t,type_real,llm,name='dtheta_rhodz')
58! Model state at current time step (RK/MLF/Euler)
59    CALL allocate_field(f_ps,field_t,type_real, name='ps')
60    CALL allocate_field(f_mass,field_t,type_real,llm,name='mass')
61    CALL allocate_field(f_u,field_u,type_real,llm,name='u')
62    CALL allocate_field(f_theta_rhodz,field_t,type_real,llm,name='theta_rhodz')
63! Model state at previous time step (RK/MLF)
64    CALL allocate_field(f_um1,field_u,type_real,llm,name='um1')
65    CALL allocate_field(f_theta_rhodzm1,field_t,type_real,llm,name='theta_rhodzm1')
66! Tracers
67    CALL allocate_field(f_q,field_t,type_real,llm,nqtot,'q')
68    CALL allocate_field(f_rhodz,field_t,type_real,llm,name='rhodz')
69! Mass fluxes
70    CALL allocate_field(f_hflux,field_u,type_real,llm)    ! instantaneous mass fluxes
71    CALL allocate_field(f_hfluxt,field_u,type_real,llm)   ! mass "fluxes" accumulated in time
72    CALL allocate_field(f_wflux,field_t,type_real,llm+1)  ! vertical mass fluxes
73    CALL allocate_field(f_dmass,field_t,type_real,llm, name='dmass')
74
75    IF(caldyn_eta == eta_mass) THEN ! eta = mass coordinate (default)
76       CALL allocate_field(f_dps,field_t,type_real,name='dps')
77       CALL allocate_field(f_psm1,field_t,type_real,name='psm1')
78       CALL allocate_field(f_wfluxt,field_t,type_real,llm+1,name='wfluxt')
79       ! the following are unused but must point to something
80!       f_massm1 => f_mass
81    ELSE ! eta = Lagrangian vertical coordinate
82       CALL allocate_field(f_massm1,field_t,type_real,llm, name='massm1')
83       ! the following are unused but must point to something
84       f_wfluxt => f_wflux
85       f_dps  => f_phis
86       f_psm1 => f_phis
87    END IF
88
89    def='runge_kutta'
90    CALL getin('scheme',def)
91
92    SELECT CASE (TRIM(def))
93      CASE('euler')
94        scheme=euler
95        nb_stage=1
96      CASE ('runge_kutta')
97        scheme=rk4
98        nb_stage=4
99      CASE ('RK2.5')
100        scheme=rk25
101        nb_stage=5
102      CASE ('leapfrog_matsuno')
103        scheme=mlf
104        matsuno_period=5
105        CALL getin('matsuno_period',matsuno_period)
106        nb_stage=matsuno_period+1
107     ! Model state 2 time steps ago (MLF)
108        CALL allocate_field(f_theta_rhodzm2,field_t,type_real,llm)
109        CALL allocate_field(f_um2,field_u,type_real,llm)
110        IF(caldyn_eta == eta_mass) THEN ! eta = mass coordinate (default)
111           CALL allocate_field(f_psm2,field_t,type_real)
112           ! the following are unused but must point to something
113           f_massm2 => f_mass
114        ELSE ! eta = Lagrangian vertical coordinate
115           CALL allocate_field(f_massm2,field_t,type_real,llm)
116           ! the following are unused but must point to something
117           f_psm2 => f_phis
118        END IF
119      CASE ('none')
120        nb_stage=0
121
122    CASE default
123       PRINT*,'Bad selector for variable scheme : <', TRIM(def),             &
124            ' > options are <euler>, <runge_kutta>, <leapfrog_matsuno>'
125       STOP
126    END SELECT
127
128    CALL init_theta2theta_rhodz
129    CALL init_dissip
130    CALL init_sponge
131    CALL init_caldyn
132    CALL init_guided
133    CALL init_advect_tracer
134    CALL init_check_conserve
135
136    CALL etat0(f_ps,f_mass,f_phis,f_theta_rhodz,f_u, f_q)
137
138    CALL transfert_request(f_phis,req_i0) 
139    CALL transfert_request(f_phis,req_i1) 
140    CALL writefield("phis",f_phis,once=.TRUE.)
141
142    CALL init_message(f_ps,req_i0,req_ps0)
143    CALL init_message(f_mass,req_i0,req_mass0)
144    CALL init_message(f_theta_rhodz,req_i0,req_theta_rhodz0)
145    CALL init_message(f_u,req_e0_vect,req_u0)
146    CALL init_message(f_q,req_i0,req_q0)
147
148  END SUBROUTINE init_timeloop
149 
150  SUBROUTINE timeloop
151  USE icosa
152  USE dissip_gcm_mod
153  USE sponge_mod
154  USE disvert_mod
155  USE caldyn_mod
156  USE caldyn_gcm_mod, ONLY : req_ps, req_mass
157  USE etat0_mod
158  USE guided_mod
159  USE advect_tracer_mod
160  USE physics_mod
161  USE mpipara
162  USE omp_para
163  USE trace
164  USE transfert_mod
165  USE check_conserve_mod
166  USE xios_mod
167  USE output_field_mod
168  USE write_etat0_mod
169  USE checksum_mod
170  IMPLICIT NONE 
171    REAL(rstd),POINTER :: q(:,:,:)
172    REAL(rstd),POINTER :: phis(:), ps(:) ,psm1(:), psm2(:), dps(:)
173    REAL(rstd),POINTER :: u(:,:) , um1(:,:), um2(:,:), du(:,:)
174    REAL(rstd),POINTER :: rhodz(:,:), mass(:,:), massm1(:,:), massm2(:,:), dmass(:,:)
175    REAL(rstd),POINTER :: theta_rhodz(:,:) , theta_rhodzm1(:,:), theta_rhodzm2(:,:), dtheta_rhodz(:,:)
176    REAL(rstd),POINTER :: hflux(:,:),wflux(:,:),hfluxt(:,:),wfluxt(:,:)
177
178    INTEGER :: ind
179    INTEGER :: it,i,j,n,  stage
180    LOGICAL :: fluxt_zero(ndomain) ! set to .TRUE. to start accumulating fluxes in time
181    LOGICAL, PARAMETER :: check=.FALSE.
182    INTEGER :: start_clock
183    INTEGER :: stop_clock
184    INTEGER :: rate_clock
185    INTEGER :: l
186    LOGICAL,SAVE :: first_physic=.TRUE.
187!$OMP THREADPRIVATE(first_physic)   
188   
189   
190    CALL switch_omp_distrib_level
191    CALL caldyn_BC(f_phis, f_wflux) ! set constant values in first/last interfaces
192   
193!$OMP BARRIER
194  DO ind=1,ndomain
195     IF (.NOT. assigned_domain(ind)) CYCLE
196     CALL swap_dimensions(ind)
197     CALL swap_geometry(ind)
198     rhodz=f_rhodz(ind); mass=f_mass(ind); ps=f_ps(ind)
199     IF(caldyn_eta==eta_mass) THEN
200        CALL compute_rhodz(.TRUE., ps, rhodz) ! save rhodz for transport scheme before dynamics update ps
201     ELSE
202       DO l=ll_begin,ll_end
203         rhodz(:,l)=mass(:,l)
204       ENDDO
205     END IF
206  END DO
207!$OMP BARRIER
208  fluxt_zero=.TRUE.
209
210!$OMP MASTER
211  CALL SYSTEM_CLOCK(start_clock)
212!$OMP END MASTER   
213
214  CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,itau0) 
215
216  CALL trace_on
217 
218  DO it=itau0+1,itau0+itaumax
219     
220    CALL check_conserve_detailed('detailed_budget 0', &
221         f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it)
222
223    IF (xios_output) CALL xios_update_calendar(it)
224    IF (it==itau0+1 .OR. MOD(it,itau_sync)==0) THEN
225      CALL send_message(f_ps,req_ps0)
226      CALL wait_message(req_ps0)
227      CALL send_message(f_mass,req_mass0)
228      CALL wait_message(req_mass0)
229      CALL send_message(f_theta_rhodz,req_theta_rhodz0) 
230      CALL wait_message(req_theta_rhodz0) 
231      CALL send_message(f_u,req_u0)
232      CALL wait_message(req_u0)
233      CALL send_message(f_q,req_q0) 
234      CALL wait_message(req_q0) 
235
236    ENDIF
237
238    IF (is_master) PRINT *,"It No :",It,"   t :",dt*It
239
240    IF (mod(it,itau_out)==0 ) THEN
241      CALL update_time_counter(dt*it)
242      CALL output_field("q",f_q)
243    ENDIF
244   
245    CALL guided(it*dt,f_ps,f_theta_rhodz,f_u,f_q)
246
247    DO stage=1,nb_stage
248       CALL caldyn((stage==1) .AND. (MOD(it,itau_out)==0), &
249            f_phis,f_ps,f_mass,f_theta_rhodz,f_u, f_q, &
250            f_hflux, f_wflux, f_dps, f_dmass, f_dtheta_rhodz, f_du)
251       SELECT CASE (scheme)
252       CASE(euler)
253          CALL euler_scheme(.TRUE.)
254       CASE (rk4)
255          CALL rk_scheme(stage, coef_rk4)
256       CASE (rk25)
257          CALL rk_scheme(stage, coef_rk25)
258       CASE (mlf)
259          CALL  leapfrog_matsuno_scheme(stage)
260       CASE DEFAULT
261          STOP
262       END SELECT
263    END DO
264
265    CALL check_conserve_detailed('detailed_budget 1', &
266         f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it)
267
268    IF (MOD(it,itau_dissip)==0) THEN
269       
270       IF(caldyn_eta==eta_mass) THEN
271!ym flush ps
272!$OMP BARRIER
273          DO ind=1,ndomain
274             IF (.NOT. assigned_domain(ind)) CYCLE
275             CALL swap_dimensions(ind)
276             CALL swap_geometry(ind)
277             mass=f_mass(ind); ps=f_ps(ind);
278             CALL compute_rhodz(.TRUE., ps, mass)
279          END DO
280       ENDIF
281
282       CALL dissip(f_u,f_du,f_mass,f_phis, f_theta_rhodz,f_dtheta_rhodz)
283
284       CALL euler_scheme(.FALSE.)  ! update only u, theta
285       IF (iflag_sponge > 0) THEN
286         CALL sponge(f_u,f_du,f_theta_rhodz,f_dtheta_rhodz)
287         CALL euler_scheme(.FALSE.)  ! update only u, theta
288       ENDIF
289    END IF
290
291    CALL check_conserve_detailed('detailed_budget 2', &
292         f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it)
293
294    IF(MOD(it,itau_adv)==0) THEN
295
296       CALL advect_tracer(f_hfluxt,f_wfluxt,f_u, f_q,f_rhodz)  ! update q and rhodz after RK step
297       fluxt_zero=.TRUE.
298
299       ! FIXME : check that rhodz is consistent with ps
300       IF (check) THEN
301         DO ind=1,ndomain
302            IF (.NOT. assigned_domain(ind)) CYCLE
303            CALL swap_dimensions(ind)
304            CALL swap_geometry(ind)
305            rhodz=f_rhodz(ind); ps=f_ps(ind);
306            CALL compute_rhodz(.FALSE., ps, rhodz)   
307         END DO
308       ENDIF
309
310    END IF
311
312
313    IF (MOD(it,itau_check_conserv)==0) THEN
314      CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) 
315    ENDIF
316     
317    IF (MOD(it,itau_physics)==0) THEN
318      CALL physics(it,f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q)
319
320!$OMP MASTER
321   IF (first_physic) CALL SYSTEM_CLOCK(start_clock)
322!$OMP END MASTER   
323      first_physic=.FALSE.
324    ENDIF
325   
326  ENDDO
327
328  CALL write_etat0(itau0+itaumax,f_ps, f_phis,f_theta_rhodz,f_u,f_q) 
329
330  CALL check_conserve(f_ps,f_dps,f_u,f_theta_rhodz,f_phis,it) 
331
332!$OMP MASTER
333  CALL SYSTEM_CLOCK(stop_clock)
334  CALL SYSTEM_CLOCK(count_rate=rate_clock)
335   
336  IF (mpi_rank==0) THEN
337    PRINT *,"Time elapsed : ",(stop_clock-start_clock)*1./rate_clock 
338  ENDIF 
339!$OMP END MASTER
340 
341 CONTAINS
342
343    SUBROUTINE Euler_scheme(with_dps)
344    IMPLICIT NONE
345    LOGICAL :: with_dps
346    INTEGER :: ind
347    INTEGER :: i,j,ij,l
348    CALL trace_start("Euler_scheme") 
349
350    DO ind=1,ndomain
351       IF (.NOT. assigned_domain(ind)) CYCLE
352       CALL swap_dimensions(ind)
353       CALL swap_geometry(ind)
354
355       IF(with_dps) THEN ! update ps/mass
356          IF(caldyn_eta==eta_mass) THEN ! update ps
357             ps=f_ps(ind) ; dps=f_dps(ind) ;             
358             IF (is_omp_first_level) THEN
359!$SIMD
360                DO ij=ij_begin,ij_end
361                    ps(ij)=ps(ij)+dt*dps(ij)
362                ENDDO
363             ENDIF
364          ELSE ! update mass
365             mass=f_mass(ind) ; dmass=f_dmass(ind) ;             
366             DO l=ll_begin,ll_end
367!$SIMD
368                DO ij=ij_begin,ij_end
369                    mass(ij,l)=mass(ij,l)+dt*dmass(ij,l)
370                ENDDO
371             END DO
372          END IF
373
374          hflux=f_hflux(ind);     hfluxt=f_hfluxt(ind)
375          wflux=f_wflux(ind);     wfluxt=f_wfluxt(ind)
376          CALL accumulate_fluxes(hflux,wflux,hfluxt,wfluxt,dt,fluxt_zero(ind))
377       END IF ! update ps/mass
378       
379       u=f_u(ind) ; theta_rhodz=f_theta_rhodz(ind)
380       du=f_du(ind) ; dtheta_rhodz=f_dtheta_rhodz(ind)
381
382       DO l=ll_begin,ll_end
383!$SIMD
384         DO ij=ij_begin,ij_end
385             u(ij+u_right,l)=u(ij+u_right,l)+dt*du(ij+u_right,l)
386             u(ij+u_lup,l)=u(ij+u_lup,l)+dt*du(ij+u_lup,l)
387             u(ij+u_ldown,l)=u(ij+u_ldown,l)+dt*du(ij+u_ldown,l)
388             theta_rhodz(ij,l)=theta_rhodz(ij,l)+dt*dtheta_rhodz(ij,l)
389         ENDDO
390       ENDDO
391    ENDDO
392
393    CALL trace_end("Euler_scheme") 
394
395    END SUBROUTINE Euler_scheme
396
397    SUBROUTINE RK_scheme(stage,coef)
398      IMPLICIT NONE
399      INTEGER :: ind, stage
400      REAL(rstd), INTENT(IN) :: coef(:)
401      REAL(rstd) :: tau
402      INTEGER :: i,j,ij,l
403 
404      CALL trace_start("RK_scheme") 
405
406      tau = dt*coef(stage)
407
408      ! if mass coordinate, deal with ps first on one core
409      IF(caldyn_eta==eta_mass) THEN
410         IF (is_omp_first_level) THEN
411
412            DO ind=1,ndomain
413               IF (.NOT. assigned_domain(ind)) CYCLE
414               CALL swap_dimensions(ind)
415               CALL swap_geometry(ind)
416               ps=f_ps(ind) ; psm1=f_psm1(ind) ; dps=f_dps(ind) 
417               
418               IF (stage==1) THEN ! first stage : save model state in XXm1
419!$SIMD
420                 DO ij=ij_begin,ij_end
421                   psm1(ij)=ps(ij)
422                 ENDDO
423               ENDIF
424           
425               ! updates are of the form x1 := x0 + tau*f(x1)
426!$SIMD
427               DO ij=ij_begin,ij_end
428                   ps(ij)=psm1(ij)+tau*dps(ij)
429               ENDDO
430            ENDDO
431         ENDIF
432!         CALL send_message(f_ps,req_ps)
433!ym no overlap for now         
434!         CALL wait_message(req_ps) 
435     
436      ELSE ! Lagrangian coordinate, deal with mass
437         DO ind=1,ndomain
438            IF (.NOT. assigned_domain(ind)) CYCLE
439            CALL swap_dimensions(ind)
440            CALL swap_geometry(ind)
441            mass=f_mass(ind); dmass=f_dmass(ind); massm1=f_massm1(ind)
442
443            IF (stage==1) THEN ! first stage : save model state in XXm1
444               DO l=ll_begin,ll_end
445!$SIMD
446                 DO ij=ij_begin,ij_end
447                    massm1(ij,l)=mass(ij,l)
448                 ENDDO
449               ENDDO
450            END IF
451
452            ! updates are of the form x1 := x0 + tau*f(x1)
453            DO l=ll_begin,ll_end
454!$SIMD
455               DO ij=ij_begin,ij_end
456                   mass(ij,l)=massm1(ij,l)+tau*dmass(ij,l)
457               ENDDO
458            ENDDO
459         END DO
460!         CALL send_message(f_mass,req_mass)
461!ym no overlap for now         
462!         CALL wait_message(req_mass) 
463
464      END IF
465
466      ! now deal with other prognostic variables
467      DO ind=1,ndomain
468         IF (.NOT. assigned_domain(ind)) CYCLE
469         CALL swap_dimensions(ind)
470         CALL swap_geometry(ind)
471         u=f_u(ind)      ; du=f_du(ind)      ; um1=f_um1(ind) 
472         theta_rhodz=f_theta_rhodz(ind)
473         theta_rhodzm1=f_theta_rhodzm1(ind)
474         dtheta_rhodz=f_dtheta_rhodz(ind)
475         
476         IF (stage==1) THEN ! first stage : save model state in XXm1
477           DO l=ll_begin,ll_end
478!$SIMD
479                DO ij=ij_begin,ij_end
480                 um1(ij+u_right,l)=u(ij+u_right,l)
481                 um1(ij+u_lup,l)=u(ij+u_lup,l)
482                 um1(ij+u_ldown,l)=u(ij+u_ldown,l)
483                 theta_rhodzm1(ij,l)=theta_rhodz(ij,l)
484             ENDDO
485           ENDDO
486         END IF       
487
488         DO l=ll_begin,ll_end
489!$SIMD
490           DO ij=ij_begin,ij_end
491               u(ij+u_right,l)=um1(ij+u_right,l)+tau*du(ij+u_right,l)
492               u(ij+u_lup,l)=um1(ij+u_lup,l)+tau*du(ij+u_lup,l)
493               u(ij+u_ldown,l)=um1(ij+u_ldown,l)+tau*du(ij+u_ldown,l)
494               theta_rhodz(ij,l)=theta_rhodzm1(ij,l)+tau*dtheta_rhodz(ij,l)
495           ENDDO
496         ENDDO
497
498         IF(stage==nb_stage) THEN ! accumulate mass fluxes at last stage
499            hflux=f_hflux(ind);     hfluxt=f_hfluxt(ind)
500            wflux=f_wflux(ind);     wfluxt=f_wfluxt(ind)
501            CALL accumulate_fluxes(hflux,wflux, hfluxt,wfluxt, dt,fluxt_zero(ind))
502         END IF
503      END DO
504     
505      CALL trace_end("RK_scheme")
506     
507    END SUBROUTINE RK_scheme
508
509    SUBROUTINE leapfrog_matsuno_scheme(stage)
510    IMPLICIT NONE
511    INTEGER :: ind, stage
512    REAL :: tau
513
514      CALL trace_start("leapfrog_matsuno_scheme")
515   
516      tau = dt/nb_stage
517      DO ind=1,ndomain
518        IF (.NOT. assigned_domain(ind)) CYCLE
519        CALL swap_dimensions(ind)
520        CALL swap_geometry(ind)
521
522        ps=f_ps(ind)   ; u=f_u(ind)   ; theta_rhodz=f_theta_rhodz(ind)
523        psm1=f_psm1(ind) ; um1=f_um1(ind) ; theta_rhodzm1=f_theta_rhodzm1(ind)
524        psm2=f_psm2(ind) ; um2=f_um2(ind) ; theta_rhodzm2=f_theta_rhodzm2(ind)
525        dps=f_dps(ind) ; du=f_du(ind) ; dtheta_rhodz=f_dtheta_rhodz(ind)
526
527       
528        IF (stage==1) THEN
529          psm1(:)=ps(:) ; um1(:,:)=u(:,:) ; theta_rhodzm1(:,:)=theta_rhodz(:,:)
530          ps(:)=psm1(:)+tau*dps(:)
531          u(:,:)=um1(:,:)+tau*du(:,:)
532          theta_rhodz(:,:)=theta_rhodzm1(:,:)+tau*dtheta_rhodz(:,:)
533
534        ELSE IF (stage==2) THEN
535
536          ps(:)=psm1(:)+tau*dps(:)
537          u(:,:)=um1(:,:)+tau*du(:,:)
538          theta_rhodz(:,:)=theta_rhodzm1(:,:)+tau*dtheta_rhodz(:,:)
539
540          psm2(:)=psm1(:) ; theta_rhodzm2(:,:)=theta_rhodzm1(:,:) ; um2(:,:)=um1(:,:) 
541          psm1(:)=ps(:) ; theta_rhodzm1(:,:)=theta_rhodz(:,:) ; um1(:,:)=u(:,:) 
542
543        ELSE
544
545          ps(:)=psm2(:)+2*tau*dps(:)
546          u(:,:)=um2(:,:)+2*tau*du(:,:)
547          theta_rhodz(:,:)=theta_rhodzm2(:,:)+2*tau*dtheta_rhodz(:,:)
548
549          psm2(:)=psm1(:) ; theta_rhodzm2(:,:)=theta_rhodzm1(:,:) ; um2(:,:)=um1(:,:) 
550          psm1(:)=ps(:) ; theta_rhodzm1(:,:)=theta_rhodz(:,:) ; um1(:,:)=u(:,:) 
551
552        ENDIF
553     
554      ENDDO
555      CALL trace_end("leapfrog_matsuno_scheme")
556     
557    END SUBROUTINE leapfrog_matsuno_scheme 
558         
559  END SUBROUTINE timeloop   
560
561 SUBROUTINE accumulate_fluxes(hflux,wflux, hfluxt,wfluxt, tau,fluxt_zero)
562    USE icosa
563    USE omp_para
564    USE disvert_mod
565    IMPLICIT NONE
566    REAL(rstd), INTENT(IN)    :: hflux(3*iim*jjm,llm), wflux(iim*jjm,llm+1)
567    REAL(rstd), INTENT(INOUT) :: hfluxt(3*iim*jjm,llm), wfluxt(iim*jjm,llm+1)
568    REAL(rstd), INTENT(IN) :: tau
569    LOGICAL, INTENT(INOUT) :: fluxt_zero
570    INTEGER :: l,i,j,ij
571
572    IF(fluxt_zero) THEN
573
574       fluxt_zero=.FALSE.
575
576       DO l=ll_begin,ll_end
577!$SIMD
578         DO ij=ij_begin_ext,ij_end_ext
579             hfluxt(ij+u_right,l) = tau*hflux(ij+u_right,l)
580             hfluxt(ij+u_lup,l) = tau*hflux(ij+u_lup,l)
581             hfluxt(ij+u_ldown,l) = tau*hflux(ij+u_ldown,l)
582         ENDDO
583       ENDDO
584
585       IF(caldyn_eta==eta_mass) THEN ! no need for vertical fluxes if eta_lag
586          DO l=ll_begin,ll_endp1
587!$SIMD
588             DO ij=ij_begin,ij_end
589                wfluxt(ij,l) = tau*wflux(ij,l)
590             ENDDO
591          ENDDO
592       END IF
593
594    ELSE
595
596       DO l=ll_begin,ll_end
597!$SIMD
598         DO ij=ij_begin_ext,ij_end_ext
599             hfluxt(ij+u_right,l) = hfluxt(ij+u_right,l)+tau*hflux(ij+u_right,l)
600             hfluxt(ij+u_lup,l) = hfluxt(ij+u_lup,l)+tau*hflux(ij+u_lup,l)
601             hfluxt(ij+u_ldown,l) = hfluxt(ij+u_ldown,l)+tau*hflux(ij+u_ldown,l)
602         ENDDO
603       ENDDO
604
605       IF(caldyn_eta==eta_mass) THEN ! no need for vertical fluxes if eta_lag
606          DO l=ll_begin,ll_endp1
607!$SIMD
608             DO ij=ij_begin,ij_end
609                   wfluxt(ij,l) = wfluxt(ij,l)+tau*wflux(ij,l)
610             ENDDO
611          ENDDO
612       END IF
613
614   END IF
615
616  END SUBROUTINE accumulate_fluxes
617 
618END MODULE timeloop_gcm_mod
Note: See TracBrowser for help on using the repository browser.