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 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  //- Helper class for list
41  template<class Type>
43  {
44  public:
45  void operator()(List<Type>& x, const List<Type> y) const
46  {
47  if (y.size())
48  {
49  if (x.size())
50  {
51  label sz = x.size();
52  x.setSize(sz + y.size());
53  forAll(y, i)
54  {
55  x[sz++] = y[i];
56  }
57  }
58  else
59  {
60  x = y;
61  }
62  }
63  }
64  };
65 }
66 
67 
68 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69 
70 template<class Type>
71 void Foam::meshToMesh::add
72 (
73  UList<Type>& fld,
74  const label offset
75 ) const
76 {
77  forAll(fld, i)
78  {
79  fld[i] += offset;
80  }
81 }
82 
83 
84 template<class Type, class CombineOp>
86 (
87  const UList<Type>& srcField,
88  const CombineOp& cop,
89  List<Type>& result
90 ) const
91 {
92  if (result.size() != tgtToSrcCellAddr_.size())
93  {
95  << "Supplied field size is not equal to target mesh size" << nl
96  << " source mesh = " << srcToTgtCellAddr_.size() << nl
97  << " target mesh = " << tgtToSrcCellAddr_.size() << nl
98  << " supplied field = " << result.size()
99  << abort(FatalError);
100  }
101 
103 
104  if (singleMeshProc_ == -1)
105  {
106  const mapDistribute& map = srcMapPtr_();
107 
108  List<Type> work(srcField);
109  map.distribute(work);
110 
111  forAll(result, celli)
112  {
113  const labelList& srcAddress = tgtToSrcCellAddr_[celli];
114  const scalarList& srcWeight = tgtToSrcCellWght_[celli];
115 
116  if (srcAddress.size())
117  {
118 // result[celli] = Zero;
119  result[celli] *= (1.0 - sum(srcWeight));
120  forAll(srcAddress, i)
121  {
122  label srcI = srcAddress[i];
123  scalar w = srcWeight[i];
124  cbop(result[celli], celli, work[srcI], w);
125  }
126  }
127  }
128  }
129  else
130  {
131  forAll(result, celli)
132  {
133  const labelList& srcAddress = tgtToSrcCellAddr_[celli];
134  const scalarList& srcWeight = tgtToSrcCellWght_[celli];
135 
136  if (srcAddress.size())
137  {
138 // result[celli] = Zero;
139  result[celli] *= (1.0 - sum(srcWeight));
140  forAll(srcAddress, i)
141  {
142  label srcI = srcAddress[i];
143  scalar w = srcWeight[i];
144  cbop(result[celli], celli, srcField[srcI], w);
145  }
146  }
147  }
148  }
149 }
150 
151 
152 template<class Type, class CombineOp>
154 (
155  const UList<Type>& srcField,
156  const UList<typename outerProduct<vector, Type>::type>& srcGradField,
157  const CombineOp& cop,
158  List<Type>& result
159 ) const
160 {
161  if (result.size() != tgtToSrcCellAddr_.size())
162  {
164  << "Supplied field size is not equal to target mesh size" << nl
165  << " source mesh = " << srcToTgtCellAddr_.size() << nl
166  << " target mesh = " << tgtToSrcCellAddr_.size() << nl
167  << " supplied field = " << result.size()
168  << abort(FatalError);
169  }
170 
172 
173  if (singleMeshProc_ == -1)
174  {
175  if (returnReduce(tgtToSrcCellVec_.size(), sumOp<label>()) == 0)
176  {
177  // No correction vectors calculated. Fall back to first order.
178  mapSrcToTgt(srcField, cop, result);
179  return;
180  }
181 
182  const mapDistribute& map = srcMapPtr_();
183 
184  List<Type> work(srcField);
185  map.distribute(work);
186 
188  (
189  srcGradField
190  );
191  map.distribute(workGrad);
192 
193  forAll(result, cellI)
194  {
195  const labelList& srcAddress = tgtToSrcCellAddr_[cellI];
196  const scalarList& srcWeight = tgtToSrcCellWght_[cellI];
197  const pointList& srcVec = tgtToSrcCellVec_[cellI];
198 
199  if (srcAddress.size())
200  {
201  result[cellI] *= (1.0 - sum(srcWeight));
202  forAll(srcAddress, i)
203  {
204  label srcI = srcAddress[i];
205  scalar w = srcWeight[i];
206  const vector& v = srcVec[i];
207  const Type srcVal = work[srcI]+(workGrad[srcI]&v);
208  cbop(result[cellI], cellI, srcVal, w);
209  }
210  }
211  }
212  }
213  else
214  {
215  if (tgtToSrcCellVec_.empty())
216  {
217  // No correction vectors calculated. Fall back to first order.
218  mapSrcToTgt(srcField, cop, result);
219  return;
220  }
221 
222  forAll(result, cellI)
223  {
224  const labelList& srcAddress = tgtToSrcCellAddr_[cellI];
225  const scalarList& srcWeight = tgtToSrcCellWght_[cellI];
226  const pointList& srcVec = tgtToSrcCellVec_[cellI];
227 
228  if (srcAddress.size())
229  {
230  // Do non-conservative interpolation
231  result[cellI] *= (1.0 - sum(srcWeight));
232  forAll(srcAddress, i)
233  {
234  label srcI = srcAddress[i];
235  scalar w = srcWeight[i];
236  const vector& v = srcVec[i];
237  const Type srcVal = srcField[srcI]+(srcGradField[srcI]&v);
238  cbop(result[cellI], cellI, srcVal, w);
239  }
240  }
241  }
242  }
243 }
244 
245 
246 template<class Type, class CombineOp>
248 (
249  const Field<Type>& srcField,
250  const CombineOp& cop
251 ) const
252 {
253  tmp<Field<Type>> tresult
254  (
255  new Field<Type>
256  (
257  tgtToSrcCellAddr_.size(),
258  Zero
259  )
260  );
261 
262  mapSrcToTgt(srcField, cop, tresult.ref());
263 
264  return tresult;
265 }
266 
267 
268 template<class Type, class CombineOp>
270 (
271  const tmp<Field<Type>>& tsrcField,
272  const CombineOp& cop
273 ) const
274 {
275  return mapSrcToTgt(tsrcField(), cop);
276 }
277 
278 
279 template<class Type>
281 (
282  const Field<Type>& srcField
283 ) const
284 {
285  return mapSrcToTgt(srcField, plusEqOp<Type>());
286 }
287 
288 
289 template<class Type>
291 (
292  const tmp<Field<Type>>& tsrcField
293 ) const
294 {
295  return mapSrcToTgt(tsrcField());
296 }
297 
298 
299 template<class Type, class CombineOp>
301 (
302  const UList<Type>& tgtField,
303  const CombineOp& cop,
304  List<Type>& result
305 ) const
306 {
307  if (result.size() != srcToTgtCellAddr_.size())
308  {
310  << "Supplied field size is not equal to source mesh size" << nl
311  << " source mesh = " << srcToTgtCellAddr_.size() << nl
312  << " target mesh = " << tgtToSrcCellAddr_.size() << nl
313  << " supplied field = " << result.size()
314  << abort(FatalError);
315  }
316 
318 
319  if (singleMeshProc_ == -1)
320  {
321  const mapDistribute& map = tgtMapPtr_();
322 
323  List<Type> work(tgtField);
324  map.distribute(work);
325 
326  forAll(result, celli)
327  {
328  const labelList& tgtAddress = srcToTgtCellAddr_[celli];
329  const scalarList& tgtWeight = srcToTgtCellWght_[celli];
330 
331  if (tgtAddress.size())
332  {
333  result[celli] *= (1.0 - sum(tgtWeight));
334  forAll(tgtAddress, i)
335  {
336  label tgtI = tgtAddress[i];
337  scalar w = tgtWeight[i];
338  cbop(result[celli], celli, work[tgtI], w);
339  }
340  }
341  }
342  }
343  else
344  {
345  forAll(result, celli)
346  {
347  const labelList& tgtAddress = srcToTgtCellAddr_[celli];
348  const scalarList& tgtWeight = srcToTgtCellWght_[celli];
349 
350  if (tgtAddress.size())
351  {
352  result[celli] *= (1.0 - sum(tgtWeight));
353  forAll(tgtAddress, i)
354  {
355  label tgtI = tgtAddress[i];
356  scalar w = tgtWeight[i];
357  cbop(result[celli], celli, tgtField[tgtI], w);
358  }
359  }
360  }
361  }
362 }
363 
364 
365 template<class Type, class CombineOp>
367 (
368  const UList<Type>& tgtField,
369  const UList<typename outerProduct<vector, Type>::type>& tgtGradField,
370  const CombineOp& cop,
371  List<Type>& result
372 ) const
373 {
374  if (result.size() != srcToTgtCellAddr_.size())
375  {
377  << "Supplied field size is not equal to source mesh size" << nl
378  << " source mesh = " << srcToTgtCellAddr_.size() << nl
379  << " target mesh = " << tgtToSrcCellAddr_.size() << nl
380  << " supplied field = " << result.size()
381  << abort(FatalError);
382  }
383 
385 
386  if (singleMeshProc_ == -1)
387  {
388  if (returnReduce(srcToTgtCellVec_.size(), sumOp<label>()) == 0)
389  {
390  // No correction vectors calculated. Fall back to first order.
391  mapTgtToSrc(tgtField, cop, result);
392  return;
393  }
394 
395  const mapDistribute& map = tgtMapPtr_();
396 
397  List<Type> work(tgtField);
398  map.distribute(work);
399 
401  (
402  tgtGradField
403  );
404  map.distribute(workGrad);
405 
406  forAll(result, cellI)
407  {
408  const labelList& tgtAddress = srcToTgtCellAddr_[cellI];
409  const scalarList& tgtWeight = srcToTgtCellWght_[cellI];
410  const pointList& tgtVec = srcToTgtCellVec_[cellI];
411 
412  if (tgtAddress.size())
413  {
414  result[cellI] *= (1.0 - sum(tgtWeight));
415  forAll(tgtAddress, i)
416  {
417  label tgtI = tgtAddress[i];
418  scalar w = tgtWeight[i];
419  const vector& v = tgtVec[i];
420  const Type tgtVal = work[tgtI]+(workGrad[tgtI]&v);
421  cbop(result[cellI], cellI, tgtVal, w);
422  }
423  }
424  }
425  }
426  else
427  {
428  forAll(result, cellI)
429  {
430  const labelList& tgtAddress = srcToTgtCellAddr_[cellI];
431  const scalarList& tgtWeight = srcToTgtCellWght_[cellI];
432  const pointList& tgtVec = srcToTgtCellVec_[cellI];
433 
434  if (tgtAddress.size())
435  {
436  result[cellI] *= (1.0 - sum(tgtWeight));
437  forAll(tgtAddress, i)
438  {
439  label tgtI = tgtAddress[i];
440  scalar w = tgtWeight[i];
441  const vector& v = tgtVec[i];
442  const Type tgtVal = tgtField[tgtI]+(tgtGradField[tgtI]&v);
443  cbop(result[cellI], cellI, tgtVal, w);
444  }
445  }
446  }
447  }
448 }
449 
450 
451 template<class Type, class CombineOp>
453 (
454  const Field<Type>& tgtField,
455  const CombineOp& cop
456 ) const
457 {
458  tmp<Field<Type>> tresult
459  (
460  new Field<Type>
461  (
462  srcToTgtCellAddr_.size(),
463  Zero
464  )
465  );
466 
467  mapTgtToSrc(tgtField, cop, tresult.ref());
468 
469  return tresult;
470 }
471 
472 
473 template<class Type, class CombineOp>
475 (
476  const tmp<Field<Type>>& ttgtField,
477  const CombineOp& cop
478 ) const
479 {
480  return mapTgtToSrc(ttgtField(), cop);
481 }
482 
483 
484 template<class Type>
486 (
487  const Field<Type>& tgtField
488 ) const
489 {
490  return mapTgtToSrc(tgtField, plusEqOp<Type>());
491 }
492 
493 
494 template<class Type>
496 (
497  const tmp<Field<Type>>& ttgtField
498 ) const
499 {
500  return mapTgtToSrc(ttgtField(), plusEqOp<Type>());
501 }
502 
503 
504 template<class Type, class CombineOp>
505 void Foam::meshToMesh::mapInternalSrcToTgt
506 (
508  const CombineOp& cop,
510  const bool secondOrder
511 ) const
512 {
513  if (secondOrder && returnReduce(tgtToSrcCellVec_.size(), sumOp<label>()))
514  {
515  mapSrcToTgt
516  (
517  field,
518  fvc::grad(field)().primitiveField(),
519  cop,
520  result.primitiveFieldRef()
521  );
522  }
523  else
524  {
525  mapSrcToTgt(field, cop, result.primitiveFieldRef());
526  }
527 }
528 
529 
530 template<class Type, class CombineOp>
531 void Foam::meshToMesh::mapAndOpSrcToTgt
532 (
533  const AMIPatchToPatchInterpolation& AMI,
534  const Field<Type>& srcField,
535  Field<Type>& tgtField,
536  const CombineOp& cop
537 ) const
538 {
539  tgtField = Type(Zero);
540 
541  AMI.interpolateToTarget
542  (
543  srcField,
544  multiplyWeightedOp<Type, CombineOp>(cop),
545  tgtField,
547  );
548 }
549 
550 
551 template<class Type, class CombineOp>
553 (
555  const CombineOp& cop,
557  const bool secondOrder
558 ) const
559 {
560  mapInternalSrcToTgt(field, cop, result, secondOrder);
561 
562  const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
563 
565  Boundary& resultBf = result.boundaryFieldRef();
566 
567  forAll(AMIList, i)
568  {
569  label srcPatchi = srcPatchID_[i];
570  label tgtPatchi = tgtPatchID_[i];
571 
572  const fvPatchField<Type>& srcField = field.boundaryField()[srcPatchi];
573  fvPatchField<Type>& tgtField = resultBf[tgtPatchi];
574 
575  // Clone and map (since rmap does not do general mapping)
576  tmp<fvPatchField<Type>> tnewTgt
577  (
579  (
580  srcField,
581  tgtField.patch(),
582  result(),
584  (
585  AMIList[i].singlePatchProc(),
586  (
587  AMIList[i].singlePatchProc() == -1
588  ? &AMIList[i].srcMap()
589  : nullptr
590  ),
591  AMIList[i].tgtAddress(),
592  AMIList[i].tgtWeights()
593  )
594  )
595  );
596 
597  // Transfer all mapped quantities (value and e.g. gradient) onto
598  // tgtField. Value will get overwritten below.
599  tgtField.rmap(tnewTgt(), identity(tgtField.size()));
600 
601  // Override value to account for CombineOp (note: is dummy template
602  // specialisation for plusEqOp)
603  mapAndOpSrcToTgt(AMIList[i], srcField, tgtField, cop);
604  }
605 
606  forAll(cuttingPatches_, i)
607  {
608  label patchi = cuttingPatches_[i];
609  fvPatchField<Type>& pf = resultBf[patchi];
610  pf == pf.patchInternalField();
611  }
612 }
613 
614 
615 template<class Type, class CombineOp>
618 (
620  const CombineOp& cop,
621  const bool secondOrder
622 ) const
623 {
625 
626  const fvMesh& tgtMesh = static_cast<const fvMesh&>(tgtRegion_);
627 
628  const fvBoundaryMesh& tgtBm = tgtMesh.boundary();
629  const typename fieldType::Boundary& srcBfld =
630  field.boundaryField();
631 
632  PtrList<fvPatchField<Type>> tgtPatchFields(tgtBm.size());
633 
634  // construct tgt boundary patch types as copy of 'field' boundary types
635  // note: this will provide place holders for fields with additional
636  // entries, but these values will need to be reset
637  forAll(tgtPatchID_, i)
638  {
639  label srcPatchi = srcPatchID_[i];
640  label tgtPatchi = tgtPatchID_[i];
641 
642  if (!tgtPatchFields.set(tgtPatchi))
643  {
644  tgtPatchFields.set
645  (
646  tgtPatchi,
648  (
649  srcBfld[srcPatchi],
650  tgtMesh.boundary()[tgtPatchi],
653  (
654  labelList(tgtMesh.boundary()[tgtPatchi].size(), -1)
655  )
656  )
657  );
658  }
659  }
660 
661  // Any unset tgtPatchFields become calculated
662  forAll(tgtPatchFields, tgtPatchi)
663  {
664  if (!tgtPatchFields.set(tgtPatchi))
665  {
666  // Note: use factory New method instead of direct generation of
667  // calculated so we keep constraints
668  tgtPatchFields.set
669  (
670  tgtPatchi,
672  (
674  tgtMesh.boundary()[tgtPatchi],
676  )
677  );
678  }
679  }
680 
681  tmp<fieldType> tresult
682  (
683  new fieldType
684  (
685  IOobject
686  (
687  type() + ":interpolate(" + field.name() + ")",
688  tgtMesh.time().timeName(),
689  tgtMesh,
692  ),
693  tgtMesh,
694  field.dimensions(),
695  Field<Type>(tgtMesh.nCells(), Zero),
696  tgtPatchFields
697  )
698  );
699 
700  mapSrcToTgt(field, cop, tresult.ref(), secondOrder);
701 
702  return tresult;
703 }
704 
705 
706 template<class Type, class CombineOp>
709 (
711  const CombineOp& cop,
712  const bool secondOrder
713 ) const
714 {
715  return mapSrcToTgt(tfield(), cop, secondOrder);
716 }
717 
718 
719 template<class Type>
722 (
724  const bool secondOrder
725 ) const
726 {
727  return mapSrcToTgt(field, plusEqOp<Type>(), secondOrder);
728 }
729 
730 
731 template<class Type>
734 (
736  const bool secondOrder
737 ) const
738 {
739  return mapSrcToTgt(tfield(), plusEqOp<Type>(), secondOrder);
740 }
741 
742 
743 template<class Type, class CombineOp>
744 void Foam::meshToMesh::mapInternalTgtToSrc
745 (
747  const CombineOp& cop,
749  const bool secondOrder
750 ) const
751 {
752  if (secondOrder && returnReduce(srcToTgtCellVec_.size(), sumOp<label>()))
753  {
754  mapTgtToSrc
755  (
756  field,
757  fvc::grad(field)().primitiveField(),
758  cop,
759  result.primitiveFieldRef()
760  );
761  }
762  else
763  {
764  mapTgtToSrc(field, cop, result.primitiveFieldRef());
765  }
766 }
767 
768 
769 template<class Type, class CombineOp>
770 void Foam::meshToMesh::mapAndOpTgtToSrc
771 (
772  const AMIPatchToPatchInterpolation& AMI,
773  Field<Type>& srcField,
774  const Field<Type>& tgtField,
775  const CombineOp& cop
776 ) const
777 {
778  srcField = Type(Zero);
779 
780  AMI.interpolateToSource
781  (
782  tgtField,
783  multiplyWeightedOp<Type, CombineOp>(cop),
784  srcField,
786  );
787 }
788 
789 
790 template<class Type, class CombineOp>
792 (
794  const CombineOp& cop,
796  const bool secondOrder
797 ) const
798 {
799  mapInternalTgtToSrc(field, cop, result, secondOrder);
800 
801  const PtrList<AMIPatchToPatchInterpolation>& AMIList = patchAMIs();
802 
803  forAll(AMIList, i)
804  {
805  label srcPatchi = srcPatchID_[i];
806  label tgtPatchi = tgtPatchID_[i];
807 
808  fvPatchField<Type>& srcField = result.boundaryFieldRef()[srcPatchi];
809  const fvPatchField<Type>& tgtField = field.boundaryField()[tgtPatchi];
810 
811 
812  // Clone and map (since rmap does not do general mapping)
813  tmp<fvPatchField<Type>> tnewSrc
814  (
816  (
817  tgtField,
818  srcField.patch(),
819  result(),
821  (
822  AMIList[i].singlePatchProc(),
823  (
824  AMIList[i].singlePatchProc() == -1
825  ? &AMIList[i].tgtMap()
826  : nullptr
827  ),
828  AMIList[i].srcAddress(),
829  AMIList[i].srcWeights()
830  )
831  )
832  );
833  // Transfer all mapped quantities (value and e.g. gradient) onto
834  // srcField. Value will get overwritten below
835  srcField.rmap(tnewSrc(), identity(srcField.size()));
836 
837  // Override value to account for CombineOp (could be dummy for
838  // plusEqOp)
839  mapAndOpTgtToSrc(AMIList[i], srcField, tgtField, cop);
840  }
841 
842  forAll(cuttingPatches_, i)
843  {
844  label patchi = cuttingPatches_[i];
845  fvPatchField<Type>& pf = result.boundaryFieldRef()[patchi];
846  pf == pf.patchInternalField();
847  }
848 }
849 
850 
851 template<class Type, class CombineOp>
854 (
856  const CombineOp& cop,
857  const bool secondOrder
858 ) const
859 {
861 
862  const fvMesh& srcMesh = static_cast<const fvMesh&>(srcRegion_);
863 
864  const fvBoundaryMesh& srcBm = srcMesh.boundary();
865  const typename fieldType::Boundary& tgtBfld =
866  field.boundaryField();
867 
868  PtrList<fvPatchField<Type>> srcPatchFields(srcBm.size());
869 
870  // construct src boundary patch types as copy of 'field' boundary types
871  // note: this will provide place holders for fields with additional
872  // entries, but these values will need to be reset
873  forAll(srcPatchID_, i)
874  {
875  label srcPatchi = srcPatchID_[i];
876  label tgtPatchi = tgtPatchID_[i];
877 
878  if (!srcPatchFields.set(tgtPatchi))
879  {
880  srcPatchFields.set
881  (
882  srcPatchi,
884  (
885  tgtBfld[srcPatchi],
886  srcMesh.boundary()[tgtPatchi],
889  (
890  labelList(srcMesh.boundary()[srcPatchi].size(), -1)
891  )
892  )
893  );
894  }
895  }
896 
897  // Any unset srcPatchFields become calculated
898  forAll(srcPatchFields, srcPatchi)
899  {
900  if (!srcPatchFields.set(srcPatchi))
901  {
902  // Note: use factory New method instead of direct generation of
903  // calculated so we keep constraints
904  srcPatchFields.set
905  (
906  srcPatchi,
908  (
910  srcMesh.boundary()[srcPatchi],
912  )
913  );
914  }
915  }
916 
917  tmp<fieldType> tresult
918  (
919  new fieldType
920  (
921  IOobject
922  (
923  type() + ":interpolate(" + field.name() + ")",
924  srcMesh.time().timeName(),
925  srcMesh,
928  ),
929  srcMesh,
930  field.dimensions(),
931  Field<Type>(srcMesh.nCells(), Zero),
932  srcPatchFields
933  )
934  );
935 
936  mapTgtToSrc(field, cop, tresult.ref(), secondOrder);
937 
938  return tresult;
939 }
940 
941 
942 template<class Type, class CombineOp>
945 (
947  const CombineOp& cop,
948  const bool secondOrder
949 ) const
950 {
951  return mapTgtToSrc(tfield(), cop, secondOrder);
952 }
953 
954 
955 template<class Type>
958 (
960  const bool secondOrder
961 ) const
962 {
963  return mapTgtToSrc(field, plusEqOp<Type>(), secondOrder);
964 }
965 
966 
967 template<class Type>
970 (
972  const bool secondOrder
973 ) const
974 {
975  return mapTgtToSrc(tfield(), plusEqOp<Type>(), secondOrder);
976 }
977 
978 
979 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
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:301
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:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::ListPlusEqOp
Helper class for list.
Definition: meshToMeshTemplates.C:42
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:764
calculatedFvPatchField.H
Foam::ListPlusEqOp::operator()
void operator()(List< Type > &x, const List< Type > y) const
Definition: meshToMeshTemplates.C:45
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:86
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
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:225
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:65
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:303
Foam::FatalError
error FatalError
Foam::AMIPatchToPatchInterpolation
AMIInterpolation< PrimitivePatch< face, SubList, const pointField & >, PrimitivePatch< face, SubList, const pointField & > > AMIPatchToPatchInterpolation
Definition: AMIPatchToPatchInterpolation.H:45
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
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
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::multiplyWeightedOp
Definition: ops.H:250
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
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::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
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::directFvPatchFieldMapper
direct fvPatchFieldMapper
Definition: directFvPatchFieldMapper.H:48
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
distributedWeightedFvPatchFieldMapper.H
Foam::plusEqOp
Definition: ops.H:72
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
Foam::UList< Type >::null
static const UList< Type > & null()
Return a UList reference to a nullObject.
Definition: UListI.H:54
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::distributedWeightedFvPatchFieldMapper
FieldMapper with weighted mapping from (optionally remote) quantities.
Definition: distributedWeightedFvPatchFieldMapper.H:49
y
scalar y
Definition: LISASMDCalcMethod1.H:14