initialize_type_crs Subroutine

private subroutine initialize_type_crs(self, domain)

Type Bound

type_crs

Arguments

Type IntentOptional Attributes Name
class(type_crs), intent(inout) :: self
type(type_domain), intent(inout) :: domain

Calls

proc~~initialize_type_crs~~CallsGraph proc~initialize_type_crs type_crs%initialize_type_crs interface~allocate_array allocate_array proc~initialize_type_crs->interface~allocate_array interface~deallocate_array deallocate_array proc~initialize_type_crs->interface~deallocate_array none~to_reordered type_reordering%to_reordered proc~initialize_type_crs->none~to_reordered proc~destory_coo type_coo%destory_coo proc~initialize_type_crs->proc~destory_coo proc~get_num_nodes type_domain%get_num_nodes proc~initialize_type_crs->proc~get_num_nodes proc~initialize_type_coo type_coo%initialize_type_coo proc~initialize_type_crs->proc~initialize_type_coo proc~allocate_rank1_int16 allocate_rank1_int16 interface~allocate_array->proc~allocate_rank1_int16 proc~allocate_rank1_int32 allocate_rank1_int32 interface~allocate_array->proc~allocate_rank1_int32 proc~allocate_rank1_int64 allocate_rank1_int64 interface~allocate_array->proc~allocate_rank1_int64 proc~allocate_rank1_int8 allocate_rank1_int8 interface~allocate_array->proc~allocate_rank1_int8 proc~allocate_rank1_logical1 allocate_rank1_logical1 interface~allocate_array->proc~allocate_rank1_logical1 proc~allocate_rank1_logical4 allocate_rank1_logical4 interface~allocate_array->proc~allocate_rank1_logical4 proc~allocate_rank1_logical8 allocate_rank1_logical8 interface~allocate_array->proc~allocate_rank1_logical8 proc~allocate_rank1_real128 allocate_rank1_real128 interface~allocate_array->proc~allocate_rank1_real128 proc~allocate_rank1_real32 allocate_rank1_real32 interface~allocate_array->proc~allocate_rank1_real32 proc~allocate_rank1_real64 allocate_rank1_real64 interface~allocate_array->proc~allocate_rank1_real64 proc~allocate_rank2_int16 allocate_rank2_int16 interface~allocate_array->proc~allocate_rank2_int16 proc~allocate_rank2_int32 allocate_rank2_int32 interface~allocate_array->proc~allocate_rank2_int32 proc~allocate_rank2_int64 allocate_rank2_int64 interface~allocate_array->proc~allocate_rank2_int64 proc~allocate_rank2_int8 allocate_rank2_int8 interface~allocate_array->proc~allocate_rank2_int8 proc~allocate_rank2_logical1 allocate_rank2_logical1 interface~allocate_array->proc~allocate_rank2_logical1 proc~allocate_rank2_logical4 allocate_rank2_logical4 interface~allocate_array->proc~allocate_rank2_logical4 proc~allocate_rank2_logical8 allocate_rank2_logical8 interface~allocate_array->proc~allocate_rank2_logical8 proc~allocate_rank2_real128 allocate_rank2_real128 interface~allocate_array->proc~allocate_rank2_real128 proc~allocate_rank2_real32 allocate_rank2_real32 interface~allocate_array->proc~allocate_rank2_real32 proc~allocate_rank2_real64 allocate_rank2_real64 interface~allocate_array->proc~allocate_rank2_real64 proc~deallocate_rank1_int32 deallocate_rank1_int32 interface~deallocate_array->proc~deallocate_rank1_int32 proc~deallocate_rank1_int64 deallocate_rank1_int64 interface~deallocate_array->proc~deallocate_rank1_int64 proc~deallocate_rank1_int8 deallocate_rank1_int8 interface~deallocate_array->proc~deallocate_rank1_int8 proc~deallocate_rank1_logical1 deallocate_rank1_logical1 interface~deallocate_array->proc~deallocate_rank1_logical1 proc~deallocate_rank1_logical4 deallocate_rank1_logical4 interface~deallocate_array->proc~deallocate_rank1_logical4 proc~deallocate_rank1_logical8 deallocate_rank1_logical8 interface~deallocate_array->proc~deallocate_rank1_logical8 proc~deallocate_rank1_real128 deallocate_rank1_real128 interface~deallocate_array->proc~deallocate_rank1_real128 proc~deallocate_rank1_real32 deallocate_rank1_real32 interface~deallocate_array->proc~deallocate_rank1_real32 proc~deallocate_rank1_real64 deallocate_rank1_real64 interface~deallocate_array->proc~deallocate_rank1_real64 proc~deallocate_rank2_int32 deallocate_rank2_int32 interface~deallocate_array->proc~deallocate_rank2_int32 proc~deallocate_rank2_int64 deallocate_rank2_int64 interface~deallocate_array->proc~deallocate_rank2_int64 proc~deallocate_rank2_int8 deallocate_rank2_int8 interface~deallocate_array->proc~deallocate_rank2_int8 proc~deallocate_rank2_logical1 deallocate_rank2_logical1 interface~deallocate_array->proc~deallocate_rank2_logical1 proc~deallocate_rank2_logical4 deallocate_rank2_logical4 interface~deallocate_array->proc~deallocate_rank2_logical4 proc~deallocate_rank2_logical8 deallocate_rank2_logical8 interface~deallocate_array->proc~deallocate_rank2_logical8 proc~deallocate_rank2_real128 deallocate_rank2_real128 interface~deallocate_array->proc~deallocate_rank2_real128 proc~deallocate_rank2_real32 deallocate_rank2_real32 interface~deallocate_array->proc~deallocate_rank2_real32 proc~deallocate_rank2_real64 deallocate_rank2_real64 interface~deallocate_array->proc~deallocate_rank2_real64 interface~to_reordered_index type_reordering%to_reordered_index none~to_reordered->interface~to_reordered_index interface~to_reordered_indices type_reordering%to_reordered_indices none~to_reordered->interface~to_reordered_indices proc~destory_coo->interface~deallocate_array proc~initialize_type_coo->interface~allocate_array proc~initialize_type_coo->interface~deallocate_array proc~initialize_type_coo->proc~get_num_nodes interface~unique unique proc~initialize_type_coo->interface~unique proc~get_computation_dimension type_domain%get_computation_dimension proc~initialize_type_coo->proc~get_computation_dimension proc~get_num_elements type_domain%get_num_elements proc~initialize_type_coo->proc~get_num_elements proc~get_num_sides type_domain%get_num_sides proc~initialize_type_coo->proc~get_num_sides proc~to_reordered_index to_reordered_index interface~to_reordered_index->proc~to_reordered_index proc~to_reordered_indices to_reordered_indices interface~to_reordered_indices->proc~to_reordered_indices proc~unique_int16 unique_int16 interface~unique->proc~unique_int16 proc~unique_int32 unique_int32 interface~unique->proc~unique_int32 proc~unique_int64 unique_int64 interface~unique->proc~unique_int64 proc~unique_int8 unique_int8 interface~unique->proc~unique_int8 proc~error_message error_message proc~allocate_rank1_int16->proc~error_message proc~allocate_rank1_int32->proc~error_message proc~allocate_rank1_int64->proc~error_message proc~allocate_rank1_int8->proc~error_message proc~allocate_rank1_logical1->proc~error_message proc~allocate_rank1_logical4->proc~error_message proc~allocate_rank1_logical8->proc~error_message proc~allocate_rank1_real128->proc~error_message proc~allocate_rank1_real32->proc~error_message proc~allocate_rank1_real64->proc~error_message proc~allocate_rank2_int16->proc~error_message proc~allocate_rank2_int32->proc~error_message proc~allocate_rank2_int64->proc~error_message proc~allocate_rank2_int8->proc~error_message proc~allocate_rank2_logical1->proc~error_message proc~allocate_rank2_logical4->proc~error_message proc~allocate_rank2_logical8->proc~error_message proc~allocate_rank2_real128->proc~error_message proc~allocate_rank2_real32->proc~error_message proc~allocate_rank2_real64->proc~error_message proc~deallocate_rank1_int32->proc~error_message proc~deallocate_rank1_int64->proc~error_message proc~deallocate_rank1_int8->proc~error_message proc~deallocate_rank1_logical1->proc~error_message proc~deallocate_rank1_logical4->proc~error_message proc~deallocate_rank1_logical8->proc~error_message proc~deallocate_rank1_real128->proc~error_message proc~deallocate_rank1_real32->proc~error_message proc~deallocate_rank1_real64->proc~error_message proc~deallocate_rank2_int32->proc~error_message proc~deallocate_rank2_int64->proc~error_message proc~deallocate_rank2_int8->proc~error_message proc~deallocate_rank2_logical1->proc~error_message proc~deallocate_rank2_logical4->proc~error_message proc~deallocate_rank2_logical8->proc~error_message proc~deallocate_rank2_real128->proc~error_message proc~deallocate_rank2_real32->proc~error_message proc~deallocate_rank2_real64->proc~error_message log_error log_error proc~error_message->log_error proc~unique_int16->interface~allocate_array sort sort proc~unique_int16->sort proc~unique_int32->interface~allocate_array proc~unique_int32->sort proc~unique_int64->interface~allocate_array proc~unique_int64->sort proc~unique_int8->interface~allocate_array proc~unique_int8->sort

