61 const std::shared_ptr<GlobalLinSys> &plinsys,
87 ASSERTL0(0,
"Unsupported solver type");
97 std::shared_ptr<MultiRegions::ExpList> expList =
98 ((
m_linsys.lock())->GetLocMat()).lock();
104 int i, j, n, cnt, gid1, gid2;
106 int nGlobal = asmMap->GetNumGlobalCoeffs();
107 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
108 int nInt = nGlobal - nDir;
115 int nElmt = expList->GetNumElmts();
116 for (n = cnt = 0; n < nElmt; ++n)
118 loc_mat = (
m_linsys.lock())->GetBlock(n);
119 loc_row = loc_mat->GetRows();
121 for (i = 0; i < loc_row; ++i)
123 gid1 = asmMap->GetLocalToGlobalMap(cnt + i) - nDir;
124 sign1 = asmMap->GetLocalToGlobalSign(cnt + i);
128 for (j = 0; j < loc_row; ++j)
130 gid2 = asmMap->GetLocalToGlobalMap(cnt + j) - nDir;
131 sign2 = asmMap->GetLocalToGlobalSign(cnt + j);
138 value = vOutput[gid1 + nDir] +
139 sign1 * sign2 * (*loc_mat)(i, j);
140 vOutput[gid1 + nDir] = value;
149 asmMap->UniversalAssemble(vOutput);
164 int nGlobalBnd = asmMap->GetNumGlobalBndCoeffs();
165 int nDirBnd = asmMap->GetNumGlobalDirBndCoeffs();
166 int rows = nGlobalBnd - nDirBnd;
172 for (
unsigned int i = 0; i < rows; ++i)
174 vOutput[nDirBnd + i] = diagonals[i];
178 asmMap->UniversalAssembleBnd(vOutput);
197 int nGlobal = (isFull) ? asmMap->GetNumGlobalCoeffs()
198 : asmMap->GetNumGlobalBndCoeffs();
199 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
200 int nNonDir = nGlobal - nDir;
205 (isFull) ? asmMap->Assemble(pInput, wk)
206 : asmMap->AssembleBnd(pInput, wk);
208 wk.data() + nDir, 1);
210 (isFull) ? asmMap->GlobalToLocal(wk, pOutput)
211 : asmMap->GlobalToLocalBnd(wk, pOutput);
230 const std::shared_ptr<GlobalLinSys> &plinsys,
258 boost::ignore_unused(isLocal);
276 const std::shared_ptr<GlobalLinSys> &plinsys,
290 auto expList = ((
m_linsys.lock())->GetLocMat()).lock();
291 std::shared_ptr<LibUtilities::SessionReader> session =
292 expList->GetSession();
296 if (session->DefinesGlobalSysSolnInfo(var,
"JacobiIterations"))
298 m_niter = boost::lexical_cast<int>(
299 session->GetGlobalSysSolnInfo(var,
"JacobiIterations").c_str());
319 ? asmMap->GetNumGlobalCoeffs()
320 : asmMap->GetNumGlobalBndCoeffs();
321 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
322 int nNonDir = nGlobal - nDir;
328 ? asmMap->GetNumLocalCoeffs()
329 : asmMap->GetNumLocalBndCoeffs();
332 asmMap->Assemble(pInput, wk);
334 wk.data() + nDir, 1);
337 for (
int n = 1; n <
m_niter; ++n)
339 asmMap->GlobalToLocal(wk, pOutput);
342 std::dynamic_pointer_cast<GlobalLinSysIterative>(
m_linsys.lock())
343 ->DoMatrixMultiply(pOutput, wk1);
347 asmMap->Assemble(wk1, pOutput);
349 1, wk.data() + nDir, 1, wk.data() + nDir, 1);
352 asmMap->GlobalToLocal(wk, pOutput);
358 wk.data() + nDir, 1);
361 for (
int n = 1; n <
m_niter; ++n)
364 std::dynamic_pointer_cast<GlobalLinSysIterative>(
m_linsys.lock())
365 ->DoMatrixMultiply(wk, wk1);
368 Vmath::Vsub(nNonDir, pInput.data(), 1, wk1.data() + nDir, 1,
369 wk1.data() + nDir, 1);
373 wk.data() + nDir, 1, wk.data() + nDir, 1);
376 Vmath::Vcopy(nNonDir, wk.data() + nDir, 1, pOutput.data(), 1);
#define ASSERTL0(condition, msg)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
virtual void v_BuildPreconditioner() override
void StaticCondDiagonalPreconditionerSum(void)
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
Array< OneD, NekDouble > m_diagonals
virtual void v_InitObject() override
void DiagonalPreconditionerSum(void)
PreconditionerDiagonal(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &IsLocal=false) override
Apply a preconditioner to the conjugate gradient method.
static std::string className
Name of class.
const std::weak_ptr< GlobalLinSys > m_linsys
Array< OneD, NekDouble > AssembleStaticCondGlobalDiagonals()
Performs global assembly of diagonal entries to global Schur complement matrix.
std::weak_ptr< AssemblyMap > m_locToGloMap
virtual void v_InitObject() override
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
PreconditionerJacobi(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
virtual void v_BuildPreconditioner() override
static std::string className
Name of class.
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &IsLocal=false) override
Apply a preconditioner to the conjugate gradient method.
virtual void v_BuildPreconditioner() override
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &isLocal=false) override
Apply a preconditioner to the conjugate gradient method.
static std::string className
Name of class.
virtual void v_InitObject() override
PreconditionerNull(const std::shared_ptr< GlobalLinSys > &plinsys, const AssemblyMapSharedPtr &pLocToGloMap)
static PreconditionerSharedPtr create(const std::shared_ptr< GlobalLinSys > &plinsys, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
std::shared_ptr< Expansion > ExpansionSharedPtr
@ eIterativeMultiLevelStaticCond
@ ePETScMultiLevelStaticCond
PreconFactory & GetPreconFactory()
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
void Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.