Nektar++
Loading...
Searching...
No Matches
CompressibleFlowSystemImplicit.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: CompressibleFlowSystemImplicit.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: Compressible flow system base class with auxiliary functions
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
39
40namespace Nektar
41{
42
45 : UnsteadySystem(pSession, pGraph), CompressibleFlowSystem(pSession, pGraph)
46{
47}
48
49/**
50 * @brief Initialization object for CFSImplicit class.
51 */
52void CFSImplicit::v_InitObject(bool DeclareFields)
53{
55 m_explicitAdvection = false;
56 m_explicitDiffusion = false;
57
58 // Initialise implicit parameters
59 m_session->MatchSolverInfo("FLAGIMPLICITITSSTATISTICS", "True",
61
62 m_session->LoadParameter("JacobiFreeEps", m_jacobiFreeEps,
64
65 m_session->LoadParameter("nPadding", m_nPadding, 4);
66
67 int ntmp;
68 m_session->LoadParameter("AdvectionJacFlag", ntmp, 1);
69 m_advectionJacFlag = (ntmp != 0);
70
71 m_session->LoadParameter("ViscousJacFlag", ntmp, 1);
72 m_viscousJacFlag = (ntmp != 0);
73
74 // Initialise implicit functors
77
79}
80
82{
83 int nvariables = m_fields.size();
84 int ntotal = m_fields[0]->GetNcoeffs() * nvariables;
85
86 std::string SolverType = "Newton";
87 if (m_session->DefinesSolverInfo("NonlinSysIterSolver"))
88 {
89 SolverType = m_session->GetSolverInfo("NonlinSysIterSolver");
90 }
93 "NekNonlinSys '" + SolverType + "' is not defined.\n");
94
95 // Create the key to hold settings for nonlin solver
97 // Load required LinSys parameters:
98 m_session->LoadParameter("NekLinSysMaxIterations",
100 m_session->LoadParameter("LinSysMaxStorage", key.m_LinSysMaxStorage, 30);
101 m_session->LoadParameter("LinSysRelativeTolInNonlin",
102 key.m_NekLinSysTolerance, 5.0E-2);
103 m_session->LoadParameter("GMRESMaxHessMatBand", key.m_KrylovMaxHessMatBand,
104 31);
105 m_session->MatchSolverInfo("GMRESLeftPrecon", "True",
106 key.m_NekLinSysLeftPrecon, false);
107 m_session->MatchSolverInfo("GMRESRightPrecon", "True",
108 key.m_NekLinSysRightPrecon, true);
109 int GMRESCentralDifference = 0;
110 m_session->LoadParameter("GMRESCentralDifference", GMRESCentralDifference,
111 0);
112 key.m_GMRESCentralDifference = (bool)GMRESCentralDifference;
113 // Load required NonLinSys parameters:
114 m_session->LoadParameter("NekNonlinSysMaxIterations",
116 m_session->LoadParameter("NewtonRelativeIteTol",
117 key.m_NekNonLinSysTolerance, 1.0E-12);
118 WARNINGL0(!m_session->DefinesParameter("NewtonAbsoluteIteTol"),
119 "Please specify NewtonRelativeIteTol instead of "
120 "NewtonAbsoluteIteTol in XML session file");
121 m_session->LoadParameter("NonlinIterTolRelativeL2",
122 key.m_NonlinIterTolRelativeL2, 1.0E-3);
123 m_session->LoadSolverInfo("LinSysIterSolverTypeInNonlin",
124 key.m_LinSysIterSolverTypeInNonlin, "GMRES");
125
126 // Initialize operator
130 this);
132
133 // Initialize trace
134 const auto locTraceToTraceMap = m_fields[0]->GetLocTraceToTraceMap();
135 locTraceToTraceMap->CalcLocTracePhysToTraceIDMap(m_fields[0]->GetTrace(),
136 m_spacedim);
137 for (int i = 1; i < nvariables; i++)
138 {
139 m_fields[i]->GetLocTraceToTraceMap()->SetLocTracePhysToTraceIDMap(
140 locTraceToTraceMap->GetLocTracephysToTraceIDMap());
141 }
142
143 // Initialize non-linear system
145 "Newton", m_session, m_comm->GetRowComm(), ntotal, key);
146 m_nonlinsol->SetSysOperators(nekSysOp);
147
148 // Initialize preconditioner
149 NekPreconCfsOperators preconOp;
151 this);
154 m_preconCfs->SetOperators(preconOp);
155}
156
164
166 const NekDouble cpuTime)
167{
169
170 if (m_infosteps && m_session->GetComm()->GetSpaceComm()->GetRank() == 0 &&
171 !((step + 1) % m_infosteps) && m_flagImplicitItsStatistics)
172 {
173 std::cout << " &&"
174 << " TotImpStages= " << m_TotImpStages
175 << " TotNewtonIts= " << m_TotNewtonIts
176 << " TotLinearIts = " << m_TotLinIts << std::endl;
177 }
178}
179
181{
183
184 if (m_session->GetComm()->GetRank() == 0 && m_flagImplicitItsStatistics)
185 {
186 std::cout << "-------------------------------------------" << std::endl
187 << "Total Implicit Stages: " << m_TotImpStages << std::endl
188 << "Total Newton Its : " << m_TotNewtonIts << std::endl
189 << "Total Linear Its : " << m_TotLinIts << std::endl
190 << "-------------------------------------------" << std::endl;
191 }
192}
193
196 const bool &flag)
197{
199 unsigned int nvariables = m_fields.size();
200 unsigned int ncoeffs = m_fields[0]->GetNcoeffs();
201 Array<OneD, Array<OneD, NekDouble>> in2D(nvariables);
202 Array<OneD, Array<OneD, NekDouble>> out2D(nvariables);
203 for (int i = 0; i < nvariables; ++i)
204 {
205 int offset = i * ncoeffs;
206 in2D[i] = inarray + offset;
207 out2D[i] = out + offset;
208 }
209
210 timer.Start();
211 NonlinSysEvaluatorCoeff(in2D, out2D, flag);
212 timer.Stop();
213 timer.AccumulateRegion("CFSImplicit::NonlinSysEvaluatorCoeff1D");
214}
215
217 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
218 Array<OneD, Array<OneD, NekDouble>> &out, const bool &flag)
219{
221 unsigned int nvariable = inarray.size();
222 unsigned int ncoeffs = m_fields[0]->GetNcoeffs();
223 unsigned int npoints = m_fields[0]->GetNpoints();
224
225 Array<OneD, Array<OneD, NekDouble>> inpnts(nvariable);
226
227 for (int i = 0; i < nvariable; ++i)
228 {
229 inpnts[i] = Array<OneD, NekDouble>(npoints, 0.0);
230 m_fields[i]->BwdTrans(inarray[i], inpnts[i]);
231 }
232
233 timer.Start();
234 DoOdeProjection(inpnts, inpnts, m_bndEvaluateTime);
235 timer.Stop();
236 timer.AccumulateRegion("CompressibleFlowSystem::DoOdeProjection", 1);
237
238 timer.Start();
239 DoOdeRhsCoeff(inpnts, out, m_bndEvaluateTime);
240 timer.Stop();
241 timer.AccumulateRegion("CFSImplicit::DoOdeRhsCoeff", 1);
242
243 for (int i = 0; i < nvariable; ++i)
244 {
245 Vmath::Svtvp(ncoeffs, -m_TimeIntegLambda, out[i], 1, inarray[i], 1,
246 out[i], 1);
247 if (flag)
248 {
249 Vmath::Vsub(ncoeffs, out[i], 1,
250 m_nonlinsol->GetRefSourceVec() + i * ncoeffs, 1, out[i],
251 1);
252 }
253 }
254}
255
257 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
258 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
259{
260 int nvariables = inarray.size();
261 int ncoeffs = m_fields[0]->GetNcoeffs();
262
263 Array<OneD, Array<OneD, NekDouble>> tmpOut(nvariables);
264 for (int i = 0; i < nvariables; ++i)
265 {
266 tmpOut[i] = Array<OneD, NekDouble>(ncoeffs);
267 }
268
269 DoOdeRhsCoeff(inarray, tmpOut, time);
270
271 for (int i = 0; i < nvariables; ++i)
272 {
273 m_fields[i]->BwdTrans(tmpOut[i], outarray[i]);
274 }
275}
276
277/**
278 * @brief Compute the right-hand side.
279 */
281 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
282 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
283{
284 ASSERTL0(
286 "Do not use Local Time-Stepping with implicit time-discretization");
287
289
290 int nvariables = inarray.size();
291 int nTracePts = GetTraceTotPoints();
292 int ncoeffs = GetNcoeffs();
293
294 m_bndEvaluateTime = time;
295
296 // Store forwards/backwards space along trace space
297 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
298 Array<OneD, Array<OneD, NekDouble>> Bwd(nvariables);
299
301 {
304 }
305 else
306 {
307 for (int i = 0; i < nvariables; ++i)
308 {
309 Fwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
310 Bwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
311 m_fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
312 }
313 m_fields[0]->PeriodicBwdRot(Bwd);
314 if (m_fields[0]->GetGraph()->GetMovement()->GetSectorRotateFlag())
315 {
316 m_fields[0]->RotLocalBwdTrace(Bwd);
317 }
318 }
319
320 timer.Start();
321 DoAdvectionCoeff(inarray, outarray, time, Fwd, Bwd);
322 timer.Stop();
323 timer.AccumulateRegion("CFSImplicit::DoAdvectionCoeff", 2);
324
325 // Negate results
326 for (int i = 0; i < nvariables; ++i)
327 {
328 Vmath::Neg(ncoeffs, outarray[i], 1);
329 }
330
331 // Add diffusion terms
332 timer.Start();
333 DoDiffusionCoeff(inarray, outarray, Fwd, Bwd);
334 timer.Stop();
335 timer.AccumulateRegion("CFSImplicit::DoDiffusionCoeff", 2);
336
337 // Add forcing terms
338 for (auto &x : m_forcing)
339 {
340 x->ApplyCoeff(m_fields, inarray, outarray, time);
341 }
342}
343
344/**
345 * @brief Compute the advection terms for the right-hand side
346 */
348 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
349 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
350 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
351 const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
352{
353 int nvariables = inarray.size();
355
356 auto advWeakDGObject =
357 std::dynamic_pointer_cast<SolverUtils::AdvectionWeakDG>(m_advObject);
358 ASSERTL0(advWeakDGObject,
359 "Use WeakDG for implicit compressible flow solver!");
360 advWeakDGObject->AdvectCoeffs(nvariables, m_fields, advVel, inarray,
361 outarray, time, pFwd, pBwd);
362}
363
364/**
365 * @brief Add the diffusions terms to the right-hand side
366 * Similar to DoDiffusion() but with outarray in coefficient space
367 */
369 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
371 const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
372 const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
373{
374 v_DoDiffusionCoeff(inarray, outarray, pFwd, pBwd);
375}
376
378 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
379 Array<OneD, Array<OneD, NekDouble>> &outpnt, const NekDouble time,
380 const NekDouble lambda)
381{
382 unsigned int nvariables = inpnts.size();
383 unsigned int ncoeffs = m_fields[0]->GetNcoeffs();
384 unsigned int ntotal = nvariables * ncoeffs;
385
386 Array<OneD, NekDouble> inarray(ntotal);
387 Array<OneD, NekDouble> outarray(ntotal);
388 Array<OneD, NekDouble> tmpArray;
389 Array<OneD, Array<OneD, NekDouble>> tmpIn(nvariables);
390 Array<OneD, Array<OneD, NekDouble>> tmpOut(nvariables);
391 Array<OneD, Array<OneD, NekDouble>> tmpoutarray(nvariables);
392
393 // Switch flag to make sure the physical shock capturing AV is updated
395
396 if (m_meshDistorted)
397 {
399 }
400 else
401 {
402 tmpIn = inpnts;
403 }
404
405 for (int i = 0; i < nvariables; ++i)
406 {
407 int noffset = i * ncoeffs;
408 tmpArray = inarray + noffset;
409 m_fields[i]->FwdTrans(tmpIn[i], tmpArray);
410 }
411
412 DoImplicitSolveCoeff(tmpIn, inarray, outarray, time, lambda);
413
414 if (m_meshDistorted)
415 {
416
417 for (int i = 0; i < nvariables; ++i)
418 {
419 tmpOut[i] = Array<OneD, NekDouble>(tmpIn[0].size(), 0.0);
420 tmpoutarray[i] = Array<OneD, NekDouble>(ncoeffs, 0.0);
421 }
422 for (int i = 0; i < nvariables; ++i)
423 {
424 int noffset = i * ncoeffs;
425 tmpArray = outarray + noffset;
426 m_fields[i]->BwdTrans(tmpArray, tmpOut[i]);
427 }
429 for (int i = 0; i < nvariables; ++i)
430 {
431 m_fields[i]->FwdTrans(tmpOut[i], tmpoutarray[i]);
432 m_fields[i]->GeneralMatrixOp(mkey, tmpoutarray[i], outpnt[i]);
433 }
434 }
435 else
436 {
437 for (int i = 0; i < nvariables; ++i)
438 {
439 int noffset = i * ncoeffs;
440 tmpArray = outarray + noffset;
441 m_fields[i]->BwdTrans(tmpArray, outpnt[i]);
442 }
443 }
444}
445
447 const Array<OneD, const Array<OneD, NekDouble>> &inpnts,
449 const NekDouble time, const NekDouble lambda)
450{
451 m_TimeIntegLambda = lambda;
452 m_bndEvaluateTime = time;
453 m_solutionPhys = inpnts;
454 unsigned int ntotal = inarray.size();
455
456 if (m_inArrayNorm < 0.0)
457 {
458 CalcRefValues(inarray);
459
460 m_nonlinsol->SetRhsMagnitude(m_inArrayNorm);
461 }
462
463 m_TotNewtonIts += m_nonlinsol->SolveSystem(ntotal, inarray, out);
464
465 m_TotLinIts += m_nonlinsol->GetNtotLinSysIts();
466
468}
469
471{
472 unsigned int nvariables = m_fields.size();
473 unsigned int ntotal = inarray.size();
474 unsigned int npoints = ntotal / nvariables;
475
476 unsigned int nTotalGlobal = ntotal;
477 m_comm->GetSpaceComm()->AllReduce(nTotalGlobal,
479 unsigned int nTotalDOF = nTotalGlobal / nvariables;
480 NekDouble invTotalDOF = 1.0 / nTotalDOF;
481
482 m_inArrayNorm = 0.0;
483 m_magnitdEstimat = Array<OneD, NekDouble>(nvariables, 0.0);
484
485 for (int i = 0; i < nvariables; ++i)
486 {
487 int offset = i * npoints;
489 Vmath::Dot(npoints, inarray + offset, inarray + offset);
490 }
491 m_comm->GetSpaceComm()->AllReduce(m_magnitdEstimat,
493
494 for (int i = 0; i < nvariables; ++i)
495 {
497 }
498
499 for (int i = 2; i < nvariables - 1; ++i)
500 {
502 }
503
504 for (int i = 2; i < nvariables - 1; ++i)
505 {
507 }
508
509 for (int i = 0; i < nvariables; ++i)
510 {
511 m_magnitdEstimat[i] = sqrt(m_magnitdEstimat[i] * invTotalDOF);
512 }
513 if (m_comm->GetRank() == 0 && m_verbose)
514 {
515 for (int i = 0; i < nvariables; ++i)
516 {
517 std::cout << "m_magnitdEstimat[" << i
518 << "] = " << m_magnitdEstimat[i] << std::endl;
519 }
520 std::cout << "m_inArrayNorm = " << m_inArrayNorm << std::endl;
521 }
522}
523
526 const bool &centralDifferenceFlag)
527{
528 const Array<OneD, const NekDouble> solref = m_nonlinsol->GetRefSolution();
529
530 unsigned int ntotal = inarray.size();
531 NekDouble magninarray = Vmath::Dot(ntotal, inarray, inarray);
532 m_comm->GetSpaceComm()->AllReduce(magninarray,
534 NekDouble eps =
535 m_jacobiFreeEps * sqrt((sqrt(m_inArrayNorm) + 1.0) / magninarray);
536
537 Array<OneD, NekDouble> solplus{ntotal};
538 Array<OneD, NekDouble> resplus{ntotal};
539 Vmath::Svtvp(ntotal, eps, inarray, 1, solref, 1, solplus, 1);
540 NonlinSysEvaluatorCoeff1D(solplus, resplus, !centralDifferenceFlag);
541
542 if (centralDifferenceFlag)
543 {
544 Array<OneD, NekDouble> solminus{ntotal};
545 Array<OneD, NekDouble> resminus{ntotal};
546 Vmath::Svtvp(ntotal, -1.0 * eps, inarray, 1, solref, 1, solminus, 1);
547 NonlinSysEvaluatorCoeff1D(solminus, resminus, false);
548 Vmath::Vsub(ntotal, resplus, 1, resminus, 1, out, 1);
549 Vmath::Smul(ntotal, 0.5 / eps, out, 1, out, 1);
550 }
551 else
552 {
553 const Array<OneD, const NekDouble> resref =
554 m_nonlinsol->GetRefResidual();
555
556 Vmath::Vsub(ntotal, resplus, 1, resref, 1, out, 1);
557 Vmath::Smul(ntotal, 1.0 / eps, out, 1, out, 1);
558 }
559}
560
562 Array<OneD, NekDouble> &outarray,
563 const bool &flag)
564{
565 LibUtilities::Timer timer, Gtimer;
566
567 Gtimer.Start();
568 if (m_preconCfs->UpdatePreconMatCheck(NullNekDouble1DArray,
570 {
571 int nvariables = m_solutionPhys.size();
572 int nphspnt = m_solutionPhys[nvariables - 1].size();
573 Array<OneD, Array<OneD, NekDouble>> intmp(nvariables);
574 for (int i = 0; i < nvariables; i++)
575 {
576 intmp[i] = Array<OneD, NekDouble>(nphspnt, 0.0);
577 }
578
579 timer.Start();
581 timer.Stop();
582 timer.AccumulateRegion("CompressibleFlowSystem::DoOdeProjection", 1);
583
584 timer.Start();
585 m_preconCfs->BuildPreconCfs(m_fields, intmp, m_bndEvaluateTime,
587 timer.Stop();
588 timer.AccumulateRegion("PreconCfsOp::BuildPreconCfs", 1);
589 }
590
591 timer.Start();
592 m_preconCfs->DoPreconCfs(m_fields, inarray, outarray, flag);
593 timer.Stop();
594 timer.AccumulateRegion("PreconCfsOp::DoPreconCfs", 1);
595
596 Gtimer.Stop();
597 Gtimer.AccumulateRegion("CFSImplicit::PreconCoeff");
598}
599
600template <typename DataType, typename TypeNekBlkMatSharedPtr>
602 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
603 const Array<OneD, const TensorOfArray2D<NekDouble>> &qfield,
605 TensorOfArray4D<DataType> &StdMatDataDBB,
606 TensorOfArray5D<DataType> &StdMatDataDBDB)
607{
608 if (StdMatDataDBB.size() == 0)
609 {
610 CalcVolJacStdMat(StdMatDataDBB, StdMatDataDBDB);
611 }
612
613 int nSpaceDim = m_graph->GetSpaceDimension();
614 int nvariable = inarray.size();
615 int npoints = m_fields[0]->GetTotPoints();
616 int nVar2 = nvariable * nvariable;
617 std::shared_ptr<LocalRegions::ExpansionVector> expvect =
618 m_fields[0]->GetExp();
619 int nTotElmt = (*expvect).size();
620
621 Array<OneD, NekDouble> mu(npoints, 0.0);
622 Array<OneD, NekDouble> DmuDT(npoints, 0.0);
624 {
625 CalcMuDmuDT(inarray, mu, DmuDT);
626 }
627
630 for (int i = 0; i < 3; i++)
631 {
632 normal3D[i] = Array<OneD, NekDouble>(3, 0.0);
633 }
634 normal3D[0][0] = 1.0;
635 normal3D[1][1] = 1.0;
636 normal3D[2][2] = 1.0;
638
639 DNekMatSharedPtr wspMat =
640 MemoryManager<DNekMat>::AllocateSharedPtr(nvariable, nvariable, 0.0);
642 nvariable - 1, nvariable, 0.0);
643
644 Array<OneD, DataType> GmatxData;
645 Array<OneD, DataType> MatData;
646
649 m_spacedim); // Nvar*Nvar*Ndir*Nelmt*Npnt
650 TensorOfArray3D<DataType> PntJacConsStd(
651 m_spacedim); // Nvar*Nvar*Ndir*Nelmt*Npnt
655 m_spacedim); // Nvar*Nvar*Ndir*Nelmt*Npnt
656 TensorOfArray4D<DataType> PntJacDervStd(
657 m_spacedim); // Nvar*Nvar*Ndir*Nelmt*Npnt
659 m_spacedim); // Nvar*Nvar*Ndir*Nelmt*Npnt
661 m_spacedim); // Nvar*Nvar*Ndir*Nelmt*Npnt
662 for (int ndir = 0; ndir < m_spacedim; ndir++)
663 {
664 PntJacDerv[ndir] = TensorOfArray3D<NekDouble>(m_spacedim);
665 PntJacDervStd[ndir] = TensorOfArray3D<DataType>(m_spacedim);
668 }
669
671 Array<OneD, NekDouble> locDmuDT;
672 Array<OneD, Array<OneD, NekDouble>> locVars(nvariable);
674 for (int ndir = 0; ndir < m_spacedim; ndir++)
675 {
676 locDerv[ndir] = Array<OneD, Array<OneD, NekDouble>>(nvariable);
677 }
678
679 int nElmtCoefOld = -1;
680 for (int ne = 0; ne < nTotElmt; ne++)
681 {
682 int nElmtCoef = (*expvect)[ne]->GetNcoeffs();
683 int nElmtCoef2 = nElmtCoef * nElmtCoef;
684 int nElmtPnt = (*expvect)[ne]->GetTotPoints();
685
686 int nQuot = nElmtCoef2 / m_nPadding;
687 int nRemd = nElmtCoef2 - nQuot * m_nPadding;
688 int nQuotPlus = nQuot;
689 if (nRemd > 0)
690 {
691 nQuotPlus++;
692 }
693 int nElmtCoef2Paded = nQuotPlus * m_nPadding;
694
695 if (nElmtPnt > PntJacCons[0].size() || nElmtCoef > nElmtCoefOld)
696 {
697 nElmtCoefOld = nElmtCoef;
698 for (int ndir = 0; ndir < 3; ndir++)
699 {
700 normalPnt[ndir] = Array<OneD, NekDouble>(npoints, 0.0);
701 }
702 tmppnts = Array<OneD, NekDouble>(nElmtPnt);
703 MatData = Array<OneD, DataType>(nElmtCoef2Paded * nVar2);
704 for (int ndir = 0; ndir < m_spacedim; ndir++)
705 {
706 ConsCurv[ndir] = Array<OneD, NekDouble>(nElmtPnt);
707 ConsStdd[ndir] = Array<OneD, NekDouble>(nElmtPnt);
708 PntJacCons[ndir] =
710 PntJacConsStd[ndir] =
712 for (int i = 0; i < nElmtPnt; i++)
713 {
714 PntJacCons[ndir][i] = Array<OneD, NekDouble>(nVar2);
715 PntJacConsStd[ndir][i] = Array<OneD, DataType>(nVar2);
716 }
717
718 for (int ndir1 = 0; ndir1 < m_spacedim; ndir1++)
719 {
720 PntJacDerv[ndir][ndir1] =
722 PntJacDervStd[ndir][ndir1] =
724 DervStdd[ndir][ndir1] = Array<OneD, NekDouble>(nElmtPnt);
725 DervCurv[ndir][ndir1] = Array<OneD, NekDouble>(nElmtPnt);
726 for (int i = 0; i < nElmtPnt; i++)
727 {
728 PntJacDerv[ndir][ndir1][i] =
730 PntJacDervStd[ndir][ndir1][i] =
732 }
733 }
734 }
735 }
736
737 int noffset = GetPhys_Offset(ne);
738 for (int j = 0; j < nvariable; j++)
739 {
740 locVars[j] = inarray[j] + noffset;
741 }
742
744 {
745 for (int nFluxDir = 0; nFluxDir < nSpaceDim; nFluxDir++)
746 {
747 normals = normal3D[nFluxDir];
748 GetFluxVectorJacDirElmt(nvariable, nElmtPnt, locVars, normals,
749 wspMat, PntJacCons[nFluxDir]);
750 }
751 }
752
754 {
755 for (int j = 0; j < nSpaceDim; j++)
756 {
757 for (int k = 0; k < nvariable; k++)
758 {
759 locDerv[j][k] = qfield[j][k] + noffset;
760 }
761 }
762 locmu = mu + noffset;
763 locDmuDT = DmuDT + noffset;
764 for (int nFluxDir = 0; nFluxDir < nSpaceDim; nFluxDir++)
765 {
766 normals = normal3D[nFluxDir];
767 MinusDiffusionFluxJacPoint(nvariable, nElmtPnt, locVars,
768 locDerv, locmu, locDmuDT, normals,
769 wspMatDrv, PntJacCons[nFluxDir]);
770 }
771 }
772
774 {
775 locmu = mu + noffset;
776 for (int nFluxDir = 0; nFluxDir < nSpaceDim; nFluxDir++)
777 {
778 Vmath::Fill(npoints, 1.0, normalPnt[nFluxDir], 1);
779 for (int nDervDir = 0; nDervDir < nSpaceDim; nDervDir++)
780 {
782 nvariable, nElmtPnt, nDervDir, locVars, locmu,
783 normalPnt, wspMatDrv, PntJacDerv[nFluxDir][nDervDir]);
784 }
785 Vmath::Fill(npoints, 0.0, normalPnt[nFluxDir], 1);
786 }
787 }
788
789 for (int n = 0; n < nvariable; n++)
790 {
791 for (int m = 0; m < nvariable; m++)
792 {
793 int nVarOffset = m + n * nvariable;
794 GmatxData = gmtxarray[m][n]->GetBlock(ne, ne)->GetPtr();
795
796 for (int ndStd0 = 0; ndStd0 < m_spacedim; ndStd0++)
797 {
798 Vmath::Zero(nElmtPnt, ConsStdd[ndStd0], 1);
799 }
800 for (int ndir = 0; ndir < m_spacedim; ndir++)
801 {
802 for (int i = 0; i < nElmtPnt; i++)
803 {
804 tmppnts[i] = PntJacCons[ndir][i][nVarOffset];
805 }
806 (*expvect)[ne]->AlignVectorToCollapsedDir(ndir, tmppnts,
807 ConsCurv);
808 for (int nd = 0; nd < m_spacedim; nd++)
809 {
810 Vmath::Vadd(nElmtPnt, ConsCurv[nd], 1, ConsStdd[nd], 1,
811 ConsStdd[nd], 1);
812 }
813 }
814
815 for (int ndir = 0; ndir < m_spacedim; ndir++)
816 {
817 (*expvect)[ne]->MultiplyByQuadratureMetric(
818 ConsStdd[ndir], ConsStdd[ndir]); // weight with metric
819 for (int i = 0; i < nElmtPnt; i++)
820 {
821 PntJacConsStd[ndir][i][nVarOffset] =
822 DataType(ConsStdd[ndir][i]);
823 }
824 }
825 }
826 }
827
829 {
830 for (int m = 0; m < nvariable; m++)
831 {
832 for (int n = 0; n < nvariable; n++)
833 {
834 int nVarOffset = m + n * nvariable;
835 for (int ndStd0 = 0; ndStd0 < m_spacedim; ndStd0++)
836 {
837 for (int ndStd1 = 0; ndStd1 < m_spacedim; ndStd1++)
838 {
839 Vmath::Zero(nElmtPnt, DervStdd[ndStd0][ndStd1], 1);
840 }
841 }
842 for (int nd0 = 0; nd0 < m_spacedim; nd0++)
843 {
844 for (int nd1 = 0; nd1 < m_spacedim; nd1++)
845 {
846 for (int i = 0; i < nElmtPnt; i++)
847 {
848 tmppnts[i] =
849 PntJacDerv[nd0][nd1][i][nVarOffset];
850 }
851
852 (*expvect)[ne]->AlignVectorToCollapsedDir(
853 nd0, tmppnts, ConsCurv);
854 for (int nd = 0; nd < m_spacedim; nd++)
855 {
856 (*expvect)[ne]->AlignVectorToCollapsedDir(
857 nd1, ConsCurv[nd], DervCurv[nd]);
858 }
859
860 for (int ndStd0 = 0; ndStd0 < m_spacedim; ndStd0++)
861 {
862 for (int ndStd1 = 0; ndStd1 < m_spacedim;
863 ndStd1++)
864 {
865 Vmath::Vadd(nElmtPnt,
866 DervCurv[ndStd0][ndStd1], 1,
867 DervStdd[ndStd0][ndStd1], 1,
868 DervStdd[ndStd0][ndStd1], 1);
869 }
870 }
871 }
872 }
873 for (int nd0 = 0; nd0 < m_spacedim; nd0++)
874 {
875 for (int nd1 = 0; nd1 < m_spacedim; nd1++)
876 {
877 (*expvect)[ne]->MultiplyByQuadratureMetric(
878 DervStdd[nd0][nd1],
879 DervStdd[nd0][nd1]); // weight with metric
880 for (int i = 0; i < nElmtPnt; i++)
881 {
882 PntJacDervStd[nd0][nd1][i][nVarOffset] =
883 -DataType(DervStdd[nd0][nd1][i]);
884 }
885 }
886 }
887 }
888 }
889 }
890
891 Vmath::Zero(nElmtCoef2Paded * nVar2, MatData, 1);
892 DataType one = 1.0;
893 for (int ndir = 0; ndir < m_spacedim; ndir++)
894 {
895 for (int i = 0; i < nElmtPnt; i++)
896 {
897 Blas::Ger(nElmtCoef2Paded, nVar2, one,
898 &StdMatDataDBB[ne][ndir][i][0], 1,
899 &PntJacConsStd[ndir][i][0], 1, &MatData[0],
900 nElmtCoef2Paded);
901 }
902 }
903
905 {
906 for (int nd0 = 0; nd0 < m_spacedim; nd0++)
907 {
908 for (int nd1 = 0; nd1 < m_spacedim; nd1++)
909 {
910 for (int i = 0; i < nElmtPnt; i++)
911 {
912 Blas::Ger(nElmtCoef2Paded, nVar2, one,
913 &StdMatDataDBDB[ne][nd0][nd1][i][0], 1,
914 &PntJacDervStd[nd0][nd1][i][0], 1,
915 &MatData[0], nElmtCoef2Paded);
916 }
917 }
918 }
919 }
920
922
923 for (int n = 0; n < nvariable; n++)
924 {
925 for (int m = 0; m < nvariable; m++)
926 {
927 int nVarOffset = m + n * nvariable;
928 GmatxData = gmtxarray[m][n]->GetBlock(ne, ne)->GetPtr();
929 Vmath::Vcopy(nElmtCoef2,
930 tmpA = MatData + nVarOffset * nElmtCoef2Paded, 1,
931 GmatxData, 1);
932 }
933 }
934 }
935}
936
937template <typename DataType>
939 TensorOfArray5D<DataType> &StdMatDataDBDB)
940{
941 std::shared_ptr<LocalRegions::ExpansionVector> expvect =
942 m_fields[0]->GetExp();
943 int nTotElmt = (*expvect).size();
944
945 StdMatDataDBB = TensorOfArray4D<DataType>(nTotElmt);
946 StdMatDataDBDB = TensorOfArray5D<DataType>(nTotElmt);
947
948 std::vector<DNekMatSharedPtr> VectStdDerivBase0;
949 std::vector<TensorOfArray3D<DataType>> VectStdDerivBase_Base;
950 std::vector<TensorOfArray4D<DataType>> VectStdDervBase_DervBase;
951 DNekMatSharedPtr MatStdDerivBase0;
954 for (int ne = 0; ne < nTotElmt; ne++)
955 {
957 stdExp = (*expvect)[ne]->GetStdExp();
959 stdExp->DetShapeType(), *stdExp);
960 MatStdDerivBase0 = stdExp->GetStdMatrix(matkey);
961
962 int nTotStdExp = VectStdDerivBase0.size();
963 int nFoundStdExp = -1;
964 for (int i = 0; i < nTotStdExp; i++)
965 {
966 if ((*VectStdDerivBase0[i]) == (*MatStdDerivBase0))
967 {
968 nFoundStdExp = i;
969 }
970 }
971 if (nFoundStdExp >= 0)
972 {
973 StdMatDataDBB[ne] = VectStdDerivBase_Base[nFoundStdExp];
974 StdMatDataDBDB[ne] = VectStdDervBase_DervBase[nFoundStdExp];
975 }
976 else
977 {
978 int nElmtCoef = (*expvect)[ne]->GetNcoeffs();
979 int nElmtCoef2 = nElmtCoef * nElmtCoef;
980 int nElmtPnt = (*expvect)[ne]->GetTotPoints();
981
982 int nQuot = nElmtCoef2 / m_nPadding;
983 int nRemd = nElmtCoef2 - nQuot * m_nPadding;
984 int nQuotPlus = nQuot;
985 if (nRemd > 0)
986 {
987 nQuotPlus++;
988 }
989 int nPaded = nQuotPlus * m_nPadding;
990
991 ArrayStdMat[0] = MatStdDerivBase0;
992 if (m_spacedim > 1)
993 {
995 StdRegions::eDerivBase1, stdExp->DetShapeType(), *stdExp);
996 ArrayStdMat[1] = stdExp->GetStdMatrix(matkey);
997
998 if (m_spacedim > 2)
999 {
1001 stdExp->DetShapeType(),
1002 *stdExp);
1003 ArrayStdMat[2] = stdExp->GetStdMatrix(matkey);
1004 }
1005 }
1006 for (int nd0 = 0; nd0 < m_spacedim; nd0++)
1007 {
1008 ArrayStdMatData[nd0] = ArrayStdMat[nd0]->GetPtr();
1009 }
1010
1012 stdExp->DetShapeType(), *stdExp);
1013 DNekMatSharedPtr BwdMat = stdExp->GetStdMatrix(matkey);
1014 Array<OneD, NekDouble> BwdMatData = BwdMat->GetPtr();
1015
1018
1019 for (int nd0 = 0; nd0 < m_spacedim; nd0++)
1020 {
1021 tmpStdDBB[nd0] = Array<OneD, Array<OneD, DataType>>(nElmtPnt);
1022 for (int i = 0; i < nElmtPnt; i++)
1023 {
1024 tmpStdDBB[nd0][i] = Array<OneD, DataType>(nPaded, 0.0);
1025 for (int nc1 = 0; nc1 < nElmtCoef; nc1++)
1026 {
1027 int noffset = nc1 * nElmtCoef;
1028 for (int nc0 = 0; nc0 < nElmtCoef; nc0++)
1029 {
1030 tmpStdDBB[nd0][i][nc0 + noffset] = DataType(
1031 ArrayStdMatData[nd0][i * nElmtCoef + nc0] *
1032 BwdMatData[i * nElmtCoef + nc1]);
1033 }
1034 }
1035 }
1036
1037 tmpStdDBDB[nd0] = TensorOfArray3D<DataType>(m_spacedim);
1038 for (int nd1 = 0; nd1 < m_spacedim; nd1++)
1039 {
1040 tmpStdDBDB[nd0][nd1] =
1042 for (int i = 0; i < nElmtPnt; i++)
1043 {
1044 tmpStdDBDB[nd0][nd1][i] =
1045 Array<OneD, DataType>(nPaded, 0.0);
1046 for (int nc1 = 0; nc1 < nElmtCoef; nc1++)
1047 {
1048 int noffset = nc1 * nElmtCoef;
1049 for (int nc0 = 0; nc0 < nElmtCoef; nc0++)
1050 {
1051 tmpStdDBDB[nd0][nd1][i][nc0 + noffset] =
1052 DataType(
1053 ArrayStdMatData[nd0]
1054 [i * nElmtCoef + nc0] *
1055 ArrayStdMatData[nd1]
1056 [i * nElmtCoef + nc1]);
1057 }
1058 }
1059 }
1060 }
1061 }
1062 VectStdDerivBase0.push_back(MatStdDerivBase0);
1063 VectStdDerivBase_Base.push_back(tmpStdDBB);
1064 VectStdDervBase_DervBase.push_back(tmpStdDBDB);
1065
1066 StdMatDataDBB[ne] = tmpStdDBB;
1067 StdMatDataDBDB[ne] = tmpStdDBDB;
1068 }
1069 }
1070}
1071
1072template <typename DataType, typename TypeNekBlkMatSharedPtr>
1074 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1079 Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
1080 TensorOfArray5D<DataType> &TraceIPSymJacArray)
1081{
1082 int nvariables = inarray.size();
1083
1084 LibUtilities::Timer timer;
1085 timer.Start();
1086 GetTraceJac(inarray, qfield, TraceJac, TraceJacDeriv, TraceJacDerivSign,
1087 TraceIPSymJacArray);
1088 timer.Stop();
1089 timer.AccumulateRegion("CFSImplicit::GetTraceJac", 10);
1090
1093
1094 timer.Start();
1095 m_advObject->AddTraceJacToMat(nvariables, m_spacedim, m_fields, TraceJac,
1096 gmtxarray, tmpJac, tmpSign);
1097 timer.Stop();
1098 timer.AccumulateRegion("Advection::AddTraceJacToMap", 10);
1099}
1100
1101template <typename DataType, typename TypeNekBlkMatSharedPtr>
1104 TypeNekBlkMatSharedPtr &gmtVar,
1105 [[maybe_unused]] const DataType &tmpDataType)
1106{
1107 int n1d = gmtxarray.size();
1108 int n2d = gmtxarray[0].size();
1109 int nConvectiveFields = n1d;
1110
1111 ASSERTL0(n1d == n2d, "ElmtVarInvMtrx requires n1d==n2d");
1112
1115
1116 gmtxarray[0][0]->GetBlockSizes(rowSizes, colSizes);
1117 int nTotElmt = rowSizes.size();
1118 int nElmtCoef = rowSizes[0] - 1;
1119 int nElmtCoef0 = -1;
1120 int blocksize = -1;
1121
1122 Array<OneD, unsigned int> tmprow(1);
1123 TypeNekBlkMatSharedPtr tmpGmtx;
1124
1125 Array<OneD, DataType> GMatData, ElmtMatData;
1126 Array<OneD, DataType> tmpArray1, tmpArray2;
1127
1128 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
1129 {
1130 int nrows = gmtxarray[0][0]->GetBlock(nelmt, nelmt)->GetRows();
1131 int ncols = gmtxarray[0][0]->GetBlock(nelmt, nelmt)->GetColumns();
1132 ASSERTL0(nrows == ncols, "ElmtVarInvMtrx requires nrows==ncols");
1133
1134 nElmtCoef = nrows;
1135
1136 if (nElmtCoef0 != nElmtCoef)
1137 {
1138 nElmtCoef0 = nElmtCoef;
1139 int nElmtCoefVar = nElmtCoef0 * nConvectiveFields;
1140 blocksize = nElmtCoefVar * nElmtCoefVar;
1141 tmprow[0] = nElmtCoefVar;
1142 AllocateNekBlkMatDig(tmpGmtx, tmprow, tmprow);
1143 GMatData = tmpGmtx->GetBlock(0, 0)->GetPtr();
1144 }
1145
1146 for (int n = 0; n < nConvectiveFields; n++)
1147 {
1148 for (int m = 0; m < nConvectiveFields; m++)
1149 {
1150 ElmtMatData = gmtxarray[m][n]->GetBlock(nelmt, nelmt)->GetPtr();
1151
1152 for (int ncl = 0; ncl < nElmtCoef; ncl++)
1153 {
1154 int Goffset =
1155 (n * nElmtCoef + ncl) * nConvectiveFields * nElmtCoef +
1156 m * nElmtCoef;
1157 int Eoffset = ncl * nElmtCoef;
1158
1159 Vmath::Vcopy(nElmtCoef, tmpArray1 = ElmtMatData + Eoffset,
1160 1, tmpArray2 = GMatData + Goffset, 1);
1161 }
1162 }
1163 }
1164
1165 tmpGmtx->GetBlock(0, 0)->Invert();
1166
1167 for (int m = 0; m < nConvectiveFields; m++)
1168 {
1169 for (int n = 0; n < nConvectiveFields; n++)
1170 {
1171 ElmtMatData = gmtxarray[m][n]->GetBlock(nelmt, nelmt)->GetPtr();
1172
1173 for (int ncl = 0; ncl < nElmtCoef; ncl++)
1174 {
1175 int Goffset =
1176 (n * nElmtCoef + ncl) * nConvectiveFields * nElmtCoef +
1177 m * nElmtCoef;
1178 int Eoffset = ncl * nElmtCoef;
1179
1180 Vmath::Vcopy(nElmtCoef, tmpArray1 = GMatData + Goffset, 1,
1181 tmpArray2 = ElmtMatData + Eoffset, 1);
1182 }
1183 }
1184 }
1185 ElmtMatData = gmtVar->GetBlock(nelmt, nelmt)->GetPtr();
1186 Vmath::Vcopy(blocksize, GMatData, 1, ElmtMatData, 1);
1187 }
1188 return;
1189}
1190
1191template <typename DataType, typename TypeNekBlkMatSharedPtr>
1193 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1196 [[maybe_unused]] Array<OneD, TypeNekBlkMatSharedPtr> &TraceJacDeriv,
1197 [[maybe_unused]] Array<OneD, Array<OneD, DataType>> &TraceJacDerivSign,
1198 TensorOfArray5D<DataType> &TraceIPSymJacArray)
1199{
1200 int nvariables = inarray.size();
1201 int nTracePts = GetTraceTotPoints();
1202
1203 // Store forwards/backwards space along trace space
1204 Array<OneD, Array<OneD, NekDouble>> Fwd(nvariables);
1205 Array<OneD, Array<OneD, NekDouble>> Bwd(nvariables);
1206
1207 TypeNekBlkMatSharedPtr FJac, BJac;
1208 Array<OneD, unsigned int> n_blks1(nTracePts, nvariables);
1209
1210 if (TraceJac.size() > 0)
1211 {
1212 FJac = TraceJac[0];
1213 BJac = TraceJac[1];
1214 }
1215 else
1216 {
1218
1219 AllocateNekBlkMatDig(FJac, n_blks1, n_blks1);
1220 AllocateNekBlkMatDig(BJac, n_blks1, n_blks1);
1221 }
1222
1224 {
1227 }
1228 else
1229 {
1230 for (int i = 0; i < nvariables; ++i)
1231 {
1232 Fwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1233 Bwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1234 m_fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
1235 }
1236 }
1237
1239
1240 NumCalcRiemFluxJac(nvariables, m_fields, AdvVel, inarray, qfield,
1241 m_bndEvaluateTime, Fwd, Bwd, FJac, BJac,
1242 TraceIPSymJacArray);
1243
1244 TraceJac[0] = FJac;
1245 TraceJac[1] = BJac;
1246}
1247
1248template <typename DataType, typename TypeNekBlkMatSharedPtr>
1250 const int nConvectiveFields,
1252 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
1253 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1254 TensorOfArray3D<NekDouble> &qfield, const NekDouble &time,
1255 const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
1256 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
1257 TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac,
1258 [[maybe_unused]] TensorOfArray5D<DataType> &TraceIPSymJacArray)
1259{
1260 const NekDouble PenaltyFactor2 = 0.0;
1261 int nvariables = nConvectiveFields;
1262 int npoints = GetNpoints();
1263 int nTracePts = GetTraceTotPoints();
1264 int nDim = m_spacedim;
1265
1266 Array<OneD, int> nonZeroIndex;
1267
1268 Array<OneD, Array<OneD, NekDouble>> tmpinarry(nvariables);
1269 for (int i = 0; i < nvariables; i++)
1270 {
1271 tmpinarry[i] = Array<OneD, NekDouble>(npoints, 0.0);
1272 Vmath::Vcopy(npoints, inarray[i], 1, tmpinarry[i], 1);
1273 }
1274
1275 // DmuDT of artificial diffusion is neglected
1276 // TODO: to consider the Jacobian of AV seperately
1279
1280 Array<OneD, Array<OneD, NekDouble>> numflux(nvariables);
1281 for (int i = 0; i < nvariables; ++i)
1282 {
1283 numflux[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1284 }
1285
1287 fields[0]->GetTraceMap();
1288 const MultiRegions::InterfaceMapDGSharedPtr InterfaceMap =
1289 fields[0]->GetInterfaceMap();
1290 TensorOfArray3D<NekDouble> qBwd(nDim);
1291 TensorOfArray3D<NekDouble> qFwd(nDim);
1292 if (m_viscousJacFlag)
1293 {
1294 for (int nd = 0; nd < nDim; ++nd)
1295 {
1296 qBwd[nd] = Array<OneD, Array<OneD, NekDouble>>(nConvectiveFields);
1297 qFwd[nd] = Array<OneD, Array<OneD, NekDouble>>(nConvectiveFields);
1298 for (int i = 0; i < nConvectiveFields; ++i)
1299 {
1300 qBwd[nd][i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1301 qFwd[nd][i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1302
1303 fields[i]->GetFwdBwdTracePhys(qfield[nd][i], qFwd[nd][i],
1304 qBwd[nd][i], true, true, false);
1305 TraceMap->GetAssemblyCommDG()->PerformExchange(qFwd[nd][i],
1306 qBwd[nd][i]);
1307 InterfaceMap->ExchangeTrace(qFwd[nd][i], qBwd[nd][i]);
1308 }
1309 }
1310 }
1311
1312 CalcTraceNumericalFlux(nConvectiveFields, nDim, npoints, nTracePts,
1313 PenaltyFactor2, fields, AdvVel, inarray, time,
1314 qfield, Fwd, Bwd, qFwd, qBwd, MuVarTrace,
1315 nonZeroIndex, numflux);
1316
1317 int nFields = nvariables;
1318 Array<OneD, Array<OneD, NekDouble>> plusFwd(nFields), plusBwd(nFields);
1319 Array<OneD, Array<OneD, NekDouble>> Jacvect(nFields);
1320 Array<OneD, Array<OneD, NekDouble>> FwdBnd(nFields);
1321 Array<OneD, Array<OneD, NekDouble>> plusflux(nFields);
1322 for (int i = 0; i < nFields; i++)
1323 {
1324 Jacvect[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1325 plusFwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1326 plusBwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1327 plusflux[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1328 FwdBnd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1329 }
1330
1331 for (int i = 0; i < nFields; i++)
1332 {
1333 Vmath::Vcopy(nTracePts, Fwd[i], 1, plusFwd[i], 1);
1334 Vmath::Vcopy(nTracePts, Bwd[i], 1, plusBwd[i], 1);
1335 }
1336
1337 NekDouble eps = 1.0E-6;
1338
1339 Array<OneD, DataType> tmpMatData;
1340 // Fwd Jacobian
1341 for (int i = 0; i < nFields; i++)
1342 {
1343 NekDouble epsvar = eps * m_magnitdEstimat[i];
1344 NekDouble oepsvar = 1.0 / epsvar;
1345 Vmath::Sadd(nTracePts, epsvar, Fwd[i], 1, plusFwd[i], 1);
1346
1347 if (m_bndConds.size())
1348 {
1349 for (int i = 0; i < nFields; i++)
1350 {
1351 Vmath::Vcopy(nTracePts, plusFwd[i], 1, FwdBnd[i], 1);
1352 }
1353 // Loop over user-defined boundary conditions
1354 for (auto &x : m_bndConds)
1355 {
1356 x->Apply(FwdBnd, tmpinarry, time);
1357 }
1358 }
1359
1360 for (int j = 0; j < nFields; j++)
1361 {
1362 m_fields[j]->FillBwdWithBoundCond(plusFwd[j], plusBwd[j]);
1363 }
1364
1365 CalcTraceNumericalFlux(nConvectiveFields, nDim, npoints, nTracePts,
1366 PenaltyFactor2, fields, AdvVel, inarray, time,
1367 qfield, plusFwd, plusBwd, qFwd, qBwd, MuVarTrace,
1368 nonZeroIndex, plusflux);
1369
1370 for (int n = 0; n < nFields; n++)
1371 {
1372 Vmath::Vsub(nTracePts, plusflux[n], 1, numflux[n], 1, Jacvect[n],
1373 1);
1374 Vmath::Smul(nTracePts, oepsvar, Jacvect[n], 1, Jacvect[n], 1);
1375 }
1376 for (int j = 0; j < nTracePts; j++)
1377 {
1378 tmpMatData = FJac->GetBlock(j, j)->GetPtr();
1379 for (int n = 0; n < nFields; n++)
1380 {
1381 tmpMatData[n + i * nFields] = DataType(Jacvect[n][j]);
1382 }
1383 }
1384
1385 Vmath::Vcopy(nTracePts, Fwd[i], 1, plusFwd[i], 1);
1386 }
1387
1388 // Reset the boundary conditions
1389 if (m_bndConds.size())
1390 {
1391 for (int i = 0; i < nFields; i++)
1392 {
1393 Vmath::Vcopy(nTracePts, Fwd[i], 1, FwdBnd[i], 1);
1394 }
1395 // Loop over user-defined boundary conditions
1396 for (auto &x : m_bndConds)
1397 {
1398 x->Apply(FwdBnd, tmpinarry, time);
1399 }
1400 }
1401
1402 for (int i = 0; i < nFields; i++)
1403 {
1404 Vmath::Vcopy(nTracePts, Bwd[i], 1, plusBwd[i], 1);
1405 }
1406
1407 for (int i = 0; i < nFields; i++)
1408 {
1409 NekDouble epsvar = eps * m_magnitdEstimat[i];
1410 NekDouble oepsvar = 1.0 / epsvar;
1411
1412 Vmath::Sadd(nTracePts, epsvar, Bwd[i], 1, plusBwd[i], 1);
1413
1414 for (int j = 0; j < nFields; j++)
1415 {
1416 m_fields[j]->FillBwdWithBoundCond(Fwd[j], plusBwd[j]);
1417 }
1418
1419 CalcTraceNumericalFlux(nConvectiveFields, nDim, npoints, nTracePts,
1420 PenaltyFactor2, fields, AdvVel, inarray, time,
1421 qfield, Fwd, plusBwd, qFwd, qBwd, MuVarTrace,
1422 nonZeroIndex, plusflux);
1423
1424 for (int n = 0; n < nFields; n++)
1425 {
1426 Vmath::Vsub(nTracePts, plusflux[n], 1, numflux[n], 1, Jacvect[n],
1427 1);
1428 Vmath::Smul(nTracePts, oepsvar, Jacvect[n], 1, Jacvect[n], 1);
1429 }
1430 for (int j = 0; j < nTracePts; j++)
1431 {
1432 tmpMatData = BJac->GetBlock(j, j)->GetPtr();
1433 for (int n = 0; n < nFields; n++)
1434 {
1435 tmpMatData[n + i * nFields] = DataType(Jacvect[n][j]);
1436 }
1437 }
1438
1439 Vmath::Vcopy(nTracePts, Bwd[i], 1, plusBwd[i], 1);
1440 }
1441}
1442
1444 const int nConvectiveFields, [[maybe_unused]] const int nDim,
1445 [[maybe_unused]] const int nPts, const int nTracePts,
1446 [[maybe_unused]] const NekDouble PenaltyFactor2,
1448 const Array<OneD, const Array<OneD, NekDouble>> &AdvVel,
1449 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1450 [[maybe_unused]] const NekDouble time, TensorOfArray3D<NekDouble> &qfield,
1451 const Array<OneD, const Array<OneD, NekDouble>> &vFwd,
1452 const Array<OneD, const Array<OneD, NekDouble>> &vBwd,
1453 [[maybe_unused]] const Array<OneD, const TensorOfArray2D<NekDouble>> &qFwd,
1454 [[maybe_unused]] const Array<OneD, const TensorOfArray2D<NekDouble>> &qBwd,
1455 [[maybe_unused]] const Array<OneD, NekDouble> &MuVarTrace,
1456 Array<OneD, int> &nonZeroIndex,
1457 Array<OneD, Array<OneD, NekDouble>> &traceflux)
1458{
1460 {
1461 auto advWeakDGObject =
1462 std::dynamic_pointer_cast<SolverUtils::AdvectionWeakDG>(
1463 m_advObject);
1464 ASSERTL0(advWeakDGObject,
1465 "Use WeakDG for implicit compressible flow solver!");
1466 advWeakDGObject->AdvectTraceFlux(nConvectiveFields, m_fields, AdvVel,
1467 inarray, traceflux, m_bndEvaluateTime,
1468 vFwd, vBwd);
1469 }
1470 else
1471 {
1472 for (int i = 0; i < nConvectiveFields; i++)
1473 {
1474 traceflux[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1475 }
1476 }
1477
1478 if (m_viscousJacFlag)
1479 {
1480 Array<OneD, Array<OneD, NekDouble>> visflux(nConvectiveFields);
1481 for (int i = 0; i < nConvectiveFields; i++)
1482 {
1483 visflux[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
1484 }
1485
1486 std::string diffName;
1487 m_session->LoadSolverInfo("DiffusionType", diffName, "InteriorPenalty");
1488 if (diffName == "InteriorPenalty")
1489 {
1490 m_diffusion->DiffuseTraceFlux(fields, inarray, qfield,
1492 vFwd, vBwd, nonZeroIndex);
1493 }
1494 else
1495 {
1496 ASSERTL1(false, "LDGNS not yet validated for implicit compressible "
1497 "flow solver");
1498 // For LDGNS, the array size should be nConvectiveFields - 1
1499 Array<OneD, Array<OneD, NekDouble>> inBwd(nConvectiveFields - 1);
1500 Array<OneD, Array<OneD, NekDouble>> inFwd(nConvectiveFields - 1);
1501 for (int i = 0; i < nConvectiveFields - 1; ++i)
1502 {
1503 inBwd[i] = vBwd[i];
1504 inFwd[i] = vFwd[i];
1505 }
1506 m_diffusion->DiffuseTraceFlux(fields, inarray, qfield,
1508 inFwd, inBwd, nonZeroIndex);
1509 }
1510 for (int i = 0; i < nConvectiveFields; i++)
1511 {
1512 Vmath::Vsub(nTracePts, traceflux[i], 1, visflux[i], 1, traceflux[i],
1513 1);
1514 }
1515 }
1516}
1517
1519 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1522 Array<OneD, SNekBlkMatSharedPtr> &TraceJacDeriv,
1523 Array<OneD, Array<OneD, NekSingle>> &TraceJacDerivSign,
1524 TensorOfArray4D<NekSingle> &TraceJacArray,
1525 [[maybe_unused]] TensorOfArray4D<NekSingle> &TraceJacDerivArray,
1526 TensorOfArray5D<NekSingle> &TraceIPSymJacArray)
1527{
1529
1530 if (m_viscousJacFlag)
1531 {
1532 CalcPhysDeriv(inarray, qfield);
1533 }
1534
1535 NekSingle zero = 0.0;
1536 Fill2DArrayOfBlkDiagonalMat(gmtxarray, zero);
1537
1538 LibUtilities::Timer timer;
1539 timer.Start();
1540 AddMatNSBlkDiagVol(inarray, qfield, gmtxarray, m_stdSMatDataDBB,
1542 timer.Stop();
1543 timer.AccumulateRegion("CFSImplicit::AddMatNSBlkDiagVol", 2);
1544
1545 timer.Start();
1546 AddMatNSBlkDiagBnd(inarray, qfield, gmtxarray, TraceJac, TraceJacDeriv,
1547 TraceJacDerivSign, TraceIPSymJacArray);
1548 timer.Stop();
1549 timer.AccumulateRegion("CFSImplicit::AddMatNSBlkDiagBnd", 2);
1550
1551 MultiplyElmtInvMassPlusSource<NekSingle>(gmtxarray, m_TimeIntegLambda);
1552
1553 timer.Start();
1554 ElmtVarInvMtrx(gmtxarray, gmtVar, zero);
1555 timer.Stop();
1556 timer.AccumulateRegion("CFSImplicit::ElmtVarInvMtrx", 2);
1557
1558 TransTraceJacMatToArray(TraceJac, TraceJacArray);
1559}
1560
1562 [[maybe_unused]] const int nConvectiveFields,
1563 [[maybe_unused]] const int nElmtPnt,
1564 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &locVars,
1565 [[maybe_unused]] const TensorOfArray3D<NekDouble> &locDerv,
1566 [[maybe_unused]] const Array<OneD, NekDouble> &locmu,
1567 [[maybe_unused]] const Array<OneD, NekDouble> &locDmuDT,
1568 [[maybe_unused]] const Array<OneD, NekDouble> &normals,
1569 [[maybe_unused]] DNekMatSharedPtr &wspMat,
1570 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
1571{
1572 // Do nothing by default
1573}
1574
1575template <typename DataType, typename TypeNekBlkMatSharedPtr>
1578 const NekDouble dtlamda)
1579{
1581 std::shared_ptr<LocalRegions::ExpansionVector> pexp = explist->GetExp();
1582 int nTotElmt = (*pexp).size();
1583 int nConvectiveFields = m_fields.size();
1584
1585 NekDouble Negdtlamda = -dtlamda;
1586
1587 Array<OneD, NekDouble> pseudotimefactor(nTotElmt, 0.0);
1588 Vmath::Fill(nTotElmt, Negdtlamda, pseudotimefactor, 1);
1589
1590 Array<OneD, DataType> GMatData;
1591 for (int m = 0; m < nConvectiveFields; m++)
1592 {
1593 for (int n = 0; n < nConvectiveFields; n++)
1594 {
1595 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
1596 {
1597 GMatData = gmtxarray[m][n]->GetBlock(nelmt, nelmt)->GetPtr();
1598 DataType factor = DataType(pseudotimefactor[nelmt]);
1599
1600 Vmath::Smul(GMatData.size(), factor, GMatData, 1, GMatData, 1);
1601 }
1602 }
1603 }
1604
1605 DNekMatSharedPtr MassMat;
1606 Array<OneD, NekDouble> BwdMatData, MassMatData, tmp;
1608 Array<OneD, DataType> MassMatDataDataType;
1610
1611 for (int nelmt = 0; nelmt < nTotElmt; nelmt++)
1612 {
1613 int nelmtcoef = GetNcoeffs(nelmt);
1614 int nelmtpnts = GetTotPoints(nelmt);
1615 LibUtilities::ShapeType ElmtTypeNow =
1616 explist->GetExp(nelmt)->DetShapeType();
1617
1618 if (tmp.size() != nelmtcoef || (ElmtTypeNow != ElmtTypePrevious))
1619 {
1621 stdExp = explist->GetExp(nelmt)->GetStdExp();
1623 stdExp->DetShapeType(), *stdExp);
1624
1625 DNekMatSharedPtr BwdMat = stdExp->GetStdMatrix(matkey);
1626 BwdMatData = BwdMat->GetPtr();
1627
1628 if (nelmtcoef != tmp.size())
1629 {
1630 tmp = Array<OneD, NekDouble>(nelmtcoef, 0.0);
1632 nelmtcoef, nelmtcoef, 0.0);
1633 MassMatData = MassMat->GetPtr();
1634 MassMatDataDataType =
1635 Array<OneD, DataType>(nelmtcoef * nelmtcoef);
1636 }
1637
1638 ElmtTypePrevious = ElmtTypeNow;
1639 }
1640
1641 for (int np = 0; np < nelmtcoef; np++)
1642 {
1643 explist->GetExp(nelmt)->IProductWRTBase(BwdMatData + np * nelmtpnts,
1644 tmp);
1645 Vmath::Vcopy(nelmtcoef, tmp, 1, tmp2 = MassMatData + np * nelmtcoef,
1646 1);
1647 }
1648 for (int i = 0; i < MassMatData.size(); i++)
1649 {
1650 MassMatDataDataType[i] = DataType(MassMatData[i]);
1651 }
1652
1653 for (int m = 0; m < nConvectiveFields; m++)
1654 {
1655 GMatData = gmtxarray[m][m]->GetBlock(nelmt, nelmt)->GetPtr();
1656 Vmath::Vadd(MassMatData.size(), MassMatDataDataType, 1, GMatData, 1,
1657 GMatData, 1);
1658 }
1659 }
1660 return;
1661}
1662
1663template <typename DataType, typename TypeNekBlkMatSharedPtr>
1666 TensorOfArray4D<DataType> &TraceJacArray)
1667{
1668 int nFwdBwd, nDiagBlks, nvar0Jac, nvar1Jac;
1669
1672 nFwdBwd = TraceJac.size();
1673 TraceJac[0]->GetBlockSizes(rowSizes, colSizes);
1674 nDiagBlks = rowSizes.size();
1675 nvar0Jac = rowSizes[1] - rowSizes[0];
1676 nvar1Jac = colSizes[1] - colSizes[0];
1677
1678 if (0 == TraceJacArray.size())
1679 {
1680 TraceJacArray = TensorOfArray4D<DataType>(nFwdBwd);
1681 for (int nlr = 0; nlr < nFwdBwd; nlr++)
1682 {
1683 TraceJacArray[nlr] = TensorOfArray3D<DataType>(nvar0Jac);
1684 for (int m = 0; m < nvar0Jac; m++)
1685 {
1686 TraceJacArray[nlr][m] =
1688 for (int n = 0; n < nvar1Jac; n++)
1689 {
1690 TraceJacArray[nlr][m][n] = Array<OneD, DataType>(nDiagBlks);
1691 }
1692 }
1693 }
1694 }
1695
1696 for (int nlr = 0; nlr < nFwdBwd; nlr++)
1697 {
1698 const TypeNekBlkMatSharedPtr tmpMat = TraceJac[nlr];
1699 TensorOfArray3D<DataType> tmpaa = TraceJacArray[nlr];
1700 TranSamesizeBlkDiagMatIntoArray(tmpMat, tmpaa);
1701 }
1702
1703 return;
1704}
1705
1707 [[maybe_unused]] const MultiRegions::ExpListSharedPtr &explist,
1708 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &normals,
1709 [[maybe_unused]] const int nDervDir,
1710 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1711 [[maybe_unused]] TensorOfArray5D<NekDouble> &ElmtJacArray,
1712 [[maybe_unused]] const int nFluxDir)
1713{
1714 NEKERROR(ErrorUtil::efatal, "v_GetFluxDerivJacDirctn not coded");
1715}
1716
1718 [[maybe_unused]] const int nConvectiveFields,
1719 [[maybe_unused]] const int nElmtPnt, [[maybe_unused]] const int nDervDir,
1720 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &locVars,
1721 [[maybe_unused]] const Array<OneD, NekDouble> &locmu,
1722 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &locnormal,
1723 [[maybe_unused]] DNekMatSharedPtr &wspMat,
1724 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
1725{
1726 NEKERROR(ErrorUtil::efatal, "v_GetFluxDerivJacDirctn not coded");
1727}
1728
1730 [[maybe_unused]] const MultiRegions::ExpListSharedPtr &explist,
1731 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &normals,
1732 [[maybe_unused]] const int nDervDir,
1733 [[maybe_unused]] const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1734 [[maybe_unused]] TensorOfArray2D<DNekMatSharedPtr> &ElmtJac)
1735{
1736}
1737
1739 const int nConvectiveFields, const int nElmtPnt,
1740 const Array<OneD, const Array<OneD, NekDouble>> &locVars,
1741 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &wspMat,
1742 Array<OneD, Array<OneD, NekDouble>> &PntJacArray)
1743{
1744 Array<OneD, NekDouble> wspMatData = wspMat->GetPtr();
1745
1746 int matsize = nConvectiveFields * nConvectiveFields;
1747
1748 Array<OneD, NekDouble> pointVar(nConvectiveFields);
1749
1750 for (int npnt = 0; npnt < nElmtPnt; npnt++)
1751 {
1752 for (int j = 0; j < nConvectiveFields; j++)
1753 {
1754 pointVar[j] = locVars[j][npnt];
1755 }
1756
1757 GetFluxVectorJacPoint(nConvectiveFields, pointVar, normals, wspMat);
1758
1759 Vmath::Vcopy(matsize, wspMatData, 1, PntJacArray[npnt], 1);
1760 }
1761 return;
1762}
1763
1765 const int nConvectiveFields, const Array<OneD, NekDouble> &conservVar,
1766 const Array<OneD, NekDouble> &normals, DNekMatSharedPtr &fluxJac)
1767{
1768 int nvariables = conservVar.size();
1769 const int nvariables3D = 5;
1770 int expDim = m_spacedim;
1771
1772 NekDouble fsw, efix_StegerWarming;
1773 efix_StegerWarming = 0.0;
1774 fsw = 0.0; // exact flux Jacobian if fsw=0.0
1775 if (nvariables > expDim + 2)
1776 {
1777 NEKERROR(ErrorUtil::efatal, "nvariables > expDim+2 case not coded")
1778 }
1779
1780 Array<OneD, NekDouble> fluxJacData;
1781 ;
1782 fluxJacData = fluxJac->GetPtr();
1783
1784 if (nConvectiveFields == nvariables3D)
1785 {
1786 PointFluxJacobianPoint(conservVar, normals, fluxJac, efix_StegerWarming,
1787 fsw);
1788 }
1789 else
1790 {
1791 DNekMatSharedPtr PointFJac3D =
1793 nvariables3D, 0.0);
1794
1795 Array<OneD, NekDouble> PointFJac3DData;
1796 PointFJac3DData = PointFJac3D->GetPtr();
1797
1798 Array<OneD, NekDouble> PointFwd(nvariables3D, 0.0);
1799
1800 Array<OneD, unsigned int> index(nvariables);
1801
1802 index[nvariables - 1] = 4;
1803 for (int i = 0; i < nvariables - 1; i++)
1804 {
1805 index[i] = i;
1806 }
1807
1808 int nj = 0;
1809 int nk = 0;
1810 for (int j = 0; j < nvariables; j++)
1811 {
1812 nj = index[j];
1813 PointFwd[nj] = conservVar[j];
1814 }
1815
1816 PointFluxJacobianPoint(PointFwd, normals, PointFJac3D,
1817 efix_StegerWarming, fsw);
1818
1819 for (int j = 0; j < nvariables; j++)
1820 {
1821 nj = index[j];
1822 for (int k = 0; k < nvariables; k++)
1823 {
1824 nk = index[k];
1825 fluxJacData[j + k * nConvectiveFields] =
1826 PointFJac3DData[nj + nk * nvariables3D];
1827 }
1828 }
1829 }
1830}
1831
1832template <typename DataType, typename TypeNekBlkMatSharedPtr>
1834 const TypeNekBlkMatSharedPtr &BlkMat, TensorOfArray3D<DataType> &MatArray)
1835{
1838 BlkMat->GetBlockSizes(rowSizes, colSizes);
1839 int nDiagBlks = rowSizes.size();
1840 int nvar0 = rowSizes[1] - rowSizes[0];
1841 int nvar1 = colSizes[1] - colSizes[0];
1842
1843 Array<OneD, DataType> ElmtMatData;
1844
1845 for (int i = 0; i < nDiagBlks; i++)
1846 {
1847 ElmtMatData = BlkMat->GetBlock(i, i)->GetPtr();
1848 for (int n = 0; n < nvar1; n++)
1849 {
1850 int noffset = n * nvar0;
1851 for (int m = 0; m < nvar0; m++)
1852 {
1853 MatArray[m][n][i] = ElmtMatData[m + noffset];
1854 }
1855 }
1856 }
1857}
1858
1859template <typename DataType, typename TypeNekBlkMatSharedPtr>
1862 const DataType valu)
1863{
1864 int n1d = gmtxarray.size();
1865
1866 for (int n1 = 0; n1 < n1d; ++n1)
1867 {
1868 Fill1DArrayOfBlkDiagonalMat(gmtxarray[n1], valu);
1869 }
1870}
1871
1872template <typename DataType, typename TypeNekBlkMatSharedPtr>
1874 Array<OneD, TypeNekBlkMatSharedPtr> &gmtxarray, const DataType valu)
1875{
1876 int n1d = gmtxarray.size();
1877
1880
1881 Array<OneD, DataType> loc_mat_arr;
1882
1883 for (int n1 = 0; n1 < n1d; ++n1)
1884 {
1885 gmtxarray[n1]->GetBlockSizes(rowSizes, colSizes);
1886 int nelmts = rowSizes.size();
1887
1888 for (int i = 0; i < nelmts; ++i)
1889 {
1890 loc_mat_arr = gmtxarray[n1]->GetBlock(i, i)->GetPtr();
1891
1892 int nrows = gmtxarray[n1]->GetBlock(i, i)->GetRows();
1893 int ncols = gmtxarray[n1]->GetBlock(i, i)->GetColumns();
1894
1895 Vmath::Fill(nrows * ncols, valu, loc_mat_arr, 1);
1896 }
1897 }
1898}
1899
1900// Currently duplacate in compressibleFlowSys
1901// if fsw=+-1 calculate the steger-Warming flux vector splitting flux Jacobian
1902// if fsw=0 calculate the Jacobian of the exact flux
1903// efix is the numerical flux entropy fix parameter
1905 const Array<OneD, NekDouble> &normals,
1906 DNekMatSharedPtr &FJac,
1907 const NekDouble efix,
1908 const NekDouble fsw)
1909{
1910 Array<OneD, NekDouble> FJacData = FJac->GetPtr();
1911 const int nvariables3D = 5;
1912
1913 NekDouble ro, vx, vy, vz, ps, gama, ae;
1914 NekDouble a, a2, h, h0, v2, vn, eps, eps2;
1915 NekDouble nx, ny, nz;
1916 NekDouble sn, osn, nxa, nya, nza, vna;
1917 NekDouble l1, l4, l5, al1, al4, al5, x1, x2, x3, y1;
1918 NekDouble c1, d1, c2, d2, c3, d3, c4, d4, c5, d5;
1919 NekDouble sml_ssf = 1.0E-12;
1920
1921 NekDouble fExactorSplt = 2.0 - abs(fsw); // if fsw=+-1 calculate
1922
1923 NekDouble rhoL = Fwd[0];
1924 NekDouble rhouL = Fwd[1];
1925 NekDouble rhovL = Fwd[2];
1926 NekDouble rhowL = Fwd[3];
1927 NekDouble EL = Fwd[4];
1928
1929 ro = rhoL;
1930 vx = rhouL / rhoL;
1931 vy = rhovL / rhoL;
1932 vz = rhowL / rhoL;
1933
1934 // Internal energy (per unit mass)
1935 NekDouble eL = (EL - 0.5 * (rhouL * vx + rhovL * vy + rhowL * vz)) / rhoL;
1936
1937 ps = m_varConv->Geteos()->GetPressure(rhoL, eL);
1938 gama = m_gamma;
1939
1940 ae = gama - 1.0;
1941 v2 = vx * vx + vy * vy + vz * vz;
1942 a2 = gama * ps / ro;
1943 h = a2 / ae;
1944
1945 h0 = h + 0.5 * v2;
1946 a = sqrt(a2);
1947
1948 nx = normals[0];
1949 ny = normals[1];
1950 nz = normals[2];
1951 vn = nx * vx + ny * vy + nz * vz;
1952 sn = std::max(sqrt(nx * nx + ny * ny + nz * nz), sml_ssf);
1953 osn = 1.0 / sn;
1954
1955 nxa = nx * osn;
1956 nya = ny * osn;
1957 nza = nz * osn;
1958 vna = vn * osn;
1959 l1 = vn;
1960 l4 = vn + sn * a;
1961 l5 = vn - sn * a;
1962
1963 eps = efix * sn;
1964 eps2 = eps * eps;
1965
1966 al1 = sqrt(l1 * l1 + eps2);
1967 al4 = sqrt(l4 * l4 + eps2);
1968 al5 = sqrt(l5 * l5 + eps2);
1969
1970 l1 = 0.5 * (fExactorSplt * l1 + fsw * al1);
1971 l4 = 0.5 * (fExactorSplt * l4 + fsw * al4);
1972 l5 = 0.5 * (fExactorSplt * l5 + fsw * al5);
1973
1974 x1 = 0.5 * (l4 + l5);
1975 x2 = 0.5 * (l4 - l5);
1976 x3 = x1 - l1;
1977 y1 = 0.5 * v2;
1978 c1 = ae * x3 / a2;
1979 d1 = x2 / a;
1980
1981 int nVar0 = 0;
1982 int nVar1 = nvariables3D;
1983 int nVar2 = 2 * nvariables3D;
1984 int nVar3 = 3 * nvariables3D;
1985 int nVar4 = 4 * nvariables3D;
1986 FJacData[nVar0] = c1 * y1 - d1 * vna + l1;
1987 FJacData[nVar1] = -c1 * vx + d1 * nxa;
1988 FJacData[nVar2] = -c1 * vy + d1 * nya;
1989 FJacData[nVar3] = -c1 * vz + d1 * nza;
1990 FJacData[nVar4] = c1;
1991 c2 = c1 * vx + d1 * nxa * ae;
1992 d2 = x3 * nxa + d1 * vx;
1993 FJacData[1 + nVar0] = c2 * y1 - d2 * vna;
1994 FJacData[1 + nVar1] = -c2 * vx + d2 * nxa + l1;
1995 FJacData[1 + nVar2] = -c2 * vy + d2 * nya;
1996 FJacData[1 + nVar3] = -c2 * vz + d2 * nza;
1997 FJacData[1 + nVar4] = c2;
1998 c3 = c1 * vy + d1 * nya * ae;
1999 d3 = x3 * nya + d1 * vy;
2000 FJacData[2 + nVar0] = c3 * y1 - d3 * vna;
2001 FJacData[2 + nVar1] = -c3 * vx + d3 * nxa;
2002 FJacData[2 + nVar2] = -c3 * vy + d3 * nya + l1;
2003 FJacData[2 + nVar3] = -c3 * vz + d3 * nza;
2004 FJacData[2 + nVar4] = c3;
2005 c4 = c1 * vz + d1 * nza * ae;
2006 d4 = x3 * nza + d1 * vz;
2007 FJacData[3 + nVar0] = c4 * y1 - d4 * vna;
2008 FJacData[3 + nVar1] = -c4 * vx + d4 * nxa;
2009 FJacData[3 + nVar2] = -c4 * vy + d4 * nya;
2010 FJacData[3 + nVar3] = -c4 * vz + d4 * nza + l1;
2011 FJacData[3 + nVar4] = c4;
2012 c5 = c1 * h0 + d1 * vna * ae;
2013 d5 = x3 * vna + d1 * h0;
2014 FJacData[4 + nVar0] = c5 * y1 - d5 * vna;
2015 FJacData[4 + nVar1] = -c5 * vx + d5 * nxa;
2016 FJacData[4 + nVar2] = -c5 * vy + d5 * nya;
2017 FJacData[4 + nVar3] = -c5 * vz + d5 * nza;
2018 FJacData[4 + nVar4] = c5 + l1;
2019}
2020
2022{
2023 bool flag =
2024 (m_time + m_timestep > m_fintime && m_fintime > 0.0) ||
2026 return flag || m_preconCfs->UpdatePreconMatCheck(NullNekDouble1DArray,
2028}
2029
2031 int spaceDim, Array<OneD, MultiRegions::ExpListSharedPtr> &fields)
2032{
2033 m_implicitALESolver = true;
2034
2035 if (fields[0]->GetGraph() != nullptr)
2036 {
2037 fields[0]->GetGraph()->GetMovement()->SetImplicitALEFlag(
2039 fields[0]->GetGraph()->GetMovement()->SetMeshDistortedFlag(
2041 }
2042
2043 m_spaceDim = spaceDim;
2044 m_fieldsALE = fields;
2045
2046 // Initialise grid velocities as 0s
2049 for (int i = 0; i < spaceDim; ++i)
2050 {
2051 m_gridVelocity[i] =
2052 Array<OneD, NekDouble>(fields[0]->GetTotPoints(), 0.0);
2054 Array<OneD, NekDouble>(fields[0]->GetTrace()->GetTotPoints(), 0.0);
2055 }
2056 ALEHelper::InitObject(spaceDim, fields);
2057}
2058
2059} // namespace Nektar
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define WARNINGL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
void CalcVolJacStdMat(TensorOfArray4D< DataType > &StdMatDataDBB, TensorOfArray5D< DataType > &StdMatDataDBDB)
void NonlinSysEvaluatorCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &out, const bool &flag)
void NumCalcRiemFluxJac(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &AdvVel, const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, const NekDouble &time, const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, TypeNekBlkMatSharedPtr &FJac, TypeNekBlkMatSharedPtr &BJac, TensorOfArray5D< DataType > &TraceIPSymJacArray)
void DoOdeImplicitRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void AddMatNSBlkDiagVol(const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const TensorOfArray2D< NekDouble > > &qfield, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, TensorOfArray4D< DataType > &StdMatDataDBB, TensorOfArray5D< DataType > &StdMatDataDBDB)
void GetFluxDerivJacDirctnElmt(const int nConvectiveFields, const int nElmtPnt, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const Array< OneD, NekDouble > &locmu, const Array< OneD, const Array< OneD, NekDouble > > &locnormal, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void CalcTraceNumericalFlux(const int nConvectiveFields, const int nDim, const int nPts, const int nTracePts, const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &AdvVel, const Array< OneD, const Array< OneD, NekDouble > > &inarray, const NekDouble time, TensorOfArray3D< NekDouble > &qfield, const Array< OneD, const Array< OneD, NekDouble > > &vFwd, const Array< OneD, const Array< OneD, NekDouble > > &vBwd, const Array< OneD, const TensorOfArray2D< NekDouble > > &qFwd, const Array< OneD, const TensorOfArray2D< NekDouble > > &qBwd, const Array< OneD, NekDouble > &MuVarTrace, Array< OneD, int > &nonZeroIndex, Array< OneD, Array< OneD, NekDouble > > &traceflux)
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inpnts, Array< OneD, Array< OneD, NekDouble > > &outpnt, const NekDouble time, const NekDouble lambda)
void DoImplicitSolveCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inpnts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const NekDouble time, const NekDouble lambda)
void GetFluxVectorJacDirElmt(const int nConvectiveFields, const int nElmtPnt, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void ElmtVarInvMtrx(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, TypeNekBlkMatSharedPtr &gmtVar, const DataType &tmpDatatype)
void CalcPreconMatBRJCoeff(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)
virtual void v_MinusDiffusionFluxJacPoint(const int nConvectiveFields, const int nElmtPnt, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const TensorOfArray3D< NekDouble > &locDerv, const Array< OneD, NekDouble > &locmu, const Array< OneD, NekDouble > &locDmuDT, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void PointFluxJacobianPoint(const Array< OneD, NekDouble > &Fwd, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &FJac, const NekDouble efix, const NekDouble fsw)
void v_DoSolve() override
Virtual function for solve implementation.
CFSImplicit(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void GetTraceJac(const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJac, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJacDeriv, Array< OneD, Array< OneD, DataType > > &TraceJacDerivSign, TensorOfArray5D< DataType > &TraceIPSymJacArray)
void PreconCoeff(const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const bool &flag)
void MultiplyElmtInvMassPlusSource(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const NekDouble dtlamda)
void NonlinSysEvaluatorCoeff1D(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &flag)
void CalcRefValues(const Array< OneD, const NekDouble > &inarray)
void TransTraceJacMatToArray(const Array< OneD, TypeNekBlkMatSharedPtr > &TraceJac, TensorOfArray4D< DataType > &TraceJacDerivArray)
void v_PrintStatusInformation(const int step, const NekDouble cpuTime) override
Print Status Information.
void CalcPhysDeriv(const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield)
void MinusDiffusionFluxJacPoint(const int nConvectiveFields, const int nElmtPnt, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const TensorOfArray3D< NekDouble > &locDerv, const Array< OneD, NekDouble > &locmu, const Array< OneD, NekDouble > &locDmuDT, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void CalcMuDmuDT(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &DmuDT)
Array< OneD, Array< OneD, NekDouble > > m_solutionPhys
void DoOdeRhsCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side.
void DoDiffusionCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, const Array< OneD, NekDouble > > &pFwd, const Array< OneD, const Array< OneD, NekDouble > > &pBwd)
Add the diffusions terms to the right-hand side Similar to DoDiffusion() but with outarray in coeffic...
void AllocateNekBlkMatDig(SNekBlkMatSharedPtr &mat, const Array< OneD, unsigned int > nrow, const Array< OneD, unsigned int > ncol)
virtual void v_GetFluxDerivJacDirctnElmt(const int nConvectiveFields, const int nElmtPnt, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &locVars, const Array< OneD, NekDouble > &locmu, const Array< OneD, const Array< OneD, NekDouble > > &locnormal, DNekMatSharedPtr &wspMat, Array< OneD, Array< OneD, NekDouble > > &PntJacArray)
void v_ALEInitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
TensorOfArray4D< NekSingle > m_stdSMatDataDBB
void TranSamesizeBlkDiagMatIntoArray(const TypeNekBlkMatSharedPtr &BlkMat, TensorOfArray3D< DataType > &MatArray)
virtual void v_GetFluxDerivJacDirctn(const MultiRegions::ExpListSharedPtr &explist, const Array< OneD, const Array< OneD, NekDouble > > &normals, const int nDervDir, const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray5D< NekDouble > &ElmtJacArray, const int nFluxDir)
void v_PrintSummaryStatistics(const NekDouble intTime) override
Print Summary Statistics.
virtual void v_DoDiffusionCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, const Array< OneD, NekDouble > > &pFwd, const Array< OneD, const Array< OneD, NekDouble > > &pBwd)
void GetFluxVectorJacPoint(const int nConvectiveFields, const Array< OneD, NekDouble > &conservVar, const Array< OneD, NekDouble > &normals, DNekMatSharedPtr &fluxJac)
TensorOfArray5D< NekSingle > m_stdSMatDataDBDB
NekDouble m_TimeIntegLambda
coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(us...
void Fill2DArrayOfBlkDiagonalMat(Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const DataType valu)
Array< OneD, NekDouble > m_magnitdEstimat
Estimate the magnitude of each conserved varibles.
void AddMatNSBlkDiagBnd(const Array< OneD, const Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray2D< TypeNekBlkMatSharedPtr > &gmtxarray, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJac, Array< OneD, TypeNekBlkMatSharedPtr > &TraceJacDeriv, Array< OneD, Array< OneD, DataType > > &TraceJacDerivSign, TensorOfArray5D< DataType > &TraceIPSymJacArray)
void MatrixMultiplyMatrixFreeCoeff(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out, const bool &centralDifferenceFlag)
LibUtilities::NekNonlinSysIterSharedPtr m_nonlinsol
void Fill1DArrayOfBlkDiagonalMat(Array< OneD, TypeNekBlkMatSharedPtr > &gmtxarray, const DataType valu)
void v_InitObject(bool DeclareFields=true) override
Initialization object for CFSImplicit class.
void DoAdvectionCoeff(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const Array< OneD, const Array< OneD, NekDouble > > &pFwd, const Array< OneD, const Array< OneD, NekDouble > > &pBwd)
Compute the advection terms for the right-hand side.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection and call the method for imposing the boundary conditions in case of discontinu...
SolverUtils::DiffusionSharedPtr m_diffusion
void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.
std::vector< CFSBndCondSharedPtr > m_bndConds
VariableConverterSharedPtr m_varConv
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
std::string m_LinSysIterSolverTypeInNonlin
Definition NekSys.h:230
void DefineNekSysResEval(FuncPointerT func, ObjectPointerT obj)
Definition NekSys.h:100
void DefineNekSysLhsEval(FuncPointerT func, ObjectPointerT obj)
Definition NekSys.h:107
void DefineNekSysPrecon(FuncPointerT func, ObjectPointerT obj)
Definition NekSys.h:114
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition Timer.cpp:70
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Describes a matrix with ordering defined by a local to global map.
void DefineCalcPreconMatBRJCoeff(FuncPointerT func, ObjectPointerT obj)
Definition PreconCfsOp.h:88
Array< OneD, MultiRegions::ExpListSharedPtr > m_fieldsALE
Definition ALEHelper.h:140
SOLVER_UTILS_EXPORT void InitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Definition ALEHelper.cpp:48
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Definition ALEHelper.h:142
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
Definition ALEHelper.h:141
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
int m_spacedim
Spatial dimension (>= expansion dim).
SpatialDomains::MeshGraphSharedPtr m_graph
Pointer to graph defining mesh.
LibUtilities::CommSharedPtr m_comm
Communicator.
NekDouble m_timestep
Time step size.
int m_infosteps
Number of time steps between outputting status information.
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNpoints()
NekDouble m_fintime
Finish time of the simulation.
SOLVER_UTILS_EXPORT int GetNcoeffs()
NekDouble m_checktime
Time between checkpoints.
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTotPoints()
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics(const NekDouble intTime)
Print Summary Statistics.
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation(const int step, const NekDouble cpuTime)
Print Status Information.
SOLVER_UTILS_EXPORT void v_DoSolve() override
Solves an unsteady problem.
static void Ger(const int &m, const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy, double *a, const int &lda)
BLAS level 2: Matrix vector multiply A = alpha*x*y**T + A where A[m x n].
Definition Blas.hpp:335
NekNonlinSysIterFactory & GetNekNonlinSysIterFactory()
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
std::shared_ptr< InterfaceMapDG > InterfaceMapDGSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static const NekDouble kNekMachineEpsilon
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition MeshGraph.h:217
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
PreconCfsFactory & GetPreconCfsFactory()
Declaration of the boundary condition factory singleton.
Definition PreconCfs.cpp:41
std::shared_ptr< SNekBlkMat > SNekBlkMatSharedPtr
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
static Array< OneD, Array< OneD, Array< OneD, NekDouble > > > NullNekDoubleTensorOfArray3D
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
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 Neg(int n, T *x, const int incx)
Negate x = -x.
Definition Vmath.hpp:292
T Dot(int n, const T *w, const T *x)
dot product
Definition Vmath.hpp:761
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 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 Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition Vmath.hpp:194
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825
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:295
scalarT< T > sqrt(scalarT< T > in)
Definition scalar.hpp:290