MODULE tl_diag_mod ! !git $Id$ !svn $Id: tl_diag.F 1180 2023-07-13 02:42:10Z arango $ !================================================== Hernan G. Arango === ! Copyright (c) 2002-2023 The ROMS/TOMS Group ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This routine computes various tangent linear diagnostic fields. ! ! ! !======================================================================= ! implicit none ! PRIVATE PUBLIC :: tl_diag ! CONTAINS ! !*********************************************************************** SUBROUTINE tl_diag (ng, tile) !*********************************************************************** ! USE mod_param USE mod_grid USE mod_ocean USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & "ROMS/Tangent/tl_diag.F" ! integer :: IminS, ImaxS, JminS, JmaxS integer :: LBi, UBi, LBj, UBj, LBij, UBij ! ! Set horizontal starting and ending indices for automatic private ! storage arrays. ! IminS=BOUNDS(ng)%Istr(tile)-3 ImaxS=BOUNDS(ng)%Iend(tile)+3 JminS=BOUNDS(ng)%Jstr(tile)-3 JmaxS=BOUNDS(ng)%Jend(tile)+3 ! ! Determine array lower and upper bounds in the I- and J-directions. ! LBi=BOUNDS(ng)%LBi(tile) UBi=BOUNDS(ng)%UBi(tile) LBj=BOUNDS(ng)%LBj(tile) UBj=BOUNDS(ng)%UBj(tile) ! ! Set array lower and upper bounds for MIN(I,J) directions and ! MAX(I,J) directions. ! LBij=BOUNDS(ng)%LBij UBij=BOUNDS(ng)%UBij ! CALL wclock_on (ng, iTLM, 7, 45, MyFile) CALL tl_diag_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nstp(ng), kstp(ng), & & GRID(ng) % h, & & GRID(ng) % omn, & & GRID(ng) % tl_Hz, & & GRID(ng) % tl_z_r, & & GRID(ng) % tl_z_w, & & OCEAN(ng) % tl_rho, & & OCEAN(ng) % tl_u, & & OCEAN(ng) % tl_v, & & OCEAN(ng) % tl_ubar, & & OCEAN(ng) % tl_vbar, & & OCEAN(ng) % tl_zeta) CALL wclock_off (ng, iTLM, 7, 65, MyFile) ! RETURN END SUBROUTINE tl_diag ! !*********************************************************************** SUBROUTINE tl_diag_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nstp, kstp, & & h, omn, & & tl_Hz, tl_z_r, tl_z_w, & & tl_rho, tl_u, tl_v, & & tl_ubar, tl_vbar, tl_zeta) !*********************************************************************** ! USE mod_param USE mod_parallel USE mod_iounits USE mod_scalars ! USE distribute_mod, ONLY : mp_reduce ! implicit none ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: nstp, kstp ! real(r8), intent(in) :: h(LBi:,LBj:) real(r8), intent(in) :: omn(LBi:,LBj:) real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:) real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:) real(r8), intent(in) :: tl_z_w(LBi:,LBj:,0:) real(r8), intent(in) :: tl_rho(LBi:,LBj:,:) real(r8), intent(in) :: tl_u(LBi:,LBj:,:,:) real(r8), intent(in) :: tl_v(LBi:,LBj:,:,:) real(r8), intent(in) :: tl_ubar(LBi:,LBj:,:) real(r8), intent(in) :: tl_vbar(LBi:,LBj:,:) real(r8), intent(in) :: tl_zeta(LBi:,LBj:,:) ! ! Local variable declarations. ! integer :: NSUB, i, j, k, trd integer :: idia, istep real(r8), dimension(3) :: rbuffer character (len=3), dimension(3) :: op_handle ! real(r8) :: cff, my_avgke, my_avgpe, my_volume real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ke2d real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: pe2d ! character (len=8 ) :: kechar, pechar character (len=22) :: DateTime ! !----------------------------------------------------------------------- ! Set lower and upper tile bounds and staggered variables bounds for ! this horizontal domain partition. Notice that if tile=-1, it will ! set the values for the global grid. !----------------------------------------------------------------------- ! integer :: Istr, IstrB, IstrP, IstrR, IstrT, IstrM, IstrU integer :: Iend, IendB, IendP, IendR, IendT integer :: Jstr, JstrB, JstrP, JstrR, JstrT, JstrM, JstrV integer :: Jend, JendB, JendP, JendR, JendT integer :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1 integer :: Iendp1, Iendp2, Iendp2i, Iendp3 integer :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1 integer :: Jendp1, Jendp2, Jendp2i, Jendp3 ! Istr =BOUNDS(ng) % Istr (tile) IstrB =BOUNDS(ng) % IstrB (tile) IstrM =BOUNDS(ng) % IstrM (tile) IstrP =BOUNDS(ng) % IstrP (tile) IstrR =BOUNDS(ng) % IstrR (tile) IstrT =BOUNDS(ng) % IstrT (tile) IstrU =BOUNDS(ng) % IstrU (tile) Iend =BOUNDS(ng) % Iend (tile) IendB =BOUNDS(ng) % IendB (tile) IendP =BOUNDS(ng) % IendP (tile) IendR =BOUNDS(ng) % IendR (tile) IendT =BOUNDS(ng) % IendT (tile) Jstr =BOUNDS(ng) % Jstr (tile) JstrB =BOUNDS(ng) % JstrB (tile) JstrM =BOUNDS(ng) % JstrM (tile) JstrP =BOUNDS(ng) % JstrP (tile) JstrR =BOUNDS(ng) % JstrR (tile) JstrT =BOUNDS(ng) % JstrT (tile) JstrV =BOUNDS(ng) % JstrV (tile) Jend =BOUNDS(ng) % Jend (tile) JendB =BOUNDS(ng) % JendB (tile) JendP =BOUNDS(ng) % JendP (tile) JendR =BOUNDS(ng) % JendR (tile) JendT =BOUNDS(ng) % JendT (tile) ! Istrm3 =BOUNDS(ng) % Istrm3 (tile) ! Istr-3 Istrm2 =BOUNDS(ng) % Istrm2 (tile) ! Istr-2 Istrm1 =BOUNDS(ng) % Istrm1 (tile) ! Istr-1 IstrUm2=BOUNDS(ng) % IstrUm2(tile) ! IstrU-2 IstrUm1=BOUNDS(ng) % IstrUm1(tile) ! IstrU-1 Iendp1 =BOUNDS(ng) % Iendp1 (tile) ! Iend+1 Iendp2 =BOUNDS(ng) % Iendp2 (tile) ! Iend+2 Iendp2i=BOUNDS(ng) % Iendp2i(tile) ! Iend+2 interior Iendp3 =BOUNDS(ng) % Iendp3 (tile) ! Iend+3 Jstrm3 =BOUNDS(ng) % Jstrm3 (tile) ! Jstr-3 Jstrm2 =BOUNDS(ng) % Jstrm2 (tile) ! Jstr-2 Jstrm1 =BOUNDS(ng) % Jstrm1 (tile) ! Jstr-1 JstrVm2=BOUNDS(ng) % JstrVm2(tile) ! JstrV-2 JstrVm1=BOUNDS(ng) % JstrVm1(tile) ! JstrV-1 Jendp1 =BOUNDS(ng) % Jendp1 (tile) ! Jend+1 Jendp2 =BOUNDS(ng) % Jendp2 (tile) ! Jend+2 Jendp2i=BOUNDS(ng) % Jendp2i(tile) ! Jend+2 interior Jendp3 =BOUNDS(ng) % Jendp3 (tile) ! Jend+3 ! !----------------------------------------------------------------------- ! Compute and report out volume averaged kinetic, potential total ! energy, and volume of tangent linear fields. Notice that the ! proper adjoint of these quantities are not coded. !----------------------------------------------------------------------- ! ! Set time timestep counter and time level indices to process. Restart ! counter after 10 billion steps. ! istep=INT(MOD(REAL(iic(ng)-1,r8),1.0E+10_r8)) idia=nstp DateTime=time_code(ng) ! ! Compute kinetic and potential energy. ! IF (MOD(istep,ninfo(ng)).eq.0) THEN DO j=Jstr,Jend DO i=Istr,Iend ke2d(i,j)=0.0_r8 pe2d(i,j)=0.5_r8*g*tl_z_w(i,j,N(ng))*tl_z_w(i,j,N(ng)) END DO cff=g/rho0 DO k=N(ng),1,-1 DO i=Istr,Iend ke2d(i,j)=ke2d(i,j)+ & & tl_Hz(i,j,k)* & & 0.25_r8*(tl_u(i ,j,k,idia)*tl_u(i ,j,k,idia)+ & & tl_u(i+1,j,k,idia)*tl_u(i+1,j,k,idia)+ & & tl_v(i,j ,k,idia)*tl_v(i,j ,k,idia)+ & & tl_v(i,j+1,k,idia)*tl_v(i,j+1,k,idia)) pe2d(i,j)=pe2d(i,j)+ & & cff*tl_Hz(i,j,k)*(tl_rho(i,j,k)+1000.0_r8)* & & (tl_z_r(i,j,k)-tl_z_w(i,j,0)) END DO END DO END DO ! ! Integrate horizontally within one tile. In order to reduce the ! round-off errors, the summation is performed in two stages. First, ! the index j is collapsed and then the accumulation is carried out ! along index i. In this order, the partial sums consist on much ! fewer number of terms than in a straightforward two-dimensional ! summation. Thus, adding numbers which are orders of magnitude ! apart is avoided. ! DO i=Istr,Iend pe2d(i,Jend+1)=0.0_r8 pe2d(i,Jstr-1)=0.0_r8 ke2d(i,Jstr-1)=0.0_r8 END DO DO j=Jstr,Jend DO i=Istr,Iend pe2d(i,Jend+1)=pe2d(i,Jend+1)+omn(i,j)*h(i,j) pe2d(i,Jstr-1)=pe2d(i,Jstr-1)+omn(i,j)*pe2d(i,j) ke2d(i,Jstr-1)=ke2d(i,Jstr-1)+omn(i,j)*ke2d(i,j) END DO END DO my_volume=0.0_r8 my_avgpe=0.0_r8 my_avgke=0.0_r8 DO i=Istr,Iend my_volume=my_volume+pe2d(i,Jend+1) my_avgpe =my_avgpe +pe2d(i,Jstr-1) my_avgke =my_avgke +ke2d(i,Jstr-1) END DO ! ! Perform global summation: whoever gets first to the critical region ! resets global sums before global summation starts; after the global ! summation is completed, thread, which is the last one to enter the ! critical region, finalizes the computation of diagnostics and prints ! them out. ! NSUB=1 ! distributed-memory !$OMP CRITICAL (TL_DIAGNOSTICS) IF (tile_count.eq.0) THEN volume=0.0_r8 avgke=0.0_r8 avgpe=0.0_r8 END IF volume=volume+my_volume avgke=avgke+my_avgke avgpe=avgpe+my_avgpe tile_count=tile_count+1 IF (tile_count.eq.NSUB) THEN tile_count=0 rbuffer(1)=volume rbuffer(2)=avgke rbuffer(3)=avgpe op_handle(1)='SUM' op_handle(2)='SUM' op_handle(3)='SUM' CALL mp_reduce (ng, iTLM, 3, rbuffer, op_handle) volume=rbuffer(1) avgke=rbuffer(2) avgpe=rbuffer(3) trd=MyMaster avgke=avgke/volume avgpe=avgpe/volume avgkp=avgke+avgpe ! ! Report global run diagnotics values for the tangent linear kernel. ! IF (first_time(ng).eq.0) THEN first_time(ng)=1 IF (Master.and.(ng.eq.1)) THEN WRITE (stdout,10) 'TIME-STEP', 'YYYY-MM-DD hh:mm:ss.ss', & & 'KINETIC_ENRG', 'POTEN_ENRG', & & 'TOTAL_ENRG', 'NET_VOLUME' 10 FORMAT (/,1x,a,1x,a,2x,a,3x,a,4x,a,4x,a) END IF END IF ! IF (Master) THEN WRITE (stdout,20) istep, DateTime, & & avgke, avgpe, avgkp, volume 20 FORMAT (i10,1x,a,4(1pe14.6)) CALL my_flush (stdout) END IF ! ! If blowing-up, set exit_flag to stop computations. ! WRITE (kechar,'(1pe8.1)') avgke WRITE (pechar,'(1pe8.1)') avgpe DO i=1,8 IF ((kechar(i:i).eq.'N').or.(pechar(i:i).eq.'N').or. & & (kechar(i:i).eq.'n').or.(pechar(i:i).eq.'n').or. & & (kechar(i:i).eq.'*').or.(pechar(i:i).eq.'*')) THEN exit_flag=1 END IF END DO END IF !$OMP END CRITICAL (TL_DIAGNOSTICS) END IF ! RETURN END SUBROUTINE tl_diag_tile END MODULE tl_diag_mod