analyse_solver Subroutine

public subroutine analyse_solver(mat)

Arguments

Type IntentOptional AttributesName
type(MAT_SOLV), intent(inout), target:: mat

high level system type


Calls

proc~~analyse_solver~~CallsGraph proc~analyse_solver analyse_solver ma48_analyse ma48_analyse proc~analyse_solver->ma48_analyse dmumps dmumps proc~analyse_solver->dmumps proc~prep_superlu prep_superlu proc~analyse_solver->proc~prep_superlu proc~s_umfpack_di_symbolic s_umfpack_di_symbolic proc~analyse_solver->proc~s_umfpack_di_symbolic interface~dcreate_compcol_matrix dCreate_CompCol_Matrix proc~prep_superlu->interface~dcreate_compcol_matrix interface~dcreate_dense_matrix dCreate_Dense_Matrix proc~prep_superlu->interface~dcreate_dense_matrix proc~umfpack_di_symbolic umfpack_di_symbolic proc~s_umfpack_di_symbolic->proc~umfpack_di_symbolic interface~c_umfpack_di_symbolic c_umfpack_di_symbolic proc~umfpack_di_symbolic->interface~c_umfpack_di_symbolic

Called by

proc~~analyse_solver~~CalledByGraph proc~analyse_solver analyse_solver proc~solve_syst solve_syst proc~solve_syst->proc~analyse_solver program~test_solvers test_solvers program~test_solvers->proc~solve_syst

Contents

Source Code


Source Code

   subroutine analyse_solver(mat)
   implicit none
   type(MAT_SOLV), intent(inout), target :: mat   !! *high level system type*

      select case(mat%slv_t)

         case(MA48)
            mat%matma48%zmat%row => mat%irow
            mat%matma48%zmat%col => mat%jcol
            mat%matma48%zmat%val => mat%a_elt

            mat%matma48%zmat%n   = mat%nn
            mat%matma48%zmat%m   = mat%nn
            mat%matma48%zmat%ne  = mat%nz

            call ma48_analyse(matrix = mat%matma48%zmat, &
                             factors = mat%matma48%fact, &
                             control = mat%matma48%ctrl, &
                               ainfo = mat%matma48%ainf, &
                               finfo = mat%matma48%finf)
            if (mat%matma48%ainf%flag < 0) then
               write(OPU,*) 'Failure of ma48_analyse with ainfop%flag = ', mat%matma48%ainf%flag
               stop
            endif

         case(MUMP)
            mat%matmump%eltptr   => mat%eltptr
            mat%matmump%eltvar   => mat%eltvar
            mat%matmump%a_elt    => mat%a_elt
            mat%matmump%rhs      => mat%b

            mat%matmump%n    = mat%nn
            mat%matmump%nelt = mat%ne

            mat%matmump%job = 1 ! performs the analysis

            call dmumps(mat%matmump)

         case(SULU)
            mat%matsulu%irow   => mat%irow
            mat%matsulu%jptr   => mat%jptr
            mat%matsulu%a_elt  => mat%a_elt
            mat%matsulu%b      => mat%b

            mat%matsulu%n     = mat%nn
            mat%matsulu%nz    = mat%nz
            mat%matsulu%nrhs  = 1
            mat%matsulu%first = .true.

            call prep_superlu(sulu = mat%matsulu)

         case(UMFP)
            call s_umfpack_di_symbolic(n_row    = mat%nn,                  &
                                       n_col    = mat%nn,                  &
                                       Ap       = mat%jptr,                &
                                       Ai       = mat%irow,                &
                                       Ax       = mat%a_elt,               &
                                       Symbolic = mat%matumfp%c_symbolic,  &
                                       Control  = mat%matumfp%c_control,   &
                                       Info     = mat%matumfp%c_info)

         case default
            stop 'ANALYSE_SOLVER : unknown solver type'

      endselect

   return
   endsubroutine analyse_solver