check_serial_singlelevel_faceDesc Subroutine

subroutine check_serial_singlelevel_faceDesc()

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

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

Arguments

None

Calls

proc~~check_serial_singlelevel_facedesc~~CallsGraph proc~check_serial_singlelevel_facedesc check_serial_singlelevel_faceDesc proc~load_env load_env proc~check_serial_singlelevel_facedesc->proc~load_env interface~positionofval~5 positionofval proc~check_serial_singlelevel_facedesc->interface~positionofval~5 proc~tem_finalize tem_finalize proc~check_serial_singlelevel_facedesc->proc~tem_finalize proc~tem_build_face_info tem_build_face_info proc~check_serial_singlelevel_facedesc->proc~tem_build_face_info proc~tem_abort tem_abort proc~check_serial_singlelevel_facedesc->proc~tem_abort proc~tem_start tem_start proc~load_env->proc~tem_start 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~load_tem load_tem proc~load_env->proc~load_tem proc~open_config_chunk open_config_chunk proc~load_env->proc~open_config_chunk proc~tem_load_general tem_load_general proc~load_env->proc~tem_load_general proc~load_tem_bc_prop load_tem_BC_prop proc~load_env->proc~load_tem_bc_prop proc~posofval_label posofval_label interface~positionofval~5->proc~posofval_label proc~tem_finalize->proc~tem_abort proc~tem_timer_dump_glob tem_timer_dump_glob proc~tem_finalize->proc~tem_timer_dump_glob proc~tem_status_run_terminate tem_status_run_terminate proc~tem_finalize->proc~tem_status_run_terminate proc~tem_gettimerval tem_getTimerVal proc~tem_finalize->proc~tem_gettimerval proc~fin_env fin_env proc~tem_finalize->proc~fin_env proc~tem_gettimername tem_getTimerName proc~tem_finalize->proc~tem_gettimername proc~print_self_status print_self_status proc~tem_finalize->proc~print_self_status proc~tem_build_facebuffers tem_build_faceBuffers proc~tem_build_face_info->proc~tem_build_facebuffers proc~tem_dimbydim_construction tem_dimByDim_construction proc~tem_build_face_info->proc~tem_dimbydim_construction 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

Called by

proc~~check_serial_singlelevel_facedesc~~CalledByGraph proc~check_serial_singlelevel_facedesc check_serial_singlelevel_faceDesc program~tem_face_test~2 tem_face_test program~tem_face_test~2->proc~check_serial_singlelevel_facedesc

Contents


Source Code

  subroutine check_serial_singlelevel_faceDesc()
!    character(len=*) :: filename
    integer i
    integer :: pos
    type(treelmesh_type) :: tree
    type(tem_bc_prop_type) :: boundary 
    type(tem_general_type) :: general
    type(tem_face_type), allocatable :: faces(:)
    integer :: level
 !   type(tem_face_type) :: faces(tree%global%minLevel:tree%global%maxLevel)
      

    ! Load the environment.
    call load_env(tree = tree, boundary = boundary, general = general ) 
    level = tree%global%minLevel

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

    write(*,*) tree%nElems

    do i=1, tree%nElems
       write(*,*) tree%treeID(i)
    end do

    write(*,*) tree%part_Last(1)
    write(*,*) tree%elemOffset


    ! Start our checks:
    
    !!!!!!!!!!!!!! Bottom Plane z=0 !!!!!!!!!!!!!!!!!!!!!!!!!!!!
    ! We only have single level, so iterate over all these faces.
    ! First, we start with the faces aligned in x direction.
    
    ! Let us have a look if the number of faces is 8.
    write(logUnit(1),*)'X-direction!'
    if( size(faces(level)%faces(1)%computeFace%leftPos).ne.8 ) then
      write(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'expected 8 faces, but face description is wrong ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'number of compute faces is 8.'
    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(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 2-1 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'found face 2-1 in face list.'
      if(faces(level)%faces(1)%faceList%rightElemId%val(pos).ne.1_long_k) then
        write(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 2-1 in face list (wrong ' // &
                  & 'right elem id) ...'
        call tem_abort()
      else 
        write(logUnit(1),*)'SUCCESS in check_serial_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 8.
    write(logUnit(1),*)'Y-direction'
    if( size(faces(level)%faces(2)%computeFace%leftPos).ne.8 ) then
      write(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'expected 8 faces, but face description is wrong ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'number of compute faces is 8.'
    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(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 3-1 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'found face 3-1 in face list.'

      if(faces(level)%faces(2)%faceList%rightElemId%val(pos).ne.1_long_k) then
        write(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 3-1 in face list (wrong ' // &
                  & 'right elem id) ...'
        call tem_abort()
      else 
        write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'found face 3-1 in face list, correct right elem id.'
      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( size(faces(level)%faces(3)%computeFace%leftPos).ne.8 ) then
      write(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'expected 8 faces, but face description is wrong ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'number of compute faces 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_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 1-5 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_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_serial_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_serial_singlelevel_faceDesc : ' // &
                  & 'found face 1-5 in face list, correct right elem id.'
      end if
    end if



    ! Checking in z-direction    
    ! 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_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 4-8 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_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_serial_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_serial_singlelevel_faceDesc : ' // &
                  & 'found face 4-8 in face list, correct right elem id.'
      end if
    end if


    !!!!!!!!!!!!!!!!! Top plane Z=1 !!!!!!!!!!!!!!!!!!!!!!!!!!!
    ! We only have single level, so iterate over all these faces.
    ! First, we start with the faces aligned in x direction.
    write(logUnit(1),*)'X-direction'
     
    ! 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(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 6-5 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'found face 6-5 in face list.'
      if(faces(level)%faces(1)%faceList%rightElemId%val(pos).ne.5_long_k) then
        write(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 6-5 in face list (wrong ' // &
                  & 'right elem id) ...'
        call tem_abort()
      else 
        write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'found face 6-5 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 8.
    write(logUnit(1),*)'Y-direction'
    if( size(faces(level)%faces(2)%computeFace%leftPos).ne.8 ) then
      write(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'expected 8 faces, but face description is wrong ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'number of compute faces is 8.'
    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(logUnit(1),*)'ERROR in check_serial_singlelevel_faceDesc : ' // &
                  & 'not able to locate face 7-5 in face list ...'
      call tem_abort()
    else 
      write(logUnit(1),*)'SUCCESS in check_serial_singlelevel_faceDesc : ' // &
                  & 'found face 7-5 in face list.'

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

    call tem_finalize(general)

    ! When we reach this point, all tests have passed...
    write(logUnit(1),*)'PASSED'

  end subroutine check_serial_singlelevel_faceDesc