Nektar++
PreconCfsBRJ.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PreconCfsBRJ.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// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: PreconCfsBRJ definition
33//
34///////////////////////////////////////////////////////////////////////////////
35
38
39using namespace std;
40
41namespace Nektar
42{
43/**
44 * @class PreconCfsBRJ
45 *
46 * Solves a linear system using iterative methods.
47 */
48std::string PreconCfsBRJ::className =
50 "PreconCfsBRJ", PreconCfsBRJ::create,
51 "Block Relaxed Jacobi Preconditioner for CFS.");
52
56 const LibUtilities::CommSharedPtr &vComm)
57 : PreconCfs(pFields, pSession, vComm)
58{
59 pSession->LoadParameter("PreconItsStep", m_PreconItsStep, 7);
60 pSession->LoadParameter("BRJRelaxParam", m_BRJRelaxParam, 1.0);
61
62 size_t nvariables = pFields.size();
63
66 for (size_t i = 0; i < nvariables; i++)
67 {
69 }
71
73}
74
75/**
76 *
77 */
79{
80}
81
82/**
83 *
84 */
87 const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray,
88 [[maybe_unused]] const bool &flag)
89{
90 size_t nBRJIterTot = m_PreconItsStep;
91 if (nBRJIterTot == 0)
92 {
93 DoNullPrecon(inarray, outarray, flag);
94 }
95 else
96 {
97 const NekDouble BRJParam = m_BRJRelaxParam;
98 const NekDouble OmBRJParam = 1.0 - BRJParam;
99
100 size_t nvariables = pFields.size();
101 size_t npoints = pFields[0]->GetNcoeffs();
102 size_t ntotpnt = inarray.size();
103
104 ASSERTL0(nvariables * npoints == ntotpnt,
105 "nvariables*npoints!=ntotpnt in PreconCoeff");
106
107 Array<OneD, NekDouble> rhs(ntotpnt);
108
109 Array<OneD, NekDouble> outN(ntotpnt);
110 Array<OneD, NekDouble> outTmp(ntotpnt);
111 Array<OneD, Array<OneD, NekDouble>> rhs2d(nvariables);
112 Array<OneD, Array<OneD, NekDouble>> out_2d(nvariables);
113 Array<OneD, Array<OneD, NekDouble>> outTmp_2d(nvariables);
114 for (size_t m = 0; m < nvariables; m++)
115 {
116 size_t moffset = m * npoints;
117 rhs2d[m] = rhs + moffset;
118 out_2d[m] = outarray + moffset;
119 outTmp_2d[m] = outTmp + moffset;
120 pFields[m]->MultiplyByMassMatrix(inarray + moffset, rhs2d[m]);
121 }
122
123 size_t nphysic = pFields[0]->GetNpoints();
124 size_t nTracePts = pFields[0]->GetTrace()->GetNpoints();
126 for (int i = 0; i < m_spacedim; i++)
127 {
128 qfield[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
129 for (size_t j = 0; j < nvariables; j++)
130 {
131 qfield[i][j] = Array<OneD, NekDouble>(nphysic);
132 }
133 }
134 size_t ntmpTrace = 4 + 2 * m_spacedim;
135 TensorOfArray3D<NekDouble> tmpTrace(ntmpTrace);
136 for (size_t i = 0; i < ntmpTrace; i++)
137 {
138 tmpTrace[i] = Array<OneD, Array<OneD, NekDouble>>(nvariables);
139 for (size_t j = 0; j < nvariables; j++)
140 {
141 tmpTrace[i][j] = Array<OneD, NekDouble>(nTracePts);
142 }
143 }
144 Array<OneD, Array<OneD, NekDouble>> FwdFluxDeriv(nvariables);
145 Array<OneD, Array<OneD, NekDouble>> BwdFluxDeriv(nvariables);
146 for (size_t j = 0; j < nvariables; j++)
147 {
148 FwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
149 BwdFluxDeriv[j] = Array<OneD, NekDouble>(nTracePts);
150 }
151
152 const size_t nwspTraceDataType = nvariables + 1;
153 Array<OneD, Array<OneD, NekSingle>> wspTraceDataType(nwspTraceDataType);
154 for (size_t m = 0; m < nwspTraceDataType; m++)
155 {
156 wspTraceDataType[m] = Array<OneD, NekSingle>(nTracePts);
157 }
158
160 timer.Start();
161 PreconBlkDiag(pFields, rhs, outarray);
162 timer.Stop();
163 timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
164
165 for (size_t nrelax = 0; nrelax < nBRJIterTot - 1; nrelax++)
166 {
167 Vmath::Smul(ntotpnt, OmBRJParam, outarray, 1, outN, 1);
168
169 timer.Start();
170 MinusOffDiag2Rhs(pFields, nvariables, npoints, rhs2d, out_2d,
171 tmpTrace, wspTraceDataType, m_TraceJacArraySingle);
172 timer.Stop();
173 timer.AccumulateRegion("PreconCfsBRJ::MinusOffDiag2Rhs", 2);
174
175 timer.Start();
176 PreconBlkDiag(pFields, outarray, outTmp);
177 timer.Stop();
178 timer.AccumulateRegion("PreconCfsBRJ::PreconBlkDiag", 2);
179
180 Vmath::Svtvp(ntotpnt, BRJParam, outTmp, 1, outN, 1, outarray, 1);
181 }
182 }
183}
184
185/**
186 *
187 */
189 [[maybe_unused]] const Array<OneD, MultiRegions::ExpListSharedPtr> &pFields,
190 const Array<OneD, const Array<OneD, NekDouble>> &intmp,
191 [[maybe_unused]] const NekDouble time, const NekDouble lambda)
192{
193 if (m_PreconItsStep > 0)
194 {
195 SNekBlkMatSharedPtr PreconMatSingle;
196 using vec_t = simd<NekSingle>;
197 int nvariables = pFields.size();
198 int nelmts = pFields[0]->GetNumElmts();
199 Array<OneD, unsigned int> matdim(nelmts);
200 for (int i = 0; i < nelmts; i++)
201 {
202 matdim[i] = pFields[0]->GetExp(i)->GetNcoeffs() * nvariables;
203 }
204 AllocateNekBlkMatDig(PreconMatSingle, matdim, matdim);
205
207 intmp, m_PreconMatVarsSingle, PreconMatSingle, m_TraceJacSingle,
211
212 if (m_verbose && m_Comm->GetRank() == 0)
213 {
214 cout << " ## CalcuPreconMat " << endl;
215 }
216
217 // copy matrix to simd layout
218 // load matrix
219 int cnt = 0;
220 const auto vecwidth = vec_t::width;
221
222 alignas(vec_t::alignment) std::array<NekSingle, vec_t::width> tmp;
223
224 for (int ne = 0; ne < nelmts; ne++)
225 {
226 const auto nElmtDof = matdim[ne];
227 const auto nblocks = nElmtDof / vecwidth;
228
229 const NekSingle *mmat =
230 PreconMatSingle->GetBlockPtr(ne, ne)->GetRawPtr();
231 /// Copy array into column major blocks of vector width
232 for (int i1 = 0; i1 < nblocks; ++i1)
233 {
234 for (int j = 0; j < nElmtDof; ++j)
235 {
236 for (int i = 0; i < vecwidth; ++i)
237 {
238 tmp[i] = mmat[j + (i1 * vecwidth + i) * nElmtDof];
239 }
240 // store contiguous vec_t array.
241 m_sBlkDiagMat[cnt++].load(tmp.data());
242 }
243 }
244
245 const auto endwidth = nElmtDof - nblocks * vecwidth;
246
247 // end rows that do not fit into vector widths
248 if (endwidth)
249 {
250 for (int j = 0; j < nElmtDof; ++j)
251 {
252 for (int i = 0; i < endwidth; ++i)
253 {
254 tmp[i] = mmat[j + (nblocks * vecwidth + i) * nElmtDof];
255 }
256
257 for (int i = endwidth; i < vecwidth; ++i)
258 {
259 tmp[i] = 0.0;
260 }
261 m_sBlkDiagMat[cnt++].load(tmp.data());
262 }
263 }
264 }
265 }
266
267 m_DtLambdaPreconMat = lambda;
268 m_CalcPreconMatFlag = false;
270}
271
272/**
273 *
274 */
276 [[maybe_unused]] const Array<OneD, const NekDouble> &res,
277 const NekDouble dtLambda)
278{
279 bool flag = (m_CalcPreconMatFlag || m_DtLambdaPreconMat != dtLambda ||
281 m_CalcPreconMatFlag = flag;
282 return flag;
283}
284
285/**
286 *
287 */
289 Array<OneD, NekDouble> &pOutput,
290 [[maybe_unused]] const bool &flag)
291{
292 Vmath::Vcopy(pInput.size(), pInput, 1, pOutput, 1);
293}
294
295/**
296 *
297 */
300 const Array<OneD, NekDouble> &inarray, Array<OneD, NekDouble> &outarray)
301{
302 unsigned int nvariables = pFields.size();
303
304 int nTotElmt = pFields[0]->GetNumElmts();
305
306 using vec_t = simd<NekSingle>;
307 const auto vecwidth = vec_t::width;
308
309 // vectorized matrix multiply
310 std::vector<vec_t, tinysimd::allocator<vec_t>> Sinarray(m_max_nblocks);
311 std::vector<vec_t, tinysimd::allocator<vec_t>> Soutarray(m_max_nElmtDof);
312
313 alignas(vec_t::alignment) std::array<NekSingle, vec_t::width> tmp;
314
315 for (int ne = 0, cnt = 0, icnt = 0, icnt1 = 0; ne < nTotElmt; ne++)
316 {
317 const auto nElmtCoef = pFields[0]->GetNcoeffs(ne);
318 const auto nElmtDof = nElmtCoef * nvariables;
319 const auto nblocks = (nElmtDof % vecwidth) ? nElmtDof / vecwidth + 1
320 : nElmtDof / vecwidth;
321
322 // gather data into blocks - could probably be done with a
323 // gather call? can be replaced with a gather op when working
324 for (int j = 0; j < nblocks; ++j, icnt += vecwidth)
325 {
326 for (int i = 0; i < vecwidth; ++i)
327 {
328 tmp[i] = inarray[m_inputIdx[icnt + i]];
329 }
330
331 Sinarray[j].load(tmp.data());
332 }
333
334 // Do matrix multiply
335 // first block just needs multiplying
336 vec_t in = Sinarray[0];
337 for (int i = 0; i < nElmtDof; ++i)
338 {
339 Soutarray[i] = m_sBlkDiagMat[cnt++] * in;
340 }
341
342 // rest of blocks are multiply add operations;
343 for (int n = 1; n < nblocks; ++n)
344 {
345 in = Sinarray[n];
346 for (int i = 0; i < nElmtDof; ++i)
347 {
348 Soutarray[i].fma(m_sBlkDiagMat[cnt++], in);
349 }
350 }
351
352 // get block aligned index for this expansion
353 NekSingle val;
354 for (int i = 0; i < nElmtDof; ++i)
355 {
356 // Get hold of data
357 Soutarray[i].store(tmp.data());
358
359 // Sum vector width
360 val = tmp[0];
361 for (int j = 1; j < vecwidth; ++j)
362 {
363 val += tmp[j];
364 }
365 // put data into outarray
366 outarray[m_inputIdx[icnt1 + i]] = NekDouble(val);
367 }
368
369 icnt1 += nblocks * vecwidth;
370 }
371}
372
373/**
374 *
375 */
376template <typename DataType>
379 const size_t nvariables, const size_t nCoeffs,
380 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
383 Array<OneD, Array<OneD, DataType>> &wspTraceDataType,
384 const TensorOfArray4D<DataType> &TraceJacArray)
385{
386 size_t nTracePts = pFields[0]->GetTrace()->GetNpoints();
387 size_t npoints = pFields[0]->GetNpoints();
388 size_t nDim = m_spacedim;
389
390 Array<OneD, Array<OneD, NekDouble>> outpnts(nvariables);
391 for (size_t i = 0; i < nvariables; i++)
392 {
393 outpnts[i] = Array<OneD, NekDouble>(npoints, 0.0);
394 pFields[i]->BwdTrans(outarray[i], outpnts[i]);
395 }
396
397 // Store forwards/backwards space along trace space
402 TensorOfArray3D<NekDouble> numDerivBwd(nDim);
403 TensorOfArray3D<NekDouble> numDerivFwd(nDim);
404 size_t indexwspTrace = 0;
405 Fwd = wspTrace[indexwspTrace], indexwspTrace++;
406 Bwd = wspTrace[indexwspTrace], indexwspTrace++;
407 FwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
408 BwdFlux = wspTrace[indexwspTrace], indexwspTrace++;
409
411 for (size_t i = 0; i < nvariables; ++i)
412 {
413 timer.Start();
414 pFields[i]->GetFwdBwdTracePhys(outpnts[i], Fwd[i], Bwd[i]);
415 timer.Stop();
416 timer.AccumulateRegion("ExpList::GetFwdBwdTracePhys", 10);
417 }
418
419 size_t indexwspTraceDataType = 0;
420 Array<OneD, Array<OneD, DataType>> Fwdarray(nvariables);
421 for (size_t m = 0; m < nvariables; ++m)
422 {
423 Fwdarray[m] = wspTraceDataType[indexwspTraceDataType],
424 indexwspTraceDataType++;
425 }
426 Array<OneD, DataType> Fwdreslt;
427 Fwdreslt = wspTraceDataType[indexwspTraceDataType], indexwspTraceDataType++;
428
429 for (size_t m = 0; m < nvariables; ++m)
430 {
431 for (size_t i = 0; i < nTracePts; ++i)
432 {
433 Fwdarray[m][i] = DataType(Fwd[m][i]);
434 }
435 }
436 for (size_t m = 0; m < nvariables; ++m)
437 {
438 Vmath::Zero(nTracePts, Fwdreslt, 1);
439 for (size_t n = 0; n < nvariables; ++n)
440 {
441 for (size_t p = 0; p < nTracePts; ++p)
442 {
443 Fwdreslt[p] += TraceJacArray[0][m][n][p] * Fwdarray[n][p];
444 }
445 }
446 for (size_t i = 0; i < nTracePts; ++i)
447 {
448 FwdFlux[m][i] = NekDouble(Fwdreslt[i]);
449 }
450 }
451
452 for (size_t m = 0; m < nvariables; ++m)
453 {
454 for (size_t i = 0; i < nTracePts; ++i)
455 {
456 Fwdarray[m][i] = DataType(Bwd[m][i]);
457 }
458 }
459 for (size_t m = 0; m < nvariables; ++m)
460 {
461 Vmath::Zero(nTracePts, Fwdreslt, 1);
462 for (size_t n = 0; n < nvariables; ++n)
463 {
464 for (size_t p = 0; p < nTracePts; ++p)
465 {
466 Fwdreslt[p] += TraceJacArray[1][m][n][p] * Fwdarray[n][p];
467 }
468 }
469 for (size_t i = 0; i < nTracePts; ++i)
470 {
471 BwdFlux[m][i] = NekDouble(Fwdreslt[i]);
472 }
473 }
474
475 for (size_t i = 0; i < nvariables; ++i)
476 {
477 Vmath::Fill(nCoeffs, 0.0, outarray[i], 1);
478 timer.Start();
479 pFields[i]->AddTraceIntegralToOffDiag(FwdFlux[i], BwdFlux[i],
480 outarray[i]);
481 timer.Stop();
482 timer.AccumulateRegion("ExpList::AddTraceIntegralToOffDiag", 10);
483 }
484
485 for (size_t i = 0; i < nvariables; ++i)
486 {
487 for (size_t p = 0; p < nCoeffs; ++p)
488 {
489 outarray[i][p] =
490 -m_DtLambdaPreconMat * outarray[i][p] + inarray[i][p];
491 }
492 }
493}
494
495/**
496 *
497 */
498template <typename TypeNekBlkMatSharedPtr>
502 const int &nscale)
503{
504
505 size_t nvars = pFields.size();
506 size_t nelmts = pFields[0]->GetNumElmts();
507 size_t nelmtcoef;
508 Array<OneD, unsigned int> nelmtmatdim(nelmts);
509 for (size_t i = 0; i < nelmts; i++)
510 {
511 nelmtcoef = pFields[0]->GetExp(i)->GetNcoeffs();
512 nelmtmatdim[i] = nelmtcoef * nscale;
513 }
514
515 for (size_t i = 0; i < nvars; i++)
516 {
517 for (size_t j = 0; j < nvars; j++)
518 {
519 AllocateNekBlkMatDig(gmtxarray[i][j], nelmtmatdim, nelmtmatdim);
520 }
521 }
522}
523
524} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
void DoCalcPreconMatBRJCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr > > &gmtxarray, SNekBlkMatSharedPtr &gmtVar, Array< OneD, SNekBlkMatSharedPtr > &TraceJac, Array< OneD, SNekBlkMatSharedPtr > &TraceJacDeriv, Array< OneD, Array< OneD, NekSingle > > &TraceJacDerivSign, TensorOfArray4D< NekSingle > &TraceJacArray, TensorOfArray4D< NekSingle > &TraceJacDerivArray, TensorOfArray5D< NekSingle > &TraceIPSymJacArray)
Definition: PreconCfsOp.h:97
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
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)
NekPreconCfsOperators m_operator
Definition: PreconCfs.h:91
int m_PreconTimesCounter
Definition: PreconCfs.h:93
int m_PreconMatFreezNumb
Definition: PreconCfs.h:92
NekDouble m_DtLambdaPreconMat
Definition: PreconCfs.h:94
bool m_CalcPreconMatFlag
Definition: PreconCfs.h:95
LibUtilities::CommSharedPtr m_Comm
Definition: PreconCfs.h:88
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
PreconCfsFactory & GetPreconCfsFactory()
Declaration of the boundary condition factory singleton.
Definition: PreconCfs.cpp:42
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
tinysimd::simd< NekDouble > vec_t
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
typename abi< ScalarType, width >::type simd
Definition: tinysimd.hpp:80