lfo.F90 Source File


This file depends on

sourcefile~~lfo.f90~~EfferentGraph sourcefile~lfo.f90 lfo.F90 sourcefile~multicoloring.f90 multicoloring.F90 sourcefile~lfo.f90->sourcefile~multicoloring.f90 sourcefile~adjacency_element.f90 adjacency_element.F90 sourcefile~multicoloring.f90->sourcefile~adjacency_element.f90 sourcefile~core.f90 core.F90 sourcefile~multicoloring.f90->sourcefile~core.f90 sourcefile~adjacency_element.f90->sourcefile~core.f90 sourcefile~element.f90 element.F90 sourcefile~adjacency_element.f90->sourcefile~element.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~allocate.f90->sourcefile~error.f90 sourcefile~deallocate.f90->sourcefile~error.f90 sourcefile~element.f90->sourcefile~core.f90 sourcefile~input.f90 input.F90 sourcefile~element.f90->sourcefile~input.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~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~array.f90->sourcefile~allocate.f90 sourcefile~array.f90->sourcefile~deallocate.f90 sourcefile~input_interface.f90 input_interface.F90 sourcefile~input.f90->sourcefile~input_interface.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 sourcefile~input_interface.f90->sourcefile~core.f90 sourcefile~project_settings.f90 project_settings.F90 sourcefile~input_interface.f90->sourcefile~project_settings.f90 sourcefile~project_settings.f90->sourcefile~core.f90

Source Code

submodule(domain_multicoloring) multicoloring_lfo
contains
! ==============================================================================
! LFO (Largest First Ordering) アルゴリズムを実行し、彩色結果をselfに格納する
! ==============================================================================
    module subroutine coloring_lfo(self, graph)
        implicit none
        class(type_coloring), intent(inout) :: self
        type(type_crs_adjacency_element), intent(in) :: graph

        integer(int32) :: num_elements
        integer(int32) :: i, j, k, current_node, current_color
        integer(int32), allocatable :: degrees(:), sorted_indices(:)
        integer(int32), allocatable :: neighbor_colors(:), neighbors(:)

        ! --- 1. 初期化 ---
        num_elements = graph%get_num_elements()
        if (allocated(self%color)) call deallocate_array(self%color)
        call allocate_array(self%color, length=num_elements)
        self%color = 0

        call allocate_array(degrees, length=num_elements)
        call allocate_array(sorted_indices, length=num_elements)

        ! --- 2. 次数の計算と順序付け ---
        do i = 1, num_elements
            degrees(i) = graph%get_degree(i)
            sorted_indices(i) = i ! 初期インデックスを格納
        end do

        ! 次数に基づいてインデックスを降順ソート (単純なバブルソートで実装)
        do i = 1, num_elements - 1
            do j = i + 1, num_elements
                if (degrees(sorted_indices(i)) < degrees(sorted_indices(j))) then
                    k = sorted_indices(i)
                    sorted_indices(i) = sorted_indices(j)
                    sorted_indices(j) = k
                end if
            end do
        end do

        deallocate (degrees)

        ! --- 3. ソートされた順序で彩色 ---
        do i = 1, num_elements
            current_node = sorted_indices(i)

            ! 利用可能な最小の色を決定
            neighbors = graph%get_neighbors(current_node)
            if (size(neighbors) > 0) then
                allocate (neighbor_colors(size(neighbors)))
                k = 0
                do j = 1, size(neighbors)
                    if (self%color(neighbors(j)) /= 0) then
                        k = k + 1
                        neighbor_colors(k) = self%color(neighbors(j))
                    end if
                end do
            else
                k = 0
            end if

            current_color = 1
            if (k > 0) then
                do
                    if (count(neighbor_colors(1:k) == current_color) == 0) exit
                    current_color = current_color + 1
                end do
            end if
            if (allocated(neighbor_colors)) deallocate (neighbor_colors)

            deallocate (neighbors)

            ! 彩色
            self%color(current_node) = current_color
        end do

        call self%populate()

        deallocate (sorted_indices)

        ! --- 4. 結果の集計 ---

    end subroutine coloring_lfo

end submodule multicoloring_lfo