hcf_base.F90 Source File


This file depends on

sourcefile~~hcf_base.f90~~EfferentGraph sourcefile~hcf_base.f90 hcf_base.F90 sourcefile~hcf_interface.f90 hcf_interface.F90 sourcefile~hcf_base.f90->sourcefile~hcf_interface.f90 sourcefile~core.f90 core.F90 sourcefile~hcf_interface.f90->sourcefile~core.f90 sourcefile~input.f90 input.F90 sourcefile~hcf_interface.f90->sourcefile~input.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~input_interface.f90 input_interface.F90 sourcefile~input.f90->sourcefile~input_interface.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~input_interface.f90->sourcefile~core.f90 sourcefile~project_settings.f90 project_settings.F90 sourcefile~input_interface.f90->sourcefile~project_settings.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~project_settings.f90->sourcefile~core.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

Source Code

submodule(calculate_hcf) calculate_hcf_base
    implicit none
contains

    module subroutine initialize_holder_hcfs(self, input, material_id)
        implicit none
        class(holder_hcfs), intent(inout) :: self
        type(type_input), intent(in) :: input
        integer(int32), intent(in) :: material_id

        select case (input%basic%materials(material_id)%hydraulic%model_number)
        case (1)
            self%p = create_type_hcf_base(input, material_id)
        case (2)
            self%p = create_type_hcf_impedance(input, material_id)
        case (3)
            self%p = create_type_hcf_viscosity(input, material_id)
        case (4)
            self%p = create_type_hcf_base_impedance(input, material_id)
        case (5)
            self%p = create_type_hcf_base_viscosity(input, material_id)
        case (6)
            self%p = create_type_hcf_impedance_viscosity(input, material_id)
        case (7)
            self%p = create_type_hcf_base_impedance_viscosity(input, material_id)
        end select

    end subroutine initialize_holder_hcfs

    module pure elemental function calc_kflh_base(self, state) result(kflh)
        implicit none
        class(type_hcf_base), intent(in) :: self
        type(type_state), intent(in) :: state
        real(real64) :: kflh

        kflh = self%k_s * self%base%calc_kr(state%pressure)

    end function calc_kflh_base

    module pure elemental function calc_kflh_impedance(self, state) result(kflh)
        implicit none
        class(type_hcf_impedance), intent(in) :: self
        type(type_state), intent(in) :: state
        real(real64) :: kflh

        kflh = self%k_s * self%impedance%calc_impedance(state%ice_content)

    end function calc_kflh_impedance

    module pure elemental function calc_kflh_viscosity(self, state) result(kflh)
        implicit none
        class(type_hcf_viscosity), intent(in) :: self
        type(type_state), intent(in) :: state
        real(real64) :: kflh

        kflh = self%k_s * self%viscosity%calc_viscosity(state%temperature)

    end function calc_kflh_viscosity

    module pure elemental function calc_kflh_base_impedance(self, state) result(kflh)
        implicit none
        class(type_hcf_base_impedance), intent(in) :: self
        type(type_state), intent(in) :: state
        real(real64) :: kflh

        kflh = self%k_s * self%base%calc_kr(state%pressure) & !&
                        * self%impedance%calc_impedance(state%ice_content)

    end function calc_kflh_base_impedance

    module pure elemental function calc_kflh_base_viscosity(self, state) result(kflh)
        implicit none
        class(type_hcf_base_viscosity), intent(in) :: self
        type(type_state), intent(in) :: state
        real(real64) :: kflh

        kflh = self%k_s * self%base%calc_kr(state%pressure) & !&
                        * self%viscosity%calc_viscosity(state%temperature)

    end function calc_kflh_base_viscosity

    module pure elemental function calc_kflh_impedance_viscosity(self, state) result(kflh)
        implicit none
        class(type_hcf_impedance_viscosity), intent(in) :: self
        type(type_state), intent(in) :: state
        real(real64) :: kflh

        kflh = self%k_s * self%impedance%calc_impedance(state%ice_content) & !&
                        * self%viscosity%calc_viscosity(state%temperature)

    end function calc_kflh_impedance_viscosity

    module pure elemental function calc_kflh_base_impedance_viscosity(self, state) result(kflh)
        implicit none
        class(type_hcf_base_impedance_viscosity), intent(in) :: self
        type(type_state), intent(in) :: state
        real(real64) :: kflh

        kflh = self%k_s * self%base%calc_kr(state%pressure) & !&
                        * self%impedance%calc_impedance(state%ice_content) & !&
                        * self%viscosity%calc_viscosity(state%temperature)

    end function calc_kflh_base_impedance_viscosity

end submodule calculate_hcf_base