* Constructor.
*
* The @p exec_type defines the paradigm where the solution is computed.
+ * It is a string and the choices are "omp" , "reference" or "cuda".
+ * The respective strings create the respective executors as given below.
+ *
* Ginkgo currently supports three different executor types:
*
* + OmpExecutor specifies that the data should be stored and the
* The @p solver_control object is the same as for other
* deal.II iterative solvers.
*/
- SolverBase(SolverControl &solver_control, std::string exec_type);
+ SolverBase(SolverControl &solver_control, const std::string &exec_type);
/**
* Destructor.
std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>> system_matrix;
/**
- * The execution paradigm to be set by the user. The choices are between
- * `omp`, `cuda` and `reference`
- * and more details can be found in Ginkgo's documentation.
+ * The execution paradigm as a string to be set by the user . The choices
+ * are between `omp`, `cuda` and `reference` and more details can be found
+ * in Ginkgo's documentation.
*/
- std::string exec_type;
+ const std::string exec_type;
};
* @p exec_type The execution paradigm for the CG solver.
*/
SolverCG(SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data = AdditionalData());
/**
*
* @p preconditioner The preconditioner for the solver.
*/
- SolverCG(SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data = AdditionalData());
+ SolverCG(SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData & data = AdditionalData());
protected:
/**
/**
- * An implementation of the solver interface using the Ginkgo BICGSTAB solver.
+ * An implementation of the solver interface using the Ginkgo Bicgstab solver.
*
* @ingroup GinkgoWrappers
*/
template <typename ValueType = double, typename IndexType = int32_t>
- class SolverBICGSTAB : public SolverBase<ValueType, IndexType>
+ class SolverBicgstab : public SolverBase<ValueType, IndexType>
{
public:
/**
* Constructor.
*
* @p solver_control The solver control object is then used to set the
- * parameters and setup the BICGSTAB solver from the BICGSTAB factory which
+ * parameters and setup the Bicgstab solver from the Bicgstab factory which
* solves the linear system.
*
- * @p exec_type The execution paradigm for the BICGSTAB solver.
+ * @p exec_type The execution paradigm for the Bicgstab solver.
*/
- SolverBICGSTAB(SolverControl & solver_control,
- std::string exec_type,
+ SolverBicgstab(SolverControl & solver_control,
+ const std::string & exec_type,
const AdditionalData &data = AdditionalData());
/**
* Constructor.
*
* @p solver_control The solver control object is then used to set the
- * parameters and setup the BICGSTAB solver from the BICGSTAB factory which
+ * parameters and setup the Bicgstab solver from the Bicgstab factory which
* solves the linear system.
*
- * @p exec_type The execution paradigm for the BICGSTAB solver.
+ * @p exec_type The execution paradigm for the Bicgstab solver.
*
* @p preconditioner The preconditioner for the solver.
*/
- SolverBICGSTAB(SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data = AdditionalData());
+ SolverBicgstab(SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData &data = AdditionalData());
protected:
/**
* @p exec_type The execution paradigm for the CGS solver.
*/
SolverCGS(SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data = AdditionalData());
/**
*
* @p preconditioner The preconditioner for the solver.
*/
- SolverCGS(SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data = AdditionalData());
+ SolverCGS(SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData &data = AdditionalData());
protected:
/**
* @p exec_type The execution paradigm for the FCG solver.
*/
SolverFCG(SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data = AdditionalData());
/**
*
* @p preconditioner The preconditioner for the solver.
*/
- SolverFCG(SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data = AdditionalData());
+ SolverFCG(SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData &data = AdditionalData());
protected:
/**
* @p exec_type The execution paradigm for the GMRES solver.
*/
SolverGMRES(SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data = AdditionalData());
/**
*
* @p preconditioner The preconditioner for the solver.
*/
- SolverGMRES(SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data = AdditionalData());
+ SolverGMRES(SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData &data = AdditionalData());
protected:
/**
* @p exec_type The execution paradigm for the IR solver.
*/
SolverIR(SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data = AdditionalData());
/**
*
* @p inner_solver The Inner solver for the IR solver.
*/
- SolverIR(SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> inner_solver,
- const AdditionalData & data = AdditionalData());
+ SolverIR(SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &inner_solver,
+ const AdditionalData & data = AdditionalData());
protected:
/**
{
template <typename ValueType, typename IndexType>
SolverBase<ValueType, IndexType>::SolverBase(SolverControl &solver_control,
- std::string exec_type)
+ const std::string &exec_type)
: solver_control(solver_control)
, exec_type(exec_type)
{
else
{
std::cerr
- << "exec_type needs to be one of the three strings: reference, cuda or omp, but provided with"
+ << "exec_type needs to be one of the three strings: \"reference\", \"cuda\" or \"omp\", but provided with"
<< exec_type << std::endl;
std::exit(-1);
}
template <typename ValueType, typename IndexType>
SolverCG<ValueType, IndexType>::SolverCG(SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
template <typename ValueType, typename IndexType>
SolverCG<ValueType, IndexType>::SolverCG(
- SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data)
+ SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData & data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
{
}
- /* ---------------------- SolverBICGSTAB ------------------------ */
+ /* ---------------------- SolverBicgstab ------------------------ */
template <typename ValueType, typename IndexType>
- SolverBICGSTAB<ValueType, IndexType>::SolverBICGSTAB(
+ SolverBicgstab<ValueType, IndexType>::SolverBicgstab(
SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
}
template <typename ValueType, typename IndexType>
- SolverBICGSTAB<ValueType, IndexType>::SolverBICGSTAB(
- SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data)
+ SolverBicgstab<ValueType, IndexType>::SolverBicgstab(
+ SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData & data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
{
/* ---------------------- SolverCGS ------------------------ */
template <typename ValueType, typename IndexType>
- SolverCGS<ValueType, IndexType>::SolverCGS(SolverControl &solver_control,
- std::string exec_type,
+ SolverCGS<ValueType, IndexType>::SolverCGS(SolverControl & solver_control,
+ const std::string &exec_type,
const AdditionalData &data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
template <typename ValueType, typename IndexType>
SolverCGS<ValueType, IndexType>::SolverCGS(
- SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data)
+ SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData & data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
{
/* ---------------------- SolverFCG ------------------------ */
template <typename ValueType, typename IndexType>
- SolverFCG<ValueType, IndexType>::SolverFCG(SolverControl &solver_control,
- std::string exec_type,
+ SolverFCG<ValueType, IndexType>::SolverFCG(SolverControl & solver_control,
+ const std::string &exec_type,
const AdditionalData &data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
template <typename ValueType, typename IndexType>
SolverFCG<ValueType, IndexType>::SolverFCG(
- SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data)
+ SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData & data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
{
template <typename ValueType, typename IndexType>
SolverGMRES<ValueType, IndexType>::SolverGMRES(SolverControl &solver_control,
- std::string exec_type,
+ const std::string &exec_type,
const AdditionalData &data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
template <typename ValueType, typename IndexType>
SolverGMRES<ValueType, IndexType>::SolverGMRES(
- SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> preconditioner,
- const AdditionalData & data)
+ SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &preconditioner,
+ const AdditionalData & data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
{
template <typename ValueType, typename IndexType>
SolverIR<ValueType, IndexType>::SolverIR(SolverControl & solver_control,
- std::string exec_type,
+ const std::string & exec_type,
const AdditionalData &data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
template <typename ValueType, typename IndexType>
SolverIR<ValueType, IndexType>::SolverIR(
- SolverControl & solver_control,
- std::string exec_type,
- std::shared_ptr<gko::LinOpFactory> inner_solver,
- const AdditionalData & data)
+ SolverControl & solver_control,
+ const std::string & exec_type,
+ const std::shared_ptr<gko::LinOpFactory> &inner_solver,
+ const AdditionalData & data)
: SolverBase<ValueType, IndexType>(solver_control, exec_type)
, additional_data(data)
{
DEALII_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(DECLARE_SOLVER_CG)
# undef DECLARE_SOLVER_CG
-# define DECLARE_SOLVER_BICGSTAB(ValueType, IndexType) \
- class SolverBICGSTAB<ValueType, IndexType>
- DEALII_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(DECLARE_SOLVER_BICGSTAB)
-# undef DECLARE_SOLVER_BICGSTAB
+# define DECLARE_SOLVER_Bicgstab(ValueType, IndexType) \
+ class SolverBicgstab<ValueType, IndexType>
+ DEALII_INSTANTIATE_FOR_EACH_VALUE_AND_INDEX_TYPE(DECLARE_SOLVER_Bicgstab)
+# undef DECLARE_SOLVER_Bicgstab
# define DECLARE_SOLVER_CGS(ValueType, IndexType) \
class SolverCGS<ValueType, IndexType>