get_absolute_resnorm Subroutine

private subroutine get_absolute_resnorm(residue, F, G, H, control, dims)

Get absolute residual for current process

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), 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(controltype), intent(in) :: control

Control parameters: number of variables

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

extent of the 3D domain


Called by

proc~~get_absolute_resnorm~~CalledByGraph proc~get_absolute_resnorm get_absolute_resnorm proc~find_resnorm find_resnorm proc~find_resnorm->proc~get_absolute_resnorm proc~iterate_one_more_time_step iterate_one_more_time_step proc~iterate_one_more_time_step->proc~find_resnorm program~main main program~main->proc~iterate_one_more_time_step

Contents

Source Code


Source Code

    subroutine get_absolute_resnorm(residue, F,G,H, control, dims)
      !< Get absolute residual for current process
      implicit none
      type(controltype), intent(in) :: control
      !< Control parameters: number of variables
      type(extent), intent(in) :: dims
      !< extent of the 3D domain
      real(wp), dimension(:, :, :, :), intent(in)  :: 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
      do i=1,control%n_var
        Res_abs(i) =(sum(Residue(:,:,:,i)**2)/Res_scale(i)**2)
      end do
      merror = (                                     &
               sum(F(  1,1:dims%jmx-1,1:dims%kmx-1,1)) &
              -sum(F(dims%imx,1:dims%jmx-1,1:dims%kmx-1,1)) &
              +sum(G(1:dims%imx-1,  1,1:dims%kmx-1,1)) &
              -sum(G(1:dims%imx-1,dims%jmx,1:dims%kmx-1,1)) &
              +sum(H(1:dims%imx-1,1:dims%jmx-1,  1,1)) &
              -sum(H(1:dims%imx-1,1:dims%jmx-1,dims%kmx,1)) &
              )
      Res_abs(0) = (merror/Res_scale(0))
    end subroutine get_absolute_resnorm