foodie_integrator_factory Subroutine

public subroutine foodie_integrator_factory(scheme, integrator, stages, tolerance, nu, alpha, iterations, autoupdate, U)

Return a concrete instance of integrator_object given a scheme selection.

This is the FOODIE integrators factory.

Arguments

Type IntentOptional AttributesName
character(len=*), intent(in) :: scheme

Selected integrator given.

class(integrator_object), intent(out), allocatable:: integrator

The FOODIE integrator.

integer(kind=I_P), intent(in), optional :: stages

Stages of multi-stage methods.

real(kind=R_P), intent(in), optional :: tolerance

Tolerance on the local truncation error.

real(kind=R_P), intent(in), optional :: nu

Williams-Robert-Asselin filter coefficient.

real(kind=R_P), intent(in), optional :: alpha

Robert-Asselin filter coefficient.

integer(kind=I_P), intent(in), optional :: iterations

Implicit iterations.

logical, intent(in), optional :: autoupdate

Enable cyclic autoupdate for multistep.

class(integrand_object), intent(in), optional :: U

Integrand molding prototype.

Called By

proc~~foodie_integrator_factory~~CalledByGraph proc~foodie_integrator_factory foodie_integrator_factory proc~integrate~15 integrate proc~integrate~15->proc~foodie_integrator_factory proc~integrate~14 integrate proc~integrate~14->proc~foodie_integrator_factory proc~integrate~16 integrate proc~integrate~16->proc~foodie_integrator_factory proc~test~2 test proc~test~2->proc~integrate~15 proc~test test proc~test->proc~integrate~14 proc~execute~3 execute proc~execute~3->proc~integrate~16
Help


