grans_obj_PrepareSolverBCS Subroutine

subroutine grans_obj_PrepareSolverBCS(this)

Prepares the boundary conditions used by GRANS. Note that the only boundary conditions required by GRANS are those in connection with the volume fraction Poisson equations. This subroutine reads pre-existing BCS file on disk. It adds additional boundary conditions for ibF, ibN, and ibS based on the specified boundary conditions for ibVF. If no Dirichlet conditions are given for ibVF, it marks the volume fraction Poisson equation as singular.

Arguments

Type IntentOptional Attributes Name
class(grans_obj), intent(inout) :: this

GRANS solver


Calls

proc~~grans_obj_preparesolverbcs~~CallsGraph proc~grans_obj_preparesolverbcs grans_obj_PrepareSolverBCS none~get~4 parser_obj%Get proc~grans_obj_preparesolverbcs->none~get~4 none~lor parallel_obj%LOR proc~grans_obj_preparesolverbcs->none~lor proc~bc_set_checkbcexists bc_set%bc_set_CheckBCExists proc~grans_obj_preparesolverbcs->proc~bc_set_checkbcexists proc~bc_set_getbcpointer bc_set%bc_set_GetBCPointer proc~grans_obj_preparesolverbcs->proc~bc_set_getbcpointer proc~bc_set_getbctype bc_set%bc_set_GetBCType proc~grans_obj_preparesolverbcs->proc~bc_set_getbctype proc~bc_set_init bc_set%bc_set_Init proc~grans_obj_preparesolverbcs->proc~bc_set_init proc~bc_set_setbc bc_set%bc_set_SetBC proc~grans_obj_preparesolverbcs->proc~bc_set_setbc proc~parser_obj_read0d parser_obj%parser_obj_read0D none~get~4->proc~parser_obj_read0d proc~parser_obj_read1d parser_obj%parser_obj_read1D none~get~4->proc~parser_obj_read1d proc~parallel_obj_lor_0d parallel_obj%parallel_obj_Lor_0d none~lor->proc~parallel_obj_lor_0d proc~parallel_obj_lor_1d parallel_obj%parallel_obj_Lor_1d none~lor->proc~parallel_obj_lor_1d proc~bc_set_getregionindex bc_set%bc_set_GetRegionIndex proc~bc_set_checkbcexists->proc~bc_set_getregionindex proc~region_obj_getbcindex region_obj%region_obj_GetBCIndex proc~bc_set_checkbcexists->proc~region_obj_getbcindex proc~bc_set_getbcpointer->proc~bc_set_getregionindex proc~bc_set_getbcpointer->proc~region_obj_getbcindex proc~bc_set_getbctype->proc~bc_set_getregionindex proc~bc_set_getbctype->proc~region_obj_getbcindex proc~hashtbl_obj_init hashtbl_obj%hashtbl_obj_Init proc~bc_set_init->proc~hashtbl_obj_init proc~bc_set_getextents bc_set%bc_set_GetExtents proc~bc_set_setbc->proc~bc_set_getextents proc~bc_set_setbc->proc~bc_set_getregionindex proc~region_obj_add region_obj%region_obj_Add proc~bc_set_setbc->proc~region_obj_add proc~bc_set_getextents->proc~bc_set_getregionindex none~get~3 hashtbl_obj%Get proc~bc_set_getregionindex->none~get~3 proc~hashtbl_obj_hashstring hashtbl_obj%hashtbl_obj_HashString proc~bc_set_getregionindex->proc~hashtbl_obj_hashstring mpi_allreduce mpi_allreduce proc~parallel_obj_lor_0d->mpi_allreduce proc~parallel_obj_lor_1d->mpi_allreduce none~assigndefault parser_obj%AssignDefault proc~parser_obj_read0d->none~assigndefault proc~parser_obj_fetchlabelid parser_obj%parser_obj_FetchLabelID proc~parser_obj_read0d->proc~parser_obj_fetchlabelid proc~parser_obj_read1d->none~assigndefault proc~parser_obj_read1d->proc~parser_obj_fetchlabelid proc~region_obj_add->proc~hashtbl_obj_hashstring proc~hashtbl_obj_put hashtbl_obj%hashtbl_obj_Put proc~region_obj_add->proc~hashtbl_obj_put proc~region_obj_expand region_obj%region_obj_Expand proc~region_obj_add->proc~region_obj_expand proc~region_obj_getbcindex->none~get~3 proc~region_obj_getbcindex->proc~hashtbl_obj_hashstring proc~parser_obj_assigndefault0d parser_obj%parser_obj_AssignDefault0D none~assigndefault->proc~parser_obj_assigndefault0d proc~parser_obj_assigndefault1d parser_obj%parser_obj_AssignDefault1D none~assigndefault->proc~parser_obj_assigndefault1d proc~hashtbl_obj_get_int4 hashtbl_obj%hashtbl_obj_Get_int4 none~get~3->proc~hashtbl_obj_get_int4 proc~hashtbl_obj_get_int8 hashtbl_obj%hashtbl_obj_Get_int8 none~get~3->proc~hashtbl_obj_get_int8 proc~hashtbl_obj_get_real_dp hashtbl_obj%hashtbl_obj_Get_real_dp none~get~3->proc~hashtbl_obj_get_real_dp proc~hashtbl_obj_get_real_sp hashtbl_obj%hashtbl_obj_Get_real_sp none~get~3->proc~hashtbl_obj_get_real_sp proc~sllist_obj_put sllist_obj%sllist_obj_Put proc~hashtbl_obj_put->proc~sllist_obj_put none~get~2 sllist_obj%Get proc~hashtbl_obj_get_int4->none~get~2 proc~hashtbl_obj_get_int8->none~get~2 proc~hashtbl_obj_get_real_dp->none~get~2 proc~hashtbl_obj_get_real_sp->none~get~2 proc~sllist_obj_put->proc~sllist_obj_put proc~sllist_obj_get_int4 sllist_obj%sllist_obj_Get_int4 none~get~2->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8 sllist_obj%sllist_obj_Get_int8 none~get~2->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp sllist_obj%sllist_obj_Get_real_dp none~get~2->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp sllist_obj%sllist_obj_Get_real_sp none~get~2->proc~sllist_obj_get_real_sp proc~sllist_obj_get_int4->proc~sllist_obj_get_int4 proc~sllist_obj_get_int8->proc~sllist_obj_get_int8 proc~sllist_obj_get_real_dp->proc~sllist_obj_get_real_dp proc~sllist_obj_get_real_sp->proc~sllist_obj_get_real_sp

