element.F90 Source File


This file depends on

sourcefile~~element.f90~~EfferentGraph sourcefile~element.f90 element.F90 sourcefile~core.f90 core.F90 sourcefile~element.f90->sourcefile~core.f90 sourcefile~input.f90 input.F90 sourcefile~element.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~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~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

Files dependent on this one

sourcefile~~element.f90~~AfferentGraph sourcefile~element.f90 element.F90 sourcefile~adjacency_element.f90 adjacency_element.F90 sourcefile~adjacency_element.f90->sourcefile~element.f90 sourcefile~domain.f90 domain.F90 sourcefile~domain.f90->sourcefile~element.f90 sourcefile~domain_manager.f90 domain_manager.F90 sourcefile~domain.f90->sourcefile~domain_manager.f90 sourcefile~element_factory.f90 element_factory.F90 sourcefile~domain.f90->sourcefile~element_factory.f90 sourcefile~reordering.f90 reordering.F90 sourcefile~domain.f90->sourcefile~reordering.f90 sourcefile~adjacency.f90 adjacency.F90 sourcefile~domain.f90->sourcefile~adjacency.f90 sourcefile~multicoloring.f90 multicoloring.F90 sourcefile~domain.f90->sourcefile~multicoloring.f90 sourcefile~domain_manager.f90->sourcefile~element.f90 sourcefile~domain_manager.f90->sourcefile~element_factory.f90 sourcefile~domain_manager.f90->sourcefile~reordering.f90 sourcefile~domain_manager.f90->sourcefile~adjacency.f90 sourcefile~domain_manager.f90->sourcefile~multicoloring.f90 sourcefile~element_factory.f90->sourcefile~element.f90 sourcefile~element_square_first.f90 element_square_first.F90 sourcefile~element_square_first.f90->sourcefile~element.f90 sourcefile~element_square_second.f90 element_square_second.F90 sourcefile~element_square_second.f90->sourcefile~element.f90 sourcefile~element_triangle_first.f90 element_triangle_first.F90 sourcefile~element_triangle_first.f90->sourcefile~element.f90 sourcefile~element_triangle_second.f90 element_triangle_second.F90 sourcefile~element_triangle_second.f90->sourcefile~element.f90 sourcefile~reordering.f90->sourcefile~element.f90 sourcefile~adjacency.f90->sourcefile~adjacency_element.f90 sourcefile~boundary_interface.f90 boundary_interface.F90 sourcefile~boundary_interface.f90->sourcefile~domain.f90 sourcefile~matrix.f90 matrix.F90 sourcefile~boundary_interface.f90->sourcefile~matrix.f90 sourcefile~boundary_manager.f90 boundary_manager.F90 sourcefile~boundary_manager.f90->sourcefile~domain.f90 sourcefile~boundary_manager.f90->sourcefile~boundary_interface.f90 sourcefile~boundary_manager.f90->sourcefile~matrix.f90 sourcefile~ftdss.f90 ftdss.F90 sourcefile~ftdss.f90->sourcefile~domain.f90 sourcefile~thermal_interface.f90 thermal_interface.F90 sourcefile~ftdss.f90->sourcefile~thermal_interface.f90 sourcefile~boundary.f90 boundary.F90 sourcefile~ftdss.f90->sourcefile~boundary.f90 sourcefile~initial.f90 initial.F90 sourcefile~ftdss.f90->sourcefile~initial.f90 sourcefile~output.f90 output.F90 sourcefile~ftdss.f90->sourcefile~output.f90 sourcefile~initial_interface.f90 initial_interface.F90 sourcefile~initial_interface.f90->sourcefile~domain.f90 sourcefile~initial_interface.f90->sourcefile~boundary.f90 sourcefile~initial_manager.f90 initial_manager.F90 sourcefile~initial_manager.f90->sourcefile~domain.f90 sourcefile~initial_manager.f90->sourcefile~initial_interface.f90 sourcefile~matrix_base.f90 matrix_base.F90 sourcefile~matrix_base.f90->sourcefile~domain.f90 sourcefile~matrix_coo.f90 matrix_coo.F90 sourcefile~matrix_coo.f90->sourcefile~domain.f90 sourcefile~matrix_coo.f90->sourcefile~matrix_base.f90 sourcefile~matrix_crs.f90 matrix_crs.F90 sourcefile~matrix_crs.f90->sourcefile~domain.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_base.f90 sourcefile~matrix_crs.f90->sourcefile~matrix_coo.f90 sourcefile~matrix_dense.f90 matrix_dense.F90 sourcefile~matrix_dense.f90->sourcefile~domain.f90 sourcefile~matrix_dense.f90->sourcefile~matrix_base.f90 sourcefile~methods.f90 methods.F90 sourcefile~methods.f90->sourcefile~reordering.f90 sourcefile~multicoloring.f90->sourcefile~adjacency_element.f90 sourcefile~output_interface.f90 output_interface.F90 sourcefile~output_interface.f90->sourcefile~domain.f90 sourcefile~output_interface.f90->sourcefile~matrix.f90 sourcefile~thermal_interface.f90->sourcefile~domain.f90 sourcefile~thermal_interface.f90->sourcefile~boundary.f90 sourcefile~thermal_interface.f90->sourcefile~matrix.f90 sourcefile~solver.f90 solver.F90 sourcefile~thermal_interface.f90->sourcefile~solver.f90 sourcefile~to_original.f90 to_original.F90 sourcefile~to_original.f90->sourcefile~reordering.f90 sourcefile~to_reordered.f90 to_reordered.F90 sourcefile~to_reordered.f90->sourcefile~reordering.f90 sourcefile~boundary.f90->sourcefile~boundary_interface.f90 sourcefile~boundary.f90->sourcefile~boundary_manager.f90 sourcefile~boundary_adiabatic.f90 boundary_adiabatic.F90 sourcefile~boundary_adiabatic.f90->sourcefile~boundary_interface.f90 sourcefile~boundary_base.f90 boundary_base.F90 sourcefile~boundary_base.f90->sourcefile~boundary_interface.f90 sourcefile~boundary_dirichlet.f90 boundary_dirichlet.F90 sourcefile~boundary_dirichlet.f90->sourcefile~boundary_interface.f90 sourcefile~dsatur.f90 dsatur.F90 sourcefile~dsatur.f90->sourcefile~multicoloring.f90 sourcefile~initial.f90->sourcefile~initial_interface.f90 sourcefile~initial.f90->sourcefile~initial_manager.f90 sourcefile~initial_laplace.f90 initial_laplace.F90 sourcefile~initial_laplace.f90->sourcefile~initial_interface.f90 sourcefile~initial_uniform.f90 initial_uniform.F90 sourcefile~initial_uniform.f90->sourcefile~initial_interface.f90 sourcefile~lfo.f90 lfo.F90 sourcefile~lfo.f90->sourcefile~multicoloring.f90 sourcefile~matrix.f90->sourcefile~matrix_base.f90 sourcefile~matrix.f90->sourcefile~matrix_coo.f90 sourcefile~matrix.f90->sourcefile~matrix_crs.f90 sourcefile~matrix.f90->sourcefile~matrix_dense.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 thermal.F90 sourcefile~thermal.f90->sourcefile~thermal_interface.f90 sourcefile~thermal_3phase.f90 thermal_3phase.F90 sourcefile~thermal_3phase.f90->sourcefile~thermal_interface.f90 sourcefile~welch_powell.f90 welch_powell.F90 sourcefile~welch_powell.f90->sourcefile~multicoloring.f90 sourcefile~solver_factory.f90 solver_factory.F90 sourcefile~solver_factory.f90->sourcefile~matrix.f90 sourcefile~solver.f90->sourcefile~solver_factory.f90

