Loading...
Searching...
No Matches
surfaceFieldValue.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2017-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "surfaceFieldValue.H"
30#include "fvMesh.H"
31#include "emptyPolyPatch.H"
32#include "coupledPolyPatch.H"
33#include "sampledSurface.H"
35#include "PatchTools.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40// Max number of warnings
41static constexpr const unsigned maxWarnings = 10u;
43namespace Foam
45namespace functionObjects
46{
47namespace fieldValues
48{
49 defineTypeNameAndDebug(surfaceFieldValue, 0);
52}
53}
54}
55
56
57const Foam::Enum
58<
59 Foam::functionObjects::fieldValues::surfaceFieldValue::regionTypes
60>
61Foam::functionObjects::fieldValues::surfaceFieldValue::regionTypeNames_
62({
63 { regionTypes::stFaceZone, "faceZone" },
64 { regionTypes::stPatch, "patch" },
65 { regionTypes::stObject, "functionObjectSurface" },
66 { regionTypes::stSampled, "sampledSurface" },
67});
68
69
70const Foam::Enum
71<
72 Foam::functionObjects::fieldValues::surfaceFieldValue::operationType
73>
74Foam::functionObjects::fieldValues::surfaceFieldValue::operationTypeNames_
75({
76 // Normal operations
77 { operationType::opNone, "none" },
78 { operationType::opMin, "min" },
79 { operationType::opMax, "max" },
80 { operationType::opSum, "sum" },
81 { operationType::opSumMag, "sumMag" },
82 { operationType::opSumDirection, "sumDirection" },
83 { operationType::opSumDirectionBalance, "sumDirectionBalance" },
84 { operationType::opAverage, "average" },
85 { operationType::opAreaAverage, "areaAverage" },
86 { operationType::opAreaIntegrate, "areaIntegrate" },
87 { operationType::opCoV, "CoV" },
88 { operationType::opAreaNormalAverage, "areaNormalAverage" },
89 { operationType::opAreaNormalIntegrate, "areaNormalIntegrate" },
90 { operationType::opUniformity, "uniformity" },
91
92 // Using weighting
93 { operationType::opWeightedSum, "weightedSum" },
94 { operationType::opWeightedAverage, "weightedAverage" },
95 { operationType::opWeightedAreaAverage, "weightedAreaAverage" },
96 { operationType::opWeightedAreaIntegrate, "weightedAreaIntegrate" },
97 { operationType::opWeightedUniformity, "weightedUniformity" },
98
99 // Using absolute weighting
100 { operationType::opAbsWeightedSum, "absWeightedSum" },
101 { operationType::opAbsWeightedAverage, "absWeightedAverage" },
102 { operationType::opAbsWeightedAreaAverage, "absWeightedAreaAverage" },
103 { operationType::opAbsWeightedAreaIntegrate, "absWeightedAreaIntegrate" },
104 { operationType::opAbsWeightedUniformity, "absWeightedUniformity" },
105});
106
107const Foam::Enum
108<
109 Foam::functionObjects::fieldValues::surfaceFieldValue::postOperationType
110>
111Foam::functionObjects::fieldValues::surfaceFieldValue::postOperationTypeNames_
112({
113 { postOperationType::postOpNone, "none" },
114 { postOperationType::postOpMag, "mag" },
115 { postOperationType::postOpSqrt, "sqrt" },
116});
117
118
119// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
120
122Foam::functionObjects::fieldValues::surfaceFieldValue::obr() const
123{
124 if (stObject == regionType_)
125 {
126 return storedObjects().lookupObject<polySurface>(regionName_);
127 }
128
129 return mesh_;
130}
131
132
133void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
134{
135 // Indices for all matches, already sorted
136 const labelList zoneIds
137 (
138 mesh_.faceZones().indices(selectionNames_)
139 );
140
141 // Total number of faces that could be selected (before patch filtering)
142 label numFaces = 0;
143 for (const label zoneId : zoneIds)
144 {
145 numFaces += mesh_.faceZones()[zoneId].size();
146 }
147
148 faceId_.resize_nocopy(numFaces);
149 facePatchId_.resize_nocopy(numFaces);
150 faceFlip_.resize_nocopy(numFaces);
151
152 numFaces = 0;
153
154 for (const label zoneId : zoneIds)
155 {
156 const faceZone& fZone = mesh_.faceZones()[zoneId];
157
158 forAll(fZone, i)
159 {
160 const label meshFacei = fZone[i];
161 const bool isFlip = fZone.flipMap()[i];
162
163 // Internal faces
164 label faceId = meshFacei;
165 label facePatchId = -1;
166
167 // Boundary faces
168 if (!mesh_.isInternalFace(meshFacei))
169 {
170 facePatchId = mesh_.boundaryMesh().whichPatch(meshFacei);
171 const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
172
173 if (isA<emptyPolyPatch>(pp))
174 {
175 continue; // Ignore empty patch
176 }
177
178 const auto* cpp = isA<coupledPolyPatch>(pp);
179
180 if (cpp && !cpp->owner())
181 {
182 continue; // Ignore neighbour side
183 }
184
185 faceId = pp.whichFace(meshFacei);
186 }
187
188 if (faceId >= 0)
189 {
190 faceId_[numFaces] = faceId;
191 facePatchId_[numFaces] = facePatchId;
192 faceFlip_[numFaces] = isFlip;
193
194 ++numFaces;
195 }
196 }
197 }
198
199 // Shrink to size used
200 faceId_.resize(numFaces);
201 facePatchId_.resize(numFaces);
202 faceFlip_.resize(numFaces);
203 nFaces_ = returnReduce(numFaces, sumOp<label>());
204
205 if (!nFaces_)
206 {
207 // Raise warning or error
208 refPtr<OSstream> os;
209 bool fatal = false;
210
211 ++nWarnings_; // Always increment (even if ignore etc)
212
213 switch (emptySurfaceError_)
214 {
215 case error::handlerTypes::IGNORE:
216 {
217 break;
218 }
219
220 case error::handlerTypes::WARN:
221 {
222 if (nWarnings_ <= maxWarnings)
223 {
225 }
226 break;
227 }
228
229 // STRICT / DEFAULT
230 default:
231 {
233 fatal = true;
234 break;
235 }
236 }
237
238 if (os)
239 {
240 os.ref()
241 << type() << ' ' << name() << ": "
242 << regionTypeNames_[regionType_]
243 << '(' << regionName_ << "):" << nl;
244
245 if (zoneIds.empty())
246 {
247 os.ref()
248 << " No matching face zones: "
249 << flatOutput(selectionNames_) << nl
250 << " Known face zones: "
251 << flatOutput(mesh_.faceZones().names()) << nl;
252 }
253 else
254 {
255 os.ref()
256 << " The face zones: "
257 << flatOutput(selectionNames_) << nl
258 << " resulted in 0 faces" << nl;
259 }
260
261 if (fatal)
262 {
263 FatalError<< exit(FatalError);
264 }
265 else if (nWarnings_ == maxWarnings)
266 {
267 os.ref()
268 << "... suppressing further warnings." << nl;
269 }
270 }
271 }
272}
273
274
275void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchFaces()
276{
277 // Patch indices for all matches
279
280 // Total number of faces selected
281 label numFaces = 0;
282
283 labelList selected
284 (
285 mesh_.boundaryMesh().indices(selectionNames_, true) // useGroup
286 );
287
288 DynamicList<label> bad;
289 for (const label patchi : selected)
290 {
291 const polyPatch& pp = mesh_.boundaryMesh()[patchi];
292
293 if (isA<emptyPolyPatch>(pp))
294 {
295 bad.append(patchi);
296 }
297 else
298 {
299 numFaces += pp.size();
300 }
301 }
302
303 if (bad.size())
304 {
305 label nGood = (selected.size() - bad.size());
306
307 auto& os = (nGood > 0 ? WarningInFunction : FatalErrorInFunction);
308
309 os << "Cannot sample an empty patch" << nl;
310
311 for (const label patchi : bad)
312 {
313 os << " "
314 << mesh_.boundaryMesh()[patchi].name() << nl;
315 }
316
317 if (nGood)
318 {
319 os << "No non-empty patches selected" << endl
320 << exit(FatalError);
321 }
322 else
323 {
324 os << "Selected " << nGood << " non-empty patches" << nl;
325 }
326
327 patchIds.resize(nGood);
328 nGood = 0;
329 for (const label patchi : selected)
330 {
331 if (!bad.contains(patchi))
332 {
333 patchIds[nGood] = patchi;
334 ++nGood;
335 }
336 }
337 }
338 else
339 {
340 patchIds = std::move(selected);
341 }
342
343 faceId_.resize_nocopy(numFaces);
344 facePatchId_.resize_nocopy(numFaces);
345 faceFlip_.resize_nocopy(numFaces);
346 nFaces_ = returnReduce(numFaces, sumOp<label>());
347
348 numFaces = 0;
349 for (const label patchi : patchIds)
350 {
351 const polyPatch& pp = mesh_.boundaryMesh()[patchi];
352 const label len = pp.size();
353
354 SubList<label>(faceId_, len, numFaces) = identity(len);
355 SubList<label>(facePatchId_, len, numFaces) = patchi;
356 SubList<bool>(faceFlip_, len, numFaces) = false;
357
358 numFaces += len;
359 }
360
361 if (!nFaces_)
362 {
363 // Raise warning or error
364 refPtr<OSstream> os;
365 bool fatal = false;
366
367 ++nWarnings_; // Always increment (even if ignore etc)
368
369 switch (emptySurfaceError_)
370 {
371 case error::handlerTypes::IGNORE:
372 {
373 break;
374 }
375
376 case error::handlerTypes::WARN:
377 {
378 if (nWarnings_ <= maxWarnings)
379 {
381 }
382 break;
383 }
384
385 // STRICT / DEFAULT
386 default:
387 {
389 fatal = true;
390 break;
391 }
392 }
393
394 if (os)
395 {
396 os.ref()
397 << type() << ' ' << name() << ": "
398 << regionTypeNames_[regionType_]
399 << '(' << regionName_ << "):" << nl;
400
401 if (patchIds.empty())
402 {
403 os.ref()
404 << " No matching patches: "
405 << flatOutput(selectionNames_) << nl
406 << " Known patch names:" << nl
407 << mesh_.boundaryMesh().names() << nl;
408 }
409 else
410 {
411 os.ref()
412 << " The patches: "
413 << flatOutput(selectionNames_) << nl
414 << " resulted in 0 faces" << nl;
415 }
416
417 if (fatal)
418 {
419 FatalError<< exit(FatalError);
420 }
421 else if (nWarnings_ == maxWarnings)
422 {
423 os.ref()
424 << "... suppressing further warnings." << nl;
425 }
426 }
427 }
428}
429
430
431void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
432(
433 faceList& faces,
434 pointField& points
435) const
436{
437 labelList whichFaces(faceId_);
438
439 // Remap patch-face ids to mesh face ids
440 forAll(whichFaces, i)
441 {
442 const label patchi = facePatchId_[i];
443 if (patchi != -1)
444 {
445 whichFaces[i] += mesh_.boundaryMesh()[patchi].start();
446 }
447 }
448
450 (
451 IndirectList<face>(mesh_.faces(), std::move(whichFaces)),
452 mesh_.points()
453 );
454
455
456 if (Pstream::parRun())
457 {
458 // Topological merge
459 labelList pointToGlobal;
460 labelList uniqueMeshPointLabels;
461 autoPtr<globalIndex> globalPoints;
462 autoPtr<globalIndex> globalFaces;
464 (
465 mesh_,
466 pp.localFaces(),
467 pp.meshPoints(),
468 pp.meshPointMap(),
469
470 pointToGlobal,
471 uniqueMeshPointLabels,
472 globalPoints,
473 globalFaces,
474
475 faces,
476 points
477 );
478 }
479 else
480 {
481 faces = pp.localFaces();
482 points = pp.localPoints();
483 }
484}
485
486
487void Foam::functionObjects::fieldValues::surfaceFieldValue::
488combineSurfaceGeometry
489(
490 faceList& faces,
491 pointField& points
492) const
493{
494 if (stObject == regionType_)
495 {
496 const auto& s = refCast<const polySurface>(obr());
497
498 if (Pstream::parRun())
499 {
500 // Dimension as fraction of surface
501 const scalar mergeDim = 1e-10*boundBox(s.points(), true).mag();
502
504 (
505 mergeDim,
506 primitivePatch(SubList<face>(s.faces()), s.points()),
507 points,
508 faces
509 );
510 }
511 else
512 {
513 faces = s.faces();
514 points = s.points();
515 }
516 }
517 else if (sampledPtr_)
518 {
519 const sampledSurface& s = *sampledPtr_;
520
521 if (Pstream::parRun())
522 {
523 // Dimension as fraction of mesh bounding box
524 const scalar mergeDim = 1e-10*mesh_.bounds().mag();
525
527 (
528 mergeDim,
529 primitivePatch(SubList<face>(s.faces()), s.points()),
530 points,
531 faces
532 );
533 }
534 else
535 {
536 faces = s.faces();
537 points = s.points();
538 }
539 }
540}
541
542
543Foam::scalar
544Foam::functionObjects::fieldValues::surfaceFieldValue::totalArea() const
545{
546 scalar totalArea = 0;
547
548 if (stObject == regionType_)
549 {
550 const auto& s = refCast<const polySurface>(obr());
551
552 totalArea = gSum(s.magSf());
553 }
554 else if (sampledPtr_)
555 {
556 totalArea = gSum(sampledPtr_->magSf());
557 }
558 else
559 {
560 totalArea = gSum(filterField(mesh_.magSf()));
561 }
562
563 return totalArea;
564}
565
566
567// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
568
569bool Foam::functionObjects::fieldValues::surfaceFieldValue::usesSf()
570const noexcept
571{
572 // Few operations do not require the Sf field
573 switch (operation_)
574 {
575 case opNone:
576 case opMin:
577 case opMax:
578 case opSum:
579 case opSumMag:
580 case opAverage:
581 return false;
582
583 default:
584 return true;
585 }
586}
587
588
589bool Foam::functionObjects::fieldValues::surfaceFieldValue::update()
590{
591 if (sampledPtr_)
592 {
593 sampledPtr_->update();
594 }
595
596 if (!needsUpdate_)
597 {
598 return false;
599 }
600
601 // Reset some values
602 totalArea_ = 0;
603 nFaces_ = 0;
604 bool checkEmptyFaces = true;
605
606 switch (regionType_)
607 {
608 case stFaceZone:
609 {
610 // Raises warning or error internally, don't check again
611 setFaceZoneFaces();
612 checkEmptyFaces = false;
613 break;
614 }
615 case stPatch:
616 {
617 // Raises warning or error internally, don't check again
618 setPatchFaces();
619 checkEmptyFaces = false;
620 break;
621 }
622 case stObject:
623 {
624 // TBD: special handling of cast errors?
625 const auto& s = refCast<const polySurface>(obr());
626 nFaces_ = returnReduce(s.size(), sumOp<label>());
627 break;
628 }
629 case stSampled:
630 {
631 nFaces_ = returnReduce(sampledPtr_->faces().size(), sumOp<label>());
632 break;
633 }
634
635 // Compiler warning if we forgot an enumeration
636 }
637
638 if (nFaces_)
639 {
640 // Appears to be successful
641 needsUpdate_ = false;
642 totalArea_ = totalArea(); // Update the area
643 nWarnings_ = 0u; // Reset the warnings counter
644 }
645 else if (checkEmptyFaces)
646 {
647 // Raise warning or error
648 refPtr<OSstream> os;
649 bool fatal = false;
650
651 ++nWarnings_; // Always increment (even if ignore etc)
652
653 switch (emptySurfaceError_)
654 {
655 case error::handlerTypes::IGNORE:
656 {
657 break;
658 }
659
660 case error::handlerTypes::WARN:
661 {
662 if (nWarnings_ <= maxWarnings)
663 {
665 }
666 break;
667 }
668
669 // STRICT / DEFAULT
670 default:
671 {
673 fatal = true;
674 break;
675 }
676 }
677
678 if (os)
679 {
680 os.ref()
681 << type() << ' ' << name() << ": "
682 << regionTypeNames_[regionType_]
683 << '(' << regionName_ << "):" << nl
684 << " Region has no faces" << endl;
685
686 if (fatal)
687 {
688 FatalError<< exit(FatalError);
689 }
690 else if (nWarnings_ == maxWarnings)
691 {
692 os.ref()
693 << "... suppressing further warnings." << nl;
694 }
695 }
696 }
697
698 Log << " total faces = " << nFaces_ << nl
699 << " total area = " << totalArea_ << nl
700 << endl;
701
702 // Emit file header on success or change of state
703 if (nWarnings_ <= 1)
704 {
705 writeFileHeader(file());
706 }
707
708 return true;
709}
710
711
712void Foam::functionObjects::fieldValues::surfaceFieldValue::writeFileHeader
713(
714 Ostream& os
715)
716{
717 if (canWriteHeader() && (operation_ != opNone))
718 {
719 writeCommented(os, "Region type : ");
720 os << regionTypeNames_[regionType_] << ' ' << regionName_ << nl;
721
722 writeHeaderValue(os, "Faces", nFaces_);
723 writeHeaderValue(os, "Area", totalArea_);
724 writeHeaderValue(os, "Scale factor", scaleFactor_);
725
726 if (weightFieldNames_.size())
727 {
728 writeHeaderValue
729 (
730 os,
731 "Weight field",
732 flatOutput(weightFieldNames_, FlatOutput::BareComma{})
733 );
734 }
735
736 writeCommented(os, "Time");
737 if (writeArea_)
738 {
739 os << tab << "Area";
740 }
741
742 // TBD: add in postOperation information?
743
744 for (const word& fieldName : fields_)
745 {
746 os << tab << operationTypeNames_[operation_]
747 << '(' << fieldName << ')';
748 }
749
750 os << endl;
751 }
752
753 writtenHeader_ = true;
754}
755
756
757template<>
758Foam::scalar
759Foam::functionObjects::fieldValues::surfaceFieldValue::processValues
760(
761 const Field<scalar>& values,
762 const vectorField& Sf,
763 const scalarField& weightField
764) const
765{
766 switch (operation_)
767 {
768 case opSumDirection:
769 {
770 const vector n(dict_.get<vector>("direction"));
771 return gSum(pos0(values*(Sf & n))*mag(values));
772 }
773 case opSumDirectionBalance:
774 {
775 const vector n(dict_.get<vector>("direction"));
776 const scalarField nv(values*(Sf & n));
777
778 return gSum(pos0(nv)*mag(values) - neg(nv)*mag(values));
779 }
780
781 case opUniformity:
782 case opWeightedUniformity:
783 case opAbsWeightedUniformity:
784 {
785 const scalar areaTotal = gSum(mag(Sf));
786 tmp<scalarField> areaVal(values * mag(Sf));
787
788 scalar mean, numer;
789
790 if (is_weightedOp() && canWeight(weightField))
791 {
792 // Weighted quantity = (Weight * phi * dA)
793
794 tmp<scalarField> weight
795 (
796 weightingFactor(weightField, is_magOp())
797 );
798
799 // Mean weighted value (area-averaged)
800 mean = gSum(weight()*areaVal()) / areaTotal;
801
802 // Abs. deviation from weighted mean value
803 numer = gSum(mag(weight*areaVal - (mean * mag(Sf))));
804 }
805 else
806 {
807 // Unweighted quantity = (1 * phi * dA)
808
809 // Mean value (area-averaged)
810 mean = gSum(areaVal()) / areaTotal;
811
812 // Abs. deviation from mean value
813 numer = gSum(mag(areaVal - (mean * mag(Sf))));
814 }
815
816 // Uniformity index
817 const scalar ui = 1 - numer/(2*mag(mean*areaTotal) + ROOTVSMALL);
818
819 return clamp(ui, zero_one{});
820 }
821
822 default:
823 {
824 // Fall through to other operations
825 return processSameTypeValues(values, Sf, weightField);
826 }
827 }
828}
829
830
831template<>
833Foam::functionObjects::fieldValues::surfaceFieldValue::processValues
834(
835 const Field<vector>& values,
836 const vectorField& Sf,
837 const scalarField& weightField
838) const
839{
840 switch (operation_)
841 {
842 case opSumDirection:
843 {
844 const vector n(dict_.get<vector>("direction").normalise());
845
846 const scalarField nv(n & values);
847 return gSum(pos0(nv)*n*(nv));
848 }
849 case opSumDirectionBalance:
850 {
851 const vector n(dict_.get<vector>("direction").normalise());
852
853 const scalarField nv(n & values);
854 return gSum(pos0(nv)*n*(nv));
855 }
856 case opAreaNormalAverage:
857 {
858 const scalar val = gSum(values & Sf)/gSum(mag(Sf));
859 return vector(val, 0, 0);
860 }
861 case opAreaNormalIntegrate:
862 {
863 const scalar val = gSum(values & Sf);
864 return vector(val, 0, 0);
865 }
866
867 case opUniformity:
868 case opWeightedUniformity:
869 case opAbsWeightedUniformity:
870 {
871 const scalar areaTotal = gSum(mag(Sf));
872 tmp<scalarField> areaVal(values & Sf);
873
874 scalar mean, numer;
875
876 if (is_weightedOp() && canWeight(weightField))
877 {
878 // Weighted quantity = (Weight * phi . dA)
879
880 tmp<scalarField> weight
881 (
882 weightingFactor(weightField, is_magOp())
883 );
884
885 // Mean weighted value (area-averaged)
886 mean = gSum(weight()*areaVal()) / areaTotal;
887
888 // Abs. deviation from weighted mean value
889 numer = gSum(mag(weight*areaVal - (mean * mag(Sf))));
890 }
891 else
892 {
893 // Unweighted quantity = (1 * phi . dA)
894
895 // Mean value (area-averaged)
896 mean = gSum(areaVal()) / areaTotal;
897
898 // Abs. deviation from mean value
899 numer = gSum(mag(areaVal - (mean * mag(Sf))));
900 }
901
902 // Uniformity index
903 const scalar ui = 1 - numer/(2*mag(mean*areaTotal) + ROOTVSMALL);
904
905 return vector(clamp(ui, zero_one{}), 0, 0);
906 }
907
908 default:
909 {
910 // Fall through to other operations
911 return processSameTypeValues(values, Sf, weightField);
912 }
913 }
914}
915
916
917// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
918
919template<>
921Foam::functionObjects::fieldValues::surfaceFieldValue::weightingFactor
922(
923 const Field<scalar>& weightField,
924 const bool useMag
925)
926{
927 if (useMag)
928 {
929 return mag(weightField);
930 }
931
932 // pass through
933 return weightField;
934}
935
936
937template<>
939Foam::functionObjects::fieldValues::surfaceFieldValue::weightingFactor
940(
941 const Field<scalar>& weightField,
942 const vectorField& Sf,
943 const bool useMag
944)
945{
946 // scalar * unit-normal
947
948 // Can skip this check - already used canWeight()
954
955 if (useMag)
956 {
957 return mag(weightField);
958 }
959
960 // pass through
961 return weightField;
962}
963
964
965template<>
967Foam::functionObjects::fieldValues::surfaceFieldValue::areaWeightingFactor
968(
969 const Field<scalar>& weightField,
970 const vectorField& Sf,
971 const bool useMag
972)
973{
974 // scalar * Area
975
976 // Can skip this check - already used canWeight()
982
983 if (useMag)
984 {
985 return mag(weightField * mag(Sf));
986 }
987
988 return (weightField * mag(Sf));
989}
990
991
992template<>
994Foam::functionObjects::fieldValues::surfaceFieldValue::weightingFactor
995(
996 const Field<vector>& weightField,
997 const vectorField& Sf,
998 const bool useMag
999)
1000{
1001 // vector (dot) unit-normal
1002
1003 // Can skip this check - already used canWeight()
1009
1010 const label len = weightField.size();
1011
1012 auto tresult = tmp<scalarField>::New(weightField.size());
1013 auto& result = tresult.ref();
1014
1015 for (label facei=0; facei < len; ++facei)
1016 {
1017 const vector unitNormal(normalised(Sf[facei]));
1018 result[facei] = (weightField[facei] & unitNormal);
1019 }
1020
1021 if (useMag)
1022 {
1023 for (scalar& val : result)
1024 {
1025 val = mag(val);
1026 }
1027 }
1028
1029 return tresult;
1030}
1031
1032
1033template<>
1035Foam::functionObjects::fieldValues::surfaceFieldValue::areaWeightingFactor
1036(
1037 const Field<vector>& weightField,
1038 const vectorField& Sf,
1039 const bool useMag
1040)
1041{
1042 // vector (dot) Area
1043
1044 // Can skip this check - already used canWeight()
1050
1051 if (useMag)
1052 {
1053 return mag(weightField & Sf);
1054 }
1055
1056 return (weightField & Sf);
1057}
1058
1059
1060// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1061
1062Foam::functionObjects::fieldValues::surfaceFieldValue::surfaceFieldValue
1063(
1064 const word& name,
1065 const Time& runTime,
1066 const dictionary& dict
1067)
1068:
1069 fieldValue(name, runTime, dict, typeName),
1070 regionType_(regionTypeNames_.get("regionType", dict)),
1071 operation_(operationTypeNames_.get("operation", dict)),
1072 postOperation_
1073 (
1074 postOperationTypeNames_.getOrDefault
1075 (
1076 "postOperation",
1077 dict,
1078 postOperationType::postOpNone,
1079 true // Failsafe behaviour
1080 )
1081 ),
1082 needsUpdate_(true),
1083 writeArea_(false),
1084 emptySurfaceError_(error::handlerTypes::DEFAULT),
1085 selectionNames_(),
1086 weightFieldNames_(),
1087 totalArea_(0),
1088 nFaces_(0),
1089 nWarnings_(0),
1090 faceId_(),
1091 facePatchId_(),
1092 faceFlip_()
1093{
1094 read(dict);
1095}
1096
1097
1098Foam::functionObjects::fieldValues::surfaceFieldValue::surfaceFieldValue
1099(
1100 const word& name,
1101 const objectRegistry& obr,
1102 const dictionary& dict
1103)
1104:
1105 fieldValue(name, obr, dict, typeName),
1106 regionType_(regionTypeNames_.get("regionType", dict)),
1107 operation_(operationTypeNames_.get("operation", dict)),
1108 postOperation_
1109 (
1110 postOperationTypeNames_.getOrDefault
1111 (
1112 "postOperation",
1113 dict,
1114 postOperationType::postOpNone,
1115 true // Failsafe behaviour
1116 )
1117 ),
1118 needsUpdate_(true),
1119 writeArea_(false),
1120 emptySurfaceError_(error::handlerTypes::DEFAULT),
1121 selectionNames_(),
1122 weightFieldNames_(),
1123 totalArea_(0),
1124 nFaces_(0),
1125 nWarnings_(0),
1126 faceId_(),
1127 facePatchId_(),
1128 faceFlip_()
1129{
1130 read(dict);
1131}
1132
1133
1134// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1135
1136// Needs completed sampledSurface, surfaceWriter
1137Foam::functionObjects::fieldValues::surfaceFieldValue::~surfaceFieldValue()
1138{}
1139
1140
1141// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1142
1143bool Foam::functionObjects::fieldValues::surfaceFieldValue::read
1144(
1145 const dictionary& dict
1146)
1147{
1149
1150 needsUpdate_ = true;
1151 writeArea_ = dict.getOrDefault("writeArea", false);
1152 emptySurfaceError_ = error::handlerNames.getOrDefault
1153 (
1154 "empty-surface",
1155 dict,
1157 true // Failsafe behaviour
1158 );
1159
1160 weightFieldNames_.clear();
1161 // future?
1162 // sampleFaceScheme_ = dict.getOrDefault<word>("sampleScheme", "cell");
1163
1164 totalArea_ = 0;
1165 nFaces_ = 0;
1166 nWarnings_ = 0;
1167 faceId_.clear();
1168 facePatchId_.clear();
1169 faceFlip_.clear();
1170 sampledPtr_.reset(nullptr);
1171 surfaceWriterPtr_.reset(nullptr);
1172
1173 // Can have "name" (word) and/or "names" (wordRes)
1174 //
1175 // If "names" exists AND contains a literal (non-regex) that can be used
1176 // as a suitable value for "name", the "name" entry becomes optional.
1177
1178 regionName_.clear();
1179 selectionNames_.clear();
1180
1181 {
1182 dict.readIfPresent("names", selectionNames_);
1183
1184 for (const auto& item : selectionNames_)
1185 {
1186 if (item.isLiteral())
1187 {
1188 regionName_ = item;
1189 break;
1190 }
1191 }
1192
1193 // The "name" entry
1194 // - mandatory if we didn't pick up a value from selectionNames_
1195 dict.readEntry
1196 (
1197 "name",
1198 regionName_,
1200 (
1201 regionName_.empty()
1203 )
1204 );
1205
1206 // Ensure there is always content for selectionNames_
1207 if (selectionNames_.empty())
1208 {
1209 selectionNames_.resize(1);
1210 selectionNames_.first() = regionName_;
1211 }
1212 }
1213
1214
1215 // Create sampled surface, but leave 'expired' (ie, no update) since it
1216 // may depend on fields or data that do not yet exist
1217 if (stSampled == regionType_)
1218 {
1219 sampledPtr_ = sampledSurface::New
1220 (
1221 name(),
1222 mesh_,
1223 dict.subDict("sampledSurfaceDict")
1224 );
1225
1226 // Internal consistency. Want face values, never point values!
1227 sampledPtr_->isPointData(false);
1228 }
1229
1230 Info<< type() << ' ' << name() << ':' << nl
1231 << " operation = ";
1232
1233 if (postOperation_ != postOpNone)
1234 {
1235 Info<< postOperationTypeNames_[postOperation_] << '('
1236 << operationTypeNames_[operation_] << ')' << nl;
1237 }
1238 else
1239 {
1240 Info<< operationTypeNames_[operation_] << nl;
1241 }
1242
1243 if (is_weightedOp())
1244 {
1245 // Can have "weightFields" or "weightField"
1246
1247 bool missing = true;
1248 if (dict.readIfPresent("weightFields", weightFieldNames_))
1249 {
1250 missing = false;
1251 }
1252 else
1253 {
1254 weightFieldNames_.resize(1);
1255
1256 if (dict.readIfPresent("weightField", weightFieldNames_.first()))
1257 {
1258 missing = false;
1259 if ("none" == weightFieldNames_.first())
1260 {
1261 // "none" == no weighting
1262 weightFieldNames_.clear();
1263 }
1264 }
1265 }
1266
1267 if (missing)
1268 {
1269 // Suggest possible alternative unweighted operation?
1271 << "The '" << operationTypeNames_[operation_]
1272 << "' operation is missing a weightField." << nl
1273 << "Either provide the weightField, "
1274 << "use weightField 'none' to suppress weighting," << nl
1275 << "or use a different operation."
1276 << exit(FatalIOError);
1277 }
1278
1279 Info<< " weight field = ";
1280 if (weightFieldNames_.empty())
1281 {
1282 Info<< "none" << nl;
1283 }
1284 else
1285 {
1286 Info<< flatOutput(weightFieldNames_) << nl;
1287 }
1288 }
1289
1290 if (stSampled == regionType_ && sampledPtr_)
1291 {
1292 Info<< " sampled surface: ";
1293 sampledPtr_->print(Info, 0);
1294 Info<< nl;
1295 }
1296
1297 if (writeFields_)
1298 {
1299 const word writerType = dict.get<word>("surfaceFormat");
1300
1301 surfaceWriterPtr_.reset
1302 (
1304 (
1305 writerType,
1307 )
1308 );
1309
1310 // Propagate field counts (per surface)
1311 surfaceWriterPtr_->nFields(fields_.size());
1312
1313 if (debug)
1314 {
1315 surfaceWriterPtr_->verbose(true);
1316 }
1317
1318 if (surfaceWriterPtr_->enabled())
1319 {
1320 Info<< " surfaceFormat = " << writerType << nl;
1321 }
1322 else
1323 {
1324 surfaceWriterPtr_->clear();
1325 }
1326 }
1327
1328 Info<< nl << endl;
1329
1330 return true;
1331}
1332
1333
1334bool Foam::functionObjects::fieldValues::surfaceFieldValue::write()
1335{
1336 if (needsUpdate_ || operation_ != opNone)
1337 {
1339 }
1340
1341 update();
1342
1343 if (operation_ != opNone)
1344 {
1345 writeCurrentTime(file());
1346 }
1347
1348 // Handle ignore/warn about empty-surface
1349 if (!nFaces_)
1350 {
1351 totalArea_ = 0; // Update the area (safety)
1352
1353 if (operation_ != opNone)
1354 {
1355 if (emptySurfaceError_ == error::handlerTypes::WARN)
1356 {
1357 if (writeArea_)
1358 {
1359 Log << " total area = " << totalArea_ << endl;
1360 file() << tab << totalArea_;
1361 }
1362
1363 file() << tab << "NaN";
1364 Log << endl;
1365 }
1366
1367 file() << endl;
1368 }
1369
1370 // Early exit on error
1371 return true;
1372 }
1373
1374 if (writeArea_)
1375 {
1376 // Update the area
1377 totalArea_ = totalArea();
1378 Log << " total area = " << totalArea_ << endl;
1379
1380 if (operation_ != opNone && UPstream::master())
1381 {
1382 file() << tab << totalArea_;
1383 }
1384 }
1385
1386 // Many operations use the Sf field
1387 vectorField Sf;
1388 if (usesSf())
1389 {
1390 if (stObject == regionType_)
1391 {
1392 const auto& s = refCast<const polySurface>(obr());
1393 Sf = s.Sf();
1394 }
1395 else if (sampledPtr_)
1396 {
1397 Sf = sampledPtr_->Sf();
1398 }
1399 else
1400 {
1401 Sf = filterField(mesh_.Sf());
1402 }
1403 }
1404
1405 // Faces and points for surface format (if specified)
1406 faceList faces;
1408
1409 if (surfaceWriterPtr_)
1410 {
1411 if (withTopologicalMerge())
1412 {
1413 combineMeshGeometry(faces, points);
1414 }
1415 else
1416 {
1417 combineSurfaceGeometry(faces, points);
1418 }
1419 }
1420
1421
1422 // Check availability and type of weight field
1423 // Only support a few weight types:
1424 // scalar: 0-N fields
1425 // vector: 0-1 fields
1426
1427 // Default is a zero-size scalar weight field (ie, weight = 1)
1428 scalarField scalarWeights;
1429 vectorField vectorWeights;
1430
1431 for (const word& weightName : weightFieldNames_)
1432 {
1433 if (validField<scalar>(weightName))
1434 {
1435 tmp<scalarField> tfld = getFieldValues<scalar>(weightName, true);
1436
1437 if (scalarWeights.empty())
1438 {
1439 scalarWeights = tfld;
1440 }
1441 else
1442 {
1443 scalarWeights *= tfld;
1444 }
1445 }
1446 else if (validField<vector>(weightName))
1447 {
1448 tmp<vectorField> tfld = getFieldValues<vector>(weightName, true);
1449
1450 if (vectorWeights.empty())
1451 {
1452 vectorWeights = tfld;
1453 }
1454 else
1455 {
1457 << "weightField " << weightName
1458 << " - only one vector weight field allowed. " << nl
1459 << "weights: " << flatOutput(weightFieldNames_) << nl
1460 << abort(FatalError);
1461 }
1462 }
1463 else if (weightName != "none")
1464 {
1465 // Silently ignore "none", flag everything else as an error
1466
1467 // TBD: treat missing "rho" like incompressible with rho = 1
1468 // and/or provided rhoRef value
1469
1471 << "weightField " << weightName
1472 << " not found or an unsupported type" << nl
1473 << abort(FatalError);
1474 }
1475 }
1476
1477
1478 // Process the fields
1479 if (returnReduceOr(!vectorWeights.empty()))
1480 {
1481 if (scalarWeights.size())
1482 {
1483 vectorWeights *= scalarWeights;
1484 }
1485
1486 writeAll(Sf, vectorWeights, points, faces);
1487 }
1488 else
1489 {
1490 writeAll(Sf, scalarWeights, points, faces);
1491 }
1492
1493
1494 if (operation_ != opNone)
1495 {
1496 file() << endl;
1497 Log << endl;
1498 }
1499
1500
1501 return true;
1502}
1503
1504
1505void Foam::functionObjects::fieldValues::surfaceFieldValue::updateMesh
1506(
1507 const mapPolyMesh& mpm
1508)
1509{
1510 needsUpdate_ = true;
1511}
1512
1513
1514void Foam::functionObjects::fieldValues::surfaceFieldValue::movePoints
1515(
1516 const polyMesh& mesh
1517)
1518{
1519 needsUpdate_ = true;
1520}
1521
1522
1523// ************************************************************************* //
#define Log
Definition PDRblock.C:28
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
EnumType getOrDefault(const word &key, const dictionary &dict, const EnumType deflt, const bool warnOnly=false) const
Find the key in the dictionary and return the corresponding enumeration element based on its name.
Definition Enum.C:173
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &pp, Field< typename PrimitivePatch< FaceList, PointField >::point_type > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::face_type > &mergedFaces, globalIndex &pointAddr, globalIndex &faceAddr, labelList &pointMergeMap=const_cast< labelList & >(labelList::null()), const bool useLocal=false)
Gather points and faces onto master and merge into single patch.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
@ WARN
Warn on errors/problems.
Definition error.H:113
@ DEFAULT
Default behaviour (local meaning).
Definition error.H:111
static const Enum< handlerTypes > handlerNames
Names of the error handler types.
Definition error.H:120
Abstract base-class for Time/database function objects.
Intermediate class for handling field value-based function objects.
Definition fieldValue.H:113
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition fieldValue.C:82
virtual bool write()
Write the function-object results.
Definition fieldValue.C:108
@ LITERAL
String literal.
Definition keyType.H:82
Registry of regIOobjects.
static autoPtr< sampledSurface > New(const word &name, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected surface.
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
A class for managing temporary objects.
Definition tmp.H:75
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
mesh update()
labelList patchIds
dynamicFvMesh & mesh
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
const pointField & points
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))
label faceId(-1)
return returnReduce(nRefine-oldNRefine, sumOp< label >())
#define WarningInFunction
Report a warning using Foam::Warning.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
type
Types of root.
Definition Roots.H:53
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Type gSum(const FieldField< Field, Type > &f)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedScalar pos0(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
messageStream Info
Information stream (stdout output on master, null elsewhere).
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
List< face > faceList
List of faces.
Definition faceListFwd.H:41
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
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
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
static constexpr const unsigned maxWarnings