materials_manager.F90 Source File


This file depends on

sourcefile~~materials_manager.f90~~EfferentGraph sourcefile~materials_manager.f90 materials_manager.F90 sourcefile~calculate.f90 calculate.F90 sourcefile~materials_manager.f90->sourcefile~calculate.f90 sourcefile~input.f90 input.F90 sourcefile~materials_manager.f90->sourcefile~input.f90 sourcefile~density_interface.f90 density_interface.F90 sourcefile~calculate.f90->sourcefile~density_interface.f90 sourcefile~gcc_interface.f90 gcc_interface.F90 sourcefile~calculate.f90->sourcefile~gcc_interface.f90 sourcefile~hcf_interface.f90 hcf_interface.F90 sourcefile~calculate.f90->sourcefile~hcf_interface.f90 sourcefile~heat_capacity_interface.f90 heat_capacity_interface.F90 sourcefile~calculate.f90->sourcefile~heat_capacity_interface.f90 sourcefile~linalg.f90 linalg.F90 sourcefile~calculate.f90->sourcefile~linalg.f90 sourcefile~specific_heat_interface.f90 specific_heat_interface.F90 sourcefile~calculate.f90->sourcefile~specific_heat_interface.f90 sourcefile~thermal_conductivity_interface.f90 thermal_conductivity_interface.F90 sourcefile~calculate.f90->sourcefile~thermal_conductivity_interface.f90 sourcefile~input_interface.f90 input_interface.F90 sourcefile~input.f90->sourcefile~input_interface.f90 sourcefile~density_interface.f90->sourcefile~input.f90 sourcefile~core.f90 core.F90 sourcefile~density_interface.f90->sourcefile~core.f90 sourcefile~gcc_interface.f90->sourcefile~input.f90 sourcefile~gcc_interface.f90->sourcefile~core.f90 sourcefile~hcf_interface.f90->sourcefile~input.f90 sourcefile~hcf_interface.f90->sourcefile~core.f90 sourcefile~heat_capacity_interface.f90->sourcefile~density_interface.f90 sourcefile~heat_capacity_interface.f90->sourcefile~input_interface.f90 sourcefile~heat_capacity_interface.f90->sourcefile~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~matrix_ops.f90 matrix_ops.F90 sourcefile~linalg.f90->sourcefile~matrix_ops.f90 sourcefile~matvec.f90 matvec.F90 sourcefile~linalg.f90->sourcefile~matvec.f90 sourcefile~vector_ops.f90 vector_ops.F90 sourcefile~linalg.f90->sourcefile~vector_ops.f90 sourcefile~specific_heat_interface.f90->sourcefile~input.f90 sourcefile~specific_heat_interface.f90->sourcefile~core.f90 sourcefile~thermal_conductivity_interface.f90->sourcefile~input.f90 sourcefile~thermal_conductivity_interface.f90->sourcefile~core.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~matrix_ops.f90->sourcefile~core.f90 sourcefile~matvec.f90->sourcefile~core.f90 sourcefile~project_settings.f90->sourcefile~core.f90 sourcefile~allocate.f90->sourcefile~error.f90 sourcefile~deallocate.f90->sourcefile~error.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~matrix.f90 matrix.F90 sourcefile~types.f90->sourcefile~matrix.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~types.f90 sourcefile~vtk.f90->sourcefile~unique.f90 sourcefile~vtk.f90->sourcefile~vtk_constants.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~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~matrix_interface.f90 matrix_interface.F90 sourcefile~matrix.f90->sourcefile~matrix_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~matrix_coo.f90->sourcefile~allocate.f90 sourcefile~matrix_coo.f90->sourcefile~deallocate.f90 sourcefile~matrix_coo.f90->sourcefile~matrix_interface.f90 sourcefile~matrix_crs.f90->sourcefile~allocate.f90 sourcefile~matrix_crs.f90->sourcefile~deallocate.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_interface.f90 sourcefile~matrix_dense.f90->sourcefile~allocate.f90 sourcefile~matrix_dense.f90->sourcefile~deallocate.f90 sourcefile~matrix_dense.f90->sourcefile~matrix_interface.f90

Files dependent on this one