Source Code

  subroutine foodie_integrator_factory(scheme, integrator, stages, tolerance, nu, alpha, iterations, autoupdate, U)
  !< Return a concrete instance of [[integrator_object]] given a scheme selection.
  !<
  !< This is the FOODIE integrators factory.
  !<
  !< @note If an error occurs the error status of [[integrator_object]] is updated.
  character(*),                          intent(in)  :: scheme                      !< Selected integrator given.
  class(integrator_object), allocatable, intent(out) :: integrator                  !< The FOODIE integrator.
  integer(I_P),            optional,     intent(in)  :: stages                      !< Stages of multi-stage methods.
  real(R_P),               optional,     intent(in)  :: tolerance                   !< Tolerance on the local truncation error.
  real(R_P),               optional,     intent(in)  :: nu                          !< Williams-Robert-Asselin filter coefficient.
  real(R_P),               optional,     intent(in)  :: alpha                       !< Robert-Asselin filter coefficient.
  integer(I_P),            optional,     intent(in)  :: iterations                  !< Implicit iterations.
  logical,                 optional,     intent(in)  :: autoupdate                  !< Enable cyclic autoupdate for multistep.
  class(integrand_object), optional,     intent(in)  :: U                           !< Integrand molding prototype.
  type(integrator_adams_bashforth)                   :: int_adams_bashforth         !< Integrator Adams Bashforth.
  type(integrator_adams_bashforth_moulton)           :: int_adams_bashforth_moulton !< Integrator Adams Bashforth Moulton.
  type(integrator_adams_moulton)                     :: int_adams_moulton           !< Integrator Adams Moulton.
  type(integrator_back_df)                           :: int_back_df                 !< Integrator back differentiation formula.
  type(integrator_euler_explicit)                    :: int_euler_explicit          !< Integrator euler explicit.
  type(integrator_leapfrog)                          :: int_leapfrog                !< Integrator leapfrog.
  type(integrator_lmm_ssp)                           :: int_lmm_ssp                 !< Integrator LMM SSP.
  type(integrator_lmm_ssp_vss)                       :: int_lmm_ssp_vss             !< Integrator LMM SSP VSS.
  type(integrator_ms_runge_kutta_ssp)                :: int_ms_runge_kutta_ssp      !< Integrator multistep Runge Kutta ssp.
  type(integrator_runge_kutta_emd)                   :: int_runge_kutta_emd         !< Integrator Runge Kutta_embdedded.
  type(integrator_runge_kutta_ls)                    :: int_runge_kutta_ls          !< Integrator Runge Kutta low storage.
  type(integrator_runge_kutta_lssp)                  :: int_runge_kutta_lssp        !< Integrator linear Runge Kutta SSP.
  type(integrator_runge_kutta_ssp)                   :: int_runge_kutta_ssp         !< Integrator Runge Kutta SSP.

  if     (index(trim(adjustl(scheme)), trim(int_adams_bashforth_moulton%class_name())) > 0) then
    allocate(integrator_adams_bashforth_moulton :: integrator)
    select type(integrator)
    type is(integrator_adams_bashforth_moulton)
      call integrator%initialize(scheme=scheme, iterations=iterations, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_adams_bashforth%class_name())) > 0) then
    allocate(integrator_adams_bashforth :: integrator)
    select type(integrator)
    type is(integrator_adams_bashforth)
      call integrator%initialize(scheme=scheme, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_adams_moulton%class_name())) > 0) then
    allocate(integrator_adams_moulton :: integrator)
    select type(integrator)
    type is(integrator_adams_moulton)
      call integrator%initialize(scheme=scheme, iterations=iterations, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_back_df%class_name())) > 0) then
    allocate(integrator_back_df :: integrator)
    select type(integrator)
    type is(integrator_back_df)
      call integrator%initialize(scheme=scheme, iterations=iterations, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_euler_explicit%class_name())) > 0) then
    allocate(integrator_euler_explicit :: integrator)
    select type(integrator)
    type is(integrator_euler_explicit)
      call integrator%initialize(scheme=scheme, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_leapfrog%class_name())) > 0) then
    allocate(integrator_leapfrog :: integrator)
    select type(integrator)
    type is(integrator_leapfrog)
      call integrator%initialize(scheme=scheme, nu=nu, alpha=alpha, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_lmm_ssp_vss%class_name())) > 0) then
    allocate(integrator_lmm_ssp_vss :: integrator)
    select type(integrator)
    type is(integrator_lmm_ssp_vss)
      call integrator%initialize(scheme=scheme, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_lmm_ssp%class_name())) > 0) then
    allocate(integrator_lmm_ssp :: integrator)
    select type(integrator)
    type is(integrator_lmm_ssp)
      call integrator%initialize(scheme=scheme, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_ms_runge_kutta_ssp%class_name())) > 0) then
    allocate(integrator_ms_runge_kutta_ssp :: integrator)
    select type(integrator)
    type is(integrator_ms_runge_kutta_ssp)
      call integrator%initialize(scheme=scheme, iterations=iterations, autoupdate=autoupdate, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_runge_kutta_emd%class_name())) > 0) then
    allocate(integrator_runge_kutta_emd :: integrator)
    select type(integrator)
    type is(integrator_runge_kutta_emd)
      call integrator%initialize(scheme=scheme, tolerance=tolerance, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_runge_kutta_lssp%class_name())) > 0) then
    allocate(integrator_runge_kutta_lssp :: integrator)
    select type(integrator)
    type is(integrator_runge_kutta_lssp)
      call integrator%initialize(scheme=scheme, stages=stages, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_runge_kutta_ls%class_name())) > 0) then
    allocate(integrator_runge_kutta_ls :: integrator)
    select type(integrator)
    type is(integrator_runge_kutta_ls)
      call integrator%initialize(scheme=scheme, U=U)
    endselect
  elseif (index(trim(adjustl(scheme)), trim(int_runge_kutta_ssp%class_name())) > 0) then
    allocate(integrator_runge_kutta_ssp :: integrator)
    select type(integrator)
    type is(integrator_runge_kutta_ssp)
      call integrator%initialize(scheme=scheme, U=U)
    endselect
  else
    write(stderr, '(A)')'error: "'//trim(adjustl(scheme))//'" scheme is unknown!'
    stop
  endif
  endsubroutine foodie_integrator_factory


