normalize_face_normals Subroutine

private subroutine normalize_face_normals(Ifaces, Jfaces, Kfaces, bc)

Normalize the face normal vectors computed to get unit vectors

Arguments

Type IntentOptional AttributesName
type(facetype), intent(inout), dimension(-2:imx+3,-2:jmx+2,-2:kmx+2):: Ifaces

Store face quantites for I faces: normal and area

type(facetype), intent(inout), dimension(-2:imx+2,-2:jmx+3,-2:kmx+2):: Jfaces

Store face quantites for J faces: normal and area

type(facetype), intent(inout), dimension(-2:imx+2,-2:jmx+2,-2:kmx+3):: Kfaces

Store face quantites for K faces: normal and area

type(boundarytype), intent(in) :: bc

Called by

proc~~normalize_face_normals~~CalledByGraph proc~normalize_face_normals normalize_face_normals proc~compute_face_areas_and_normals compute_face_areas_and_normals proc~compute_face_areas_and_normals->proc~normalize_face_normals proc~setup_geometry setup_geometry proc~setup_geometry->proc~compute_face_areas_and_normals proc~setup_solver setup_solver proc~setup_solver->proc~setup_geometry proc~start_run start_run proc~start_run->proc~setup_solver program~main main program~main->proc~start_run

Contents


