MODULE check_conserve_mod USE icosa IMPLICIT NONE PRIVATE TYPE(t_field),POINTER,SAVE :: f_pk(:), f_pks(:), f_p(:) TYPE(t_field),POINTER,SAVE :: f_vort(:) TYPE(t_field),POINTER,SAVE :: f_rhodz(:) PUBLIC init_check_conserve, check_conserve REAL(rstd),SAVE:: mtot0,ztot0,etot0,ang0,stot0,rmsv0 !$OMP THREADPRIVATE(mtot0,ztot0,etot0,ang0,stot0,rmsv0) REAL(rstd),SAVE:: etot,ang,stot,rmsv !$OMP THREADPRIVATE(etot,ang,stot,rmsv) REAL(rstd),SAVE:: ztot !$OMP THREADPRIVATE(ztot) CONTAINS !--------------------------------------------------------------------- SUBROUTINE init_check_conserve USE icosa IMPLICIT NONE 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) CALL allocate_field(f_rhodz,field_t,type_real,llm) CALL allocate_field(f_vort,field_z,type_real,llm) END SUBROUTINE init_check_conserve SUBROUTINE check_conserve(f_ps,f_dps,f_ue,f_theta_rhodz,f_phis,it) USE icosa USE pression_mod USE vorticity_mod USE caldyn_gcm_mod USE exner_mod USE mpipara, ONLY : is_mpi_root, comm_icosa IMPLICIT NONE TYPE(t_field),POINTER :: f_ps(:) TYPE(t_field),POINTER :: f_dps(:) TYPE(t_field),POINTER :: f_ue(:) TYPE(t_field),POINTER :: f_theta_rhodz(:) TYPE(t_field),POINTER :: f_phis(:) INTEGER::it REAL(rstd),POINTER :: p(:,:),rhodz(:,:) INTEGER::ind,ierr REAL(rstd) :: mtot, rmsdpdt etot=0.0; ang=0.0;stot=0.0;rmsv=0.0 ztot = 0.0 CALL pression(f_ps,f_p) DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) p=f_p(ind) rhodz=f_rhodz(ind) CALL compute_rhodz(p,rhodz) END DO CALL vorticity(f_ue,f_vort) CALL check_mass_conserve(f_ps,f_dps,mtot,rmsdpdt) CALL check_PV CALL exner(f_ps,f_p,f_pks,f_pk) CALL check_EN(f_ue,f_theta_rhodz,f_phis) IF (is_mpi_root) THEN !$OMP MASTER IF ( it == itau0 ) Then ztot0 = ztot mtot0 = mtot etot0 = etot ang0 = ang stot0 = stot END IF rmsv=SQRT(rmsv/mtot) ztot=ztot/ztot0-1. ; mtot=mtot/mtot0-1. etot=etot/etot0-1. ; ang=ang/ang0-1. ; stot=stot/stot0-1. rmsdpdt= daysec*1.e-2*sqrt(rmsdpdt/ncell_glo) OPEN(134,file="checkconsicosa.txt",position='append') WRITE(134,4000)mtot,rmsdpdt,etot,ztot,stot,rmsv,ang WRITE(134,*)mtot,rmsdpdt,etot,ztot,stot,rmsv,ang WRITE(134,*)"==================================================" WRITE(*,4000)mtot,rmsdpdt,etot,ztot,stot,rmsv,ang 4000 FORMAT(10x,'masse',5x,'rmsdpdt',5x,'energie',5x,'enstrophie' & ,5x,'entropie',5x,'rmsv',5x,'mt.ang',/,'GLOB ' & ,e10.3,e13.6,5e10.3/) close(134) !$OMP END MASTER END IF END SUBROUTINE check_conserve !--------------------------------------------------------------------- SUBROUTINE check_mass_conserve(f_ps,f_dps,mtot,rmsdpdt) USE mpi_mod USE mpipara USE transfert_omp_mod USE icosa IMPLICIT NONE TYPE(t_field),POINTER :: f_ps(:) TYPE(t_field),POINTER :: f_dps(:) REAL(rstd),POINTER :: ps(:),dps(:) REAL(rstd), INTENT(OUT) :: mtot, rmsdpdt INTEGER :: ind,i,j,ij REAL :: mloc, rmsloc REAL :: mloc_mpi, rmsloc_mpi mloc=0.0; rmsloc=0.0 DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) dps=f_dps(ind) ps=f_ps(ind) DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i IF (domain(ind)%own(i,j)) THEN mloc=mloc+ps(ij)*Ai(ij) rmsloc=rmsloc+dps(ij)*dps(ij) ENDIF ENDDO ENDDO ENDDO IF (using_mpi) THEN CALL reduce_sum_omp(mloc, mloc_mpi) CALL reduce_sum_omp(rmsloc, rmsloc_mpi) !$OMP MASTER CALL MPI_REDUCE(mloc_mpi, mtot, 1, MPI_REAL8, MPI_SUM, 0, comm_icosa, ierr) CALL MPI_REDUCE(rmsloc_mpi, rmsdpdt, 1, MPI_REAL8, MPI_SUM, 0, comm_icosa, ierr) !$OMP END MASTER ELSE ENDIF END SUBROUTINE check_mass_conserve !--------------------------------------------------------------------- SUBROUTINE check_en(f_ue,f_theta_rhodz,f_phis) USE icosa USE pression_mod USE vorticity_mod IMPLICIT NONE TYPE(t_field), POINTER :: f_ue(:) TYPE(t_field), POINTER :: f_theta_rhodz(:) TYPE(t_field), POINTER :: f_phis(:) REAL(rstd), POINTER :: ue(:,:) REAL(rstd), POINTER :: p(:,:) REAL(rstd), POINTER :: theta_rhodz(:,:) REAL(rstd), POINTER :: pk(:,:) REAL(rstd), POINTER :: phis(:) REAL(rstd), POINTER :: rhodz(:,:) INTEGER :: ind DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) ue=f_ue(ind) p=f_p(ind) rhodz=f_rhodz(ind) theta_rhodz=f_theta_rhodz(ind) pk=f_pk(ind) phis=f_phis(ind) CALL compute_en(ind,ue,p,rhodz,theta_rhodz,pk,phis) END DO END SUBROUTINE check_en SUBROUTINE compute_en(ind,u,p,rhodz,theta_rhodz,pk,phis) USE icosa USE disvert_mod USE wind_mod IMPLICIT NONE INTEGER,INTENT(IN)::ind REAL(rstd),INTENT(IN) :: u(3*iim*jjm,llm) REAL(rstd),INTENT(IN) :: p(iim*jjm,llm+1) REAL(rstd),INTENT(IN) :: rhodz(iim*jjm,llm) REAL(rstd),INTENT(IN) :: theta_rhodz(iim*jjm,llm) REAL(rstd),INTENT(IN) :: pk(iim*jjm,llm) REAL(rstd),INTENT(IN) :: phis(iim*jjm) REAL(rstd):: theta(iim*jjm,llm) REAL(rstd)::KE(iim*jjm,llm) REAL(rstd) :: ucenter(iim*jjm,3,llm) REAL(rstd) :: ulon(iim*jjm,llm) REAL(rstd) :: ulat(iim*jjm,llm) REAL(rstd)::masse(iim*jjm,llm) REAL(rstd)::rad,radomeg,lat,lon REAL(rstd) ::xyz(3) INTEGER :: i,j,ij,l,ij2 DO l = 1, llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i masse(ij,l) = rhodz(ij,l)*Ai(ij) theta(ij,l) = theta_rhodz(ij,l)/rhodz(ij,l) IF (domain(ind)%own(i,j)) THEN stot = stot + Ai(ij)*theta_rhodz(ij,l) END IF END DO END DO END DO DO l = 1, llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i IF (domain(ind)%own(i,j)) THEN KE(ij,l)= 1/(4*Ai(ij))*(le(ij+u_right)*de(ij+u_right)*u(ij+u_right,l)**2 + & le(ij+u_rup)*de(ij+u_rup)*u(ij+u_rup,l)**2 + & le(ij+u_lup)*de(ij+u_lup)*u(ij+u_lup,l)**2 + & le(ij+u_left)*de(ij+u_left)*u(ij+u_left,l)**2 + & le(ij+u_ldown)*de(ij+u_ldown)*u(ij+u_ldown,l)**2 + & le(ij+u_rdown)*de(ij+u_rdown)*u(ij+u_rdown,l)**2 ) rmsv = rmsv + 2*masse(ij,l)*KE(ij,l) END IF END DO END DO END DO DO l = 1, llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i IF (domain(ind)%own(i,j)) THEN etot = etot + masse(ij,l)*(phis(ij)+theta(ij,l)*pk(ij,l)+KE(ij,l)) END IF END DO END DO END DO !------------------------------ ANGULAR VELOCITY CALL compute_wind_centered(u,ucenter) CALL compute_wind_centered_lonlat_compound(ucenter,ulon,ulat) rad = radius radomeg = rad*omega DO l = 1, llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i IF (domain(ind)%own(i,j)) THEN lat=lat_i(ij) ang=ang+rad*cos(lat)*masse(ij,l)*(ulon(ij,l)+ radomeg*cos(lat)) END IF END DO END DO END DO END SUBROUTINE compute_en !--------------------------------------------------------------------- SUBROUTINE check_PV USE icosa IMPLICIT NONE REAL(rstd), POINTER :: vort(:,:) REAL(rstd), POINTER :: rhodz(:,:) INTEGER :: ind DO ind=1,ndomain IF (.NOT. assigned_domain(ind)) CYCLE CALL swap_dimensions(ind) CALL swap_geometry(ind) vort=f_vort(ind) rhodz=f_rhodz(ind) CALL compute_PV(vort,rhodz) ENDDO END SUBROUTINE check_PV SUBROUTINE compute_PV(vort,rhodz) USE icosa USE disvert_mod IMPLICIT NONE REAL(rstd),INTENT(IN) :: vort(2*iim*jjm,llm) REAL(rstd),INTENT(IN) :: rhodz(iim*jjm,llm) REAL(rstd)::qv1,qv2 REAL(rstd)::hv1,hv2 INTEGER :: i,j,ij,l,ij2 hv1 = 0.0 ; hv2 = 0.0 DO l = 1,llm DO j=jj_begin+1,jj_end-1 DO i=ii_begin+1,ii_end-1 ij=(j-1)*iim+i hv1 = Riv2(ij,vup) * rhodz(ij,l) & + Riv2(ij+t_rup,vldown) * rhodz(ij+t_rup,l) & + Riv2(ij+t_lup,vrdown) * rhodz(ij+t_lup,l) qv1 = ( vort(ij+z_up,l)+fv(ij+z_up) )/hv1 hv2 = Riv2(ij,vdown) * rhodz(ij,l) & + Riv2(ij+t_ldown,vrup) * rhodz(ij+t_ldown,l) & + Riv2(ij+t_rdown,vlup) * rhodz(ij+t_rdown,l) qv2 =( vort(ij+z_down,l)+fv(ij+z_down) )/hv2 ztot = ztot + (qv1*qv1*hv1 + qv2*qv2*hv2) ENDDO ENDDO ENDDO END SUBROUTINE compute_PV !--------------------------------------------------------------------- SUBROUTINE compute_rhodz(p,rhodz) USE icosa IMPLICIT NONE REAL(rstd),INTENT(IN) :: p(iim*jjm,llm+1) REAL(rstd),INTENT(OUT):: rhodz(iim*jjm,llm) REAL(rstd) ::mass(iim*jjm,llm) INTEGER :: i,j,ij,l DO l = 1, llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end ij=(j-1)*iim+i mass(ij,l) = ( p(ij,l) - p(ij,l+1) ) * Ai(ij)/g rhodz(ij,l) = mass(ij,l) / Ai(ij) ENDDO ENDDO ENDDO END SUBROUTINE compute_rhodz !==================================================================== END MODULE check_conserve_mod