hypre_obj_SolveIJ Subroutine

private impure subroutine hypre_obj_SolveIJ(this, sol)

Solves the system Ax=b and return the solution, IJ interface.

Type Bound

hypre_obj

Arguments

Type IntentOptional Attributes Name
class(hypre_obj), intent(inout) :: this

Hypre machinery

type(eulerian_obj_r), intent(inout) :: sol

Solution vector


Calls

proc~~hypre_obj_solveij~~CallsGraph proc~hypre_obj_solveij hypre_obj%hypre_obj_SolveIJ hypre_boomeramggetfinalreltvres hypre_boomeramggetfinalreltvres proc~hypre_obj_solveij->hypre_boomeramggetfinalreltvres hypre_boomeramggetnumiterations hypre_boomeramggetnumiterations proc~hypre_obj_solveij->hypre_boomeramggetnumiterations hypre_boomeramgsolve hypre_boomeramgsolve proc~hypre_obj_solveij->hypre_boomeramgsolve hypre_ijvectorgetvalues hypre_ijvectorgetvalues proc~hypre_obj_solveij->hypre_ijvectorgetvalues hypre_parcsrpcggetfinalrelative hypre_parcsrpcggetfinalrelative proc~hypre_obj_solveij->hypre_parcsrpcggetfinalrelative hypre_parcsrpcggetnumiterations hypre_parcsrpcggetnumiterations proc~hypre_obj_solveij->hypre_parcsrpcggetnumiterations hypre_parcsrpcgsolve hypre_parcsrpcgsolve proc~hypre_obj_solveij->hypre_parcsrpcgsolve proc~eulerian_obj_updateghostcells eulerian_obj_base%eulerian_obj_UpdateGhostCells proc~hypre_obj_solveij->proc~eulerian_obj_updateghostcells proc~eulerian_obj_updateghostcells_x eulerian_obj_base%eulerian_obj_UpdateGhostCells_x proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_x proc~eulerian_obj_updateghostcells_y eulerian_obj_base%eulerian_obj_UpdateGhostCells_y proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_y proc~eulerian_obj_updateghostcells_z eulerian_obj_base%eulerian_obj_UpdateGhostCells_z proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_z cell cell proc~eulerian_obj_updateghostcells_x->cell mpi_irecv mpi_irecv proc~eulerian_obj_updateghostcells_x->mpi_irecv mpi_isend mpi_isend proc~eulerian_obj_updateghostcells_x->mpi_isend mpi_waitall mpi_waitall proc~eulerian_obj_updateghostcells_x->mpi_waitall proc~eulerian_obj_updateghostcells_y->cell proc~eulerian_obj_updateghostcells_y->mpi_irecv proc~eulerian_obj_updateghostcells_y->mpi_isend proc~eulerian_obj_updateghostcells_y->mpi_waitall proc~eulerian_obj_updateghostcells_z->cell proc~eulerian_obj_updateghostcells_z->mpi_irecv proc~eulerian_obj_updateghostcells_z->mpi_isend proc~eulerian_obj_updateghostcells_z->mpi_waitall

Called by

proc~~hypre_obj_solveij~~CalledByGraph proc~hypre_obj_solveij hypre_obj%hypre_obj_SolveIJ proc~hypre_obj_solve hypre_obj%hypre_obj_Solve proc~hypre_obj_solve->proc~hypre_obj_solveij proc~cdifs_obj_advancesolutioncorrector cdifs_obj_AdvanceSolutionCorrector proc~cdifs_obj_advancesolutioncorrector->proc~hypre_obj_solve proc~cdifs_obj_computesolidvf cdifs_obj%cdifs_obj_ComputeSolidVF proc~cdifs_obj_computesolidvf->proc~hypre_obj_solve proc~grans_obj_computesolidvf grans_obj%grans_obj_ComputeSolidVF proc~grans_obj_computesolidvf->proc~hypre_obj_solve proc~marker_set_computesolidvolfrac marker_set%marker_set_ComputeSolidVolFrac proc~marker_set_computesolidvolfrac->proc~hypre_obj_solve proc~cdifs_obj_advancesolution cdifs_obj_AdvanceSolution proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutioncorrector proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionrp proc~cdifs_obj_advancesolutionrp->proc~cdifs_obj_computesolidvf proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->proc~cdifs_obj_computesolidvf proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->proc~grans_obj_computesolidvf interface~grans_obj_writeoutputdata grans_obj%grans_obj_WriteOutputData proc~grans_obj_preparesolver->interface~grans_obj_writeoutputdata proc~grans_obj_writeoutputdata grans_obj_WriteOutputData proc~grans_obj_writeoutputdata->proc~grans_obj_computesolidvf interface~cdifs_obj_advancesolution cdifs_obj%cdifs_obj_AdvanceSolution interface~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolution interface~cdifs_obj_preparesolver cdifs_obj%cdifs_obj_PrepareSolver interface~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolver interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver interface~grans_obj_writeoutputdata->proc~grans_obj_writeoutputdata

Source Code

    impure subroutine hypre_obj_SolveIJ(this,sol)
      !> Solves the system Ax=b and return the solution, IJ interface.
      implicit none
      class(hypre_obj),     intent(inout) :: this                              !! Hypre machinery
      type(eulerian_obj_r), intent(inout) :: sol                               !! Solution vector
      ! Work variables
      integer :: ierr
      integer :: i,j,k,n

      select case (this%solver_name)
      case (HYPRE_SOL_IJ_AMG_NONE)
        call HYPRE_BoomerAMGSolve           (this%solver, this%par_matrix, this%par_rhs, this%par_sol, ierr)
        call HYPRE_BoomerAMGGetNumIterations(this%solver, this%it,  ierr)
        call HYPRE_BoomerAMGGetFinalReltvRes(this%solver, this%rel, ierr)
      case (HYPRE_SOL_IJ_PCG_AMG)
        call HYPRE_ParCSRPCGSolve           (this%solver, this%par_matrix, this%par_rhs, this%par_sol, ierr)
        call HYPRE_ParCSRPCGGetNumIterations(this%solver, this%it,  ierr)
        call HYPRE_ParCSRPCGGetFinalRelative(this%solver, this%rel, ierr)
      case (HYPRE_SOL_IJ_PCG_DS)
        call HYPRE_ParCSRPCGSolve           (this%solver, this%par_matrix, this%par_rhs, this%par_sol, ierr)
        call HYPRE_ParCSRPCGGetNumIterations(this%solver, this%it,  ierr)
        call HYPRE_ParCSRPCGGetFinalRelative(this%solver, this%rel, ierr)
      case default
        call this%parallel%Stop("Unknown HYPRE solver")
      end select

      ! - Transfer solution
      associate(lo=>this%block%lo, hi=>this%block%hi, irow=>this%irow, &
        irow_hi=>this%irow_hi,irow_lo =>this%irow_lo)

        call HYPRE_IJVectorGetValues(this%sol, irow_hi-irow_lo+1,this%rows,this%sol_values, ierr)

        do k=lo(3),hi(3)
          do j=lo(2),hi(2)
            do i=lo(1),hi(1)
              n = irow%cell(i,j,k) - irow_lo  + 1
              sol%cell(i,j,k) = this%sol_values(n)
            end do
          end do
        end do

        call sol%UpdateGhostCells
      end associate

      return
    end subroutine hypre_obj_SolveIJ