1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217
|
#include <iostream>
#include <vector>
#include <string>
#include <complex>
#include <boost/program_options.hpp>
#include <boost/property_tree/ptree.hpp>
#include <boost/property_tree/json_parser.hpp>
#include <amgcl/backend/builtin.hpp>
#include <amgcl/value_type/complex.hpp>
#include <amgcl/adapter/crs_tuple.hpp>
#include <amgcl/mpi/util.hpp>
#include <amgcl/mpi/make_solver.hpp>
#include <amgcl/mpi/preconditioner.hpp>
#include <amgcl/mpi/solver/runtime.hpp>
#include <amgcl/io/mm.hpp>
#include <amgcl/io/binary.hpp>
#include <amgcl/profiler.hpp>
namespace amgcl {
profiler<> prof;
}
namespace math = amgcl::math;
//---------------------------------------------------------------------------
ptrdiff_t assemble_poisson3d(amgcl::mpi::communicator comm,
ptrdiff_t n, int block_size,
std::vector<ptrdiff_t> &ptr,
std::vector<ptrdiff_t> &col,
std::vector<std::complex<double>> &val,
std::vector<std::complex<double>> &rhs)
{
ptrdiff_t n3 = n * n * n;
ptrdiff_t chunk = (n3 + comm.size - 1) / comm.size;
if (chunk % block_size != 0) {
chunk += block_size - chunk % block_size;
}
ptrdiff_t row_beg = std::min(n3, chunk * comm.rank);
ptrdiff_t row_end = std::min(n3, row_beg + chunk);
chunk = row_end - row_beg;
ptr.clear(); ptr.reserve(chunk + 1);
col.clear(); col.reserve(chunk * 7);
val.clear(); val.reserve(chunk * 7);
rhs.resize(chunk);
std::fill(rhs.begin(), rhs.end(), 1.0);
const double h2i = (n - 1) * (n - 1);
ptr.push_back(0);
for (ptrdiff_t idx = row_beg; idx < row_end; ++idx) {
ptrdiff_t k = idx / (n * n);
ptrdiff_t j = (idx / n) % n;
ptrdiff_t i = idx % n;
if (k > 0) {
col.push_back(idx - n * n);
val.push_back(-h2i);
}
if (j > 0) {
col.push_back(idx - n);
val.push_back(-h2i);
}
if (i > 0) {
col.push_back(idx - 1);
val.push_back(-h2i);
}
col.push_back(idx);
val.push_back(6 * h2i);
if (i + 1 < n) {
col.push_back(idx + 1);
val.push_back(-h2i);
}
if (j + 1 < n) {
col.push_back(idx + n);
val.push_back(-h2i);
}
if (k + 1 < n) {
col.push_back(idx + n * n);
val.push_back(-h2i);
}
ptr.push_back( col.size() );
}
return chunk;
}
//---------------------------------------------------------------------------
void solve_scalar(
amgcl::mpi::communicator comm,
ptrdiff_t chunk,
const std::vector<ptrdiff_t> &ptr,
const std::vector<ptrdiff_t> &col,
const std::vector<std::complex<double>> &val,
const boost::property_tree::ptree &prm,
const std::vector<std::complex<double>> &rhs
)
{
typedef amgcl::backend::builtin<std::complex<double>> Backend;
typedef
amgcl::mpi::make_solver<
amgcl::runtime::mpi::preconditioner<Backend>,
amgcl::runtime::mpi::solver::wrapper<Backend>
>
Solver;
using amgcl::prof;
prof.tic("setup");
Solver solve(comm, std::tie(chunk, ptr, col, val), prm);
prof.toc("setup");
if (comm.rank == 0) {
std::cout << solve << std::endl;
}
std::vector<std::complex<double>> x(chunk);
int iters;
double error;
prof.tic("solve");
std::tie(iters, error) = solve(rhs, x);
prof.toc("solve");
if (comm.rank == 0) {
std::cout
<< "Iterations: " << iters << std::endl
<< "Error: " << error << std::endl
<< prof << std::endl;
}
}
//---------------------------------------------------------------------------
int main(int argc, char *argv[]) {
amgcl::mpi::init_thread mpi(&argc, &argv);
amgcl::mpi::communicator comm(MPI_COMM_WORLD);
if (comm.rank == 0)
std::cout << "World size: " << comm.size << std::endl;
using amgcl::prof;
// Read configuration from command line
namespace po = boost::program_options;
po::options_description desc("Options");
desc.add_options()
("help,h", "show help")
(
"size,n",
po::value<ptrdiff_t>()->default_value(128),
"domain size"
)
("prm-file,P",
po::value<std::string>(),
"Parameter file in json format. "
)
(
"prm,p",
po::value< std::vector<std::string> >()->multitoken(),
"Parameters specified as name=value pairs. "
"May be provided multiple times. Examples:\n"
" -p solver.tol=1e-3\n"
" -p precond.coarse_enough=300"
)
;
po::positional_options_description p;
p.add("prm", -1);
po::variables_map vm;
po::store(po::command_line_parser(argc, argv).options(desc).positional(p).run(), vm);
po::notify(vm);
if (vm.count("help")) {
if (comm.rank == 0) std::cout << desc << std::endl;
return 0;
}
boost::property_tree::ptree prm;
if (vm.count("prm-file")) {
read_json(vm["prm-file"].as<std::string>(), prm);
}
if (vm.count("prm")) {
for(const std::string &v : vm["prm"].as<std::vector<std::string> >()) {
amgcl::put(prm, v);
}
}
ptrdiff_t n;
std::vector<ptrdiff_t> ptr;
std::vector<ptrdiff_t> col;
std::vector<std::complex<double>> val;
std::vector<std::complex<double>> rhs;
prof.tic("assemble");
n = assemble_poisson3d(comm, vm["size"].as<ptrdiff_t>(), 1, ptr, col, val, rhs);
prof.toc("assemble");
solve_scalar(comm, n, ptr, col, val, prm, rhs);
}
|