MODULE dissip_gcm_mod USE icosa PRIVATE TYPE(t_field),POINTER,SAVE :: f_due_diss1(:) TYPE(t_field),POINTER,SAVE :: f_due_diss2(:) TYPE(t_field),POINTER,SAVE :: f_theta(:), f_phi(:), f_pk(:), f_pks(:), f_p(:) TYPE(t_field),POINTER,SAVE :: f_dtheta_diss(:) TYPE(t_field),POINTER,SAVE :: f_dtheta_rhodz_diss(:) TYPE(t_message),SAVE :: req_due, req_dtheta INTEGER,SAVE :: nitergdiv=1 !$OMP THREADPRIVATE(nitergdiv) INTEGER,SAVE :: nitergrot=1 !$OMP THREADPRIVATE(nitergrot) INTEGER,SAVE :: niterdivgrad=1 !$OMP THREADPRIVATE(niterdivgrad) REAL,ALLOCATABLE,SAVE :: tau_graddiv(:) !$OMP THREADPRIVATE(tau_graddiv) REAL,ALLOCATABLE,SAVE :: tau_gradrot(:) !$OMP THREADPRIVATE(tau_gradrot) REAL,ALLOCATABLE,SAVE :: tau_divgrad(:) !$OMP THREADPRIVATE(tau_divgrad) REAL,SAVE :: cgraddiv !$OMP THREADPRIVATE(cgraddiv) REAL,SAVE :: cgradrot !$OMP THREADPRIVATE(cgradrot) REAL,SAVE :: cdivgrad !$OMP THREADPRIVATE(cdivgrad) INTEGER, SAVE :: rayleigh_friction_type, rayleigh_shear !$OMP THREADPRIVATE(rayleigh_friction_type) REAL, SAVE :: rayleigh_tau !$OMP THREADPRIVATE(rayleigh_shear) REAL,SAVE :: dtdissip !$OMP THREADPRIVATE(dtdissip) PUBLIC init_dissip, dissip CONTAINS SUBROUTINE allocate_dissip USE icosa IMPLICIT NONE CALL allocate_field(f_due_diss1,field_u,type_real,llm) CALL allocate_field(f_due_diss2,field_u,type_real,llm) CALL allocate_field(f_theta,field_t,type_real,llm) CALL allocate_field(f_dtheta_diss,field_t,type_real,llm) CALL allocate_field(f_dtheta_rhodz_diss,field_t,type_real,llm) CALL allocate_field(f_phi,field_t,type_real,llm) CALL allocate_field(f_pk,field_t,type_real,llm) CALL allocate_field(f_p,field_t,type_real,llm+1) CALL allocate_field(f_pks,field_t,type_real) ALLOCATE(tau_graddiv(llm)) ALLOCATE(tau_gradrot(llm)) ALLOCATE(tau_divgrad(llm)) END SUBROUTINE allocate_dissip SUBROUTINE init_dissip USE icosa USE disvert_mod USE mpi_mod USE mpipara USE transfert_mod USE time_mod USE transfert_omp_mod IMPLICIT NONE TYPE(t_field),POINTER,SAVE :: f_u(:) TYPE(t_field),POINTER,SAVE :: f_du(:) REAL(rstd),POINTER :: u(:) REAL(rstd),POINTER :: du(:) TYPE(t_field),POINTER,SAVE :: f_theta(:) TYPE(t_field),POINTER ,SAVE :: f_dtheta(:) REAL(rstd),POINTER :: theta(:) REAL(rstd),POINTER :: dtheta(:) REAL(rstd) :: dumax,dumax1 REAL(rstd) :: dthetamax,dthetamax1 REAL(rstd) :: r REAL(rstd) :: tau REAL(rstd) :: zz, zvert, fact INTEGER :: l CHARACTER(len=255) :: rayleigh_friction_key REAL(rstd) :: mintau INTEGER :: seed_size INTEGER,ALLOCATABLE :: seed(:) INTEGER :: i,j,ij,ind,it,iter rayleigh_friction_key='none' CALL getin("rayleigh_friction_type",rayleigh_friction_key) SELECT CASE(TRIM(rayleigh_friction_key)) CASE('none') rayleigh_friction_type=0 IF (is_mpi_root) PRINT *, 'No Rayleigh friction' CASE('dcmip2_schaer_noshear') rayleigh_friction_type=1 rayleigh_shear=0 IF (is_mpi_root) PRINT *, 'Rayleigh friction : Schaer-like mountain without shear DCMIP2.1' CASE('dcmip2_schaer_shear') rayleigh_shear=1 rayleigh_friction_type=2 IF (is_mpi_root) PRINT *, 'Rayleigh friction : Schaer-like mountain with shear DCMIP2.2' CASE DEFAULT IF (is_mpi_root) PRINT *, 'Bad selector : rayleigh_friction_type =', TRIM(rayleigh_friction_key), ' in dissip_gcm.f90/init_dissip' STOP END SELECT IF(rayleigh_friction_type>0) THEN rayleigh_tau=0. CALL getin("rayleigh_friction_tau",rayleigh_tau) rayleigh_tau = rayleigh_tau / scale_factor IF(rayleigh_tau<=0) THEN IF (is_mpi_root) PRINT *, 'Forbidden : negative value for rayleigh_friction_tau =',rayleigh_tau STOP END IF END IF CALL allocate_dissip CALL allocate_field(f_u,field_u,type_real) CALL allocate_field(f_du,field_u,type_real) CALL allocate_field(f_theta,field_t,type_real) CALL allocate_field(f_dtheta,field_t,type_real) CALL init_message(f_due_diss1,req_e1_vect,req_due) CALL init_message(f_dtheta_diss,req_i1,req_dtheta) tau_graddiv(:)=5000 CALL getin("tau_graddiv",tau) tau_graddiv(:)=tau/scale_factor CALL getin("nitergdiv",nitergdiv) tau_gradrot(:)=5000 CALL getin("tau_gradrot",tau) tau_gradrot(:)=tau/scale_factor CALL getin("nitergrot",nitergrot) tau_divgrad(:)=5000 CALL getin("tau_divgrad",tau) tau_divgrad(:)=tau/scale_factor CALL getin("niterdivgrad",niterdivgrad) cgraddiv=-1 cdivgrad=-1 cgradrot=-1 DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i CALL RANDOM_NUMBER(r) u(ij+u_right)=r-0.5 CALL RANDOM_NUMBER(r) u(ij+u_lup)=r-0.5 CALL RANDOM_NUMBER(r) u(ij+u_ldown)=r-0.5 ENDDO ENDDO ENDDO DO it=1,20 dumax=0 DO iter=1,nitergdiv CALL transfert_request(f_u,req_e1_vect) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) CALL compute_gradiv(u,du,1,1) u=du ENDDO ENDDO CALL transfert_request(f_du,req_e1_vect) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i if (le(ij+u_right)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_right))) if (le(ij+u_lup)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_lup))) if (le(ij+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_ldown))) ENDDO ENDDO ENDDO IF (using_mpi) THEN CALL reduce_sum_omp(dumax,dumax1) !$OMP MASTER CALL MPI_ALLREDUCE(dumax1,dumax,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr) !$OMP END MASTER CALL bcast_omp(dumax) ELSE CALL allreduce_sum_omp(dumax,dumax1) dumax=dumax1 ENDIF DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) u=du/dumax ENDDO IF (is_mpi_root) PRINT *,"gradiv : it :",it ,": dumax",dumax ENDDO IF (is_mpi_root) PRINT *,"gradiv : dumax",dumax IF (is_mpi_root) PRINT *, 'mean T-cell edge size (km)', 1.45*radius/iim_glo/1000., & 'effective T-cell half-edge size (km)', dumax**(-.5/nitergdiv)/1000 IF (is_mpi_root) PRINT *, 'Max. time step assuming c=340 m/s and Courant number=2.8 :', & 2.8/340.*dumax**(-.5/nitergdiv) cgraddiv=dumax**(-1./nitergdiv) IF (is_mpi_root) PRINT *,"cgraddiv : ",cgraddiv DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i CALL RANDOM_NUMBER(r) u(ij+u_right)=r-0.5 CALL RANDOM_NUMBER(r) u(ij+u_lup)=r-0.5 CALL RANDOM_NUMBER(r) u(ij+u_ldown)=r-0.5 ENDDO ENDDO ENDDO DO it=1,20 dumax=0 DO iter=1,nitergrot CALL transfert_request(f_u,req_e1_vect) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) CALL compute_gradrot(u,du,1,1) u=du ENDDO ENDDO CALL transfert_request(f_du,req_e1_vect) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i if (le(ij+u_right)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_right))) if (le(ij+u_lup)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_lup))) if (le(ij+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(ij+u_ldown))) ENDDO ENDDO ENDDO IF (using_mpi) THEN CALL reduce_sum_omp(dumax,dumax1) !$OMP MASTER CALL MPI_ALLREDUCE(dumax1,dumax,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr) !$OMP END MASTER CALL bcast_omp(dumax) ELSE CALL allreduce_sum_omp(dumax,dumax1) dumax=dumax1 ENDIF DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) u=du/dumax ENDDO IF (is_mpi_root) PRINT *,"gradrot : it :",it ,": dumax",dumax ENDDO IF (is_mpi_root) PRINT *,"gradrot : dumax",dumax cgradrot=dumax**(-1./nitergrot) IF (is_mpi_root) PRINT *,"cgradrot : ",cgradrot DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i CALL RANDOM_NUMBER(r) theta(ij)=r-0.5 ENDDO ENDDO ENDDO DO it=1,20 dthetamax=0 DO iter=1,niterdivgrad CALL transfert_request(f_theta,req_i1) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) dtheta=f_dtheta(ind) CALL compute_divgrad(theta,dtheta,1,1) theta=dtheta ENDDO ENDDO CALL transfert_request(f_dtheta,req_i1) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) dtheta=f_dtheta(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i dthetamax=MAX(dthetamax,ABS(dtheta(ij))) ENDDO ENDDO ENDDO IF (using_mpi) THEN CALL reduce_sum_omp(dthetamax ,dthetamax1) !$OMP MASTER CALL MPI_ALLREDUCE(dthetamax1,dthetamax,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr) !$OMP END MASTER CALL bcast_omp(dthetamax) ELSE CALL allreduce_sum_omp(dthetamax,dthetamax1) dumax=dumax1 ENDIF IF (is_mpi_root) PRINT *,"divgrad : it :",it ,": dthetamax",dthetamax DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) dtheta=f_dtheta(ind) theta=dtheta/dthetamax ENDDO ENDDO IF (is_mpi_root) PRINT *,"divgrad : divgrad",dthetamax cdivgrad=dthetamax**(-1./niterdivgrad) IF (is_mpi_root) PRINT *,"cdivgrad : ",cdivgrad fact=2 DO l=1,llm IF(ap_bp_present) THEN ! height-dependent dissipation zz= 1. - preff/presnivs(l) ELSE zz = 0. END IF zvert=fact-(fact-1)/(1+zz*zz) tau_graddiv(l) = tau_graddiv(l)/zvert tau_gradrot(l) = tau_gradrot(l)/zvert tau_divgrad(l) = tau_divgrad(l)/zvert ENDDO mintau=tau_graddiv(1) DO l=1,llm mintau=MIN(mintau,tau_graddiv(l)) mintau=MIN(mintau,tau_gradrot(l)) mintau=MIN(mintau,tau_divgrad(l)) ENDDO IF(mintau>0) THEN itau_dissip=INT(mintau/dt) dtdissip=itau_dissip*dt ELSE IF (is_mpi_root) PRINT *,"No dissipation time set, setting itau_dissip to 1000000000" itau_dissip=100000000 END IF itau_dissip=MAX(1,itau_dissip) IF (is_mpi_root) PRINT *,"mintau ",mintau,"itau_dissip",itau_dissip," dtdissip ",dtdissip END SUBROUTINE init_dissip SUBROUTINE dissip(f_ue,f_due,f_mass,f_phis,f_theta_rhodz,f_dtheta_rhodz) USE icosa USE theta2theta_rhodz_mod USE pression_mod USE exner_mod USE geopotential_mod USE trace USE time_mod USE omp_para IMPLICIT NONE TYPE(t_field),POINTER :: f_ue(:) TYPE(t_field),POINTER :: f_due(:) TYPE(t_field),POINTER :: f_mass(:), f_phis(:) TYPE(t_field),POINTER :: f_theta_rhodz(:) TYPE(t_field),POINTER :: f_dtheta_rhodz(:) REAL(rstd),POINTER :: due(:,:) REAL(rstd),POINTER :: phi(:,:), ue(:,:) REAL(rstd),POINTER :: due_diss1(:,:) REAL(rstd),POINTER :: due_diss2(:,:) REAL(rstd),POINTER :: dtheta_rhodz(:,:) REAL(rstd),POINTER :: dtheta_rhodz_diss(:,:) INTEGER :: ind, shear INTEGER :: l,ij !$OMP BARRIER CALL trace_start("dissip") CALL gradiv(f_ue,f_due_diss1) CALL gradrot(f_ue,f_due_diss2) CALL divgrad_theta_rhodz(f_mass,f_theta_rhodz,f_dtheta_rhodz_diss) ! later for openmp ! IF(rayleigh_friction_type>0) THEN ! CALL pression(f_ps, f_p) ! CALL exner(f_ps, f_p, f_pks, f_pk) ! CALL geopotential(f_phis,f_pks,f_pk,f_theta,f_phi) ! END IF DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) due=f_due(ind) due_diss1=f_due_diss1(ind) due_diss2=f_due_diss2(ind) dtheta_rhodz=f_dtheta_rhodz(ind) dtheta_rhodz_diss=f_dtheta_rhodz_diss(ind) DO l=ll_begin,ll_end !$SIMD DO ij=ij_begin,ij_end 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 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 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 dtheta_rhodz(ij,l) = -0.5*dtheta_rhodz_diss(ij,l)/tau_divgrad(l)*itau_dissip ENDDO ENDDO ! dtheta_rhodz=0 ! due=0 ! later for openmp ! IF(rayleigh_friction_type>0) THEN ! phi=f_phi(ind) ! ue=f_ue(ind) ! DO l=1,llm ! DO j=jj_begin,jj_end ! DO i=ii_begin,ii_end ! n=(j-1)*iim+i ! CALL relax(t_right, u_right) ! CALL relax(t_lup, u_lup) ! CALL relax(t_ldown, u_ldown) ! ENDDO ! ENDDO ! END DO ! END IF END DO CALL trace_end("dissip") !$OMP BARRIER CONTAINS SUBROUTINE relax(shift_t, shift_u) USE dcmip_initial_conditions_test_1_2_3 REAL(rstd) :: z, ulon,ulat, lon,lat, & ! input to test2_schaer_mountain p,hyam,hybm,w,t,phis,ps,rho,q, & ! unused input/output to test2_schaer_mountain fz, u3d(3), uref REAL(rstd), PARAMETER :: zh=2e4,ztop=3e4 ! DCMIP values LOGICAL :: hybrid_eta INTEGER :: shift_u, shift_t, zcoords, nn z = (phi(ij,l)+phi(ij+shift_t,l))/(2.*g) IF(z>zh) THEN ! relax only in the sponge layer z>zh nn = ij+shift_u CALL xyz2lonlat(xyz_e(nn,:),lon,lat) zcoords = 1 ; hybrid_eta = .FALSE. ! use z instead of p or hyam/hybm CALL test2_schaer_mountain(lon,lat,p,z,zcoords,hybrid_eta, & hyam,hybm,shear,ulon,ulat,w,t,phis,ps,rho,q) ! u3d = ulon*elon_e(nn,:) + ulat*elat_e(nn,:) u3d = ulon*elon_e(nn,:) ! ulat=0 uref = sum(u3d*ep_e(nn,:)) fz = sin((pi/2)*(z-zh)/(ztop-zh)) fz = fz*fz/rayleigh_tau ! fz = 1/rayleigh_tau due(nn,l) = due(nn,l) - fz*(ue(nn,l)-uref) ! due(nn,l) = due(nn,l) - fz*ue(nn,l) END IF END SUBROUTINE relax END SUBROUTINE dissip SUBROUTINE gradiv(f_ue,f_due) USE icosa USE trace USE omp_para IMPLICIT NONE TYPE(t_field),POINTER :: f_ue(:) TYPE(t_field),POINTER :: f_due(:) REAL(rstd),POINTER :: ue(:,:) REAL(rstd),POINTER :: due(:,:) INTEGER :: ind INTEGER :: it,l,ij CALL trace_start("gradiv") DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) ue=f_ue(ind) due=f_due(ind) DO l = ll_begin, ll_end !$SIMD DO ij=ij_begin,ij_end due(ij+u_right,l)=ue(ij+u_right,l) due(ij+u_lup,l)=ue(ij+u_lup,l) due(ij+u_ldown,l)=ue(ij+u_ldown,l) ENDDO ENDDO ENDDO DO it=1,nitergdiv CALL send_message(f_due,req_due) CALL wait_message(req_due) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) due=f_due(ind) CALL compute_gradiv(due,due,ll_begin,ll_end) ENDDO ENDDO CALL trace_end("gradiv") END SUBROUTINE gradiv SUBROUTINE gradrot(f_ue,f_due) USE icosa USE trace USE omp_para IMPLICIT NONE TYPE(t_field),POINTER :: f_ue(:) TYPE(t_field),POINTER :: f_due(:) REAL(rstd),POINTER :: ue(:,:) REAL(rstd),POINTER :: due(:,:) INTEGER :: ind INTEGER :: it,l,ij CALL trace_start("gradrot") DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) ue=f_ue(ind) due=f_due(ind) DO l = ll_begin, ll_end !$SIMD DO ij=ij_begin,ij_end due(ij+u_right,l)=ue(ij+u_right,l) due(ij+u_lup,l)=ue(ij+u_lup,l) due(ij+u_ldown,l)=ue(ij+u_ldown,l) ENDDO ENDDO ENDDO DO it=1,nitergrot CALL send_message(f_due,req_due) CALL wait_message(req_due) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) due=f_due(ind) CALL compute_gradrot(due,due,ll_begin,ll_end) ENDDO ENDDO CALL trace_end("gradrot") END SUBROUTINE gradrot SUBROUTINE divgrad(f_theta,f_dtheta) USE icosa USE trace USE omp_para IMPLICIT NONE TYPE(t_field),POINTER :: f_theta(:) TYPE(t_field),POINTER :: f_dtheta(:) REAL(rstd),POINTER :: theta(:,:) REAL(rstd),POINTER :: dtheta(:,:) INTEGER :: ind INTEGER :: it CALL trace_start("divgrad") DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) dtheta=f_dtheta(ind) dtheta=theta ENDDO DO it=1,niterdivgrad CALL transfert_request(f_dtheta,req_i1) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) dtheta=f_dtheta(ind) CALL compute_divgrad(dtheta,dtheta,ll_begin,ll_end) ENDDO ENDDO CALL trace_end("divgrad") END SUBROUTINE divgrad SUBROUTINE divgrad_theta_rhodz(f_mass,f_theta_rhodz,f_dtheta_rhodz) USE icosa USE trace USE omp_para IMPLICIT NONE TYPE(t_field),POINTER :: f_mass(:) TYPE(t_field),POINTER :: f_theta_rhodz(:) TYPE(t_field),POINTER :: f_dtheta_rhodz(:) REAL(rstd),POINTER :: mass(:,:) REAL(rstd),POINTER :: theta_rhodz(:,:) REAL(rstd),POINTER :: dtheta_rhodz(:,:) INTEGER :: ind INTEGER :: it,l,ij CALL trace_start("divgrad") DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) mass=f_mass(ind) theta_rhodz=f_theta_rhodz(ind) dtheta_rhodz=f_dtheta_rhodz(ind) DO l = ll_begin, ll_end !$SIMD DO ij=ij_begin,ij_end dtheta_rhodz(ij,l) = theta_rhodz(ij,l) / mass(ij,l) ENDDO ENDDO ENDDO DO it=1,niterdivgrad CALL send_message(f_dtheta_rhodz,req_dtheta) CALL wait_message(req_dtheta) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) dtheta_rhodz=f_dtheta_rhodz(ind) CALL compute_divgrad(dtheta_rhodz,dtheta_rhodz,ll_begin,ll_end) ENDDO ENDDO DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) dtheta_rhodz=f_dtheta_rhodz(ind) mass=f_mass(ind) DO l = ll_begin, ll_end !$SIMD DO ij=ij_begin,ij_end dtheta_rhodz(ij,l) = dtheta_rhodz(ij,l) * mass(ij,l) ENDDO ENDDO ENDDO CALL trace_end("divgrad") END SUBROUTINE divgrad_theta_rhodz SUBROUTINE compute_gradiv(ue,gradivu_e,llb,lle) USE icosa IMPLICIT NONE INTEGER,INTENT(IN) :: llb INTEGER,INTENT(IN) :: lle REAL(rstd),INTENT(IN) :: ue(iim*3*jjm,llm) REAL(rstd),INTENT(OUT) :: gradivu_e(iim*3*jjm,llm) REAL(rstd) :: divu_i(iim*jjm,llb:lle) INTEGER :: ij,l DO l=llb,lle !$SIMD DO ij=ij_begin,ij_end divu_i(ij,l)=1./Ai(ij)*(ne(ij,right)*ue(ij+u_right,l)*le(ij+u_right) + & ne(ij,rup)*ue(ij+u_rup,l)*le(ij+u_rup) + & ne(ij,lup)*ue(ij+u_lup,l)*le(ij+u_lup) + & ne(ij,left)*ue(ij+u_left,l)*le(ij+u_left) + & ne(ij,ldown)*ue(ij+u_ldown,l)*le(ij+u_ldown) + & ne(ij,rdown)*ue(ij+u_rdown,l)*le(ij+u_rdown)) ENDDO ENDDO DO l=llb,lle !$SIMD DO ij=ij_begin,ij_end 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) ) 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)) 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) ) ENDDO ENDDO DO l=llb,lle !$SIMD DO ij=ij_begin,ij_end gradivu_e(ij+u_right,l)=-gradivu_e(ij+u_right,l)*cgraddiv gradivu_e(ij+u_lup,l)=-gradivu_e(ij+u_lup,l)*cgraddiv gradivu_e(ij+u_ldown,l)=-gradivu_e(ij+u_ldown,l)*cgraddiv ENDDO ENDDO END SUBROUTINE compute_gradiv SUBROUTINE compute_divgrad(theta,divgrad_i,llb,lle) USE icosa IMPLICIT NONE INTEGER,INTENT(IN) :: llb INTEGER,INTENT(IN) :: lle REAL(rstd),INTENT(IN) :: theta(iim*jjm,llm) REAL(rstd),INTENT(OUT) :: divgrad_i(iim*jjm,llm) REAL(rstd) :: grad_e(3*iim*jjm,llb:lle) INTEGER :: ij,l DO l=llb,lle !$SIMD DO ij=ij_begin_ext,ij_end_ext 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) ) 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 )) 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) ) ENDDO ENDDO DO l=llb,lle !$SIMD DO ij=ij_begin,ij_end divgrad_i(ij,l)=1./Ai(ij)*(ne(ij,right)*grad_e(ij+u_right,l)*le(ij+u_right) + & ne(ij,rup)*grad_e(ij+u_rup,l)*le(ij+u_rup) + & ne(ij,lup)*grad_e(ij+u_lup,l)*le(ij+u_lup) + & ne(ij,left)*grad_e(ij+u_left,l)*le(ij+u_left) + & ne(ij,ldown)*grad_e(ij+u_ldown,l)*le(ij+u_ldown) + & ne(ij,rdown)*grad_e(ij+u_rdown,l)*le(ij+u_rdown)) ENDDO ENDDO DO l=llb,lle DO ij=ij_begin,ij_end divgrad_i(ij,l)=-divgrad_i(ij,l)*cdivgrad ENDDO ENDDO END SUBROUTINE compute_divgrad SUBROUTINE compute_gradrot(ue,gradrot_e,llb,lle) USE icosa IMPLICIT NONE INTEGER,INTENT(IN) :: llb INTEGER,INTENT(IN) :: lle REAL(rstd),INTENT(IN) :: ue(iim*3*jjm,llm) REAL(rstd),INTENT(OUT) :: gradrot_e(iim*3*jjm,llm) REAL(rstd) :: rot_v(2*iim*jjm,llb:lle) INTEGER :: ij,l DO l=llb,lle !$SIMD DO ij=ij_begin_ext,ij_end_ext rot_v(ij+z_up,l)= 1./Av(ij+z_up)*( ne(ij,rup)*ue(ij+u_rup,l)*de(ij+u_rup) & + ne(ij+t_rup,left)*ue(ij+t_rup+u_left,l)*de(ij+t_rup+u_left) & - ne(ij,lup)*ue(ij+u_lup,l)*de(ij+u_lup) ) rot_v(ij+z_down,l) = 1./Av(ij+z_down)*( ne(ij,ldown)*ue(ij+u_ldown,l)*de(ij+u_ldown) & + ne(ij+t_ldown,right)*ue(ij+t_ldown+u_right,l)*de(ij+t_ldown+u_right) & - ne(ij,rdown)*ue(ij+u_rdown,l)*de(ij+u_rdown) ) ENDDO ENDDO DO l=llb,lle !$SIMD DO ij=ij_begin,ij_end 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)) 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)) 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)) ENDDO ENDDO DO l=llb,lle !$SIMD DO ij=ij_begin,ij_end gradrot_e(ij+u_right,l)=-gradrot_e(ij+u_right,l)*cgradrot gradrot_e(ij+u_lup,l)=-gradrot_e(ij+u_lup,l)*cgradrot gradrot_e(ij+u_ldown,l)=-gradrot_e(ij+u_ldown,l)*cgradrot ENDDO ENDDO END SUBROUTINE compute_gradrot END MODULE dissip_gcm_mod