solution_solver Subroutine

public subroutine solution_solver(mat)

Arguments

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

high level system type


Calls

proc~~solution_solver~~CallsGraph proc~solution_solver solution_solver ma48_solve ma48_solve proc~solution_solver->ma48_solve proc~solv_superlu solv_superlu proc~solution_solver->proc~solv_superlu mpi_finalize mpi_finalize proc~solution_solver->mpi_finalize dmumps dmumps proc~solution_solver->dmumps proc~s_umfpack_di_solve s_umfpack_di_solve proc~solution_solver->proc~s_umfpack_di_solve interface~dgssvx dgssvx proc~solv_superlu->interface~dgssvx interface~destroy_supernode_matrix Destroy_SuperNode_Matrix proc~solv_superlu->interface~destroy_supernode_matrix interface~statfree StatFree proc~solv_superlu->interface~statfree proc~prep_superlu prep_superlu proc~solv_superlu->proc~prep_superlu interface~destroy_compcol_matrix Destroy_CompCol_Matrix proc~solv_superlu->interface~destroy_compcol_matrix interface~statinit StatInit proc~solv_superlu->interface~statinit proc~umfpack_di_solve umfpack_di_solve proc~s_umfpack_di_solve->proc~umfpack_di_solve interface~c_umfpack_di_solve c_umfpack_di_solve proc~umfpack_di_solve->interface~c_umfpack_di_solve 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

Called by

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

Contents

Source Code


Source Code

   subroutine solution_solver(mat)
   implicit none
   type(MAT_SOLV), target, intent(inout) :: mat   !! *high level system type*
      integer(kind=I4) :: ierr

      select case(mat%slv_t)

         case(MA48)
            call ma48_solve(matrix = mat%matma48%zmat,   &
                           factors = mat%matma48%fact,   &
                               rhs = mat%b,              &
                                 x = mat%x,              &
                           control = mat%matma48%ctrl,   &
                             sinfo = mat%matma48%sinf,   &
                             resid = mat%matma48%resid,  &
                             error = mat%error)

         case(MUMP)
            mat%matmump%job = 3
            call dmumps(mat%matmump)
            if (mat%matmump%infog(1) < 0) then
               write(OPU,'(a,a,i6,a,i9)') ' error return: ',                              &
                                          '  mumps_par%infog(1)= ', mat%matmump%infog(1), &
                                          '  mumps_par%infog(2)= ', mat%matmump%infog(2)
               call mpi_finalize(ierr)
               stop 'Error in SOLUTION_SOLVER'
            endif
            ! solution has been assembled on the host
            if ( mat%matmump%myid == 0 ) then
               mat%x(1:mat%nn) = mat%matmump%rhs(1:mat%nn)
            endif

         case(SULU)
            call solv_superlu(sol_x = mat%x,        &
                              sulu  = mat%matsulu,  &
                           verbose  = (mat%matsulu%options%PrintStat==1))
            mat%matsulu%first = .false.

         case(UMFP)
            ! Numeric factors must exist
            if (.not.C_ASSOCIATED(mat%matumfp%c_numeric)) call solve_syst(mat, 'fac')
            mat%matumfp%c_control(UMFPACK_IRSTEP) = 0 ! solve ax=b, without iterative refinement

            ! If you want to evaluate the required RAM (Go)
            ! write(*,*) mat%matumfp%c_info(UMFPACK_PEAK_MEMORY_ESTIMATE)/mat%matumfp%c_info(UMFPACK_SIZE_OF_UNIT)/1e9
            ! write(*,*) sizeof(mat%a_elt)/1e9

            call s_umfpack_di_solve(sys = UMFPACK_A,              &
                                      x = mat%x,                  &
                                      b = mat%b,                  &
                                numeric = mat%matumfp%c_numeric,  &
                                control = mat%matumfp%c_control,  &
                                   info = mat%matumfp%c_info)

            if (mat%matumfp%c_info(UMFPACK_STATUS) < 0) then
               write(OPU, *) 'error occurred in umfpack_di_solve: ', mat%matumfp%c_info(UMFPACK_STATUS)
               stop 'Error in SOLUTION_SOLVER'
            endif

         case default
            stop 'Unknown solver type, SOLUTION_SOLVER'

      endselect

   return
   endsubroutine solution_solver