solver_factory.F90 Source File


This file depends on

sourcefile~~solver_factory.f90~~EfferentGraph sourcefile~solver_factory.f90 solver_factory.F90 sourcefile~input.f90 input.F90 sourcefile~solver_factory.f90->sourcefile~input.f90 sourcefile~matrix.f90 matrix.F90 sourcefile~solver_factory.f90->sourcefile~matrix.f90 sourcefile~input_interface.f90 input_interface.F90 sourcefile~input.f90->sourcefile~input_interface.f90 sourcefile~matrix_base.f90 matrix_base.F90 sourcefile~matrix.f90->sourcefile~matrix_base.f90 sourcefile~matrix_coo.f90 matrix_coo.F90 sourcefile~matrix.f90->sourcefile~matrix_coo.f90 sourcefile~matrix_crs.f90 matrix_crs.F90 sourcefile~matrix.f90->sourcefile~matrix_crs.f90 sourcefile~matrix_dense.f90 matrix_dense.F90 sourcefile~matrix.f90->sourcefile~matrix_dense.f90 sourcefile~core.f90 core.F90 sourcefile~input_interface.f90->sourcefile~core.f90 sourcefile~project_settings.f90 project_settings.F90 sourcefile~input_interface.f90->sourcefile~project_settings.f90 sourcefile~domain.f90 domain.F90 sourcefile~matrix_base.f90->sourcefile~domain.f90 sourcefile~matrix_coo.f90->sourcefile~matrix_base.f90 sourcefile~matrix_coo.f90->sourcefile~core.f90 sourcefile~matrix_coo.f90->sourcefile~domain.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_base.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_coo.f90 sourcefile~matrix_crs.f90->sourcefile~core.f90 sourcefile~matrix_crs.f90->sourcefile~domain.f90 sourcefile~matrix_dense.f90->sourcefile~matrix_base.f90 sourcefile~matrix_dense.f90->sourcefile~core.f90 sourcefile~matrix_dense.f90->sourcefile~domain.f90 sourcefile~allocate.f90 allocate.F90 sourcefile~core.f90->sourcefile~allocate.f90 sourcefile~check_range.f90 check_range.F90 sourcefile~core.f90->sourcefile~check_range.f90 sourcefile~deallocate.f90 deallocate.F90 sourcefile~core.f90->sourcefile~deallocate.f90 sourcefile~error.f90 error.F90 sourcefile~core.f90->sourcefile~error.f90 sourcefile~fortran_utils.f90 fortran_utils.F90 sourcefile~core.f90->sourcefile~fortran_utils.f90 sourcefile~string_utils.f90 string_utils.F90 sourcefile~core.f90->sourcefile~string_utils.f90 sourcefile~types.f90 types.F90 sourcefile~core.f90->sourcefile~types.f90 sourcefile~unique.f90 unique.F90 sourcefile~core.f90->sourcefile~unique.f90 sourcefile~vtk.f90 vtk.F90 sourcefile~core.f90->sourcefile~vtk.f90 sourcefile~vtk_constants.f90 vtk_constants.F90 sourcefile~core.f90->sourcefile~vtk_constants.f90 sourcefile~adjacency.f90 adjacency.F90 sourcefile~domain.f90->sourcefile~adjacency.f90 sourcefile~domain_manager.f90 domain_manager.F90 sourcefile~domain.f90->sourcefile~domain_manager.f90 sourcefile~element.f90 element.F90 sourcefile~domain.f90->sourcefile~element.f90 sourcefile~element_factory.f90 element_factory.F90 sourcefile~domain.f90->sourcefile~element_factory.f90 sourcefile~multicoloring.f90 multicoloring.F90 sourcefile~domain.f90->sourcefile~multicoloring.f90 sourcefile~reordering.f90 reordering.F90 sourcefile~domain.f90->sourcefile~reordering.f90 sourcefile~side.f90 side.F90 sourcefile~domain.f90->sourcefile~side.f90 sourcefile~side_factory.f90 side_factory.F90 sourcefile~domain.f90->sourcefile~side_factory.f90 sourcefile~project_settings.f90->sourcefile~core.f90 sourcefile~adjacency_element.f90 adjacency_element.F90 sourcefile~adjacency.f90->sourcefile~adjacency_element.f90 sourcefile~adjacency_node.f90 adjacency_node.F90 sourcefile~adjacency.f90->sourcefile~adjacency_node.f90 sourcefile~allocate.f90->sourcefile~error.f90 sourcefile~deallocate.f90->sourcefile~error.f90 sourcefile~domain_manager.f90->sourcefile~input.f90 sourcefile~domain_manager.f90->sourcefile~core.f90 sourcefile~domain_manager.f90->sourcefile~adjacency.f90 sourcefile~domain_manager.f90->sourcefile~element.f90 sourcefile~domain_manager.f90->sourcefile~element_factory.f90 sourcefile~domain_manager.f90->sourcefile~multicoloring.f90 sourcefile~domain_manager.f90->sourcefile~reordering.f90 sourcefile~domain_manager.f90->sourcefile~side.f90 sourcefile~domain_manager.f90->sourcefile~side_factory.f90 sourcefile~element.f90->sourcefile~input.f90 sourcefile~element.f90->sourcefile~core.f90 sourcefile~element_factory.f90->sourcefile~input.f90 sourcefile~element_factory.f90->sourcefile~core.f90 sourcefile~element_factory.f90->sourcefile~element.f90 sourcefile~memory_stats_wrapper.f90 memory_stats_wrapper.F90 sourcefile~fortran_utils.f90->sourcefile~memory_stats_wrapper.f90 sourcefile~signal_flag_wrapper.f90 signal_flag_wrapper.F90 sourcefile~fortran_utils.f90->sourcefile~signal_flag_wrapper.f90 sourcefile~system_info_wrapper.f90 system_info_wrapper.F90 sourcefile~fortran_utils.f90->sourcefile~system_info_wrapper.f90 sourcefile~multicoloring.f90->sourcefile~core.f90 sourcefile~multicoloring.f90->sourcefile~adjacency_element.f90 sourcefile~reordering.f90->sourcefile~core.f90 sourcefile~reordering.f90->sourcefile~element.f90 sourcefile~reordering.f90->sourcefile~adjacency_node.f90 sourcefile~side.f90->sourcefile~input.f90 sourcefile~side.f90->sourcefile~core.f90 sourcefile~side_factory.f90->sourcefile~input.f90 sourcefile~side_factory.f90->sourcefile~core.f90 sourcefile~side_factory.f90->sourcefile~side.f90 sourcefile~string_utils.f90->sourcefile~allocate.f90 sourcefile~array.f90 array.F90 sourcefile~types.f90->sourcefile~array.f90 sourcefile~gauss.f90 gauss.F90 sourcefile~types.f90->sourcefile~gauss.f90 sourcefile~pointer.f90 pointer.F90 sourcefile~types.f90->sourcefile~pointer.f90 sourcefile~variable.f90 variable.F90 sourcefile~types.f90->sourcefile~variable.f90 sourcefile~vector.f90 vector.F90 sourcefile~types.f90->sourcefile~vector.f90 sourcefile~unique.f90->sourcefile~allocate.f90 sourcefile~vtk.f90->sourcefile~allocate.f90 sourcefile~vtk.f90->sourcefile~deallocate.f90 sourcefile~vtk.f90->sourcefile~unique.f90 sourcefile~vtk.f90->sourcefile~vtk_constants.f90 sourcefile~vtk.f90->sourcefile~array.f90 sourcefile~vtk_wrapper.f90 vtk_wrapper.F90 sourcefile~vtk.f90->sourcefile~vtk_wrapper.f90 sourcefile~vtu_wrapper.f90 vtu_wrapper.F90 sourcefile~vtk.f90->sourcefile~vtu_wrapper.f90 sourcefile~adjacency_element.f90->sourcefile~core.f90 sourcefile~adjacency_element.f90->sourcefile~element.f90 sourcefile~adjacency_node.f90->sourcefile~core.f90 sourcefile~array.f90->sourcefile~allocate.f90 sourcefile~array.f90->sourcefile~deallocate.f90 sourcefile~c_utils.f90 c_utils.F90 sourcefile~memory_stats_wrapper.f90->sourcefile~c_utils.f90 sourcefile~signal_flag.f90 signal_flag.F90 sourcefile~signal_flag_wrapper.f90->sourcefile~signal_flag.f90 sourcefile~system_info_wrapper.f90->sourcefile~c_utils.f90 sourcefile~variable.f90->sourcefile~allocate.f90 sourcefile~c_utils.f90->sourcefile~signal_flag.f90 sourcefile~memory_stats.f90 memory_stats.F90 sourcefile~c_utils.f90->sourcefile~memory_stats.f90 sourcefile~system_info.f90 system_info.F90 sourcefile~c_utils.f90->sourcefile~system_info.f90

