output_overall_vtu_fields Subroutine

subroutine output_overall_vtu_fields(self, file_counts, domain, porosity, temperature, si, pressure, water_flux)

Arguments

Type IntentOptional Attributes Name
class(type_output_overall), intent(inout) :: self
integer(kind=int32), intent(in) :: file_counts
type(type_domain), intent(in) :: domain
real(kind=real64), intent(in), optional :: porosity(:)
real(kind=real64), intent(in), optional :: temperature(:)
real(kind=real64), intent(in), optional :: si(:)
real(kind=real64), intent(in), optional :: pressure(:)
type(type_dp_3d), intent(in), optional :: water_flux

Calls

proc~~output_overall_vtu_fields~~CallsGraph proc~output_overall_vtu_fields output_overall_vtu_fields finalize finalize proc~output_overall_vtu_fields->finalize initialize initialize proc~output_overall_vtu_fields->initialize interface~allocate_array allocate_array proc~output_overall_vtu_fields->interface~allocate_array interface~deallocate_array deallocate_array proc~output_overall_vtu_fields->interface~deallocate_array none~to_original_value type_reordering%to_original_value proc~output_overall_vtu_fields->none~to_original_value write_connectivity write_connectivity proc~output_overall_vtu_fields->write_connectivity write_dataarray write_dataarray proc~output_overall_vtu_fields->write_dataarray write_geo write_geo proc~output_overall_vtu_fields->write_geo write_piece write_piece proc~output_overall_vtu_fields->write_piece 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~deallocate_rank1_int32 deallocate_rank1_int32 interface~deallocate_array->proc~deallocate_rank1_int32 proc~deallocate_rank1_int64 deallocate_rank1_int64 interface~deallocate_array->proc~deallocate_rank1_int64 proc~deallocate_rank1_int8 deallocate_rank1_int8 interface~deallocate_array->proc~deallocate_rank1_int8 proc~deallocate_rank1_logical1 deallocate_rank1_logical1 interface~deallocate_array->proc~deallocate_rank1_logical1 proc~deallocate_rank1_logical4 deallocate_rank1_logical4 interface~deallocate_array->proc~deallocate_rank1_logical4 proc~deallocate_rank1_logical8 deallocate_rank1_logical8 interface~deallocate_array->proc~deallocate_rank1_logical8 proc~deallocate_rank1_real128 deallocate_rank1_real128 interface~deallocate_array->proc~deallocate_rank1_real128 proc~deallocate_rank1_real32 deallocate_rank1_real32 interface~deallocate_array->proc~deallocate_rank1_real32 proc~deallocate_rank1_real64 deallocate_rank1_real64 interface~deallocate_array->proc~deallocate_rank1_real64 proc~deallocate_rank2_int32 deallocate_rank2_int32 interface~deallocate_array->proc~deallocate_rank2_int32 proc~deallocate_rank2_int64 deallocate_rank2_int64 interface~deallocate_array->proc~deallocate_rank2_int64 proc~deallocate_rank2_int8 deallocate_rank2_int8 interface~deallocate_array->proc~deallocate_rank2_int8 proc~deallocate_rank2_logical1 deallocate_rank2_logical1 interface~deallocate_array->proc~deallocate_rank2_logical1 proc~deallocate_rank2_logical4 deallocate_rank2_logical4 interface~deallocate_array->proc~deallocate_rank2_logical4 proc~deallocate_rank2_logical8 deallocate_rank2_logical8 interface~deallocate_array->proc~deallocate_rank2_logical8 proc~deallocate_rank2_real128 deallocate_rank2_real128 interface~deallocate_array->proc~deallocate_rank2_real128 proc~deallocate_rank2_real32 deallocate_rank2_real32 interface~deallocate_array->proc~deallocate_rank2_real32 proc~deallocate_rank2_real64 deallocate_rank2_real64 interface~deallocate_array->proc~deallocate_rank2_real64 interface~to_original_values_int32 type_reordering%to_original_values_int32 none~to_original_value->interface~to_original_values_int32 interface~to_original_values_real64 type_reordering%to_original_values_real64 none~to_original_value->interface~to_original_values_real64 proc~to_original_values_int32 to_original_values_int32 interface~to_original_values_int32->proc~to_original_values_int32 proc~to_original_values_real64 to_original_values_real64 interface~to_original_values_real64->proc~to_original_values_real64 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 proc~deallocate_rank1_int32->proc~error_message proc~deallocate_rank1_int64->proc~error_message proc~deallocate_rank1_int8->proc~error_message proc~deallocate_rank1_logical1->proc~error_message proc~deallocate_rank1_logical4->proc~error_message proc~deallocate_rank1_logical8->proc~error_message proc~deallocate_rank1_real128->proc~error_message proc~deallocate_rank1_real32->proc~error_message proc~deallocate_rank1_real64->proc~error_message proc~deallocate_rank2_int32->proc~error_message proc~deallocate_rank2_int64->proc~error_message proc~deallocate_rank2_int8->proc~error_message proc~deallocate_rank2_logical1->proc~error_message proc~deallocate_rank2_logical4->proc~error_message proc~deallocate_rank2_logical8->proc~error_message proc~deallocate_rank2_real128->proc~error_message proc~deallocate_rank2_real32->proc~error_message proc~deallocate_rank2_real64->proc~error_message log_error log_error proc~error_message->log_error

