45 namespace MultiRegions
50 string PreconditionerDiagonal::className =
52 "Diagonal", PreconditionerDiagonal::create,
"Diagonal Preconditioning");
60 PreconditionerDiagonal::PreconditionerDiagonal(
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);
195 ? asmMap->GetNumGlobalCoeffs()
196 : asmMap->GetNumGlobalBndCoeffs();
197 int nDir = asmMap->GetNumGlobalDirBndCoeffs();
198 int nNonDir = nGlobal - nDir;
213 const std::shared_ptr<GlobalLinSys> &plinsys,
#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)
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Apply a preconditioner to the conjugate gradient method.
Array< OneD, NekDouble > m_diagonals
virtual void v_InitObject() override
void DiagonalPreconditionerSum(void)
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_BuildPreconditioner() override
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.
virtual void v_DoPreconditioner(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput) override
Apply a preconditioner to the conjugate gradient method.
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 Sdiv(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha/y.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)