bc_set_UpdateBoundaryScalar Subroutine

private impure subroutine bc_set_UpdateBoundaryScalar(this, var)

Imposes boundary conditions for a scalar variable.

Type Bound

bc_set

Arguments

Type IntentOptional Attributes Name
class(bc_set), intent(in) :: this

Boundary conditions utility

type(eulerian_obj_r), intent(inout) :: var

Eulerian variable


Calls

proc~~bc_set_updateboundaryscalar~~CallsGraph proc~bc_set_updateboundaryscalar bc_set%bc_set_UpdateBoundaryScalar proc~bc_set_checkbcexists bc_set%bc_set_CheckBCExists proc~bc_set_updateboundaryscalar->proc~bc_set_checkbcexists proc~bc_set_getbctype bc_set%bc_set_GetBCType proc~bc_set_updateboundaryscalar->proc~bc_set_getbctype proc~bc_set_updateboundarydirichlet bc_set%bc_set_UpdateBoundaryDirichlet proc~bc_set_updateboundaryscalar->proc~bc_set_updateboundarydirichlet proc~bc_set_updateboundaryneumann bc_set%bc_set_UpdateBoundaryNeumann proc~bc_set_updateboundaryscalar->proc~bc_set_updateboundaryneumann proc~bc_set_updateboundarysymmetryplus bc_set%bc_set_UpdateBoundarySymmetryPlus proc~bc_set_updateboundaryscalar->proc~bc_set_updateboundarysymmetryplus proc~eulerian_obj_updateghostcells eulerian_obj_base%eulerian_obj_UpdateGhostCells proc~bc_set_updateboundaryscalar->proc~eulerian_obj_updateghostcells 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_getbctype->proc~bc_set_getregionindex proc~bc_set_getbctype->proc~region_obj_getbcindex proc~bc_set_getbcpointer bc_set%bc_set_GetBCPointer proc~bc_set_updateboundarydirichlet->proc~bc_set_getbcpointer proc~bc_set_getextents bc_set%bc_set_GetExtents proc~bc_set_updateboundarydirichlet->proc~bc_set_getextents proc~bc_set_getsidedirbyregion bc_set%bc_set_GetSideDirByRegion proc~bc_set_updateboundarydirichlet->proc~bc_set_getsidedirbyregion proc~bc_set_updateboundaryneumann->proc~bc_set_getextents proc~bc_set_updateboundaryneumann->proc~bc_set_getsidedirbyregion proc~bc_set_updateboundarysymmetryplus->proc~bc_set_getextents proc~bc_set_updateboundarysymmetryplus->proc~bc_set_getsidedirbyregion proc~eulerian_obj_updateghostcells_x eulerian_obj_base%eulerian_obj_UpdateGhostCells_x proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_x proc~eulerian_obj_updateghostcells_y eulerian_obj_base%eulerian_obj_UpdateGhostCells_y proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_y proc~eulerian_obj_updateghostcells_z eulerian_obj_base%eulerian_obj_UpdateGhostCells_z proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_z proc~bc_set_getbcpointer->proc~bc_set_getregionindex proc~bc_set_getbcpointer->proc~region_obj_getbcindex 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 proc~bc_set_getsidedirbyregion->proc~bc_set_getregionindex cell cell proc~eulerian_obj_updateghostcells_x->cell mpi_irecv mpi_irecv proc~eulerian_obj_updateghostcells_x->mpi_irecv mpi_isend mpi_isend proc~eulerian_obj_updateghostcells_x->mpi_isend mpi_waitall mpi_waitall proc~eulerian_obj_updateghostcells_x->mpi_waitall proc~eulerian_obj_updateghostcells_y->cell proc~eulerian_obj_updateghostcells_y->mpi_irecv proc~eulerian_obj_updateghostcells_y->mpi_isend proc~eulerian_obj_updateghostcells_y->mpi_waitall proc~eulerian_obj_updateghostcells_z->cell proc~eulerian_obj_updateghostcells_z->mpi_irecv proc~eulerian_obj_updateghostcells_z->mpi_isend proc~eulerian_obj_updateghostcells_z->mpi_waitall proc~region_obj_getbcindex->none~get~3 proc~region_obj_getbcindex->proc~hashtbl_obj_hashstring 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 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_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~~bc_set_updateboundaryscalar~~CalledByGraph proc~bc_set_updateboundaryscalar bc_set%bc_set_UpdateBoundaryScalar none~updateboundary bc_set%UpdateBoundary none~updateboundary->proc~bc_set_updateboundaryscalar proc~cdifs_obj_advancesolutioncorrector cdifs_obj_AdvanceSolutionCorrector proc~cdifs_obj_advancesolutioncorrector->none~updateboundary proc~cdifs_obj_advancesolutionib cdifs_obj_AdvanceSolutionIB proc~cdifs_obj_advancesolutionib->none~updateboundary proc~cdifs_obj_advancesolutionpredictor cdifs_obj_AdvanceSolutionPredictor proc~cdifs_obj_advancesolutionpredictor->none~updateboundary proc~cdifs_obj_computesolidvf cdifs_obj%cdifs_obj_ComputeSolidVF proc~cdifs_obj_computesolidvf->none~updateboundary proc~respart_set_updatenormals ResPart_set%ResPart_set_UpdateNormals proc~cdifs_obj_computesolidvf->proc~respart_set_updatenormals proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->none~updateboundary proc~cdifs_obj_preparesolver->proc~cdifs_obj_computesolidvf proc~cdifs_obj_preparesolverbcs cdifs_obj_PrepareSolverBCS proc~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolverbcs proc~cdifs_obj_preparesolverbcs->none~updateboundary proc~grans_obj_computesolidvf grans_obj%grans_obj_ComputeSolidVF proc~grans_obj_computesolidvf->none~updateboundary proc~grans_obj_computesolidvf->proc~respart_set_updatenormals proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->none~updateboundary proc~grans_obj_preparesolver->proc~grans_obj_computesolidvf interface~grans_obj_writeoutputdata grans_obj%grans_obj_WriteOutputData proc~grans_obj_preparesolver->interface~grans_obj_writeoutputdata proc~grans_obj_writeoutputdata grans_obj_WriteOutputData proc~grans_obj_writeoutputdata->none~updateboundary proc~grans_obj_writeoutputdata->proc~grans_obj_computesolidvf proc~marker_set_updatenormals marker_set%marker_set_UpdateNormals proc~marker_set_updatenormals->none~updateboundary interface~cdifs_obj_preparesolver cdifs_obj%cdifs_obj_PrepareSolver interface~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolver interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver interface~grans_obj_writeoutputdata->proc~grans_obj_writeoutputdata proc~cdifs_obj_advancesolution cdifs_obj_AdvanceSolution proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutioncorrector proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionib proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionpredictor proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionrp proc~cdifs_obj_advancesolutionrp->proc~cdifs_obj_computesolidvf proc~respart_set_updatenormals->proc~marker_set_updatenormals interface~cdifs_obj_advancesolution cdifs_obj%cdifs_obj_AdvanceSolution interface~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolution

