stop.F90 7.16 KB
Newer Older
1 2 3 4 5 6
!--------------------------------------------------------------------------------
! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
! This file is part of FLEUR and available as free software under the conditions
! of the MIT license as expressed in the LICENSE file in more detail.
!--------------------------------------------------------------------------------

7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63
MODULE m_juDFT_stop
  !-----------------------------------------------
  !    module to terminate Calculation, should be used instead
  !    of a simple STOP
  !
  !    error(message,calledby,hint,no,warning)
  !         message  : message string
  !         calledby : subroutine in which error occurs(optional)
  !         hint     : string with more information (optional)
  !         no       : error number (optional)
  !         warning  : logical indicating a warning message (optional)
  !
  !    warn(message,calledby,hint,no)
  !         shortcut for calling error with warning=.true.
  !
  !    juDFT_end(message)
  !         call this to terminate without error
  !
  !   IF the file "JUDFT_WARN_ONLY" is not present, warnings will lead to errors.
  !
  !   If the file "JUDFT_TRACE" is present, a stacktrace will be generated
  !   on some compilers
  !
  !
  !                 Daniel Wortmann (2010)
  !-----------------------------------------------
  USE m_judft_time
  USE m_judft_sysinfo
  IMPLICIT NONE
  PRIVATE
  PUBLIC juDFT_error,juDFT_warn,juDFT_end,judft_file_readable
