Loading...
Searching...
No Matches
IntegralScaleBox.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "IntegralScaleBox.H"
29#include "cartesianCS.H"
30#include "OBJstream.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34template<class Type>
35const Foam::Enum
36<
37 typename Foam::turbulence::IntegralScaleBox<Type>::kernelType
38>
39Foam::turbulence::IntegralScaleBox<Type>::kernelTypeNames
40({
41 { kernelType::GAUSSIAN, "Gaussian" },
42 { kernelType::EXPONENTIAL , "exponential" }
43});
44
45
46template<class Type>
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
52template<class Type>
54Foam::turbulence::IntegralScaleBox<Type>::calcCoordinateSystem
55(
56 const dictionary& dict
57) const
58{
59 return coordinateSystem::NewIfPresent(dict);
60}
61
62
63template<class Type>
64void Foam::turbulence::IntegralScaleBox<Type>::calcCoordinateSystem()
65{
66 // Get patch normal direction into the domain
67 const vector nf((-gAverage(p_.nf())).normalise());
68
69 // Find the second local coordinate direction
70 direction minCmpt = 0;
71 scalar minMag = mag(nf[minCmpt]);
72 for (direction cmpt = 1; cmpt < pTraits<vector>::nComponents; ++cmpt)
73 {
74 const scalar s = mag(nf[cmpt]);
75 if (s < minMag)
76 {
77 minMag = s;
78 minCmpt = cmpt;
79 }
80 }
81
82 // Create the second local coordinate direction
83 vector e2(Zero);
84 e2[minCmpt] = 1;
85
86 // Remove normal component
87 e2 -= (nf&e2)*nf;
88
89 // Create the local coordinate system - default e3-e1 order
90 csysPtr_.reset
91 (
92 new coordSystem::cartesian
93 (
94 Zero, // origin
95 nf^e2, // e3
96 nf // e1
97 )
98 );
99}
100
101
102template<class Type>
104Foam::turbulence::IntegralScaleBox<Type>::calcBoundBox() const
105{
106 // Convert patch points into local coordinate system
107 const pointField localPos
108 (
109 csysPtr_->localPosition
110 (
112 (
113 p_.patch().points(),
114 p_.patch().meshPoints()
115 )
116 )
117 );
118
119 // Calculate bounding-box span and min
120 const bool globalReduction = true;
121 const boundBox bb(localPos, globalReduction);
122
123 return Vector2D<vector>(bb.span(), bb.min());
124}
125
126
127template<class Type>
129Foam::turbulence::IntegralScaleBox<Type>::calcDelta() const
130{
131 return Vector2D<scalar>
132 (
133 boundingBoxSpan_[1]/n_.x(),
134 boundingBoxSpan_[2]/n_.y()
135 );
136}
137
138
139template<class Type>
140Foam::labelList Foam::turbulence::IntegralScaleBox<Type>::calcSpans() const
141{
142 if (!Pstream::master())
143 {
144 return labelList();
145 }
146
147 labelList spans(pTraits<TypeL>::nComponents, label(1));
148 const Vector<label> slice(label(1), n_.x(), n_.y());
149 const TypeL L(convert(L_));
150
151 label j = 0;
152 if (fsm_)
153 {
154 j = pTraits<Type>::nComponents;
155 }
156
157 for (label i = j; i < pTraits<TypeL>::nComponents; ++i)
158 {
159 const label slicei = label(i/pTraits<Type>::nComponents);
160
161 const label n = ceil(L[i]);
162 const label twiceN = 4*n;
163
164 spans[i] = slice[slicei] + twiceN;
165 }
166
167 return spans;
168}
169
170
171template<class Type>
173Foam::turbulence::IntegralScaleBox<Type>::calcKernel() const
174{
175 if (!Pstream::master())
176 {
177 return scalarListList();
178 }
179
180 scalarListList kernel
181 (
182 pTraits<TypeL>::nComponents,
183 scalarList(1, scalar(1))
184 );
185
186 const TypeL L(convert(L_));
187
188 label i = 0;
189 if (fsm_)
190 {
191 i = pTraits<Type>::nComponents;
192 }
193
194 for (direction dir = i; dir < pTraits<TypeL>::nComponents; ++dir)
195 {
196 // The smallest kernel width larger than integral scale
197 // (KSJ:'n' in Eq. 15)
198 const label n = ceil(L[dir]);
199
200 // Full kernel-width (except mid-zero) according to the condition
201 // (KSJ:Eq. 15 whereat N is minimum = 2n)
202 const label twiceN = 4*n;
203
204 // Initialise kernel-coeffs containers with full kernel-width size
205 // Extra elem is mid-zero within [-N, N]
206 kernel[dir] = scalarList(twiceN + 1, Zero);
207
208 // First element: -N within [-N, N]
209 const scalar initElem = -2*scalar(n);
210
211 // Container initialised with [-N, N] (KSJ:p. 658, item-b)
212 std::iota
213 (
214 kernel[dir].begin(),
215 kernel[dir].end(),
216 initElem
217 );
218
219 // Compute kernel coefficients (KSJ:Eq. 14 (Gaussian))
220 scalarList kTemp(kernel[dir]);
221 scalar kSum = 0;
222
223 // Model constant shaping the autocorrelation function (KSJ:Eq. 14)
224 const scalar C = -0.5*constant::mathematical::pi;
225
226 if (kernelType_)
227 {
228 const scalar nSqr = n*n;
229
230 kTemp = sqr(Foam::exp(C*sqr(kTemp)/nSqr));
231 kSum = Foam::sqrt(sum(kTemp));
232
233 kernel[dir] = Foam::exp(C*sqr(kernel[dir])/nSqr)/kSum;
234 }
235 else
236 {
237 kTemp = sqr(Foam::exp(C*mag(kTemp)/n));
238 kSum = Foam::sqrt(sum(kTemp));
239
240 kernel[dir] = Foam::exp(C*mag(kernel[dir])/n)/kSum;
241 }
242 }
243
244 return kernel;
245}
246
247
248template<class Type>
249Foam::scalarListList Foam::turbulence::IntegralScaleBox<Type>::calcBox()
250{
251 if (!Pstream::master())
252 {
253 return scalarListList();
254 }
255
256 scalarListList box(pTraits<Type>::nComponents, scalarList());
257
258 // Initialise: Remaining convenience factors for (e1 e2 e3)
259 for (direction dir = 0; dir < pTraits<Type>::nComponents; ++dir)
260 {
261 scalarList& randomSet = box[dir];
262
263 randomSet = scalarList
264 (
265 spans_[dir]
266 *spans_[dir+pTraits<TypeL>::nComponents/3]
267 *spans_[dir+2*(pTraits<TypeL>::nComponents/3)]
268 );
269
270 if (randomSet.size() > 1e8)
271 {
273 << "Size of random-number set is relatively high:" << nl
274 << " size = " << randomSet.size() << nl
275 << " Please consider to use the forward-stepwise method."
276 << endl;
277 }
278
279 // Initialise: Integral-scale box content with random-number
280 // sets obeying the standard normal distribution
281 std::generate
282 (
283 randomSet.begin(),
284 randomSet.end(),
285 [&]{ return rndGen_.GaussNormal<scalar>(); }
286 );
287 }
288
289 return box;
290}
291
292
293template<class Type>
295Foam::turbulence::IntegralScaleBox<Type>::calcPatchPoints() const
296{
297 if (!Pstream::master())
298 {
299 return pointField();
300 }
301
302 // List of vertex points of the virtual patch in local coordinate system
303 const label nx = n_.x();
304 const label ny = n_.y();
305 const label nPoints = (nx + 1)*(ny + 1);
307
308 label pointi = 0;
309 for (label j = 0; j <= ny; ++j)
310 {
311 for (label i = 0; i <= nx; ++i)
312 {
313 const point p
314 (
315 boundingBoxMin_[0],
316 boundingBoxMin_[1] + i*delta_.x(),
317 boundingBoxMin_[2] + j*delta_.y()
318 );
319 points[pointi] = p;
320 ++pointi;
321 }
322 }
323
324 points = csysPtr_->globalPosition(points);
325
326 return points;
327}
328
329
330template<class Type>
331Foam::faceList Foam::turbulence::IntegralScaleBox<Type>::calcPatchFaces() const
332{
333 if (!Pstream::master())
334 {
335 return faceList();
336 }
337
338 // List of faces of the virtual patch
339 const label nx = n_.x();
340 const label ny = n_.y();
341 const label nFaces = nx*ny;
342 faceList faces(nFaces);
343
344 label m = 0;
345 for (label j = 0; j < ny; ++j)
346 {
347 for (label i = 0; i < nx; ++i)
348 {
349 const label k = j*(nx+1) + i;
350 faces[m] = face({k, k+(nx+1), k+(nx+2), k+1});
351 ++m;
352 }
353 }
354
355 return faces;
356}
357
358
359template<class Type>
360void Foam::turbulence::IntegralScaleBox<Type>::calcPatch()
361{
362 if (debug && Pstream::master())
363 {
364 const auto& tm = p_.patch().boundaryMesh().mesh().time();
365 OBJstream os(tm.path()/"patch.obj");
366 os.write(patchFaces_, patchPoints_, false);
367 }
368
369 if (!patchPtr_)
370 {
371 patchPtr_.reset
372 (
373 new primitivePatch
374 (
375 SubList<face>(patchFaces_),
376 patchPoints_
377 )
378 );
379 }
380}
381
382
383template<class Type>
384typename Foam::turbulence::IntegralScaleBox<Type>::TypeL
385Foam::turbulence::IntegralScaleBox<Type>::convert
386(
387 const typename Foam::turbulence::IntegralScaleBox<Type>::TypeL& L
388) const
389{
390 TypeL Ls(L);
391
392 const scalar deltaT =
393 p_.patch().boundaryMesh().mesh().time().deltaTValue();
394
395 for (direction dir = 0; dir < pTraits<Type>::nComponents; ++dir)
396 {
397 // (KSJ:Eq. 13)
398 // Integral time scales
399 Ls[dir] /= deltaT;
400 // Integral length scales
401 Ls[dir + Ls.size()/3] /= delta_.x();
402 Ls[dir + 2*(Ls.size()/3)] /= delta_.y();
403 }
404
405 return Ls;
406}
407
408
409template<class Type>
410Foam::scalar Foam::turbulence::IntegralScaleBox<Type>::calcC1
411(
412 const vector& L
413) const
414{
415 constexpr scalar c1 = -0.25*constant::mathematical::pi;
416 return Foam::exp(c1/L.x());
417}
418
419
420template<class Type>
421Foam::vector Foam::turbulence::IntegralScaleBox<Type>::calcC1
422(
423 const tensor& L
424) const
425{
426 constexpr scalar c1 = -0.25*constant::mathematical::pi;
427
428 vector C1(Zero);
429 forAll(C1, i)
430 {
431 C1[i] = Foam::exp(c1/L.x()[i]);
432 }
433
434 return C1;
435}
436
437
438template<class Type>
439Foam::scalar Foam::turbulence::IntegralScaleBox<Type>::calcC2
440(
441 const vector& L
442) const
443{
444 constexpr scalar c2 = -0.5*constant::mathematical::pi;
445 return Foam::sqrt(scalar(1) - Foam::exp(c2/L.x()));
446}
447
448
449template<class Type>
450Foam::vector Foam::turbulence::IntegralScaleBox<Type>::calcC2
451(
452 const tensor& L
453) const
454{
455 constexpr scalar c2 = -0.5*constant::mathematical::pi;
456
457 vector C2(Zero);
458 forAll(C2, i)
459 {
460 C2[i] = Foam::sqrt(scalar(1) - Foam::exp(c2/L.x()[i]));
461 }
462
463 return C2;
464}
465
466
467template<class Type>
468void Foam::turbulence::IntegralScaleBox<Type>::updateC1C2()
469{
470 if (p_.patch().boundaryMesh().mesh().time().isAdjustTimeStep())
471 {
472 C1_ = calcC1(convert(L_));
473 C2_ = calcC2(convert(L_));
475}
476
477
478// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
479
480template<class Type>
482(
483 const fvPatch& p
484)
485:
486 p_(p),
487 patchPtr_(nullptr),
488 csysPtr_(nullptr),
489 kernelType_(kernelType::GAUSSIAN),
490 rndGen_(0),
491 n_(Zero),
492 delta_(Zero),
493 boundingBoxSpan_(Zero),
494 boundingBoxMin_(Zero),
495 L_(Zero),
496 spans_(Zero),
497 box_(Zero),
498 kernel_(Zero),
499 patchPoints_(Zero),
500 patchFaces_(Zero),
501 fsm_(false),
502 C1_(Zero),
503 C2_(Zero),
504 slice_(Zero)
505{}
506
507
508template<class Type>
510(
511 const fvPatch& p,
512 const IntegralScaleBox& b
513)
514:
515 p_(p),
516 patchPtr_(nullptr),
517 csysPtr_(b.csysPtr_.clone()),
518 kernelType_(b.kernelType_),
519 rndGen_(b.rndGen_),
520 n_(b.n_),
521 delta_(b.delta_),
522 boundingBoxSpan_(b.boundingBoxSpan_),
523 boundingBoxMin_(b.boundingBoxMin_),
524 L_(b.L_),
525 spans_(b.spans_),
526 box_(b.box_),
527 kernel_(b.kernel_),
528 patchPoints_(b.patchPoints_),
529 patchFaces_(b.patchFaces_),
530 fsm_(b.fsm_),
531 C1_(b.C1_),
532 C2_(b.C2_),
533 slice_(b.slice_)
534{}
535
536
537template<class Type>
539(
540 const fvPatch& p,
541 const dictionary& dict
542)
543:
544 p_(p),
545 patchPtr_(nullptr),
546 csysPtr_(calcCoordinateSystem(dict)),
547 kernelType_
548 (
549 kernelTypeNames.getOrDefault
550 (
551 "kernelType",
552 dict,
553 kernelType::GAUSSIAN
554 )
555 ),
556 rndGen_(time(0)),
557 n_(dict.get<Vector2D<label>>("n")),
558 delta_(Zero),
559 boundingBoxSpan_(Zero),
560 boundingBoxMin_(Zero),
561 L_(dict.get<TypeL>("L")),
562 spans_(Zero),
563 box_(Zero),
564 kernel_(Zero),
565 patchPoints_(Zero),
566 patchFaces_(Zero),
567 fsm_(dict.getOrDefault("fsm", false)),
568 C1_(Zero),
569 C2_(Zero),
570 slice_(Zero)
571{
572 if (cmptMin(L_) < ROOTVSMALL)
573 {
575 << "Integral scale set contains a very small input" << nl
576 << " L = " << L_
577 << exit(FatalIOError);
578 }
579
580 if (min(n_.x(), n_.y()) <= 0)
581 {
582 FatalIOErrorInFunction(dict)
583 << "Number of faces on box inlet plane has non-positive input"
584 << " n = " << n_
585 << exit(FatalIOError);
586 }
587}
588
589
590template<class Type>
592(
593 const IntegralScaleBox& b
594)
595:
596 p_(b.p_),
597 patchPtr_(nullptr),
598 csysPtr_(b.csysPtr_.clone()),
599 kernelType_(b.kernelType_),
600 rndGen_(b.rndGen_),
601 n_(b.n_),
602 delta_(b.delta_),
603 boundingBoxSpan_(b.boundingBoxSpan_),
604 boundingBoxMin_(b.boundingBoxMin_),
605 L_(b.L_),
606 spans_(b.spans_),
607 box_(b.box_),
608 kernel_(b.kernel_),
609 patchPoints_(b.patchPoints_),
610 patchFaces_(b.patchFaces_),
611 fsm_(b.fsm_),
612 C1_(b.C1_),
613 C2_(b.C2_),
614 slice_(b.slice_)
615{}
616
617
618// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
619
620template<class Type>
622{
623 if (!csysPtr_)
624 {
625 calcCoordinateSystem();
626 }
627
628 if (debug && csysPtr_)
629 {
630 Info<< "Local coordinate system:" << nl
631 << " - origin = " << csysPtr_->origin() << nl
632 << " - e1-axis = " << csysPtr_->e1() << nl
633 << " - e2-axis = " << csysPtr_->e2() << nl
634 << " - e3-axis = " << csysPtr_->e3() << nl << endl;
635 }
636
637 {
638 const Vector2D<vector> bb(calcBoundBox());
639
640 boundingBoxSpan_ = bb.x();
641
642 boundingBoxMin_ = bb.y();
643 }
644
645 delta_ = calcDelta();
646
647 spans_ = calcSpans();
648
649 kernel_ = calcKernel();
650
651 box_ = calcBox();
652
653 patchPoints_ = calcPatchPoints();
654
655 patchFaces_ = calcPatchFaces();
656
657 calcPatch();
658
659 if (fsm_)
660 {
661 C1_ = calcC1(convert(L_));
662
663 C2_ = calcC2(convert(L_));
665 slice_ = Field<Type>(p_.size(), Zero);
666 }
667}
668
669
670template<class Type>
672{
673 for (direction dir = 0; dir < pTraits<Type>::nComponents; ++dir)
674 {
675 scalarList& slice = box_[dir];
676
677 // Slice span: span of each inlet-normal slice of integral-scale box
678 // e.g. for U: (Lyu*Lzu, Lyv*Lzv, Lyw*Lzw)
679 const label sliceSpan =
680 spans_[dir + pTraits<TypeL>::nComponents/3]
681 *spans_[dir + 2*(pTraits<TypeL>::nComponents/3)];
682
683 // Shift forward from the back to the front
684 inplaceRotateList(slice, sliceSpan);
685 }
686}
687
688
689template<class Type>
691{
692 for (direction dir = 0; dir < pTraits<Type>::nComponents; ++dir)
693 {
694 scalarList& slice = box_[dir];
695
696 const label sliceSpan =
697 spans_[dir + pTraits<TypeL>::nComponents/3]
698 *spans_[dir + 2*(pTraits<TypeL>::nComponents/3)];
699
700 // Refill the back with a new random-number set
701 for (label i = 0; i < sliceSpan; ++i)
702 {
703 slice[i] = rndGen_.GaussNormal<scalar>();
705 }
706}
707
708
709template<class Type>
712{
713 Field<Type> outFld(n_.x()*n_.y(), Zero);
714
715 for (direction dir = 0; dir < pTraits<Type>::nComponents; ++dir)
716 {
717 const scalarList& in = box_[dir];
718
719 Field<scalar> out(n_.x()*n_.y(), Zero);
720
721 const scalarList& kernel1 = kernel_[dir];
722 const scalarList& kernel2 =
723 kernel_[dir + pTraits<TypeL>::nComponents/3];
724 const scalarList& kernel3 =
725 kernel_[dir + 2*(pTraits<TypeL>::nComponents/3)];
726
727 const label szkernel1 = kernel1.size();
728 const label szkernel2 = kernel2.size();
729 const label szkernel3 = kernel3.size();
730
731 const label sz1 = spans_[dir];
732 const label sz2 = spans_[dir + pTraits<TypeL>::nComponents/3];
733 const label sz3 = spans_[dir + 2*(pTraits<TypeL>::nComponents/3)];
734 const label sz23 = sz2*sz3;
735 const label sz123 = sz1*sz23;
736
737 const label validSlice2 = sz2 - (szkernel2 - 1);
738 const label validSlice3 = sz3 - (szkernel3 - 1);
739
740 // Convolution summation - Along 1st direction
741 scalarField tmp(sz123, Zero);
742 {
743 const label filterCentre = label(szkernel2/label(2));
744 const label endIndex = sz2 - filterCentre;
745 label i0 = 0;
746 label i1 = 0;
747
748 for (label i = 0; i < sz1; ++i)
749 {
750 for (label j = 0; j < sz3; ++j)
751 {
752 i1 += filterCentre;
753
754 for (label k = filterCentre; k < endIndex; ++k)
755 {
756 label q = 0;
757
758 for (label p = szkernel2 - 1; p >= 0; --p, ++q)
759 {
760 tmp[i1] += in[i0 + q]*kernel2[p];
761 }
762 ++i0;
763 ++i1;
764 }
765 i0 += 2*filterCentre;
766 i1 += filterCentre;
767 }
768 }
769 }
770
771 // Convolution summation - Along 2nd direction
772 {
773 const scalarField tmp2(tmp);
774 const label filterCentre = label(szkernel3/label(2));
775 const label endIndex = sz3 - filterCentre;
776 label i1 = 0;
777 label i2 = 0;
778
779 for (label i = 0; i < sz1; ++i)
780 {
781 const label sl = i*sz23;
782
783 for (label j = 0; j < sz2; ++j)
784 {
785 i1 = j + sl;
786 i2 = i1;
787
788 for (label k = 0; k < endIndex - filterCentre; ++k)
789 {
790 tmp[i1] = 0;
791
792 for (label p = szkernel3 - 1, q = 0; p >= 0; --p, ++q)
793 {
794 tmp[i1] += tmp2[i2 + q*sz2]*kernel3[p];
795 }
796 i1 += sz2;
797 i2 += sz2;
798 }
799 i1 += (sz2 + filterCentre);
800 i2 += (sz2 + filterCentre);
801 }
802 }
803 }
804
805 // Convolution summation - Along 3rd direction
806 {
807 label i1 = (szkernel2 - label(1))/label(2);
808 label i2 = (szkernel2 - label(1))/label(2);
809 label i3 = 0;
810
811 for (label i = 0; i < validSlice3; ++i)
812 {
813 for (label j = 0; j < validSlice2; ++j)
814 {
815 scalar sum = 0;
816 i1 = i2 + j;
817
818 for (label k = szkernel1 - 1; k >= 0; --k)
819 {
820 sum += tmp[i1]*kernel1[k];
821 i1 += sz23;
822 }
823 out[i3] = sum;
824 ++i3;
825 }
826 i2 += sz2;
827 }
828 }
829
830 outFld.replace(dir, out);
832
833 return outFld;
834}
835
836
837template<class Type>
839(
841)
842{
843 updateC1C2();
844
845 fld *= C2_;
846 fld += slice_*C1_;
848 // Store current field for the next time-step
849 slice_ = fld;
850}
851
852
853template<class Type>
855(
857)
858{
859 updateC1C2();
860
861 for (direction dir = 0; dir < pTraits<vector>::nComponents; ++dir)
862 {
863 fld.replace
864 (
865 dir,
866 slice_.component(dir)*C1_[dir] + fld.component(dir)*C2_[dir]
867 );
868 }
870 // Store current field for the next time-step
871 slice_ = fld;
872}
873
874
875template<class Type>
877(
878 Ostream& os
879) const
880{
881 os.writeEntryIfDifferent<bool>("fsm", false, fsm_);
882 os.writeEntry("n", n_);
883 os.writeEntry("L", L_);
884 os.writeEntry("kernelType", kernelTypeNames[kernelType_]);
885 if (csysPtr_)
886 {
887 csysPtr_->writeEntry(os);
888 }
889}
890
891
892// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
label k
label n
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition Field.C:619
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Templated 2D Vector derived from VectorSpace adding construction from 2 components,...
Definition Vector2D.H:54
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector2D.H:127
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector2D.H:132
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
Tensor of scalars, i.e. Tensor<scalar>.
A class for managing temporary objects.
Definition tmp.H:75
Class to describe the integral-scale container being used in the turbulentDigitalFilterInletFvPatchFi...
void initialise()
Initialise integral-scale box properties.
void write(Ostream &) const
Write integral-scale box settings.
void refill()
Add a new integral-scale box slice to the rear of the box.
void shift()
Discard current time-step integral-scale box slice (the closest to the patch) by shifting from the ba...
IntegralScaleBox(const fvPatch &p)
Construct from patch.
void correlate(scalarField &fld)
Apply forward-stepwise correlation for scalar fields.
static int debug
Flag to activate debug statements.
Field< Type > convolve() const
Embed two-point correlations, i.e. L.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
label nPoints
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
const char * end
Definition SVGTools.H:223
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Namespace for handling debugging switches.
Definition debug.C:45
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
dimensionedScalar exp(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
void inplaceRotateList(ListType< DataType > &list, label n)
Inplace reversal of a list using the Reversal Block Swapping algorithm.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
const vector L(dict.get< vector >("L"))