59 "Iterative solver for full matrix system.");
70 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
75 "This routine should only be used when using an Iterative "
76 "conjugate gradient matrix solve.");
110 std::shared_ptr<MultiRegions::ExpList> expList =
m_expList.lock();
114 if (std::dynamic_pointer_cast<AssemblyMapCG>(pLocToGloMap))
118 else if (std::dynamic_pointer_cast<AssemblyMapDG>(pLocToGloMap))
127 bool dirForcCalculated = (bool)pDirForcing.size();
128 int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
129 int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
130 int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
132 int nDirTotal = nDirDofs;
133 expList->GetComm()->GetRowComm()->AllReduce(nDirTotal,
142 if (dirForcCalculated)
144 Vmath::Vsub(nLocDofs, pLocInput, 1, pDirForcing, 1, tmp1, 1);
150 expList->GeneralMatrixOp(
m_linSysKey, pLocOutput, tmp);
160 int offset = expList->GetCoeff_Offset(n);
164 for (rBC = r.second; rBC; rBC = rBC->next)
166 vExp->AddRobinTraceContribution(
167 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs,
168 pLocOutput + offset, tmploc = tmp + offset);
172 Vmath::Vsub(nLocDofs, pLocInput, 1, tmp, 1, tmp1, 1);
180 Vmath::Vadd(nLocDofs, tmp, 1, pLocOutput, 1, pLocOutput, 1);
184 ASSERTL0(
false,
"Need DG solve if using Dir BCs");
202 std::shared_ptr<MultiRegions::ExpList> expList =
m_expList.lock();
206 int ncoeffs = expList->GetNcoeffs();
215 expList->GeneralMatrixOp(
m_linSysKey, InputLoc, OutputLoc);
222 asmMap->GlobalToLocal(pInput, InputLoc);
224 expList->GeneralMatrixOp(
m_linSysKey, InputLoc, OutputLoc);
235 int offset = expList->GetCoeff_Offset(n);
239 for (rBC = r.second; rBC; rBC = rBC->next)
241 vExp->AddRobinTraceContribution(
242 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, InputLoc + offset,
243 tmp = OutputLoc + offset);
248 if (isLocal ==
false)
250 asmMap->Assemble(OutputLoc, pOutput);
274 m_expList.lock()->GetComm()->GetRowComm();
282 "' is not defined.\n");
304 plocToGloMap->Assemble(pInput, gloIn);
306 plocToGloMap->GlobalToLocal(gloOut, pOutput);
312 int nLocDofs = plocToGloMap->GetNumLocalCoeffs();
320 plocToGloMap->Assemble(pInput, gloIn);
322 plocToGloMap->GlobalToLocal(gloOut, pOutput);
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
void SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir=0)
Solve the linear system for given input and output vectors.
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
PreconditionerSharedPtr CreatePrecon(AssemblyMapSharedPtr asmMap)
Create a preconditioner object from the parameters defined in the supplied assembly map.
virtual ~GlobalLinSysIterativeFull()
virtual void v_UniqueMap() override
static GlobalLinSysSharedPtr create(const GlobalLinSysKey &pLinSysKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
GlobalLinSysIterativeFull(const GlobalLinSysKey &pLinSysKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Constructor for full direct matrix solve.
static std::string className
Name of class.
virtual void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
std::weak_ptr< AssemblyMap > m_locToGloMap
virtual void v_Solve(const Array< OneD, const NekDouble > &in, Array< OneD, NekDouble > &out, const AssemblyMapSharedPtr &locToGloMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray) override
Solve the linear system for given input and output vectors using a specified local to global map.
virtual void v_SolveLinearSystem(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const AssemblyMapSharedPtr &locToGloMap, const int pNumDir) override
Solve the matrix system.
LibUtilities::NekLinSysIterSharedPtr m_linsol
PreconditionerSharedPtr m_precon
void DoProjection(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir, const NekDouble tol, const bool isAconjugate)
projection technique
bool m_useProjection
Whether to apply projection technique.
NekDouble m_tolerance
Tolerance of iterative solver.
Array< OneD, int > m_map
Global to universal unique map.
NekDouble m_rhs_magnitude
dot product of rhs to normalise stopping criterion
LibUtilities::NekSysOperators m_NekSysOp
std::string m_linSysIterSolver
Iterative solver: Conjugate Gradient, GMRES.
Describe a linear system.
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
NekLinSysIterFactory & GetNekLinSysIterFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
The above copyright notice and this permission notice shall be included.
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.