changes towards supporting new Marc2017 input file format

still not working
This commit is contained in:
Franz Roters 2018-01-10 17:13:25 +01:00
parent 515056e9ff
commit 7149f9599f
2 changed files with 637 additions and 525 deletions

View File

@ -50,6 +50,7 @@ module IO
IO_skipChunks, &
IO_extractValue, &
IO_countDataLines, &
IO_countNumericalDataLines, &
IO_countContinuousIntValues, &
IO_continuousIntValues, &
IO_error, &
@ -1247,6 +1248,37 @@ integer(pInt) function IO_countDataLines(fileUnit)
end function IO_countDataLines
!--------------------------------------------------------------------------------------------------
!> @brief count lines containig data up to next *keyword
!--------------------------------------------------------------------------------------------------
integer(pInt) function IO_countNumericalDataLines(fileUnit)
implicit none
integer(pInt), intent(in) :: fileUnit !< file handle
integer(pInt), allocatable, dimension(:) :: chunkPos
character(len=65536) :: line, &
tmp
IO_countNumericalDataLines = 0_pInt
line = ''
do while (trim(line) /= IO_EOF)
line = IO_read(fileUnit)
chunkPos = IO_stringPos(line)
tmp = IO_lc(IO_stringValue(line,chunkPos,1_pInt))
if (scan(tmp,"abcdefghijklmnopqrstuvwxyz")/=0) then ! found keyword
line = IO_read(fileUnit, .true.) ! reset IO_read
exit
else
IO_countNumericalDataLines = IO_countNumericalDataLines + 1_pInt
endif
enddo
backspace(fileUnit)
end function IO_countNumericalDataLines
!--------------------------------------------------------------------------------------------------
!> @brief count items in consecutive lines depending on lines
!> @details Marc: ints concatenated by "c" as last char or range of values a "to" b
@ -1313,7 +1345,7 @@ end function IO_countContinuousIntValues
!--------------------------------------------------------------------------------------------------
!> @brief return integer list corrsponding to items in consecutive lines.
!> @brief return integer list corresponding to items in consecutive lines.
!! First integer in array is counter
!> @details Marc: ints concatenated by "c" as last char, range of a "to" b, or named set
!! Abaqus: triplet of start,stop,inc or named set

View File