Called by

proc~~initialize_type_crs~~CalledByGraph proc~initialize_type_crs type_crs%initialize_type_crs proc~construct_type_thermal_3phase_2d construct_type_thermal_3phase_2d proc~construct_type_thermal_3phase_2d->proc~initialize_type_crs interface~construct_type_thermal_3phase_2d construct_type_thermal_3phase_2d interface~construct_type_thermal_3phase_2d->proc~construct_type_thermal_3phase_2d interface~type_thermal_3phase_2d type_thermal_3phase_2d interface~type_thermal_3phase_2d->interface~construct_type_thermal_3phase_2d

Source Code

    subroutine initialize_type_crs(self, domain)
        implicit none
        class(type_crs), intent(inout) :: self
        type(type_domain), intent(inout) :: domain

        ! --- ローカル変数 ---
        type(type_coo) :: coo
        integer(int32) :: i, r, pos
        integer(int32), allocatable :: rcm_row(:), rcm_col(:)
        integer(int32), allocatable :: write_pos(:)

        ! =================================================================
        ! Step 1: coo%initializeを呼び出し、COO行列を作成する
        ! =================================================================
        call coo%initialize(domain)

        self%nnz = coo%nnz
        self%num_row = domain%get_num_nodes()
        self%num_ptr = self%num_row + 1

        if (self%nnz == 0) then
            call allocate_array(self%ptr, self%num_ptr)
            self%ptr = 1
            return
        end if

        ! =================================================================
        ! Step 2: COOの行・列インデックスをRCM番号に変換
        ! =================================================================
        call allocate_array(rcm_row, self%nnz)
        call allocate_array(rcm_col, self%nnz)

        call domain%reordering%to_reordered(coo%row, rcm_row)
        call domain%reordering%to_reordered(coo%col, rcm_col)

        ! =================================================================
        ! Step 3: RCM適用後のデータから直接CRS行列を構築
        ! =================================================================
        call allocate_array(self%ptr, self%num_ptr)
        call allocate_array(self%ind, self%nnz)
        call allocate_array(self%val, self%nnz)
        self%ptr = 0

        ! -- Pass 1: 各行の非ゼロ要素数をカウント
        do i = 1, self%nnz
            self%ptr(rcm_row(i) + 1) = self%ptr(rcm_row(i) + 1) + 1
        end do

        ! -- Pass 2: 累積和を計算して、各行の開始ポインタを決定
        self%ptr(1) = 1
        do i = 2, self%num_ptr
            self%ptr(i) = self%ptr(i) + self%ptr(i - 1)
        end do

        ! -- Pass 3: indとvalを正しい位置に配置
        call allocate_array(write_pos, self%num_row)
        write_pos(:) = self%ptr(1:self%num_row)

        do i = 1, self%nnz
            r = rcm_row(i)
            pos = write_pos(r)

            self%ind(pos) = rcm_col(i)
            self%val(pos) = 0.0d0

            write_pos(r) = pos + 1
        end do

        call deallocate_array(write_pos)
        call deallocate_array(rcm_row)
        call deallocate_array(rcm_col)

        call coo%destory()

    end subroutine initialize_type_crs