a a a a a add_burgers add_euler add_euler add_euler add_lorenz allocate_integrand_members allocate_integrand_members allocate_integrand_members allocate_integrand_members allocate_integrand_members amplitude_phase assign_abstract assign_integrand assign_integrand assign_integrand assign_multistage assign_multistage_multistep assign_multistep assign_real assign_real assign_real average_solution burgers_assign_burgers burgers_assign_real burgers_local_error burgers_multiply_burgers burgers_multiply_real check_error check_scheme_has_fast_mode check_scheme_has_fast_mode check_scheme_has_fast_mode class_name class_name class_name class_name class_name class_name class_name class_name class_name class_name class_name class_name class_name compute_dt compute_dt compute_dt compute_dt compute_dt compute_dt compute_dt compute_dx compute_inter_states compute_inter_states compute_inter_states compute_inter_states compute_inter_states conservative2primitive conservative2primitive conservative2primitive conservative2primitive conservative2primitive d2Burgers_dx2 dBurgers_dt dBurgers_dx description description description description destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy destroy_abstract destroy_multistage destroy_multistage_multistep destroy_multistep destroy_rk destroy_rk dEuler_dt dEuler_dt dEuler_dt dEuler_dt dEuler_dt dLorenz_dt dt_ratio dU_dt dU_dt dU_dt E E E E E error error error euler_assign_euler euler_assign_euler euler_assign_euler euler_assign_euler euler_assign_euler euler_assign_real euler_assign_real euler_assign_real euler_local_error euler_local_error euler_local_error euler_multiply_euler euler_multiply_euler euler_multiply_euler euler_multiply_real euler_multiply_real euler_multiply_real exact_solution exact_solution exact_solution execute execute execute export_tecplot export_tecplot export_tecplot foodie_integrator_class_names foodie_integrator_factory foodie_integrator_schemes H H H H H has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode has_fast_mode impose_boundary_conditions impose_boundary_conditions impose_boundary_conditions impose_boundary_conditions impose_boundary_conditions impose_boundary_conditions init init init init init init init init init init init init init init_rk init_rk initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize_order_s initialize_order_s_1 integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integr_assign_integr integrand_add_integrand integrand_add_integrand integrand_add_integrand integrand_add_integrand_fast integrand_add_integrand_fast integrand_add_integrand_fast integrand_add_real integrand_add_real integrand_add_real integrand_dimension integrand_dimension integrand_dimension integrand_multiply_integrand integrand_multiply_integrand integrand_multiply_integrand integrand_multiply_integrand_fast integrand_multiply_integrand_fast integrand_multiply_integrand_fast integrand_multiply_real integrand_multiply_real integrand_multiply_real integrand_multiply_real_scalar integrand_multiply_real_scalar integrand_multiply_real_scalar integrand_multiply_real_scalar_fast integrand_multiply_real_scalar_fast integrand_multiply_real_scalar_fast integrand_sub_integrand integrand_sub_integrand integrand_sub_integrand integrand_sub_real integrand_sub_real integrand_sub_real integrand_subtract_integrand_fast integrand_subtract_integrand_fast integrand_subtract_integrand_fast integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_fast integrate_order_2 integrate_order_2_fast integrate_order_3 integrate_order_3_fast integrate_order_s integrate_order_s_1 integrate_order_s_1_fast integrate_order_s_fast integrate_rk integrate_rk is_admissible is_available is_class_available is_multistage is_multistage is_multistage is_multistep is_multistep is_multistep is_scheme_available is_supported is_supported is_supported is_supported is_supported is_supported is_supported is_supported is_supported is_supported is_supported is_supported is_supported local_error local_error local_error lorenz_assign_lorenz lorenz_assign_real lorenz_local_error lorenz_multiply_lorenz lorenz_multiply_real new_Dt observed_order observed_order observed_order omega output output output output output output output output output output p p p p p parse_cli parse_cli parse_cli previous_step previous_step previous_step primitive2conservative primitive2conservative primitive2conservative primitive2conservative primitive2conservative r r r r r real_add_integrand real_add_integrand real_add_integrand real_multiply_burgers real_multiply_euler real_multiply_euler real_multiply_euler real_multiply_integrand real_multiply_integrand real_multiply_integrand real_multiply_lorenz real_scalar_multiply_integrand real_scalar_multiply_integrand real_scalar_multiply_integrand real_sub_integrand real_sub_integrand real_sub_integrand reconstruct_interfaces reconstruct_interfaces_states reconstruct_interfaces_states reconstruct_interfaces_states reconstruct_interfaces_states reconstruct_interfaces_states registers_number riemann_solver riemann_solver riemann_solver riemann_solver riemann_solver save_results save_results save_results save_results save_results save_results save_results save_results save_results save_time_serie save_time_serie save_time_serie save_time_serie save_time_serie scheme_number set_cli set_cli set_cli set_sin_wave_initial_state set_square_wave_initial_state stages_number stages_number stages_number steps_number steps_number steps_number sub_burgers sub_euler sub_euler sub_euler sub_lorenz supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes supported_schemes synchronize synchronize synchronize synchronize t_fast t_fast t_fast test test test_ab test_ab test_ab test_euler test_euler test_euler test_leapfrog test_leapfrog test_leapfrog test_ls_rk test_ls_rk test_ls_rk test_tvd_rk test_tvd_rk test_tvd_rk tokenize trigger_error update_previous update_previous update_previous_steps update_previous_steps update_previous_steps