KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
amgcl_mpi_solver.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Denis Demidov
11 // Riccardo Rossi
12 //
13 
14 #if !defined (KRATOS_AMGCL_MPI_SOLVER_H_INCLUDED)
15 #define KRATOS_AMGCL_MPI_SOLVER_H_INCLUDED
16 
17 #ifndef AMGCL_PARAM_UNKNOWN
18 #include "input_output/logger.h"
19 # define AMGCL_PARAM_UNKNOWN(name) \
20  Kratos::Logger("AMGCL") << KRATOS_CODE_LOCATION << Kratos::Logger::Severity::WARNING << "Unknown parameter " << name << std::endl
21 #endif
22 
23 // System includes
24 #include <iostream>
25 #include <fstream>
26 #include <utility>
27 
28 // Project includes
29 #include "includes/define.h"
33 
34 namespace Kratos
35 {
38 
49 template< class TSparseSpaceType, class TDenseSpaceType,
50  class TReordererType = Reorderer<TSparseSpaceType, TDenseSpaceType> >
51 class AmgclMPISolver : public AMGCLSolver< TSparseSpaceType,TDenseSpaceType, TReordererType >
52 {
53 public:
56 
59 
62 
65 
68 
71 
73  typedef std::size_t IndexType;
74 
76  typedef std::size_t SizeType;
77 
81 
86  AmgclMPISolver(Parameters ThisParameters = Parameters(R"({})"))
87  :
88  AMGCLSolver< TSparseSpaceType,TDenseSpaceType, TReordererType>(ThisParameters) { }
89 
91  AmgclMPISolver(const AmgclMPISolver& Other) = delete;
92 
94  ~AmgclMPISolver() override = default;
95 
99 
101  AmgclMPISolver& operator=(const AmgclMPISolver& Other) = delete;
102 
106 
108  void PrintInfo(std::ostream& rOStream) const override
109  {
110  rOStream << "AMGCL-MPI-Solver";
111  }
112 
114 
115 }; // Class AmgclMPISolver
116 
118 template<class TSparseSpaceType, class TDenseSpaceType, class TReordererType>
119 inline std::ostream& operator << (std::ostream& rOStream,
120  const AmgclMPISolver<TSparseSpaceType,
121  TDenseSpaceType, TReordererType>& rThis)
122 {
123  rThis.PrintInfo(rOStream);
124  rOStream << std::endl;
125  rThis.PrintData(rOStream);
126 
127  return rOStream;
128 }
129 
130 } // namespace Kratos.
131 
132 #endif // KRATOS_AMGCL_MPI_SOLVER_H_INCLUDED defined
This is a multigrid solver based on the AMGCL library.
Definition: amgcl_solver.h:112
This is a multigrid solver based on the AMGCL library.
Definition: amgcl_mpi_solver.h:52
TSparseSpaceType::MatrixType SparseMatrixType
The sparse matric type.
Definition: amgcl_mpi_solver.h:61
AmgclMPISolver(const AmgclMPISolver &Other)=delete
Copy constructor.
TDenseSpaceType::MatrixType DenseMatrixType
Dense matrix type.
Definition: amgcl_mpi_solver.h:67
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: amgcl_mpi_solver.h:108
AmgclMPISolver(Parameters ThisParameters=Parameters(R"({})"))
This is the default constructor.
Definition: amgcl_mpi_solver.h:86
~AmgclMPISolver() override=default
Destructor.
std::size_t IndexType
The index type definition.
Definition: amgcl_mpi_solver.h:73
KRATOS_CLASS_POINTER_DEFINITION(AmgclMPISolver)
Pointer definition of AmgclMPISolver.
TSparseSpaceType::VectorType VectorType
Vector type definition.
Definition: amgcl_mpi_solver.h:64
AmgclMPISolver & operator=(const AmgclMPISolver &Other)=delete
Assignment operator.
ModelPart::DofsArrayType DofsArrayType
DofArray type.
Definition: amgcl_mpi_solver.h:70
std::size_t SizeType
The size type definition.
Definition: amgcl_mpi_solver.h:76
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
Vector VectorType
Definition: geometrical_transformation_utilities.h:56
Matrix MatrixType
Definition: geometrical_transformation_utilities.h:55
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432