sourcefile~~materials_manager.f90~~AfferentGraph sourcefile~materials_manager.f90 materials_manager.F90 sourcefile~properties.f90 properties.F90 sourcefile~properties.f90->sourcefile~materials_manager.f90 sourcefile~properties_manager.f90 properties_manager.F90 sourcefile~properties.f90->sourcefile~properties_manager.f90 sourcefile~properties_manager.f90->sourcefile~materials_manager.f90 sourcefile~ftdss.f90 ftdss.F90 sourcefile~ftdss.f90->sourcefile~properties.f90 sourcefile~hydraulic.f90 hydraulic.F90 sourcefile~ftdss.f90->sourcefile~hydraulic.f90 sourcefile~output.f90 output.F90 sourcefile~ftdss.f90->sourcefile~output.f90 sourcefile~thermal.f90 thermal.F90 sourcefile~ftdss.f90->sourcefile~thermal.f90 sourcefile~hydraulic_assemble.f90 hydraulic_assemble.F90 sourcefile~hydraulic_assemble.f90->sourcefile~properties.f90 sourcefile~hydraulic_interface.f90 hydraulic_interface.F90 sourcefile~hydraulic_interface.f90->sourcefile~properties.f90 sourcefile~hydraulic_interface.f90->sourcefile~hydraulic_assemble.f90 sourcefile~output_interface.f90 output_interface.F90 sourcefile~output_interface.f90->sourcefile~properties.f90 sourcefile~thermal_interface.f90 thermal_interface.F90 sourcefile~thermal_interface.f90->sourcefile~properties.f90 sourcefile~hydraulic.f90->sourcefile~hydraulic_interface.f90 sourcefile~hydraulic_crs.f90 hydraulic_crs.F90 sourcefile~hydraulic_crs.f90->sourcefile~hydraulic_interface.f90 sourcefile~output.f90->sourcefile~output_interface.f90 sourcefile~output_base.f90 output_base.F90 sourcefile~output_base.f90->sourcefile~output_interface.f90 sourcefile~output_observation.f90 output_observation.F90 sourcefile~output_observation.f90->sourcefile~output_interface.f90 sourcefile~output_overall_base.f90 output_overall_base.F90 sourcefile~output_overall_base.f90->sourcefile~output_interface.f90 sourcefile~output_overall_vtk.f90 output_overall_vtk.F90 sourcefile~output_overall_vtk.f90->sourcefile~output_interface.f90 sourcefile~output_overall_vtu.f90 output_overall_vtu.F90 sourcefile~output_overall_vtu.f90->sourcefile~output_interface.f90 sourcefile~output_system_logger.f90 output_system_logger.F90 sourcefile~output_system_logger.f90->sourcefile~output_interface.f90 sourcefile~thermal.f90->sourcefile~thermal_interface.f90 sourcefile~thermal_crs.f90 thermal_crs.F90 sourcefile~thermal_crs.f90->sourcefile~thermal_interface.f90

Source Code

module properties_material_manager
    use, intrinsic :: iso_fortran_env, only: int32, real64
    use :: module_input, only:type_input
    use :: module_calculate, only: &
        holder_gccs, holder_wrfs, holder_dens, holder_sphs, holder_vhcs, holder_thcs, holder_hcfs, &
        abst_gcc, abst_wrf, abst_den, abst_sph, abst_vhc, abst_thc, abst_hcf
    implicit none

    public :: type_material_manager

    type :: type_material_manager
        private
        type(holder_thcs), allocatable :: thc(:)
        type(holder_dens), allocatable :: den(:)
        type(holder_sphs), allocatable :: sph(:)
        type(holder_vhcs), allocatable :: vhc(:)
        type(holder_gccs), allocatable :: gcc(:)
        type(holder_wrfs), allocatable :: wrf(:)
        type(holder_hcfs), allocatable :: hcf(:)

        integer(int32), allocatable :: region_id_map(:)
    contains
        procedure, public, pass(self) :: initialize => initialize_type_material_manager

        procedure, public, pass(self) :: get_thc => get_thc_ptr
        procedure, public, pass(self) :: get_den => get_den_ptr
        procedure, public, pass(self) :: get_sph => get_sph_ptr
        procedure, public, pass(self) :: get_vhc => get_vhc_ptr
        procedure, public, pass(self) :: get_gcc => get_gcc_ptr
        procedure, public, pass(self) :: get_wrf => get_wrf_ptr
        procedure, public, pass(self) :: get_hcf => get_hcf_ptr

    end type type_material_manager

