check_parallel_singlelevel_faceDesc Subroutine

subroutine check_parallel_singlelevel_faceDesc()

!!!!!!!!!!!! Bottom Plane z=0 !!!!!!!!!!!!!!!!!!!!!!!!!!!!

!!!!!!!!!!!!!!! Top plane Z=1 !!!!!!!!!!!!!!!!!!!!!!!!!!!

ComputeFace checks for

Arguments

None

Calls

proc~~check_parallel_singlelevel_facedesc~~CallsGraph proc~check_parallel_singlelevel_facedesc check_parallel_singlelevel_faceDesc proc~load_env load_env proc~check_parallel_singlelevel_facedesc->proc~load_env interface~positionofval~5 positionofval proc~check_parallel_singlelevel_facedesc->interface~positionofval~5 proc~tem_build_face_info tem_build_face_info proc~check_parallel_singlelevel_facedesc->proc~tem_build_face_info proc~tem_abort tem_abort proc~check_parallel_singlelevel_facedesc->proc~tem_abort proc~computeface_to_treeidcheck ComputeFace_to_TreeIDcheck proc~check_parallel_singlelevel_facedesc->proc~computeface_to_treeidcheck proc~tem_load_general tem_load_general proc~load_env->proc~tem_load_general proc~open_config_chunk open_config_chunk proc~load_env->proc~open_config_chunk proc~load_tem load_tem proc~load_env->proc~load_tem proc~tem_start tem_start proc~load_env->proc~tem_start proc~load_tem_bc_prop load_tem_BC_prop proc~load_env->proc~load_tem_bc_prop proc~tem_logging_load_primary tem_logging_load_primary proc~load_env->proc~tem_logging_load_primary proc~close_config close_config proc~load_env->proc~close_config proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~tem_dimbydim_construction tem_dimByDim_construction proc~tem_build_face_info->proc~tem_dimbydim_construction proc~tem_build_facebuffers tem_build_faceBuffers proc~tem_build_face_info->proc~tem_build_facebuffers proc~tem_facedep_vertical tem_faceDep_vertical proc~tem_build_face_info->proc~tem_facedep_vertical proc~tem_build_facelists tem_build_faceLists proc~tem_build_face_info->proc~tem_build_facelists proc~tem_collect_faces tem_collect_faces proc~tem_build_face_info->proc~tem_collect_faces proc~tem_extend_remoteprp tem_extend_remotePrp proc~tem_build_face_info->proc~tem_extend_remoteprp mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~computeface_to_treeidcheck->proc~tem_abort

Called by

proc~~check_parallel_singlelevel_facedesc~~CalledByGraph proc~check_parallel_singlelevel_facedesc check_parallel_singlelevel_faceDesc program~tem_face_test tem_face_test program~tem_face_test->proc~check_parallel_singlelevel_facedesc

Contents


