Nektar++
GlobalLinSysXxtStaticCond.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: GlobalLinSysXxtStaticCond.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: GlobalLinSysXxtStaticCond definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
38using namespace std;
39
41{
42/**
43 * @class GlobalLinSysIterativeStaticCond
44 *
45 * Solves a linear system iteratively using single- or multi-level
46 * static condensation.
47 */
48
49/**
50 * Registers the class with the Factory.
51 */
54 "XxtStaticCond", GlobalLinSysXxtStaticCond::create,
55 "Iterative static condensation.");
56
59 "XxtMultiLevelStaticCond", GlobalLinSysXxtStaticCond::create,
60 "Xxt multi-level static condensation.");
61
62/**
63 * For a matrix system of the form @f[
64 * \left[ \begin{array}{cc}
65 * \boldsymbol{A} & \boldsymbol{B}\\
66 * \boldsymbol{C} & \boldsymbol{D}
67 * \end{array} \right]
68 * \left[ \begin{array}{c} \boldsymbol{x_1}\\ \boldsymbol{x_2}
69 * \end{array}\right]
70 * = \left[ \begin{array}{c} \boldsymbol{y_1}\\ \boldsymbol{y_2}
71 * \end{array}\right],
72 * @f]
73 * where @f$\boldsymbol{D}@f$ and
74 * @f$(\boldsymbol{A-BD^{-1}C})@f$ are invertible, store and assemble
75 * a static condensation system, according to a given local to global
76 * mapping. #m_linSys is constructed by AssembleSchurComplement().
77 * @param mKey Associated matrix key.
78 * @param pLocMatSys LocalMatrixSystem
79 * @param locToGloMap Local to global mapping.
80 */
82 const GlobalLinSysKey &pKey, const std::weak_ptr<ExpList> &pExpList,
83 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
84 : GlobalLinSys(pKey, pExpList, pLocToGloMap),
85 GlobalLinSysXxt(pKey, pExpList, pLocToGloMap),
86 GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
87{
90 "This constructor is only valid when using static "
91 "condensation");
93 pLocToGloMap->GetGlobalSysSolnType(),
94 "The local to global map is not set up for the requested "
95 "solution type");
96}
97
98/**
99 *
100 */
102 const GlobalLinSysKey &pKey, const std::weak_ptr<ExpList> &pExpList,
103 const DNekScalBlkMatSharedPtr pSchurCompl,
105 const DNekScalBlkMatSharedPtr pInvD,
106 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
107 : GlobalLinSys(pKey, pExpList, pLocToGloMap),
108 GlobalLinSysXxt(pKey, pExpList, pLocToGloMap),
109 GlobalLinSysStaticCond(pKey, pExpList, pLocToGloMap)
110{
111 m_schurCompl = pSchurCompl;
112 m_BinvD = pBinvD;
113 m_C = pC;
114 m_invD = pInvD;
115 m_locToGloMap = pLocToGloMap;
116}
117
118/**
119 *
120 */
122{
123}
124
125/**
126 * Construct the local matrix row index, column index and value index
127 * arrays and initialize the XXT data structure with this information.
128 * @param locToGloMap Local to global mapping information.
129 */
131 std::shared_ptr<AssemblyMap> pLocToGloMap)
132{
133 ExpListSharedPtr vExp = m_expList.lock();
134 unsigned int nElmt = m_schurCompl->GetNumberOfBlockRows();
135 DNekScalMatSharedPtr loc_mat;
136 unsigned int iCount = 0;
137 unsigned int rCount = 0;
138 unsigned int nRows = 0;
139 unsigned int nEntries = 0;
140 unsigned int numDirBnd = pLocToGloMap->GetNumGlobalDirBndCoeffs();
141 unsigned int nLocal = pLocToGloMap->GetNumLocalBndCoeffs();
142 const Array<OneD, NekDouble> &vMapSign =
143 pLocToGloMap->GetLocalToGlobalBndSign();
144 bool doSign = pLocToGloMap->GetSignChange();
145 unsigned int i = 0, j = 0, k = 0, n = 0;
146 int gid1;
147 Array<OneD, unsigned int> vSizes(nElmt);
148
149 // First construct a map of the number of local DOFs in each block
150 // and the number of matrix entries for each block
151 for (n = 0; n < nElmt; ++n)
152 {
153 loc_mat = m_schurCompl->GetBlock(n, n);
154 vSizes[n] = loc_mat->GetRows();
155 nEntries += vSizes[n] * vSizes[n];
156 }
157
158 // Set up i-index, j-index and value arrays
161 m_Ar = Array<OneD, double>(nEntries, 0.0);
162
163 // Set up the universal ID array for XXT
164 Array<OneD, unsigned long> vId(nLocal);
165
166 // Loop over each elemental block, extract matrix indices and value
167 // and set the universal ID array
168 for (n = iCount = 0; n < nElmt; ++n)
169 {
170 loc_mat = m_schurCompl->GetBlock(n, n);
171 nRows = loc_mat->GetRows();
172
173 for (i = 0; i < nRows; ++i)
174 {
175 gid1 = pLocToGloMap->GetLocalToGlobalBndMap(iCount + i);
176 for (j = 0; j < nRows; ++j)
177 {
178 k = rCount + i * vSizes[n] + j;
179 m_Ai[k] = iCount + i;
180 m_Aj[k] = iCount + j;
181 m_Ar[k] = (*loc_mat)(i, j);
182 if (doSign)
183 {
184 m_Ar[k] *= vMapSign[iCount + i] * vMapSign[iCount + j];
185 }
186 }
187
188 // Dirichlet DOFs are not included in the solve, so we set
189 // these to the special XXT id=0.
190 if (gid1 < numDirBnd)
191 {
192 vId[iCount + i] = 0;
193 }
194 else
195 {
196 vId[iCount + i] =
197 pLocToGloMap->GetGlobalToUniversalBndMap()[gid1];
198 }
199 }
200 iCount += vSizes[n];
201 rCount += vSizes[n] * vSizes[n];
202 }
203
204 // Set up XXT and output some stats
205 LibUtilities::CommSharedPtr vComm = pLocToGloMap->GetComm()->GetRowComm();
206 m_crsData = Xxt::Init(nLocal, vId, m_Ai, m_Aj, m_Ar, vComm);
207 if (m_verbose)
208 {
210 }
211}
212
214 const GlobalLinSysKey &mkey, const std::weak_ptr<ExpList> &pExpList,
215 const DNekScalBlkMatSharedPtr pSchurCompl,
217 const DNekScalBlkMatSharedPtr pInvD,
218 const std::shared_ptr<AssemblyMap> &l2gMap)
219{
222 mkey, pExpList, pSchurCompl, pBinvD, pC, pInvD, l2gMap);
223 sys->Initialise(l2gMap);
224 return sys;
225}
226
227/// Solve the linear system for given input and output vectors.
229 [[maybe_unused]] const int pNumRows,
231 const AssemblyMapSharedPtr &pLocToGloMap,
232 [[maybe_unused]] const int pNumDir)
233{
234 int nLocal = pLocToGloMap->GetLocalToGlobalBndSign().size();
235 Vmath::Zero(nLocal, pOutput, 1);
236
237 if (pLocToGloMap->GetSignChange())
238 {
239 Array<OneD, NekDouble> vlocal(nLocal);
240 Vmath::Vmul(nLocal, pLocToGloMap->GetLocalToGlobalBndSign(), 1, pInput,
241 1, vlocal, 1);
242
243 Xxt::Solve(pOutput, m_crsData, vlocal);
244
245 Vmath::Vmul(nLocal, pLocToGloMap->GetLocalToGlobalBndSign(), 1, pOutput,
246 1, pOutput, 1);
247 }
248 else
249 {
250 Xxt::Solve(pOutput, m_crsData, pInput);
251 }
252}
253
254} // namespace Nektar::MultiRegions
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
A global linear system.
Definition: GlobalLinSys.h:70
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:122
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
DNekScalBlkMatSharedPtr m_schurCompl
Block Schur complement matrix.
std::weak_ptr< AssemblyMap > m_locToGloMap
Local to global map.
DNekScalBlkMatSharedPtr m_BinvD
Block matrix.
DNekScalBlkMatSharedPtr m_C
Block matrix.
DNekScalBlkMatSharedPtr m_invD
Block matrix.
Array< OneD, unsigned int > m_Aj
Array< OneD, unsigned int > m_Ai
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 linear system for given input and output vectors.
GlobalLinSysStaticCondSharedPtr v_Recurse(const GlobalLinSysKey &mkey, const std::weak_ptr< ExpList > &pExpList, const DNekScalBlkMatSharedPtr pSchurCompl, const DNekScalBlkMatSharedPtr pBinvD, const DNekScalBlkMatSharedPtr pC, const DNekScalBlkMatSharedPtr pInvD, const std::shared_ptr< AssemblyMap > &locToGloMap) override
GlobalLinSysXxtStaticCond(const GlobalLinSysKey &mkey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &locToGloMap)
Constructor for full direct matrix solve.
void v_AssembleSchurComplement(std::shared_ptr< AssemblyMap > locToGloMap) override
Assemble the Schur complement matrix.
static GlobalLinSysSharedPtr create(const GlobalLinSysKey &pLinSysKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< GlobalLinSysXxtStaticCond > GlobalLinSysXxtStaticCondSharedPtr
std::shared_ptr< GlobalLinSysStaticCond > GlobalLinSysStaticCondSharedPtr
GlobalLinSysFactory & GetGlobalLinSysFactory()
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
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.
Definition: Vmath.hpp:72
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void nektar_crs_stats(struct crs_data *data)
static struct crs_data * Init(unsigned int pRank, const Nektar::Array< OneD, unsigned long > pId, const Nektar::Array< OneD, unsigned int > pAi, const Nektar::Array< OneD, unsigned int > pAj, const Nektar::Array< OneD, NekDouble > pAr, const LibUtilities::CommSharedPtr &pComm)
Initialise the matrix-solve.
Definition: Xxt.hpp:158
static void Solve(Nektar::Array< OneD, NekDouble > pX, struct crs_data *pCrs, Nektar::Array< OneD, NekDouble > pB)
Solve the matrix system for a given input vector b.
Definition: Xxt.hpp:186