Files dependent on this one

sourcefile~~solver_factory.f90~~AfferentGraph sourcefile~solver_factory.f90 solver_factory.F90 sourcefile~solver.f90 solver.F90 sourcefile~solver.f90->sourcefile~solver_factory.f90 sourcefile~thermal_interface.f90 thermal_interface.F90 sourcefile~thermal_interface.f90->sourcefile~solver.f90 sourcefile~ftdss.f90 ftdss.F90 sourcefile~ftdss.f90->sourcefile~thermal_interface.f90 sourcefile~thermal.f90 thermal.F90 sourcefile~thermal.f90->sourcefile~thermal_interface.f90 sourcefile~thermal_3phase.f90 thermal_3phase.F90 sourcefile~thermal_3phase.f90->sourcefile~thermal_interface.f90

Source Code

module solver_solver_factory
    use :: module_input
    use :: module_matrix, only:type_crs
    use :: solver_solve
    implicit none

    public :: create_solver

    interface create_solver
        module procedure :: create_solver_crs
    end interface

contains
    function create_solver_crs(input, target_solver, target_matrix, num_node) result(solver)
        implicit none
        type(type_input), intent(in) :: input
        class(abst_solver), allocatable :: solver
        character(*), intent(in) :: target_solver
        type(type_crs), intent(in) :: target_matrix
        integer(int32), intent(in) :: num_node

        if (allocated(solver)) deallocate (solver)

        select case (trim(adjustl(target_solver)))
        case ('thermal')
            select case (input%basic%solver_settings%linear_solver%thermal%method)
            case ('direct')
                solver = type_solver_sparse_crs_lu(N=num_node, &
                                                   MAXFCT=1, &
                                                   MNUM=1, &
                                                   MTYPE=1, &
                                                   PHASE=13, &
                                                   NRHS=1, &
                                                   MSGVLV=0, &
                                                   a=target_matrix)
            case ('iterative')
                select case (input%basic%solver_settings%linear_solver%thermal%iterative_solver%solver_type)
                case (4)
                    solver = type_solver_sparse_crs_bicgstab( &
                             N=num_node, &
                             tolerance=input%basic%solver_settings%linear_solver%thermal%iterative_solver%tolerance, &
                             max_iterations=input%basic%solver_settings%linear_solver%thermal%iterative_solver%max_iterations, &
                             preconditioner=input%basic%solver_settings%linear_solver%thermal%iterative_solver%preconditioner_type)
                end select
            end select
        case ('hydraulic')
            select case (input%basic%solver_settings%linear_solver%hydraulic%method)
            case ('direct')
                solver = type_solver_sparse_crs_lu(N=num_node, &
                                                   MAXFCT=1, &
                                                   MNUM=1, &
                                                   MTYPE=1, &
                                                   PHASE=13, &
                                                   NRHS=1, &
                                                   MSGVLV=0, &
                                                   a=target_matrix)
            case ('iterative')
                select case (input%basic%solver_settings%linear_solver%hydraulic%iterative_solver%solver_type)
                case (4)
                    solver = type_solver_sparse_crs_bicgstab( &
                             N=num_node, &
                             tolerance=input%basic%solver_settings%linear_solver%hydraulic%iterative_solver%tolerance, &
                             max_iterations=input%basic%solver_settings%linear_solver%hydraulic%iterative_solver%max_iterations, &
                             preconditioner=input%basic%solver_settings%linear_solver%hydraulic%iterative_solver%preconditioner_type)
                end select
            end select
        end select

    end function create_solver_crs
end module solver_solver_factory