@ -73,8 +73,10 @@ module mesh
#ifdef Marc4DAMASK
integer(pInt), private :: &
MarcVersion, & !< Version of input file format (Marc only)
hypoelasticTableStyle, & !< Table style (Marc only)
initialcondTableStyle !< Table style (Marc only)
initialcondTableStyle, & !< Table style (Marc only)
Marc_matNumber !< Material number of hypoelastic material (Marc only)
#endif
integer(pInt), dimension(2), private :: &
@ -428,7 +430,9 @@ module mesh
mesh_spectral_build_elements, &
mesh_spectral_build_ipNeighborhood, &
#elif defined Marc4DAMASK
mesh_marc_get_fileFormat, &
mesh_marc_get_tableStyles, &
mesh_marc_get_matNumber, &
mesh_marc_count_nodesAndElements, &
mesh_marc_count_elementSets, &
mesh_marc_map_elementSets, &
@ -579,8 +583,14 @@ subroutine mesh_init(ip,el)
#elif defined Marc4DAMASK
call IO_open_inputFile(FILEUNIT,modelName) ! parse info from input file...
if (myDebug) write(6,'(a)') ' Opened input file'; flush(6)
call mesh_marc_get_fileFormat(FILEUNIT)
if (myDebug) write(6,'(a)') ' Got input file format'; flush(6)
call mesh_marc_get_tableStyles(FILEUNIT)
if (myDebug) write(6,'(a)') ' Got table styles'; flush(6)
if (MarcVersion > 12) then
call mesh_marc_get_matNumber(FILEUNIT)
if (myDebug) write(6,'(a)') ' Got hypoleastic material number'; flush(6)
endif
call mesh_marc_count_nodesAndElements(FILEUNIT)
if (myDebug) write(6,'(a)') ' Counted nodes/elements'; flush(6)
call mesh_marc_count_elementSets(FILEUNIT)
@ -1534,6 +1544,37 @@ end function mesh_nodesAroundCentres
#endif
#ifdef Marc4DAMASK
!--------------------------------------------------------------------------------------------------
!> @brief Figures out version of Marc input file format and stores ist as MarcVersion
!--------------------------------------------------------------------------------------------------
subroutine mesh_marc_get_fileFormat(fileUnit)
use IO, only: &
IO_lc, &
IO_intValue, &
IO_stringValue, &
IO_stringPos
implicit none
integer(pInt), intent(in) :: fileUnit
integer(pInt), allocatable, dimension(:) :: chunkPos
character(len=300) line
610 FORMAT(A300)
rewind(fileUnit)
do
read (fileUnit,610,END=620) line
chunkPos = IO_stringPos(line)
if ( IO_lc(IO_stringValue(line,chunkPos,1_pInt)) == 'version') then
MarcVersion = IO_intValue(line,chunkPos,2_pInt)
exit
endif
enddo
620 end subroutine mesh_marc_get_fileFormat
!--------------------------------------------------------------------------------------------------
!> @brief Figures out table styles (Marc only) and stores to 'initialcondTableStyle' and
!! 'hypoelasticTableStyle'
@ -1570,6 +1611,41 @@ subroutine mesh_marc_get_tableStyles(fileUnit)
620 end subroutine mesh_marc_get_tableStyles
!--------------------------------------------------------------------------------------------------
!> @brief Figures out material number of hypoelastic material and sores it as Marc_matNumber
!--------------------------------------------------------------------------------------------------
subroutine mesh_marc_get_matNumber(fileUnit)
use IO, only: &
IO_lc, &
IO_intValue, &
IO_stringValue, &
IO_stringPos
implicit none
integer(pInt), intent(in) :: fileUnit
integer(pInt), allocatable, dimension(:) :: chunkPos
integer(pInt) :: i
character(len=300) line
610 FORMAT(A300)
rewind(fileUnit)
do
read (fileUnit,610,END=620) line
chunkPos = IO_stringPos(line)
if ( IO_lc(IO_stringValue(line,chunkPos,1_pInt)) == 'hypoelastic') then
do i=1_pInt,1_pInt+hypoelasticTableStyle ! Skip 1 or 2 lines
read (fileUnit,610,END=620) line
enddo
Marc_matNumber = IO_intValue(line,chunkPos,1_pInt)
exit
endif
enddo
620 end subroutine mesh_marc_get_matNumber
!--------------------------------------------------------------------------------------------------
!> @brief Count overall number of nodes and elements in mesh and stores the numbers in
@ -1699,13 +1775,14 @@ subroutine mesh_marc_count_cpElements(fileUnit)
IO_stringPos, &
IO_countContinuousIntValues, &
IO_error, &
IO_intValue
IO_intValue, &
IO_countNumericalDataLines
implicit none
integer(pInt), intent(in) :: fileUnit
integer(pInt), allocatable, dimension(:) :: chunkPos
integer(pInt) :: i, version
integer(pInt) :: i
character(len=300):: line
mesh_NcpElems = 0_pInt
@ -1713,13 +1790,7 @@ subroutine mesh_marc_count_cpElements(fileUnit)
610 FORMAT(A300)
rewind(fileUnit)
do
read (fileUnit,610,END=620) line
chunkPos = IO_stringPos(line)
if ( IO_lc(IO_stringValue(line,chunkPos,1_pInt)) == 'version') then
version = IO_intValue(line,chunkPos,2_pInt)
if (version < 13) then ! Marc 2016 or earlier
rewind(fileUnit)
if (MarcVersion < 13) then ! Marc 2016 or earlier
do
read (fileUnit,610,END=620) line
chunkPos = IO_stringPos(line)
@ -1727,15 +1798,25 @@ subroutine mesh_marc_count_cpElements(fileUnit)
do i=1_pInt,3_pInt+hypoelasticTableStyle ! Skip 3 or 4 lines
read (fileUnit,610,END=620) line
enddo
mesh_NcpElems = mesh_NcpElems + IO_countContinuousIntValues(fileUnit) ! why not simply mesh_NcpElems = IO_countContinuousIntValues(fileUnit)? keyword hypoelastic might appear several times
mesh_NcpElems = mesh_NcpElems + IO_countContinuousIntValues(fileUnit) ! why not simply mesh_NcpElems = IO_countContinuousIntValues(fileUnit)? keyword hypoelastic might appear several times so the next line probably should not be there
exit
endif
enddo
else ! Marc2017 and later
call IO_error(error_ID=701_pInt)
end if
end if
do
read (fileUnit,610,END=620) line
chunkPos = IO_stringPos(line)
if ( IO_lc(IO_stringValue(line,chunkPos,1_pInt)) == 'connectivity') then
read (fileUnit,610,END=620) line
chunkPos = IO_stringPos(line)
if (IO_intValue(line,chunkPos,6_pInt)==Marc_matNumber) then
mesh_NcpElems = mesh_NcpElems + IO_countNumericalDataLines(fileUnit)
exit
endif
endif
enddo
end if
620 end subroutine mesh_marc_count_cpElements
@ -2938,7 +3019,6 @@ subroutine mesh_build_sharedElems
allocate(node_seen(maxval(FE_NmatchingNodes)))
node_count = 0_pInt
do e = 1_pInt,mesh_NcpElems