Nektar++
FilterLagrangianPoints.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: FilterLagrangianPoints.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: Evaluate physics values at a set of points.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <iomanip>
36using namespace std;
37
42#include <boost/format.hpp>
43
44namespace Nektar::SolverUtils
45{
46
48 GetFilterFactory().RegisterCreatorFunction("LagrangianPoints",
50
52{
53 for (int d = 0; d < m_dim; ++d)
54 {
55 gcoords[d] = m_coords[0][d][id];
56 }
57}
58
60 const Array<OneD, NekDouble> &gcoords)
61{
62 for (int d = 0; d < m_dim; ++d)
63 {
64 m_coords[0][d][id] = gcoords[d];
65 }
66}
67
69{
70 for (int d = 0; d < m_dim; ++d)
71 {
72 data[d] = m_velocity[0][d][id];
73 }
74}
75
77 const Array<OneD, NekDouble> &data)
78{
79 for (int d = 0; d < m_dim; ++d)
80 {
81 m_velocity[0][d][id] = data[d];
82 }
83 for (size_t d = 0; d < m_extraPhysVars.size(); ++d)
84 {
85 m_extraPhysics[d][id] = data[d + m_dim];
86 }
87}
88
90{
91 int n = data.size();
92 if (n <= 1)
93 {
94 return;
95 }
96 Array<OneD, Array<OneD, NekDouble>> res = data[n - 1];
97 for (int i = n - 1; i > 0; --i)
98 {
99 data[i] = data[i - 1];
100 }
101 data[0] = res;
102}
103
105{
106 order = min(order, m_intOrder);
107 if (order == 1)
108 {
109 for (int d = 0; d < m_dim; ++d)
110 {
112 m_coords[0][d], 1);
113 }
115 }
116 else if (order == 2)
117 {
118 for (int d = 0; d < m_dim; ++d)
119 {
121 Vmath::Svtvm(m_totPts, 3., m_velocity[0][d], 1, m_velocity[1][d], 1,
122 res, 1);
123 Vmath::Svtvp(m_totPts, 0.5 * m_dt, res, 1, m_coords[0][d], 1, res,
124 1);
125 }
128 }
129 else
130 {
131 ASSERTL0(false, "Integration order not supported.");
132 }
133}
134
136{
137 if (Np > m_totPts)
138 {
139 for (int i = 0; i < m_coords.size(); ++i)
140 {
141 for (int d = 0; d < m_dim; ++d)
142 {
145 }
146 }
147 }
148 m_totPts = Np;
149 m_localIDToGlobal.clear();
150 m_globalIDToLocal.clear();
151}
152
153static int BinaryWrite(std::ofstream &ofile, std::string str)
154{
155 int tmp = 0;
156 for (size_t i = 0; i < str.size(); ++i)
157 {
158 tmp = str[i];
159 ofile.write((char *)&tmp, 4);
160 }
161 tmp = 0;
162 ofile.write((char *)&tmp, 4);
163 return 4;
164}
165
166static int OutputTec360_binary(const std::string filename,
167 const std::vector<std::string> &variables,
168 const std::vector<int> &rawN,
169 std::vector<Array<OneD, NekDouble>> &data,
170 int isdouble)
171{
172 std::vector<int> N = rawN;
173 for (int i = N.size(); i < 3; ++i)
174 {
175 N.push_back(1);
176 }
177 std::ofstream odata;
178 odata.open(filename, std::ios::binary);
179 if (!odata.is_open())
180 {
181 printf("error unable to open file %s\n", filename.c_str());
182 return -1;
183 }
184 char tecplotversion[] = "#!TDV112";
185 odata.write((char *)tecplotversion, 8);
186 int value1 = 1;
187 odata.write((char *)&value1, 4);
188 int filetype = 0;
189 odata.write((char *)&filetype, 4);
190 // read file title and variable names
191 std::string filetitle = "";
192 BinaryWrite(odata, filetitle);
193 int nvar = variables.size();
194 odata.write((char *)&nvar, 4); // number of variables
195 std::vector<std::string> vartitle;
196 for (int i = 0; i < nvar; ++i)
197 {
198 BinaryWrite(odata, variables[i]);
199 }
200 float zonemarker299I = 299.0f;
201 odata.write((char *)&zonemarker299I, 4);
202 // zone title
203 std::string zonetitle("ZONE 0");
204 BinaryWrite(odata, zonetitle);
205 int parentzone = -1;
206 odata.write((char *)&parentzone, 4);
207 int strandid = -1;
208 odata.write((char *)&strandid, 4);
209 double soltime = 0.0;
210 odata.write((char *)&soltime, 8);
211 int unused = -1;
212 odata.write((char *)&unused, 4);
213 int zonetype = 0;
214 odata.write((char *)&zonetype, 4);
215 int zero = 0;
216 odata.write((char *)&zero, 4);
217 odata.write((char *)&zero, 4);
218 odata.write((char *)&zero, 4);
219 for (int i = 0; i < 3; ++i)
220 {
221 int tmp = N[i];
222 odata.write((char *)&tmp, 4);
223 }
224
225 odata.write((char *)&zero, 4);
226 float eohmarker357 = 357.0f;
227 odata.write((char *)&eohmarker357, 4);
228 float zonemarker299II = 299.0f;
229 odata.write((char *)&zonemarker299II, 4);
230 std::vector<int> binarydatatype(nvar, 1 + (isdouble > 0));
231 odata.write((char *)binarydatatype.data(), 4 * nvar);
232 odata.write((char *)&zero, 4);
233 odata.write((char *)&zero, 4);
234 int minus1 = -1;
235 odata.write((char *)&minus1, 4);
236
237 int datanumber, datasize;
238 datanumber = N[0] * N[1] * N[2];
239 datasize = N[0] * N[1] * N[2] * 4;
240 for (int i = 0; i < nvar; ++i)
241 {
242 double minv = 0., maxv = 1.;
243 for (int p = 0; p < datanumber; ++p)
244 {
245 if (maxv < data[i][p])
246 {
247 maxv = data[i][p];
248 }
249 if (minv > data[i][p])
250 {
251 minv = data[i][p];
252 }
253 }
254 odata.write((char *)&minv, 8);
255 odata.write((char *)&maxv, 8);
256 }
257
258 std::vector<float> vardata(datanumber);
259 for (int i = 0; i < nvar; ++i)
260 {
261 if (isdouble)
262 {
263 odata.write((char *)data[i].data(), datasize * 2);
264 }
265 else
266 {
267 std::vector<float> fdata(datanumber);
268 for (int j = 0; j < datanumber; ++j)
269 {
270 fdata[j] = data[i][j];
271 }
272 odata.write((char *)fdata.data(), datasize);
273 }
274 }
275 odata.close();
276 return 0;
277}
278
280 std::string filename, bool verbose,
281 std::map<std::string, NekDouble> &params)
282{
283 std::vector<std::string> COORDS = {"x", "y", "z"};
284 std::vector<std::string> VELOCI = {"u", "v", "w"};
285 std::vector<std::string> variables;
286 std::vector<Array<OneD, NekDouble>> data;
287 for (int d = 0; d < m_dim; ++d)
288 {
289 variables.push_back(COORDS[d]);
290 }
291 for (int d = 0; d < m_dim; ++d)
292 {
293 variables.push_back(VELOCI[d]);
294 }
295 if (params.size() > 0)
296 {
297 for (int d = 0; d < m_dim; ++d)
298 {
300 Vmath::Sadd(m_totPts, params[COORDS[d]], m_coords[0][d], 1, x, 1);
301 data.push_back(x);
302 }
303 for (int d = 0; d < m_dim; ++d)
304 {
306 Vmath::Sadd(m_totPts, params[VELOCI[d]], m_velocity[0][d], 1, x, 1);
307 data.push_back(x);
308 }
309 }
310 else
311 {
312 for (int d = 0; d < m_dim; ++d)
313 {
314 data.push_back(m_coords[0][d]);
315 }
316 for (int d = 0; d < m_dim; ++d)
317 {
318 data.push_back(m_velocity[0][d]);
319 }
320 }
321 for (size_t d = 0; d < m_extraPhysics.size(); ++d)
322 {
323 data.push_back(m_extraPhysics[d]);
324 variables.push_back(m_extraPhysVars[d]);
325 }
326 OutputTec360_binary(filename, variables, m_N, data, 1);
327 if (verbose)
328 {
329 int Np = m_coords[0][0].size();
330 for (int d = 0; d < m_dim; ++d)
331 {
332 NekDouble value = 0.;
333 for (int i = 0; i < Np; ++i)
334 {
335 value += m_coords[0][d][i] * m_coords[0][d][i];
336 }
337 value = sqrt(value / Np);
338 cout << "L 2 error (variable L" << COORDS[d] << ") : " << value
339 << endl;
340 }
341 for (int d = 0; d < m_dim; ++d)
342 {
343 NekDouble value = 0.;
344 for (int i = 0; i < Np; ++i)
345 {
346 value += m_velocity[0][d][i] * m_velocity[0][d][i];
347 }
348 value = sqrt(value / Np);
349 cout << "L 2 error (variable L" << VELOCI[d] << ") : " << value
350 << endl;
351 }
352 }
353 cout << "Write file " << filename << endl;
354}
355
357 const Array<OneD, NekDouble> &gcoords)
358{
359 StationaryPoints::v_AssignPoint(id, pid, gcoords);
360 for (int d = 0; d < m_dim; ++d)
361 {
362 m_coords[0][d][id] = gcoords[d];
363 }
364}
365
366StatLagrangianPoints::StatLagrangianPoints([[maybe_unused]] int rank, int dim,
367 int intOrder, int idOffset,
368 NekDouble dt,
369 const std::vector<int> &Np,
370 const std::vector<NekDouble> &Box,
371 std::vector<std::string> extraVars)
372{
373 m_idOffset = idOffset;
374 m_dim = dim;
375 if (intOrder < 0)
376 {
377 m_intOrder = -intOrder;
378 m_dt = -dt;
379 }
380 else
381 {
382 m_intOrder = intOrder;
383 m_dt = dt;
384 }
385 if (Np.size() < m_dim || Box.size() < m_dim * 2)
386 {
387 throw ErrorUtil::NekError("Not enough input for StationaryPoints.");
388 }
389 Array<OneD, NekDouble> delta(3, 1.);
390 m_N.resize(3, 1);
391 m_totPts = 1;
392 for (int d = 0; d < m_dim; ++d)
393 {
394 m_N[d] = Np[d];
395 m_totPts *= Np[d];
396 delta[d] = Np[d] < 2 ? 1. : (Box[d * 2 + 1] - Box[d * 2]) / (Np[d] - 1);
397 }
400 for (int i = 0; i < m_intOrder; ++i)
401 {
404 for (int d = 0; d < m_dim; ++d)
405 {
408 }
409 }
411 for (size_t i = 0; i < extraVars.size(); ++i)
412 {
413 m_extraPhysVars.push_back(extraVars[i]);
415 }
416 // initialise points
417 int count = 0;
418 for (int k = 0; k < m_N[2]; ++k)
419 {
420 for (int j = 0; j < m_N[1]; ++j)
421 {
422 for (int i = 0; i < m_N[0]; ++i)
423 {
424 if (m_dim > 0)
425 {
426 m_coords[0][0][count] = Box[0] + delta[0] * i;
427 }
428 if (m_dim > 1)
429 {
430 m_coords[0][1][count] = Box[2] + delta[1] * j;
431 }
432 if (m_dim > 2)
433 {
434 m_coords[0][2][count] = Box[4] + delta[2] * k;
435 }
436 m_localIDToGlobal[count] = m_idOffset + count;
437 m_globalIDToLocal[m_idOffset + count] = count;
438 ++count;
439 }
440 }
441 }
442}
443
444/**
445 *
446 */
449 const std::weak_ptr<EquationSystem> &pEquation,
450 const std::map<std::string, std::string> &pParams)
451 : Filter(pSession, pEquation), v_params(pParams)
452{
453}
454
455/**
456 *
457 */
459{
461}
462
463/**
464 *
465 */
468 const NekDouble &time)
469{
470 auto it = v_params.find("DefaultValues");
471 if (it != v_params.end())
472 {
474 }
476 NekDouble dt;
477 m_session->LoadParameter("TimeStep", dt);
478
479 it = v_params.find("FrameVelocity");
480 if (it != v_params.end())
481 {
482 std::vector<std::string> strs;
483 ParseUtils::GenerateVector(it->second, strs);
484 for (int i = strs.size(); i < m_spacedim; ++i)
485 {
486 strs.push_back("0");
487 }
488 for (int i = 0; i < m_spacedim; ++i)
489 {
492 pFields[0]->GetSession()->GetInterpreter(), strs[i]);
493 }
494 }
495
496 it = v_params.find("FrameDisplacement");
497 if (it != v_params.end())
498 {
499 std::vector<std::string> strs;
500 ParseUtils::GenerateVector(it->second, strs);
501 for (int i = strs.size(); i < m_spacedim; ++i)
502 {
503 strs.push_back("0");
504 }
505 for (int i = 0; i < m_spacedim; ++i)
506 {
509 pFields[0]->GetSession()->GetInterpreter(), strs[i]);
510 }
511 }
512
513 int intOrder = 1;
514 it = v_params.find("IntOrder");
515 if (it != v_params.end())
516 {
517 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
518 intOrder = round(equ.Evaluate());
519 }
520 if (intOrder == 0)
521 {
522 intOrder = 1;
523 m_isMovablePoints = false;
524 }
525 else
526 {
527 m_isMovablePoints = true;
528 }
529
530 it = v_params.find("RootOutputL2Norm");
531 if (it != v_params.end())
532 {
533 std::string value = it->second;
534 m_outputL2Norm = m_rank == 0 && (value == "1" || value == "Yes");
535 }
536 else
537 {
538 m_outputL2Norm = false;
539 }
540
541 it = v_params.find("OutputFile");
542 std::string outputFilename;
543 if (it == v_params.end())
544 {
545 outputFilename = m_session->GetSessionName();
546 }
547 else
548 {
549 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
550 outputFilename = it->second;
551 }
552 m_outputFile = outputFilename + "_R%05d_T%010.6lf.plt";
553
554 // OutputFrequency
555 it = v_params.find("OutputFrequency");
556 if (it == v_params.end())
557 {
559 }
560 else
561 {
562 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
563 m_outputFrequency = round(equ.Evaluate());
564 }
565
566 it = v_params.find("OutputSampleFrequency");
567 if (it == v_params.end())
568 {
570 }
571 else
572 {
573 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
574 m_outputSampleFrequency = round(equ.Evaluate());
575 }
576
577 it = v_params.find("Box_" + std::to_string(m_rank));
578 if (it != v_params.end())
579 {
580 std::vector<NekDouble> values;
581 ParseUtils::GenerateVector(it->second, values);
582 if (values.size() < m_spacedim * 3 + 1)
583 {
585 "Box formate error "
586 "(offset,Nx,Ny,Nx,xmin,xmax,ymin,ymax,zmin,zmax): " +
587 it->second);
588 }
589 m_frame.Update(time);
590 std::vector<int> Np(m_spacedim, 1);
591 m_box.resize(2 * m_spacedim, 0.);
592 int idOffset = std::round(values[0]);
593 for (int d = 0; d < m_spacedim; ++d)
594 {
595 Np[d] = std::round(values[d + 1]);
596 }
597 for (int d = 0; d < 2 * m_spacedim; ++d)
598 {
599 m_box[d] = values[1 + d + m_spacedim] - m_frame.m_frameDisp[d / 2];
600 }
601 std::vector<std::string> extraVars;
602 ExtraPhysicsVars(extraVars);
604 m_rank, m_spacedim, intOrder, idOffset, dt, Np, m_box, extraVars);
606
607 it = v_params.find("OutputSamplePoints");
608 std::vector<unsigned int> samplepts;
609 if (it != v_params.end())
610 {
611 ParseUtils::GenerateSeqVector(it->second, samplepts);
612 }
613 std::set<int> sampleIds(samplepts.begin(), samplepts.end());
614 it = v_params.find("OutputSamplePointsCondition");
615 LibUtilities::EquationSharedPtr sampleCondition;
616 if (it != v_params.end())
617 {
618 sampleCondition =
620 pFields[0]->GetSession()->GetInterpreter(), it->second);
621 }
622
623 for (int p = 0; p < m_staticPts->GetTotPoints(); ++p)
624 {
626 m_staticPts->GetCoords(p, data);
627 int globalId = m_staticPts->LocalToGlobal(p);
628 NekDouble x = data[0] + m_frame.m_frameDisp[0],
629 y = data[1] + m_frame.m_frameDisp[1];
630 NekDouble z =
631 m_spacedim > 2 ? data[2] : 0. + m_frame.m_frameDisp[2];
632 bool issample =
633 (sampleCondition == nullptr)
634 ? false
635 : (sampleCondition->Evaluate(x, y, z, time) > 0);
636 if (sampleIds.find(globalId) != sampleIds.end() || issample)
637 {
638 m_samplePointIDs.insert(globalId);
639 }
640 }
641 }
642 if (!m_samplePointIDs.empty())
643 {
644 std::string sampleFilename = outputFilename + "_Sample_R%05d.dat";
645 boost::format filename(sampleFilename);
646 filename % m_rank;
647 m_ofstreamSamplePoints.open(filename.str(), std::ofstream::out);
648 std::vector<std::string> COORDS = {"x", "y", "z"};
649 std::vector<std::string> VELOCI = {"u", "v", "w"};
650 m_ofstreamSamplePoints << "variables = pid, time ";
651 for (int d = 0; d < m_spacedim; ++d)
652 {
653 m_ofstreamSamplePoints << COORDS[d] << " ";
654 }
655 for (int d = 0; d < m_spacedim; ++d)
656 {
657 m_ofstreamSamplePoints << VELOCI[d] << " ";
658 }
660 }
662 v_Update(pFields, time);
663}
664
666 [[maybe_unused]] Array<OneD, NekDouble> gcoords,
667 [[maybe_unused]] NekDouble time, Array<OneD, NekDouble> vel)
668{
669 if (m_frame.m_frameVelFunction.size() > 0)
670 {
671 for (int i = 0; i < m_spacedim; ++i)
672 {
673 vel[i] -= m_frame.m_frameVel[i];
674 }
675 }
676}
677
679 std::vector<std::string> &extraVars)
680{
681 if (m_spacedim == 2)
682 {
683 extraVars.push_back("u_x");
684 extraVars.push_back("u_y");
685 extraVars.push_back("v_x");
686 extraVars.push_back("v_y");
687 extraVars.push_back("LW_z");
688 }
689 else if (m_spacedim == 3)
690 {
691 extraVars.push_back("u_x");
692 extraVars.push_back("u_y");
693 extraVars.push_back("u_z");
694 extraVars.push_back("v_x");
695 extraVars.push_back("v_y");
696 extraVars.push_back("v_z");
697 extraVars.push_back("w_x");
698 extraVars.push_back("w_y");
699 extraVars.push_back("w_z");
700 extraVars.push_back("LW_x");
701 extraVars.push_back("LW_y");
702 extraVars.push_back("LW_z");
703 }
704}
705
708 std::vector<Array<OneD, NekDouble>> &PhysicsData)
709{
710 int npts = pFields[0]->GetTotPoints();
711 for (int i = 0; i < m_spacedim; ++i)
712 {
713 PhysicsData.push_back(pFields[i]->UpdatePhys());
714 }
715 int offset = PhysicsData.size();
716 for (int d = 0; d < m_spacedim; ++d)
717 {
719 for (int i = 0; i < m_spacedim; ++i)
720 {
721 data[i] = Array<OneD, NekDouble>(npts, 0.);
722 }
723 if (m_spacedim == 2)
724 {
725 pFields[d]->PhysDeriv(pFields[d]->UpdatePhys(), data[0], data[1]);
726 }
727 else if (m_spacedim == 3)
728 {
729 pFields[d]->PhysDeriv(pFields[d]->UpdatePhys(), data[0], data[1],
730 data[2]);
731 }
732 for (int i = 0; i < m_spacedim; ++i)
733 {
734 PhysicsData.push_back(data[i]);
735 }
736 }
737 int vortexdim = m_spacedim == 3 ? 3 : 1;
738 Array<OneD, Array<OneD, NekDouble>> vorticity(vortexdim);
739 Array<OneD, Array<OneD, NekDouble>> Ldata(vortexdim);
741 Array<OneD, NekDouble> temps(npts, 0.);
742 for (int d = 0; d < vortexdim; ++d)
743 {
744 Ldata[d] = Array<OneD, NekDouble>(npts, 0.);
745 vorticity[d] = Array<OneD, NekDouble>(npts, 0.);
746 }
747 for (int d = 0; d < m_spacedim; ++d)
748 {
749 temp[d] = Array<OneD, NekDouble>(npts, 0.);
750 }
751 if (m_spacedim == 2)
752 {
753 // 0 ux, 1 uy, 2 vx, 3 vy
754 Vmath::Vsub(npts, PhysicsData[offset + 2], 1, PhysicsData[offset + 1],
755 1, vorticity[0], 1); // vx-uy
756 }
757 else if (m_spacedim == 3)
758 {
759 // 0 ux, 1 uy, 2 uz, 3 vx, 4 vy, 5 vz, 6 wx, 7 wy, 8 wz
760 Vmath::Vsub(npts, PhysicsData[offset + 7], 1, PhysicsData[offset + 5],
761 1, vorticity[0], 1); // wy-vz
762 Vmath::Vsub(npts, PhysicsData[offset + 2], 1, PhysicsData[offset + 6],
763 1, vorticity[1], 1); // uz-wx
764 Vmath::Vsub(npts, PhysicsData[offset + 3], 1, PhysicsData[offset + 1],
765 1, vorticity[2], 1); // vx-uy
766 }
767 for (int d = 0; d < vortexdim; ++d)
768 {
769 if (m_spacedim == 2)
770 {
771 pFields[0]->PhysDeriv(vorticity[d], temp[0], temp[1]);
772 }
773 else if (m_spacedim == 3)
774 {
775 pFields[0]->PhysDeriv(vorticity[d], temp[0], temp[1], temp[2]);
776 }
777 for (int j = 0; j < m_spacedim; ++j)
778 {
779 pFields[0]->PhysDeriv(j, temp[j], temps);
780 Vmath::Vadd(npts, temps, 1, Ldata[d], 1, Ldata[d], 1);
781 }
782 PhysicsData.push_back(Ldata[d]);
783 }
784}
785
788 const NekDouble &time)
789{
790 if (m_isMovablePoints || m_index == 0)
791 {
792 PartitionMobilePts(pFields);
793 }
794 std::vector<Array<OneD, NekDouble>> PhysicsData;
795 GetPhysicsData(pFields, PhysicsData);
796 m_frame.Update(time);
797 EvaluateMobilePhys(pFields[0], PhysicsData, time);
798 std::map<int, std::set<int>> callbackUpdateMobCoords;
799 PassMobilePhysToStatic(callbackUpdateMobCoords);
800 // output Lagrangian coordinates
801 if (m_index % m_outputFrequency == 0 && m_staticPts != nullptr)
802 {
803 OutputStatPoints(time);
804 }
805 // output sample points
806 if (m_index % m_outputSampleFrequency == 0 && !m_samplePointIDs.empty())
807 {
808 OutputSamplePoints(time);
809 }
810 // update Lagrangian coordinates
812 {
813 if (m_staticPts != nullptr)
814 {
815 m_staticPts->TimeAdvance(m_index + 1);
816 }
817 // return the updated global coordinates to mobile points
818 PassStaticCoordsToMobile(callbackUpdateMobCoords);
819 }
820 ++m_index;
821}
822
824{
825 m_ofstreamSamplePoints << "ZONE T=\"" << time
826 << "\", I=" << m_samplePointIDs.size() << endl;
827 for (auto &pid : m_samplePointIDs)
828 {
829 m_ofstreamSamplePoints << pid << " " << time << " ";
831 m_staticPts->GetCoordsByPID(pid, data);
832 for (int d = 0; d < m_spacedim; ++d)
833 {
834 m_ofstreamSamplePoints << data[d] + m_frame.m_frameDisp[d] << " ";
835 }
836 m_staticPts->GetPhysicsByPID(pid, data);
837 for (int d = 0; d < m_spacedim; ++d)
838 {
839 m_ofstreamSamplePoints << data[d] + m_frame.m_frameVel[d] << " ";
840 }
842 }
843}
844
846{
847 boost::format filename(m_outputFile);
848 filename % m_rank % time;
849 std::map<std::string, NekDouble> params;
850 if (m_frame.m_frameVelFunction.size())
851 {
852 std::vector<std::string> COORDS = {"x", "y", "z"};
853 std::vector<std::string> VELOCI = {"u", "v", "w"};
854 for (auto it : m_frame.m_frameVelFunction)
855 {
856 params[VELOCI[it.first]] = m_frame.m_frameVel[it.first];
857 }
858 for (auto it : m_frame.m_frameDispFunction)
859 {
860 params[COORDS[it.first]] = m_frame.m_frameDisp[it.first];
861 }
862 }
863 m_staticPts->OutputData(filename.str(), m_outputL2Norm, params);
864}
865
868 &pFields,
869 [[maybe_unused]] const NekDouble &time)
870{
871}
872
874{
875 return true;
876}
877
878} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
Nektar::ErrorUtil::NekError NekError
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector.
Definition: ParseUtils.cpp:104
SOLVER_UTILS_EXPORT void EvaluateMobilePhys(const MultiRegions::ExpListSharedPtr &pField, std::vector< Array< OneD, NekDouble > > &PhysicsData, NekDouble time)
SOLVER_UTILS_EXPORT void Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time, std::vector< std::string > &defaultValues)
SOLVER_UTILS_EXPORT void PassMobilePhysToStatic(std::map< int, std::set< int > > &callbackUpdateMobCoords)
SOLVER_UTILS_EXPORT void PartitionMobilePts(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
SOLVER_UTILS_EXPORT void SetUpCommInfo()
SOLVER_UTILS_EXPORT void CopyStaticPtsToMobile()
StationaryPointsSharedPtr m_staticPts
SOLVER_UTILS_EXPORT void PassStaticCoordsToMobile(std::map< int, std::set< int > > &callbackUpdateMobCoords)
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
void ExtraPhysicsVars(std::vector< std::string > &extraVars)
SOLVER_UTILS_EXPORT FilterLagrangianPoints(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
SOLVER_UTILS_EXPORT ~FilterLagrangianPoints() override
struct Nektar::SolverUtils::FilterLagrangianPoints::MovingFrame m_frame
std::map< std::string, std::string > v_params
void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
void v_ModifyVelocity(Array< OneD, NekDouble > gcoords, NekDouble time, Array< OneD, NekDouble > vel) override
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
void GetPhysicsData(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &PhysicsData)
StatLagrangianPoints(int rank, int dim, int intOrder, int idOffset, NekDouble dt, const std::vector< int > &Np, const std::vector< NekDouble > &Box, std::vector< std::string > extraVars)
void v_GetCoords(int pid, Array< OneD, NekDouble > &gcoords) override
void v_SetPhysics(int pid, const Array< OneD, NekDouble > &data) override
void v_OutputData(std::string filename, bool verbose, std::map< std::string, NekDouble > &params) override
void v_SetCoords(int pid, const Array< OneD, NekDouble > &gcoords) override
void v_GetPhysics(int pid, Array< OneD, NekDouble > &data) override
Array< OneD, Array< OneD, NekDouble > > m_extraPhysics
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_velocity
void v_AssignPoint(int id, int pid, const Array< OneD, NekDouble > &gcoords) override
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_coords
virtual void v_AssignPoint(int id, int pid, const Array< OneD, NekDouble > &gcoords)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
static void RollOver(Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &data)
static int OutputTec360_binary(const std::string filename, const std::vector< std::string > &variables, const std::vector< int > &rawN, std::vector< Array< OneD, NekDouble > > &data, int isdouble)
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:39
static int BinaryWrite(std::ofstream &ofile, std::string str)
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
void Svtvm(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvm (scalar times vector minus vector): z = alpha*x - y.
Definition: Vmath.hpp:424
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 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 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 > sqrt(scalarT< T > in)
Definition: scalar.hpp:294
std::map< int, LibUtilities::EquationSharedPtr > m_frameVelFunction
std::map< int, LibUtilities::EquationSharedPtr > m_frameDispFunction