initialize_abst_mesh Subroutine

private subroutine initialize_abst_mesh(self, id, type, group, dimension, order, num_nodes, connectivity, num_gauss, weight, gauss, global_coordinate)

Type Bound

abst_mesh

Arguments

Type IntentOptional Attributes Name
class(abst_mesh), intent(inout) :: self
integer(kind=int32), intent(in) :: id
integer(kind=int32), intent(in) :: type
integer(kind=int32), intent(in) :: group
integer(kind=int32), intent(in) :: dimension
integer(kind=int32), intent(in) :: order
integer(kind=int32), intent(in) :: num_nodes
integer(kind=int32), intent(in) :: connectivity(:)
integer(kind=int32), intent(in) :: num_gauss
real(kind=real64), intent(in) :: weight(:)
real(kind=real64), intent(in) :: gauss(:,:)
type(type_dp_3d), intent(in), pointer :: global_coordinate

Calls

proc~~initialize_abst_mesh~~CallsGraph proc~initialize_abst_mesh abst_mesh%initialize_abst_mesh interface~allocate_array allocate_array proc~initialize_abst_mesh->interface~allocate_array none~set~3 type_dp_vector_3d%set proc~initialize_abst_mesh->none~set~3 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~set_dp_vector_3d type_dp_vector_3d%set_dp_vector_3d none~set~3->proc~set_dp_vector_3d proc~set_dp_vector_3d_array type_dp_vector_3d%set_dp_vector_3d_array none~set~3->proc~set_dp_vector_3d_array 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 log_error log_error proc~error_message->log_error

Called by

proc~~initialize_abst_mesh~~CalledByGraph proc~initialize_abst_mesh abst_mesh%initialize_abst_mesh proc~construct_side_first construct_side_first proc~construct_side_first->proc~initialize_abst_mesh proc~construct_side_second construct_side_second proc~construct_side_second->proc~initialize_abst_mesh proc~construct_square_first construct_square_first proc~construct_square_first->proc~initialize_abst_mesh proc~construct_square_second construct_square_second proc~construct_square_second->proc~initialize_abst_mesh proc~construct_triangle_first construct_triangle_first proc~construct_triangle_first->proc~initialize_abst_mesh proc~construct_triangle_second construct_triangle_second proc~construct_triangle_second->proc~initialize_abst_mesh interface~construct_side_first construct_side_first interface~construct_side_first->proc~construct_side_first interface~construct_side_second construct_side_second interface~construct_side_second->proc~construct_side_second interface~construct_square_first construct_square_first interface~construct_square_first->proc~construct_square_first interface~construct_square_second construct_square_second interface~construct_square_second->proc~construct_square_second interface~construct_triangle_first construct_triangle_first interface~construct_triangle_first->proc~construct_triangle_first interface~construct_triangle_second construct_triangle_second interface~construct_triangle_second->proc~construct_triangle_second interface~type_side_first type_side_first interface~type_side_first->interface~construct_side_first interface~type_side_second type_side_second interface~type_side_second->interface~construct_side_second interface~type_square_first type_square_first interface~type_square_first->interface~construct_square_first interface~type_square_second type_square_second interface~type_square_second->interface~construct_square_second interface~type_triangle_first type_triangle_first interface~type_triangle_first->interface~construct_triangle_first interface~type_triangle_second type_triangle_second interface~type_triangle_second->interface~construct_triangle_second

Source Code

    subroutine initialize_abst_mesh(self, id, type, group, dimension, order, num_nodes, connectivity, &
                                    num_gauss, weight, gauss, global_coordinate)
        implicit none
        class(abst_mesh), intent(inout) :: self
        integer(int32), intent(in) :: id
        integer(int32), intent(in) :: type
        integer(int32), intent(in) :: group
        integer(int32), intent(in) :: dimension
        integer(int32), intent(in) :: order
        integer(int32), intent(in) :: num_nodes
        integer(int32), intent(in) :: connectivity(:)
        integer(int32), intent(in) :: num_gauss
        real(real64), intent(in) :: weight(:)
        real(real64), intent(in) :: gauss(:, :)
        type(type_dp_3d), pointer, intent(in) :: global_coordinate

        integer(int32) :: i

        self%id = id
        self%type = type
        self%group = group
        self%dimension = dimension
        self%order = order
        self%num_nodes = num_nodes
        self%num_gauss = num_gauss

        call allocate_array(self%connectivity, self%num_nodes)
        do i = 1, self%num_nodes
            self%connectivity(i) = connectivity(i)
        end do

        call allocate_array(self%weight, self%num_gauss)
        do i = 1, self%num_gauss
            self%weight(i) = weight(i)
        end do

        allocate (self%gauss(self%num_gauss))
        do i = 1, self%num_gauss
            call self%gauss(i)%set(gauss(1, i), gauss(2, i), gauss(3, i))
        end do

        allocate (self%x(self%num_nodes))
        allocate (self%y(self%num_nodes))
        allocate (self%z(self%num_nodes))
        do i = 1, self%num_nodes
            nullify (self%x(i)%val)
            nullify (self%y(i)%val)
            nullify (self%z(i)%val)
            self%x(i)%val => global_coordinate%x(self%connectivity(i))
            self%y(i)%val => global_coordinate%y(self%connectivity(i))
            self%z(i)%val => global_coordinate%z(self%connectivity(i))
        end do

    end subroutine initialize_abst_mesh