Load the treelmesh
Depending on the configuration given in conf, this will either load a a \ref tem_distributed_octree "tree mesh" from disk with the given prefix or \ref tem_load_internal "create the mesh internally".
to calculate levelWeight, global has to be loaded fisrt, maybe load_global shoule be moved out of this routine
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(treelmesh_type), | intent(out) | :: | me |
Structure to load the mesh to |
||
type(flu_State) | :: | conf |
Directory containing the mesh informations |
|||
integer, | intent(in) | :: | myPart |
Partition to use on the calling process (= MPI Rank in comm) |
||
integer, | intent(in) | :: | nParts |
Number of partitions, the mesh is partitioned into (= Number of MPI processes in comm). |
||
integer, | intent(in) | :: | comm |
MPI Communicator to use now defaults to the one given in the tree%global%comm |
||
real(kind=rk), | intent(in), | optional | :: | levelWeight(globalMaxLevels) |
Balancing weight for elements on different levels. If these weights are present, the mesh is partitioned initially according to the given weights. This can not be combined with given offsets. If offset is also provided they will overwrite the levelWeight distribution. |
|
character(len=*), | intent(out), | optional | :: | meshDir |
output Mesh directory name |
|
integer, | intent(in), | optional | :: | parent |
optional parent handle |
subroutine load_tem( me, conf, myPart, nParts, comm, levelWeight, meshDir, & & parent) ! -------------------------------------------------------------------- ! !> Structure to load the mesh to type(treelmesh_type), intent(out) :: me !> Directory containing the mesh informations type(flu_State) :: conf !> Partition to use on the calling process (= MPI Rank in comm) integer, intent(in) :: myPart !> Number of partitions, the mesh is partitioned into (= Number of MPI !! processes in comm). integer, intent(in) :: nParts !> MPI Communicator to use !! now defaults to the one given in the tree%global%comm integer, intent(in) :: comm !> Balancing weight for elements on different levels. !! If these weights are present, the mesh is partitioned initially according !! to the given weights. !! This can not be combined with given offsets. If offset is also provided !! they will overwrite the levelWeight distribution. real(kind=rk), intent(in), optional :: levelWeight(globalMaxLevels) !> output Mesh directory name character(len=*), intent(out), optional :: meshDir !> optional parent handle integer, intent(in), optional :: parent ! -------------------------------------------------------------------- ! integer :: tem_handle character(len=pathLen) :: dirname character(len=4) :: EndianSuffix integer :: mesh_error, weights_error ! The communicator which is used integer :: commLocal real(kind=rk), allocatable :: weights(:) integer(kind=long_k) :: chunksize, remainder logical :: found_weights type( tem_sparta_type ) :: sparta ! -------------------------------------------------------------------- ! ! Use the incoming communicator commLocal = comm call tem_horizontalSpacer(fUnit = logUnit(1)) EndianSuffix = tem_create_EndianSuffix() ! Make sure the dirname has some sane value. dirname = '' call aot_table_open(L = conf, & & parent = parent, & & thandle = tem_handle, & & key = 'mesh' ) if (tem_handle /= 0) then ! key mesh is a table ! The mesh is actually given as a table, parse it and ! generate a mesh internally. write(logUnit(1),*) 'Generating an internally defined mesh' call tem_load_internal(me, conf, tem_handle, myPart, nParts, commLocal) else ! tem_handle == 0, mesh is not a table ! The mesh is not a table, try to interpret it as a string. call aot_get_val( L = conf, & & thandle = parent, & & key = 'mesh', & & val = dirname, & & ErrCode = mesh_error, & & default = 'mesh/' ) if (btest(mesh_error, aoterr_Fatal)) then write(logUnit(0),*) 'FATAL Error occured, while retrieving mesh folder!' if (btest(mesh_error, aoterr_NonExistent))then write(logUnit(0),*) 'Variable not existent!' write(logUnit(0),*) 'Using the default value mesh folder: mesh/' else if (btest(mesh_error, aoterr_WrongType))then write(logUnit(0),*) 'Variable has wrong type!' write(logUnit(0),*) 'STOPPING' call tem_abort() end if end if ! Load global information ! @todo: this routine may already be called outside if (associated(me%global%property)) deallocate(me%global%property) call load_tem_global(me%global, dirname, myPart, nParts, commLocal) ! Homogeneous distribution chunksize = me%global%nElems / nParts remainder = mod(me%global%nElems, int(nParts, kind=long_k)) me%nElems = int(chunksize) if (myPart < int(remainder)) then me%nElems = me%nElems + 1 end if me%elemOffset = chunksize * int(myPart, kind=long_k) & & + min(int(myPart, kind=long_k), remainder) call load_treelmesh( me = me, & & nParts = nParts ) end if ! If requested, provide the directory name to the caller. if ( present(meshDir) ) then meshDir = trim( dirname ) end if ! Make sure each process has at least one element! if (me%global%nElems/me%global%nParts == 0) then write(logunit(1),*) '' write(logunit(1),*) 'FATAL Error occured:' write(logunit(1),*) 'Less elements then processes, this is not' write(logunit(1),*) 'supported! And you should reconfigure your' write(logunit(1),*) 'setup.' write(logunit(1),*) '' write(logunit(1),*) 'Restart this simulation with at most ', & & me%global%nElems, ' processes!' call tem_abort() end if me%levelWeight = 0.0_rk if (present(levelWeight)) then me%levelWeight = levelWeight end if ! Get the name for the weights file. ! This defaults to the mesh prefix with weights appended. call aot_get_val( L = conf, & & key = 'weights', & & thandle = parent, & & val = me%weights_file, & & ErrCode = weights_error, & & default = '' ) if (.not. btest(weights_error, aoterr_NonExistent)) then if (myPart == 0 .and. trim(me%weights_file) /= '') then inquire(file=trim(me%weights_file)//EndianSuffix, exist=found_weights) if (.not. found_weights) then write(logUnit(1),*) '!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!' write(logUnit(1),*) 'Weight file ' & & // trim(me%weights_file)//EndianSuffix & & // ' can NOT be found!' write(logUnit(1),*) 'May proceed WITHOUT BALANCING!' write(logUnit(1),*) '!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!' end if end if end if call aot_get_val( L = conf, & & key = 'write_weights', & & thandle = parent, & & val = me%write_weights, & & ErrCode = weights_error, & & default = me%weights_file ) if (trim(me%write_weights) /= '') then write(logunit(1),*) 'The computational weight of each element will be' write(logunit(1),*) 'written at the end of the simulation into file' write(logunit(1),*) trim(me%write_weights)//EndianSuffix write(logunit(1),*) 'If this is not desired, set write_weights explicitly' write(logunit(1),*) 'to an empty string.' end if ! Do sparta balancing, but only if there are actually multiple processes. if (me%global%nParts > 1) then allocate( weights(me%nElems) ) call tem_load_weights( me = me, & & weights = weights, & & success = found_weights ) if ( found_weights ) then call tem_init_sparta( sparta, nParts ) call tem_balance_Sparta(weight = weights, & & myPart = myPart, & & nParts = nParts, & & comm = commLocal, & & myElems = me%nElems, & & sparta = sparta, & & offset = me%elemOffset ) call exchange_elements( me, sparta ) call tem_destroy_sparta( sparta ) end if deallocate( weights ) end if ! Output basic mesh infos to the screen call dump_meshHeader( me ) call aot_table_close(L=conf, thandle=tem_handle) call tem_horizontalSpacer(fUnit = logUnit(1)) end subroutine load_tem