Nektar++
GlobalLinSysXxtFull.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: GlobalLinSysXxtFull.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: GlobalLinSysXxtFull definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
38
39using namespace std;
40
42{
43/**
44 * @class GlobalLinSysXxtFull
45 */
46
47/**
48 * Registers the class with the Factory.
49 */
52 "XxtFull", GlobalLinSysXxtFull::create, "Xxt Full Matrix.");
53
54/// Constructor for full direct matrix solve.
56 const GlobalLinSysKey &pLinSysKey, const std::weak_ptr<ExpList> &pExp,
57 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
58 : GlobalLinSys(pLinSysKey, pExp, pLocToGloMap),
59 GlobalLinSysXxt(pLinSysKey, pExp, pLocToGloMap)
60{
61
63 "This routine should only be used when using a Full XXT"
64 " matrix solve");
65
66 AssembleMatrixArrays(pLocToGloMap);
67}
68
69/**
70 * Construct the local matrix row index, column index and value index
71 * arrays and initialize the XXT data structure with this information.
72 * @param locToGloMap Local to global mapping information.
73 */
75 const std::shared_ptr<AssemblyMap> &pLocToGloMap)
76{
77 ExpListSharedPtr vExp = m_expList.lock();
78 unsigned int nElmt = vExp->GetNumElmts();
80 unsigned int iCount = 0;
81 unsigned int rCount = 0;
82 unsigned int nRows = 0;
83 unsigned int nEntries = 0;
84 unsigned int numDirBnd = pLocToGloMap->GetNumGlobalDirBndCoeffs();
85 unsigned int nLocal = pLocToGloMap->GetNumLocalCoeffs();
86 const Array<OneD, NekDouble> &vMapSign =
87 pLocToGloMap->GetLocalToGlobalSign();
88 bool doSign = pLocToGloMap->GetSignChange();
89 unsigned int i = 0, j = 0, k = 0, n = 0;
90 int gid1;
91 Array<OneD, unsigned int> vSizes(nElmt);
92
93 // First construct a map of the number of local DOFs in each block
94 // and the number of matrix entries for each block
95
96 // Dimension of matrix is just the linear vertex space
99 {
100 for (n = 0; n < nElmt; ++n)
101 {
102 vSizes[n] = vExp->GetExp(n)->GetNverts();
103 nEntries += vSizes[n] * vSizes[n];
104 }
105 }
106 else
107 {
108 for (n = 0; n < nElmt; ++n)
109 {
110 vSizes[n] = vExp->GetExp(n)->GetNcoeffs();
111 nEntries += vSizes[n] * vSizes[n];
112 }
113 }
114
115 // Set up i-index, j-index and value arrays
118 m_Ar = Array<OneD, double>(nEntries, 0.0);
119
120 // Set up the universal ID array for XXT
121 Array<OneD, unsigned long> vId(nLocal);
122
123 // Loop over each elemental block, extract matrix indices and value
124 // and set the universal ID array
125 for (n = iCount = 0; n < nElmt; ++n)
126 {
127 loc_mat = GetBlock(n);
128 nRows = loc_mat->GetRows();
129
130 for (i = 0; i < nRows; ++i)
131 {
132 gid1 = pLocToGloMap->GetLocalToGlobalMap(iCount + i);
133 for (j = 0; j < nRows; ++j)
134 {
135 k = rCount + i * vSizes[n] + j;
136 m_Ai[k] = iCount + i;
137 m_Aj[k] = iCount + j;
138 m_Ar[k] = (*loc_mat)(i, j);
139 if (doSign)
140 {
141 m_Ar[k] *= vMapSign[iCount + i] * vMapSign[iCount + j];
142 }
143 }
144
145 // Dirichlet DOFs are not included in the solve, so we set
146 // these to the special XXT id=0.
147 if (gid1 < numDirBnd)
148 {
149 vId[iCount + i] = 0;
150 }
151 else
152 {
153 vId[iCount + i] = pLocToGloMap->GetGlobalToUniversalMap(gid1);
154 }
155 }
156 iCount += vSizes[n];
157 rCount += vSizes[n] * vSizes[n];
158 }
159
160 // Set up XXT and output some stats
161 LibUtilities::CommSharedPtr vComm = pLocToGloMap->GetComm();
162 m_crsData = Xxt::Init(nLocal, vId, m_Ai, m_Aj, m_Ar, vComm);
163 if (m_verbose)
164 {
166 }
167}
168
169/**
170 * Solve the linear system using a full global matrix system.
171 */
173 const Array<OneD, const NekDouble> &pLocInput,
174 Array<OneD, NekDouble> &pLocOutput,
175 const AssemblyMapSharedPtr &pLocToGloMap,
176 const Array<OneD, const NekDouble> &pDirForcing)
177{
178 bool dirForcCalculated = (bool)pDirForcing.size();
179 int nDirDofs = pLocToGloMap->GetNumGlobalDirBndCoeffs();
180 int nLocDofs = pLocToGloMap->GetNumLocalCoeffs();
181
182 if (nDirDofs)
183 {
184 std::shared_ptr<MultiRegions::ExpList> expList = m_expList.lock();
185 Array<OneD, NekDouble> rhs(nLocDofs);
186
187 // Calculate the Dirichlet forcing
188 if (dirForcCalculated)
189 {
190 // Assume pDirForcing is in local space
191 ASSERTL0(
192 pDirForcing.size() >= nLocDofs,
193 "DirForcing is not of sufficient size. Is it in local space?");
194 Vmath::Vsub(nLocDofs, pLocInput, 1, pDirForcing, 1, rhs, 1);
195 }
196 else
197 {
198 // Calculate initial condition and Dirichlet forcing and subtract it
199 // from the rhs
200 expList->GeneralMatrixOp(m_linSysKey, pLocOutput, rhs);
201
202 // Iterate over all the elements computing Robin BCs where
203 // necessary
204 for (auto &r : m_robinBCInfo) // add robin mass matrix
205 {
208
209 int n = r.first;
210 int offset = expList->GetCoeff_Offset(n);
211
212 LocalRegions::ExpansionSharedPtr vExp = expList->GetExp(n);
213 // Add local matrix contribution
214 for (rBC = r.second; rBC; rBC = rBC->next)
215 {
216 vExp->AddRobinTraceContribution(
217 rBC->m_robinID, rBC->m_robinPrimitiveCoeffs,
218 pLocOutput + offset, rhsloc = rhs + offset);
219 }
220 }
221 Vmath::Vsub(nLocDofs, pLocInput, 1, rhs, 1, rhs, 1);
222 }
223
224 Array<OneD, NekDouble> diff(nLocDofs);
225
226 // Solve for perturbation from initial guess in pOutput
227 SolveLinearSystem(nLocDofs, rhs, diff, pLocToGloMap);
228
229 // Add back initial and boundary condition
230 Vmath::Vadd(nLocDofs, diff, 1, pLocOutput, 1, pLocOutput, 1);
231 }
232 else
233 {
234 SolveLinearSystem(nLocDofs, pLocInput, pLocOutput, pLocToGloMap);
235 }
236}
237
238/// Solve the linear system for given input and output vectors.
240 [[maybe_unused]] const int pNumRows,
242 const AssemblyMapSharedPtr &pLocToGloMap,
243 [[maybe_unused]] const int pNumDir)
244{
245 int nLocal = pNumRows;
246
247 Vmath::Zero(nLocal, pOutput, 1);
248
249 // Set Output into correct sign
250 if (pLocToGloMap->GetSignChange())
251 {
252 Array<OneD, NekDouble> vlocal(nLocal);
253 Vmath::Vmul(nLocal, pLocToGloMap->GetLocalToGlobalSign(), 1, pInput, 1,
254 vlocal, 1);
255
256 Xxt::Solve(pOutput, m_crsData, vlocal);
257
258 Vmath::Vmul(nLocal, pLocToGloMap->GetLocalToGlobalSign(), 1, pOutput, 1,
259 pOutput, 1);
260 }
261 else
262 {
263 Xxt::Solve(pOutput, m_crsData, pInput);
264 }
265}
266
267} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#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.
A global linear system.
Definition: GlobalLinSys.h:70
const std::weak_ptr< ExpList > m_expList
Local Matrix System.
Definition: GlobalLinSys.h:122
const std::map< int, RobinBCInfoSharedPtr > m_robinBCInfo
Robin boundary info.
Definition: GlobalLinSys.h:124
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.
Definition: GlobalLinSys.h:190
const GlobalLinSysKey m_linSysKey
Key associated with this linear system.
Definition: GlobalLinSys.h:120
DNekScalMatSharedPtr GetBlock(unsigned int n)
Definition: GlobalLinSys.h:209
GlobalSysSolnType GetGlobalSysSolnType() const
Return the associated solution type.
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=0) override
Solve the linear system for given input and output vectors.
static GlobalLinSysSharedPtr create(const GlobalLinSysKey &pLinSysKey, const std::weak_ptr< ExpList > &pExpList, const std::shared_ptr< AssemblyMap > &pLocToGloMap)
Creates an instance of this class.
void AssembleMatrixArrays(const std::shared_ptr< AssemblyMap > &pLocToGloMap)
GlobalLinSysXxtFull(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.
Array< OneD, unsigned int > m_Aj
Array< OneD, unsigned int > m_Ai
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< RobinBCInfo > RobinBCInfoSharedPtr
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
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 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.
Definition: Vmath.hpp:180
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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.
Definition: Vmath.hpp:220
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
STL namespace.