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 -------------------------------------------------------------------------------
11 License
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"
32 #include "calculatedFvPatchField.H"
33 #include "fvcGrad.H"
35 
36 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
37 
38 template<class Type>
39 void 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 
52 template<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 
70  multiplyWeightedOp<Type, CombineOp> cbop(cop);
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  }
114  }
115  }
116  }
117 }
118 
119 
120 template<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 
139  multiplyWeightedOp<Type, CombineOp> cbop(cop);
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  }
208  }
209  }
210  }
211 }
212 
213 
214 template<class Type, class CombineOp>
216 (
217  const Field<Type>& srcField,
218  const CombineOp& cop
219 ) const
220 {
221  tmp<Field<Type>> tresult
222  (
223  new Field<Type>
224  (
225  tgtToSrcCellAddr_.size(),
226  Zero
227  )
228  );
229 
230  mapSrcToTgt(srcField, cop, tresult.ref());
232  return tresult;
233 }
234 
235 
236 template<class Type, class CombineOp>
238 (
239  const tmp<Field<Type>>& tsrcField,
240  const CombineOp& cop
241 ) const
242 {
243  return mapSrcToTgt(tsrcField(), cop);
244 }
245 
246 
247 template<class Type>
249 (
250  const Field<Type>& srcField
251 ) const
252 {
253  return mapSrcToTgt(srcField, plusEqOp<Type>());
254 }
255 
256 
257 template<class Type>
259 (
260  const tmp<Field<Type>>& tsrcField
261 ) const
262 {
263  return mapSrcToTgt(tsrcField());
264 }
265 
266 
267 template<class Type, class CombineOp>
269 (
270  const UList<Type>& tgtField,
271  const CombineOp& cop,
272  List<Type>& result
273 ) const
274 {
275  if (result.size() != srcToTgtCellAddr_.size())
276  {
278  << "Supplied field size is not equal to source mesh size" << nl
279  << " source mesh = " << srcToTgtCellAddr_.size() << nl
280  << " target mesh = " << tgtToSrcCellAddr_.size() << nl
281  << " supplied field = " << result.size()
282  << abort(FatalError);
283  }
284 
285  multiplyWeightedOp<Type, CombineOp> cbop(cop);
286 
287  if (distributed())
288  {
289  const mapDistribute& map = tgtMapPtr_();
290 
291  List<Type> work(tgtField);
292  map.distribute(work);
293 
294  forAll(result, celli)
295  {
296  const labelList& tgtAddress = srcToTgtCellAddr_[celli];
297  const scalarList& tgtWeight = srcToTgtCellWght_[celli];
298 
299  if (tgtAddress.size())
300  {
301  result[celli] *= (1.0 - sum(tgtWeight));
302  forAll(tgtAddress, i)
303  {
304  label tgtI = tgtAddress[i];
305  scalar w = tgtWeight[i];
306  cbop(result[celli], celli, work[tgtI], w);
307  }
308  }
309  }
310  }
311  else
312  {
313  forAll(result, celli)
314  {
315  const labelList& tgtAddress = srcToTgtCellAddr_[celli];
316  const scalarList& tgtWeight = srcToTgtCellWght_[celli];
317 
318  if (tgtAddress.size())
319  {
320  result[celli] *= (1.0 - sum(tgtWeight));
321  forAll(tgtAddress, i)
322  {
323  label tgtI = tgtAddress[i];
324  scalar w = tgtWeight[i];
325  cbop(result[celli], celli, tgtField[tgtI], w);
326  }
327  }
328  }
329  }
330 }
331 
332 
333 template<class Type, class CombineOp>
335 (
336  const UList<Type>& tgtField,
337  const UList<typename outerProduct<vector, Type>::type>& tgtGradField,
338  const CombineOp& cop,
339  List<Type>& result
340 ) const
341 {
342  if (result.size() != srcToTgtCellAddr_.size())
343  {
345  << "Supplied field size is not equal to source mesh size" << nl
346  << " source mesh = " << srcToTgtCellAddr_.size() << nl
347  << " target mesh = " << tgtToSrcCellAddr_.size() << nl
348  << " supplied field = " << result.size()
349  << abort(FatalError);
350  }
351 
352  multiplyWeightedOp<Type, CombineOp> cbop(cop);
353 
354  if (distributed())
355  {
356  if (returnReduceAnd(srcToTgtCellVec_.empty()))
357  {
358  // No correction vectors calculated. Fall back to first order.
359  mapTgtToSrc(tgtField, cop, result);
360  return;
361  }
362 
363  const mapDistribute& map = tgtMapPtr_();
364 
365  List<Type> work(tgtField);
366  map.distribute(work);
367 
369  (
370  tgtGradField
371  );
372  map.distribute(workGrad);
373 
374  forAll(result, cellI)
375  {
376  const labelList& tgtAddress = srcToTgtCellAddr_[cellI];
377  const scalarList& tgtWeight = srcToTgtCellWght_[cellI];
378  const pointList& tgtVec = srcToTgtCellVec_[cellI];
379 
380  if (tgtAddress.size())
381  {
382  result[cellI] *= (1.0 - sum(tgtWeight));
383  forAll(tgtAddress, i)
384  {
385  label tgtI = tgtAddress[i];
386  scalar w = tgtWeight[i];
387  const vector& v = tgtVec[i];
388  const Type tgtVal = work[tgtI]+(workGrad[tgtI]&v);
389  cbop(result[cellI], cellI, tgtVal, w);
390  }
391  }
392  }
393  }
394  else
395  {
396  forAll(result, cellI)
397  {
398  const labelList& tgtAddress = srcToTgtCellAddr_[cellI];
399  const scalarList& tgtWeight = srcToTgtCellWght_[cellI];
400  const pointList& tgtVec = srcToTgtCellVec_[cellI];
401 
402  if (tgtAddress.size())
403  {
404  result[cellI] *= (1.0 - sum(tgtWeight));
405  forAll(tgtAddress, i)
406  {
407  label tgtI = tgtAddress[i];
408  scalar w = tgtWeight[i];
409  const vector& v = tgtVec[i];
410  const Type tgtVal = tgtField[tgtI]+(tgtGradField[tgtI]&v);
411  cbop(result[cellI], cellI, tgtVal, w);
412  }
413  }
414  }
415  }
416 }
417 
418 
419 template<class Type, class CombineOp>
421 (
422  const Field<Type>& tgtField,
423  const CombineOp& cop
424 ) const
425 {
426  tmp<Field<Type>> tresult
427  (
428  new Field<Type>
429  (
430  srcToTgtCellAddr_.size(),
431  Zero
432  )
433  );
434 
435  mapTgtToSrc(tgtField, cop, tresult.ref());
437  return tresult;
438 }
439 
440 
441 template<class Type, class CombineOp>
443 (
444  const tmp<Field<Type>>& ttgtField,
445  const CombineOp& cop
446 ) const
447 {
448  return mapTgtToSrc(ttgtField(), cop);
449 }
450 
451 
452 template<class Type>
454 (
455  const Field<Type>& tgtField
456 ) const
457 {
458  return mapTgtToSrc(tgtField, plusEqOp<Type>());
459 }
460 
461 
462 template<class Type>
464 (
465  const tmp<Field<Type>>& ttgtField
466 ) const
467 {
468  return mapTgtToSrc(ttgtField(), plusEqOp<Type>());
469 }
470 
471 
472 template<class Type, class CombineOp>
473 void Foam::meshToMesh::mapInternalSrcToTgt
474 (
475  const VolumeField<Type>& field,
476  const CombineOp& cop,
477  VolumeField<Type>& result,
478  const bool secondOrder
479 ) const
480 {
481  if (secondOrder && returnReduceOr(tgtToSrcCellVec_.size()))
482  {
483  mapSrcToTgt
484  (
485  field,
486  fvc::grad(field)().primitiveField(),
487  cop,
488  result.primitiveFieldRef()
489  );
490  }
491  else
492  {
493  mapSrcToTgt(field, cop, result.primitiveFieldRef());
494  }
495 }
496 
497 
498 template<class Type, class CombineOp>
499 void Foam::meshToMesh::mapAndOpSrcToTgt
500 (
501  const AMIPatchToPatchInterpolation& AMI,
502  const Field<Type>& srcField,
503  Field<Type>& tgtField,
504  const CombineOp& cop
505 ) const
506 {
507  tgtField = Type(Zero);
508 
509  AMI.interpolateToTarget
510  (
511  srcField,
512  multiplyWeightedOp<Type, CombineOp>(cop),
513  tgtField,
515  );
516 }
517 
518 
519 template<class Type, class CombineOp>
521 (
522  const VolumeField<Type>& field,
523  const CombineOp& cop,
524  VolumeField<Type>& result,
525  const bool secondOrder
526 ) const
527 {
528  mapInternalSrcToTgt(field, cop, result, secondOrder);
529 
530  const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
531 
532  auto& resultBf = result.boundaryFieldRef();
533 
534  forAll(AMIList, i)
535  {
536  label srcPatchi = srcPatchID_[i];
537  label tgtPatchi = tgtPatchID_[i];
538 
539  const fvPatchField<Type>& srcField = field.boundaryField()[srcPatchi];
540  fvPatchField<Type>& tgtField = resultBf[tgtPatchi];
541 
542  // Clone and map (since rmap does not do general mapping)
543  tmp<fvPatchField<Type>> tnewTgt
544  (
546  (
547  srcField,
548  tgtField.patch(),
549  result(),
551  (
552  AMIList[i].singlePatchProc(),
553  (
554  AMIList[i].distributed()
555  ? AMIList[i].hasSrcMap() // pointer to map
556  : nullptr
557  ),
558  AMIList[i].tgtAddress(),
559  AMIList[i].tgtWeights()
560  )
561  )
562  );
563 
564  // Transfer all mapped quantities (value and e.g. gradient) onto
565  // tgtField. Value will get overwritten below.
566  tgtField.rmap(tnewTgt(), identity(tgtField.size()));
567 
568  // Override value to account for CombineOp (note: is dummy template
569  // specialisation for plusEqOp)
570  mapAndOpSrcToTgt(AMIList[i], srcField, tgtField, cop);
571  }
572 
573  forAll(cuttingPatches_, i)
574  {
575  label patchi = cuttingPatches_[i];
576  fvPatchField<Type>& pf = resultBf[patchi];
577  pf == pf.patchInternalField();
578  }
579 }
580 
581 
582 template<class Type, class CombineOp>
585 (
586  const VolumeField<Type>& field,
587  const CombineOp& cop,
588  const bool secondOrder
589 ) const
590 {
591  const fvMesh& tgtMesh = static_cast<const fvMesh&>(tgtRegion_);
592 
593  const fvBoundaryMesh& tgtBm = tgtMesh.boundary();
594  const auto& srcBfld = field.boundaryField();
595 
596  PtrList<fvPatchField<Type>> tgtPatchFields(tgtBm.size());
597 
598  // construct tgt boundary patch types as copy of 'field' boundary types
599  // note: this will provide place holders for fields with additional
600  // entries, but these values will need to be reset
601  forAll(tgtPatchID_, i)
602  {
603  label srcPatchi = srcPatchID_[i];
604  label tgtPatchi = tgtPatchID_[i];
605 
606  if (!tgtPatchFields.set(tgtPatchi))
607  {
608  tgtPatchFields.set
609  (
610  tgtPatchi,
612  (
613  srcBfld[srcPatchi],
614  tgtMesh.boundary()[tgtPatchi],
617  (
618  labelList(tgtMesh.boundary()[tgtPatchi].size(), -1)
619  )
620  )
621  );
622  }
623  }
624 
625  // Any unset tgtPatchFields become calculated
626  forAll(tgtPatchFields, tgtPatchi)
627  {
628  if (!tgtPatchFields.set(tgtPatchi))
629  {
630  // Note: use factory New method instead of direct generation of
631  // calculated so we keep constraints
632  tgtPatchFields.set
633  (
634  tgtPatchi,
636  (
638  tgtMesh.boundary()[tgtPatchi],
640  )
641  );
642  }
643  }
644 
645  auto tresult =
646  tmp<VolumeField<Type>>::New
647  (
648  IOobject
649  (
650  type() + ":interpolate(" + field.name() + ")",
651  tgtMesh.time().timeName(),
652  tgtMesh,
655  ),
656  tgtMesh,
657  field.dimensions(),
658  Field<Type>(tgtMesh.nCells(), Zero),
659  tgtPatchFields
660  );
661 
662  mapSrcToTgt(field, cop, tresult.ref(), secondOrder);
663 
664  return tresult;
665 }
666 
667 
668 template<class Type, class CombineOp>
671 (
672  const tmp<VolumeField<Type>>& tfield,
673  const CombineOp& cop,
674  const bool secondOrder
675 ) const
676 {
677  return mapSrcToTgt(tfield(), cop, secondOrder);
678 }
679 
680 
681 template<class Type>
684 (
685  const VolumeField<Type>& field,
686  const bool secondOrder
687 ) const
688 {
689  return mapSrcToTgt(field, plusEqOp<Type>(), secondOrder);
690 }
691 
692 
693 template<class Type>
696 (
697  const tmp<VolumeField<Type>>& tfield,
698  const bool secondOrder
699 ) const
700 {
701  return mapSrcToTgt(tfield(), plusEqOp<Type>(), secondOrder);
702 }
703 
704 
705 template<class Type, class CombineOp>
706 void Foam::meshToMesh::mapInternalTgtToSrc
707 (
708  const VolumeField<Type>& field,
709  const CombineOp& cop,
710  VolumeField<Type>& result,
711  const bool secondOrder
712 ) const
713 {
714  if (secondOrder && returnReduceOr(srcToTgtCellVec_.size()))
715  {
716  mapTgtToSrc
717  (
718  field,
719  fvc::grad(field)().primitiveField(),
720  cop,
721  result.primitiveFieldRef()
722  );
723  }
724  else
725  {
726  mapTgtToSrc(field, cop, result.primitiveFieldRef());
727  }
728 }
729 
730 
731 template<class Type, class CombineOp>
732 void Foam::meshToMesh::mapAndOpTgtToSrc
733 (
734  const AMIPatchToPatchInterpolation& AMI,
735  Field<Type>& srcField,
736  const Field<Type>& tgtField,
737  const CombineOp& cop
738 ) const
739 {
740  srcField = Type(Zero);
741 
742  AMI.interpolateToSource
743  (
744  tgtField,
745  multiplyWeightedOp<Type, CombineOp>(cop),
746  srcField,
748  );
749 }
750 
751 
752 template<class Type, class CombineOp>
754 (
755  const VolumeField<Type>& field,
756  const CombineOp& cop,
757  VolumeField<Type>& result,
758  const bool secondOrder
759 ) const
760 {
761  mapInternalTgtToSrc(field, cop, result, secondOrder);
762 
763  const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
764 
765  forAll(AMIList, i)
766  {
767  label srcPatchi = srcPatchID_[i];
768  label tgtPatchi = tgtPatchID_[i];
769 
770  fvPatchField<Type>& srcField = result.boundaryFieldRef()[srcPatchi];
771  const fvPatchField<Type>& tgtField = field.boundaryField()[tgtPatchi];
772 
773 
774  // Clone and map (since rmap does not do general mapping)
775  tmp<fvPatchField<Type>> tnewSrc
776  (
778  (
779  tgtField,
780  srcField.patch(),
781  result(),
783  (
784  AMIList[i].singlePatchProc(),
785  (
786  AMIList[i].distributed()
787  ? AMIList[i].hasTgtMap() // pointer to map
788  : nullptr
789  ),
790  AMIList[i].srcAddress(),
791  AMIList[i].srcWeights()
792  )
793  )
794  );
795  // Transfer all mapped quantities (value and e.g. gradient) onto
796  // srcField. Value will get overwritten below
797  srcField.rmap(tnewSrc(), identity(srcField.size()));
798 
799  // Override value to account for CombineOp (could be dummy for
800  // plusEqOp)
801  mapAndOpTgtToSrc(AMIList[i], srcField, tgtField, cop);
802  }
803 
804  forAll(cuttingPatches_, i)
805  {
806  label patchi = cuttingPatches_[i];
807  fvPatchField<Type>& pf = result.boundaryFieldRef()[patchi];
808  pf == pf.patchInternalField();
809  }
810 }
811 
812 
813 template<class Type, class CombineOp>
816 (
817  const VolumeField<Type>& field,
818  const CombineOp& cop,
819  const bool secondOrder
820 ) const
821 {
822  const fvMesh& srcMesh = static_cast<const fvMesh&>(srcRegion_);
823 
824  const fvBoundaryMesh& srcBm = srcMesh.boundary();
825  const auto& tgtBfld = field.boundaryField();
826 
827  PtrList<fvPatchField<Type>> srcPatchFields(srcBm.size());
828 
829  // construct src boundary patch types as copy of 'field' boundary types
830  // note: this will provide place holders for fields with additional
831  // entries, but these values will need to be reset
832  forAll(srcPatchID_, i)
833  {
834  label srcPatchi = srcPatchID_[i];
835  label tgtPatchi = tgtPatchID_[i];
836 
837  if (!srcPatchFields.set(srcPatchi))
838  {
839  srcPatchFields.set
840  (
841  srcPatchi,
843  (
844  tgtBfld[tgtPatchi],
845  srcMesh.boundary()[srcPatchi],
848  (
849  labelList(srcMesh.boundary()[srcPatchi].size(), -1)
850  )
851  )
852  );
853  }
854  }
855 
856  // Any unset srcPatchFields become calculated
857  forAll(srcPatchFields, srcPatchi)
858  {
859  if (!srcPatchFields.set(srcPatchi))
860  {
861  // Note: use factory New method instead of direct generation of
862  // calculated so we keep constraints
863  srcPatchFields.set
864  (
865  srcPatchi,
867  (
869  srcMesh.boundary()[srcPatchi],
871  )
872  );
873  }
874  }
875 
876  auto tresult =
877  tmp<VolumeField<Type>>::New
878  (
879  IOobject
880  (
881  type() + ":interpolate(" + field.name() + ")",
882  srcMesh.time().timeName(),
883  srcMesh,
886  ),
887  srcMesh,
888  field.dimensions(),
889  Field<Type>(srcMesh.nCells(), Zero),
890  srcPatchFields
891  );
892 
893  mapTgtToSrc(field, cop, tresult.ref(), secondOrder);
894 
895  return tresult;
896 }
897 
898 
899 template<class Type, class CombineOp>
902 (
903  const tmp<VolumeField<Type>>& tfield,
904  const CombineOp& cop,
905  const bool secondOrder
906 ) const
907 {
908  return mapTgtToSrc(tfield(), cop, secondOrder);
909 }
910 
911 
912 template<class Type>
915 (
916  const VolumeField<Type>& field,
917  const bool secondOrder
918 ) const
919 {
920  return mapTgtToSrc(field, plusEqOp<Type>(), secondOrder);
921 }
922 
923 
924 template<class Type>
927 (
928  const tmp<VolumeField<Type>>& tfield,
929  const bool secondOrder
930 ) const
931 {
932  return mapTgtToSrc(tfield(), plusEqOp<Type>(), secondOrder);
933 }
934 
935 
936 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
type
Types of root.
Definition: Roots.H:52
rDeltaTY field()
FieldMapper with weighted mapping from (optionally remote) quantities.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
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.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:118
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
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...
Definition: labelLists.C:44
Generic templated field type.
Definition: Field.H:62
Calculate the gradient of the given field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:299
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
List< point > pointList
List of point.
Definition: pointList.H:32
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void mapSrcToTgt(const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
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.
void mapTgtToSrc(const UList< Type > &tgtFld, const CombineOp &cop, List< Type > &result) const
Map field from tgt to src mesh with defined operation.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127