source: codes/icosagcm/trunk/src/dissip_gcm.f90 @ 288

Last change on this file since 288 was 286, checked in by dubos, 10 years ago

Partial etat0 cleanup (removed calls to xyz2lonlat)

File size: 25.7 KB
RevLine 
[15]1MODULE dissip_gcm_mod
[19]2  USE icosa
[15]3
[26]4  PRIVATE
5
6  TYPE(t_field),POINTER,SAVE :: f_due_diss1(:)
7  TYPE(t_field),POINTER,SAVE :: f_due_diss2(:)
8
[109]9  TYPE(t_field),POINTER,SAVE :: f_theta(:), f_phi(:), f_pk(:), f_pks(:), f_p(:)
[26]10  TYPE(t_field),POINTER,SAVE :: f_dtheta_diss(:)
11  TYPE(t_field),POINTER,SAVE :: f_dtheta_rhodz_diss(:)
[186]12  TYPE(t_message),SAVE :: req_due, req_dtheta 
[15]13 
[26]14  INTEGER,SAVE :: nitergdiv=1
[186]15!$OMP THREADPRIVATE(nitergdiv)
[26]16  INTEGER,SAVE :: nitergrot=1
[186]17!$OMP THREADPRIVATE(nitergrot)
[26]18  INTEGER,SAVE :: niterdivgrad=1
[186]19!$OMP THREADPRIVATE(niterdivgrad)
[15]20
21  REAL,ALLOCATABLE,SAVE :: tau_graddiv(:)
[186]22!$OMP THREADPRIVATE(tau_graddiv)
[15]23  REAL,ALLOCATABLE,SAVE :: tau_gradrot(:)
[186]24!$OMP THREADPRIVATE(tau_gradrot)
[15]25  REAL,ALLOCATABLE,SAVE :: tau_divgrad(:)
[186]26!$OMP THREADPRIVATE(tau_divgrad)
[15]27 
28  REAL,SAVE :: cgraddiv
[186]29!$OMP THREADPRIVATE(cgraddiv)
[15]30  REAL,SAVE :: cgradrot
[186]31!$OMP THREADPRIVATE(cgradrot)
[15]32  REAL,SAVE :: cdivgrad
[186]33!$OMP THREADPRIVATE(cdivgrad)
[109]34
35  INTEGER, SAVE :: rayleigh_friction_type, rayleigh_shear
[186]36!$OMP THREADPRIVATE(rayleigh_friction_type)
37  REAL, SAVE    :: rayleigh_tau
38!$OMP THREADPRIVATE(rayleigh_shear)
[15]39 
[26]40  REAL,SAVE    :: dtdissip
[186]41!$OMP THREADPRIVATE(dtdissip)
[151]42 
[26]43  PUBLIC init_dissip, dissip
[15]44CONTAINS
45
46  SUBROUTINE allocate_dissip
[19]47  USE icosa
[15]48  IMPLICIT NONE 
[26]49    CALL allocate_field(f_due_diss1,field_u,type_real,llm)
50    CALL allocate_field(f_due_diss2,field_u,type_real,llm)
51    CALL allocate_field(f_theta,field_t,type_real,llm)
52    CALL allocate_field(f_dtheta_diss,field_t,type_real,llm)
53    CALL allocate_field(f_dtheta_rhodz_diss,field_t,type_real,llm)
[109]54
55    CALL allocate_field(f_phi,field_t,type_real,llm)
56    CALL allocate_field(f_pk,field_t,type_real,llm)
57    CALL allocate_field(f_p,field_t,type_real,llm+1)
58    CALL allocate_field(f_pks,field_t,type_real)
[26]59   
[15]60    ALLOCATE(tau_graddiv(llm))
61    ALLOCATE(tau_gradrot(llm))   
62    ALLOCATE(tau_divgrad(llm))
63  END SUBROUTINE allocate_dissip
64 
[98]65  SUBROUTINE init_dissip
[19]66  USE icosa
[15]67  USE disvert_mod
[26]68  USE mpi_mod
69  USE mpipara
[151]70  USE transfert_mod
[148]71  USE time_mod
[186]72  USE transfert_omp_mod
[15]73  IMPLICIT NONE
74 
[186]75  TYPE(t_field),POINTER,SAVE  :: f_u(:)
76  TYPE(t_field),POINTER,SAVE  :: f_du(:)
77  REAL(rstd),POINTER     :: u(:)
78  REAL(rstd),POINTER     :: du(:)
79  TYPE(t_field),POINTER,SAVE  :: f_theta(:)
80  TYPE(t_field),POINTER ,SAVE :: f_dtheta(:)
[15]81  REAL(rstd),POINTER    :: theta(:)
82  REAL(rstd),POINTER    :: dtheta(:)
[59]83  REAL(rstd)            :: dumax,dumax1
84  REAL(rstd)            :: dthetamax,dthetamax1
[15]85  REAL(rstd)            :: r
86  REAL(rstd)            :: tau
87  REAL(rstd)            :: zz, zvert, fact
88  INTEGER               :: l
[109]89  CHARACTER(len=255)    :: rayleigh_friction_key
[148]90  REAL(rstd)            :: mintau
[174]91  INTEGER               :: seed_size
92  INTEGER,ALLOCATABLE   :: seed(:)
[15]93 
94           
[174]95  INTEGER :: i,j,ij,ind,it,iter
[15]96
[109]97   rayleigh_friction_key='none'
98   CALL getin("rayleigh_friction_type",rayleigh_friction_key)
99   SELECT CASE(TRIM(rayleigh_friction_key))
100   CASE('none')
101      rayleigh_friction_type=0
[131]102      IF (is_mpi_root) PRINT *, 'No Rayleigh friction'
[109]103   CASE('dcmip2_schaer_noshear')
104      rayleigh_friction_type=1
105      rayleigh_shear=0
[131]106      IF (is_mpi_root) PRINT *, 'Rayleigh friction : Schaer-like mountain without shear DCMIP2.1'
[109]107   CASE('dcmip2_schaer_shear')
108      rayleigh_shear=1
109      rayleigh_friction_type=2
[131]110      IF (is_mpi_root) PRINT *, 'Rayleigh friction : Schaer-like mountain with shear DCMIP2.2'
[109]111   CASE DEFAULT
[151]112      IF (is_mpi_root) PRINT *, 'Bad selector : rayleigh_friction_type =', TRIM(rayleigh_friction_key), ' in dissip_gcm.f90/init_dissip'
[109]113      STOP
114   END SELECT
115
116   IF(rayleigh_friction_type>0) THEN
117      rayleigh_tau=0.
118      CALL getin("rayleigh_friction_tau",rayleigh_tau)
119      rayleigh_tau = rayleigh_tau / scale_factor
120      IF(rayleigh_tau<=0) THEN
[131]121         IF (is_mpi_root) PRINT *, 'Forbidden : negative value for rayleigh_friction_tau =',rayleigh_tau
[109]122         STOP
123      END IF
124   END IF
125
[15]126    CALL allocate_dissip
127    CALL allocate_field(f_u,field_u,type_real)
128    CALL allocate_field(f_du,field_u,type_real)
129    CALL allocate_field(f_theta,field_t,type_real)
130    CALL allocate_field(f_dtheta,field_t,type_real)
[151]131   
132    CALL init_message(f_due_diss1,req_e1_vect,req_due)
133    CALL init_message(f_dtheta_diss,req_i1,req_dtheta)
[15]134
135    tau_graddiv(:)=5000
136    CALL getin("tau_graddiv",tau)
[32]137    tau_graddiv(:)=tau/scale_factor
[26]138
139    CALL getin("nitergdiv",nitergdiv)
[15]140   
141    tau_gradrot(:)=5000
[186]142    CALL getin("tau_gradrot",tau)
[32]143    tau_gradrot(:)=tau/scale_factor
[26]144
145    CALL getin("nitergrot",nitergrot)
[15]146   
147
148    tau_divgrad(:)=5000
149    CALL getin("tau_divgrad",tau)
[32]150    tau_divgrad(:)=tau/scale_factor
[26]151
152    CALL getin("niterdivgrad",niterdivgrad)
[15]153
154
155    cgraddiv=-1
156    cdivgrad=-1
157    cgradrot=-1
158   
159    DO ind=1,ndomain
[186]160      IF (.NOT. assigned_domain(ind)) CYCLE
[15]161      CALL swap_dimensions(ind)
162      CALL swap_geometry(ind)
163      u=f_u(ind)
164
165      DO j=jj_begin,jj_end
166        DO i=ii_begin,ii_end
[174]167          ij=(j-1)*iim+i   
[15]168          CALL RANDOM_NUMBER(r)
[174]169          u(ij+u_right)=r-0.5
[15]170          CALL RANDOM_NUMBER(r)
[174]171          u(ij+u_lup)=r-0.5
[15]172          CALL RANDOM_NUMBER(r)
[174]173          u(ij+u_ldown)=r-0.5
174        ENDDO
175      ENDDO       
[15]176    ENDDO
177 
178
179
[29]180    DO it=1,20
[26]181     
[15]182      dumax=0
[26]183      DO iter=1,nitergdiv
[146]184        CALL transfert_request(f_u,req_e1_vect)
[26]185        DO ind=1,ndomain
[186]186          IF (.NOT. assigned_domain(ind)) CYCLE
[26]187          CALL swap_dimensions(ind)
188          CALL swap_geometry(ind)
189          u=f_u(ind)
190          du=f_du(ind)
[151]191          CALL compute_gradiv(u,du,1,1)
[26]192          u=du
193        ENDDO
[15]194      ENDDO
[26]195
[146]196      CALL transfert_request(f_du,req_e1_vect)
[15]197     
198      DO ind=1,ndomain
[186]199        IF (.NOT. assigned_domain(ind)) CYCLE
[15]200        CALL swap_dimensions(ind)
201        CALL swap_geometry(ind)
202        u=f_u(ind)
203        du=f_du(ind)
[26]204         
[15]205        DO j=jj_begin,jj_end
206          DO i=ii_begin,ii_end
[174]207            ij=(j-1)*iim+i   
208            if (le(ij+u_right)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_right)))
209            if (le(ij+u_lup)>1e-100)   dumax=MAX(dumax,ABS(du(ij+u_lup)))
210            if (le(ij+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_ldown)))
[15]211          ENDDO
212        ENDDO
213      ENDDO
214
[59]215      IF (using_mpi) THEN
[186]216        CALL reduce_sum_omp(dumax,dumax1)
217!$OMP MASTER       
218        CALL MPI_ALLREDUCE(dumax1,dumax,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr)
219!$OMP END MASTER     
220        CALL bcast_omp(dumax) 
221      ELSE
222        CALL allreduce_sum_omp(dumax,dumax1)
[59]223        dumax=dumax1
[186]224      ENDIF 
[59]225                       
[15]226      DO ind=1,ndomain
[186]227        IF (.NOT. assigned_domain(ind)) CYCLE
[15]228        CALL swap_dimensions(ind)
229        CALL swap_geometry(ind)
230        u=f_u(ind)
231        du=f_du(ind)
232        u=du/dumax
233      ENDDO
[131]234      IF (is_mpi_root) PRINT *,"gradiv : it :",it ,": dumax",dumax
[15]235
236    ENDDO 
[131]237    IF (is_mpi_root) PRINT *,"gradiv : dumax",dumax
238    IF (is_mpi_root) PRINT *, 'mean T-cell edge size (km)', 1.45*radius/iim_glo/1000., &
239                              'effective T-cell half-edge size (km)', dumax**(-.5/nitergdiv)/1000
240    IF (is_mpi_root)  PRINT *, 'Max. time step assuming c=340 m/s and Courant number=2.8 :', &
241                               2.8/340.*dumax**(-.5/nitergdiv)
[129]242
[15]243    cgraddiv=dumax**(-1./nitergdiv)
[131]244    IF (is_mpi_root) PRINT *,"cgraddiv : ",cgraddiv
[15]245
246    DO ind=1,ndomain
[186]247      IF (.NOT. assigned_domain(ind)) CYCLE
[15]248      CALL swap_dimensions(ind)
249      CALL swap_geometry(ind)
250      u=f_u(ind)
251
[174]252      DO j=jj_begin,jj_end
[15]253        DO i=ii_begin,ii_end
[174]254          ij=(j-1)*iim+i   
[15]255          CALL RANDOM_NUMBER(r)
[174]256          u(ij+u_right)=r-0.5
[15]257          CALL RANDOM_NUMBER(r)
[174]258          u(ij+u_lup)=r-0.5
[15]259          CALL RANDOM_NUMBER(r)
[174]260          u(ij+u_ldown)=r-0.5
261        ENDDO
262      ENDDO       
[15]263    ENDDO
264
265
[29]266    DO it=1,20
[15]267 
268      dumax=0
[26]269      DO iter=1,nitergrot
[146]270        CALL transfert_request(f_u,req_e1_vect)
[26]271        DO ind=1,ndomain
[186]272          IF (.NOT. assigned_domain(ind)) CYCLE
[26]273          CALL swap_dimensions(ind)
274          CALL swap_geometry(ind)
275          u=f_u(ind)
276          du=f_du(ind)
[151]277          CALL compute_gradrot(u,du,1,1)
[26]278          u=du
279        ENDDO
[15]280      ENDDO
281
[146]282      CALL transfert_request(f_du,req_e1_vect)
[15]283     
284      DO ind=1,ndomain
[186]285        IF (.NOT. assigned_domain(ind)) CYCLE
[15]286        CALL swap_dimensions(ind)
287        CALL swap_geometry(ind)
288        u=f_u(ind)
289        du=f_du(ind)
290       
291        DO j=jj_begin,jj_end
292          DO i=ii_begin,ii_end
[174]293            ij=(j-1)*iim+i   
294            if (le(ij+u_right)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_right)))
295            if (le(ij+u_lup)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_lup)))
296            if (le(ij+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_ldown)))
[15]297          ENDDO
298        ENDDO
299      ENDDO
[26]300
[186]301      IF (using_mpi) THEN
302        CALL reduce_sum_omp(dumax,dumax1)
303!$OMP MASTER       
304        CALL MPI_ALLREDUCE(dumax1,dumax,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr)
305!$OMP END MASTER     
306        CALL bcast_omp(dumax) 
307      ELSE
308        CALL allreduce_sum_omp(dumax,dumax1)
[59]309        dumax=dumax1
[186]310      ENDIF 
311
[15]312     
313      DO ind=1,ndomain
[186]314        IF (.NOT. assigned_domain(ind)) CYCLE
[15]315        CALL swap_dimensions(ind)
316        CALL swap_geometry(ind)
317        u=f_u(ind)
318        du=f_du(ind)
319        u=du/dumax
320      ENDDO
321     
[131]322      IF (is_mpi_root) PRINT *,"gradrot : it :",it ,": dumax",dumax
[15]323
324    ENDDO 
[131]325    IF (is_mpi_root) PRINT *,"gradrot : dumax",dumax
[15]326 
[26]327    cgradrot=dumax**(-1./nitergrot) 
[131]328    IF (is_mpi_root) PRINT *,"cgradrot : ",cgradrot
[15]329   
330
331
332    DO ind=1,ndomain
[186]333      IF (.NOT. assigned_domain(ind)) CYCLE
[15]334      CALL swap_dimensions(ind)
335      CALL swap_geometry(ind)
336      theta=f_theta(ind)
337 
[174]338      DO j=jj_begin,jj_end
339        DO i=ii_begin,ii_end
340          ij=(j-1)*iim+i   
341          CALL RANDOM_NUMBER(r)
342          theta(ij)=r-0.5
[15]343        ENDDO
[174]344      ENDDO 
[15]345    ENDDO
346
[29]347    DO it=1,20
[15]348
349      dthetamax=0
[26]350      DO iter=1,niterdivgrad
351        CALL transfert_request(f_theta,req_i1)
352        DO ind=1,ndomain
[186]353          IF (.NOT. assigned_domain(ind)) CYCLE
[26]354          CALL swap_dimensions(ind)
355          CALL swap_geometry(ind)
356          theta=f_theta(ind)
357          dtheta=f_dtheta(ind)
[151]358          CALL compute_divgrad(theta,dtheta,1,1)
[26]359          theta=dtheta
360        ENDDO
[15]361      ENDDO
362
363      CALL transfert_request(f_dtheta,req_i1)
364     
365      DO ind=1,ndomain
[186]366        IF (.NOT. assigned_domain(ind)) CYCLE
[15]367        CALL swap_dimensions(ind)
368        CALL swap_geometry(ind)
369        theta=f_theta(ind)
370        dtheta=f_dtheta(ind)
371       
372        DO j=jj_begin,jj_end
373          DO i=ii_begin,ii_end
[174]374            ij=(j-1)*iim+i   
375            dthetamax=MAX(dthetamax,ABS(dtheta(ij)))
[15]376          ENDDO
377        ENDDO
378      ENDDO
[186]379
380      IF (using_mpi) THEN
381        CALL reduce_sum_omp(dthetamax ,dthetamax1)
382!$OMP MASTER       
383        CALL MPI_ALLREDUCE(dthetamax1,dthetamax,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr)
384!$OMP END MASTER     
385        CALL bcast_omp(dthetamax) 
386      ELSE
387        CALL allreduce_sum_omp(dthetamax,dthetamax1)
388        dumax=dumax1
389      ENDIF 
[174]390     
[131]391      IF (is_mpi_root) PRINT *,"divgrad : it :",it ,": dthetamax",dthetamax
[15]392
393      DO ind=1,ndomain
[186]394        IF (.NOT. assigned_domain(ind)) CYCLE
[15]395        CALL swap_dimensions(ind)
396        CALL swap_geometry(ind)
397        theta=f_theta(ind)
398        dtheta=f_dtheta(ind)
399        theta=dtheta/dthetamax
400      ENDDO
401    ENDDO 
402
[131]403    IF (is_mpi_root) PRINT *,"divgrad : divgrad",dthetamax
[15]404 
[26]405    cdivgrad=dthetamax**(-1./niterdivgrad) 
[131]406    IF (is_mpi_root) PRINT *,"cdivgrad : ",cdivgrad
[15]407
408     
409    fact=2
410    DO l=1,llm
[167]411       IF(ap_bp_present) THEN ! height-dependent dissipation
412          zz= 1. - preff/presnivs(l)
413       ELSE
414          zz = 0.
415       END IF
[15]416       zvert=fact-(fact-1)/(1+zz*zz)
417       tau_graddiv(l) = tau_graddiv(l)/zvert
418       tau_gradrot(l) = tau_gradrot(l)/zvert
419       tau_divgrad(l) = tau_divgrad(l)/zvert
420    ENDDO
[148]421
422    mintau=tau_graddiv(1)
423    DO l=1,llm
424      mintau=MIN(mintau,tau_graddiv(l))
425      mintau=MIN(mintau,tau_gradrot(l))
426      mintau=MIN(mintau,tau_divgrad(l))
427    ENDDO
[15]428       
[250]429    IF(mintau>0) THEN
430       itau_dissip=INT(mintau/dt)
431       dtdissip=itau_dissip*dt
432    ELSE
433       IF (is_mpi_root) PRINT *,"No dissipation time set, setting itau_dissip to 1000000000"
434       itau_dissip=100000000
435    END IF
[148]436    itau_dissip=MAX(1,itau_dissip)
437    IF (is_mpi_root) PRINT *,"mintau ",mintau,"itau_dissip",itau_dissip," dtdissip ",dtdissip
[15]438
439  END SUBROUTINE init_dissip 
440 
441 
[167]442  SUBROUTINE dissip(f_ue,f_due,f_mass,f_phis,f_theta_rhodz,f_dtheta_rhodz)
[19]443  USE icosa
[15]444  USE theta2theta_rhodz_mod
[109]445  USE pression_mod
446  USE exner_mod
447  USE geopotential_mod
[145]448  USE trace
[148]449  USE time_mod
[151]450  USE omp_para
[15]451  IMPLICIT NONE
452    TYPE(t_field),POINTER :: f_ue(:)
453    TYPE(t_field),POINTER :: f_due(:)
[167]454    TYPE(t_field),POINTER :: f_mass(:), f_phis(:)
[15]455    TYPE(t_field),POINTER :: f_theta_rhodz(:)
456    TYPE(t_field),POINTER :: f_dtheta_rhodz(:)
[26]457
[15]458    REAL(rstd),POINTER         :: due(:,:)
[109]459    REAL(rstd),POINTER         :: phi(:,:), ue(:,:)
[26]460    REAL(rstd),POINTER         :: due_diss1(:,:)
461    REAL(rstd),POINTER         :: due_diss2(:,:)
[15]462    REAL(rstd),POINTER         :: dtheta_rhodz(:,:)
[26]463    REAL(rstd),POINTER         :: dtheta_rhodz_diss(:,:)
[15]464
[109]465    INTEGER :: ind, shear
[174]466    INTEGER :: l,ij
[151]467
468!$OMP BARRIER
[15]469   
[145]470    CALL trace_start("dissip")
[26]471    CALL gradiv(f_ue,f_due_diss1)
472    CALL gradrot(f_ue,f_due_diss2)
[167]473
474    CALL divgrad_theta_rhodz(f_mass,f_theta_rhodz,f_dtheta_rhodz_diss)
[26]475   
[151]476! later for openmp   
477!    IF(rayleigh_friction_type>0) THEN
478!       CALL pression(f_ps, f_p)
479!       CALL exner(f_ps, f_p, f_pks, f_pk)
480!       CALL geopotential(f_phis,f_pks,f_pk,f_theta,f_phi)
481!    END IF
[109]482
[15]483    DO ind=1,ndomain
[186]484      IF (.NOT. assigned_domain(ind)) CYCLE
[15]485      CALL swap_dimensions(ind)
486      CALL swap_geometry(ind)
487      due=f_due(ind) 
[26]488      due_diss1=f_due_diss1(ind)
489      due_diss2=f_due_diss2(ind)
490      dtheta_rhodz=f_dtheta_rhodz(ind)
491      dtheta_rhodz_diss=f_dtheta_rhodz_diss(ind)
[151]492
493      DO l=ll_begin,ll_end
[174]494!$SIMD
495        DO ij=ij_begin,ij_end
[15]496
[174]497            due(ij+u_right,l) = -0.5*( due_diss1(ij+u_right,l)/tau_graddiv(l) + due_diss2(ij+u_right,l)/tau_gradrot(l))*itau_dissip 
498            due(ij+u_lup,l)   = -0.5*( due_diss1(ij+u_lup,l)  /tau_graddiv(l) + due_diss2(ij+u_lup,l)  /tau_gradrot(l))*itau_dissip
499            due(ij+u_ldown,l) = -0.5*( due_diss1(ij+u_ldown,l)/tau_graddiv(l) + due_diss2(ij+u_ldown,l)/tau_gradrot(l))*itau_dissip
[26]500
[174]501            dtheta_rhodz(ij,l) = -0.5*dtheta_rhodz_diss(ij,l)/tau_divgrad(l)*itau_dissip
[26]502        ENDDO
503      ENDDO
[109]504
[151]505!      dtheta_rhodz=0
506!      due=0
507
508! later for openmp 
509!      IF(rayleigh_friction_type>0) THEN
510!         phi=f_phi(ind)
511!         ue=f_ue(ind)
512!         DO l=1,llm
513!            DO j=jj_begin,jj_end
514!               DO i=ii_begin,ii_end
515!                  n=(j-1)*iim+i
516!                  CALL relax(t_right, u_right)
517!                  CALL relax(t_lup,   u_lup)
518!                  CALL relax(t_ldown, u_ldown)
519!               ENDDO
520!            ENDDO
521!         END DO
522!      END IF
[109]523   END DO
524
[145]525   CALL trace_end("dissip")
526
[151]527!$OMP BARRIER
528
[109]529    CONTAINS
530      SUBROUTINE relax(shift_t, shift_u)
531        USE dcmip_initial_conditions_test_1_2_3
532        REAL(rstd) :: z, ulon,ulat, lon,lat, & ! input to test2_schaer_mountain
533             p,hyam,hybm,w,t,phis,ps,rho,q, &   ! unused input/output to test2_schaer_mountain
534             fz, u3d(3), uref
535        REAL(rstd), PARAMETER :: zh=2e4,ztop=3e4  ! DCMIP values
536        LOGICAL :: hybrid_eta
537        INTEGER :: shift_u, shift_t, zcoords, nn
[174]538        z = (phi(ij,l)+phi(ij+shift_t,l))/(2.*g)
[109]539        IF(z>zh) THEN  ! relax only in the sponge layer z>zh
[151]540
[174]541           nn = ij+shift_u
[109]542           zcoords = 1 ; hybrid_eta = .FALSE. ! use z instead of p or hyam/hybm
[286]543           CALL test2_schaer_mountain(lon_e(nn),lat_e(nn),p,z,zcoords,hybrid_eta, &
[109]544                hyam,hybm,shear,ulon,ulat,w,t,phis,ps,rho,q)
545!           u3d = ulon*elon_e(nn,:) + ulat*elat_e(nn,:)
546           u3d = ulon*elon_e(nn,:) ! ulat=0
547           uref = sum(u3d*ep_e(nn,:))
548
549           fz = sin((pi/2)*(z-zh)/(ztop-zh))
550           fz = fz*fz/rayleigh_tau
551!           fz = 1/rayleigh_tau
552           due(nn,l) = due(nn,l) - fz*(ue(nn,l)-uref)
553!           due(nn,l) = due(nn,l) - fz*ue(nn,l)
554         END IF
555       END SUBROUTINE relax
[15]556     
557  END SUBROUTINE dissip
558
[26]559  SUBROUTINE gradiv(f_ue,f_due)
[19]560  USE icosa
[145]561  USE trace
[151]562  USE omp_para
[15]563  IMPLICIT NONE
[26]564    TYPE(t_field),POINTER :: f_ue(:)
565    TYPE(t_field),POINTER :: f_due(:)
566    REAL(rstd),POINTER    :: ue(:,:)
567    REAL(rstd),POINTER    :: due(:,:)
568    INTEGER :: ind
[174]569    INTEGER :: it,l,ij
[26]570       
[145]571    CALL trace_start("gradiv")
572
[26]573    DO ind=1,ndomain
[186]574      IF (.NOT. assigned_domain(ind)) CYCLE
[26]575      CALL swap_dimensions(ind)
576      CALL swap_geometry(ind)
577      ue=f_ue(ind)
578      due=f_due(ind) 
[151]579      DO  l = ll_begin, ll_end
[174]580!$SIMD
581        DO ij=ij_begin,ij_end
[151]582             due(ij+u_right,l)=ue(ij+u_right,l)
583             due(ij+u_lup,l)=ue(ij+u_lup,l)
584             due(ij+u_ldown,l)=ue(ij+u_ldown,l)
585        ENDDO
586      ENDDO
[26]587    ENDDO
[15]588
[26]589    DO it=1,nitergdiv
590       
[151]591      CALL send_message(f_due,req_due)
592      CALL wait_message(req_due)
[26]593       
594      DO ind=1,ndomain
[186]595        IF (.NOT. assigned_domain(ind)) CYCLE
[26]596        CALL swap_dimensions(ind)
597        CALL swap_geometry(ind)
598        due=f_due(ind) 
[151]599        CALL compute_gradiv(due,due,ll_begin,ll_end)
[26]600      ENDDO
601    ENDDO
[15]602
[145]603   CALL trace_end("gradiv")
604
[26]605  END SUBROUTINE gradiv
606 
607
608  SUBROUTINE gradrot(f_ue,f_due)
609  USE icosa
[145]610  USE trace
[151]611  USE omp_para
[26]612  IMPLICIT NONE
613    TYPE(t_field),POINTER :: f_ue(:)
614    TYPE(t_field),POINTER :: f_due(:)
615    REAL(rstd),POINTER    :: ue(:,:)
616    REAL(rstd),POINTER    :: due(:,:)
[15]617    INTEGER :: ind
[174]618    INTEGER :: it,l,ij
[26]619       
[145]620    CALL trace_start("gradrot")
621
[26]622    DO ind=1,ndomain
[186]623      IF (.NOT. assigned_domain(ind)) CYCLE
[26]624      CALL swap_dimensions(ind)
625      CALL swap_geometry(ind)
626      ue=f_ue(ind)
627      due=f_due(ind) 
[151]628      DO  l = ll_begin, ll_end
[174]629!$SIMD
630        DO ij=ij_begin,ij_end
[151]631             due(ij+u_right,l)=ue(ij+u_right,l)
632             due(ij+u_lup,l)=ue(ij+u_lup,l)
633             due(ij+u_ldown,l)=ue(ij+u_ldown,l)
634        ENDDO
635      ENDDO
[26]636    ENDDO
[15]637
[26]638    DO it=1,nitergrot
639       
[151]640      CALL send_message(f_due,req_due)
641      CALL wait_message(req_due)
[26]642       
643      DO ind=1,ndomain
[186]644        IF (.NOT. assigned_domain(ind)) CYCLE
[26]645        CALL swap_dimensions(ind)
646        CALL swap_geometry(ind)
647        due=f_due(ind) 
[151]648        CALL compute_gradrot(due,due,ll_begin,ll_end)
[15]649      ENDDO
650
[26]651    ENDDO
[15]652
[145]653    CALL trace_end("gradrot")
654
[26]655  END SUBROUTINE gradrot
656 
657  SUBROUTINE divgrad(f_theta,f_dtheta)
658  USE icosa
[145]659  USE trace
[151]660  USE omp_para
[26]661  IMPLICIT NONE
662    TYPE(t_field),POINTER :: f_theta(:)
663    TYPE(t_field),POINTER :: f_dtheta(:)
664    REAL(rstd),POINTER    :: theta(:,:)
665    REAL(rstd),POINTER    :: dtheta(:,:)
666    INTEGER :: ind
667    INTEGER :: it
[145]668
669    CALL trace_start("divgrad")
[26]670       
671    DO ind=1,ndomain
[186]672      IF (.NOT. assigned_domain(ind)) CYCLE
[26]673      CALL swap_dimensions(ind)
674      CALL swap_geometry(ind)
675      theta=f_theta(ind)
676      dtheta=f_dtheta(ind) 
677      dtheta=theta
678    ENDDO
[15]679
[26]680    DO it=1,niterdivgrad
681       
682      CALL transfert_request(f_dtheta,req_i1)
683       
684      DO ind=1,ndomain
[186]685        IF (.NOT. assigned_domain(ind)) CYCLE
[26]686        CALL swap_dimensions(ind)
687        CALL swap_geometry(ind)
688        dtheta=f_dtheta(ind) 
[151]689        CALL compute_divgrad(dtheta,dtheta,ll_begin,ll_end)
[26]690      ENDDO
[15]691
[26]692    ENDDO
693
[145]694    CALL trace_end("divgrad")
695
[26]696  END SUBROUTINE divgrad
697   
[167]698  SUBROUTINE divgrad_theta_rhodz(f_mass,f_theta_rhodz,f_dtheta_rhodz)
[151]699  USE icosa
700  USE trace
701  USE omp_para
702  IMPLICIT NONE
[167]703    TYPE(t_field),POINTER :: f_mass(:)
[151]704    TYPE(t_field),POINTER :: f_theta_rhodz(:)
705    TYPE(t_field),POINTER :: f_dtheta_rhodz(:)
706   
[167]707    REAL(rstd),POINTER :: mass(:,:)
[151]708    REAL(rstd),POINTER :: theta_rhodz(:,:)
709    REAL(rstd),POINTER :: dtheta_rhodz(:,:)
[26]710
[151]711    INTEGER :: ind
[174]712    INTEGER :: it,l,ij
[151]713
714    CALL trace_start("divgrad")
715       
716    DO ind=1,ndomain
[186]717      IF (.NOT. assigned_domain(ind)) CYCLE
[151]718      CALL swap_dimensions(ind)
719      CALL swap_geometry(ind)
[167]720      mass=f_mass(ind)
[151]721      theta_rhodz=f_theta_rhodz(ind)
722      dtheta_rhodz=f_dtheta_rhodz(ind) 
723      DO  l = ll_begin, ll_end
[174]724!$SIMD
725        DO ij=ij_begin,ij_end
[167]726            dtheta_rhodz(ij,l) = theta_rhodz(ij,l) / mass(ij,l)
[151]727        ENDDO
728      ENDDO
729    ENDDO
730
731    DO it=1,niterdivgrad
732       
733      CALL send_message(f_dtheta_rhodz,req_dtheta)
734      CALL wait_message(req_dtheta)
735      DO ind=1,ndomain
[186]736        IF (.NOT. assigned_domain(ind)) CYCLE
[151]737        CALL swap_dimensions(ind)
738        CALL swap_geometry(ind)
739        dtheta_rhodz=f_dtheta_rhodz(ind) 
740        CALL compute_divgrad(dtheta_rhodz,dtheta_rhodz,ll_begin,ll_end)
741      ENDDO
742
743    ENDDO
744
745    DO ind=1,ndomain
[186]746      IF (.NOT. assigned_domain(ind)) CYCLE
[151]747      CALL swap_dimensions(ind)
748      CALL swap_geometry(ind)
749      dtheta_rhodz=f_dtheta_rhodz(ind) 
[175]750      mass=f_mass(ind)
[151]751   
752      DO  l = ll_begin, ll_end
[174]753!$SIMD
754        DO ij=ij_begin,ij_end
[167]755            dtheta_rhodz(ij,l) = dtheta_rhodz(ij,l) * mass(ij,l)
[151]756        ENDDO
757      ENDDO
758    ENDDO
759
760
761    CALL trace_end("divgrad")
762
763  END SUBROUTINE divgrad_theta_rhodz 
[15]764 
[151]765  SUBROUTINE compute_gradiv(ue,gradivu_e,llb,lle)
[19]766  USE icosa
[15]767  IMPLICIT NONE
[151]768    INTEGER,INTENT(IN)     :: llb
769    INTEGER,INTENT(IN)     :: lle
770    REAL(rstd),INTENT(IN)  :: ue(iim*3*jjm,llm)
771    REAL(rstd),INTENT(OUT) :: gradivu_e(iim*3*jjm,llm)
772    REAL(rstd) :: divu_i(iim*jjm,llb:lle)
[15]773   
[174]774    INTEGER :: ij,l
[15]775
[151]776    DO l=llb,lle
[174]777!$SIMD
778      DO ij=ij_begin,ij_end
779          divu_i(ij,l)=1./Ai(ij)*(ne(ij,right)*ue(ij+u_right,l)*le(ij+u_right)  +  &
780                             ne(ij,rup)*ue(ij+u_rup,l)*le(ij+u_rup)        +  & 
781                             ne(ij,lup)*ue(ij+u_lup,l)*le(ij+u_lup)        +  & 
782                             ne(ij,left)*ue(ij+u_left,l)*le(ij+u_left)     +  & 
783                             ne(ij,ldown)*ue(ij+u_ldown,l)*le(ij+u_ldown)  +  & 
784                             ne(ij,rdown)*ue(ij+u_rdown,l)*le(ij+u_rdown))
[15]785      ENDDO
786    ENDDO
787   
[151]788    DO l=llb,lle
[174]789!$SIMD
790      DO ij=ij_begin,ij_end
[15]791 
[174]792          gradivu_e(ij+u_right,l)=-1/de(ij+u_right)*(ne(ij,right)*divu_i(ij,l)+ ne(ij+t_right,left)*divu_i(ij+t_right,l) )       
[15]793
[174]794          gradivu_e(ij+u_lup,l)=-1/de(ij+u_lup)*(ne(ij,lup)*divu_i(ij,l)+ ne(ij+t_lup,rdown)*divu_i(ij+t_lup,l))       
[15]795   
[174]796          gradivu_e(ij+u_ldown,l)=-1/de(ij+u_ldown)*(ne(ij,ldown)*divu_i(ij,l)+ne(ij+t_ldown,rup)*divu_i(ij+t_ldown,l) )
[15]797
798      ENDDO
799    ENDDO
800
[151]801    DO l=llb,lle
[174]802!$SIMD
803      DO ij=ij_begin,ij_end
804          gradivu_e(ij+u_right,l)=-gradivu_e(ij+u_right,l)*cgraddiv
805          gradivu_e(ij+u_lup,l)=-gradivu_e(ij+u_lup,l)*cgraddiv
806          gradivu_e(ij+u_ldown,l)=-gradivu_e(ij+u_ldown,l)*cgraddiv
[15]807      ENDDO
808    ENDDO
809
[148]810   
[26]811  END SUBROUTINE compute_gradiv
[15]812 
[151]813  SUBROUTINE compute_divgrad(theta,divgrad_i,llb,lle)
[19]814  USE icosa
[15]815  IMPLICIT NONE
[151]816    INTEGER,INTENT(IN)     :: llb
817    INTEGER,INTENT(IN)     :: lle
818    REAL(rstd),INTENT(IN)  :: theta(iim*jjm,llm)
819    REAL(rstd),INTENT(OUT) :: divgrad_i(iim*jjm,llm)
820    REAL(rstd)  :: grad_e(3*iim*jjm,llb:lle)
[15]821
[174]822    INTEGER :: ij,l
[15]823
[148]824       
[151]825    DO l=llb,lle
[174]826!$SIMD
827      DO ij=ij_begin_ext,ij_end_ext
[15]828 
[174]829          grad_e(ij+u_right,l)=-1/de(ij+u_right)*(ne(ij,right)*theta(ij,l)+ ne(ij+t_right,left)*theta(ij+t_right,l) )       
[15]830
[174]831          grad_e(ij+u_lup,l)=-1/de(ij+u_lup)*(ne(ij,lup)*theta(ij,l)+ ne(ij+t_lup,rdown)*theta(ij+t_lup,l ))       
[15]832   
[174]833          grad_e(ij+u_ldown,l)=-1/de(ij+u_ldown)*(ne(ij,ldown)*theta(ij,l)+ne(ij+t_ldown,rup)*theta(ij+t_ldown,l) )
[15]834
835      ENDDO
836    ENDDO
837   
838   
[151]839    DO l=llb,lle
[174]840!$SIMD
[175]841      DO ij=ij_begin,ij_end
[174]842
843          divgrad_i(ij,l)=1./Ai(ij)*(ne(ij,right)*grad_e(ij+u_right,l)*le(ij+u_right)  +  &
844                             ne(ij,rup)*grad_e(ij+u_rup,l)*le(ij+u_rup)              +  & 
845                             ne(ij,lup)*grad_e(ij+u_lup,l)*le(ij+u_lup)              +  & 
846                             ne(ij,left)*grad_e(ij+u_left,l)*le(ij+u_left)           +  & 
847                             ne(ij,ldown)*grad_e(ij+u_ldown,l)*le(ij+u_ldown)        +  & 
848                             ne(ij,rdown)*grad_e(ij+u_rdown,l)*le(ij+u_rdown))
[15]849      ENDDO
850    ENDDO
851   
[151]852    DO l=llb,lle
[174]853      DO ij=ij_begin,ij_end
854          divgrad_i(ij,l)=-divgrad_i(ij,l)*cdivgrad
[15]855      ENDDO
856    ENDDO
857
[26]858  END SUBROUTINE compute_divgrad
[15]859
860   
[151]861  SUBROUTINE compute_gradrot(ue,gradrot_e,llb,lle)
[19]862  USE icosa
[15]863  IMPLICIT NONE
[151]864    INTEGER,INTENT(IN)     :: llb
865    INTEGER,INTENT(IN)     :: lle
866    REAL(rstd),INTENT(IN)  :: ue(iim*3*jjm,llm)
867    REAL(rstd),INTENT(OUT) :: gradrot_e(iim*3*jjm,llm)
868    REAL(rstd) :: rot_v(2*iim*jjm,llb:lle)
[15]869
[174]870    INTEGER :: ij,l
[15]871     
[151]872    DO l=llb,lle
[174]873!$SIMD
874      DO ij=ij_begin_ext,ij_end_ext
[15]875       
[174]876          rot_v(ij+z_up,l)= 1./Av(ij+z_up)*(  ne(ij,rup)*ue(ij+u_rup,l)*de(ij+u_rup)                     &
877                                + ne(ij+t_rup,left)*ue(ij+t_rup+u_left,l)*de(ij+t_rup+u_left)          &
878                                - ne(ij,lup)*ue(ij+u_lup,l)*de(ij+u_lup) ) 
[15]879                             
[174]880          rot_v(ij+z_down,l) = 1./Av(ij+z_down)*( ne(ij,ldown)*ue(ij+u_ldown,l)*de(ij+u_ldown)                 &
881                                     + ne(ij+t_ldown,right)*ue(ij+t_ldown+u_right,l)*de(ij+t_ldown+u_right)  &
882                                     - ne(ij,rdown)*ue(ij+u_rdown,l)*de(ij+u_rdown) )
[15]883 
884      ENDDO
885    ENDDO                             
886 
[151]887    DO l=llb,lle
[174]888!$SIMD
889      DO ij=ij_begin,ij_end
[15]890 
[174]891          gradrot_e(ij+u_right,l)=1/le(ij+u_right)*ne(ij,right)*(rot_v(ij+z_rdown,l)-rot_v(ij+z_rup,l)) 
[15]892         
[174]893          gradrot_e(ij+u_lup,l)=1/le(ij+u_lup)*ne(ij,lup)*(rot_v(ij+z_up,l)-rot_v(ij+z_lup,l)) 
[15]894       
[174]895          gradrot_e(ij+u_ldown,l)=1/le(ij+u_ldown)*ne(ij,ldown)*(rot_v(ij+z_ldown,l)-rot_v(ij+z_down,l))
[15]896       
897      ENDDO
898    ENDDO
899
[151]900    DO l=llb,lle
[174]901!$SIMD
902      DO ij=ij_begin,ij_end
[15]903   
[174]904          gradrot_e(ij+u_right,l)=-gradrot_e(ij+u_right,l)*cgradrot       
905          gradrot_e(ij+u_lup,l)=-gradrot_e(ij+u_lup,l)*cgradrot
906          gradrot_e(ij+u_ldown,l)=-gradrot_e(ij+u_ldown,l)*cgradrot
[15]907       
908      ENDDO
909    ENDDO 
910   
[26]911  END SUBROUTINE compute_gradrot
[15]912
913
914END MODULE dissip_gcm_mod
915       
Note: See TracBrowser for help on using the repository browser.