Mesh and boundary regions
Mesh regions and boundary regions provide a functionality to identify parts of a mesh or parts of an analytically specified boundary, respectively. These structures are commonly used for the creation of boundary conditions.
Mesh regions
The structure t_meshRegion
from meshregion.f90
allows to identify parts of the mesh like vertices, edges, faces and elements. A mesh region is very often created to identify parts on the boundary where to impose boundary conditions. The following subroutines define the basic maintainance for mesh regions:
Subroutine | Description |
---|---|
mshreg_createFromNodalProp | Create a mesh region based on the nodal property |
mshreg_createFromHitTest | Create a mesh region based on a "hit test" callback subroutine |
mshreg_createFromExpression | Create a mesh region based on a "hit test" expression |
mshreg_done | Destroys a mesh region after use |
Mesh regions based on nodal properties
Using the nodal property allows to quickly collect all vertices, edges and faces on a boundary component. Note that the nodal property of vertices/edges/faces is =0 for inner entities, while a value >0 identifies the boundary component the entity is associated to. A mesh region for boundary entities can quickly be created using mshreg_createFromNodalProp
. The parameter cidxCalc
specifies what to collect; a value MSHREG_IDX_ALL
collects all entities, so vertices, edges and faces.
Example: The following example collects all vertices, edges and faces on the complete boundary:
type(t_triangulation) :: rtriangulation
type(t_meshRegion) :: rmeshRegion
...
call mshreg_createFromNodalProp(rmeshRegion, rtriangulation, MSHREG_IDX_ALL)
...
call msgreg_done (rmeshRegion)
Example: The following example collects all vertices, edges and faces on the 2nd boundary component:
type(t_triangulation) :: rtriangulation
type(t_meshRegion) :: rmeshRegion
...
call mshreg_createFromNodalProp(&
rmeshRegion, rtriangulation, MSHREG_IDX_ALL, (/ 2 /))
...
call msgreg_done (rmeshRegion)
Example: The following example collects only the vertices on the 1st boundary component:
type(t_triangulation) :: rtriangulation
type(t_meshRegion) :: rmeshRegion
...
call mshreg_createFromNodalProp( &
rmeshRegion, rtriangulation, MSHREG_IDX_VERTEX, (/ 1 /))
...
call msgreg_done (rmeshRegion)
Mesh regions based on expression hit tests
Expressions provide a very flexible way to fetch entities from the mesh into a mesh region. The subroutine mshreg_createFromExpression
provides this functionality. An expression must be specified, which has to return "1" if a point (identified by an x/y/z coordinate) should be put into the mesh region. Edges, faces and elements are identified via their midpoint. Via an additional parameter bonlyBndry
it can be specified whether ot not only entities on the boundary should be fetched.
Example: The following example collects all boundary vertices, edges and faces which have an X-coordinate <=0 and an Y-coordinate <= 0
type(t_triangulation) :: rtriangulation
type(t_meshRegion) :: rmeshRegion
...
call mshreg_createFromExpression( &
rmeshRegion, rtriangulation, MSHREG_IDX_ALL, .true., &
"IF (X<=0 AND Y<=0; 1 ; 0)" )
...
call msgreg_done (rmeshRegion)
Mesh regions based on general hit tests
For complicated hit tests, the subroutine mshreg_createFromHitTest
provides a fast functionality to let a callback subroutine decide on whether a point is taken into a mesh region or not. The specified callback routine is called with the x/y/z coordinate of a set of points to check; edges, faces and elements are identified via their midpoint. The routine has to return whether for every point if it is "in" the mesh region or not. Via an additional parameter bonlyBndry
it can be specified whether ot not only entities on the boundary should be fetched.
Example: The following 2D example collects all boundary vertices and edges which have a Y-coordinate <=0 and an X-coordinate <= 0.
subroutine fmshregHitTest(inumCells,Dcoords,Ihit,rcollection)
integer, intent(in) :: inumCells
real(DP), dimension(:,:), intent(in) :: Dcoords
integer, dimension(:), intent(out) :: Ihit
type(t_collection), intent(inout), optional :: rcollection
integer :: i
real(DP) :: dx,dy
! Check every point
do i = 1,inumCells
dx = Dcoords(1,i)
dy = Dcoords(1,i)
if ( (dx .le. 0.0_DP) .and. (dy .le. 0.0_DP)) then
Ihit(i) = 1
else
Ihit(i) = 0
end do
end do
end subroutine
type(t_triangulation) :: rtriangulation
type(t_meshRegion) :: rmeshRegion
...
call mshreg_createFromHitTest( &
rmeshRegion, rtriangulation, MSHREG_IDX_ALL, .true., &
fmshregHitTest )
...
call msgreg_done (rmeshRegion)
Boundary regions
The structure t_boundaryRegion
from boundary.f90
allows to identify regions on the analytical boundary via a start parameter value and an end parameter value. A boundary region is very often created to identify parts on the boundary where to impose boundary conditions. The structure has basically the following entries:
Variable | Description |
---|---|
(t_boundaryRegion)%iboundCompIdx | Boundary component |
(t_boundaryRegion)%dminParam | Minimum parameter value |
(t_boundaryRegion)%dmaxParam | Maximum parameter value |
(t_boundaryRegion)%cparType | Type of parametrisation; 0-1 is the default |
(t_boundaryRegion)%iproperties | Specifies whether start/endpoint are included or not. By default, the starting point is included, the endpoint not. |
The structure can be set up manually or using the subroutines from boundary.f90
. The following subroutines support boundary regions:
Subroutine | Description |
---|---|
boundary_createRegion | Creates a boundary region for a boundary component or a boundary segment |
boundary_getRegionLength | Calculates the length of a boundary region |
boundary_isInRegion | Checks if a point is in a boundary region |
Boundary regions do not have to be destroyed, they do not contain allocatable data.
By default, a boundary region includes the starting point but excludes the endpoint. The variable (t_boundaryRegion)%iproperties
decides on this. The following values are allowed:
Value | Description |
---|---|
0 | t_boundaryRegion describes (dminParam,dmaxParam) |
BDR_PROP_WITHSTART | t_boundaryRegion describes [dminParam,dmaxParam) |
BDR_PROP_WITHEND | t_boundaryRegion describes (dminParam,dmaxParam] |
BDR_PROP_WITHSTART + BDR_PROP_WITHEND |
t_boundaryRegion describes [dminParam,dmaxParam] |
Example: The following piece of code creates a boundary region that covers the complete boundary component 1:
type(t_boundary) :: rboundary
type(t_boundaryRegion) :: rregion
...
call boundary_createRegion (rboundary, 1, 0, rregion)
Example: The following piece of code creates a boundary region that covers the 3rd segment of boundary component 1 including the starting point but excluding the endpoint:
type(t_boundary) :: rboundary
type(t_boundaryRegion) :: rregion
...
call boundary_createRegion (rboundary, 1, 3, rregion)
Example: The following piece of code creates a boundary region that covers the parameter set [1,3] on boundary component 1 including starting point and endpoint:
type(t_boundary) :: rboundary
type(t_boundaryRegion) :: rregion
...
rregion%iboundCompIdx = 1
rregion%dminParam = 1.0_DP
rregion%dminParam = 3.0_DP
rregion%iproperties = BDR_PROP_WITHSTART + BDR_PROP_WITHEND