Source Code

module domain_element
    use, intrinsic :: iso_fortran_env, only: int32, real64
    use :: stdlib_logger
    use :: module_core, only:type_dp_3d, type_dp_pointer, type_vtk_cell, allocate_array, deallocate_array
    use :: module_input, only:type_geometry_settings
    implicit none
    private

    public :: abst_element
    public :: type_triangle_first
    public :: type_triangle_second
    public :: type_square_first
    public :: type_square_second
    public :: holder_elements

    public :: interpolate_reordered_triangle_first
    public :: interpolate_reordered_triangle_second
    public :: interpolate_reordered_square_first
    public :: interpolate_reordered_square_second

    public :: interpolate
    public :: interpolate_reordered

    public :: get_connectivity
    public :: get_connectivity_reordered

    !--------------------------------------------------------------------------------------
    ! Holder for polymorphic element objects
    !--------------------------------------------------------------------------------------
    type :: holder_elements
        class(abst_element), allocatable :: e
    end type holder_elements

    !--------------------------------------------------------------------------------------
    !   Abstract base type for 2D elements
    !--------------------------------------------------------------------------------------
    type, abstract :: abst_element
        integer(int32), private :: id !! Element ID
        integer(int32), private :: type !! Element type (5: triangle 1st, 9: square 1st)
        integer(int32), private :: num_nodes !! Number of nodes in the element
        integer(int32), private :: group !! Element group number
        integer(int32), private :: dimension
        integer(int32), private :: order
        integer(int32), allocatable :: connectivity(:) !! connectivity information
        integer(int32), allocatable :: connectivity_reordered(:) !! reordered connectivity information
        type(type_dp_pointer), allocatable :: x(:) !! X coordinate
        type(type_dp_pointer), allocatable :: y(:) !! Y coordinate
        type(type_dp_pointer), allocatable :: z(:) !! Z coordinate
        !----------------------------------------------------------------------------------
        ! Gauss Quadrature points and weights
        !  - Gauss Quadrature points are defined in the local coordinate system
        !  - The number of Gauss points is determined by the element type
        !  - The weights are used for numerical integration over the element
        !  - The Gauss points are used to evaluate the shape functions and their derivatives
        !----------------------------------------------------------------------------------
        integer(int32) :: num_gauss !! Number of Gauss Quadrature points
        real(real64), allocatable :: weight(:) !! Gauss weight
        real(real64), allocatable :: gauss(:, :) !! Gauss Quadrature points Coordinate
        !----------------------------------------------------------------------------------
        ! Interpolation functions
        !----------------------------------------------------------------------------------
        procedure(abst_interpolate),      pass(self), pointer :: interpolate => null() !&
        procedure(abst_get_connectivity), pass(self), pointer :: get_connectivity => null() !&
    contains
        procedure(abst_get_id),        pass(self), deferred :: get_id !&
        procedure(abst_get_type),      pass(self), deferred :: get_type !&
        procedure(abst_get_num_nodes), pass(self), deferred :: get_num_nodes !&
        procedure(abst_get_group),     pass(self), deferred :: get_group !&
        procedure(abst_get_order),     pass(self), deferred :: get_order !&
        procedure(abst_get_dimension), pass(self), deferred :: get_dimension !&
        procedure(abst_get_num_gauss), pass(self), deferred :: get_num_gauss !&
        !----------------------------------------------------------------------------------
        procedure(abst_psi),           pass(self), deferred :: psi !&
        procedure(abst_dpsi_dxi),      pass(self), deferred :: dpsi_dxi !&
        procedure(abst_dpsi_deta),     pass(self), deferred :: dpsi_deta !&
        procedure(abst_jacobian),      pass(self), deferred :: jacobian !&
        procedure(abst_jacobian_det),  pass(self), deferred :: jacobian_det !&
        procedure(abst_is_inside),     pass(self), deferred :: is_inside !&
    end type abst_element

    !--------------------------------------------------------------------------------------
    !   Triangle First Order Element Type
    !--------------------------------------------------------------------------------------
    type, extends(abst_element) :: type_triangle_first
    contains
        procedure, pass(self) :: get_id        => get_id_triangle_first !&
        procedure, pass(self) :: get_type      => get_type_triangle_first !&
        procedure, pass(self) :: get_num_nodes => get_num_nodes_triangle_first !&
        procedure, pass(self) :: get_group     => get_group_triangle_first !&
        procedure, pass(self) :: get_order     => get_order_triangle_first !&
        procedure, pass(self) :: get_dimension => get_dimension_triangle_first !&
        procedure, pass(self) :: get_num_gauss => get_num_gauss_triangle_first !&
        !----------------------------------------------------------------------------------
        procedure, pass(self) :: psi           => psi_triangle_first !&
        procedure, pass(self) :: dpsi_dxi      => dpsi_dxi_triangle_first !&
        procedure, pass(self) :: dpsi_deta     => dpsi_deta_triangle_first !&
        procedure, pass(self) :: jacobian      => jacobian_triangle_first !&
        procedure, pass(self) :: jacobian_det  => jacobian_det_triangle_first !&
        procedure, pass(self) :: is_inside     => is_in_triangle_first !&
    end type type_triangle_first

    !--------------------------------------------------------------------------------------
    !   Square First Order Element Type
    !--------------------------------------------------------------------------------------
    type, extends(abst_element) :: type_square_first
    contains
        procedure, pass(self) :: get_id        => get_id_square_first !&
        procedure, pass(self) :: get_type      => get_type_square_first !&
        procedure, pass(self) :: get_num_nodes => get_num_nodes_square_first !&
        procedure, pass(self) :: get_group     => get_group_square_first !&
        procedure, pass(self) :: get_order     => get_order_square_first !&
        procedure, pass(self) :: get_dimension => get_dimension_square_first !&
        procedure, pass(self) :: get_num_gauss => get_num_gauss_square_first !&
        !----------------------------------------------------------------------------------
        procedure, pass(self) :: psi           => psi_square_first !&
        procedure, pass(self) :: dpsi_dxi      => dpsi_dxi_square_first !&
        procedure, pass(self) :: dpsi_deta     => dpsi_deta_square_first !&
        procedure, pass(self) :: jacobian      => jacobian_square_first !&
        procedure, pass(self) :: jacobian_det  => jacobian_det_square_first !&
        procedure, pass(self) :: is_inside     => is_in_square_first !&
    end type type_square_first

    !--------------------------------------------------------------------------------------
    !   Triangle Second Order Element Type
    !--------------------------------------------------------------------------------------
    type, extends(abst_element) :: type_triangle_second
    contains
        procedure, pass(self) :: get_id        => get_id_triangle_second !&
        procedure, pass(self) :: get_type      => get_type_triangle_second !&
        procedure, pass(self) :: get_num_nodes => get_num_nodes_triangle_second !&
        procedure, pass(self) :: get_group     => get_group_triangle_second !&
        procedure, pass(self) :: get_order     => get_order_triangle_second !&
        procedure, pass(self) :: get_dimension => get_dimension_triangle_second !&
        procedure, pass(self) :: get_num_gauss => get_num_gauss_triangle_second !&
        !----------------------------------------------------------------------------------
        procedure, pass(self) :: psi           => psi_triangle_second !&
        procedure, pass(self) :: dpsi_dxi      => dpsi_dxi_triangle_second !&
        procedure, pass(self) :: dpsi_deta     => dpsi_deta_triangle_second !&
        procedure, pass(self) :: jacobian      => jacobian_triangle_second !&
        procedure, pass(self) :: jacobian_det  => jacobian_det_triangle_second !&
        procedure, pass(self) :: is_inside     => is_in_triangle_second !&
    end type type_triangle_second

    !--------------------------------------------------------------------------------------
    !   Square Second Order Element Type
    !--------------------------------------------------------------------------------------
    type, extends(abst_element) :: type_square_second
    contains
        procedure, pass(self) :: get_id        => get_id_square_second !&
        procedure, pass(self) :: get_type      => get_type_square_second !&
        procedure, pass(self) :: get_num_nodes => get_num_nodes_square_second !&
        procedure, pass(self) :: get_group     => get_group_square_second !&
        procedure, pass(self) :: get_order     => get_order_square_second !&
        procedure, pass(self) :: get_dimension => get_dimension_square_second !&
        procedure, pass(self) :: get_num_gauss => get_num_gauss_square_second !&
        !----------------------------------------------------------------------------------
        procedure, pass(self) :: psi           => psi_square_second !&
        procedure, pass(self) :: dpsi_dxi      => dpsi_dxi_square_second !&
        procedure, pass(self) :: dpsi_deta     => dpsi_deta_square_second !&
        procedure, pass(self) :: jacobian      => jacobian_square_second !&
        procedure, pass(self) :: jacobian_det  => jacobian_det_square_second !&
        procedure, pass(self) :: is_inside     => is_in_square_second !&
    end type type_square_second

    !
    !----- 抽象インターフェース定義 -----
    !
    abstract interface
        function abst_get_id(self) result(id)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32) :: id
        end function abst_get_id

        function abst_get_type(self) result(type)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32) :: type
        end function abst_get_type

        function abst_get_num_nodes(self) result(num_nodes)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32) :: num_nodes

        end function abst_get_num_nodes

        function abst_get_group(self) result(group)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32) :: group
        end function abst_get_group

        function abst_get_order(self) result(order)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32) :: order
        end function abst_get_order

        function abst_get_dimension(self) result(dimension)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32) :: dimension
        end function abst_get_dimension

        function abst_get_num_gauss(self) result(num_gauss)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32) :: num_gauss

        end function abst_get_num_gauss

        function abst_get_connectivity(self, index) result(connectivity)
            import :: abst_element, int32
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32), intent(in) :: index
            integer(int32) :: connectivity

        end function abst_get_connectivity

        function abst_psi(self, i, xi, eta) result(psi)
            import :: abst_element, int32, real64
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: psi
        end function abst_psi

        function abst_dpsi_dxi(self, i, xi, eta) result(dpsi)
            import :: abst_element, int32, real64
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi
        end function abst_dpsi_dxi

        function abst_dpsi_deta(self, i, xi, eta) result(dpsi)
            import :: abst_element, int32, real64
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi
        end function abst_dpsi_deta

        function abst_jacobian(self, i, j, xi, eta) result(Jval)
            import :: abst_element, int32, real64
            implicit none
            class(abst_element), intent(in) :: self
            integer(int32), intent(in) :: i, j
            real(real64), intent(in) :: xi, eta
            real(real64) :: Jval
        end function abst_jacobian

        function abst_jacobian_det(self, xi, eta) result(J_Det)
            import :: abst_element, int32, real64
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64) :: J_Det
        end function abst_jacobian_det

        subroutine abst_is_inside(self, px, py, pxi, peta, is_in)
            import abst_element, real64
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: px, py
            real(real64), intent(inout) :: pxi, peta
            logical, intent(inout) :: is_in
        end subroutine abst_is_inside

        function abst_interpolate(self, xi, eta, value) result(interpolated_value)
            import :: abst_element, real64
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value
        end function abst_interpolate

    end interface

    !--------------------------------------------------------------------------------------
    !   三角形一次要素型 procedures interface
    !--------------------------------------------------------------------------------------
    interface
        module function construct_triangle_first(id, global_coordinate, cell_info, integration) result(element)
            implicit none
            integer(int32), intent(in) :: id
            type(type_dp_3d), pointer, intent(in) :: global_coordinate
            type(type_vtk_cell), intent(in) :: cell_info
            type(type_geometry_settings), intent(in) :: integration
            class(abst_element), allocatable :: element

        end function construct_triangle_first

        module function get_id_triangle_first(self) result(id)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32) :: id

        end function get_id_triangle_first

        module function get_type_triangle_first(self) result(type)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32) :: type

        end function get_type_triangle_first

        module function get_num_nodes_triangle_first(self) result(num_nodes)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32) :: num_nodes

        end function get_num_nodes_triangle_first

        module function get_group_triangle_first(self) result(group)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32) :: group

        end function get_group_triangle_first

        module function get_order_triangle_first(self) result(order)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32) :: order

        end function get_order_triangle_first

        module function get_dimension_triangle_first(self) result(dimension)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32) :: dimension

        end function get_dimension_triangle_first

        module function get_num_gauss_triangle_first(self) result(num_gauss)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32) :: num_gauss

        end function get_num_gauss_triangle_first

        module function psi_triangle_first(self, i, xi, eta) result(N)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: N

        end function psi_triangle_first

        module function dpsi_dxi_triangle_first(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_dxi_triangle_first

        module function dpsi_deta_triangle_first(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_deta_triangle_first

        module function jacobian_triangle_first(self, i, j, xi, eta) result(Jval)
            implicit none
            class(type_triangle_first), intent(in) :: self
            integer(int32), intent(in) :: i, j
            real(real64), intent(in) :: xi, eta
            real(real64) :: Jval

        end function jacobian_triangle_first

        module function jacobian_det_triangle_first(self, xi, eta) result(J_Det)
            implicit none
            class(type_triangle_first), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64) :: J_Det

        end function jacobian_det_triangle_first

        module subroutine is_in_triangle_first(self, px, py, pxi, peta, is_in)
            implicit none
            class(type_triangle_first), intent(in) :: self
            real(real64), intent(in) :: px, py
            real(real64), intent(inout) :: pxi, peta
            logical, intent(inout) :: is_in

        end subroutine is_in_triangle_first

        module function interpolate_triangle_first(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_triangle_first

        module function interpolate_reordered_triangle_first(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_reordered_triangle_first
    end interface

    !--------------------------------------------------------------------------------------
    !   四角形一次要素型 procedures interface
    !--------------------------------------------------------------------------------------
    interface
        module function construct_square_first(id, global_coordinate, cell_info, integration) result(element)
            implicit none
            integer(int32), intent(in) :: id
            type(type_dp_3d), pointer, intent(in) :: global_coordinate
            type(type_vtk_cell), intent(in) :: cell_info
            type(type_geometry_settings), intent(in) :: integration
            class(abst_element), allocatable :: element

        end function construct_square_first

        module function get_id_square_first(self) result(id)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32) :: id

        end function get_id_square_first

        module function get_type_square_first(self) result(type)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32) :: type

        end function get_type_square_first

        module function get_num_nodes_square_first(self) result(num_nodes)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32) :: num_nodes

        end function get_num_nodes_square_first

        module function get_group_square_first(self) result(group)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32) :: group

        end function get_group_square_first

        module function get_order_square_first(self) result(order)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32) :: order

        end function get_order_square_first

        module function get_dimension_square_first(self) result(dimension)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32) :: dimension

        end function get_dimension_square_first

        module function get_num_gauss_square_first(self) result(num_gauss)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32) :: num_gauss

        end function get_num_gauss_square_first

        module function psi_square_first(self, i, xi, eta) result(psi)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: psi

        end function psi_square_first

        module function dpsi_dxi_square_first(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_dxi_square_first

        module function dpsi_deta_square_first(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_deta_square_first

        module function jacobian_square_first(self, i, j, xi, eta) result(Jval)
            implicit none
            class(type_square_first), intent(in) :: self
            integer(int32), intent(in) :: i, j
            real(real64), intent(in) :: xi, eta
            real(real64) :: Jval

        end function jacobian_square_first

        module function jacobian_det_square_first(self, xi, eta) result(J_Det)
            implicit none
            class(type_square_first), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64) :: J_Det

        end function jacobian_det_square_first

        module subroutine is_in_square_first(self, px, py, pxi, peta, is_in)
            implicit none
            class(type_square_first), intent(in) :: self
            real(real64), intent(in) :: px, py
            real(real64), intent(inout) :: pxi, peta
            logical, intent(inout) :: is_in

        end subroutine is_in_square_first

        module function interpolate_square_first(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_square_first

        module function interpolate_reordered_square_first(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_reordered_square_first
    end interface

    !--------------------------------------------------------------------------------------
    !   三角形二次要素型 procedures interface
    !--------------------------------------------------------------------------------------
    interface
        module function construct_triangle_second(id, global_coordinate, cell_info, integration) result(element)
            implicit none
            integer(int32), intent(in) :: id
            type(type_dp_3d), pointer, intent(in) :: global_coordinate
            type(type_vtk_cell), intent(in) :: cell_info
            type(type_geometry_settings), intent(in) :: integration
            class(abst_element), allocatable :: element

        end function construct_triangle_second

        module function get_id_triangle_second(self) result(id)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32) :: id

        end function get_id_triangle_second

        module function get_type_triangle_second(self) result(type)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32) :: type

        end function get_type_triangle_second

        module function get_num_nodes_triangle_second(self) result(num_nodes)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32) :: num_nodes

        end function get_num_nodes_triangle_second

        module function get_group_triangle_second(self) result(group)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32) :: group

        end function get_group_triangle_second

        module function get_order_triangle_second(self) result(order)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32) :: order

        end function get_order_triangle_second

        module function get_dimension_triangle_second(self) result(dimension)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32) :: dimension

        end function get_dimension_triangle_second

        module function get_num_gauss_triangle_second(self) result(num_gauss)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32) :: num_gauss

        end function get_num_gauss_triangle_second

        module function psi_triangle_second(self, i, xi, eta) result(N)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: N

        end function psi_triangle_second

        module function dpsi_dxi_triangle_second(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_dxi_triangle_second

        module function dpsi_deta_triangle_second(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_deta_triangle_second

        module function jacobian_triangle_second(self, i, j, xi, eta) result(Jval)
            implicit none
            class(type_triangle_second), intent(in) :: self
            integer(int32), intent(in) :: i, j
            real(real64), intent(in) :: xi, eta
            real(real64) :: Jval

        end function jacobian_triangle_second

        module function jacobian_det_triangle_second(self, xi, eta) result(J_Det)
            implicit none
            class(type_triangle_second), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64) :: J_Det

        end function jacobian_det_triangle_second

        module subroutine is_in_triangle_second(self, px, py, pxi, peta, is_in)
            implicit none
            class(type_triangle_second), intent(in) :: self
            real(real64), intent(in) :: px, py
            real(real64), intent(inout) :: pxi, peta
            logical, intent(inout) :: is_in

        end subroutine is_in_triangle_second

        module function interpolate_triangle_second(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_triangle_second

        module function interpolate_reordered_triangle_second(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_reordered_triangle_second
    end interface

    !--------------------------------------------------------------------------------------
    !   四角形二次要素型 procedures interface
    !--------------------------------------------------------------------------------------
    interface
        module function construct_square_second(id, global_coordinate, cell_info, integration) result(element)
            implicit none
            integer(int32), intent(in) :: id
            type(type_dp_3d), pointer, intent(in) :: global_coordinate
            type(type_vtk_cell), intent(in) :: cell_info
            type(type_geometry_settings), intent(in) :: integration
            class(abst_element), allocatable :: element

        end function construct_square_second

        module function get_id_square_second(self) result(id)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32) :: id

        end function get_id_square_second

        module function get_type_square_second(self) result(type)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32) :: type

        end function get_type_square_second

        module function get_num_nodes_square_second(self) result(num_nodes)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32) :: num_nodes

        end function get_num_nodes_square_second

        module function get_group_square_second(self) result(group)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32) :: group

        end function get_group_square_second

        module function get_order_square_second(self) result(order)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32) :: order

        end function get_order_square_second

        module function get_dimension_square_second(self) result(dimension)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32) :: dimension

        end function get_dimension_square_second

        module function get_num_gauss_square_second(self) result(num_gauss)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32) :: num_gauss

        end function get_num_gauss_square_second
        module function psi_square_second(self, i, xi, eta) result(psi)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: psi

        end function psi_square_second

        module function dpsi_dxi_square_second(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_dxi_square_second

        module function dpsi_deta_square_second(self, i, xi, eta) result(dpsi)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32), intent(in) :: i
            real(real64), intent(in) :: xi, eta
            real(real64) :: dpsi

        end function dpsi_deta_square_second

        module function jacobian_square_second(self, i, j, xi, eta) result(Jval)
            implicit none
            class(type_square_second), intent(in) :: self
            integer(int32), intent(in) :: i, j
            real(real64), intent(in) :: xi, eta
            real(real64) :: Jval

        end function jacobian_square_second

        module function jacobian_det_square_second(self, xi, eta) result(J_Det)
            implicit none
            class(type_square_second), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64) :: J_Det

        end function jacobian_det_square_second

        module subroutine is_in_square_second(self, px, py, pxi, peta, is_in)
            implicit none
            class(type_square_second), intent(in) :: self
            real(real64), intent(in) :: px, py
            real(real64), intent(inout) :: pxi, peta
            logical, intent(inout) :: is_in

        end subroutine is_in_square_second

        module function interpolate_square_second(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_square_second

        module function interpolate_reordered_square_second(self, xi, eta, value) result(interpolated_value)
            implicit none
            class(abst_element), intent(in) :: self
            real(real64), intent(in) :: xi, eta
            real(real64), intent(in) :: value(:)
            real(real64) :: interpolated_value

        end function interpolate_reordered_square_second
    end interface

    interface type_triangle_first
        module procedure :: construct_triangle_first
    end interface

    interface type_square_first
        module procedure :: construct_square_first
    end interface

    interface type_triangle_second
        module procedure :: construct_triangle_second
    end interface

    interface type_square_second
        module procedure :: construct_square_second
    end interface

contains
    function interpolate(self, xi, eta, value) result(interpolated_value)
        implicit none
        class(abst_element), intent(in) :: self
        real(real64), intent(in) :: xi, eta
        real(real64), intent(in) :: value(:)
        real(real64) :: interpolated_value
        integer(int32) :: i

        interpolated_value = 0.0d0
        do i = 1, self%num_nodes
            interpolated_value = interpolated_value + self%psi(i, xi, eta) * value(self%connectivity(i))
        end do
    end function interpolate

    function interpolate_reordered(self, xi, eta, value) result(interpolated_value)
        implicit none
        class(abst_element), intent(in) :: self
        real(real64), intent(in) :: xi, eta
        real(real64), intent(in) :: value(:)
        real(real64) :: interpolated_value
        integer(int32) :: i

        interpolated_value = 0.0d0
        do i = 1, self%num_nodes
            interpolated_value = interpolated_value + self%psi(i, xi, eta) * value(self%connectivity_reordered(i))
        end do
    end function interpolate_reordered

    function get_connectivity(self, index) result(connectivity)
        implicit none
        class(abst_element), intent(in) :: self
        integer(int32), intent(in) :: index
        integer(int32) :: connectivity

        if (index < 1 .or. index > self%num_nodes) then
            error stop 'Index out of bounds in get_connectivity'
        end if

        connectivity = self%connectivity(index)
    end function get_connectivity

    function get_connectivity_reordered(self, index) result(connectivity)
        implicit none
        class(abst_element), intent(in) :: self
        integer(int32), intent(in) :: index
        integer(int32) :: connectivity

        if (index < 1 .or. index > self%num_nodes) then
            error stop 'Index out of bounds in get_connectivity_reordered'
        end if

        connectivity = self%connectivity_reordered(index)
    end function get_connectivity_reordered

end module domain_element