Read in an arbitrary shape from a lua file defined in a single table
read a shape like for example inside a tracking table
tracking = { variable = { 'velocity' },
shape = { kind = 'canoND',
object = { origin = { 1.0, 1.0, 1.0 },
vec = { 2.0, 2.0, 2.0 },
segments = { 10, 20, 30 } }
}
} -- tracking table
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_shape_type), | intent(out) | :: | me |
shape type defined in a lua file |
||
type(flu_State) | :: | conf |
lua state |
|||
character(len=*), | intent(in), | optional | :: | key |
optional key to load from |
|
integer, | intent(in), | optional | :: | parent |
optional parent handle |
|
integer, | intent(in), | optional | :: | sub_table |
shape table handle |
|
integer, | intent(out), | optional | :: | iError |
error flag |
|
logical, | intent(in), | optional | :: | reqSegments |
Is true if use_get_point is true in output table |
subroutine tem_load_shape_single( me, conf, key, parent, sub_table, iError, & & reqSegments ) !--------------------------------------------------------------------------- !> shape type defined in a lua file type(tem_shape_type), intent(out) :: me !> lua state type(flu_state) :: conf !> optional key to load from character(len=*), optional, intent(in) :: key !> optional parent handle integer, optional, intent(in) :: parent !> shape table handle integer, optional, intent(in) :: sub_table !> error flag integer, optional, intent(out) :: iError !> Is true if use_get_point is true in output table logical, optional, intent(in) :: reqSegments !--------------------------------------------------------------------------- character(len=PathLen) :: buffer integer :: iError_loc integer :: shape_table character(len=32) :: localKey type(tem_transformation_type) :: transform !--------------------------------------------------------------------------- if( present( iError )) iError = 0 if( present( key )) then localKey = key else localKey = 'shape' endif ! open the table ! shape = {} if( present( sub_table )) then shape_table = sub_table else call aot_table_open( L = conf, & & thandle = shape_table, & & parent = parent, & & key = trim(localKey) ) endif if (shape_table > 0 ) then ! load transformation table if defined and pass the transformation ! to geom table to apply the transformation to geometry object call tem_load_transformation(transform = transform, conf = conf, & & thandle = shape_table) ! Get the entry 'kind' ! shape = {{kind='test'}} call aot_get_val( L = conf, & & thandle = shape_table, & & val = buffer, & & ErrCode = iError_loc, & & key = 'kind', & & default = 'none' ) me%kind = trim(buffer) write(logUnit(1),*)' Loading shape kind: '//trim(me%kind) ! Read in the type and configuration of tracking shape ! (point, line, plane, ...) select case( trim(buffer) ) case('global', 'all', 'globalmesh') me%shapeID = tem_global_shape case('canonicalND', 'canoND') me%shapeID = tem_geometrical_shape me%kind = 'canoND' call tem_load_canonicalND( me = me%canoND, & & conf = conf, & & thandle = shape_table, & & transform = transform, & & reqSegments = reqSegments ) case('triangle') me%shapeID = tem_geometrical_shape call tem_load_triangle( me = me%triangle, & & conf = conf, & & thandle = shape_table, & & transform = transform ) case('stl') me%shapeID = tem_geometrical_shape call tem_load_stl( stl_data = me%stl_data, & & conf = conf, & & thandle = shape_table, & & transform = transform ) case('sphere') me%shapeID = tem_geometrical_shape call tem_load_sphere( me = me%sphere, & & conf = conf, & & thandle = shape_table, & & transform = transform ) case('ellipsoid') me%shapeID = tem_geometrical_shape call tem_load_ellipsoid( me = me%ellipsoid, & & conf = conf, & & thandle = shape_table, & & transform = transform ) case('cylinder') me%shapeID = tem_geometrical_shape call tem_load_cylinder( me = me%cylinder, & & conf = conf, & & thandle = shape_table, & & transform = transform ) case('property') me%shapeID = tem_property_shape call tem_shape_load_propLabel( me%propBits, conf, shape_table ) case( 'level' ) me%shapeID = tem_level_shape call tem_shape_load_level( me%minLevel, me%maxLevel, conf, shape_table ) case('boundary') me%shapeID = tem_boundary_shape call tem_shape_load_bcLabels( me%bcLabels, conf, shape_table ) call aot_get_val( L = conf, & & thandle = shape_table, & & val = me%cutOffQVal, & & ErrCode = iError_loc, & & key = 'cutoff_qvalue', & & default = 1.0_rk ) write(logUnit(1),*)' Cutoff qValue: ', me%cutOffQVal ! If find unknown shape, code will stop! case default write(logUnit(1),*) ' Unknown shape specified: '//trim(buffer) write(logUnit(1),*) ' Choose: ' write(logUnit(1),*) ' shape = { kind = "canoND", ' write(logUnit(1),*) ' object = {origin = {-0.5,0.,0.}, ' write(logUnit(1),*) ' vec = {1.0, 0.0, 0.0},' write(logUnit(1),*) ' segments = 1000 } } or' write(logUnit(1),*) ' shape = { kind = "all" } or' write(logUnit(1),*) ' shape = { kind = "property",' write(logUnit(1),*) ' property = {"boundary"} }' write(logUnit(1),*) ' shape = { kind = "boundary",' write(logUnit(1),*) ' boundary = {"bc1", "bc2"} }' write(logUnit(1),*) ' Stopping ...' call tem_abort() end select call aot_get_val( L = conf, & & thandle = shape_table, & & val = me%inverted, & & ErrCode = iError_loc, & & key = 'inverted', & & default = .false. ) write(logUnit(1),*)' Inverted shape: ',me%inverted if( .not. present( parent )) & & call aot_table_close( L = conf, thandle = shape_table ) else write(logUnit(2),*) 'Shape table is not defined' write(logunit(2),*) '... using global mesh' me%kind = 'all' me%shapeID = tem_global_shape if( present( iError )) iError = ibset( iError, aoterr_NonExistent ) endif end subroutine tem_load_shape_single