model_template_m.f90 Source File


Contents

Source Code


Source Code

module model_template_m
    !! Implementation of a template model, from where development of a new model can be started.
    use perf_m, only : perf_print, perf_reset
    use precision_grillix_m, only : GP
    use screen_io_m, only : get_stdout
    use parallelisation_setup_m, only : is_rank_info_writer
    use comm_handler_m, only : comm_handler_t
    use params_template_m, only : read_params_template_paths, display_params_template_paths, read_all_params_template, &
        dtau, tau_fin, tau_snaps, tau_diag, tau_perf, path_snapsdir, diagfile_path, mms_on
    use mms_template_m, only : init_mms_template, display_mms_params_template, mms_diagnostics_template
    use state_template_m, only : alloc_state_template, write_template_state_netcdf, &
        tau, tstep
    use init_template_m, only : init_template
    use timestep_template_m, only : timestep_template
    use diagnostics_template_m, only : diagnostics_template
    implicit none

contains

    subroutine model_template(comm_handler)
        !! Model template subroutine
        type(comm_handler_t), intent(in) :: comm_handler
        !! Communication handler

        real(GP) :: next_snaps_time, next_diag_time, next_perf_time

        if (is_rank_info_writer) then
            write(get_stdout(),*) new_line('a') // repeat('-',80)
            write(get_stdout(),*)'START RUNNING TEMPLATE MODEL'
        endif

        call read_params_template_paths()
        call display_params_template_paths()
        call read_all_params_template()

        call alloc_state_template(comm_handler)

        if (mms_on) then
            call init_mms_template()
            call display_mms_params_template()
        endif

        call init_template(comm_handler)
        
        ! For new simulations, started from scratch, write out initial state
        if (tstep == 0) then
            call write_template_state_netcdf(comm_handler, dirname=trim(path_snapsdir))
        endif

        next_snaps_time = tau + tau_snaps
        next_diag_time = tau + tau_diag
        next_perf_time = tau + tau_perf

        timeloop: do
            if (tau > (tau_fin + 0.01_GP * dtau) ) then
                exit timeloop
            endif

            ! Write snapshot
            if (tau >= next_snaps_time - 0.01_GP * dtau) then
                if (is_rank_info_writer) then
                    write(get_stdout(),*)new_line('a') // repeat('-',80)
                    write(get_stdout(),503)tau, tstep
                    503 FORMAT('Writing snapshot at tau = ',ES20.12E3,X,'; tstep = ',I10)
                endif
                call write_template_state_netcdf(comm_handler, dirname=trim(path_snapsdir))
                next_snaps_time = next_snaps_time + tau_snaps
                if (is_rank_info_writer) then
                    write(get_stdout(),*) repeat('-',80) // new_line('a')
                endif

                if (mms_on) then
                    call mms_diagnostics_template(comm_handler, tau)
                endif
            endif

            ! Perform diagnostics
            if (tau >= next_diag_time - 0.01_GP * dtau) then
                if (is_rank_info_writer) then
                    write(get_stdout(),*)new_line('a') // repeat('-',80)
                    write(get_stdout(),504)tau, tstep
                    504 FORMAT('Performing diagnostics at tau = ',ES20.12E3,X,'; tstep = ',I10)
                endif
                call diagnostics_template(comm_handler, diagfile_path)
                next_diag_time = next_diag_time + tau_diag
                if (is_rank_info_writer) then
                    write(get_stdout(),*)repeat('-',80) // new_line('a')
                endif
            endif

            ! Print performance statistics and reset performance measurements
            if (tau >= next_perf_time - 0.01_GP * dtau) then
                call perf_print()
                call perf_reset()
                next_perf_time = next_perf_time + tau_perf
            endif

            call timestep_template(comm_handler, tau)

            ! Advance time
            tstep = tstep + 1
            tau = tau + dtau

        enddo timeloop

        if (is_rank_info_writer) then
            write(get_stdout(),*)'TEMPLATE MODEL RUN SUCCESFULLY'
            write(get_stdout(),*) repeat('-',80) // new_line('a') 
        endif

    end subroutine

end module