Source Code

    impure subroutine bc_set_UpdateBoundaryScalar(this,var)
      !> Imposes boundary conditions for a scalar variable.
      class(bc_set),        intent(in)    :: this                              !! Boundary conditions utility
      type(eulerian_obj_r), intent(inout) :: var                               !! Eulerian variable
      ! Work variables
      logical :: found
      integer :: n

      if (.not.allocated(this%region)) return

      do n=1,this%count

        ! Check Whether we have a BC for this variable on this region
        found = this%CheckBCExists(this%region(n)%name,var%name)
        if (.not.found) cycle

        select case (this%GetBCType(this%region(n)%name,var%name))
        case (BC_DIRICHLET)
          call this%UpdateBoundaryDirichlet(this%region(n)%name,var)
        case (BC_NEUMANN  )
          call this%UpdateBoundaryNeumann  (this%region(n)%name,var)
        case (BC_SYMMETRY )
          call this%UpdateBoundarySymmetryPLUS(this%region(n)%name,var)
        end select
      end do

      if (this%count.ne.0.and. (.not.found)) &
        call this%parallel%Stop('Unable to find any boundary conditions for variable '//trim(adjustl(var%name)))

      call var%UpdateGhostCells()
      return
    end subroutine bc_set_UpdateBoundaryScalar