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(:) INTEGER,SAVE :: nitergdiv=1 INTEGER,SAVE :: nitergrot=1 INTEGER,SAVE :: niterdivgrad=1 REAL,ALLOCATABLE,SAVE :: tau_graddiv(:) REAL,ALLOCATABLE,SAVE :: tau_gradrot(:) REAL,ALLOCATABLE,SAVE :: tau_divgrad(:) REAL,SAVE :: cgraddiv REAL,SAVE :: cgradrot REAL,SAVE :: cdivgrad INTEGER, SAVE :: rayleigh_friction_type, rayleigh_shear REAL, save :: rayleigh_tau INTEGER,SAVE :: idissip REAL,SAVE :: 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 IMPLICIT NONE TYPE(t_field),POINTER :: f_u(:) TYPE(t_field),POINTER :: f_du(:) REAL(rstd),POINTER :: u(:) REAL(rstd),POINTER :: du(:) TYPE(t_field),POINTER :: f_theta(:) TYPE(t_field),POINTER :: 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 INTEGER :: i,j,n,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) 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_gradrot) 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) ! CALL create_request(field_u,req_dissip) ! DO ind=1,ndomain ! DO i=ii_begin,ii_end ! CALL request_add_point(ind,i,jj_begin-1,req_dissip,rup) ! CALL request_add_point(ind,i+1,jj_begin-1,req_dissip,lup) ! ENDDO ! DO j=jj_begin,jj_end ! CALL request_add_point(ind,ii_end+1,j,req_dissip,left) ! CALL request_add_point(ind,ii_end+1,j-1,req_dissip,lup) ! ENDDO ! ! DO i=ii_begin,ii_end ! CALL request_add_point(ind,i,jj_end+1,req_dissip,ldown) ! CALL request_add_point(ind,i-1,jj_end+1,req_dissip,rdown) ! ENDDO ! DO j=jj_begin,jj_end ! CALL request_add_point(ind,ii_begin-1,j,req_dissip,right) ! CALL request_add_point(ind,ii_begin-1,j+1,req_dissip,rdown) ! ENDDO ! DO i=ii_begin+1,ii_end-1 ! CALL request_add_point(ind,i,jj_begin,req_dissip,right) ! CALL request_add_point(ind,i,jj_end,req_dissip,right) ! ENDDO ! ! DO j=jj_begin+1,jj_end-1 ! CALL request_add_point(ind,ii_begin,j,req_dissip,rup) ! CALL request_add_point(ind,ii_end,j,req_dissip,rup) ! ENDDO ! CALL request_add_point(ind,ii_begin+1,jj_begin,req_dissip,left) ! CALL request_add_point(ind,ii_begin,jj_begin+1,req_dissip,ldown) ! CALL request_add_point(ind,ii_begin+1,jj_end,req_dissip,left) ! CALL request_add_point(ind,ii_end,jj_begin+1,req_dissip,ldown) ! ! ENDDO cgraddiv=-1 cdivgrad=-1 cgradrot=-1 CALL RANDOM_SEED() DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i CALL RANDOM_NUMBER(r) u(n+u_right)=r-0.5 CALL RANDOM_NUMBER(r) u(n+u_lup)=r-0.5 CALL RANDOM_NUMBER(r) u(n+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) DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) CALL compute_gradiv(u,du,1) u=du ENDDO ENDDO CALL transfert_request(f_du,req_e1) DO ind=1,ndomain 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 n=(j-1)*iim+i if (le(n+u_right)>1e-100) dumax=MAX(dumax,ABS(du(n+u_right))) if (le(n+u_lup)>1e-100) dumax=MAX(dumax,ABS(du(n+u_lup))) if (le(n+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(n+u_ldown))) ENDDO ENDDO ENDDO IF (using_mpi) THEN CALL MPI_ALLREDUCE(dumax,dumax1,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr) dumax=dumax1 ENDIF DO ind=1,ndomain 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 CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i CALL RANDOM_NUMBER(r) u(n+u_right)=r-0.5 CALL RANDOM_NUMBER(r) u(n+u_lup)=r-0.5 CALL RANDOM_NUMBER(r) u(n+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) DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) u=f_u(ind) du=f_du(ind) CALL compute_gradrot(u,du,1) u=du ENDDO ENDDO CALL transfert_request(f_du,req_e1) DO ind=1,ndomain 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 n=(j-1)*iim+i if (le(n+u_right)>1e-100) dumax=MAX(dumax,ABS(du(n+u_right))) if (le(n+u_lup)>1e-100) dumax=MAX(dumax,ABS(du(n+u_lup))) if (le(n+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(n+u_ldown))) ENDDO ENDDO ENDDO IF (using_mpi) THEN CALL MPI_ALLREDUCE(dumax,dumax1,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr) dumax=dumax1 ENDIF DO ind=1,ndomain 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 CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i CALL RANDOM_NUMBER(r) theta(n)=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 CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) dtheta=f_dtheta(ind) CALL compute_divgrad(theta,dtheta,1) theta=dtheta ENDDO ENDDO ! CALL writefield("divgrad",f_dtheta) CALL transfert_request(f_dtheta,req_i1) DO ind=1,ndomain 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 n=(j-1)*iim+i dthetamax=MAX(dthetamax,ABS(dtheta(n))) ENDDO ENDDO ENDDO IF (using_mpi) THEN CALL MPI_ALLREDUCE(dthetamax,dthetamax1,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr) dthetamax=dthetamax1 ENDIF IF (is_mpi_root) PRINT *,"divgrad : it :",it ,": dthetamax",dthetamax DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) theta=f_theta(ind) dtheta=f_dtheta(ind) theta=dtheta/dthetamax ENDDO ENDDO ! CALL writefield("divgrad",f_dtheta) 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 zz= 1. - preff/presnivs(l) 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 ! idissip=INT(MIN(tetagdiv,tetagrot)/(2.*dt)) ! idissip=MAX(1,idissip) ! dtdissip=idissip*dt ! PRINT *,"idissip",idissip," dtdissip ",dtdissip END SUBROUTINE init_dissip SUBROUTINE dissip(f_ue,f_due,f_ps,f_phis,f_theta_rhodz,f_dtheta_rhodz) USE icosa USE theta2theta_rhodz_mod USE pression_mod USE exner_mod USE geopotential_mod IMPLICIT NONE TYPE(t_field),POINTER :: f_ue(:) TYPE(t_field),POINTER :: f_due(:) TYPE(t_field),POINTER :: f_ps(:), 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,i,j,n CALL gradiv(f_ue,f_due_diss1) CALL gradrot(f_ue,f_due_diss2) CALL theta_rhodz2theta(f_ps,f_theta_rhodz,f_theta) CALL divgrad(f_theta,f_dtheta_diss) CALL theta2theta_rhodz(f_ps,f_dtheta_diss,f_dtheta_rhodz_diss) 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 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=1,llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i due(n+u_right,l) = -0.5*( due_diss1(n+u_right,l)/tau_graddiv(l) + due_diss2(n+u_right,l)/tau_gradrot(l) ) due(n+u_lup,l) = -0.5*( due_diss1(n+u_lup,l) /tau_graddiv(l) + due_diss2(n+u_lup,l) /tau_gradrot(l) ) due(n+u_ldown,l) = -0.5*( due_diss1(n+u_ldown,l)/tau_graddiv(l) + due_diss2(n+u_ldown,l)/tau_gradrot(l) ) dtheta_rhodz(n,l) = -0.5*dtheta_rhodz_diss(n,l)/tau_divgrad(l) ENDDO ENDDO ENDDO 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 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(n,l)+phi(n+shift_t,l))/(2.*g) IF(z>zh) THEN ! relax only in the sponge layer z>zh ! PRINT *, 'z>zh : z,zh,l',z,zh,l ! STOP nn = n+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 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 DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) ue=f_ue(ind) due=f_due(ind) due=ue ENDDO DO it=1,nitergdiv CALL transfert_request(f_due,req_e1) DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) due=f_due(ind) CALL compute_gradiv(due,due,llm) ENDDO ENDDO END SUBROUTINE gradiv SUBROUTINE gradrot(f_ue,f_due) USE icosa 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 DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) ue=f_ue(ind) due=f_due(ind) due=ue ENDDO DO it=1,nitergrot CALL transfert_request(f_due,req_e1) DO ind=1,ndomain CALL swap_dimensions(ind) CALL swap_geometry(ind) due=f_due(ind) CALL compute_gradrot(due,due,llm) ENDDO ENDDO END SUBROUTINE gradrot SUBROUTINE divgrad(f_theta,f_dtheta) USE icosa 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 DO ind=1,ndomain 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 CALL swap_dimensions(ind) CALL swap_geometry(ind) dtheta=f_dtheta(ind) CALL compute_divgrad(dtheta,dtheta,llm) ENDDO ENDDO END SUBROUTINE divgrad ! SUBROUTINE compute_dissip(ue,due,ps,theta_rhodz,dtheta_rhodz) ! USE icosa ! USE theta2theta_rhodz_mod ! IMPLICIT NONE ! REAL(rstd) :: ue(3*iim*jjm,llm) ! REAL(rstd) :: due(3*iim*jjm,llm) ! REAL(rstd) :: ps(iim*jjm) ! REAL(rstd) :: theta_rhodz(iim*jjm,llm) ! REAL(rstd) :: dtheta_rhodz(iim*jjm,llm) ! ! REAL(rstd),SAVE,ALLOCATABLE :: theta(:,:) ! REAL(rstd),SAVE,ALLOCATABLE :: du_dissip(:,:) ! REAL(rstd),SAVE,ALLOCATABLE :: dtheta_dissip(:,:) ! REAL(rstd),SAVE,ALLOCATABLE :: dtheta_rhodz_dissip(:,:) ! ! INTEGER :: ind ! INTEGER :: l,i,j,n ! !!$OMP BARRIER !!$OMP MASTER ! ALLOCATE(theta(iim*jjm,llm)) ! ALLOCATE(du_dissip(3*iim*jjm,llm)) ! ALLOCATE(dtheta_dissip(iim*jjm,llm)) ! ALLOCATE(dtheta_rhodz_dissip(iim*jjm,llm)) !!$OMP END MASTER !!$OMP BARRIER ! ! CALL gradiv(ue,du_dissip,llm) ! DO l=1,llm !!$OMP DO ! DO j=jj_begin,jj_end ! DO i=ii_begin,ii_end ! n=(j-1)*iim+i ! due(n+u_right,l)=due(n+u_right,l)+du_dissip(n+u_right,l)/tau_graddiv(l)*0.5 ! due(n+u_lup,l)=due(n+u_lup,l)+du_dissip(n+u_lup,l)/tau_graddiv(l)*0.5 ! due(n+u_ldown,l)=due(n+u_ldown,l)+du_dissip(n+u_ldown,l)/tau_graddiv(l)*0.5 ! ENDDO ! ENDDO ! ENDDO ! ! CALL gradrot(ue,du_dissip,llm) ! ! DO l=1,llm !!$OMP DO ! DO j=jj_begin,jj_end ! DO i=ii_begin,ii_end ! n=(j-1)*iim+i ! due(n+u_right,l)=due(n+u_right,l)+du_dissip(n+u_right,l)/tau_gradrot(l)*0.5 ! due(n+u_lup,l)=due(n+u_lup,l)+du_dissip(n+u_lup,l)/tau_gradrot(l)*0.5 ! due(n+u_ldown,l)=due(n+u_ldown,l)+du_dissip(n+u_ldown,l)/tau_gradrot(l)*0.5 ! ENDDO ! ENDDO ! ENDDO ! ! CALL compute_theta_rhodz2theta(ps,theta_rhodz,theta,1) ! CALL divgrad(theta,dtheta_dissip,llm) ! CALL compute_theta2theta_rhodz(ps,dtheta_dissip,dtheta_rhodz_dissip,0) ! ! DO l=1,llm !!$OMP DO ! DO j=jj_begin,jj_end ! DO i=ii_begin,ii_end ! n=(j-1)*iim+i ! dtheta_rhodz(n,l)=dtheta_rhodz(n,l)+dtheta_rhodz_dissip(n,l)/tau_divgrad(l)*0.5 ! ENDDO ! ENDDO ! ENDDO ! !!$OMP BARRIER !!$OMP MASTER ! DEALLOCATE(theta) ! DEALLOCATE(du_dissip) ! DEALLOCATE(dtheta_dissip) ! DEALLOCATE(dtheta_rhodz_dissip) !!$OMP END MASTER !!$OMP BARRIER ! ! END SUBROUTINE compute_dissip SUBROUTINE compute_gradiv(ue,gradivu_e,ll) USE icosa IMPLICIT NONE INTEGER,INTENT(IN) :: ll REAL(rstd),INTENT(IN) :: ue(iim*3*jjm,ll) REAL(rstd),INTENT(OUT) :: gradivu_e(iim*3*jjm,ll) REAL(rstd),SAVE,ALLOCATABLE :: divu_i(:,:) INTEGER :: i,j,n,l !$OMP BARRIER !$OMP MASTER ALLOCATE(divu_i(iim*jjm,ll)) !$OMP END MASTER !$OMP BARRIER DO l=1,ll !$OMP DO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i divu_i(n,l)=1./Ai(n)*(ne(n,right)*ue(n+u_right,l)*le(n+u_right) + & ne(n,rup)*ue(n+u_rup,l)*le(n+u_rup) + & ne(n,lup)*ue(n+u_lup,l)*le(n+u_lup) + & ne(n,left)*ue(n+u_left,l)*le(n+u_left) + & ne(n,ldown)*ue(n+u_ldown,l)*le(n+u_ldown) + & ne(n,rdown)*ue(n+u_rdown,l)*le(n+u_rdown)) ENDDO ENDDO ENDDO DO l=1,ll !$OMP DO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i gradivu_e(n+u_right,l)=-1/de(n+u_right)*(ne(n,right)*divu_i(n,l)+ ne(n+t_right,left)*divu_i(n+t_right,l) ) gradivu_e(n+u_lup,l)=-1/de(n+u_lup)*(ne(n,lup)*divu_i(n,l)+ ne(n+t_lup,rdown)*divu_i(n+t_lup,l)) gradivu_e(n+u_ldown,l)=-1/de(n+u_ldown)*(ne(n,ldown)*divu_i(n,l)+ne(n+t_ldown,rup)*divu_i(n+t_ldown,l) ) ENDDO ENDDO ENDDO DO l=1,ll !$OMP DO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i gradivu_e(n+u_right,l)=-gradivu_e(n+u_right,l)*cgraddiv gradivu_e(n+u_lup,l)=-gradivu_e(n+u_lup,l)*cgraddiv gradivu_e(n+u_ldown,l)=-gradivu_e(n+u_ldown,l)*cgraddiv ENDDO ENDDO ENDDO !$OMP BARRIER !$OMP MASTER DEALLOCATE(divu_i) !$OMP END MASTER !$OMP BARRIER END SUBROUTINE compute_gradiv SUBROUTINE compute_divgrad(theta,divgrad_i,ll) USE icosa IMPLICIT NONE INTEGER,INTENT(IN) :: ll REAL(rstd),INTENT(IN) :: theta(iim*jjm,ll) REAL(rstd),INTENT(OUT) :: divgrad_i(iim*jjm,ll) REAL(rstd),SAVE,ALLOCATABLE :: grad_e(:,:) INTEGER :: i,j,n,l !$OMP BARRIER !$OMP MASTER ALLOCATE(grad_e(3*iim*jjm,ll)) !$OMP END MASTER !$OMP BARRIER DO l=1,ll !$OMP DO DO j=jj_begin-1,jj_end+1 DO i=ii_begin-1,ii_end+1 n=(j-1)*iim+i grad_e(n+u_right,l)=-1/de(n+u_right)*(ne(n,right)*theta(n,l)+ ne(n+t_right,left)*theta(n+t_right,l) ) grad_e(n+u_lup,l)=-1/de(n+u_lup)*(ne(n,lup)*theta(n,l)+ ne(n+t_lup,rdown)*theta(n+t_lup,l )) grad_e(n+u_ldown,l)=-1/de(n+u_ldown)*(ne(n,ldown)*theta(n,l)+ne(n+t_ldown,rup)*theta(n+t_ldown,l) ) ENDDO ENDDO ENDDO DO l=1,ll !$OMP DO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i divgrad_i(n,l)=1./Ai(n)*(ne(n,right)*grad_e(n+u_right,l)*le(n+u_right) + & ne(n,rup)*grad_e(n+u_rup,l)*le(n+u_rup) + & ne(n,lup)*grad_e(n+u_lup,l)*le(n+u_lup) + & ne(n,left)*grad_e(n+u_left,l)*le(n+u_left) + & ne(n,ldown)*grad_e(n+u_ldown,l)*le(n+u_ldown) + & ne(n,rdown)*grad_e(n+u_rdown,l)*le(n+u_rdown)) ENDDO ENDDO ENDDO DO l=1,ll !$OMP DO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i divgrad_i(n,l)=-divgrad_i(n,l)*cdivgrad ENDDO ENDDO ENDDO !$OMP BARRIER !$OMP MASTER DEALLOCATE(grad_e) !$OMP END MASTER !$OMP BARRIER END SUBROUTINE compute_divgrad SUBROUTINE compute_gradrot(ue,gradrot_e,ll) USE icosa IMPLICIT NONE INTEGER,INTENT(IN) :: ll REAL(rstd),INTENT(IN) :: ue(iim*3*jjm,ll) REAL(rstd),INTENT(OUT) :: gradrot_e(iim*3*jjm,ll) REAL(rstd),SAVE,ALLOCATABLE :: rot_v(:,:) INTEGER :: i,j,n,l !$OMP BARRIER !$OMP MASTER ALLOCATE(rot_v(2*iim*jjm,ll)) !$OMP END MASTER !$OMP BARRIER DO l=1,ll !$OMP DO DO j=jj_begin-1,jj_end+1 DO i=ii_begin-1,ii_end+1 n=(j-1)*iim+i rot_v(n+z_up,l)= 1./Av(n+z_up)*( ne(n,rup)*ue(n+u_rup,l)*de(n+u_rup) & + ne(n+t_rup,left)*ue(n+t_rup+u_left,l)*de(n+t_rup+u_left) & - ne(n,lup)*ue(n+u_lup,l)*de(n+u_lup) ) rot_v(n+z_down,l) = 1./Av(n+z_down)*( ne(n,ldown)*ue(n+u_ldown,l)*de(n+u_ldown) & + ne(n+t_ldown,right)*ue(n+t_ldown+u_right,l)*de(n+t_ldown+u_right) & - ne(n,rdown)*ue(n+u_rdown,l)*de(n+u_rdown) ) ENDDO ENDDO ENDDO DO l=1,ll !$OMP DO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i gradrot_e(n+u_right,l)=1/le(n+u_right)*ne(n,right)*(rot_v(n+z_rdown,l)-rot_v(n+z_rup,l)) gradrot_e(n+u_lup,l)=1/le(n+u_lup)*ne(n,lup)*(rot_v(n+z_up,l)-rot_v(n+z_lup,l)) gradrot_e(n+u_ldown,l)=1/le(n+u_ldown)*ne(n,ldown)*(rot_v(n+z_ldown,l)-rot_v(n+z_down,l)) ENDDO ENDDO ENDDO DO l=1,ll !$OMP DO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i gradrot_e(n+u_right,l)=-gradrot_e(n+u_right,l)*cgradrot gradrot_e(n+u_lup,l)=-gradrot_e(n+u_lup,l)*cgradrot gradrot_e(n+u_ldown,l)=-gradrot_e(n+u_ldown,l)*cgradrot ENDDO ENDDO ENDDO !$OMP BARRIER !$OMP MASTER DEALLOCATE(rot_v) !$OMP END MASTER !$OMP BARRIER END SUBROUTINE compute_gradrot END MODULE dissip_gcm_mod