!> A module which contains a series of high level tests used in the !! linear and nonlinear test cases, as well as a driver function !! which runs GS2. module functional_tests use unit_tests, only: announce_check, agrees_with, process_check use unit_tests, only: process_test, announce_test, print_with_stars use runtime_tests, only: verbosity implicit none private !> Check that the relative error of the growth rates as a !! function of ky with respect to the given rslt is less than !! err. rslt should be an array of length naky which contains !! the growth rates as a function of ky. !! Returns .true. for pass and .false. for fail public :: check_growth_rate !> As check_growth_rate but for the frequency instead of growth rate public :: check_frequency !> Run gs2 and then call the test_function to check the results !! corresponding to the input file provided. test_name should !! be a character(*) containing the title of the tests, and !! test_function should be a logical function which returns !! .true. for pass and false for .fail. public :: test_gs2 public :: check_growth_rates_equal_in_list public :: test_function_interface interface !> Interface for functions passed to [[test_gs2]] logical function test_function_interface() implicit none end function test_function_interface end interface contains !> FIXME : Add documentation subroutine announce_functional_test(test_name) character(*), intent(in) :: test_name if (verbosity() .gt. 0) call print_with_stars('Starting functional test: ', test_name) end subroutine announce_functional_test !> FIXME : Add documentation subroutine close_functional_test(test_name) character(*), intent(in) :: test_name if (verbosity() .gt. 0) call print_with_stars('Completed functional test: ', test_name) end subroutine close_functional_test !> FIXME : Add documentation function check_growth_rate(rslt, err) use gs2_diagnostics, only: omegaavg_old => omegaavg use kt_grids, only: ntheta0, naky use mp, only: proc0 use run_parameters, only: use_old_diagnostics #ifdef NEW_DIAG use gs2_diagnostics_new, only: gnostics #endif real, dimension(naky), intent(in) :: rslt real, intent(in) :: err logical :: check_growth_rate logical :: check_result complex, dimension(ntheta0, naky) :: omegaavg integer :: ik check_growth_rate = .true. if (proc0) then call announce_check('growth rate') if (use_old_diagnostics) then omegaavg = omegaavg_old else #ifdef NEW_DIAG omegaavg = gnostics%current_results%omega_average #endif end if check_result = agrees_with(aimag(omegaavg(1,:)), rslt, err) call process_check(check_growth_rate, check_result, 'growth rate') if (.not. check_result) then if (naky == 1) then write(*, '("Expected : ",E14.6E3," but got ",E14.6E3)') rslt, aimag(omegaavg(1,:)) else do ik = 1, naky write(*, '("Expected : ",E14.6E3," but got ",E14.6E3," for ik=",I0)') rslt(ik), aimag(omegaavg(1,ik)),ik end do end if end if end if end function check_growth_rate !> FIXME : Add documentation logical function check_frequency(rslt, err) use gs2_diagnostics, only: omegaavg_old => omegaavg use kt_grids, only: ntheta0, naky use mp, only: proc0 use run_parameters, only: use_old_diagnostics #ifdef NEW_DIAG use gs2_diagnostics_new, only: gnostics #endif real, dimension(naky), intent(in) :: rslt real, intent(in) :: err logical :: check_result complex, dimension(ntheta0, naky) :: omegaavg integer :: ik check_frequency = .true. if (proc0) then call announce_check('frequency') if (use_old_diagnostics) then omegaavg = omegaavg_old else #ifdef NEW_DIAG omegaavg = gnostics%current_results%omega_average #endif end if check_result = agrees_with(real(omegaavg(1,:)), rslt, err) call process_check(check_frequency, check_result, 'frequency') if (.not. check_result) then if (naky == 1) then write(*, '("Expected : ",E14.6E3," but got ",E14.6E3)') rslt, real(omegaavg(1,:)) else do ik = 1, naky write(*, '("Expected : ",E14.6E3," but got ",E14.6E3," for ik=",I0)') rslt(ik), real(omegaavg(1,ik)),ik end do end if end if end if end function check_frequency !> FIXME : Add documentation function check_growth_rates_equal_in_list(err) use gs2_diagnostics, only: omegaavg_old => omegaavg use kt_grids, only: ntheta0, naky use mp, only: proc0, grp0, scope, subprocs, allprocs use mp, only: iproc, sum_allreduce use kt_grids, only: aky use run_parameters, only: wstar_units use job_manage, only: njobs use run_parameters, only: use_old_diagnostics #ifdef NEW_DIAG use gs2_diagnostics_new, only: gnostics #endif implicit none real, intent(in) :: err real, dimension(naky, njobs) :: omegas logical :: check_growth_rates_equal_in_list logical :: check_result integer :: i, ijob complex, dimension(ntheta0, naky) :: omegaavg check_growth_rates_equal_in_list = .true. omegas=0. omegaavg=0. if (proc0) then if (use_old_diagnostics) then omegaavg = omegaavg_old else #ifdef NEW_DIAG omegaavg = gnostics%current_results%omega_average #endif end if !check_result = agrees_with(aimag(omegaavg(1,:)), rslt, err) !call process_check(check_growth_rate, check_result, 'growth rate') call scope(allprocs) ! work out which job in the list we are do i = 1,njobs if (iproc==grp0(i-1)) ijob = i end do ! Get the omegas from this job omegas(:,ijob) = aimag(omegaavg(1,:)) if (wstar_units) then omegas(:, ijob) = omegas(:,ijob) * aky(:) / 2.0 end if else call scope(allprocs) end if call sum_allreduce(omegas) call scope(subprocs) do i = 2,njobs call announce_check('growth rate') check_result = agrees_with(omegas(:,1), omegas(:,i), err) call process_check( & check_growth_rates_equal_in_list, check_result, 'growth rate') end do end function check_growth_rates_equal_in_list !> Run GS2 and check that `test_function` returns true subroutine test_gs2(test_name, test_function) use gs2_main, only: run_gs2, gs2_program_state_type, finish_gs2 use mp, only: init_mp, proc0, finish_mp, mp_comm implicit none character(*), intent(in) :: test_name procedure(test_function_interface) :: test_function type(gs2_program_state_type) :: state call init_mp state%mp_comm = mp_comm if (proc0) call announce_functional_test(test_name) call run_gs2(finalize=.false., state=state) call announce_test('results') call process_test(test_function(), 'results') if (proc0) call close_functional_test(test_name) call finish_gs2(state, quiet=.true.) call finish_mp end subroutine test_gs2 end module functional_tests