Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | info |
Informaion flag |
||
character(len=*), | intent(in) | :: | origin |
The name of the subroutine from which the flag originates |
||
character(len=*), | intent(in), | optional | :: | module |
The name of the module in which the call happens |
|
character(len=*), | intent(in), | optional | :: | procedure |
The name of the procedure in which the call happens |
|
character(len=*), | intent(in), | optional | :: | info_msg |
subroutine check_info(info, origin, module, procedure, info_msg) integer, intent(in) :: info !! Informaion flag character(len=*), intent(in) :: origin !! The name of the subroutine from which the flag originates character(len=*), optional, intent(in) :: module !! The name of the module in which the call happens character(len=*), optional, intent(in) :: procedure !! The name of the procedure in which the call happens character(len=*), optional, intent(in) :: info_msg character(len=128) :: str !! Optional extra message ! internals character(len=256) :: msg integer :: ierr str = optval(info_msg, '') ierr = 0 if (info == 0) then ! Successful exit --> only log on debug write (msg, '(A)') 'The subroutine "'//trim(origin)//'" returned successfully. '//trim(str) call log_debug(trim(msg), module=module, procedure=procedure) else ! ! LAPACK ! if (trim(to_lower(origin)) == 'getref') then ! GETREF if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A,I0,A,I0,A)') "U(", info, ",", info, ") is exactly zero. The factorization ", & & "has been completed but the factor U is exactly singular. ", & & "Division by zero will occur if used to solve Ax=b. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'getri') then ! GETRI if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A,I0,A)') "U(", info, ",", info, ") is exactly zero. ", & & "The matrix is singular and its inverse cannot be computed. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'geev') then ! GEEV if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A,I0,A)') "The QR alg. failed to compute all of the eigenvalues.", & & "No eigenvector has been computed. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'syev') then ! SYEV if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A)') "The QR alg. failed to compute all of the eigenvalues.", & & "No eigenvector has been computed. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'heev') then ! HEEV if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A)') "The QR alg. failed to compute all of the eigenvalues.", & & "No eigenvector has been computed. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'gels') then ! GELS if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'gees') then ! GEES if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A)') "The QR alg. failed to compute all of the eigenvalues.", & & "No eigenvector has been computed. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'trsen') then ! GEES if (info < 0) then write (msg, '(A,I0,A)') "The ", -info, "-th argument has illegal value. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else if (info == 1) then write (msg, '(A)') "The reordering of T failed because some eigenvalues are too", & & "close to separate (the problem is very ill-conditioned); ", & & "T may have been partially reordered, and WR and WI ", & & "contain the eigenvalues in the same order as in T; S and", & & "SEP (if requested) are set to zero. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if ! ! LightKrylov_Utils ! else if (trim(to_lower(origin)) == 'sqrtm') then if (info == 1) then write (msg, '(A)') 'The input matrix is singular to tolerance. The singular eigenvalues are set to zero. ' call log_warning(trim(msg), module=module, procedure=procedure) else if (info == -1) then write (msg, '(A)') "The input matrix is not positive (semi-)definite. " call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if ! ! LightKrylov_BaseKrylov ! else if (trim(to_lower(origin)) == 'orthogonalize_against_basis') then ! the regular case ! orthogonalization if (info > 0) then write (msg, '(A,I0,A)') 'Orthogonalization: The ', info, 'th input vector is numerically zero.' call log_debug(trim(msg), module=module, procedure=procedure) else if (info == -1) then write (msg, '(A)') 'The input Krylov basis is not orthonormal.' call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else if (info == -2) then write (msg, '(A)') 'Orthogonalization: The last column of the input basis is zero.' call log_debug(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'orthogonalize_against_basis_p1') then ! orthogonalization if (info > 0) then write (msg, '(A,I0,A)') 'Orthogonalization: The ', info, 'th input vector is numerically zero.' call log_debug(trim(msg), module=module, procedure=procedure) else if (info == -1) then write (msg, '(A)') 'The input Krylov basis is not orthonormal.' call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else if (info == -2) then write (msg, '(A)') 'Orthogonalization: The last column of the input basis is zero.' call log_warning(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'orthogonalize_against_basis_p2') then ! orthogonalization if (info > 0) then ! show this information only for debugging write (msg, '(A,I0,A)') 'Orthogonalization: The ', info, 'th input vector is numerically zero.' call log_debug(trim(msg), module=module, procedure=procedure) else if (info == -1) then write (msg, '(A)') 'The input Krylov basis is not orthonormal.' call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else if (info == -2) then write (msg, '(A)') 'Orthogonalization: The last column of the input basis is zero.' call log_warning(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'double_gram_schmidt_step') then ! orthogonalization if (info > 0) then write (msg, '(A,I0,A)') 'Orthogonalization: The ', info, 'th input vector is numerically zero.' call log_debug(trim(msg), module=module, procedure=procedure) else if (info == -1) then write (msg, '(A)') 'The input Krylov basis is not orthonormal.' call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 else if (info == -2) then write (msg, '(A)') 'Orthogonalization: The last column of the input basis is zero.' call log_warning(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'qr') then ! qr if (info > 0) then write (msg, '(A,I0,A)') 'QR factorization: Colinear column detected in column ', info, & & '. NOTE: Other subsequent columns may also be colinear.' call log_debug(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'qr_pivot') then ! qr_pivot if (info > 0) then write (msg, '(A,I0,A)') 'QR factorization: Invariant subspace found after ', info, ' steps.' call log_debug(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'arnoldi') then ! arnoldi if (info > 0) then write (msg, '(A,I0,A)') 'Arnoldi factorization: Invariant subspace computed after ', info, ' iterations.' call log_debug(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'bidiagonalization') then ! lanczos_bidiagonalization if (info > 0) then write (msg, '(A,I0,A)') 'Lanczos Bidiagonalisation: Invariant subspace found after ', info, ' steps.' call log_debug(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'lanczos') then ! lanczos_tridiagonalization if (info > 0) then write (msg, '(A,I0,A)') 'Lanczos Tridiagonalisation: Invariant subspace found after ', info, ' steps.' call log_debug(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if ! ! LightKrylov_IterativeSolvers ! else if (trim(to_lower(origin)) == 'eigs') then ! GMRES if (info > 0) then write (msg, '(A,I0,A)') 'eigs iteration converged after ', info, ' iterations' call log_information(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'eighs') then ! GMRES if (info > 0) then write (msg, '(A,I0,A)') 'eigs iteration converged after ', info, ' iterations' call log_information(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'svds') then ! GMRES if (info > 0) then write (msg, '(A,I0,A)') 'svds iteration converged after ', info, ' iterations' call log_information(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'gmres') then ! GMRES if (info > 0) then write (msg, '(A,I0,A)') 'GMRES iteration converged after ', info, ' iterations' call log_message(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'fgmres') then ! GMRES if (info > 0) then write (msg, '(A,I0,A)') 'FGMRES iteration converged after ', info, ' iterations' call log_message(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'cg') then ! CG if (info > 0) then write (msg, '(A,I0,A)') 'CG iteration converged after ', info, ' iterations' call log_message(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if else if (trim(to_lower(origin)) == 'linear_solver') then ! Abstract linear solver if (info > 0) then write (msg, '(A,I0,A)') 'The linear solver converged after ', info, ' iterations' call log_message(trim(msg), module=module, procedure=procedure) else write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if ! ! LightKrylov_ExpmLib ! else if (trim(to_lower(origin)) == 'kexpm') then ! Krylov Matrix Exponential if (info > 0) then write (msg, '(A,I0,A)') 'kexpm converged. Estimated error below tolerance using ', info, ' Krylov vectors.' call log_debug(trim(msg), module=module, procedure=procedure) else if (info == -2) then write (msg, '(A)') 'kexpm converged. Arnoldi iteration breakdown. Approximation is exact to arnoldi tolerance.' call log_debug(trim(msg), module=module, procedure=procedure) else if (info == -1) then write (msg, '(A)') 'kexpm did not converge. Maximum number of Krylov vectors reached.' call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 write (msg, '(A)') "Undocumented error. "//trim(str) call log_error(origin, module=module, procedure=procedure, stat=info, errmsg=trim(msg)) ierr = -1 end if ! ! stop error ! else if (trim(origin) == 'STOP_ERROR') then call log_error(trim(origin), module=module, procedure=procedure, stat=info, errmsg=trim(str)) ierr = -1 ! ! Default ! else write (msg, '(A)') 'subroutine "'//trim(origin)//'" returned with a non-zero error flag.' call log_error(trim(msg), module=module, procedure=procedure, stat=info, errmsg=trim(str)) ierr = -1 end if end if ! info /= 0 call error_handler(ierr) end subroutine check_info