Skip to content
Snippets Groups Projects
Commit b727a9ae authored by ESCOBAR MUNOZ Juan's avatar ESCOBAR MUNOZ Juan
Browse files

Juan 17/03/2023:ZSOLVER/multigrid.f90,finalize_mnh.f90, Managed output in...

Juan 17/03/2023:ZSOLVER/multigrid.f90,finalize_mnh.f90, Managed output in timing.txt , verbose=1/residual , verbose=10/timing by level, verbose=100/timing by iteration
parent d4c1b0fe
No related branches found
No related merge requests found
!MNH_LIC Copyright 2021-2021 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
!-----------------------------------------------------------------
! Author:
! P. Wautelet 06/07/2021
! Modifications:
!
!-----------------------------------------------------------------
MODULE MODE_FINALIZE_MNH
IMPLICIT NONE
CONTAINS
SUBROUTINE FINALIZE_MNH
USE MODD_CONF, only: LCHECK, NMODEL
USE MODD_IO, only: NIO_VERB, NVERB_DEBUG
USE MODD_LUNIT, only: TLUOUT0
USE MODD_LUNIT_n, only: LUNIT_MODEL
USE MODD_MNH_SURFEX_n, only: SURFEX_DEALLO_LIST
#ifdef CPLOASIS
USE MODD_SFX_OASIS, only: LOASIS
#endif
USE MODE_INIT_ll, only: END_PARA_ll
USE MODE_IO_FILE, only: IO_File_close
USE MODE_IO_MANAGE_STRUCT, only: IO_Filelist_print
#ifdef MNH_OPENACC
USE MODE_MNH_ZWORK, only: PRINT_FLATPOOL_STATS
#endif
USE multigrid, only: mg_finalise
USE MODE_MPPDB, only: MPPDB_BARRIER
USE MODE_MSG, only: MSG_STATS
#ifdef CPLOASIS
USE MODI_SFX_OASIS_END
#endif
IMPLICIT NONE
INTEGER :: IRESP
INTEGER :: JMODEL
!Print the list of all files and some statistics on them
IF ( NIO_VERB >= NVERB_DEBUG ) CALL IO_Filelist_print()
#ifdef MNH_OPENACC
CALL PRINT_FLATPOOL_STATS()
#endif
call mg_finalise()
!Print the number of printed messages via Print_msg
CALL MSG_STATS()
!Close all the opened 'output listing' files
IF ( TLUOUT0%LOPENED ) CALL IO_File_close(TLUOUT0)
DO JMODEL = 1, NMODEL
IF ( ASSOCIATED( LUNIT_MODEL(JMODEL)%TLUOUT ) ) THEN
IF ( LUNIT_MODEL(JMODEL)%TLUOUT%LOPENED ) CALL IO_File_close( LUNIT_MODEL(JMODEL)%TLUOUT)
END IF
END DO
!Finalize the parallel libraries
IF ( LCHECK ) THEN
CALL MPPDB_BARRIER()
ELSE
CALL END_PARA_ll( IRESP )
#ifdef CPLOASIS
IF ( LOASIS ) THEN
CALL SFX_OASIS_END()
END IF
#endif
END IF
!Free SURFEX structures if necessary
CALL SURFEX_DEALLO_LIST()
END SUBROUTINE FINALIZE_MNH
END MODULE MODE_FINALIZE_MNH
......@@ -132,6 +132,8 @@ private
type(time), allocatable, dimension(:,:) :: t_coarsesolve
type(time), allocatable, dimension(:,:) :: t_total
logical :: Gmg_initialized = .FALSE.
contains
!==================================================================
......@@ -172,6 +174,8 @@ contains
real , dimension(:,:,:) , pointer , contiguous :: zxu_mg_st,zxb_mg_st,zxr_mg_st
Gmg_initialized = .TRUE.
if (i_am_master_mpi) &
write(STDOUT,*) '*** Initialising multigrid ***'
! Check that cell counts are valid
......@@ -413,24 +417,28 @@ contains
character(len=80) :: s
integer :: ierr
If ( .not. Gmg_initialized ) return
if (i_am_master_mpi) &
write(STDOUT,*) '*** Finalising multigrid ***'
! Deallocate memory
level = mg_param%n_lev
m = pproc
reduced_m = .false.
!if (mg_param%verbose >= 0 )
call print_timerinfo("--- V-cycle timing results ---")
do while (level > 0)
write(s,'("level = ",I3,", m = ",I3)') level,m
call print_timerinfo(s)
call print_elapsed(t_smooth(level,m),.True.,1.0_rl)
call print_elapsed(t_restrict(level,m),.True.,1.0_rl)
call print_elapsed(t_prolongate(level,m),.True.,1.0_rl)
call print_elapsed(t_residual(level,m),.True.,1.0_rl)
call print_elapsed(t_addcorr(level,m),.True.,1.0_rl)
call print_elapsed(t_coarsesolve(level,m),.True.,1.0_rl)
call print_elapsed(t_total(level,m),.True.,1.0_rl)
if (mg_param%verbose >= 10 ) then
write(s,'("level = ",I3,", m = ",I3)') level,m
call print_timerinfo(s)
call print_elapsed(t_smooth(level,m),.True.,1.0_rl)
call print_elapsed(t_restrict(level,m),.True.,1.0_rl)
call print_elapsed(t_prolongate(level,m),.True.,1.0_rl)
call print_elapsed(t_residual(level,m),.True.,1.0_rl)
call print_elapsed(t_addcorr(level,m),.True.,1.0_rl)
call print_elapsed(t_coarsesolve(level,m),.True.,1.0_rl)
call print_elapsed(t_total(level,m),.True.,1.0_rl)
endif
if (LUseO) then
deallocate(xu_mg(level,m)%s)
deallocate(xb_mg(level,m)%s)
......@@ -468,7 +476,8 @@ contains
deallocate(t_residual)
deallocate(t_addcorr)
deallocate(t_coarsesolve)
if (i_am_master_mpi) write(STDOUT,'("")')
!if ( (i_am_master_mpi) ) write(STDOUT,'("")')
if ( (i_am_master_mpi) .and. (mg_param%verbose >= 10) ) write(STDOUT,'("")')
end subroutine mg_finalise
!==================================================================
......@@ -1763,6 +1772,7 @@ contains
end if
end if
end if
if (mg_param%verbose >= 100 ) then
call print_timerinfo("--- Main iteration timing results ---")
call print_elapsed(t_apply,.True.,1.0_rl)
call print_elapsed(t_res,.True.,1.0_rl)
......@@ -1770,6 +1780,7 @@ contains
call print_elapsed(t_l2norm,.True.,1.0_rl)
call print_elapsed(t_scalprod,.True.,1.0_rl)
call print_elapsed(t_mainloop,.True.,1.0_rl)
end if
!!$ if (i_am_master_mpi) write(STDOUT,'("")')
end subroutine mg_solve_mnh
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment