|
ADflow
v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
|
Data Types | |
| type | arr1int4 |
| type | arr3int4 |
| type | arr4int4 |
Functions/Subroutines | |
| subroutine | setupamg (inputMat, nCell, blockSize, levels, nSmooth) |
| subroutine | destroyamg |
| subroutine | applyshellpc (pc, x, y, ierr) |
| subroutine | setupshellpc (pc, ierr) |
| subroutine | destroyshellpc (pc, ierr) |
| subroutine | restrictvec (x, y, interp) |
| subroutine | prolongvec (x, y, interp) |
| recursive subroutine | mgprecon (r, y, k) |
Variables | |
| integer(kind=inttype) | amglevels |
| integer(kind=inttype) | amgouterits |
| integer(kind=inttype) | amgnsmooth |
| integer(kind=inttype) | amglocalpreconits |
| integer(kind=inttype) | amglocalpreconitsfine |
| integer(kind=inttype) | amglocalpreconitscoarse |
| integer(kind=inttype) | amgasmoverlap |
| integer(kind=inttype) | amgasmoverlapfine |
| integer(kind=inttype) | amgasmoverlapcoarse |
| integer(kind=inttype) | amgfilllevel |
| integer(kind=inttype) | amgfilllevelfine |
| integer(kind=inttype) | amgfilllevelcoarse |
| character(len=maxstringlen) | amgmatrixordering |
| type(arr1int4), dimension(:), allocatable | interps |
| type(arr3int4), dimension(:, :), allocatable, target | coarseindices |
| type(arr4int4), dimension(:, :), allocatable, target | coarseoversetindices |
| logical | amgsetup = .False. |
| integer | bs |
| subroutine amg::applyshellpc | ( | pc, | |
| x, | |||
| y, | |||
| integer(kind=inttype) | ierr | ||
| ) |
Definition at line 467 of file amg.F90.
References amglevels, amgouterits, utils::echk(), mgprecon(), constants::one, and constants::zero.
Referenced by anksolver::formjacobianank(), adjointapi::setuppetscksp(), and adjointutils::setupstandardmultigrid().


| subroutine amg::destroyamg |
Definition at line 433 of file amg.F90.
References amglevels, amgsetup, coarseindices, coarseoversetindices, utils::echk(), and interps.
Referenced by anksolver::destroyanksolver(), nksolver::destroynksolver(), and adjointutils::destroypetscvars().


| subroutine amg::destroyshellpc | ( | pc, | |
| integer(kind=inttype) | ierr | ||
| ) |
Definition at line 619 of file amg.F90.
Referenced by anksolver::formjacobianank(), adjointapi::setuppetscksp(), and adjointutils::setupstandardmultigrid().

| recursive subroutine amg::mgprecon | ( | r, | |
| y, | |||
| integer(kind=inttype), intent(in) | k | ||
| ) |
Definition at line 712 of file amg.F90.
References amglevels, utils::echk(), interps, constants::one, prolongvec(), and restrictvec().
Referenced by applyshellpc().


| subroutine amg::prolongvec | ( | x, | |
| y, | |||
| integer(kind=inttype), dimension(:), intent(in) | interp | ||
| ) |
Definition at line 672 of file amg.F90.
References bs, and utils::echk().
Referenced by mgprecon().


| subroutine amg::restrictvec | ( | x, | |
| y, | |||
| integer(kind=inttype), dimension(:), intent(in) | interp | ||
| ) |
Definition at line 630 of file amg.F90.
References bs, utils::echk(), and constants::zero.
Referenced by mgprecon().


| subroutine amg::setupamg | ( | target | inputMat, |
| integer(kind=inttype), intent(in) | nCell, | ||
| integer(kind=inttype), intent(in) | blockSize, | ||
| integer(kind=inttype), intent(in) | levels, | ||
| integer(kind=inttype), intent(in) | nSmooth | ||
| ) |
Definition at line 75 of file amg.F90.
References communication::adflow_comm_world, amglevels, amgnsmooth, amgsetup, bs, coarseindices, coarseoversetindices, communication::commpatterncell_2nd, blockpointers::ib, blockpointers::iblank, blockpointers::il, communication::internalcell_2nd, interps, blockpointers::jb, blockpointers::jl, blockpointers::kb, blockpointers::kl, communication::myid, adjointvars::ncellslocal, communication::nproc, blockpointers::nx, blockpointers::ny, blockpointers::nz, constants::one, utils::setpointers(), and haloexchange::whalo1to1intgeneric().
Referenced by adjointapi::createpetscvars(), anksolver::setupanksolver(), and nksolver::setupnksolver().


