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 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 
71 
72  if (singleMeshProc_ == -1)
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 
140 
141  if (singleMeshProc_ == -1)
142  {
143  if (returnReduce(tgtToSrcCellVec_.size(), sumOp<label>()) == 0)
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());
231 
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 
286 
287  if (singleMeshProc_ == -1)
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 
353 
354  if (singleMeshProc_ == -1)
355  {
356  if (returnReduce(srcToTgtCellVec_.size(), sumOp<label>()) == 0)
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());
436 
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 (
476  const CombineOp& cop,
478  const bool secondOrder
479 ) const
480 {
481  if (secondOrder && returnReduce(tgtToSrcCellVec_.size(), sumOp<label>()))
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 (
523  const CombineOp& cop,
525  const bool secondOrder
526 ) const
527 {
528  mapInternalSrcToTgt(field, cop, result, secondOrder);
529 
530  const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
531 
533  Boundary& resultBf = result.boundaryFieldRef();
534 
535  forAll(AMIList, i)
536  {
537  label srcPatchi = srcPatchID_[i];
538  label tgtPatchi = tgtPatchID_[i];
539 
540  const fvPatchField<Type>& srcField = field.boundaryField()[srcPatchi];
541  fvPatchField<Type>& tgtField = resultBf[tgtPatchi];
542 
543  // Clone and map (since rmap does not do general mapping)
544  tmp<fvPatchField<Type>> tnewTgt
545  (
547  (
548  srcField,
549  tgtField.patch(),
550  result(),
552  (
553  AMIList[i].singlePatchProc(),
554  (
555  AMIList[i].singlePatchProc() == -1
556  ? &AMIList[i].srcMap()
557  : nullptr
558  ),
559  AMIList[i].tgtAddress(),
560  AMIList[i].tgtWeights()
561  )
562  )
563  );
564 
565  // Transfer all mapped quantities (value and e.g. gradient) onto
566  // tgtField. Value will get overwritten below.
567  tgtField.rmap(tnewTgt(), identity(tgtField.size()));
568 
569  // Override value to account for CombineOp (note: is dummy template
570  // specialisation for plusEqOp)
571  mapAndOpSrcToTgt(AMIList[i], srcField, tgtField, cop);
572  }
573 
574  forAll(cuttingPatches_, i)
575  {
576  label patchi = cuttingPatches_[i];
577  fvPatchField<Type>& pf = resultBf[patchi];
578  pf == pf.patchInternalField();
579  }
580 }
581 
582 
583 template<class Type, class CombineOp>
586 (
588  const CombineOp& cop,
589  const bool secondOrder
590 ) const
591 {
593 
594  const fvMesh& tgtMesh = static_cast<const fvMesh&>(tgtRegion_);
595 
596  const fvBoundaryMesh& tgtBm = tgtMesh.boundary();
597  const typename fieldType::Boundary& srcBfld =
598  field.boundaryField();
599 
600  PtrList<fvPatchField<Type>> tgtPatchFields(tgtBm.size());
601 
602  // construct tgt boundary patch types as copy of 'field' boundary types
603  // note: this will provide place holders for fields with additional
604  // entries, but these values will need to be reset
605  forAll(tgtPatchID_, i)
606  {
607  label srcPatchi = srcPatchID_[i];
608  label tgtPatchi = tgtPatchID_[i];
609 
610  if (!tgtPatchFields.set(tgtPatchi))
611  {
612  tgtPatchFields.set
613  (
614  tgtPatchi,
616  (
617  srcBfld[srcPatchi],
618  tgtMesh.boundary()[tgtPatchi],
621  (
622  labelList(tgtMesh.boundary()[tgtPatchi].size(), -1)
623  )
624  )
625  );
626  }
627  }
628 
629  // Any unset tgtPatchFields become calculated
630  forAll(tgtPatchFields, tgtPatchi)
631  {
632  if (!tgtPatchFields.set(tgtPatchi))
633  {
634  // Note: use factory New method instead of direct generation of
635  // calculated so we keep constraints
636  tgtPatchFields.set
637  (
638  tgtPatchi,
640  (
642  tgtMesh.boundary()[tgtPatchi],
644  )
645  );
646  }
647  }
648 
649  tmp<fieldType> tresult
650  (
651  new fieldType
652  (
653  IOobject
654  (
655  type() + ":interpolate(" + field.name() + ")",
656  tgtMesh.time().timeName(),
657  tgtMesh,
660  ),
661  tgtMesh,
662  field.dimensions(),
663  Field<Type>(tgtMesh.nCells(), Zero),
664  tgtPatchFields
665  )
666  );
667 
668  mapSrcToTgt(field, cop, tresult.ref(), secondOrder);
669 
670  return tresult;
671 }
672 
673 
674 template<class Type, class CombineOp>
677 (
679  const CombineOp& cop,
680  const bool secondOrder
681 ) const
682 {
683  return mapSrcToTgt(tfield(), cop, secondOrder);
684 }
685 
686 
687 template<class Type>
690 (
692  const bool secondOrder
693 ) const
694 {
695  return mapSrcToTgt(field, plusEqOp<Type>(), secondOrder);
696 }
697 
698 
699 template<class Type>
702 (
704  const bool secondOrder
705 ) const
706 {
707  return mapSrcToTgt(tfield(), plusEqOp<Type>(), secondOrder);
708 }
709 
710 
711 template<class Type, class CombineOp>
712 void Foam::meshToMesh::mapInternalTgtToSrc
713 (
715  const CombineOp& cop,
717  const bool secondOrder
718 ) const
719 {
720  if (secondOrder && returnReduce(srcToTgtCellVec_.size(), sumOp<label>()))
721  {
722  mapTgtToSrc
723  (
724  field,
725  fvc::grad(field)().primitiveField(),
726  cop,
727  result.primitiveFieldRef()
728  );
729  }
730  else
731  {
732  mapTgtToSrc(field, cop, result.primitiveFieldRef());
733  }
734 }
735 
736 
737 template<class Type, class CombineOp>
738 void Foam::meshToMesh::mapAndOpTgtToSrc
739 (
740  const AMIPatchToPatchInterpolation& AMI,
741  Field<Type>& srcField,
742  const Field<Type>& tgtField,
743  const CombineOp& cop
744 ) const
745 {
746  srcField = Type(Zero);
747 
748  AMI.interpolateToSource
749  (
750  tgtField,
751  multiplyWeightedOp<Type, CombineOp>(cop),
752  srcField,
754  );
755 }
756 
757 
758 template<class Type, class CombineOp>
760 (
762  const CombineOp& cop,
764  const bool secondOrder
765 ) const
766 {
767  mapInternalTgtToSrc(field, cop, result, secondOrder);
768 
769  const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
770 
771  forAll(AMIList, i)
772  {
773  label srcPatchi = srcPatchID_[i];
774  label tgtPatchi = tgtPatchID_[i];
775 
776  fvPatchField<Type>& srcField = result.boundaryFieldRef()[srcPatchi];
777  const fvPatchField<Type>& tgtField = field.boundaryField()[tgtPatchi];
778 
779 
780  // Clone and map (since rmap does not do general mapping)
781  tmp<fvPatchField<Type>> tnewSrc
782  (
784  (
785  tgtField,
786  srcField.patch(),
787  result(),
789  (
790  AMIList[i].singlePatchProc(),
791  (
792  AMIList[i].singlePatchProc() == -1
793  ? &AMIList[i].tgtMap()
794  : nullptr
795  ),
796  AMIList[i].srcAddress(),
797  AMIList[i].srcWeights()
798  )
799  )
800  );
801  // Transfer all mapped quantities (value and e.g. gradient) onto
802  // srcField. Value will get overwritten below
803  srcField.rmap(tnewSrc(), identity(srcField.size()));
804 
805  // Override value to account for CombineOp (could be dummy for
806  // plusEqOp)
807  mapAndOpTgtToSrc(AMIList[i], srcField, tgtField, cop);
808  }
809 
810  forAll(cuttingPatches_, i)
811  {
812  label patchi = cuttingPatches_[i];
813  fvPatchField<Type>& pf = result.boundaryFieldRef()[patchi];
814  pf == pf.patchInternalField();
815  }
816 }
817 
818 
819 template<class Type, class CombineOp>
822 (
824  const CombineOp& cop,
825  const bool secondOrder
826 ) const
827 {
829 
830  const fvMesh& srcMesh = static_cast<const fvMesh&>(srcRegion_);
831 
832  const fvBoundaryMesh& srcBm = srcMesh.boundary();
833  const typename fieldType::Boundary& tgtBfld =
834  field.boundaryField();
835 
836  PtrList<fvPatchField<Type>> srcPatchFields(srcBm.size());
837 
838  // construct src boundary patch types as copy of 'field' boundary types
839  // note: this will provide place holders for fields with additional
840  // entries, but these values will need to be reset
841  forAll(srcPatchID_, i)
842  {
843  label srcPatchi = srcPatchID_[i];
844  label tgtPatchi = tgtPatchID_[i];
845 
846  if (!srcPatchFields.set(tgtPatchi))
847  {
848  srcPatchFields.set
849  (
850  srcPatchi,
852  (
853  tgtBfld[srcPatchi],
854  srcMesh.boundary()[tgtPatchi],
857  (
858  labelList(srcMesh.boundary()[srcPatchi].size(), -1)
859  )
860  )
861  );
862  }
863  }
864 
865  // Any unset srcPatchFields become calculated
866  forAll(srcPatchFields, srcPatchi)
867  {
868  if (!srcPatchFields.set(srcPatchi))
869  {
870  // Note: use factory New method instead of direct generation of
871  // calculated so we keep constraints
872  srcPatchFields.set
873  (
874  srcPatchi,
876  (
878  srcMesh.boundary()[srcPatchi],
880  )
881  );
882  }
883  }
884 
885  tmp<fieldType> tresult
886  (
887  new fieldType
888  (
889  IOobject
890  (
891  type() + ":interpolate(" + field.name() + ")",
892  srcMesh.time().timeName(),
893  srcMesh,
896  ),
897  srcMesh,
898  field.dimensions(),
899  Field<Type>(srcMesh.nCells(), Zero),
900  srcPatchFields
901  )
902  );
903 
904  mapTgtToSrc(field, cop, tresult.ref(), secondOrder);
905 
906  return tresult;
907 }
908 
909 
910 template<class Type, class CombineOp>
913 (
915  const CombineOp& cop,
916  const bool secondOrder
917 ) const
918 {
919  return mapTgtToSrc(tfield(), cop, secondOrder);
920 }
921 
922 
923 template<class Type>
926 (
928  const bool secondOrder
929 ) const
930 {
931  return mapTgtToSrc(field, plusEqOp<Type>(), secondOrder);
932 }
933 
934 
935 template<class Type>
938 (
940  const bool secondOrder
941 ) const
942 {
943  return mapTgtToSrc(tfield(), plusEqOp<Type>(), secondOrder);
944 }
945 
946 
947 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::meshToMesh::mapTgtToSrc
void mapTgtToSrc(const UList< Type > &tgtFld, const CombineOp &cop, List< Type > &result) const
Map field from tgt to src mesh with defined operation.
Definition: meshToMeshTemplates.C:269
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::outerProduct::type
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:114
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
calculatedFvPatchField.H
Foam::DimensionedField::null
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
Definition: DimensionedFieldI.H:33
Foam::meshToMesh::mapSrcToTgt
void mapSrcToTgt(const UList< Type > &srcFld, const CombineOp &cop, List< Type > &result) const
Map field from src to tgt mesh with defined operation.
Definition: meshToMeshTemplates.C:54
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
directFvPatchFieldMapper.H
fld
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))
Definition: gmvOutputLagrangian.H:23
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:152
Foam::fvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:311
Foam::FatalError
error FatalError
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam::calculatedFvPatchField
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Definition: calculatedFvPatchField.H:66
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
Foam::multiplyWeightedOp
Definition: ops.H:250
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::Vector< scalar >
Foam::List< Type >
fvcGrad.H
Calculate the gradient of the given field.
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::AMIPatchToPatchInterpolation
AMIInterpolation AMIPatchToPatchInterpolation
Definition: AMIPatchToPatchInterpolation.H:38
Foam::UList< Type >
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::directFvPatchFieldMapper
direct fvPatchFieldMapper
Definition: directFvPatchFieldMapper.H:48
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
distributedWeightedFvPatchFieldMapper.H
Foam::plusEqOp
Definition: ops.H:72
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::UList< Type >::null
static const UList< Type > & null()
Return a UList reference to a nullObject.
Definition: UListI.H:53
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::distributedWeightedFvPatchFieldMapper
FieldMapper with weighted mapping from (optionally remote) quantities.
Definition: distributedWeightedFvPatchFieldMapper.H:49