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);
71 ~PreconCfsBRJ() override = default;
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
92 void v_InitObject() override;
93
94 void v_DoPreconCfs(
96 const Array<OneD, NekDouble> &pInput, Array<OneD, NekDouble> &pOutput,
97 const bool &flag) override;
98
101 const Array<OneD, const Array<OneD, NekDouble>> &intmp,
102 const NekDouble time, const NekDouble lambda) override;
103
105 const NekDouble dtLambda) override;
106
107private:
108 void DoNullPrecon(const Array<OneD, NekDouble> &pInput,
109 Array<OneD, NekDouble> &pOutput, const bool &flag);
110
111 void PreconBlkDiag(
113 const Array<OneD, NekDouble> &inarray,
114 Array<OneD, NekDouble> &outarray);
115
116 template <typename DataType>
117 void MinusOffDiag2Rhs(
119 const size_t nvariables, const size_t nCoeffs,
120 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
123 Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
124 const TensorOfArray4D<DataType> &TraceJacArray);
125
126 template <typename TypeNekBlkMatSharedPtr>
130 const int &nscale = 1);
131
132 /// This function creates the matrix structure for the block
133 /// diagonal operator. It organizes the way that the matrix
134 // is loaded, multiplied and unpacked in order to be operated
135 /// by SIMD instructions. The degrees of freedom of each element
136 /// are reorganized, so they are placed in the correct location
137 /// to perfom SIMD instructions.
140 {
141 using vec_t = simd<NekSingle>;
142
143 int TotMatLen = 0;
144 int TotLen = 0;
145 const auto nTotElmt = pFields[0]->GetNumElmts();
146 const auto nvariables = pFields.size();
147 const auto vecwidth = vec_t::width;
148
149 m_max_nblocks = 0;
150 m_max_nElmtDof = 0;
151
152 for (int ne = 0; ne < nTotElmt; ne++)
153 {
154 const auto nElmtDof = pFields[0]->GetNcoeffs(ne) * nvariables;
155 const auto nblocks = nElmtDof / vecwidth;
156 unsigned int totblocks =
157 (nElmtDof % vecwidth) ? nblocks + 1 : nblocks;
158
160 (m_max_nblocks > totblocks) ? m_max_nblocks : totblocks;
162 (m_max_nElmtDof > nElmtDof) ? m_max_nElmtDof : nElmtDof;
163
164 TotLen += totblocks * vecwidth;
165 TotMatLen += nElmtDof * totblocks;
166 }
167
168 m_sBlkDiagMat.resize(TotMatLen);
169 m_inputIdx.resize(TotLen);
170
171 // generate a index list of vector width aware mapping from
172 // local coeff storage over all variables to elemental storage
173 // over variables
174 unsigned int ncoeffs = pFields[0]->GetNcoeffs();
175 for (int ne = 0, cnt1 = 0; ne < nTotElmt; ne++)
176 {
177 const auto nElmtCoeff = pFields[0]->GetNcoeffs(ne);
178 const auto nElmtDof = nElmtCoeff * nvariables;
179 const auto nblocks = nElmtDof / vecwidth;
180 const auto nCoefOffset = pFields[0]->GetCoeff_Offset(ne);
181 int i = 0;
182 int i0 = 0;
183 int inOffset, j;
184
185 for (int m = 0; m < nvariables; m++)
186 {
187 inOffset = m * ncoeffs + nCoefOffset;
188
189 if (m && (vecwidth - i0 < nElmtCoeff))
190 {
191 // May need to add entries from later variables to
192 // remainder of last variable if the vector width
193 // was not exact multiple of number of elemental
194 // coeffs
195 for (i = 0; i0 < vecwidth; ++i, ++i0)
196 {
197 m_inputIdx[cnt1++] = inOffset + i;
198 }
199 }
200 else
201 {
202 i = 0;
203 }
204
205 // load up other vectors in variable that fit into vector
206 // width
207 for (j = 0; (j + 1) * vecwidth < nElmtCoeff - i; ++j)
208 {
209 for (i0 = 0; i0 < vecwidth; ++i0)
210 {
211 m_inputIdx[cnt1++] = inOffset + i + j * vecwidth + i0;
212 }
213 }
214
215 // load up any residual data for this variable
216 for (i0 = 0, j = j * vecwidth; j < nElmtCoeff - i; ++j, ++i0)
217 {
218 m_inputIdx[cnt1++] = inOffset + i + j;
219 }
220 }
221
222 const auto endwidth = nElmtDof - nblocks * vecwidth;
223
224 // fill out rest of index to match vector width with last entry
225 if (endwidth)
226 {
227 for (i0 = endwidth; i0 < vecwidth; ++i0)
228 {
229 m_inputIdx[cnt1++] = inOffset + i + j - 1;
230 }
231 }
232 ASSERTL1(cnt1 <= TotLen, "m_inputIdx over extended");
233 }
234 }
235
237 const Array<OneD, unsigned int> nrow,
238 const Array<OneD, unsigned int> ncol)
239 {
240 mat =
242 SNekMatSharedPtr loc_matNvar;
243 for (size_t nelm = 0; nelm < nrow.size(); ++nelm)
244 {
245 int nrowsVars = nrow[nelm];
246 int ncolsVars = ncol[nelm];
247
249 nrowsVars, ncolsVars, 0.0);
250 mat->SetBlock(nelm, nelm, loc_matNvar);
251 }
252 }
253};
254} // namespace Nektar
255
256#endif
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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 v_BuildPreconCfs(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, const Array< OneD, NekDouble > > &intmp, const NekDouble time, const NekDouble lambda) override
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
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, TensorOfArray3D< NekDouble > &wspTrace, Array< OneD, Array< OneD, DataType > > &wspTraceDataType, const TensorOfArray4D< DataType > &TraceJacArray)
Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > m_PreconMatVarsSingle
Definition: PreconCfsBRJ.h:77
~PreconCfsBRJ() override=default
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
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:236
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:138
void PreconBlkDiag(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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:55
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
std::shared_ptr< SNekMat > SNekMatSharedPtr
std::shared_ptr< PreconCfs > PreconCfsSharedPtr
Definition: PreconCfs.h:112
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