This serves as an simple grid generation for performance or scaling analysis without being obliged to use Seeder. You have to specify the generic grid parameters in the lua file instead of the mesh folder
mesh = { predefined='line',
origin = {0.,0.,0.},
length = 10.,
refinementLevel = 6 }
\n You have to specify the shape 'line', a bounding box origin, its length and also the refinement level, which results in different amount elements in the grid.\n The result of this routine is mainly the treeID list with the additional lists for saving the properties. The generated line will be all elements along the X-Axis with periodicity in Y and Z direction.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(treelmesh_type), | intent(out) | :: | me |
Mesh to generate |
||
real(kind=rk), | intent(in) | :: | origin(3) |
Corner of the cube |
||
real(kind=rk), | intent(in) | :: | length |
Length of cube |
||
integer, | intent(in) | :: | level |
Resolution level |
||
integer, | intent(in) | :: | myPart |
Partition of the caller (starts with 0) |
||
integer, | intent(in) | :: | nParts |
Number of partitions |
||
integer, | intent(in) | :: | comm |
communicator to be used |
subroutine generate_treelm_line( me, origin, length, level, myPart, &
& nParts, comm )
! -------------------------------------------------------------------- !
!> Mesh to generate
type(treelmesh_type), intent(out) :: me
!> Corner of the cube
real(kind=rk), intent(in) :: origin(3)
!> Length of cube
real(kind=rk), intent(in) :: length
!> Resolution level
integer, intent(in) :: level
!> Partition of the caller (starts with 0)
integer, intent(in) :: myPart
!> Number of partitions
integer, intent(in) :: nParts
!> communicator to be used
integer, intent(in) :: comm
! -------------------------------------------------------------------- !
integer(kind=long_k) :: firstID, lastID
integer(kind=long_k) :: share
integer :: remainder
integer :: iPart, iElem
integer :: coord(4)
integer :: lastcoord
! -------------------------------------------------------------------- !
me%global%nParts = nParts
me%global%myPart = myPart
me%global%comm = comm
me%global%origin = origin
me%global%BoundingCubeLength = length
me%global%minLevel = level
me%global%maxLevel = level
me%global%label = 'Generic_Line'
me%global%predefined = 'line'
write(me%global%comment,'(a15,i7,a16,i2,a1)') &
& 'Generated with ', nParts, ' Parts on Level ', level, '.'
me%global%dirname = './'
! Boundary property to define periodic boundary
me%global%nProperties = 1
if (associated(me%global%property)) deallocate(me%global%property)
if (associated(me%property)) deallocate(me%property)
allocate(me%global%Property(me%global%nProperties))
allocate(me%Property(me%global%nProperties))
allocate(me%Part_First(nParts))
allocate(me%Part_Last(nParts))
! Compute the treeIDs of the mesh:
firstID = tem_firstIdAtLevel(level)
lastcoord = 2**level - 1
lastID = tem_IdOfCoord( coord = [lastcoord, 0, 0, level], offset = firstID )
! Total number of elements in this mesh
me%global%nElems = 2_long_k**level
share = me%global%nElems / int(nParts, kind=long_k)
remainder = int(mod(me%global%nElems, int(nParts, kind=long_k)))
! The first partition starts always with the firstID
me%Part_First(1) = firstID
! Up to remainder partitions have share + 1 elements
coord = 0
coord(4) = level
do iPart=2,remainder+1
coord(1) = coord(1) + int(share)
me%Part_Last(iPart-1) = tem_idofcoord(coord, offset = firstID)
coord(1) = coord(1) + 1
me%Part_First(iPart) = tem_idofcoord(coord, offset = firstID)
end do
! The remaining elements get exactly the share elements:
do iPart=remainder+2,nParts
coord(1) = coord(1) + int(share) - 1
me%Part_Last(iPart-1) = tem_idofcoord(coord, offset = firstID)
coord(1) = coord(1) + 1
me%Part_First(iPart) = tem_idofcoord(coord, offset = firstID)
end do
! The last partition ends always with the lastID
me%Part_Last(nParts) = lastID
! Local data:
if (myPart < remainder) then
me%nElems = int(share+1)
me%elemOffset = int(myPart, kind=long_k) * (share+1_long_k)
else
me%nElems = int(share)
me%elemOffset = (int(myPart, kind=long_k) * share) &
& + int(remainder, kind=long_k)
end if
! All elements have (periodic) boundaries
me%Property(1)%nElems = me%nElems
me%Property(1)%offset = me%elemOffset
allocate(me%Property(1)%ElemID(me%Property(1)%nElems))
! Please note, that the tem_bc_prop_module will set boundary conditions
! based on this label accordingly!
me%global%Property(1)%label = 'internal 1D BC'
me%global%Property(1)%bitpos = prp_hasBnd
me%global%Property(1)%nElems = me%Property(1)%nElems
allocate(me%treeID(me%nElems))
allocate(me%ElemPropertyBits(me%nElems))
! Only has boundary property:
me%ElemPropertyBits = ibset(0_long_k, prp_hasBnd)
! Filling the treeIDs:
do iElem=1,me%nElems
me%Property(1)%ElemID(iElem) = iElem
! We can only have 2**20 elements per dimension, so in this case where
! we create a line, we won't exceed the integer limit. Thus we safely can
! cast the long integer elemOffset to a normal integer.
coord(1) = int(me%elemOffset) + iElem - 1
me%treeID(iElem) = tem_idofcoord(coord, offset=firstID)
end do
end subroutine generate_treelm_line