Skip to content

Commit

Permalink
clear up function, update example
Browse files Browse the repository at this point in the history
Co-authored-by: Pratik Nayak <pratikvn@protonmail.com>
Co-authored-by: Terry Cojean <terry.cojean@kit.edu>
  • Loading branch information
3 people committed Nov 6, 2022
1 parent af76d58 commit b609166
Show file tree
Hide file tree
Showing 4 changed files with 29 additions and 26 deletions.
28 changes: 14 additions & 14 deletions core/solver/multigrid.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -160,6 +160,14 @@ void handle_list(
}


template <typename Vec>
void clear_and_reserve(Vec& vec, size_type size)
{
vec.clear();
vec.reserve(size);
}


} // namespace


Expand All @@ -178,20 +186,12 @@ void MultigridState::generate(const LinOp* system_matrix_in,
auto mg_level_list = multigrid->get_mg_level_list();
auto list_size = mg_level_list.size();
auto cycle = multigrid->get_cycle();
r_list.clear();
g_list.clear();
e_list.clear();
b_list.clear();
x_list.clear();
next_one_list.clear();
one_list.clear();
neg_one_list.clear();
r_list.reserve(list_size);
g_list.reserve(list_size);
e_list.reserve(list_size);
one_list.reserve(list_size);
next_one_list.reserve(list_size);
neg_one_list.reserve(list_size);
clear_and_reserve(r_list, list_size);
clear_and_reserve(g_list, list_size);
clear_and_reserve(e_list, list_size);
clear_and_reserve(one_list, list_size);
clear_and_reserve(next_one_list, list_size);
clear_and_reserve(neg_one_list, list_size);
// Allocate memory first such that reusing allocation in each iter.
for (int i = 0; i < mg_level_list.size(); i++) {
auto next_nrows = mg_level_list.at(i)->get_coarse_op()->get_size()[0];
Expand Down
1 change: 1 addition & 0 deletions cuda/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -110,6 +110,7 @@ target_include_directories(ginkgo_cuda
target_include_directories(ginkgo_cuda
PRIVATE ${CMAKE_CURRENT_BINARY_DIR}/..)
target_link_libraries(ginkgo_cuda PRIVATE ${CUDA_RUNTIME_LIBS} ${CUBLAS} ${CUSPARSE} ${CURAND} ${CUFFT})
target_link_libraries(ginkgo_cuda PUBLIC ginkgo_device)
target_compile_options(ginkgo_cuda
PRIVATE "$<$<COMPILE_LANGUAGE:CUDA>:${GINKGO_CUDA_ARCH_FLAGS}>")
# we handle CUDA architecture flags for now, disable CMake handling
Expand Down
24 changes: 14 additions & 10 deletions examples/mixed-multigrid-solver/mixed-multigrid-solver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -45,19 +45,19 @@ int main(int argc, char* argv[])
{
// Some shortcuts
using ValueType = double;
using ValueType2 = float;
using MixedType = float;
using IndexType = int;
using vec = gko::matrix::Dense<ValueType>;
using mtx = gko::matrix::Csr<ValueType, IndexType>;
using fcg = gko::solver::Fcg<ValueType>;
using cg = gko::solver::Cg<ValueType2>;
using cg = gko::solver::Cg<MixedType>;
using ir = gko::solver::Ir<ValueType>;
using ir2 = gko::solver::Ir<ValueType2>;
using ir2 = gko::solver::Ir<MixedType>;
using mg = gko::solver::Multigrid;
using bj = gko::preconditioner::Jacobi<ValueType, IndexType>;
using bj2 = gko::preconditioner::Jacobi<ValueType2, IndexType>;
using bj2 = gko::preconditioner::Jacobi<MixedType, IndexType>;
using pgm = gko::multigrid::Pgm<ValueType, IndexType>;
using pgm2 = gko::multigrid::Pgm<ValueType2, IndexType>;
using pgm2 = gko::multigrid::Pgm<MixedType, IndexType>;

// Print version information
std::cout << gko::version_info::get() << std::endl;
Expand Down Expand Up @@ -87,7 +87,7 @@ int main(int argc, char* argv[])
// executor where Ginkgo will perform the computation
const auto exec = exec_map.at(executor_string)(); // throws if not valid
const bool use_mixed = argc >= 5;
std::cout << "Use Mixed: " << use_mixed << std::endl;
std::cout << "Using mixed precision? " << use_mixed << std::endl;
// Read data
auto A = share(gko::read<mtx>(std::ifstream(argv[2]), exec));
// Create RHS as 1 and initial guess as 0
Expand Down Expand Up @@ -140,7 +140,7 @@ int main(int argc, char* argv[])
auto smoother_gen2 = gko::share(
ir2::build()
.with_solver(bj2::build().with_max_block_size(1u).on(exec))
.with_relaxation_factor(static_cast<ValueType2>(0.9))
.with_relaxation_factor(static_cast<MixedType>(0.9))
.with_criteria(
gko::stop::Iteration::build().with_max_iters(1u).on(exec))
.on(exec));
Expand All @@ -160,14 +160,13 @@ int main(int argc, char* argv[])
auto coarsest_solver_gen2 = gko::share(
ir2::build()
.with_solver(bj2::build().with_max_block_size(1u).on(exec))
.with_relaxation_factor(static_cast<ValueType2>(0.9))
.with_relaxation_factor(static_cast<MixedType>(0.9))
.with_criteria(
gko::stop::Iteration::build().with_max_iters(4u).on(exec))
.on(exec));
// Create multigrid factory
std::shared_ptr<gko::LinOpFactory> multigrid_gen;
if (use_mixed) {
std::cout << "USE" << std::endl;
multigrid_gen =
mg::build()
.with_max_levels(10u)
Expand All @@ -177,7 +176,12 @@ int main(int argc, char* argv[])
.with_mg_level(mg_level_gen, mg_level_gen2)
.with_level_selector([](const gko::size_type level,
const gko::LinOp*) -> gko::size_type {
std::cout << "level " << level << std::endl;
// The first (index 0) level will use the first
// mg_level_gen, smoother_gen which are the factories with
// ValueType. The rest of levels (>= 1) will use the second
// (index 1) mg_level_gen2 and smoother_gen2 which use the
// MixedType. The rest of levels will use different type
// than the normal multigrid.
return level >= 1 ? 1 : 0;
})
.with_coarsest_solver(coarsest_solver_gen2)
Expand Down
2 changes: 0 additions & 2 deletions include/ginkgo/core/solver/multigrid.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -162,8 +162,6 @@ struct MultigridState {

// current level's nrows x nrhs
std::vector<std::shared_ptr<LinOp>> r_list;
std::vector<std::shared_ptr<LinOp>> x_list;
std::vector<std::shared_ptr<LinOp>> b_list;
// next level's nrows x nrhs
std::vector<std::shared_ptr<LinOp>> g_list;
std::vector<std::shared_ptr<LinOp>> e_list;
Expand Down

0 comments on commit b609166

Please sign in to comment.