contains

    ! 初期化(holder内部のinitialize呼ぶ)
    subroutine initialize_type_material_manager(self, input, ierr)
        class(type_material_manager), intent(inout) :: self
        type(type_input), intent(in) :: input
        integer(int32), intent(inout) :: ierr

        integer(int32) :: model_idx
        integer(int32) :: num_unique_regions
        integer(int32) :: max_region_id
        integer(int32), allocatable :: unique_material_ids(:)
        integer(int32) :: current_material_id

        ierr = 0
        call input%geometry%vtk%get_active_region_info(unique_material_ids, ierr)
        if (ierr /= 0) return
        if (.not. allocated(unique_material_ids) .or. size(unique_material_ids) == 0) then
            ierr = -1 ! エラーコード
            print *, "Error: No active material regions found."
            stop 1
        end if

        ! 実際に存在するユニーク領域の数を正とする
        num_unique_regions = size(unique_material_ids)
        max_region_id = maxval(unique_material_ids)

        if (input%basic%analysis_controls%calculate_thermal) then
            allocate (self%thc(num_unique_regions))
            allocate (self%den(num_unique_regions))
            allocate (self%sph(num_unique_regions))
            allocate (self%vhc(num_unique_regions))
            allocate (self%gcc(num_unique_regions))
            allocate (self%wrf(num_unique_regions))
        end if

        if (input%basic%analysis_controls%calculate_hydraulic) then
            allocate (self%hcf(num_unique_regions))
        end if

        ! region_idの最大値でマップ配列を確保し、0(無効値)で初期化
        allocate (self%region_id_map(max_region_id), source=0)

        do model_idx = 1, num_unique_regions
            current_material_id = unique_material_ids(model_idx)
            if (input%basic%analysis_controls%calculate_thermal) then
                call self%thc(model_idx)%initialize(input, current_material_id)
                call self%den(model_idx)%initialize(input, current_material_id)
                call self%sph(model_idx)%initialize(input, current_material_id)
                call self%vhc(model_idx)%initialize(input, current_material_id)
                call self%gcc(model_idx)%initialize(input, current_material_id)
                call self%wrf(model_idx)%initialize(input, current_material_id)
            end if

            if (input%basic%analysis_controls%calculate_hydraulic) then
                call self%hcf(model_idx)%initialize(input, current_material_id)
            end if

            self%region_id_map(current_material_id) = model_idx
        end do
    end subroutine initialize_type_material_manager

    ! THC getter
    function get_thc_ptr(self, region_id) result(thc_ptr)
        implicit none
        class(type_material_manager), intent(in), target :: self
        integer(int32), intent(in) :: region_id
        class(abst_thc), pointer :: thc_ptr

        integer(int32) :: model_index

#ifdef USE_DEBUG
        if (region_id < 1 .or. region_id > size(self%region_id_map)) then
            print *, "Error: Invalid region_id in get_thc_ptr:", region_id
            nullify (thc_ptr)
            stop 1
        end if
#endif

        model_index = self%region_id_map(region_id)

#ifdef USE_DEBUG
        if (model_index == 0) then
            print *, "Error: region_id not mapped in get_thc_ptr:", region_id
            nullify (thc_ptr)
            stop 1
        end if
#endif

        thc_ptr => self%thc(model_index)%p
    end function get_thc_ptr

    ! DEN getter
    function get_den_ptr(self, region_id) result(den_ptr)
        implicit none
        class(type_material_manager), intent(in), target :: self
        integer(int32), intent(in) :: region_id
        class(abst_den), pointer :: den_ptr

        integer(int32) :: model_index

#ifdef USE_DEBUG
        if (region_id < 1 .or. region_id > size(self%region_id_map)) then
            print *, "Error: Invalid region_id in get_den_ptr:", region_id
            nullify (den_ptr)
            stop 1
        end if
#endif

        model_index = self%region_id_map(region_id)

#ifdef USE_DEBUG
        if (model_index == 0) then
            print *, "Error: region_id not mapped in get_den_ptr:", region_id
            nullify (den_ptr)
            stop 1
        end if
#endif

        den_ptr => self%den(model_index)%p
    end function get_den_ptr

    ! SPH getter
    function get_sph_ptr(self, region_id) result(sph_ptr)
        implicit none
        class(type_material_manager), intent(in), target :: self
        integer(int32), intent(in) :: region_id
        class(abst_sph), pointer :: sph_ptr

        integer(int32) :: model_index

