source: codes/icosagcm/branches/SATURN_DYNAMICO/ICOSAGCM/src/timeloop_gcm.f90 @ 303

Last change on this file since 303 was 282, checked in by milmd, 10 years ago

Sponge layer in dynamics.

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