57 "Iterative solver for full matrix system.");
68 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
73 "This routine should only be used when using an Iterative "
74 "conjugate gradient matrix solve.");
103 bool dirForcCalculated = (bool)pDirForcing.size();
104 int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
105 int nGlobDofs = pLocToGloMap->GetNumGlobalCoeffs();
106 int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
108 int nDirTotal = nDirDofs;
109 std::shared_ptr<MultiRegions::ExpList> expList =
m_expList.lock();
110 expList->GetComm()->GetRowComm()->AllReduce(nDirTotal,
118 if (dirForcCalculated)
122 pDirForcing.size() >= nLocDofs,
123 "DirForcing is not of sufficient size. Is it in local space?");
124 Vmath::Vsub(nLocDofs, pLocInput, 1, pDirForcing, 1, rhs, 1);
130 expList->GeneralMatrixOp(
m_linSysKey, pLocOutput, rhs);
140 int offset = expList->GetCoeff_Offset(n);
144 for (rBC = r.second; rBC; rBC = rBC->next)
146 vExp->AddRobinTraceContribution(
147 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs,
148 pLocOutput + offset, rhsloc = rhs + offset);
151 Vmath::Vsub(nLocDofs, pLocInput, 1, rhs, 1, rhs, 1);
154 if (std::dynamic_pointer_cast<AssemblyMapCG>(pLocToGloMap))
162 Vmath::Vadd(nLocDofs, diff, 1, pLocOutput, 1, pLocOutput, 1);
166 ASSERTL0(
false,
"Need DG solve if using Dir BCs");
184 std::shared_ptr<MultiRegions::ExpList> expList =
m_expList.lock();
188 int ncoeffs = expList->GetNcoeffs();
197 expList->GeneralMatrixOp(
m_linSysKey, InputLoc, OutputLoc);
204 asmMap->GlobalToLocal(pInput, InputLoc);
206 expList->GeneralMatrixOp(
m_linSysKey, InputLoc, OutputLoc);
217 int offset = expList->GetCoeff_Offset(n);
221 for (rBC = r.second; rBC; rBC = rBC->next)
223 vExp->AddRobinTraceContribution(
224 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs, InputLoc + offset,
225 tmp = OutputLoc + offset);
230 if (isLocal ==
false)
232 asmMap->Assemble(OutputLoc, pOutput);
256 m_expList.lock()->GetComm()->GetRowComm();
264 "' is not defined.\n");
268 string variable = plocToGloMap->GetVariable();
272 if (pSession->DefinesGlobalSysSolnInfo(variable,
273 "IterativeSolverTolerance"))
275 sysKey.m_NekLinSysTolerance = boost::lexical_cast<double>(
277 ->GetGlobalSysSolnInfo(variable,
"IterativeSolverTolerance")
282 pSession->LoadParameter(
"IterativeSolverTolerance",
283 sysKey.m_NekLinSysTolerance,
287 if (pSession->DefinesGlobalSysSolnInfo(variable,
288 "NekLinSysMaxIterations"))
290 sysKey.m_NekLinSysMaxIterations = boost::lexical_cast<int>(
292 ->GetGlobalSysSolnInfo(variable,
"NekLinSysMaxIterations")
297 pSession->LoadParameter(
"NekLinSysMaxIterations",
298 sysKey.m_NekLinSysMaxIterations, 5000);
301 if (pSession->DefinesGlobalSysSolnInfo(variable,
"LinSysMaxStorage"))
303 sysKey.m_LinSysMaxStorage = boost::lexical_cast<int>(
304 pSession->GetGlobalSysSolnInfo(variable,
"LinSysMaxStorage")
309 pSession->LoadParameter(
"LinSysMaxStorage",
310 sysKey.m_LinSysMaxStorage, 100);
313 if (pSession->DefinesGlobalSysSolnInfo(variable,
"GMRESMaxHessMatBand"))
315 sysKey.m_KrylovMaxHessMatBand = boost::lexical_cast<int>(
316 pSession->GetGlobalSysSolnInfo(variable,
"GMRESMaxHessMatBand")
321 pSession->LoadParameter(
"GMRESMaxHessMatBand",
322 sysKey.m_KrylovMaxHessMatBand,
323 sysKey.m_LinSysMaxStorage + 1);
328 pSession->MatchSolverInfo(
"GMRESLeftPrecon",
"True",
329 sysKey.m_NekLinSysLeftPrecon,
false);
330 pSession->MatchSolverInfo(
"GMRESRightPrecon",
"True",
331 sysKey.m_NekLinSysRightPrecon,
true);
353 plocToGloMap->Assemble(pInput, gloIn);
355 plocToGloMap->GlobalToLocal(gloOut, pOutput);
361 int nLocDofs = plocToGloMap->GetNumLocalCoeffs();
363 m_linsol->SolveSystem(nLocDofs, pInput, pOutput, nDir);
369 plocToGloMap->Assemble(pInput, gloIn);
370 m_linsol->SolveSystem(nGlobal, gloIn, gloOut, nDir);
371 plocToGloMap->GlobalToLocal(gloOut, pOutput);
#define ASSERTL0(condition, msg)
#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.
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.
void v_DoMatrixMultiply(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
std::weak_ptr< AssemblyMap > m_locToGloMap
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.
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
bool m_isAbsoluteTolerance
PreconditionerSharedPtr m_precon
void DoProjection(const int pNumRows, const Array< OneD, const NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const int pNumDir, const bool isAconjugate)
projection technique
bool m_useProjection
Whether to apply projection technique.
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
static const NekDouble kNekIterativeTol
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.