Identify elements matching a given shape to build a subTree.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_shape_type), | intent(in) | :: | me |
The shape to identify elements for |
||
integer, | intent(in) | :: | iShape |
Numbering of the shape (only for logging) |
||
type(treelmesh_type), | intent(in) | :: | inTree |
The tree to look in for elements that match the shape definition |
||
logical, | intent(in) | :: | storePnts |
Whether to store point values |
||
type(dyn_intarray_type), | intent(inout) | :: | map2global |
Mapping to global elements in the tree |
||
type(tem_grwPoints_type), | intent(inout) | :: | grwPnts |
Growing list of Points to be observed |
||
integer, | intent(out) | :: | countElems(globalMaxLevels) |
Number of elements on each level matching the shape |
||
integer, | intent(out) | :: | countPnts |
Number of points to be observed |
||
integer, | intent(out), | allocatable | :: | bcIDs(:) |
Field of boundary ids that are to be tracked |
|
type(tem_BC_prop_type), | intent(in), | optional | :: | bc_prop |
Boundary condition property to identify boundary elements |
|
type(tem_stencilHeader_type), | intent(in), | optional | :: | stencil |
Stencil associated with the boundary to find respective neighbors |
subroutine tem_shape2subTree( me, iShape, inTree, storePnts, map2global, & & grwPnts, countElems, countPnts, bcIDs, & & bc_prop, stencil ) ! ---------------------------------------------------------------------- ! !> The shape to identify elements for type(tem_shape_type), intent(in) :: me !> Numbering of the shape (only for logging) integer, intent(in) :: iShape !> The tree to look in for elements that match the shape definition type(treelmesh_type), intent(in) :: inTree !> Whether to store point values logical, intent(in) :: storePnts !> Mapping to global elements in the tree type(dyn_intArray_type), intent(inout) :: map2global !> Growing list of Points to be observed type(tem_grwPoints_type), intent(inout) :: grwPnts !> Number of elements on each level matching the shape integer, intent(out) :: countElems(globalMaxLevels) !> Number of points to be observed integer, intent(out) :: countPnts !> Field of boundary ids that are to be tracked integer, allocatable, intent(out) :: bcIDs(:) !> Boundary condition property to identify boundary elements type(tem_bc_prop_type), optional, intent(in) :: bc_prop !> Stencil associated with the boundary to find respective neighbors type(tem_stencilHeader_type), optional, intent(in) :: stencil ! ---------------------------------------------------------------------- ! logical :: foundAny integer :: nShapeElems(globalMaxLevels) ! ---------------------------------------------------------------------- ! foundAny = .false. nShapeElems = 0 select case( me%shapeID ) case( tem_geometrical_shape ) ! Use elements intersecting a geometrical object write(logUnit(5),*) 'iShape ', iShape, ' is a geometrical shape.' call tem_shape_subTreeFromGeomInters( me = me, & & inTree = inTree, & & countElems = nShapeElems, & & countPoints = countPnts, & & grwPnts = grwPnts, & & storePnts = storePnts, & & map2global = map2global ) case( tem_property_shape ) ! Only use elements with a certain property write(logUnit(5),*) 'iShape ', iShape, ' is a property shape.' call tem_shape_initPropElements( me%propBits, inTree, & & nShapeElems, map2global ) case( tem_boundary_shape ) ! Only use elements belong to certain boundaries write(logUnit(5),*) 'iShape ', iShape, ' is a boundary shape.' if (present(bc_prop) .and. present(stencil)) then call tem_shape_findElemByBCLabels( bcLabels = me%bcLabels, & & cutOffQVal = me%cutOffQVal, & & bc_prop = bc_prop, & & foundAny = foundAny, & & map2global = map2Global, & & inTree = inTree, & & countPoints = countPnts, & & grwPnts = grwPnts, & & storePnts = storePnts, & & bcIDs = bcIDs, & & stencil = stencil ) if (foundAny) nShapeElems = 1 else call tem_abort('In tem_shape2subTree: Stencil or bc_prop not passed!') end if case( tem_level_shape ) write(logUnit(5),*) 'iShape ', iShape, ' is a level shape.' call tem_shape_initByLevels( inTree, me%minLevel, & & me%maxLevel, & & nShapeElems, map2global ) end select ! shapeID countElems = countElems + nShapeElems end subroutine tem_shape2subTree