Nektar++
DiffusionIP.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: DiffusionIP.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: IP diffusion class.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37#include <cmath>
38#include <iomanip>
39#include <iostream>
40
42
43namespace Nektar::SolverUtils
44{
46 "InteriorPenalty", DiffusionIP::create, "Interior Penalty");
47
49{
50}
51
55{
56 m_session = pSession;
57
58 m_session->LoadParameter("IPSymmFluxCoeff", m_IPSymmFluxCoeff,
59 0.0); // SIP=1.0; NIP=-1.0; IIP=0.0
60
61 m_session->LoadParameter("IP2ndDervCoeff", m_IP2ndDervCoeff,
62 0.0); // 1.0/12.0
63
64 m_session->LoadParameter("IPPenaltyCoeff", m_IPPenaltyCoeff,
65 4.0); // 1.0/12.0
66
67 // Setting up the normals
68 size_t nDim = pFields[0]->GetCoordim(0);
69 size_t nVariable = pFields.size();
70 size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
71
73 for (int i = 0; i < nDim; ++i)
74 {
75 m_traceNormals[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
76 }
79 for (int i = 0; i < nVariable; ++i)
80 {
81 m_traceAver[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
82 m_traceJump[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
83 }
84
85 pFields[0]->GetTrace()->GetNormals(m_traceNormals);
86
87 // Compute the length of the elements in the boundary-normal direction
88 // The function "GetElmtNormalLength" returns half the length of the left
89 // and right adjacent element in the lengthFwd and lengthBwd arrays. If
90 // the element belongs to a boundary (including periodic boundaries) or
91 // a parallel interface, the Bwd array will contain zeros.
92 Array<OneD, NekDouble> lengthFwd{nTracePts, 0.0};
93 Array<OneD, NekDouble> lengthBwd{nTracePts, 0.0};
94 pFields[0]->GetTrace()->GetElmtNormalLength(lengthFwd, lengthBwd);
95
96 // Copy Fwd to Bwd on parallel interfaces
97 // TODO: Move this into GetElmtNormalLength()
98 pFields[0]->GetTraceMap()->GetAssemblyCommDG()->PerformExchange(lengthFwd,
99 lengthBwd);
100 // Copy Fwd to Bwd on periodic interfaces
101 // TODO: Move this into GetElmtNormalLength()
102 pFields[0]->PeriodicBwdCopy(lengthFwd, lengthBwd);
103 // Scale the length by 0.5 on boundaries, and copy Fwd into Bwd
104 // Notes:
105 // - It is not quite clear why we need to scale by 0.5 on the boundaries
106 // - The current implementation is not perfect, it would be nicer to call
107 // a function similar to DiscontField::v_FillBwdWithBoundCond() with
108 // PutFwdInBwdOnBCs = true. If we wouldn't do the scaling by 0.5, this
109 // function could have been used.
110 for (int i = 0; i < nTracePts; ++i)
111 {
112 if (std::abs(lengthBwd[i]) < NekConstants::kNekMachineEpsilon)
113 {
114 lengthFwd[i] *= 0.5;
115 lengthBwd[i] = lengthFwd[i];
116 }
117 }
118
119 // Compute the average element normal length along the edge
120 Array<OneD, NekDouble> lengthsum(nTracePts, 0.0);
121 Array<OneD, NekDouble> lengthmul(nTracePts, 0.0);
122 Vmath::Vadd(nTracePts, lengthBwd, 1, lengthFwd, 1, lengthsum, 1);
123 Vmath::Vmul(nTracePts, lengthBwd, 1, lengthFwd, 1, lengthmul, 1);
124 Vmath::Vdiv(nTracePts, lengthsum, 1, lengthmul, 1, lengthFwd, 1);
125 m_traceNormDirctnElmtLength = lengthsum;
129
132 pFields[0]->GetBwdWeight(m_tracBwdWeightAver, m_tracBwdWeightJump);
133 Array<OneD, NekDouble> tmpBwdWeight{nTracePts, 0.0};
134 Array<OneD, NekDouble> tmpBwdWeightJump{nTracePts, 0.0};
135 for (int i = 1; i < nVariable; ++i)
136 {
137 pFields[i]->GetBwdWeight(tmpBwdWeight, tmpBwdWeightJump);
138 Vmath::Vsub(nTracePts, tmpBwdWeight, 1, m_tracBwdWeightAver, 1,
139 tmpBwdWeight, 1);
140 Vmath::Vabs(nTracePts, tmpBwdWeight, 1, tmpBwdWeight, 1);
141 NekDouble norm = 0.0;
142 for (int j = 0; j < nTracePts; ++j)
143 {
144 norm += tmpBwdWeight[j];
145 }
146 ASSERTL0(norm < 1.0E-11,
147 "different BWD for different variable not coded yet");
148 }
149
150 // workspace for v_diffuse
151 size_t nCoeffs = pFields[0]->GetNcoeffs();
153 for (int i = 0; i < nVariable; ++i)
154 {
155 m_wspDiff[i] = Array<OneD, NekDouble>{nCoeffs, 0.0};
156 }
157
158 // workspace for callnumtraceflux
161 for (int nd = 0; nd < nDim; ++nd)
162 {
165 for (int i = 0; i < nVariable; ++i)
166 {
167 m_wspNumDerivBwd[nd][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
168 m_wspNumDerivFwd[nd][i] = Array<OneD, NekDouble>{nTracePts, 0.0};
169 }
170 }
171}
172
174 const size_t nConvectiveFields,
176 const Array<OneD, Array<OneD, NekDouble>> &inarray,
178 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
179 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
180{
181
183 timer.Start();
184 DiffusionIP::v_DiffuseCoeffs(nConvectiveFields, fields, inarray, m_wspDiff,
185 pFwd, pBwd);
186 for (int i = 0; i < nConvectiveFields; ++i)
187 {
188 fields[i]->BwdTrans(m_wspDiff[i], outarray[i]);
189 }
190 timer.Stop();
191 timer.AccumulateRegion("Diffusion IP");
192}
193
195 const size_t nConvectiveFields,
197 const Array<OneD, Array<OneD, NekDouble>> &inarray,
199 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
200 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
201{
202 if (fields[0]->GetGraph()->GetMovement()->GetMoveFlag()) // i.e. if
203 // m_ALESolver
204 {
205 fields[0]->GetTrace()->GetNormals(m_traceNormals);
206 }
207
209 timer.Start();
210
211 size_t nDim = fields[0]->GetCoordim(0);
212 size_t nPts = fields[0]->GetTotPoints();
213 size_t nCoeffs = fields[0]->GetNcoeffs();
214 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
215
216 // pre-allocate this?
217 Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
218 Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
219 TensorOfArray3D<NekDouble> elmtFlux{nDim};
220 TensorOfArray3D<NekDouble> qfield{nDim};
221
222 for (int j = 0; j < nDim; ++j)
223 {
224 qfield[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
225 elmtFlux[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
226 for (int i = 0; i < nConvectiveFields; ++i)
227 {
228 qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
229 }
230 for (int i = 0; i < nConvectiveFields; ++i)
231 {
232 elmtFlux[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
233 }
234 }
235
236 // pre-allocate this?
237 Array<OneD, Array<OneD, NekDouble>> vFwd{nConvectiveFields};
238 Array<OneD, Array<OneD, NekDouble>> vBwd{nConvectiveFields};
239 // when does this happen?
241 {
242 for (int i = 0; i < nConvectiveFields; ++i)
243 {
244 vFwd[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
245 vBwd[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
246 }
247 for (int i = 0; i < nConvectiveFields; ++i)
248 {
249 fields[i]->GetFwdBwdTracePhys(inarray[i], vFwd[i], vBwd[i]);
250 }
251 }
252 else
253 {
254 for (int i = 0; i < nConvectiveFields; ++i)
255 {
256 vFwd[i] = pFwd[i];
257 vBwd[i] = pBwd[i];
258 }
259 }
260
261 DiffuseCalcDerivative(fields, inarray, qfield, vFwd, vBwd);
262 Array<OneD, int> nonZeroIndex;
263 DiffuseVolumeFlux(fields, inarray, qfield, elmtFlux, nonZeroIndex);
264 timer.Stop();
265 timer.AccumulateRegion("Diffusion:Volumeflux", 10);
266 timer.Start();
267
268 // pre-allocate this?
269 Array<OneD, Array<OneD, NekDouble>> tmpFluxIprdct{nDim};
270 // volume intergration: the nonZeroIndex indicates which flux is nonzero
271 for (int i = 0; i < nonZeroIndex.size(); ++i)
272 {
273 int j = nonZeroIndex[i];
274 for (int k = 0; k < nDim; ++k)
275 {
276 tmpFluxIprdct[k] = elmtFlux[k][j];
277 }
278 fields[j]->IProductWRTDerivBase(tmpFluxIprdct, outarray[j]);
279 Vmath::Neg(nCoeffs, outarray[j], 1);
280 }
281 timer.Stop();
282 timer.AccumulateRegion("Diffusion:IPWRTDB", 10);
283
284 // release qfield, elmtFlux and muvar;
285 timer.Start();
286 for (int j = 0; j < nDim; ++j)
287 {
288 elmtFlux[j] = NullNekDoubleArrayOfArray;
289 }
290
291 // pre-allocate this?
292 Array<OneD, Array<OneD, NekDouble>> Traceflux{nConvectiveFields};
293 for (int j = 0; j < nConvectiveFields; ++j)
294 {
295 Traceflux[j] = Array<OneD, NekDouble>{nTracePts, 0.0};
296 }
297
298 DiffuseTraceFlux(fields, inarray, qfield, elmtFlux, Traceflux, vFwd, vBwd,
299 nonZeroIndex);
300 timer.Stop();
301 timer.AccumulateRegion("Diffusion:TraceFlux", 10);
302
303 for (int i = 0; i < nonZeroIndex.size(); ++i)
304 {
305 int j = nonZeroIndex[i];
306
307 fields[j]->AddTraceIntegral(Traceflux[j], outarray[j]);
308 fields[j]->SetPhysState(false);
309 }
310
311 AddDiffusionSymmFluxToCoeff(nConvectiveFields, fields, inarray, qfield,
312 elmtFlux, outarray, vFwd, vBwd);
313
314 timer.Start();
315
316 if (!fields[0]->GetGraph()->GetMovement()->GetMoveFlag() ||
317 fields[0]
318 ->GetGraph()
319 ->GetMovement()
320 ->GetImplicitALESolverFlag()) // i.e. if
321 // m_ALESolver
322 {
323 for (int i = 0; i < nonZeroIndex.size(); ++i)
324 {
325 int j = nonZeroIndex[i];
326
327 fields[j]->MultiplyByElmtInvMass(outarray[j], outarray[j]);
328 }
329 }
330
331 timer.Stop();
332 timer.AccumulateRegion("DiffIP:Diffusion Coeff", 10);
333}
334
337 const Array<OneD, Array<OneD, NekDouble>> &inarray,
339 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &pFwd,
340 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &pBwd)
341{
343 size_t nDim = fields[0]->GetCoordim(0);
344 for (int nd = 0; nd < 3; ++nd)
345 {
346 qtmp[nd] = NullNekDouble1DArray;
347 }
348
349 size_t nConvectiveFields = inarray.size();
350 for (int i = 0; i < nConvectiveFields; ++i)
351 {
352 for (int nd = 0; nd < nDim; ++nd)
353 {
354 qtmp[nd] = qfield[nd][i];
355 }
356 fields[i]->PhysDeriv(inarray[i], qtmp[0], qtmp[1], qtmp[2]);
357 }
358}
359
362 const Array<OneD, Array<OneD, NekDouble>> &inarray,
364 Array<OneD, int> &nonZeroIndex)
365{
366 size_t nDim = fields[0]->GetCoordim(0);
367
369
371 timer.Start();
372 m_FunctorDiffusionfluxCons(nDim, inarray, qfield, VolumeFlux, nonZeroIndex,
373 tmparray2D);
374 timer.Stop();
375 timer.AccumulateRegion("DiffIP:_FunctorDiffFluxCons", 10);
376}
377
380 const Array<OneD, Array<OneD, NekDouble>> &inarray,
382 [[maybe_unused]] TensorOfArray3D<NekDouble> &VolumeFlux,
383 Array<OneD, Array<OneD, NekDouble>> &TraceFlux,
384 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
385 const Array<OneD, Array<OneD, NekDouble>> &pBwd,
386 Array<OneD, int> &nonZeroIndex)
387{
388 TensorOfArray3D<NekDouble> traceflux3D(1);
389 traceflux3D[0] = TraceFlux;
390
391 size_t nConvectiveFields = fields.size();
393 timer.Start();
394 CalcTraceNumFlux(m_IP2ndDervCoeff, fields, inarray, qfield, pFwd, pBwd,
395 nonZeroIndex, traceflux3D, m_traceAver, m_traceJump);
396 timer.Stop();
397 timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux", 10);
398
399 ApplyFluxBndConds(nConvectiveFields, fields, TraceFlux);
400}
401
403 const std::size_t nConvectiveFields,
405 const Array<OneD, Array<OneD, NekDouble>> &inarray,
408 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
409 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
410{
411 if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
412 {
413 size_t nDim = fields[0]->GetCoordim(0);
414 size_t nPts = fields[0]->GetTotPoints();
415 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
416 TensorOfArray3D<NekDouble> traceSymflux{nDim};
417 for (int nd = 0; nd < nDim; ++nd)
418 {
419 traceSymflux[nd] =
420 Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
421 for (int j = 0; j < nConvectiveFields; ++j)
422 {
423 traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
424 }
425 }
426 Array<OneD, int> nonZeroIndex;
427 DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
428 VolumeFlux, traceSymflux, pFwd, pBwd,
429 nonZeroIndex);
430
431 AddSymmFluxIntegralToCoeff(nConvectiveFields, nDim, nPts, nTracePts,
432 fields, nonZeroIndex, traceSymflux,
433 outarray);
434 }
435}
436
438 const std::size_t nConvectiveFields,
440 const Array<OneD, Array<OneD, NekDouble>> &inarray,
443 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
444 const Array<OneD, Array<OneD, NekDouble>> &pBwd)
445{
446 if (fabs(m_IPSymmFluxCoeff) > 1.0E-12)
447 {
448 size_t nDim = fields[0]->GetCoordim(0);
449 size_t nPts = fields[0]->GetTotPoints();
450 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
451 TensorOfArray3D<NekDouble> traceSymflux{nDim};
452 for (int nd = 0; nd < nDim; ++nd)
453 {
454 traceSymflux[nd] =
455 Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
456 for (int j = 0; j < nConvectiveFields; ++j)
457 {
458 traceSymflux[nd][j] = Array<OneD, NekDouble>{nTracePts, 0.0};
459 }
460 }
461 Array<OneD, int> nonZeroIndex;
462 DiffuseTraceSymmFlux(nConvectiveFields, fields, inarray, qfield,
463 VolumeFlux, traceSymflux, pFwd, pBwd,
464 nonZeroIndex);
465
466 AddSymmFluxIntegralToPhys(nConvectiveFields, nDim, nPts, nTracePts,
467 fields, nonZeroIndex, traceSymflux, outarray);
468 }
469}
470
472 const std::size_t nConvectiveFields,
474 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &inarray,
475 [[maybe_unused]] TensorOfArray3D<NekDouble> &qfield,
476 [[maybe_unused]] TensorOfArray3D<NekDouble> &VolumeFlux,
478 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &pFwd,
479 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &pBwd,
480 Array<OneD, int> &nonZeroIndex)
481{
482 size_t nDim = fields[0]->GetCoordim(0);
483
484 CalcTraceSymFlux(nConvectiveFields, nDim, m_traceAver, m_traceJump,
485 nonZeroIndex, SymmFlux);
486}
487
489 const std::size_t nConvectiveFields, const size_t nDim,
490 const Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
491 Array<OneD, Array<OneD, NekDouble>> &solution_jump,
492 Array<OneD, int> &nonZeroIndexsymm,
493 TensorOfArray3D<NekDouble> &traceSymflux)
494{
495 size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
496
497 m_FunctorSymmetricfluxCons(nDim, solution_Aver, solution_jump, traceSymflux,
498 nonZeroIndexsymm, m_traceNormals);
499
500 for (int nd = 0; nd < nDim; ++nd)
501 {
502 for (int j = 0; j < nonZeroIndexsymm.size(); ++j)
503 {
504 int i = nonZeroIndexsymm[j];
505 Vmath::Smul(nTracePts, -0.5 * m_IPSymmFluxCoeff,
506 traceSymflux[nd][i], 1, traceSymflux[nd][i], 1);
507 }
508 }
509}
510
512 const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts,
513 [[maybe_unused]] const size_t nTracePts,
515 const Array<OneD, const int> &nonZeroIndex,
518{
519 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
520 Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
522 for (int i = 0; i < nDim; ++i)
523 {
524 tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
525 }
526 int nv = 0;
527 for (int j = 0; j < nonZeroIndex.size(); ++j)
528 {
529 nv = nonZeroIndex[j];
530 MultiRegions::ExpListSharedPtr tracelist = fields[nv]->GetTrace();
531 for (int nd = 0; nd < nDim; ++nd)
532 {
533 Vmath::Zero(nPts, tmpfield[nd], 1);
534
535 tracelist->MultiplyByQuadratureMetric(tracflux[nd][nv],
536 tracflux[nd][nv]);
537
538 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
539 tracflux[nd][nv], tmpfield[nd]);
540 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
541 }
542 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
543 Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
544 }
545}
546
548 const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts,
549 [[maybe_unused]] const size_t nTracePts,
551 const Array<OneD, const int> &nonZeroIndex,
554{
555 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
556 Array<OneD, NekDouble> tmpCoeff{nCoeffs, 0.0};
557 Array<OneD, NekDouble> tmpPhysi{nPts, 0.0};
559 for (int i = 0; i < nDim; ++i)
560 {
561 tmpfield[i] = Array<OneD, NekDouble>{nPts, 0.0};
562 }
563 for (int j = 0; j < nonZeroIndex.size(); ++j)
564 {
565 int nv = nonZeroIndex[j];
566 for (int nd = 0; nd < nDim; ++nd)
567 {
568 Vmath::Zero(nPts, tmpfield[nd], 1);
569
570 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
571 tracflux[nd][nv], tmpfield[nd]);
572 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
573 }
574 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
575 fields[nv]->BwdTrans(tmpCoeff, tmpPhysi);
576 Vmath::Vadd(nPts, tmpPhysi, 1, outarray[nv], 1, outarray[nv], 1);
577 }
578}
579
583{
584 MultiRegions::ExpListSharedPtr tracelist = fields[0]->GetTrace();
585 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
586 tracelist->GetExp();
587 size_t ntotTrac = (*traceExp).size();
588 int nTracPnt, noffset;
589
590 const MultiRegions::LocTraceToTraceMapSharedPtr locTraceToTraceMap =
591 fields[0]->GetLocTraceToTraceMap();
592
593 const Array<OneD, const Array<OneD, int>> LRAdjExpid =
594 locTraceToTraceMap->GetLeftRightAdjacentExpId();
595 const Array<OneD, const Array<OneD, bool>> LRAdjflag =
596 locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
597
598 std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
599 fields[0]->GetExp();
600
601 Array<OneD, NekDouble> factorFwdBwd{2, 0.0};
602
603 NekDouble spaceDim = NekDouble(fields[0]->GetCoordim(0));
604
605 int ntmp, numModes;
606
607 for (int ntrace = 0; ntrace < ntotTrac; ++ntrace)
608 {
609 noffset = tracelist->GetPhys_Offset(ntrace);
610 nTracPnt = tracelist->GetTotPoints(ntrace);
611
612 factorFwdBwd[0] = 0.0;
613 factorFwdBwd[1] = 0.0;
614
615 for (int nlr = 0; nlr < 2; ++nlr)
616 {
617 if (LRAdjflag[nlr][ntrace])
618 {
619 numModes = 0;
620 for (int nd = 0; nd < spaceDim; nd++)
621 {
622 ntmp = fields[0]
623 ->GetExp(LRAdjExpid[nlr][ntrace])
624 ->GetBasisNumModes(nd);
625 numModes = std::max(ntmp, numModes);
626 }
627 factorFwdBwd[nlr] = (numModes) * (numModes);
628 }
629 }
630
631 for (int np = 0; np < nTracPnt; ++np)
632 {
633 factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
634 }
635 }
636}
637
639 const std::size_t nConvectiveFields, const size_t nPts,
640 const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
641 const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
644{
645 std::vector<NekDouble> vFwdTmp(nConvectiveFields),
646 vBwdTmp(nConvectiveFields), averTmp(nConvectiveFields);
647 for (size_t p = 0; p < nPts; ++p)
648 {
649 // re-arrange data
650 for (size_t f = 0; f < nConvectiveFields; ++f)
651 {
652 vFwdTmp[f] = vFwd[f][p];
653 vBwdTmp[f] = vBwd[f][p];
654 }
655
656 ConsVarAve(nConvectiveFields, m_tracBwdWeightAver[p], vFwdTmp, vBwdTmp,
657 averTmp);
658
659 // store
660 for (size_t f = 0; f < nConvectiveFields; ++f)
661 {
662 aver[f][p] = averTmp[f];
663 }
664 }
665
666 // if this can be make function of points, the nPts loop can be lifted more
667 m_SpecialBndTreat(aver);
668
669 // note: here the jump is 2.0*(aver-vFwd)
670 // because Viscous wall use a symmetry value as the Bwd,
671 // not the target one
672
673 for (size_t f = 0; f < nConvectiveFields; ++f)
674 {
675 for (size_t p = 0; p < nPts; ++p)
676 {
678 NekDouble Fweight = 2.0 - Bweight;
679
680 NekDouble tmpF = aver[f][p] - vFwd[f][p];
681 NekDouble tmpB = vBwd[f][p] - aver[f][p];
682 jump[f][p] = tmpF * Fweight + tmpB * Bweight;
683 }
684 }
685}
686
688 const NekDouble PenaltyFactor2,
690 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &inarray,
691 const TensorOfArray3D<NekDouble> &qfield,
692 const Array<OneD, Array<OneD, NekDouble>> &vFwd,
693 const Array<OneD, Array<OneD, NekDouble>> &vBwd,
694 Array<OneD, int> &nonZeroIndexflux, TensorOfArray3D<NekDouble> &traceflux,
695 Array<OneD, Array<OneD, NekDouble>> &solution_Aver,
696 Array<OneD, Array<OneD, NekDouble>> &solution_jump)
697{
698 size_t nDim = fields[0]->GetCoordim(0);
699 size_t nPts = fields[0]->GetTotPoints();
700 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
701 size_t nConvectiveFields = fields.size();
702
704 fields[0]->GetTraceMap();
705 const MultiRegions::InterfaceMapDGSharedPtr InterfaceMap =
706 fields[0]->GetInterfaceMap();
707
709 timer.Start();
710 // with further restructuring this iniziatilization could be eliminated
711 for (int nd = 0; nd < nDim; ++nd)
712 {
713 for (int i = 0; i < nConvectiveFields; ++i)
714 {
715 Vmath::Zero(nTracePts, m_wspNumDerivBwd[nd][i], 1);
716 Vmath::Zero(nTracePts, m_wspNumDerivFwd[nd][i], 1);
717 }
718 }
719
720 // could this be pre-allocated?
721 Array<OneD, NekDouble> Fwd{nTracePts, 0.0};
722 Array<OneD, NekDouble> Bwd{nTracePts, 0.0};
723
724 timer.Stop();
725 timer.AccumulateRegion("DiffIP:_CalcTraceNumFlux_alloc", 10);
726
727 timer.Start();
728 if (fabs(PenaltyFactor2) > 1.0E-12)
729 {
730 AddSecondDerivToTrace(nConvectiveFields, nDim, nPts, nTracePts,
731 PenaltyFactor2, fields, qfield, m_wspNumDerivFwd,
733 }
734
735 timer.Stop();
736 timer.AccumulateRegion("DiffIP:_AddSecondDerivToTrace", 10);
737
738 for (int nd = 0; nd < nDim; ++nd)
739 {
740 for (int i = 0; i < nConvectiveFields; ++i)
741 {
742 // this sequence of operations is really exepensive,
743 // it should be done collectively for all fields together instead of
744 // one by one
745 timer.Start();
746 fields[i]->GetFwdBwdTracePhys(qfield[nd][i], Fwd, Bwd, true, true,
747 false);
748 timer.Stop();
749 timer.AccumulateRegion("DiffIP:_GetFwdBwdTracePhys", 10);
750
751 for (size_t p = 0; p < nTracePts; ++p)
752 {
753 m_wspNumDerivBwd[nd][i][p] += 0.5 * Bwd[p];
754 m_wspNumDerivFwd[nd][i][p] += 0.5 * Fwd[p];
755 }
756
757 timer.Start();
758 TraceMap->GetAssemblyCommDG()->PerformExchange(
759 m_wspNumDerivFwd[nd][i], m_wspNumDerivBwd[nd][i]);
760 InterfaceMap->ExchangeTrace(m_wspNumDerivFwd[nd][i],
761 m_wspNumDerivBwd[nd][i]);
762 timer.Stop();
763 timer.AccumulateRegion("DiffIP:_PerformExchange", 10);
764
765 Vmath::Vadd(nTracePts, m_wspNumDerivFwd[nd][i], 1,
766 m_wspNumDerivBwd[nd][i], 1, m_wspNumDerivFwd[nd][i], 1);
767 }
768 }
769
770 timer.Start();
771 ConsVarAveJump(nConvectiveFields, nTracePts, vFwd, vBwd, solution_Aver,
772 solution_jump);
773 timer.Stop();
774 timer.AccumulateRegion("DiffIP:_ConsVarAveJump", 10);
775
776 Array<OneD, NekDouble> penaltyCoeff(nTracePts, 0.0);
777 GetPenaltyFactor(fields, penaltyCoeff);
778 for (size_t p = 0; p < nTracePts; ++p)
779 {
780 NekDouble PenaltyFactor =
781 penaltyCoeff[p] * m_traceNormDirctnElmtLengthRecip[p]; // load 1x
782
783 for (size_t f = 0; f < nConvectiveFields; ++f)
784 {
785 NekDouble jumpTmp = solution_jump[f][p] * PenaltyFactor; // load 1x
786 for (size_t d = 0; d < nDim; ++d)
787 {
788 NekDouble tmp = m_traceNormals[d][p] * jumpTmp +
789 m_wspNumDerivFwd[d][f][p]; // load 2x
790 m_wspNumDerivFwd[d][f][p] = tmp; // store 1x
791 }
792 }
793 }
794
795 timer.Start();
796 // Calculate normal viscous flux
798 traceflux, nonZeroIndexflux,
800 timer.Stop();
801 timer.AccumulateRegion("DiffIP:_FunctorDiffFluxConsTrace", 10);
802}
803
805 const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts,
806 const size_t nTracePts, const NekDouble PenaltyFactor2,
808 const TensorOfArray3D<NekDouble> &qfield,
809 TensorOfArray3D<NekDouble> &numDerivFwd,
810 TensorOfArray3D<NekDouble> &numDerivBwd)
811{
812 Array<OneD, NekDouble> Fwd(nTracePts, 0.0);
813 Array<OneD, NekDouble> Bwd(nTracePts, 0.0);
814 std::vector<NekDouble> tmp(nTracePts);
815
816 Array<OneD, Array<OneD, NekDouble>> elmt2ndDerv{nDim};
817 for (int nd1 = 0; nd1 < nDim; ++nd1)
818 {
819 elmt2ndDerv[nd1] = Array<OneD, NekDouble>{nPts, 0.0};
820 }
821
823 for (int nd = 0; nd < 3; ++nd)
824 {
825 qtmp[nd] = NullNekDouble1DArray;
826 }
827 for (int nd2 = 0; nd2 < nDim; ++nd2)
828 {
829 qtmp[nd2] = elmt2ndDerv[nd2];
830 }
831
832 for (int i = 0; i < nTracePts; ++i)
833 {
834 tmp[i] = PenaltyFactor2 * m_traceNormDirctnElmtLength[i];
835 }
836 // the derivatives are assumed to be exchangable
837 for (int nd1 = 0; nd1 < nDim; ++nd1)
838 {
839 for (int i = 0; i < nConvectiveFields; ++i)
840 {
841 fields[i]->PhysDeriv(qfield[nd1][i], qtmp[0], qtmp[1], qtmp[2]);
842
843 for (int nd2 = nd1; nd2 < nDim; ++nd2)
844 {
845 Vmath::Zero(nTracePts, Bwd, 1);
846 fields[i]->GetFwdBwdTracePhys(elmt2ndDerv[nd2], Fwd, Bwd, true,
847 true, false);
848 for (int p = 0; p < nTracePts; ++p)
849 {
850 Bwd[p] *= tmp[p];
851 numDerivBwd[nd1][i][p] += m_traceNormals[nd2][p] * Bwd[p];
852 Fwd[p] *= tmp[p];
853 numDerivFwd[nd1][i][p] = -(m_traceNormals[nd2][p] * Fwd[p] -
854 numDerivFwd[nd1][i][p]);
855 }
856 if (nd2 != nd1)
857 {
858 for (int p = 0; p < nTracePts; ++p)
859 {
860 numDerivBwd[nd2][i][p] +=
861 m_traceNormals[nd1][p] * Bwd[p];
862 numDerivFwd[nd2][i][p] =
863 -(m_traceNormals[nd1][p] * Fwd[p] -
864 numDerivFwd[nd2][i][p]);
865 }
866 }
867 }
868 }
869 }
870}
871
872/**
873 * @brief aplly Neuman boundary conditions on flux
874 * Currently only consider WallAdiabatic
875 *
876 **/
878 const int nConvectiveFields,
881{
882 int nengy = nConvectiveFields - 1;
883 int cnt;
884 // Compute boundary conditions for Energy
885 cnt = 0;
886 size_t nBndRegions = fields[nengy]->GetBndCondExpansions().size();
887 for (size_t j = 0; j < nBndRegions; ++j)
888 {
889 if (fields[nengy]->GetBndConditions()[j]->GetBoundaryConditionType() ==
891 {
892 continue;
893 }
894
895 size_t nBndEdges =
896 fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
897 for (size_t e = 0; e < nBndEdges; ++e)
898 {
899 size_t nBndEdgePts = fields[nengy]
900 ->GetBndCondExpansions()[j]
901 ->GetExp(e)
902 ->GetTotPoints();
903
904 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
905 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
906
907 if (fields[0]->GetBndConditions()[j]->GetUserDefined() ==
908 "WallAdiabatic")
909 {
910 Vmath::Zero(nBndEdgePts, &flux[nengy][id2], 1);
911 }
912 }
913 }
914}
915
916} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
Definition: Diffusion.h:202
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
Definition: Diffusion.h:357
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion Volume FLux.
Definition: Diffusion.h:215
SOLVER_UTILS_EXPORT void DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion term Trace Flux.
Definition: Diffusion.h:226
DiffusionFluxCons m_FunctorDiffusionfluxCons
Definition: Diffusion.h:356
SpecialBndTreat m_SpecialBndTreat
Definition: Diffusion.h:358
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Definition: Diffusion.h:359
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
Definition: DiffusionIP.h:121
void v_DiffuseCoeffs(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
static DiffusionSharedPtr create(std::string diffType)
Definition: DiffusionIP.h:46
Array< OneD, NekDouble > m_tracBwdWeightAver
Definition: DiffusionIP.h:118
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
Definition: DiffusionIP.h:120
void ApplyFluxBndConds(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, Array< OneD, NekDouble > > &flux)
aplly Neuman boundary conditions on flux Currently only consider WallAdiabatic
void GetPenaltyFactor(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get IP penalty factor based on order.
void ConsVarAve(const size_t nConvectiveFields, const T &Bweight, const std::vector< T > &vFwd, const std::vector< T > &vBwd, std::vector< T > &aver)
Calculate the average of conservative variables on traces.
Definition: DiffusionIP.h:128
void AddSecondDerivToTrace(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &numDerivFwd, TensorOfArray3D< NekDouble > &numDerivBwd)
Add second derivative term to trace jump (for DDG scheme)
void CalcTraceNumFlux(const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble > > &vFwd, const Array< OneD, Array< OneD, NekDouble > > &vBwd, Array< OneD, int > &nonZeroIndexflux, TensorOfArray3D< NekDouble > &traceflux, Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump)
Calculate numerical flux on traces.
void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
Diffusion Flux, calculate the physical derivatives.
void AddDiffusionSymmFluxToPhys(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Array< OneD, Array< OneD, NekDouble > > m_traceJump
Definition: DiffusionIP.h:111
void AddSymmFluxIntegralToCoeff(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
Add symmetric flux integration to field (in coefficient space)
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
Definition: DiffusionIP.h:115
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
Definition: DiffusionIP.h:116
LibUtilities::SessionReaderSharedPtr m_session
Definition: DiffusionIP.h:122
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Definition: DiffusionIP.h:110
Array< OneD, NekDouble > m_tracBwdWeightJump
Definition: DiffusionIP.h:119
void v_Diffuse(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
void ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble > > &vFwd, const Array< OneD, const Array< OneD, NekDouble > > &vBwd, Array< OneD, Array< OneD, NekDouble > > &aver, Array< OneD, Array< OneD, NekDouble > > &jump)
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Definition: DiffusionIP.cpp:52
void v_DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex) override
Diffusion Volume Flux.
void AddDiffusionSymmFluxToCoeff(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
void DiffuseTraceSymmFlux(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex)
Calculate symmetric flux on traces interface.
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
Workspace for v_Diffusion.
Definition: DiffusionIP.h:113
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Definition: DiffusionIP.h:109
void v_DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex) override
Diffusion term Trace Flux.
void CalcTraceSymFlux(const std::size_t nConvectiveFields, const size_t nDim, const Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump, Array< OneD, int > &nonZeroIndexsymm, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &traceSymflux)
Calculate symmetric flux on traces.
void AddSymmFluxIntegralToPhys(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
Add symmetric flux integration to field (in physical space)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
Definition: AssemblyMapDG.h:46
std::shared_ptr< InterfaceMapDG > InterfaceMapDGSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr
static const NekDouble kNekMachineEpsilon
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:39
std::vector< double > d(NPUPPER *NPUPPER)
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
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 Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.hpp:352
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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 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 Vdiv(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:126
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
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298