CONTAINS

  SUBROUTINE judfT_file_readable(filename,warning)
    IMPLICIT NONE
    CHARACTER(len=*),INTENT(IN):: filename
    LOGICAL,INTENT(IN),OPTIONAL:: warning
    LOGICAL  :: l_exist

    INQUIRE(file=filename,exist=l_exist)
    IF (.not.l_exist) CALL judft_error("File not readable:"//filename,hint="FLEUR wants to read a file that is not present",warning=warning)
  END SUBROUTINE judfT_file_readable

  SUBROUTINE juDFT_error(message,calledby,hint,no,warning,file,line)
    IMPLICIT NONE
    CHARACTER*(*),INTENT(IN)          :: message
    CHARACTER*(*),OPTIONAL,INTENT(IN) :: calledby,hint
    INTEGER,OPTIONAL,INTENT(IN)       :: no
    LOGICAL,OPTIONAL,INTENT(IN)       :: warning
    CHARACTER*(*),OPTIONAL,INTENT(IN) :: file
    INTEGER,INTENT(IN),OPTIONAL       :: line

    LOGICAL :: callstop,warn
    CHARACTER(LEN=4)::PE
    !store all output in variable for single call to write in MPI case
    CHARACTER(len=300)::text(10)
    INTEGER           ::linenr,n
64
#ifdef CPP_MPI
65 66 67 68
    include 'mpif.h'
    INTEGER :: irank,e
    CALL MPI_COMM_RANK(MPI_COMM_WORLD,irank,e)
    WRITE(PE,"(i4)") irank
69
#else
70
    PE="****"
71
#endif
72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114
    warn = .FALSE.
    IF (PRESENT(warning)) warn = warning
    IF (warn) THEN
       !check if we stop nevertheless
       INQUIRE(FILE ="JUDFT_WARN_ONLY",EXIST= callstop)
       callstop  = .NOT.callstop
    ELSE
       callstop = .TRUE.
    ENDIF

    IF (.NOT.warn) THEN
       WRITE(text(1),*) PE,"**************juDFT-Error*****************"
    ELSE
       WRITE(text(1),*) PE,"************juDFT-Warning*****************"
    ENDIF
    WRITE(text(2),"(3a)") PE,"Error message:",message
    linenr=3
    IF (PRESENT(calledby)) THEN
       WRITE(text(3),"(3a)") PE,"Error occurred in subroutine:",calledby
       linenr=4
    ENDIF
    IF (PRESENT(hint)) THEN
       WRITE(text(linenr),"(3a)") PE,"Hint:",hint
       linenr=linenr+1
    ENDIF
    IF (PRESENT(no)) THEN
       WRITE(text(linenr),"(2a,i0)") PE,"Error number:",no
       linenr=linenr+1
    ENDIF
    IF (present(file)) THEN
       if (present(line)) THEN
          write(text(linenr),"(4a,i0)") PE,"Source:",file,":",line
       ELSE
          write(text(linenr),"(3a)") PE,"Source:",file
       ENDIF
       linenr=linenr+1
    ENDIF
    WRITE(text(linenr),*) PE,"*****************************************"

    if (.not.warn) CALL juDFT_time_lastlocation(PE)

    IF (callstop) THEN
       IF (warn) THEN
115
          linenr=linenr+1
116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140
          WRITE(text(linenr),'(a)')"Warnings not ignored. Touch 'JUDFT_WARN_ONLY' to make the warning nonfatal"
       ENDIF
       write(0,"(10(a,/))") (trim(text(n)),n=1,linenr)
       CALL juDFT_STOP()
    ENDIF
    write(0,"(10(a,/))") (trim(text(n)),n=1,linenr)
  END SUBROUTINE juDFT_error

  SUBROUTINE juDFT_warn(message,calledby,hint,no,file,line)
    IMPLICIT NONE
    CHARACTER*(*),INTENT(IN)          :: message
    CHARACTER*(*),OPTIONAL,INTENT(IN) :: calledby,hint
    INTEGER,OPTIONAL,INTENT(IN)       :: no
    CHARACTER*(*),OPTIONAL,INTENT(IN) :: file
    INTEGER,INTENT(IN),OPTIONAL       :: line

    CALL juDFT_error(message,calledby,hint,no,warning = .TRUE.,file=file,line=line)

  END SUBROUTINE juDFT_warn

  SUBROUTINE juDFT_END(message, irank, l_endXML)
    ! If irank is present every mpi process has to call this routine.
    ! Otherwise only a single mpi process is allowed to call the routine.
    USE m_xmlOutput
    IMPLICIT NONE
141
#ifdef CPP_MPI
142 143
    INCLUDE 'mpif.h'
    INTEGER :: ierr
144
#endif
145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171
    CHARACTER*(*), INTENT(IN)      :: message
    INTEGER, OPTIONAL, INTENT(IN)  :: irank
    LOGICAL, OPTIONAL, INTENT(IN)  :: l_endXML

    LOGICAL l_endXML_local

    l_endXML_local = .TRUE.
    IF(PRESENT(l_endXML)) THEN
       l_endXML_local = l_endXML
    END IF

    IF(l_endXML_local) THEN
       IF(PRESENT(irank)) THEN
          IF (irank.EQ.0) CALL endXMLOutput()
       ELSE
          ! It is assumed that this is the only mpi process calling this routine.
          CALL endXMLOutput()
       END IF
    END IF

    WRITE(0,*) "*****************************************"
    WRITE(0,*) "Run finished successfully"
    WRITE(0,*) "Stop message:"
    WRITE(0,*) "  ",message
    WRITE(0,*) "*****************************************"
    CALL writetimes()
    CALL print_memory_info()
172
#ifdef CPP_MPI
173 174 175 176 177 178
    IF(PRESENT(irank)) THEN
       CALL MPI_BARRIER(MPI_COMM_WORLD,ierr)
       CALL MPI_FINALIZE(ierr)
    ELSE
       CALL juDFT_STOP(0)
    END IF
179
#endif
180 181
    STOP 'OK'
  END SUBROUTINE juDFT_END
182

183 184 185
  !this is a private subroutine that stops the calculations
  !different compilers might have to be added here
  SUBROUTINE juDFT_stop(errorCode)
186
#ifdef __INTEL_COMPILER
187
    USE ifcore
188 189
#endif
#ifdef CPP_MPI
190
    INCLUDE 'mpif.h'
191
#endif
192 193 194 195
    INTEGER, OPTIONAL, INTENT(IN)  :: errorCode
    INTEGER :: error
    LOGICAL :: calltrace
    LOGICAL,ALLOCATABLE:: a(:)
196
#ifdef CPP_MPI
197
    INTEGER :: ierr
198
#endif
199 200 201 202 203 204 205 206 207 208 209 210
    error = 0

    IF(PRESENT(errorCode)) THEN
       error = errorCode
    END IF
    !try to print times
    call writelocation()
    CALL writetimes(.TRUE.)
    CALL print_memory_info()
    INQUIRE(FILE="JUDFT_TRACE",EXIST=calltrace)
    IF (error.EQ.1) calltrace = .TRUE.
    IF (calltrace) THEN
211
#ifdef __INTEL_COMPILER
212
       CALL tracebackqq(USER_EXIT_CODE=-1) !return after traceback
213
#elif (defined(CPP_AIX)&&!defined(__PGI))
214
       CALL xl__trbk()
215
#endif
216 217
       DEALLOCATE(a)!will generate an error that can be found by the compiler
    ENDIF
218 219

#if defined(CPP_MPI)
220 221 222 223 224 225 226 227
    IF(error.EQ.0) THEN
       WRITE(0,*) ""
       WRITE(0,*) "Terminating all MPI processes."
       WRITE(0,*) "Note: This is a normal procedure."
       WRITE(0,*) "      Error messages in the following lines can be ignored."
       WRITE(0,*) ""
    END IF
    CALL MPI_ABORT(MPI_COMM_WORLD,error,ierr)
228
#endif
229 230
    STOP 'juDFT-STOPPED'
  END SUBROUTINE juDFT_stop
231 232


233
END MODULE m_juDFT_stop
234 235