#ifdef USE_DEBUG
        if (region_id < 1 .or. region_id > size(self%region_id_map)) then
            print *, "Error: Invalid region_id in get_sph_ptr:", region_id
            nullify (sph_ptr)
            stop 1
        end if
#endif

        model_index = self%region_id_map(region_id)

#ifdef USE_DEBUG
        if (model_index == 0) then
            print *, "Error: region_id not mapped in get_sph_ptr:", region_id
            nullify (sph_ptr)
            stop 1
        end if
#endif

        sph_ptr => self%sph(model_index)%p
    end function get_sph_ptr

    ! VHC getter
    function get_vhc_ptr(self, region_id) result(vhc_ptr)
        implicit none
        class(type_material_manager), intent(in), target :: self
        integer(int32), intent(in) :: region_id
        class(abst_vhc), pointer :: vhc_ptr

        integer(int32) :: model_index

#ifdef USE_DEBUG
        if (region_id < 1 .or. region_id > size(self%region_id_map)) then
            print *, "Error: Invalid region_id in get_vhc_ptr:", region_id
            nullify (vhc_ptr)
            stop 1
        end if
#endif

        model_index = self%region_id_map(region_id)

#ifdef USE_DEBUG
        if (model_index == 0) then
            print *, "Error: region_id not mapped in get_vhc_ptr:", region_id
            nullify (vhc_ptr)
            stop 1
        end if
#endif

        vhc_ptr => self%vhc(model_index)%p
    end function get_vhc_ptr

    ! GCC getter
    function get_gcc_ptr(self, region_id) result(gcc_ptr)
        implicit none
        class(type_material_manager), intent(in), target :: self
        integer(int32), intent(in) :: region_id
        class(abst_gcc), pointer :: gcc_ptr

        integer(int32) :: model_index

#ifdef USE_DEBUG
        if (region_id < 1 .or. region_id > size(self%region_id_map)) then
            print *, "Error: Invalid region_id in get_gcc_ptr:", region_id
            nullify (gcc_ptr)
            stop 1
        end if
#endif

        model_index = self%region_id_map(region_id)

#ifdef USE_DEBUG
        if (model_index == 0) then
            print *, "Error: region_id not mapped in get_gcc_ptr:", region_id
            nullify (gcc_ptr)
            stop 1
        end if
#endif

        gcc_ptr => self%gcc(model_index)%p
    end function get_gcc_ptr

    ! WRF getter
    function get_wrf_ptr(self, region_id) result(wrf_ptr)
        implicit none
        class(type_material_manager), intent(in), target :: self
        integer(int32), intent(in) :: region_id
        class(abst_wrf), pointer :: wrf_ptr

        integer(int32) :: model_index

#ifdef USE_DEBUG
        if (region_id < 1 .or. region_id > size(self%region_id_map)) then
            print *, "Error: Invalid region_id in get_wrf_ptr:", region_id
            nullify (wrf_ptr)
            stop 1
        end if
#endif

        model_index = self%region_id_map(region_id)

#ifdef USE_DEBUG
        if (model_index == 0) then
            print *, "Error: region_id not mapped in get_wrf_ptr:", region_id
            nullify (wrf_ptr)
            stop 1
        end if
#endif

        wrf_ptr => self%wrf(model_index)%p
    end function get_wrf_ptr

    ! HCF getter
    function get_hcf_ptr(self, region_id) result(hcf_ptr)
        implicit none
        class(type_material_manager), intent(in), target :: self
        integer(int32), intent(in) :: region_id
        class(abst_hcf), pointer :: hcf_ptr

        integer(int32) :: model_index

#ifdef USE_DEBUG
        if (region_id < 1 .or. region_id > size(self%region_id_map)) then
            print *, "Error: Invalid region_id in get_hcf_ptr:", region_id
            nullify (hcf_ptr)
            stop 1
        end if
#endif

        model_index = self%region_id_map(region_id)

#ifdef USE_DEBUG
        if (model_index == 0) then
            print *, "Error: region_id not mapped in get_hcf_ptr:", region_id
            nullify (hcf_ptr)
            stop 1
        end if
#endif

        hcf_ptr => self%hcf(model_index)%p
    end function get_hcf_ptr

end module properties_material_manager