more flexible reporting

This commit is contained in:
Martin Diehl 2022-05-27 09:29:42 +02:00
parent cfa2be44e9
commit 6c6b3e64b3
4 changed files with 35 additions and 28 deletions

View File

@ -129,7 +129,7 @@ function IO_read(fileName) result(fileContent)
inquire(file = fileName, size=fileLength)
open(newunit=fileUnit, file=fileName, access='stream',&
status='old', position='rewind', action='read',iostat=myStat)
if (myStat /= 0) call IO_error(100,ext_msg=trim(fileName))
if (myStat /= 0) call IO_error(100,trim(fileName))
allocate(character(len=fileLength)::fileContent)
if (fileLength==0) then
close(fileUnit)
@ -137,7 +137,7 @@ function IO_read(fileName) result(fileContent)
endif
read(fileUnit,iostat=myStat) fileContent
if (myStat /= 0) call IO_error(102,ext_msg=trim(fileName))
if (myStat /= 0) call IO_error(102,trim(fileName))
close(fileUnit)
if (scan(fileContent(:index(fileContent,LF)),CR//LF) /= 0) fileContent = CRLF2LF(fileContent)
@ -206,7 +206,7 @@ function IO_stringValue(string,chunkPos,myChunk)
validChunk: if (myChunk > chunkPos(1) .or. myChunk < 1) then
IO_stringValue = ''
call IO_error(110,el=myChunk,ext_msg='IO_stringValue: "'//trim(string)//'"')
call IO_error(110,'IO_stringValue: "'//trim(string)//'"',label1='chunk',ID1=myChunk)
else validChunk
IO_stringValue = string(chunkPos(myChunk*2):chunkPos(myChunk*2+1))
endif validChunk
@ -303,10 +303,10 @@ integer function IO_stringAsInt(string)
valid: if (verify(string,VALIDCHARS) == 0) then
read(string,*,iostat=readStatus) IO_stringAsInt
if (readStatus /= 0) call IO_error(111,ext_msg=string)
if (readStatus /= 0) call IO_error(111,string)
else valid
IO_stringAsInt = 0
call IO_error(111,ext_msg=string)
call IO_error(111,string)
endif valid
end function IO_stringAsInt
@ -325,10 +325,10 @@ real(pReal) function IO_stringAsFloat(string)
valid: if (verify(string,VALIDCHARS) == 0) then
read(string,*,iostat=readStatus) IO_stringAsFloat
if (readStatus /= 0) call IO_error(112,ext_msg=string)
if (readStatus /= 0) call IO_error(112,string)
else valid
IO_stringAsFloat = 0.0_pReal
call IO_error(112,ext_msg=string)
call IO_error(112,string)
endif valid
end function IO_stringAsFloat
@ -348,32 +348,30 @@ logical function IO_stringAsBool(string)
IO_stringAsBool = .false.
else
IO_stringAsBool = .false.
call IO_error(113,ext_msg=string)
call IO_error(113,string)
endif
end function IO_stringAsBool
!--------------------------------------------------------------------------------------------------
!> @brief Write error statements to standard out and terminate the run with exit #9xxx
!> @brief Write error statements to standard out and terminate the run with exit #9xxx.
!--------------------------------------------------------------------------------------------------
subroutine IO_error(error_ID,el,ext_msg)
subroutine IO_error(error_ID,ext_msg,label1,ID1,label2,ID2)
integer, intent(in) :: error_ID
integer, optional, intent(in) :: el
character(len=*), optional, intent(in) :: ext_msg
character(len=*), optional, intent(in) :: ext_msg,label1,label2
integer, optional, intent(in) :: ID1,ID2
external :: quit
character(len=:), allocatable :: msg
character(len=pStringLen) :: formatString
select case (error_ID)
if (present(ID1) .and. .not. present(label1)) error stop 'error value without label (1)'
if (present(ID2) .and. .not. present(label2)) error stop 'error value without label (2)'
!--------------------------------------------------------------------------------------------------
! internal errors
case (0)
msg = 'internal check failed:'
select case (error_ID)
!--------------------------------------------------------------------------------------------------
! file handling errors
@ -558,8 +556,16 @@ subroutine IO_error(error_ID,el,ext_msg)
max(1,72-len_trim(ext_msg)-4),'x,a)'
write(IO_STDERR,formatString) '│ ',trim(ext_msg), '│'
endif
if (present(el)) &
write(IO_STDERR,'(a19,1x,i9,44x,a3)') ' │ at element ',el, '│'
if (present(label1)) then
write(formatString,'(a,i6.6,a,i6.6,a)') '(1x,a7,a',max(1,len_trim(label1)),',i9,',&
max(1,72-len_trim(label1)-9-7),'x,a)'
write(IO_STDERR,formatString) '│ at ',trim(label1),ID1, '│'
endif
if (present(label2)) then
write(formatString,'(a,i6.6,a,i6.6,a)') '(1x,a7,a',max(1,len_trim(label2)),',i9,',&
max(1,72-len_trim(label2)-9-7),'x,a)'
write(IO_STDERR,formatString) '│ at ',trim(label2),ID2, '│'
endif
write(IO_STDERR,'(a,69x,a)') ' │', '│'
write(IO_STDERR,'(a)') ' └'//DIVIDER//'┘'
flush(IO_STDERR)

View File

@ -80,13 +80,13 @@ subroutine discretization_Marc_init
num_commercialFEM => config_numerics%get('commercialFEM',defaultVal = emptyDict)
mesh_unitlength = num_commercialFEM%get_asFloat('unitlength',defaultVal=1.0_pReal) ! set physical extent of a length unit in mesh
if (mesh_unitlength <= 0.0_pReal) call IO_error(301,ext_msg='unitlength')
if (mesh_unitlength <= 0.0_pReal) call IO_error(301,'unitlength')
call inputRead(elem,node0_elem,connectivity_elem,materialAt)
nElems = size(connectivity_elem,2)
if (debug_e < 1 .or. debug_e > nElems) call IO_error(602,ext_msg='element')
if (debug_i < 1 .or. debug_i > elem%nIPs) call IO_error(602,ext_msg='IP')
if (debug_e < 1 .or. debug_e > nElems) call IO_error(602,'element')
if (debug_i < 1 .or. debug_i > elem%nIPs) call IO_error(602,'IP')
allocate(cellNodeDefinition(elem%nNodes-1))
allocate(connectivity_cell(elem%NcellNodesPerCell,elem%nIPs,nElems))
@ -579,7 +579,7 @@ subroutine inputRead_elemType(elem, &
character(len=*), dimension(:), intent(in) :: fileContent !< file content, separated per lines
integer, allocatable, dimension(:) :: chunkPos
integer :: i,j,t,l,remainingChunks
integer :: i,j,t,t_,l,remainingChunks
t = -1
@ -594,7 +594,8 @@ subroutine inputRead_elemType(elem, &
t = mapElemtype(IO_stringValue(fileContent(l+1+i+j),chunkPos,2))
call elem%init(t)
else
if (t /= mapElemtype(IO_stringValue(fileContent(l+1+i+j),chunkPos,2))) call IO_error(191,el=t)
t_ = mapElemtype(IO_stringValue(fileContent(l+1+i+j),chunkPos,2))
if (t /= t_) call IO_error(191,IO_stringValue(fileContent(l+1+i+j),chunkPos,2),label1='type',ID1=t)
endif
remainingChunks = elem%nNodes - (chunkPos(1) - 2)
do while(remainingChunks > 0)
@ -617,7 +618,7 @@ subroutine inputRead_elemType(elem, &
character(len=*), intent(in) :: what
select case (IO_lc(what))
select case (what)
case ( '6')
mapElemtype = 1 ! Two-dimensional Plane Strain Triangle
case ( '125') ! 155, 128 (need test)
@ -645,7 +646,7 @@ subroutine inputRead_elemType(elem, &
case ( '21')
mapElemtype = 13 ! Three-dimensional Arbitrarily Distorted quadratic hexahedral
case default
call IO_error(error_ID=190,ext_msg=IO_lc(what))
call IO_error(190,what)
end select
end function mapElemtype

View File

@ -290,7 +290,7 @@ program DAMASK_grid
if (loadCases(l)%f_restart < huge(0)) &
print'(2x,a,1x,i0)', 'f_restart:', loadCases(l)%f_restart
if (errorID > 0) call IO_error(error_ID = errorID, el = l)
if (errorID > 0) call IO_error(errorID,label1='line',ID1=l)
endif reportAndCheck
enddo

View File

@ -576,7 +576,7 @@ function lattice_C66_trans(Ntrans,C_parent66,lattice_target, &
do i = 1,6
if (abs(C_target_unrotated66(i,i))<tol_math_check) &
call IO_error(135,el=i,ext_msg='matrix diagonal "el"ement in transformation')
call IO_error(135,'matrix diagonal in transformation',label1='entry',ID1=i)
enddo
call buildTransformationSystem(Q,S,Ntrans,cOverA_trans,a_cF,a_cI)