Nektar++
LinearisedAdvection.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LinearisedAdvection.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: Evaluation of the linearised advective term
32//
33///////////////////////////////////////////////////////////////////////////////
34
37#include <boost/format.hpp>
38
39using namespace std;
40
41namespace Nektar
42{
45 "Linearised", LinearisedAdvection::create,
46 "Linearised Non-Conservative");
47
48/**
49 * Constructor. Creates ...
50 *
51 * \param
52 * \param
53 */
54
56{
57}
58
62{
63 Advection::v_InitObject(pSession, pFields);
64
65 m_session = pSession;
68 m_session, pFields[0]->GetGraph());
69 m_spacedim = pFields[0]->GetGraph()->GetSpaceDimension();
70 m_expdim = pFields[0]->GetGraph()->GetMeshDimension();
71
72 // Setting parameters for homogeneous problems
73 m_HomoDirec = 0;
74 m_useFFT = false;
76 m_singleMode = false;
77 m_halfMode = false;
78 m_multipleModes = false;
79
80 if (m_session->DefinesSolverInfo("HOMOGENEOUS"))
81 {
82 std::string HomoStr = m_session->GetSolverInfo("HOMOGENEOUS");
83 m_spacedim = 3;
84
85 if ((HomoStr == "HOMOGENEOUS1D") || (HomoStr == "Homogeneous1D") ||
86 (HomoStr == "1D") || (HomoStr == "Homo1D"))
87 {
89 m_LhomZ = m_session->GetParameter("LZ");
90 m_HomoDirec = 1;
91
92 ASSERTL0(m_session->DefinesSolverInfo("ModeType"),
93 "Need to specify ModeType as HalfMode,SingleMode or "
94 "MultipleModes");
95
96 m_session->MatchSolverInfo("ModeType", "SingleMode", m_singleMode,
97 false);
98 m_session->MatchSolverInfo("ModeType", "HalfMode", m_halfMode,
99 false);
100 m_session->MatchSolverInfo("ModeType", "MultipleModes",
101 m_multipleModes, false);
102
103 if (m_singleMode)
104 {
105 m_npointsZ = 2;
106 }
107 else if (m_halfMode)
108 {
109 m_npointsZ = 1;
110 }
111 else if (m_multipleModes)
112 {
113 m_npointsZ = m_session->GetParameter("HomModesZ");
114 }
115 }
116
117 if ((HomoStr == "HOMOGENEOUS2D") || (HomoStr == "Homogeneous2D") ||
118 (HomoStr == "2D") || (HomoStr == "Homo2D"))
119 {
121 m_session->LoadParameter("HomModesY", m_npointsY);
122 m_session->LoadParameter("LY", m_LhomY);
123 m_session->LoadParameter("HomModesZ", m_npointsZ);
124 m_session->LoadParameter("LZ", m_LhomZ);
125 m_HomoDirec = 2;
126 }
127
128 if ((HomoStr == "HOMOGENEOUS3D") || (HomoStr == "Homogeneous3D") ||
129 (HomoStr == "3D") || (HomoStr == "Homo3D"))
130 {
132 m_session->LoadParameter("HomModesX", m_npointsX);
133 m_session->LoadParameter("LX", m_LhomX);
134 m_session->LoadParameter("HomModesY", m_npointsY);
135 m_session->LoadParameter("LY", m_LhomY);
136 m_session->LoadParameter("HomModesZ", m_npointsZ);
137 m_session->LoadParameter("LZ", m_LhomZ);
138 m_HomoDirec = 3;
139 }
140
141 if (m_session->DefinesSolverInfo("USEFFT"))
142 {
143 m_useFFT = true;
144 }
145 }
146 else
147 {
148 m_npointsZ = 1; // set to default value so can use to identify 2d or 3D
149 // (homogeneous) expansions
150 }
151
152 size_t nvar = m_session->GetVariables().size();
154 for (size_t i = 0; i < nvar; ++i)
155 {
156 m_baseflow[i] = Array<OneD, NekDouble>(pFields[i]->GetTotPoints(), 0.0);
157 }
158
159 size_t nBaseDerivs = (m_halfMode || m_singleMode) ? 2 : m_spacedim;
160 m_gradBase = Array<OneD, Array<OneD, NekDouble>>(nvar * nBaseDerivs);
161 for (size_t i = 0; i < nvar; ++i)
162 {
163 for (size_t j = 0; j < nBaseDerivs; ++j)
164 {
165 m_gradBase[i * nBaseDerivs + j] =
166 Array<OneD, NekDouble>(pFields[i]->GetTotPoints(), 0.0);
167 }
168 }
169
170 ASSERTL0(m_session->DefinesFunction("BaseFlow"),
171 "Base flow must be defined for linearised forms.");
172 string file = m_session->GetFunctionFilename("BaseFlow", 0);
173
174 // Periodic base flows
175 if (m_session->DefinesParameter("N_slices"))
176 {
177 m_session->LoadParameter("N_slices", m_slices);
178 if (m_slices > 1)
179 {
180 ASSERTL0(m_session->GetFunctionType("BaseFlow", 0) ==
182 "Base flow should be a sequence of files.");
183 m_session->LoadParameter("BaseFlow_interporder", m_interporder, 0);
185 "BaseFlow_interporder should be smaller than or equal to "
186 "N_slices.");
188 m_session->LoadParameter("N_start", m_start, 0);
189 m_session->LoadParameter("N_skip", m_skip, 1);
190 DFT(file, pFields);
191 }
192 else
193 {
194 ASSERTL0(false, "Number of slices must be a positive number "
195 "greater than 1");
196 }
197 }
198 // Steady base-flow
199 else
200 {
201 m_slices = 1;
202
203 // BaseFlow from file
204 if (m_session->GetFunctionType("BaseFlow", m_session->GetVariable(0)) ==
206 {
207 ImportFldBase(file, pFields, 1);
208 }
209 // analytic base flow
210 else
211 {
212 size_t nq = pFields[0]->GetNpoints();
216
217 // get the coordinates (assuming all fields have the same
218 // discretisation)
219 pFields[0]->GetCoords(x0, x1, x2);
220
221 for (size_t i = 0; i < pFields.size(); i++)
222 {
224 m_session->GetFunction("BaseFlow", i);
225
226 ifunc->Evaluate(x0, x1, x2, m_baseflow[i]);
227 }
228 }
229 }
230
231 for (size_t i = 0; i < nvar; ++i)
232 {
233 UpdateGradBase(i, pFields[i]);
234 }
235
236 if (m_session->DefinesParameter("period"))
237 {
238 m_period = m_session->GetParameter("period");
239 }
240 else
241 {
242 m_period =
243 (m_session->GetParameter("TimeStep") * m_slices) / (m_slices - 1.);
244 }
245 if (m_session->GetComm()->GetRank() == 0)
246 {
247 cout << "baseflow info : interpolation order " << m_interporder
248 << ", period " << m_period << ", periodicity ";
249 if (m_isperiodic)
250 {
251 cout << "yes\n";
252 }
253 else
254 {
255 cout << "no\n";
256 }
257 cout << "baseflow info : files from " << m_start << " to "
258 << (m_start + (m_slices - 1) * m_skip) << " (skip " << m_skip
259 << ") with " << (m_slices - (m_interporder > 1))
260 << " time intervals" << endl;
261 }
262}
263
265{
266}
267
268// Advection function
269
271 const int nConvectiveFields,
273 const Array<OneD, Array<OneD, NekDouble>> &advVel,
274 const Array<OneD, Array<OneD, NekDouble>> &inarray,
275 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble &time,
276 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &pFwd,
277 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &pBwd)
278{
279 ASSERTL1(nConvectiveFields == inarray.size(),
280 "Number of convective fields and Inarray are not compatible");
281
282 size_t nPointsTot = fields[0]->GetNpoints();
283 size_t ndim = advVel.size();
284 size_t nBaseDerivs = (m_halfMode || m_singleMode) ? 2 : m_spacedim;
285 size_t nDerivs = (m_halfMode) ? 2 : m_spacedim;
286
288 for (size_t i = 0; i < ndim; ++i)
289 {
290 if (fields[i]->GetWaveSpace() && !m_singleMode && !m_halfMode)
291 {
292 velocity[i] = Array<OneD, NekDouble>(nPointsTot, 0.0);
293 fields[i]->HomogeneousBwdTrans(nPointsTot, advVel[i], velocity[i]);
294 }
295 else
296 {
297 velocity[i] = advVel[i];
298 }
299 }
300
302 for (size_t i = 0; i < nDerivs; ++i)
303 {
304 grad[i] = Array<OneD, NekDouble>(nPointsTot);
305 }
306
307 // Evaluation of the base flow for periodic cases
308 if (m_slices > 1)
309 {
310 for (size_t i = 0; i < ndim; ++i)
311 {
312 UpdateBase(m_interp[i], m_baseflow[i], time);
313 UpdateGradBase(i, fields[i]);
314 }
315 }
316
317 // Evaluate the linearised advection term
318 for (size_t i = 0; i < (size_t)nConvectiveFields; ++i)
319 {
320 // Calculate gradient
321 switch (nDerivs)
322 {
323 case 1:
324 {
325 fields[i]->PhysDeriv(inarray[i], grad[0]);
326 }
327 break;
328 case 2:
329 {
330 fields[i]->PhysDeriv(inarray[i], grad[0], grad[1]);
331 }
332 break;
333 case 3:
334 {
335 fields[i]->PhysDeriv(inarray[i], grad[0], grad[1], grad[2]);
336 if (m_multipleModes)
337 {
338 // transform gradients into physical Fourier space
339 fields[i]->HomogeneousBwdTrans(nPointsTot, grad[0],
340 grad[0]);
341 fields[i]->HomogeneousBwdTrans(nPointsTot, grad[1],
342 grad[1]);
343 fields[i]->HomogeneousBwdTrans(nPointsTot, grad[2],
344 grad[2]);
345 }
346 }
347 break;
348 }
349
350 // Calculate U_j du'_i/dx_j
351 Vmath::Vmul(nPointsTot, grad[0], 1, m_baseflow[0], 1, outarray[i], 1);
352 for (size_t j = 1; j < nDerivs; ++j)
353 {
354 Vmath::Vvtvp(nPointsTot, grad[j], 1, m_baseflow[j], 1, outarray[i],
355 1, outarray[i], 1);
356 }
357
358 // Add u'_j dU_i/dx_j
359 size_t lim = (m_halfMode || m_singleMode) ? 2 : ndim;
360 if (m_halfMode && i == 2)
361 {
362 lim = 0;
363 }
364 for (size_t j = 0; j < lim; ++j)
365 {
366 Vmath::Vvtvp(nPointsTot, m_gradBase[i * nBaseDerivs + j], 1,
367 velocity[j], 1, outarray[i], 1, outarray[i], 1);
368 }
369
370 if (m_multipleModes)
371 {
372 fields[i]->HomogeneousFwdTrans(nPointsTot, outarray[i],
373 outarray[i]);
374 }
375 Vmath::Neg(nPointsTot, outarray[i], 1);
376 }
377}
378
380 const Array<OneD, Array<OneD, NekDouble>> &inarray,
382{
383 if (m_session->GetSolverInfo("EqType") == "UnsteadyNavierStokes")
384 {
385 // The SFD method is only applied to the velocity variables in
386 // incompressible
387 ASSERTL1(inarray.size() == (m_baseflow.size() - 1),
388 "Number of base flow variables does not match what is "
389 "expected.");
390 }
391 else
392 {
393 ASSERTL1(
394 inarray.size() == (m_baseflow.size()),
395 "Number of base flow variables does not match what is expected.");
396 }
397
398 size_t npts = inarray[0].size();
399
400 for (size_t i = 0; i < inarray.size(); ++i)
401 {
402 ASSERTL1(npts == m_baseflow[i].size(),
403 "Size of base flow array does not match expected.");
404 Vmath::Vcopy(npts, inarray[i], 1, m_baseflow[i], 1);
405 UpdateGradBase(i, fields[i]);
406 }
407}
408
409/**
410 * Import field from infile and load into \a m_fields. This routine will
411 * also perform a \a BwdTrans to ensure data is in both the physical and
412 * coefficient storage.
413 * @param pInFile Filename to read.
414 * @param pFields Array of expansion lists
415 */
417 std::string pInfile, Array<OneD, MultiRegions::ExpListSharedPtr> &pFields,
418 int pSlice)
419{
420 std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef;
421 std::vector<std::vector<NekDouble>> FieldData;
422
423 size_t nqtot = m_baseflow[0].size();
424 Array<OneD, NekDouble> tmp_coeff(pFields[0]->GetNcoeffs(), 0.0);
425
426 size_t numexp = pFields[0]->GetExpSize();
427 Array<OneD, int> ElementGIDs(numexp);
428
429 // Define list of global element ids
430 for (size_t i = 0; i < numexp; ++i)
431 {
432 ElementGIDs[i] = pFields[0]->GetExp(i)->GetGeom()->GetGlobalID();
433 }
434
435 // Get Homogeneous
438 fld->Import(pInfile, FieldDef, FieldData,
440
441 size_t nSessionVar = m_session->GetVariables().size();
442 size_t nSessionConvVar = nSessionVar - 1;
443 size_t nFileVar = FieldDef[0]->m_fields.size();
444
445 std::unordered_map<int, int> zIdToPlane;
447 {
448 zIdToPlane[0] = 0;
449 }
450
451 for (size_t j = 0; j < nFileVar; ++j)
452 {
453 size_t k = 0;
454 for (; k < nSessionConvVar; ++k)
455 {
456 if (m_session->GetVariable(k) == FieldDef[0]->m_fields[j])
457 {
458 break;
459 }
460 }
461 if (k == nSessionConvVar)
462 {
463 continue;
464 }
465 for (size_t i = 0; i < FieldDef.size(); ++i)
466 {
467 pFields[j]->ExtractDataToCoeffs(FieldDef[i], FieldData[i],
468 FieldDef[i]->m_fields[j], tmp_coeff,
469 zIdToPlane);
470 }
471
473 {
474 pFields[j]->GetPlane(0)->BwdTrans(tmp_coeff, m_baseflow[k]);
475
476 if (m_singleMode)
477 {
478 // copy the bwd trans into the second plane for single
479 // Mode Analysis
480 int ncplane = (pFields[0]->GetNpoints()) / m_npointsZ;
481 Vmath::Vcopy(ncplane, &m_baseflow[k][0], 1,
482 &m_baseflow[k][ncplane], 1);
483 }
484 }
485 else // fully 3D base flow - put in physical space.
486 {
487 bool oldwavespace = pFields[j]->GetWaveSpace();
488 pFields[j]->SetWaveSpace(false);
489 pFields[j]->BwdTrans(tmp_coeff, m_baseflow[k]);
490 pFields[j]->SetWaveSpace(oldwavespace);
491 }
492 }
493
494 // If time-periodic, put loaded data into the slice storage.
495 if (m_slices > 1)
496 {
497 for (size_t i = 0; i < nSessionConvVar; ++i)
498 {
499 Vmath::Vcopy(nqtot, &m_baseflow[i][0], 1,
500 &m_interp[i][pSlice * nqtot], 1);
501 }
502 }
503}
504
506 const Array<OneD, const NekDouble> &inarray,
507 Array<OneD, NekDouble> &outarray, const NekDouble time)
508{
509 int npoints = m_baseflow[0].size();
510 if (m_isperiodic)
511 {
512 NekDouble BetaT = 2 * M_PI * fmod(time, m_period) / m_period;
513 NekDouble phase;
514 Array<OneD, NekDouble> auxiliary(npoints);
515
516 Vmath::Vcopy(npoints, &inarray[0], 1, &outarray[0], 1);
517 Vmath::Svtvp(npoints, cos(0.5 * m_slices * BetaT), &inarray[npoints], 1,
518 &outarray[0], 1, &outarray[0], 1);
519
520 for (int i = 2; i < m_slices; i += 2)
521 {
522 phase = (i >> 1) * BetaT;
523
524 Vmath::Svtvp(npoints, cos(phase), &inarray[i * npoints], 1,
525 &outarray[0], 1, &outarray[0], 1);
526 Vmath::Svtvp(npoints, -sin(phase), &inarray[(i + 1) * npoints], 1,
527 &outarray[0], 1, &outarray[0], 1);
528 }
529 }
530 else
531 {
532 NekDouble x = time;
533 x = x / m_period * (m_slices - 1);
534 int ix = x;
535 if (ix < 0)
536 {
537 ix = 0;
538 }
539 if (ix > m_slices - 2)
540 {
541 ix = m_slices - 2;
542 }
543 int padleft = m_interporder / 2 - 1;
544 if (padleft > ix)
545 {
546 padleft = ix;
547 }
548 int padright = m_interporder - 1 - padleft;
549 if (padright > m_slices - 1 - ix)
550 {
551 padright = m_slices - 1 - ix;
552 }
553 padleft = m_interporder - 1 - padright;
555 for (int i = 0; i < m_interporder; ++i)
556 {
557 for (int j = 0; j < m_interporder; ++j)
558 {
559 if (i != j)
560 {
561 coeff[i] *= (x - ix + padleft - (NekDouble)j) /
562 ((NekDouble)i - (NekDouble)j);
563 }
564 }
565 }
566 Vmath::Zero(npoints, &outarray[0], 1);
567 for (int i = ix - padleft; i < ix + padright + 1; ++i)
568 {
569 Vmath::Svtvp(npoints, coeff[i - ix + padleft],
570 &inarray[i * npoints], 1, &outarray[0], 1,
571 &outarray[0], 1);
572 }
573 }
574}
575
577 const int var, const MultiRegions::ExpListSharedPtr &field)
578{
579 int nBaseDerivs = (m_halfMode || m_singleMode) ? 2 : m_spacedim;
580 switch (m_spacedim)
581 {
582 case 1: // 1D
583 {
584 field->PhysDeriv(m_baseflow[var],
585 m_gradBase[var * nBaseDerivs + 0]);
586 }
587 break;
588 case 2: // 2D
589 {
590 field->PhysDeriv(m_baseflow[var], m_gradBase[var * nBaseDerivs + 0],
591 m_gradBase[var * nBaseDerivs + 1]);
592 }
593 break;
594 case 3:
595 {
596 if (m_halfMode) // can assume W = 0 and d/dz == 0
597 {
598 if (var < 2)
599 {
600 field->PhysDeriv(m_baseflow[var],
601 m_gradBase[var * nBaseDerivs + 0],
602 m_gradBase[var * nBaseDerivs + 1]);
603 }
604 }
605 else if (m_singleMode) // single mode where d/dz = 0
606 {
607 field->PhysDeriv(m_baseflow[var],
608 m_gradBase[var * nBaseDerivs + 0],
609 m_gradBase[var * nBaseDerivs + 1]);
610 }
611 else
612 {
613 // Differentiate base flow in physical space
614 bool oldwavespace = field->GetWaveSpace();
615 field->SetWaveSpace(false);
616 field->PhysDeriv(m_baseflow[var],
617 m_gradBase[var * nBaseDerivs + 0],
618 m_gradBase[var * nBaseDerivs + 1],
619 m_gradBase[var * nBaseDerivs + 2]);
620 field->SetWaveSpace(oldwavespace);
621 }
622 }
623 break;
624 }
625}
626
628 [[maybe_unused]] FloquetMatType mattype,
629 [[maybe_unused]] bool UseContCoeffs) const
630{
631 DNekMatSharedPtr loc_mat;
632 DNekBlkMatSharedPtr BlkMatrix;
633 size_t n_exp = 0;
634
635 n_exp = m_baseflow[0].size(); // will operatore on m_phys
636
637 Array<OneD, unsigned int> nrows(n_exp);
638 Array<OneD, unsigned int> ncols(n_exp);
639
640 nrows = Array<OneD, unsigned int>(n_exp, m_slices);
641 ncols = Array<OneD, unsigned int>(n_exp, m_slices);
642
643 MatrixStorage blkmatStorage = eDIAGONAL;
644 BlkMatrix = MemoryManager<DNekBlkMat>::AllocateSharedPtr(nrows, ncols,
645 blkmatStorage);
646
650 StdRegions::StdSegExp StdSeg(BK);
651
653 StdSeg.DetShapeType(), StdSeg);
654
655 loc_mat = StdSeg.GetStdMatrix(matkey);
656
657 // set up array of block matrices.
658 for (size_t i = 0; i < n_exp; ++i)
659 {
660 BlkMatrix->SetBlock(i, i, loc_mat);
661 }
662
663 return BlkMatrix;
664}
665
666// Discrete Fourier Transform for Floquet analysis
668 const string file, Array<OneD, MultiRegions::ExpListSharedPtr> &pFields)
669{
670 size_t ConvectedFields = m_baseflow.size() - 1;
671 size_t npoints = m_baseflow[0].size();
673
674 for (size_t i = 0; i < ConvectedFields; ++i)
675 {
676 m_interp[i] = Array<OneD, NekDouble>(npoints * m_slices, 0.0);
677 }
678
679 // Import the slides into auxiliary vector
680 // The base flow should be stored in the form "filename_%d.ext"
681 // A subdirectory can also be included, such as "dir/filename_%d.ext"
682 size_t found = file.find("%d");
683 ASSERTL0(found != string::npos &&
684 file.find("%d", found + 1) == string::npos,
685 "Since N_slices is specified, the filename provided for function "
686 "'BaseFlow' must include exactly one instance of the format "
687 "specifier '%d', to index the time-slices.");
688 size_t nstart = m_start;
689 for (size_t i = nstart; i < nstart + m_slices * m_skip; i += m_skip)
690 {
691 boost::format filename(file);
692 filename % i;
693 ImportFldBase(filename.str(), pFields, (i - nstart) / m_skip);
694 if (m_session->GetComm()->GetRank() == 0)
695 {
696 cout << "read base flow file " << filename.str() << endl;
697 }
698 }
699 if (!m_isperiodic)
700 {
701 return;
702 }
703
704 // Discrete Fourier Transform of the fields
705 for (size_t k = 0; k < ConvectedFields; ++k)
706 {
707#ifdef NEKTAR_USING_FFTW
708
709 // Discrete Fourier Transform using FFTW
710 Array<OneD, NekDouble> fft_in(npoints * m_slices);
711 Array<OneD, NekDouble> fft_out(npoints * m_slices);
712
715
716 // Shuffle the data
717 for (int j = 0; j < m_slices; ++j)
718 {
719 Vmath::Vcopy(npoints, &m_interp[k][j * npoints], 1, &(fft_in[j]),
720 m_slices);
721 }
722
724 m_slices);
725
726 // FFT Transform
727 for (size_t i = 0; i < npoints; i++)
728 {
729 m_FFT->FFTFwdTrans(m_tmpIN = fft_in + i * m_slices,
730 m_tmpOUT = fft_out + i * m_slices);
731 }
732
733 // Reshuffle data
734 for (int s = 0; s < m_slices; ++s)
735 {
736 Vmath::Vcopy(npoints, &fft_out[s], m_slices,
737 &m_interp[k][s * npoints], 1);
738 }
739
740 Vmath::Zero(fft_in.size(), &fft_in[0], 1);
741 Vmath::Zero(fft_out.size(), &fft_out[0], 1);
742#else
743 // Discrete Fourier Transform using MVM
744 DNekBlkMatSharedPtr blkmat;
746
747 int nrows = blkmat->GetRows();
748 int ncols = blkmat->GetColumns();
749
750 Array<OneD, NekDouble> sortedinarray(ncols);
751 Array<OneD, NekDouble> sortedoutarray(nrows);
752
753 // Shuffle the data
754 for (int j = 0; j < m_slices; ++j)
755 {
756 Vmath::Vcopy(npoints, &m_interp[k][j * npoints], 1,
757 &(sortedinarray[j]), m_slices);
758 }
759
760 // Create NekVectors from the given data arrays
761 NekVector<NekDouble> in(ncols, sortedinarray, eWrapper);
762 NekVector<NekDouble> out(nrows, sortedoutarray, eWrapper);
763
764 // Perform matrix-vector multiply.
765 out = (*blkmat) * in;
766
767 // Reshuffle data
768 for (int s = 0; s < m_slices; ++s)
769 {
770 Vmath::Vcopy(npoints, &sortedoutarray[s], m_slices,
771 &m_interp[k][s * npoints], 1);
772 }
773
774 for (size_t r = 0; r < sortedinarray.size(); ++r)
775 {
776 sortedinarray[0] = 0;
777 sortedoutarray[0] = 0;
778 }
779
780#endif
781 }
782}
783
784} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Describes the specification for a Basis.
Definition: Basis.h:45
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:224
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Defines a specification for a set of points.
Definition: Points.h:50
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
int m_npointsZ
number of points in Z direction (if homogeneous)
NekDouble m_LhomZ
physical length in Z direction (if homogeneous)
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
Initialises the advection object.
enum HomogeneousType m_HomogeneousType
void UpdateBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble time)
NekDouble m_period
period length
void v_Advect(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray) override
Advects a vector field.
int m_npointsX
number of points in X direction (if homogeneous)
int m_npointsY
number of points in Y direction (if homogeneous)
int m_HomoDirec
number of homogenous directions
static SolverUtils::AdvectionSharedPtr create(std::string)
Creates an instance of this class.
bool m_useFFT
flag to determine if use or not the FFT for transformations
void v_SetBaseFlow(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
Overrides the base flow used during linearised advection.
LibUtilities::SessionReaderSharedPtr m_session
Array< OneD, NekDouble > m_tmpIN
void UpdateGradBase(const int var, const MultiRegions::ExpListSharedPtr &field)
DNekBlkMatSharedPtr GetFloquetBlockMatrix(FloquetMatType mattype, bool UseContCoeffs=false) const
bool m_singleMode
flag to determine if use single mode or not
Array< OneD, Array< OneD, NekDouble > > m_baseflow
Storage for base flow.
Array< OneD, NekDouble > m_tmpOUT
Array< OneD, Array< OneD, NekDouble > > m_gradBase
bool m_multipleModes
flag to determine if use multiple mode or not
void DFT(const std::string file, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
LibUtilities::NektarFFTSharedPtr m_FFT
auxiliary variables
bool m_halfMode
flag to determine if use half mode or not
NekDouble m_LhomY
physical length in Y direction (if homogeneous)
void ImportFldBase(std::string pInfile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, int slice)
Import Base flow.
NekDouble m_LhomX
physical length in X direction (if homogeneous)
Array< OneD, Array< OneD, NekDouble > > m_interp
interpolation vector
static std::string className
Name of class.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
An abstract base class encapsulating the concept of advection of a vector field.
Definition: Advection.h:81
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:295
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:367
Class representing a segment element in reference space All interface of this class sits in StdExpans...
Definition: StdSegExp.h:45
NektarFFTFactory & GetNektarFFTFactory()
Definition: NektarFFT.cpp:65
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51
@ eFourierEvenlySpaced
1D Evenly-spaced points using Fourier Fit
Definition: PointsType.h:74
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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 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
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.