heat_capacity_3phase_apparent.F90 Source File


This file depends on

sourcefile~~heat_capacity_3phase_apparent.f90~~EfferentGraph sourcefile~heat_capacity_3phase_apparent.f90 heat_capacity_3phase_apparent.F90 sourcefile~heat_capacity_interface.f90 heat_capacity_interface.F90 sourcefile~heat_capacity_3phase_apparent.f90->sourcefile~heat_capacity_interface.f90 sourcefile~core.f90 core.F90 sourcefile~heat_capacity_interface.f90->sourcefile~core.f90 sourcefile~density_interface.f90 density_interface.F90 sourcefile~heat_capacity_interface.f90->sourcefile~density_interface.f90 sourcefile~input_interface.f90 input_interface.F90 sourcefile~heat_capacity_interface.f90->sourcefile~input_interface.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~density_interface.f90->sourcefile~core.f90 sourcefile~input.f90 input.F90 sourcefile~density_interface.f90->sourcefile~input.f90 sourcefile~input_interface.f90->sourcefile~core.f90 sourcefile~project_settings.f90 project_settings.F90 sourcefile~input_interface.f90->sourcefile~project_settings.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.f90->sourcefile~input_interface.f90 sourcefile~project_settings.f90->sourcefile~core.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~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

Source Code

submodule(calculate_volumetric_heat_capacity) calc_vhc_3phase_apparent
    implicit none
contains
    !----------------------------------------------------------------------------------------------------
    ! Construct each type of density
    !----------------------------------------------------------------------------------------------------
    module function construct_type_vhc_3phase_apparent(input, i_material) result(property)
        implicit none
        class(abst_vhc), allocatable :: property
        type(type_input), intent(in) :: input
        integer(int32), intent(in) :: i_material

        if (allocated(property)) deallocate (property)
        allocate (type_vhc_3phase_apparent :: property)

        property%material_id = i_material

        property%material1 = input%basic%materials(i_material)%thermal%density(1) * &
                             input%basic%materials(i_material)%thermal%specific_heat(1)
        property%material2 = input%basic%materials(i_material)%thermal%density(2) * &
                             input%basic%materials(i_material)%thermal%specific_heat(2)
        property%material3 = input%basic%materials(i_material)%thermal%density(3) * &
                             input%basic%materials(i_material)%thermal%specific_heat(3)

    end function construct_type_vhc_3phase_apparent

    module function calc_vhc_gauss_point_3phase_apparent_holder(self, state, DEN, LatentHeat, dQi_dT) result(VHC)
        implicit none
        class(type_vhc_3phase_apparent), intent(in) :: self
        type(type_gauss_point_state), intent(in) :: state
        type(holder_dens), intent(in), optional :: DEN
        real(real64), intent(in), optional :: LatentHeat
        real(real64), intent(in), optional :: dQi_dT
        real(real64) :: VHC

        real(real64) :: phi1, phi2, phi3

        phi1 = 1.0d0 - state%porosity
        phi2 = state%water_content
        phi3 = 1.0d0 - phi1 - phi2

        VHC = calc_vhc_3a(self%material1, phi1, self%material2, phi2, self%material3, phi3, &
                          LatentHeat, DEN%p%material3, dQi_dT)
    end function calc_vhc_gauss_point_3phase_apparent_holder

    module function calc_vhc_gauss_point_3phase_apparent_ptr(self, state, DEN, LatentHeat, dQi_dT) result(VHC)
        implicit none
        class(type_vhc_3phase_apparent), intent(in) :: self
        type(type_gauss_point_state), intent(in) :: state
        class(abst_den), pointer, intent(in), optional :: DEN
        real(real64), intent(in), optional :: LatentHeat
        real(real64), intent(in), optional :: dQi_dT
        real(real64) :: VHC

        real(real64) :: phi1, phi2, phi3

        phi1 = 1.0d0 - state%porosity
        phi2 = state%water_content
        phi3 = 1.0d0 - phi1 - phi2

        VHC = calc_vhc_3a(self%material1, phi1, self%material2, phi2, self%material3, phi3, &
                          LatentHeat, DEN%material3, dQi_dT)
    end function calc_vhc_gauss_point_3phase_apparent_ptr

end submodule calc_vhc_3phase_apparent