| subroutine amg::setupshellpc | ( | pc, | |
| integer(kind=inttype) | ierr | ||
| ) |
Definition at line 519 of file amg.F90.
References communication::adflow_comm_world, amgasmoverlap, amgasmoverlapcoarse, amgasmoverlapfine, amgfilllevel, amgfilllevelcoarse, amgfilllevelfine, amglevels, amglocalpreconits, amglocalpreconitscoarse, amglocalpreconitsfine, amgmatrixordering, amgnsmooth, and utils::echk().
Referenced by anksolver::formjacobianank(), adjointapi::setuppetscksp(), and adjointutils::setupstandardmultigrid().


| integer(kind=inttype) amg::amgasmoverlap |
Definition at line 44 of file amg.F90.
Referenced by setupshellpc().
| integer(kind=inttype) amg::amgasmoverlapcoarse |
Definition at line 46 of file amg.F90.
Referenced by setupshellpc(), and adjointutils::setupstandardmultigrid().
| integer(kind=inttype) amg::amgasmoverlapfine |
Definition at line 45 of file amg.F90.
Referenced by setupshellpc(), and adjointutils::setupstandardmultigrid().
| integer(kind=inttype) amg::amgfilllevel |
Definition at line 49 of file amg.F90.
Referenced by setupshellpc().
| integer(kind=inttype) amg::amgfilllevelcoarse |
Definition at line 51 of file amg.F90.
Referenced by setupshellpc(), and adjointutils::setupstandardmultigrid().
| integer(kind=inttype) amg::amgfilllevelfine |
Definition at line 50 of file amg.F90.
Referenced by setupshellpc(), and adjointutils::setupstandardmultigrid().
| integer(kind=inttype) amg::amglevels |
Definition at line 30 of file amg.F90.
Referenced by applyshellpc(), destroyamg(), mgprecon(), setblock(), setupamg(), setupshellpc(), and adjointutils::setupstateresidualmatrix().
| integer(kind=inttype) amg::amglocalpreconits |
Definition at line 39 of file amg.F90.
Referenced by setupshellpc().
| integer(kind=inttype) amg::amglocalpreconitscoarse |
Definition at line 41 of file amg.F90.
Referenced by setupshellpc(), and adjointutils::setupstandardmultigrid().
| integer(kind=inttype) amg::amglocalpreconitsfine |
Definition at line 40 of file amg.F90.
Referenced by setupshellpc(), and adjointutils::setupstandardmultigrid().
| character(len=maxstringlen) amg::amgmatrixordering |
Definition at line 54 of file amg.F90.
Referenced by setupshellpc(), and adjointutils::setupstandardmultigrid().
| integer(kind=inttype) amg::amgnsmooth |
Definition at line 36 of file amg.F90.
Referenced by setupamg(), and setupshellpc().
| integer(kind=inttype) amg::amgouterits |
Definition at line 33 of file amg.F90.
Referenced by applyshellpc(), and adjointutils::setupstandardmultigrid().
| logical amg::amgsetup = .False. |
Definition at line 71 of file amg.F90.
Referenced by destroyamg(), and setupamg().
| integer amg::bs |
Definition at line 72 of file amg.F90.
Referenced by prolongvec(), restrictvec(), and setupamg().
| type(arr3int4), dimension(:, :), allocatable, target amg::coarseindices |
Definition at line 68 of file amg.F90.
Referenced by anksolver::computetimestepmat(), destroyamg(), anksolver::formjacobianank(), setupamg(), and adjointutils::setupstateresidualmatrix().
| type(arr4int4), dimension(:, :), allocatable, target amg::coarseoversetindices |
Definition at line 69 of file amg.F90.
Referenced by destroyamg(), setupamg(), and adjointutils::setupstateresidualmatrix().
| type(arr1int4), dimension(:), allocatable amg::interps |
Definition at line 65 of file amg.F90.
Referenced by destroyamg(), mgprecon(), and setupamg().