compute_residue Subroutine

public subroutine compute_residue(residue, F, G, H, dims)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(out), dimension(:, :, :, :):: residue

Store residue at each cell-center

real(kind=wp), intent(in), dimension(:, :, :, :):: F

Store fluxes throught the I faces

real(kind=wp), intent(in), dimension(:, :, :, :):: G

Store fluxes throught the J faces

real(kind=wp), intent(in), dimension(:, :, :, :):: H

Store fluxes throught the K faces

type(extent), intent(in) :: dims

extent of the 3D domain


Calls

proc~~compute_residue~~CallsGraph proc~compute_residue compute_residue debugcall debugcall proc~compute_residue->debugcall

Called by

proc~~compute_residue~~CalledByGraph proc~compute_residue compute_residue proc~get_total_conservative_residue get_total_conservative_Residue proc~get_total_conservative_residue->proc~compute_residue proc~get_next_solution get_next_solution proc~get_next_solution->proc~get_total_conservative_residue proc~iterate_one_more_time_step iterate_one_more_time_step proc~iterate_one_more_time_step->proc~get_next_solution program~main main program~main->proc~iterate_one_more_time_step

Contents

Source Code


Source Code

        subroutine compute_residue(residue,F,G,H,dims)
            
            implicit none
            type(extent), intent(in) :: dims
            !< extent of the 3D domain
            real(wp), dimension(:, :, :, :), intent(out)  :: residue
            !< Store residue at each cell-center
            real(wp), dimension(:, :, :, :), intent(in) :: F
            !< Store fluxes throught the I faces
            real(wp), dimension(:, :, :, :), intent(in) :: G
            !< Store fluxes throught the J faces
            real(wp), dimension(:, :, :, :), intent(in) :: H
            !< Store fluxes throught the K faces
            
            integer :: i, j, k, l

            DebugCall('compute_residue')

            do l = 1, dims%n_var
             do k = 1, dims%kmx - 1
              do j = 1, dims%jmx - 1
               do i = 1, dims%imx - 1
               residue(i, j, k, l) = (F(i+1, j, k, l) - F(i, j, k, l)) &
                                   + (G(i, j+1, k, l) - G(i, j, k, l)) &
                                   + (H(i, j, k+1, l) - H(i, j, k, l))
               end do
              end do
             end do
            end do
        
        end subroutine compute_residue