MODULE advect_mod USE icosa IMPLICIT NONE CONTAINS !========================================================================== SUBROUTINE init_advect(normal,tangent) USE domain_mod USE dimensions USE geometry USE metric USE vector IMPLICIT NONE REAL(rstd),INTENT(OUT) :: normal(3*iim*jjm,3) REAL(rstd),INTENT(OUT) :: tangent(3*iim*jjm,3) INTEGER :: i,j,n DO j=jj_begin-1,jj_end+1 DO i=ii_begin-1,ii_end+1 n=(j-1)*iim+i CALL cross_product2(xyz_v(n+z_rdown,:),xyz_v(n+z_rup,:),normal(n+u_right,:)) normal(n+u_right,:)=normal(n+u_right,:)/sqrt(sum(normal(n+u_right,:)**2)+1e-50)*ne(n,right) CALL cross_product2(xyz_v(n+z_up,:),xyz_v(n+z_lup,:),normal(n+u_lup,:)) normal(n+u_lup,:)=normal(n+u_lup,:)/sqrt(sum(normal(n+u_lup,:)**2)+1e-50)*ne(n,lup) CALL cross_product2(xyz_v(n+z_ldown,:),xyz_v(n+z_down,:),normal(n+u_ldown,:)) normal(n+u_ldown,:)=normal(n+u_ldown,:)/sqrt(sum(normal(n+u_ldown,:)**2)+1e-50)*ne(n,ldown) tangent(n+u_right,:)=xyz_v(n+z_rup,:)-xyz_v(n+z_rdown,:) tangent(n+u_right,:)=tangent(n+u_right,:)/sqrt(sum(tangent(n+u_right,:)**2)+1e-50) tangent(n+u_lup,:)=xyz_v(n+z_lup,:)-xyz_v(n+z_up,:) tangent(n+u_lup,:)=tangent(n+u_lup,:)/sqrt(sum(tangent(n+u_lup,:)**2)+1e-50) tangent(n+u_ldown,:)=xyz_v(n+z_down,:)-xyz_v(n+z_ldown,:) tangent(n+u_ldown,:)=tangent(n+u_ldown,:)/sqrt(sum(tangent(n+u_ldown,:)**2)+1e-50) END DO ENDDO END SUBROUTINE init_advect !======================================================================================= SUBROUTINE compute_gradq3d(qi,gradq3d) USE domain_mod USE dimensions USE geometry USE metric IMPLICIT NONE REAL(rstd),INTENT(IN) :: qi(iim*jjm,llm) REAL(rstd),INTENT(OUT) :: gradq3d(iim*jjm,llm,3) REAL(rstd) :: maxq,minq,minq_c,maxq_c REAL(rstd) :: alphamx,alphami,alpha ,maggrd,leng REAL(rstd) :: leng1,leng2 REAL(rstd) :: arr(2*iim*jjm) REAL(rstd) :: gradtri(2*iim*jjm,llm,3) REAL(rstd) :: ar INTEGER :: i,j,n,k,ind,l !========================================================================================== GRADIENT Do l = 1,llm DO j=jj_begin-1,jj_end+1 DO i=ii_begin-1,ii_end+1 n=(j-1)*iim+i CALL gradq(n,n+t_rup,n+t_lup,n+z_up,qi(:,l),gradtri(n+z_up,l,:),arr(n+z_up)) CALL gradq(n,n+t_ldown,n+t_rdown,n+z_down,qi(:,l),gradtri(n+z_down,l,:),arr(n+z_down)) END DO END DO END DO ! Do l =1,llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i gradq3d(n,:,:) = gradtri(n+z_up,:,:) + gradtri(n+z_down,:,:) + & gradtri(n+z_rup,:,:) + gradtri(n+z_ldown,:,:) + & gradtri(n+z_lup,:,:)+ gradtri(n+z_rdown,:,:) ar = arr(n+z_up)+arr(n+z_lup)+arr(n+z_ldown)+arr(n+z_down)+arr(n+z_rdown)+arr(n+z_rup) gradq3d(n,:,:) = gradq3d(n,:,:)/(ar+1.e-50) END DO END DO ! END DO !============================================================================================= LIMITING ! GO TO 120 DO l =1,llm DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i maggrd = dot_product(gradq3d(n,l,:),gradq3d(n,l,:)) maggrd = sqrt(maggrd) leng = max(sum((xyz_v(n+z_up,:) - xyz_i(n,:))**2),sum((xyz_v(n+z_down,:) - xyz_i(n,:))**2), & sum((xyz_v(n+z_rup,:) - xyz_i(n,:))**2),sum((xyz_v(n+z_rdown,:) - xyz_i(n,:))**2), & sum((xyz_v(n+z_lup,:) - xyz_i(n,:))**2),sum((xyz_v(n+z_ldown,:) - xyz_i(n,:))**2)) maxq_c = qi(n,l) + maggrd*sqrt(leng) minq_c = qi(n,l) - maggrd*sqrt(leng) maxq = max(qi(n,l),qi(n+t_right,l),qi(n+t_lup,l),qi(n+t_rup,l),qi(n+t_left,l), & qi(n+t_rdown,l),qi(n+t_ldown,l)) minq = min(qi(n,l),qi(n+t_right,l),qi(n+t_lup,l),qi(n+t_rup,l),qi(n+t_left,l), & qi(n+t_rdown,l),qi(n+t_ldown,l)) alphamx = (maxq - qi(n,l)) ; alphamx = alphamx/(maxq_c - qi(n,l) ) alphamx = max(alphamx,0.0) alphami = (minq - qi(n,l)); alphami = alphami/(minq_c - qi(n,l)) alphami = max(alphami,0.0) alpha = min(alphamx,alphami,1.0) gradq3d(n,l,:) = alpha*gradq3d(n,l,:) END DO END DO END DO END SUBROUTINE compute_gradq3d !=================================================================================================== SUBROUTINE compute_advect_horiz(normal,tangent,qi,gradq3d,him,ue,he,bigt) USE domain_mod USE dimensions USE geometry USE metric IMPLICIT NONE REAL(rstd),INTENT(IN) :: normal(3*iim*jjm,3) REAL(rstd),INTENT(IN) :: tangent(3*iim*jjm,3) REAL(rstd),INTENT(INOUT) :: qi(iim*jjm,llm) REAL(rstd),INTENT(IN) :: gradq3d(iim*jjm,llm,3) REAL(rstd),INTENT(INOUT) :: him(iim*jjm,llm) REAL(rstd),INTENT(IN) :: ue(iim*3*jjm,llm) REAL(rstd),INTENT(IN) :: he(3*iim*jjm,llm) ! mass flux REAL(rstd),INTENT(IN) :: bigt REAL(rstd) :: dqi(iim*jjm,llm),dhi(iim*jjm,llm) REAL(rstd) :: cc(3*iim*jjm,llm,3) REAL(rstd) :: v_e(3*iim*jjm,llm,3) REAL(rstd) :: up_e REAL(rstd) :: qe(3*iim*jjm,llm) REAL(rstd) :: ed(3) REAL(rstd) :: flxx(3*iim*jjm,llm) INTEGER :: i,j,n,k,l REAL(rstd):: him_old(llm) !go to 123 DO l = 1,llm DO j=jj_begin-1,jj_end+1 DO i=ii_begin-1,ii_end+1 n=(j-1)*iim+i up_e =1/de(n+u_right)*( & wee(n+u_right,1,1)*le(n+u_rup)*ue(n+u_rup,l)+ & wee(n+u_right,2,1)*le(n+u_lup)*ue(n+u_lup,l)+ & wee(n+u_right,3,1)*le(n+u_left)*ue(n+u_left,l)+ & wee(n+u_right,4,1)*le(n+u_ldown)*ue(n+u_ldown,l)+ & wee(n+u_right,5,1)*le(n+u_rdown)*ue(n+u_rdown,l)+ & wee(n+u_right,1,2)*le(n+t_right+u_ldown)*ue(n+t_right+u_ldown,l)+ & wee(n+u_right,2,2)*le(n+t_right+u_rdown)*ue(n+t_right+u_rdown,l)+ & wee(n+u_right,3,2)*le(n+t_right+u_right)*ue(n+t_right+u_right,l)+ & wee(n+u_right,4,2)*le(n+t_right+u_rup)*ue(n+t_right+u_rup,l)+ & wee(n+u_right,5,2)*le(n+t_right+u_lup)*ue(n+t_right+u_lup,l) & ) v_e(n+u_right,l,:)= ue(n+u_right,l)*normal(n+u_right,:) + up_e*tangent(n+u_right,:) up_e=1/de(n+u_lup)*( & wee(n+u_lup,1,1)*le(n+u_left)*ue(n+u_left,l)+ & wee(n+u_lup,2,1)*le(n+u_ldown)*ue(n+u_ldown,l)+ & wee(n+u_lup,3,1)*le(n+u_rdown)*ue(n+u_rdown,l)+ & wee(n+u_lup,4,1)*le(n+u_right)*ue(n+u_right,l)+ & wee(n+u_lup,5,1)*le(n+u_rup)*ue(n+u_rup,l)+ & wee(n+u_lup,1,2)*le(n+t_lup+u_right)*ue(n+t_lup+u_right,l)+ & wee(n+u_lup,2,2)*le(n+t_lup+u_rup)*ue(n+t_lup+u_rup,l)+ & wee(n+u_lup,3,2)*le(n+t_lup+u_lup)*ue(n+t_lup+u_lup,l)+ & wee(n+u_lup,4,2)*le(n+t_lup+u_left)*ue(n+t_lup+u_left,l)+ & wee(n+u_lup,5,2)*le(n+t_lup+u_ldown)*ue(n+t_lup+u_ldown,l) & ) v_e(n+u_lup,l,:)= ue(n+u_lup,l)*normal(n+u_lup,:) + up_e*tangent(n+u_lup,:) up_e=1/de(n+u_ldown)*( & wee(n+u_ldown,1,1)*le(n+u_rdown)*ue(n+u_rdown,l)+ & wee(n+u_ldown,2,1)*le(n+u_right)*ue(n+u_right,l)+ & wee(n+u_ldown,3,1)*le(n+u_rup)*ue(n+u_rup,l)+ & wee(n+u_ldown,4,1)*le(n+u_lup)*ue(n+u_lup,l)+ & wee(n+u_ldown,5,1)*le(n+u_left)*ue(n+u_left,l)+ & wee(n+u_ldown,1,2)*le(n+t_ldown+u_lup)*ue(n+t_ldown+u_lup,l)+ & wee(n+u_ldown,2,2)*le(n+t_ldown+u_left)*ue(n+t_ldown+u_left,l)+ & wee(n+u_ldown,3,2)*le(n+t_ldown+u_ldown)*ue(n+t_ldown+u_ldown,l)+ & wee(n+u_ldown,4,2)*le(n+t_ldown+u_rdown)*ue(n+t_ldown+u_rdown,l)+ & wee(n+u_ldown,5,2)*le(n+t_ldown+u_right)*ue(n+t_ldown+u_right,l) & ) v_e(n+u_ldown,l,:)= ue(n+u_ldown,l)*normal(n+u_ldown,:) + up_e*tangent(n+u_ldown,:) cc(n+u_right,l,:) = xyz_e(n+u_right,:) - v_e(n+u_right,l,:)*0.5*bigt !! redge is mid point of edge i cc(n+u_lup,l,:) = xyz_e(n+u_lup,:) - v_e(n+u_lup,l,:)*0.5*bigt cc(n+u_ldown,l,:) = xyz_e(n+u_ldown,:) - v_e(n+u_ldown,l,:)*0.5*bigt ENDDO ENDDO END DO !123 continue !========================================================================================================== DO l = 1,llm DO j=jj_begin-1,jj_end+1 DO i=ii_begin-1,ii_end+1 n=(j-1)*iim+i if (ne(n,right)*ue(n+u_right,l)>0) then ed = cc(n+u_right,l,:) - xyz_i(n,:) qe(n+u_right,l)=qi(n,l)+sum2(gradq3d(n,l,:),ed) else ed = cc(n+u_right,l,:) - xyz_i(n+t_right,:) qe(n+u_right,l)=qi(n+t_right,l)+sum2(gradq3d(n+t_right,l,:),ed) endif if (ne(n,lup)*ue(n+u_lup,l)>0) then ed = cc(n+u_lup,l,:) - xyz_i(n,:) qe(n+u_lup,l)=qi(n,l)+sum2(gradq3d(n,l,:),ed) else ed = cc(n+u_lup,l,:) - xyz_i(n+t_lup,:) qe(n+u_lup,l)=qi(n+t_lup,l)+sum2(gradq3d(n+t_lup,l,:),ed) endif if (ne(n,ldown)*ue(n+u_ldown,l)>0) then ed = cc(n+u_ldown,l,:) - xyz_i(n,:) qe(n+u_ldown,l)=qi(n,l)+ sum2(gradq3d(n,l,:),ed) else ed = cc(n+u_ldown,l,:) - xyz_i(n+t_ldown,:) qe(n+u_ldown,l)=qi(n+t_ldown,l)+sum2(gradq3d(n+t_ldown,l,:),ed) endif end do end do END DO DO j=jj_begin-1,jj_end+1 DO i=ii_begin-1,ii_end+1 n=(j-1)*iim+i flxx(n+u_right,:) = he(n+u_right,:)*qe(n+u_right,:)*ne(n,right) flxx(n+u_lup,:) = he(n+u_lup,:)*qe(n+u_lup,:)*ne(n,lup) flxx(n+u_ldown,:) = he(n+u_ldown,:)*qe(n+u_ldown,:)*ne(n,ldown) ENDDO ENDDO DO j=jj_begin,jj_end DO i=ii_begin,ii_end n=(j-1)*iim+i dhi(n,:)= -(1/Ai(n))*(he(n+u_right,:)*ne(n,right) + he(n+u_lup,:)*ne(n,lup) & + he(n+u_ldown,:)*ne(n,ldown) + he(n+u_rup,:)*ne(n,rup) & + he(n+u_left,:)*ne(n,left) + he(n+u_rdown,:)*ne(n,rdown) ) dqi(n,:)= -(1/Ai(n))*(flxx(n+u_right,:)+flxx(n+u_lup,:) & +flxx(n+u_ldown,:) - flxx(n+u_rup,:) & - flxx(n+u_left,:) - flxx(n+u_rdown,:) ) him_old(:) = him(n,:) him(n,:) = him(n,:) + dhi(n,:) qi(n,:) = (qi(n,:)*him_old(:) + dqi(n,:))/him(n,:) END DO END DO CONTAINS !==================================================================================== REAL FUNCTION sum2(a1,a2) IMPLICIT NONE REAL,INTENT(IN):: a1(3), a2(3) sum2 = a1(1)*a2(1)+a1(2)*a2(2)+a1(3)*a2(3) END FUNCTION sum2 END SUBROUTINE compute_advect_horiz !========================================================================== SUBROUTINE gradq(n0,n1,n2,n3,q,dq,det) USE geometry USE metric USE dimensions IMPLICIT NONE INTEGER, INTENT(IN) :: n0,n1,n2,n3 REAL,INTENT(IN) :: q(iim*jjm) REAL,INTENT(OUT) :: dq(3) REAL(rstd) ::det,detx,dety,detz INTEGER :: info INTEGER :: IPIV(3) REAL(rstd) :: A(3,3) REAL(rstd) :: B(3) A(1,1)=xyz_i(n1,1) -xyz_i(n0,1); A(1,2)=xyz_i(n1,2)- xyz_i(n0,2); A(1,3)=xyz_i(n1,3) - xyz_i(n0,3) A(2,1)=xyz_i(n2,1) - xyz_i(n0,1); A(2,2)=xyz_i(n2,2) - xyz_i(n0,2); A(2,3)=xyz_i(n2,3) - xyz_i(n0,3) A(3,1)=xyz_v(n3,1); A(3,2)= xyz_v(n3,2); A(3,3)= xyz_v(n3,3) dq(1) = q(n1)-q(n0) dq(2) = q(n2)-q(n0) dq(3) = 0.0 ! CALL DGESV(3,1,A,3,IPIV,dq(:),3,info) CALL determinant(A(:,1),A(:,2),A(:,3),det) CALL determinant(dq,A(:,2),A(:,3),detx) CALL determinant(A(:,1),dq,A(:,3),dety) CALL determinant(A(:,1),A(:,2),dq,detz) dq(1) = detx dq(2) = dety dq(3) = detz END SUBROUTINE gradq !========================================================================== SUBROUTINE determinant(a1,a2,a3,det) IMPLICIT NONE REAL, DIMENSION(3) :: a1, a2,a3 REAL :: det,x1,x2,x3 x1 = a1(1) * (a2(2) * a3(3) - a2(3) * a3(2)) x2 = a1(2) * (a2(1) * a3(3) - a2(3) * a3(1)) x3 = a1(3) * (a2(1) * a3(2) - a2(2) * a3(1)) det = x1 - x2 + x3 END SUBROUTINE determinant END MODULE advect_mod