Source Code

    subroutine output_overall_vtu_fields(self, file_counts, domain, porosity, temperature, si, pressure, water_flux)
        implicit none
        class(type_output_overall), intent(inout) :: self
        integer(int32), intent(in) :: file_counts
        type(type_domain), intent(in) :: domain
        real(real64), intent(in), optional :: porosity(:)
        real(real64), intent(in), optional :: temperature(:)
        real(real64), intent(in), optional :: si(:)
        real(real64), intent(in), optional :: pressure(:)
        type(type_dp_3d), intent(in), optional :: water_flux

        type(vtk_file) :: vtu
        integer(int32) :: status
        integer(int32) :: unit_num

        real(real64), allocatable :: original(:), original_vector(:, :)

        integer(int32) :: nsize, i

        character(256) :: output_name

        ! Initialize VTK file
        write (output_name, self%format_output) trim(self%dir_output_field), "Out_", file_counts, self%file_extension

        status = vtu%initialize(format='ascii', filename=trim(output_name), mesh_topology='UnstructuredGrid')

        ! Write data
        status = vtu%xml_writer%write_piece(np=self%vtk%num_points, &
                                            nc=self%vtk%num_cells)
        status = vtu%xml_writer%write_geo(np=self%vtk%num_points, &
                                          nc=self%vtk%num_cells, &
                                          x=self%vtk%coordinate%x, &
                                          y=self%vtk%coordinate%y, &
                                          z=self%vtk%coordinate%z)
        status = vtu%xml_writer%write_connectivity(nc=self%vtk%num_cells, &
                                                   connectivity=self%VTK%connectivities, &
                                                   offset=self%vtk%offsets, &
                                                   cell_type=self%vtk%cell_types)

        ! --- データセクション ---

        do i = 1, size(self%variable_names)
            if (i == 1) status = vtu%xml_writer%write_dataarray(location='node', action='open')
            select case (self%variable_names(i))
            case ("temperature")
                if (present(temperature)) then
                    call allocate_array(original, self%vtk%num_points)
                    call domain%reordering%to_original_value(temperature, original)
                    status = vtu%xml_writer%write_dataarray(data_name='Temperature', &
                                                            x=original)
                    call deallocate_array(original)
                end if
            case ("ice_saturation")
                if (present(si)) then
                    call allocate_array(original, self%vtk%num_points)
                    call domain%reordering%to_original_value(si, original)
                    status = vtu%xml_writer%write_dataarray(data_name='Si', &
                                                            x=original)
                    call deallocate_array(original)
                end if
            case ("thermal_conductivity")
                print *, "Warning: 'thermal_conductivity' is not implemented in VTK output."
            case ("volumetric_heat_capacity")
                print *, "Warning: 'volumetric_heat_capacity' is not implemented in VTK output."
            case ("pressure")
                if (present(pressure)) then
                    call allocate_array(original, self%vtk%num_points)
                    call domain%reordering%to_original_value(pressure, original)
                    status = vtu%xml_writer%write_dataarray(data_name='Pressure', &
                                                            x=original)
                    call deallocate_array(original)
                end if
            case ("water_flux")
                if (present(water_flux)) then
                    call allocate_array(original_vector, 3_int32, self%vtk%num_points)
                    call domain%reordering%to_original_value(water_flux%x, original_vector(:, 1))
                    call domain%reordering%to_original_value(water_flux%y, original_vector(:, 2))
                    call domain%reordering%to_original_value(water_flux%z, original_vector(:, 3))
                    status = vtu%xml_writer%write_dataarray(data_name='waterFlux', &
                                                            x=original_vector(:, 1), &
                                                            y=original_vector(:, 2), &
                                                            z=original_vector(:, 3))
                    call deallocate_array(original_vector)
                end if
            case ("hydraulic_conductivity")
                print *, "Warning: 'hydraulic_conductivity' is not implemented in VTK output."
            end select

        end do
        status = vtu%xml_writer%write_dataarray(location='node', action='close')
        status = vtu%xml_writer%write_piece()

        ! Finalize VTK file
        status = vtu%finalize()

    end subroutine output_overall_vtu_fields