Called by

proc~~grans_obj_preparesolverbcs~~CalledByGraph proc~grans_obj_preparesolverbcs grans_obj_PrepareSolverBCS proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->proc~grans_obj_preparesolverbcs interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver

Source Code

    subroutine grans_obj_PrepareSolverBCS(this)
      !> Prepares the boundary conditions used by GRANS.
      ! Note that the only boundary conditions required by GRANS
      ! are those in connection with the volume fraction Poisson equations.
      ! This subroutine reads pre-existing BCS file on disk.
      ! It adds additional boundary conditions for ibF, ibN, and ibS based
      ! on the specified boundary conditions for ibVF. If no Dirichlet
      ! conditions are given for ibVF, it marks the volume fraction Poisson
      ! equation as singular.
      implicit none
      class(grans_obj), intent(inout) :: this                                  !! GRANS solver
      ! Work variables
      real(wp)  :: time
      integer   :: iter
      integer   :: n
      real(wp), pointer   :: BCVal(:,:,:)
      character(len=str64):: filename
      logical             :: ibVFDirichlet
      logical             :: buffL

      call this%parser%Get('BC file ',  filename, default= 'bcs' )

      ! Initialize and read boundary conditions
      call this%bcs%Initialize(this%block,this%parallel)
      call this%bcs%Read(iter,time,filename)

      associate (bcs => this%bcs)
        if (this%use_RP .or. this%use_IB) then

          ! Check boundary conditions for ibVF
          ibVFDirichlet = .false.
          do n=1,bcs%count
            ! Check Whether we have a BC for this variable on this region
            if (.not.bcs%CheckBCExists(bcs%region(n)%name,'ibVF')) cycle

            select case (bcs%GetBCType(bcs%region(n)%name,'ibVF'))
            case (BC_DIRICHLET)
              ! Found Dirichlet condition
              ibVFDirichlet = .true.

              call bcs%SetBC(bcs%region(n)%name,BC_DIRICHLET,'ibS' )
              call bcs%SetBC(bcs%region(n)%name,BC_DIRICHLET,'ibN1')
              call bcs%SetBC(bcs%region(n)%name,BC_DIRICHLET,'ibN2')
              call bcs%SetBC(bcs%region(n)%name,BC_DIRICHLET,'ibN3')
              call bcs%SetBC(bcs%region(n)%name,BC_DIRICHLET,'ibF1')
              call bcs%SetBC(bcs%region(n)%name,BC_DIRICHLET,'ibF2')
              call bcs%SetBC(bcs%region(n)%name,BC_DIRICHLET,'ibF3')

              call bcs%GetBCPointer(bcs%region(n)%name,'ibS',  BCVal); BCVal = 0.0_wp
              call bcs%GetBCPointer(bcs%region(n)%name,'ibN1', BCVal); BCVal = 0.0_wp
              call bcs%GetBCPointer(bcs%region(n)%name,'ibN2', BCVal); BCVal = 0.0_wp
              call bcs%GetBCPointer(bcs%region(n)%name,'ibN3', BCVal); BCVal = 0.0_wp
              call bcs%GetBCPointer(bcs%region(n)%name,'ibF1', BCVal); BCVal = 0.0_wp
              call bcs%GetBCPointer(bcs%region(n)%name,'ibF2', BCVal); BCVal = 0.0_wp
              call bcs%GetBCPointer(bcs%region(n)%name,'ibF3', BCVal); BCVal = 0.0_wp
            case (BC_NEUMANN)
            case (BC_SYMMETRY)
              call bcs%SetBC(bcs%region(n)%name,BC_SYMMETRY,'ibS' )
              call bcs%SetBC(bcs%region(n)%name,BC_SYMMETRY,'ibN1')
              call bcs%SetBC(bcs%region(n)%name,BC_SYMMETRY,'ibN2')
              call bcs%SetBC(bcs%region(n)%name,BC_SYMMETRY,'ibN3')
              call bcs%SetBC(bcs%region(n)%name,BC_SYMMETRY,'ibF1')
              call bcs%SetBC(bcs%region(n)%name,BC_SYMMETRY,'ibF2')
              call bcs%SetBC(bcs%region(n)%name,BC_SYMMETRY,'ibF3')
            end select
          end do

          call this%parallel%Lor(ibVFDirichlet,buffL); ibVFDirichlet=buffL

          ! Determine if we need to apply volume fraction corrections
          if (.not.ibVFDirichlet) this%singularVF = .true.
        end if
      end associate

      return
    end subroutine grans_obj_PrepareSolverBCS