Source Code

  subroutine check_parallel_singlelevel_faceDesc()
    integer :: pos
    type(tem_general_type) :: general
    type(treelmesh_type) :: tree
    type(tem_bc_prop_type) :: boundary 
    type(tem_face_type), allocatable :: faces(:)
    integer :: level
      
    ! Load the environment.
    call load_env(tree = tree, boundary = boundary, general = general ) 
    level = tree%global%minLevel

    ! Create face description
    write(*,*) 'Building face info'
    call tem_build_face_info( tree = tree, boundary = boundary, &
                            & commPattern = general%commPattern, &
                            & proc = general%proc, &
                            & faces = faces, nEligibleChildren = 4 )

    write(*,*) 'Rank',  general%proc%rank

   ! Start our checks:
  
   !!!!!!!!!!!!!! Bottom Plane z=0 !!!!!!!!!!!!!!!!!!!!!!!!!!!!
   ! process_0 contains plane z=0  for nprocs=2
   ! Check faces in x and y-direction 
   if(general%proc%rank.eq.0)  then
      ! Let us have a look if the number of faces is 4.
      write(*,*) 'Process : ', general%proc%rank
      write(*,*) "X-direction!"
      if( faces(level)%faces(1)%faceList%faceId%nvals.ne.4 ) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : //&
                    &expected 4 faces, but face description is wrong ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'number of nvals is 4.'
      end if
       
      ! Let us have a look if the face 2-1 (left elem id - right elem id) is
      ! in the list of faces.
      pos = PositionOfVal( me = faces(level)%faces(1)%faceList%faceId, val = 2_long_k )
      if(pos .le. 0) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'not able to locate face 2-1 in face list ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'found face 2-1 in face list.'
        if(faces(level)%faces(1)%faceList%rightElemId%val(pos).ne.1_long_k) then
          write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'not able to locate face 2-1 in face list (wrong ' // &
                    & 'right elem id) ...'
          call tem_abort()
        else 
          write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'found face 2-1 in face list, correct right elem id.'
        end if
      end if

      ! Checking in y-direction    
      ! Let us have a look if the number of faces is 4.
      write(*,*) 'Y-direction'
      if(  faces(level)%faces(2)%faceList%faceId%nvals.ne.4 ) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'expected 4 faces, but face description is wrong ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'number of nvals is 4.'
      end if
       
      ! Let us have a look if the face 3-1 (left elem id - right elem id) is
      ! in the list of faces.
      pos = PositionOfVal( me = faces(level)%faces(2)%faceList%faceId, val = 3_long_k )
      if(pos .le. 0) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'not able to locate face 3-1 in face list ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'found face 3-1 in face list.'
  
        if(faces(level)%faces(2)%faceList%rightElemId%val(pos).ne.1_long_k) then
          write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'not able to locate face 3-1 in face list (wrong ' // &
                    & 'right elem id) ...'
          call tem_abort()
        else 
          write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'found face 3-1 in face list, correct right elem id.'
        end if
      end if
    end if



    ! Checking in z-direction    
    ! Let us have a look if the number of faces is 8.
    write(logUnit(1),*)'Z-direction'
    if( faces(level)%faces(3)%faceList%faceId%nvals.ne.8 ) then
      write(logUnit(1),*)'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                  & 'expected 8 faces, but face description is wrong ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                  & 'number of nvals is 8.'
    end if
     
    ! Let us have a look if the face 1-5 (left elem id - right elem id) is
    ! in the list of faces.
    pos = PositionOfVal( me = faces(level)%faces(3)%faceList%faceId, val = 1_long_k )
    if(pos .le. 0) then
      write(logUnit(1),*)'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 1-5 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                  & 'found face 1-5 in face list.'

      if(faces(level)%faces(3)%faceList%rightElemId%val(pos).ne.5_long_k) then
        write(logUnit(1),*)'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 1-5 in face list (wrong ' // &
                  & 'right elem id) ...'
        call tem_abort()
      else 
        write(logUnit(1),*)'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                  & 'found face 1-5 in face list, correct right elem id.'
      end if
    end if

    ! Let us have a look if the face 4-8 (left elem id - right elem id) is
    ! in the list of faces.
    pos = PositionOfVal( me = faces(level)%faces(3)%faceList%faceId, val = 4_long_k )
    if(pos .le. 0) then
      write(logUnit(1),*)'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 4-8 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                  & 'found face 4-8 in face list.'

      if(faces(level)%faces(3)%faceList%rightElemId%val(pos).ne.8_long_k) then
        write(logUnit(1),*)'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 4-8 in face list (wrong ' // &
                  & 'right elem id) ...'
        call tem_abort()
      else 
        write(logUnit(1),*)'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                  & 'found face 4-8 in face list, correct right elem id.'
      end if
    end if


   !!!!!!!!!!!!!!!!! Top plane Z=1 !!!!!!!!!!!!!!!!!!!!!!!!!!!
   ! process_1 contains the plane Z=1 for nprocs=2
   ! Check faces in x and y-direction
    if(general%proc%rank.eq.1) then
      ! x-direction
      write(*,*) 'Process :', general%proc%rank
      write(*,*) 'X-direction'
      if( faces(level)%faces(1)%faceList%faceId%nvals.ne.4 ) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'expected 4 faces, but face description is wrong ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'number of nvals is 4.'
      end if

      ! Let us have a look if the face 6-5 (left elem id - right elem id) is
      ! in the list of faces.
      pos = PositionOfVal( me = faces(level)%faces(1)%faceList%faceId, val = 6_long_k )
      if(pos .le. 0) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'not able to locate face 6-5 in face list ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'found face 6-5 in face list.'
        if(faces(level)%faces(1)%faceList%rightElemId%val(pos).ne.5_long_k) then
          write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                     & 'not able to locate face 6-5 in face list (wrong ' // &
                     & 'right elem id) ...'
          call tem_abort()
        else 
          write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                     & 'found face 6-5 in face list, correct right elem id.'
        end if
      end if
  
  
      ! Checking in y-direction    
      write(*,*) 'y-direction'
      if( faces(level)%faces(2)%faceList%faceId%nvals.ne.4 ) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'expected 4 faces, but face description is wrong ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'number of nvals is 4.'
      end if
      ! Let us have a look if the face 7-5 (left elem id - right elem id) is
      ! in the list of faces.
      pos = PositionOfVal( me = faces(level)%faces(2)%faceList%faceId, val = 7_long_k )
      if(pos .le. 0) then
        write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                    & 'not able to locate face 7-5 in face list ...'
        call tem_abort()
      else 
        write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                    & 'found face 7-5 in face list.'
  
        if(faces(level)%faces(2)%faceList%rightElemId%val(pos).ne.5_long_k) then
          write(*,*) 'ERROR in check_parallel_singlelevel_faceDesc : ' // &
                     & 'not able to locate face 7-5 in face list (wrong ' // &
                     & 'right elem id) ...'
          call tem_abort()
        else 
          write(*,*) 'SUCCESS in check_parallel_singlelevel_faceDesc : ' // &
                     & 'found face 7-5 in face list, correct right elem id.'
        end if
      end if
    end if



    !! ComputeFace checks for
    ! process 0
    write(logUnit(1),*)"ComputeFace Checks in z-direction"
    call ComputeFace_to_TreeIDcheck(1,5,0, general%proc, faces=faces)
    call ComputeFace_to_TreeIDcheck(2,6,0, general%proc, faces=faces)
    call ComputeFace_to_TreeIDcheck(3,7,0, general%proc, faces=faces)
    call ComputeFace_to_TreeIDcheck(4,8,0, general%proc, faces=faces)
    ! process 1
    call ComputeFace_to_TreeIDcheck(5,1,1, general%proc, faces=faces)
    call ComputeFace_to_TreeIDcheck(6,2,1, general%proc, faces=faces)
    call ComputeFace_to_TreeIDcheck(7,3,1, general%proc, faces=faces)
    call ComputeFace_to_TreeIDcheck(8,4,1, general%proc, faces=faces)

  end subroutine check_parallel_singlelevel_faceDesc