Loading...
Searching...
No Matches
meshToMeshTemplates.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) 2012-2016 OpenFOAM Foundation
9 Copyright (C) 2015-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 "fvMesh.H"
30#include "volFields.H"
33#include "fvcGrad.H"
35
36// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
37
38template<class Type>
39void Foam::meshToMesh::add
40(
41 UList<Type>& fld,
42 const label offset
43) const
44{
45 forAll(fld, i)
46 {
47 fld[i] += offset;
48 }
49}
50
51
52template<class Type, class CombineOp>
54(
55 const UList<Type>& srcField,
56 const CombineOp& cop,
57 List<Type>& result
58) const
59{
60 if (result.size() != tgtToSrcCellAddr_.size())
61 {
63 << "Supplied field size is not equal to target mesh size" << nl
64 << " source mesh = " << srcToTgtCellAddr_.size() << nl
65 << " target mesh = " << tgtToSrcCellAddr_.size() << nl
66 << " supplied field = " << result.size()
67 << abort(FatalError);
68 }
69
71
72 if (distributed())
73 {
74 const mapDistribute& map = srcMapPtr_();
75
76 List<Type> work(srcField);
77 map.distribute(work);
78
79 forAll(result, celli)
80 {
81 const labelList& srcAddress = tgtToSrcCellAddr_[celli];
82 const scalarList& srcWeight = tgtToSrcCellWght_[celli];
83
84 if (srcAddress.size())
85 {
86// result[celli] = Zero;
87 result[celli] *= (1.0 - sum(srcWeight));
88 forAll(srcAddress, i)
89 {
90 label srcI = srcAddress[i];
91 scalar w = srcWeight[i];
92 cbop(result[celli], celli, work[srcI], w);
93 }
94 }
95 }
96 }
97 else
98 {
99 forAll(result, celli)
100 {
101 const labelList& srcAddress = tgtToSrcCellAddr_[celli];
102 const scalarList& srcWeight = tgtToSrcCellWght_[celli];
103
104 if (srcAddress.size())
105 {
106// result[celli] = Zero;
107 result[celli] *= (1.0 - sum(srcWeight));
108 forAll(srcAddress, i)
109 {
110 label srcI = srcAddress[i];
111 scalar w = srcWeight[i];
112 cbop(result[celli], celli, srcField[srcI], w);
113 }
115 }
116 }
117}
118
119
120template<class Type, class CombineOp>
122(
123 const UList<Type>& srcField,
124 const UList<typename outerProduct<vector, Type>::type>& srcGradField,
125 const CombineOp& cop,
126 List<Type>& result
127) const
128{
129 if (result.size() != tgtToSrcCellAddr_.size())
130 {
132 << "Supplied field size is not equal to target mesh size" << nl
133 << " source mesh = " << srcToTgtCellAddr_.size() << nl
134 << " target mesh = " << tgtToSrcCellAddr_.size() << nl
135 << " supplied field = " << result.size()
136 << abort(FatalError);
137 }
138
140
141 if (distributed())
142 {
143 if (returnReduceAnd(tgtToSrcCellVec_.empty()))
144 {
145 // No correction vectors calculated. Fall back to first order.
146 mapSrcToTgt(srcField, cop, result);
147 return;
148 }
149
150 const mapDistribute& map = srcMapPtr_();
151
152 List<Type> work(srcField);
153 map.distribute(work);
154
156 (
157 srcGradField
158 );
159 map.distribute(workGrad);
160
161 forAll(result, cellI)
162 {
163 const labelList& srcAddress = tgtToSrcCellAddr_[cellI];
164 const scalarList& srcWeight = tgtToSrcCellWght_[cellI];
165 const pointList& srcVec = tgtToSrcCellVec_[cellI];
166
167 if (srcAddress.size())
168 {
169 result[cellI] *= (1.0 - sum(srcWeight));
170 forAll(srcAddress, i)
171 {
172 label srcI = srcAddress[i];
173 scalar w = srcWeight[i];
174 const vector& v = srcVec[i];
175 const Type srcVal = work[srcI]+(workGrad[srcI]&v);
176 cbop(result[cellI], cellI, srcVal, w);
177 }
178 }
179 }
180 }
181 else
182 {
183 if (tgtToSrcCellVec_.empty())
184 {
185 // No correction vectors calculated. Fall back to first order.
186 mapSrcToTgt(srcField, cop, result);
187 return;
188 }
189
190 forAll(result, cellI)
191 {
192 const labelList& srcAddress = tgtToSrcCellAddr_[cellI];
193 const scalarList& srcWeight = tgtToSrcCellWght_[cellI];
194 const pointList& srcVec = tgtToSrcCellVec_[cellI];
195
196 if (srcAddress.size())
197 {
198 // Do non-conservative interpolation
199 result[cellI] *= (1.0 - sum(srcWeight));
200 forAll(srcAddress, i)
201 {
202 label srcI = srcAddress[i];
203 scalar w = srcWeight[i];
204 const vector& v = srcVec[i];
205 const Type srcVal = srcField[srcI]+(srcGradField[srcI]&v);
206 cbop(result[cellI], cellI, srcVal, w);
207 }
209 }
210 }
211}
212
213
214template<class Type, class CombineOp>
215Foam::tmp<Foam::Field<Type>> Foam::meshToMesh::mapSrcToTgt
216(
217 const Field<Type>& srcField,
218 const CombineOp& cop
219) const
220{
221 auto tresult = tmp<Field<Type>>::New(tgtToSrcCellAddr_.size(), Zero);
222
223 mapSrcToTgt(srcField, cop, tresult.ref());
224
225 return tresult;
226}
227
228
229template<class Type, class CombineOp>
231(
232 const tmp<Field<Type>>& tsrcField,
233 const CombineOp& cop
234) const
235{
236 return mapSrcToTgt(tsrcField(), cop);
237}
238
239
240template<class Type>
242(
243 const Field<Type>& srcField
244) const
245{
246 return mapSrcToTgt(srcField, plusEqOp<Type>());
247}
248
249
250template<class Type>
252(
253 const tmp<Field<Type>>& tsrcField
254) const
255{
256 return mapSrcToTgt(tsrcField());
257}
258
259
260template<class Type, class CombineOp>
262(
263 const UList<Type>& tgtField,
264 const CombineOp& cop,
265 List<Type>& result
266) const
267{
268 if (result.size() != srcToTgtCellAddr_.size())
269 {
271 << "Supplied field size is not equal to source mesh size" << nl
272 << " source mesh = " << srcToTgtCellAddr_.size() << nl
273 << " target mesh = " << tgtToSrcCellAddr_.size() << nl
274 << " supplied field = " << result.size()
275 << abort(FatalError);
276 }
277
278 multiplyWeightedOp<Type, CombineOp> cbop(cop);
279
280 if (distributed())
281 {
282 const mapDistribute& map = tgtMapPtr_();
283
284 List<Type> work(tgtField);
285 map.distribute(work);
286
287 forAll(result, celli)
288 {
289 const labelList& tgtAddress = srcToTgtCellAddr_[celli];
290 const scalarList& tgtWeight = srcToTgtCellWght_[celli];
291
292 if (tgtAddress.size())
293 {
294 result[celli] *= (1.0 - sum(tgtWeight));
295 forAll(tgtAddress, i)
296 {
297 label tgtI = tgtAddress[i];
298 scalar w = tgtWeight[i];
299 cbop(result[celli], celli, work[tgtI], w);
300 }
301 }
302 }
303 }
304 else
305 {
306 forAll(result, celli)
307 {
308 const labelList& tgtAddress = srcToTgtCellAddr_[celli];
309 const scalarList& tgtWeight = srcToTgtCellWght_[celli];
310
311 if (tgtAddress.size())
312 {
313 result[celli] *= (1.0 - sum(tgtWeight));
314 forAll(tgtAddress, i)
315 {
316 label tgtI = tgtAddress[i];
317 scalar w = tgtWeight[i];
318 cbop(result[celli], celli, tgtField[tgtI], w);
319 }
321 }
322 }
323}
324
325
326template<class Type, class CombineOp>
328(
329 const UList<Type>& tgtField,
330 const UList<typename outerProduct<vector, Type>::type>& tgtGradField,
331 const CombineOp& cop,
332 List<Type>& result
333) const
334{
335 if (result.size() != srcToTgtCellAddr_.size())
336 {
338 << "Supplied field size is not equal to source mesh size" << nl
339 << " source mesh = " << srcToTgtCellAddr_.size() << nl
340 << " target mesh = " << tgtToSrcCellAddr_.size() << nl
341 << " supplied field = " << result.size()
342 << abort(FatalError);
343 }
344
346
347 if (distributed())
348 {
349 if (returnReduceAnd(srcToTgtCellVec_.empty()))
350 {
351 // No correction vectors calculated. Fall back to first order.
352 mapTgtToSrc(tgtField, cop, result);
353 return;
354 }
355
356 const mapDistribute& map = tgtMapPtr_();
357
358 List<Type> work(tgtField);
359 map.distribute(work);
360
362 (
363 tgtGradField
364 );
365 map.distribute(workGrad);
366
367 forAll(result, cellI)
368 {
369 const labelList& tgtAddress = srcToTgtCellAddr_[cellI];
370 const scalarList& tgtWeight = srcToTgtCellWght_[cellI];
371 const pointList& tgtVec = srcToTgtCellVec_[cellI];
372
373 if (tgtAddress.size())
374 {
375 result[cellI] *= (1.0 - sum(tgtWeight));
376 forAll(tgtAddress, i)
377 {
378 label tgtI = tgtAddress[i];
379 scalar w = tgtWeight[i];
380 const vector& v = tgtVec[i];
381 const Type tgtVal = work[tgtI]+(workGrad[tgtI]&v);
382 cbop(result[cellI], cellI, tgtVal, w);
383 }
384 }
385 }
386 }
387 else
388 {
389 forAll(result, cellI)
390 {
391 const labelList& tgtAddress = srcToTgtCellAddr_[cellI];
392 const scalarList& tgtWeight = srcToTgtCellWght_[cellI];
393 const pointList& tgtVec = srcToTgtCellVec_[cellI];
394
395 if (tgtAddress.size())
396 {
397 result[cellI] *= (1.0 - sum(tgtWeight));
398 forAll(tgtAddress, i)
399 {
400 label tgtI = tgtAddress[i];
401 scalar w = tgtWeight[i];
402 const vector& v = tgtVec[i];
403 const Type tgtVal = tgtField[tgtI]+(tgtGradField[tgtI]&v);
404 cbop(result[cellI], cellI, tgtVal, w);
405 }
407 }
408 }
409}
410
411
412template<class Type, class CombineOp>
413Foam::tmp<Foam::Field<Type>> Foam::meshToMesh::mapTgtToSrc
414(
415 const Field<Type>& tgtField,
416 const CombineOp& cop
417) const
418{
419 auto tresult = tmp<Field<Type>>::New(srcToTgtCellAddr_.size(), Zero);
420
421 mapTgtToSrc(tgtField, cop, tresult.ref());
422
423 return tresult;
424}
425
426
427template<class Type, class CombineOp>
429(
430 const tmp<Field<Type>>& ttgtField,
431 const CombineOp& cop
432) const
433{
434 return mapTgtToSrc(ttgtField(), cop);
435}
436
437
438template<class Type>
440(
441 const Field<Type>& tgtField
442) const
443{
444 return mapTgtToSrc(tgtField, plusEqOp<Type>());
445}
446
447
448template<class Type>
450(
451 const tmp<Field<Type>>& ttgtField
452) const
453{
454 return mapTgtToSrc(ttgtField(), plusEqOp<Type>());
455}
456
457
458template<class Type, class CombineOp>
459void Foam::meshToMesh::mapInternalSrcToTgt
460(
462 const CombineOp& cop,
463 VolumeField<Type>& result,
464 const bool secondOrder
465) const
466{
467 if (secondOrder && returnReduceOr(tgtToSrcCellVec_.size()))
468 {
469 mapSrcToTgt
470 (
471 field,
472 fvc::grad(field)().primitiveField(),
473 cop,
474 result.primitiveFieldRef()
475 );
476 }
477 else
478 {
479 mapSrcToTgt(field, cop, result.primitiveFieldRef());
480 }
481}
482
483
484template<class Type, class CombineOp>
485void Foam::meshToMesh::mapAndOpSrcToTgt
486(
488 const Field<Type>& srcField,
489 Field<Type>& tgtField,
490 const CombineOp& cop
491) const
492{
493 tgtField = Type(Zero);
494
495 AMI.interpolateToTarget
496 (
497 srcField,
499 tgtField,
501 );
502}
503
504
505template<class Type, class CombineOp>
507(
509 const CombineOp& cop,
510 VolumeField<Type>& result,
511 const bool secondOrder
512) const
513{
514 mapInternalSrcToTgt(field, cop, result, secondOrder);
515
516 const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
517
518 auto& resultBf = result.boundaryFieldRef();
519
520 forAll(AMIList, i)
521 {
522 label srcPatchi = srcPatchID_[i];
523 label tgtPatchi = tgtPatchID_[i];
524
525 const fvPatchField<Type>& srcField = field.boundaryField()[srcPatchi];
526 fvPatchField<Type>& tgtField = resultBf[tgtPatchi];
527
528 // Clone and map (since rmap does not do general mapping)
530 (
532 (
533 srcField,
534 tgtField.patch(),
535 result(),
537 (
538 AMIList[i].singlePatchProc(),
539 AMIList[i].comm(), // for testing only
540 AMIList[i].hasSrcMap(), // pointer to map
541 AMIList[i].tgtAddress(),
542 AMIList[i].tgtWeights()
543 )
544 )
545 );
546
547 // Transfer all mapped quantities (value and e.g. gradient) onto
548 // tgtField. Value will get overwritten below.
549 tgtField.rmap(tnewTgt(), identity(tgtField.size()));
550
551 // Override value to account for CombineOp (note: is dummy template
552 // specialisation for plusEqOp)
553 mapAndOpSrcToTgt(AMIList[i], srcField, tgtField, cop);
554 }
555
556 forAll(cuttingPatches_, i)
557 {
558 label patchi = cuttingPatches_[i];
559 fvPatchField<Type>& pf = resultBf[patchi];
560 pf == pf.patchInternalField();
561 }
562}
563
564
565template<class Type, class CombineOp>
568(
570 const CombineOp& cop,
571 const bool secondOrder
572) const
573{
574 const fvMesh& tgtMesh = static_cast<const fvMesh&>(tgtRegion_);
575
576 const fvBoundaryMesh& tgtBm = tgtMesh.boundary();
577 const auto& srcBfld = field.boundaryField();
578
579 PtrList<fvPatchField<Type>> tgtPatchFields(tgtBm.size());
580
581 // construct tgt boundary patch types as copy of 'field' boundary types
582 // note: this will provide place holders for fields with additional
583 // entries, but these values will need to be reset
584 forAll(tgtPatchID_, i)
585 {
586 label srcPatchi = srcPatchID_[i];
587 label tgtPatchi = tgtPatchID_[i];
588
589 if (!tgtPatchFields.set(tgtPatchi))
590 {
591 tgtPatchFields.set
592 (
593 tgtPatchi,
595 (
596 srcBfld[srcPatchi],
597 tgtMesh.boundary()[tgtPatchi],
600 (
601 labelList(tgtMesh.boundary()[tgtPatchi].size(), -1)
602 )
603 )
604 );
605 }
606 }
607
608 // Any unset tgtPatchFields become calculated
609 forAll(tgtPatchFields, tgtPatchi)
610 {
611 if (!tgtPatchFields.set(tgtPatchi))
612 {
613 // Note: use factory New method instead of direct generation of
614 // calculated so we keep constraints
615 tgtPatchFields.set
616 (
617 tgtPatchi,
619 (
621 tgtMesh.boundary()[tgtPatchi],
623 )
624 );
625 }
626 }
627
628 auto tresult =
630 (
631 tgtMesh.newIOobject
632 (
634 (
635 type(),
636 "interpolate(" + field.name() + ")"
637 )
638 ),
639 tgtMesh,
640 field.dimensions(),
641 Field<Type>(tgtMesh.nCells(), Zero),
642 tgtPatchFields
643 );
644
645 mapSrcToTgt(field, cop, tresult.ref(), secondOrder);
647 return tresult;
648}
649
650
651template<class Type, class CombineOp>
654(
655 const tmp<VolumeField<Type>>& tfield,
656 const CombineOp& cop,
657 const bool secondOrder
658) const
660 return mapSrcToTgt(tfield(), cop, secondOrder);
661}
662
663
664template<class Type>
667(
669 const bool secondOrder
670) const
672 return mapSrcToTgt(field, plusEqOp<Type>(), secondOrder);
673}
674
675
676template<class Type>
679(
680 const tmp<VolumeField<Type>>& tfield,
681 const bool secondOrder
682) const
683{
684 return mapSrcToTgt(tfield(), plusEqOp<Type>(), secondOrder);
685}
686
687
688template<class Type, class CombineOp>
689void Foam::meshToMesh::mapInternalTgtToSrc
690(
692 const CombineOp& cop,
693 VolumeField<Type>& result,
694 const bool secondOrder
695) const
696{
697 if (secondOrder && returnReduceOr(srcToTgtCellVec_.size()))
698 {
699 mapTgtToSrc
700 (
701 field,
702 fvc::grad(field)().primitiveField(),
703 cop,
704 result.primitiveFieldRef()
705 );
706 }
707 else
708 {
709 mapTgtToSrc(field, cop, result.primitiveFieldRef());
710 }
711}
712
713
714template<class Type, class CombineOp>
715void Foam::meshToMesh::mapAndOpTgtToSrc
716(
718 Field<Type>& srcField,
719 const Field<Type>& tgtField,
720 const CombineOp& cop
721) const
722{
723 srcField = Type(Zero);
724
725 AMI.interpolateToSource
726 (
727 tgtField,
729 srcField,
731 );
732}
733
734
735template<class Type, class CombineOp>
737(
739 const CombineOp& cop,
740 VolumeField<Type>& result,
741 const bool secondOrder
742) const
743{
744 mapInternalTgtToSrc(field, cop, result, secondOrder);
745
746 const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
747
748 forAll(AMIList, i)
749 {
750 label srcPatchi = srcPatchID_[i];
751 label tgtPatchi = tgtPatchID_[i];
752
753 fvPatchField<Type>& srcField = result.boundaryFieldRef()[srcPatchi];
754 const fvPatchField<Type>& tgtField = field.boundaryField()[tgtPatchi];
755
756 // Clone and map (since rmap does not do general mapping)
758 (
760 (
761 tgtField,
762 srcField.patch(),
763 result(),
765 (
766 AMIList[i].singlePatchProc(),
767 AMIList[i].comm(), // only used for testing
768 AMIList[i].hasTgtMap(), // pointer to map
769 AMIList[i].srcAddress(),
770 AMIList[i].srcWeights()
771 )
772 )
773 );
774 // Transfer all mapped quantities (value and e.g. gradient) onto
775 // srcField. Value will get overwritten below
776 srcField.rmap(tnewSrc(), identity(srcField.size()));
777
778 // Override value to account for CombineOp (could be dummy for
779 // plusEqOp)
780 mapAndOpTgtToSrc(AMIList[i], srcField, tgtField, cop);
781 }
782
783 forAll(cuttingPatches_, i)
784 {
785 label patchi = cuttingPatches_[i];
786 fvPatchField<Type>& pf = result.boundaryFieldRef()[patchi];
787 pf == pf.patchInternalField();
788 }
789}
790
791
792template<class Type, class CombineOp>
795(
797 const CombineOp& cop,
798 const bool secondOrder
799) const
800{
801 const fvMesh& srcMesh = static_cast<const fvMesh&>(srcRegion_);
802
803 const fvBoundaryMesh& srcBm = srcMesh.boundary();
804 const auto& tgtBfld = field.boundaryField();
805
806 PtrList<fvPatchField<Type>> srcPatchFields(srcBm.size());
807
808 // construct src boundary patch types as copy of 'field' boundary types
809 // note: this will provide place holders for fields with additional
810 // entries, but these values will need to be reset
811 forAll(srcPatchID_, i)
812 {
813 label srcPatchi = srcPatchID_[i];
814 label tgtPatchi = tgtPatchID_[i];
815
816 if (!srcPatchFields.set(srcPatchi))
817 {
818 srcPatchFields.set
819 (
820 srcPatchi,
822 (
823 tgtBfld[tgtPatchi],
824 srcMesh.boundary()[srcPatchi],
827 (
828 labelList(srcMesh.boundary()[srcPatchi].size(), -1)
829 )
830 )
831 );
832 }
833 }
834
835 // Any unset srcPatchFields become calculated
836 forAll(srcPatchFields, srcPatchi)
837 {
838 if (!srcPatchFields.set(srcPatchi))
839 {
840 // Note: use factory New method instead of direct generation of
841 // calculated so we keep constraints
842 srcPatchFields.set
843 (
844 srcPatchi,
846 (
848 srcMesh.boundary()[srcPatchi],
850 )
851 );
852 }
853 }
854
855 auto tresult =
857 (
858 srcMesh.newIOobject
859 (
861 (
862 type(),
863 "interpolate(" + field.name() + ")"
864 )
865 ),
866 srcMesh,
867 field.dimensions(),
868 Field<Type>(srcMesh.nCells(), Zero),
869 srcPatchFields
870 );
871
872 mapTgtToSrc(field, cop, tresult.ref(), secondOrder);
874 return tresult;
875}
876
877
878template<class Type, class CombineOp>
881(
882 const tmp<VolumeField<Type>>& tfield,
883 const CombineOp& cop,
884 const bool secondOrder
885) const
887 return mapTgtToSrc(tfield(), cop, secondOrder);
888}
889
890
891template<class Type>
894(
896 const bool secondOrder
897) const
899 return mapTgtToSrc(field, plusEqOp<Type>(), secondOrder);
900}
901
902
903template<class Type>
906(
907 const tmp<VolumeField<Type>>& tfield,
908 const bool secondOrder
909) const
910{
911 return mapTgtToSrc(tfield(), plusEqOp<Type>(), secondOrder);
912}
913
914
915// ************************************************************************* //
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))
static const this_type & null() noexcept
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
static const UList< T > & null() noexcept
Return a null UList (reference to a nullObject). Behaves like an empty UList.
Definition UList.H:225
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
FieldMapper with weighted mapping from (optionally remote) quantities.
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
static const word & calculatedType() noexcept
The type name for calculated patch fields.
const fvPatch & patch() const noexcept
Return the patch.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Class containing processor-to-processor mapping information.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
bool distributed() const noexcept
Distributed across processors (singleMeshProc == -1).
Definition meshToMeshI.H:80
const mapDistribute * hasTgtMap() const noexcept
Pointer to the target map (if distributed). Can be checked as a bool.
Definition meshToMeshI.H:94
Foam::tmp< Foam::Field< Type > > mapTgtToSrc(const tmp< Field< Type > > &ttgtField, const CombineOp &cop) const
const PtrList< AMIPatchToPatchInterpolation > & patchAMIs() const
Return the list of AMIs between source and target patches.
const mapDistribute * hasSrcMap() const noexcept
Pointer to the source map (if distributed). Can be checked as a bool.
Definition meshToMeshI.H:87
Foam::tmp< Foam::Field< Type > > mapSrcToTgt(const tmp< Field< Type > > &tsrcField, const CombineOp &cop) const
void mapTgtToSrc(const UList< Type > &tgtFld, const CombineOp &cop, List< Type > &result) const
Map field from tgt to src mesh with defined operation.
void mapSrcToTgt(const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
IOobject newIOobject(const word &name, IOobjectOption ioOpt) const
Create an IOobject at the current time instance (timeName) with the specified options.
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
label nCells() const noexcept
Number of mesh cells.
A class for managing temporary objects.
Definition tmp.H:75
rDeltaTY field()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Calculate the gradient of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< Type, fvPatchField, volMesh > VolumeField
A volume field for a given type.
List< label > labelList
A List of labels.
Definition List.H:62
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
List< point > pointList
List of point.
Definition pointList.H:32
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...
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
errorManip< error > abort(error &err)
Definition errorManip.H:139
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299