Nektar++
PreconCfsBRJ.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PreconCfsBRJ.h
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: PreconCfsBRJ header
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_PRECONCFSBRJ
36#define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_PRECONCFSBRJ
37
39
40namespace Nektar
41{
42
43using namespace tinysimd;
44
45/**
46 * Block Relaxed(weighted) Jacobi iterative (BRJ) Preconditioner for CFS
47 *
48 */
49class PreconCfsBRJ : public PreconCfs
50{
51public:
52 friend class MemoryManager<PreconCfsBRJ>;
53
54 /// Creates an instance of this class
58 const LibUtilities::CommSharedPtr &vComm)
59 {
61 pFields, pSession, vComm);
62 return p;
63 }
64
65 /// Name of the class
66 static std::string className;
67
70 const LibUtilities::CommSharedPtr &vComm);
72
73protected:
76
78
79 unsigned int m_max_nblocks;
80 unsigned int m_max_nElmtDof;
81 std::vector<simd<NekSingle>, tinysimd::allocator<simd<NekSingle>>>
83 std::vector<int> m_inputIdx;
84
91
93
94 virtual void v_InitObject() override;
95
96 virtual void v_DoPreconCfs(
98 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
99 const bool &flag) override;
100
101 virtual void v_BuildPreconCfs(
103 const Array<OneD, const Array<OneD, NekDouble>> &intmp,
104 const NekDouble time, const NekDouble lambda) override;
105
107 const NekDouble dtLambda) override;
108
109private:
110 void DoNullPrecon(const Array<OneD, NekDouble> &pInput,
111 Array<OneD, NekDouble> &pOutput, const bool &flag);
112
113 void PreconBlkDiag(
115 const Array<OneD, NekDouble> &inarray,
116 Array<OneD, NekDouble> &outarray);
117
118 template <typename DataType>
119 void MinusOffDiag2Rhs(
121 const size_t nvariables, const size_t nCoeffs,
122 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
123 Array<OneD, Array<OneD, NekDouble>> &outarray, bool flagUpdateDervFlux,
124 Array<OneD, Array<OneD, NekDouble>> &FwdFluxDeriv,
125 Array<OneD, Array<OneD, NekDouble>> &BwdFluxDeriv,
128 Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
129 const TensorOfArray4D<DataType> &TraceJacArray,
130 const TensorOfArray4D<DataType> &TraceJacDerivArray,
131 const Array<OneD, const Array<OneD, DataType>> &TraceJacDerivSign,
132 const TensorOfArray5D<DataType> &TraceIPSymJacArray);
133
134 template <typename TypeNekBlkMatSharedPtr>
138 const int &nscale = 1);
139
140 /// This function creates the matrix structure for the block
141 /// diagonal operator. It organizes the way that the matrix
142 // is loaded, multiplied and unpacked in order to be operated
143 /// by SIMD instructions. The degrees of freedom of each element
144 /// are reorganized, so they are placed in the correct location
145 /// to perfom SIMD instructions.
148 {
149 using vec_t = simd<NekSingle>;
150
151 int TotMatLen = 0;
152 int TotLen = 0;
153 const auto nTotElmt = pFields[0]->GetNumElmts();
154 const auto nvariables = pFields.size();
155 const auto vecwidth = vec_t::width;
156
157 m_max_nblocks = 0;
158 m_max_nElmtDof = 0;
159
160 for (int ne = 0; ne < nTotElmt; ne++)
161 {
162 const auto nElmtDof = pFields[0]->GetNcoeffs(ne) * nvariables;
163 const auto nblocks = nElmtDof / vecwidth;
164 unsigned int totblocks =
165 (nElmtDof % vecwidth) ? nblocks + 1 : nblocks;
166
168 (m_max_nblocks > totblocks) ? m_max_nblocks : totblocks;
170 (m_max_nElmtDof > nElmtDof) ? m_max_nElmtDof : nElmtDof;
171
172 TotLen += totblocks * vecwidth;
173 TotMatLen += nElmtDof * totblocks;
174 }
175
176 m_sBlkDiagMat.resize(TotMatLen);
177 m_inputIdx.resize(TotLen);
178
179 // generate a index list of vector width aware mapping from
180 // local coeff storage over all variables to elemental storage
181 // over variables
182 unsigned int ncoeffs = pFields[0]->GetNcoeffs();
183 for (int ne = 0, cnt1 = 0; ne < nTotElmt; ne++)
184 {
185 const auto nElmtCoeff = pFields[0]->GetNcoeffs(ne);
186 const auto nElmtDof = nElmtCoeff * nvariables;
187 const auto nblocks = nElmtDof / vecwidth;
188 const auto nCoefOffset = pFields[0]->GetCoeff_Offset(ne);
189 int i = 0;
190 int i0 = 0;
191 int inOffset, j;
192
193 for (int m = 0; m < nvariables; m++)
194 {
195 inOffset = m * ncoeffs + nCoefOffset;
196
197 if (m && (vecwidth - i0 < nElmtCoeff))
198 {
199 // May need to add entries from later variables to
200 // remainder of last variable if the vector width
201 // was not exact multiple of number of elemental
202 // coeffs
203 for (i = 0; i0 < vecwidth; ++i, ++i0)
204 {
205 m_inputIdx[cnt1++] = inOffset + i;
206 }
207 }
208 else
209 {
210 i = 0;
211 }
212
213 // load up other vectors in variable that fit into vector
214 // width
215 for (j = 0; (j + 1) * vecwidth < nElmtCoeff - i; ++j)
216 {
217 for (i0 = 0; i0 < vecwidth; ++i0)
218 {
219 m_inputIdx[cnt1++] = inOffset + i + j * vecwidth + i0;
220 }
221 }
222
223 // load up any residual data for this variable
224 for (i0 = 0, j = j * vecwidth; j < nElmtCoeff - i; ++j, ++i0)
225 {
226 m_inputIdx[cnt1++] = inOffset + i + j;
227 }
228 }
229
230 const auto endwidth = nElmtDof - nblocks * vecwidth;
231
232 // fill out rest of index to match vector width with last entry
233 if (endwidth)
234 {
235 for (i0 = endwidth; i0 < vecwidth; ++i0)
236 {
237 m_inputIdx[cnt1++] = inOffset + i + j - 1;
238 }
239 }
240 ASSERTL1(cnt1 <= TotLen, "m_inputIdx over extended");
241 }
242 }
243
245 const Array<OneD, unsigned int> nrow,
246 const Array<OneD, unsigned int> ncol)
247 {
248 mat =
250 SNekMatSharedPtr loc_matNvar;
251 for (size_t nelm = 0; nelm < nrow.size(); ++nelm)
252 {
253 int nrowsVars = nrow[nelm];
254 int ncolsVars = ncol[nelm];
255
257 nrowsVars, ncolsVars, 0.0);
258 mat->SetBlock(nelm, nelm, loc_matNvar);
259 }
260 }
261};
262} // namespace Nektar
263
264#endif
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacSingle
Definition: PreconCfsBRJ.h:85
void MinusOffDiag2Rhs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const size_t nvariables, const size_t nCoeffs, const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, bool flagUpdateDervFlux, Array< OneD, Array< OneD, NekDouble > > &FwdFluxDeriv, Array< OneD, Array< OneD, NekDouble > > &BwdFluxDeriv, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &wspTrace, Array< OneD, Array< OneD, DataType > > &wspTraceDataType, const TensorOfArray4D< DataType > &TraceJacArray, const TensorOfArray4D< DataType > &TraceJacDerivArray, const Array< OneD, const Array< OneD, DataType > > &TraceJacDerivSign, const TensorOfArray5D< DataType > &TraceIPSymJacArray)
virtual void v_BuildPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble > > &intmp, const NekDouble time, const NekDouble lambda) override
virtual void v_DoPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag) override
TensorOfArray4D< NekSingle > m_TraceJacArraySingle
Definition: PreconCfsBRJ.h:86
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:77
static std::string className
Name of the class.
Definition: PreconCfsBRJ.h:66
Array< OneD, SNekBlkMatSharedPtr > m_TraceJacDerivSingle
Definition: PreconCfsBRJ.h:87
unsigned int m_max_nblocks
Definition: PreconCfsBRJ.h:79
Array< OneD, Array< OneD, NekSingle > > m_TraceJacDerivSignSingle
Definition: PreconCfsBRJ.h:89
static PreconCfsSharedPtr create(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
Creates an instance of this class.
Definition: PreconCfsBRJ.h:55
PrecType m_PreconMatStorage
Definition: PreconCfsBRJ.h:92
virtual void v_InitObject() override
TensorOfArray4D< NekSingle > m_TraceJacDerivArraySingle
Definition: PreconCfsBRJ.h:88
PreconCfsBRJ(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, const LibUtilities::CommSharedPtr &vComm)
unsigned int m_max_nElmtDof
Definition: PreconCfsBRJ.h:80
void AllocatePreconBlkDiagCoeff(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const int &nscale=1)
TensorOfArray5D< NekSingle > m_TraceIPSymJacArraySingle
Definition: PreconCfsBRJ.h:90
std::vector< int > m_inputIdx
Definition: PreconCfsBRJ.h:83
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
Definition: PreconCfsBRJ.h:244
void AllocateSIMDPreconBlkMatDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
This function creates the matrix structure for the block diagonal operator. It organizes the way that...
Definition: PreconCfsBRJ.h:146
void PreconBlkDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual bool v_UpdatePreconMatCheck(const Array< OneD, const NekDouble > &res, const NekDouble dtLambda) override
std::vector< simd< NekSingle >, tinysimd::allocator< simd< NekSingle > > > m_sBlkDiagMat
Definition: PreconCfsBRJ.h:82
void DoNullPrecon(const Array< OneD, NekDouble > &pInput, Array< OneD, NekDouble > &pOutput, const bool &flag)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
std::shared_ptr< SNekMat > SNekMatSharedPtr
PrecType
Definition: PreconCfs.h:47
std::shared_ptr< PreconCfs > PreconCfsSharedPtr
Definition: PreconCfs.h:126
tinysimd::simd< NekDouble > vec_t
double NekDouble
typename abi< ScalarType, width >::type simd
Definition: tinysimd.hpp:80
boost::alignment::aligned_allocator< T, T::alignment > allocator
Definition: allocator.hpp:49