Source Code

        subroutine normalize_face_normals(Ifaces, Jfaces, Kfaces, bc)
            !< Normalize the face normal vectors computed to get unit
            !< vectors
            !-----------------------------------------------------------
            
            implicit none
            type(facetype), dimension(-2:imx+3,-2:jmx+2,-2:kmx+2), intent(inout) :: Ifaces
            !< Store face quantites for I faces: normal and area
            type(facetype), dimension(-2:imx+2,-2:jmx+3,-2:kmx+2), intent(inout) :: Jfaces
            !< Store face quantites for J faces: normal and area
            type(facetype), dimension(-2:imx+2,-2:jmx+2,-2:kmx+3), intent(inout) :: Kfaces
            !< Store face quantites for K faces: normal and area
            type(boundarytype), intent(in) :: bc
            integer :: i,j,k

            do k = -2,kmx+2
              do j = -2,jmx+2
                do i = -2,imx+3
                  if(Ifaces(i,j,k)%A/=0.) then
                    Ifaces(i,j,k)%nx = Ifaces(i,j,k)%nx/Ifaces(i,j,k)%A
                    Ifaces(i,j,k)%ny = Ifaces(i,j,k)%ny/Ifaces(i,j,k)%A
                    Ifaces(i,j,k)%nz = Ifaces(i,j,k)%nz/Ifaces(i,j,k)%A
                  end if
                end do
              end do
            end do

            do k = -2,kmx+2
              do j = -2,jmx+3
                do i = -2,imx+2
                  if(Jfaces(i,j,k)%A/=0.) then
                    Jfaces(i,j,k)%nx = Jfaces(i,j,k)%nx/Jfaces(i,j,k)%A
                    Jfaces(i,j,k)%ny = Jfaces(i,j,k)%ny/Jfaces(i,j,k)%A
                    Jfaces(i,j,k)%nz = Jfaces(i,j,k)%nz/Jfaces(i,j,k)%A
                  end if
                end do
              end do
            end do

            do k = -2,kmx+3
              do j = -2,jmx+2
                do i = -2,imx+2
                  if(Kfaces(i,j,k)%A/=0.) then
                    Kfaces(i,j,k)%nx = Kfaces(i,j,k)%nx/Kfaces(i,j,k)%A
                    Kfaces(i,j,k)%ny = Kfaces(i,j,k)%ny/Kfaces(i,j,k)%A
                    Kfaces(i,j,k)%nz = Kfaces(i,j,k)%nz/Kfaces(i,j,k)%A
                  end if
                end do
              end do
            end do

            ! pole boundary condition
            if(bc%imin_id==-7) then
              Ifaces( 1,:,:)%nx=Ifaces(2,:,:)%nx
              Ifaces( 0,:,:)%nx=Ifaces(2,:,:)%nx
              Ifaces(-1,:,:)%nx=Ifaces(2,:,:)%nx
              Ifaces(-2,:,:)%nx=Ifaces(2,:,:)%nx
              Ifaces( 1,:,:)%ny=Ifaces(2,:,:)%ny
              Ifaces( 0,:,:)%ny=Ifaces(2,:,:)%ny
              Ifaces(-1,:,:)%ny=Ifaces(2,:,:)%ny
              Ifaces(-2,:,:)%ny=Ifaces(2,:,:)%ny
              Ifaces( 1,:,:)%nz=Ifaces(2,:,:)%nz
              Ifaces( 0,:,:)%nz=Ifaces(2,:,:)%nz
              Ifaces(-1,:,:)%nz=Ifaces(2,:,:)%nz
              Ifaces(-2,:,:)%nz=Ifaces(2,:,:)%nz
            end if

            if(bc%imax_id==-7) then
              Ifaces(imx+0,:,:)%nx=Ifaces(imx-1,:,:)%nx
              Ifaces(imx+1,:,:)%nx=Ifaces(imx-1,:,:)%nx
              Ifaces(imx+2,:,:)%nx=Ifaces(imx-1,:,:)%nx
              Ifaces(imx+3,:,:)%nx=Ifaces(imx-1,:,:)%nx
              Ifaces(imx+0,:,:)%ny=Ifaces(imx-1,:,:)%ny
              Ifaces(imx+1,:,:)%ny=Ifaces(imx-1,:,:)%ny
              Ifaces(imx+2,:,:)%ny=Ifaces(imx-1,:,:)%ny
              Ifaces(imx+3,:,:)%ny=Ifaces(imx-1,:,:)%ny
              Ifaces(imx+0,:,:)%nz=Ifaces(imx-1,:,:)%nz
              Ifaces(imx+1,:,:)%nz=Ifaces(imx-1,:,:)%nz
              Ifaces(imx+2,:,:)%nz=Ifaces(imx-1,:,:)%nz
              Ifaces(imx+3,:,:)%nz=Ifaces(imx-1,:,:)%nz
            end if

            if(bc%jmin_id==-7) then
              Jfaces(:, 1,:)%nx=Jfaces(:,2,:)%nx
              Jfaces(:, 0,:)%nx=Jfaces(:,2,:)%nx
              Jfaces(:,-1,:)%nx=Jfaces(:,2,:)%nx
              Jfaces(:,-2,:)%nx=Jfaces(:,2,:)%nx
              Jfaces(:, 1,:)%ny=Jfaces(:,2,:)%ny
              Jfaces(:, 0,:)%ny=Jfaces(:,2,:)%ny
              Jfaces(:,-1,:)%ny=Jfaces(:,2,:)%ny
              Jfaces(:,-2,:)%ny=Jfaces(:,2,:)%ny
              Jfaces(:, 1,:)%nz=Jfaces(:,2,:)%nz
              Jfaces(:, 0,:)%nz=Jfaces(:,2,:)%nz
              Jfaces(:,-1,:)%nz=Jfaces(:,2,:)%nz
              Jfaces(:,-2,:)%nz=Jfaces(:,2,:)%nz
            end if

            if(bc%jmax_id==-7) then
              Jfaces(:,jmx+0,:)%nx=Jfaces(:,jmx-1,:)%nx
              Jfaces(:,jmx+1,:)%nx=Jfaces(:,jmx-1,:)%nx
              Jfaces(:,jmx+2,:)%nx=Jfaces(:,jmx-1,:)%nx
              Jfaces(:,jmx+3,:)%nx=Jfaces(:,jmx-1,:)%nx
              Jfaces(:,jmx+0,:)%ny=Jfaces(:,jmx-1,:)%ny
              Jfaces(:,jmx+1,:)%ny=Jfaces(:,jmx-1,:)%ny
              Jfaces(:,jmx+2,:)%ny=Jfaces(:,jmx-1,:)%ny
              Jfaces(:,jmx+3,:)%ny=Jfaces(:,jmx-1,:)%ny
              Jfaces(:,jmx+0,:)%nz=Jfaces(:,jmx-1,:)%nz
              Jfaces(:,jmx+1,:)%nz=Jfaces(:,jmx-1,:)%nz
              Jfaces(:,jmx+2,:)%nz=Jfaces(:,jmx-1,:)%nz
              Jfaces(:,jmx+3,:)%nz=Jfaces(:,jmx-1,:)%nz
            end if

            if(bc%kmin_id==-7) then
              Kfaces(:,:, 1)%nx=Kfaces(:,:,2)%nx
              Kfaces(:,:, 0)%nx=Kfaces(:,:,2)%nx
              Kfaces(:,:,-1)%nx=Kfaces(:,:,2)%nx
              Kfaces(:,:,-2)%nx=Kfaces(:,:,2)%nx
              Kfaces(:,:, 1)%ny=Kfaces(:,:,2)%ny
              Kfaces(:,:, 0)%ny=Kfaces(:,:,2)%ny
              Kfaces(:,:,-1)%ny=Kfaces(:,:,2)%ny
              Kfaces(:,:,-2)%ny=Kfaces(:,:,2)%ny
              Kfaces(:,:, 1)%nz=Kfaces(:,:,2)%nz
              Kfaces(:,:, 0)%nz=Kfaces(:,:,2)%nz
              Kfaces(:,:,-1)%nz=Kfaces(:,:,2)%nz
              Kfaces(:,:,-2)%nz=Kfaces(:,:,2)%nz
            end if

            if(bc%kmax_id==-7) then
              Kfaces(:,:,kmx+0)%nx=Kfaces(:,:,kmx-1)%nx
              Kfaces(:,:,kmx+1)%nx=Kfaces(:,:,kmx-1)%nx
              Kfaces(:,:,kmx+2)%nx=Kfaces(:,:,kmx-1)%nx
              Kfaces(:,:,kmx+3)%nx=Kfaces(:,:,kmx-1)%nx
              Kfaces(:,:,kmx+0)%ny=Kfaces(:,:,kmx-1)%ny
              Kfaces(:,:,kmx+1)%ny=Kfaces(:,:,kmx-1)%ny
              Kfaces(:,:,kmx+2)%ny=Kfaces(:,:,kmx-1)%ny
              Kfaces(:,:,kmx+3)%ny=Kfaces(:,:,kmx-1)%ny
              Kfaces(:,:,kmx+0)%nz=Kfaces(:,:,kmx-1)%nz
              Kfaces(:,:,kmx+1)%nz=Kfaces(:,:,kmx-1)%nz
              Kfaces(:,:,kmx+2)%nz=Kfaces(:,:,kmx-1)%nz
              Kfaces(:,:,kmx+3)%nz=Kfaces(:,:,kmx-1)%nz
            end if

            
        end subroutine normalize_face_normals