mappedPatchFieldBase.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2021 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 "mappedPatchFieldBase.H"
30 #include "mappedPatchBase.H"
31 #include "interpolationCell.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const dictionary& dict,
39  const bool mandatory
40 )
41 {
42  if (mandatory)
43  {
44  return dict.get<Type>("average");
45  }
46 
47  return Zero;
48 }
49 
50 
51 template<class Type>
52 template<class T>
54 (
55  const objectRegistry& obr,
56  const word& region,
57  const word& patch,
58  const label myComm,
59  const labelListList& procToMap,
60  const word& fieldName,
61  const Field<T>& fld
62 ) const
63 {
64  // Store my data onto database
65 
66  const auto& procIDs = UPstream::procID(myComm);
67 
68  forAll(procToMap, ranki)
69  {
70  const labelList& map = procToMap[ranki];
71  const label proci = procIDs[ranki];
72 
73  if (map.size())
74  {
75  const Field<T> subFld(fld, map);
76 
77  auto& subObr = const_cast<objectRegistry&>
78  (
79  mappedPatchBase::subRegistry
80  (
81  obr,
82  mapper_.sendPath(proci)
83  / region
84  / patch
85  )
86  );
87 
89  {
90  Pout<< "*** STORING :"
91  << " field:" << fieldName
92  << " values:" << flatOutput(subFld)
93  << " as:" << subObr.objectPath() << endl;
94  }
95 
96  mappedPatchBase::storeField(subObr, fieldName, subFld);
97  }
98  }
99 }
100 
101 
102 template<class Type>
103 template<class T>
105 (
106  const bool allowUnset,
107  const objectRegistry& obr,
108  const word& region,
109  const word& patch,
110  const label myComm,
111  const labelListList& procToMap,
112  const word& fieldName,
113  Field<T>& fld
114 ) const
115 {
116  const auto& procIDs = UPstream::procID(myComm);
117 
118  bool ok = true;
119 
120  forAll(procToMap, ranki)
121  {
122  const labelList& map = procToMap[ranki];
123  const label proci = procIDs[ranki];
124 
125  if (map.size())
126  {
127  auto& subObr = const_cast<objectRegistry&>
128  (
129  mappedPatchBase::subRegistry
130  (
131  obr,
132  mapper_.receivePath(proci)
133  / region
134  / patch
135  )
136  );
137 
138  const IOField<T>* subFldPtr = subObr.getObjectPtr<IOField<T>>
139  (
140  fieldName
141  );
142  if (subFldPtr)
143  {
144  if (subFldPtr->size() != map.size())
145  {
146  // This is the dummy value inserted at start-up since the
147  // map is always non-zero size (checked above)
148  //Pout<< "*** RETRIEVED DUMMY :"
149  // << " field:" << fieldName
150  // << " subFldPtr:" << subFldPtr->size()
151  // << " map:" << map.size() << endl;
152 
153  ok = false;
154  }
155  else
156  {
157  UIndirectList<T>(fld, map) = *subFldPtr;
158 
160  {
161  Pout<< "*** RETRIEVED :"
162  << " field:" << fieldName
163  << " values:" << flatOutput(fld)
164  << " from:" << subObr.objectPath() << endl;
165  }
166  }
167  }
168  else if (allowUnset)
169  {
171  {
172  WarningInFunction << "Not found"
173  << " field:" << fieldName
174  << " in:" << subObr.objectPath() << endl;
175  }
176 
177  // Store dummy value so the database has something on it.
178  // Note that size 0 should never occur naturally so we can
179  // detect it if necessary.
180  const Field<T> dummyFld(0);
181 
182  mappedPatchBase::storeField(subObr, fieldName, dummyFld);
183 
184  ok = false;
185  }
186  else
187  {
188  // Not found. Make it fail
189  (void)subObr.lookupObject<IOField<T>>(fieldName);
190  ok = false;
191  }
192  }
193  }
194  return ok;
195 }
196 
197 
198 template<class Type>
199 template<class T>
201 (
202  const objectRegistry& obr,
203  const word& region,
204  const word& patch,
205  const labelListList& map,
206  const word& fieldName,
207  const Field<T>& fld
208 ) const
209 {
210  // Old code. Likely not quite correct...
211 
212  // Store my data onto database
213  const label nProcs = Pstream::nProcs(0); // comm_
214 
215  for (label domain = 0; domain < nProcs; domain++)
216  {
217  const labelList& constructMap = map[domain];
218  if (constructMap.size())
219  {
220  auto& subObr = const_cast<objectRegistry&>
221  (
222  mappedPatchBase::subRegistry
223  (
224  obr,
225  mapper_.receivePath(domain)
226  / region
227  / patch
228  )
229  );
230 
231  const Field<T> receiveFld(fld, constructMap);
232 
234  {
235  Pout<< "*** STORING INITIAL :"
236  << " field:" << fieldName << " values:"
237  << flatOutput(receiveFld)
238  << " from:" << flatOutput(fld)
239  << " constructMap:" << flatOutput(constructMap)
240  << " as:" << subObr.objectPath() << endl;
241  }
242 
243  mappedPatchBase::storeField(subObr, fieldName, receiveFld);
244  }
245  }
246 }
247 
248 
249 template<class Type>
250 template<class T>
252 (
253  const word& fieldName,
254  const label myComm,
255  const labelListList& subMap,
256  const label constructSize,
257  const labelListList& constructMap,
258  const labelListList& address,
259  const scalarListList& weights,
260  Field<T>& fld
261 ) const
262 {
263  storeField
264  (
265  patchField_.internalField().time(),
266  patchField_.patch().boundaryMesh().mesh().name(),
267  patchField_.patch().name(),
268  myComm,
269  subMap,
270  fieldName,
271  fld
272  );
273 
274  Field<T> work(constructSize);
275  const bool ok = retrieveField
276  (
277  true, // allow unset
278  patchField_.internalField().time(),
279  mapper_.sampleRegion(),
280  mapper_.samplePatch(),
281  myComm,
282  constructMap,
283  fieldName,
284  work
285  );
286 
287  if (ok)
288  {
289  // Do interpolation
290 
291  fld.setSize(address.size());
292  fld = Zero;
293 
294  const plusEqOp<T> cop;
295  const multiplyWeightedOp<T, plusEqOp<T>> mop(cop);
296 
297  forAll(address, facei)
298  {
299  const labelList& slots = address[facei];
300  const scalarList& w = weights[facei];
301 
302  forAll(slots, i)
303  {
304  mop(fld[facei], facei, work[slots[i]], w[i]);
305  }
306  }
307  }
308  else
309  {
310  // Leave fld intact
311  }
312 
313  return ok;
314 }
315 
316 
317 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
318 
319 template<class Type>
321 (
322  const mappedPatchBase& mapper,
323  const fvPatchField<Type>& patchField,
324  const word& fieldName,
325  const bool setAverage,
326  const Type average,
327  const word& interpolationScheme
328 )
329 :
330  mapper_(mapper),
331  patchField_(patchField),
332  fieldName_(fieldName),
333  setAverage_(setAverage),
334  average_(average),
335  interpolationScheme_(interpolationScheme)
336 {}
337 
338 
339 template<class Type>
341 (
342  const mappedPatchBase& mapper,
343  const fvPatchField<Type>& patchField,
344  const dictionary& dict
345 )
346 :
347  mapper_(mapper),
348  patchField_(patchField),
349  fieldName_
350  (
351  dict.template getOrDefault<word>
352  (
353  "field",
354  patchField_.internalField().name()
355  )
356  ),
357  setAverage_(dict.getOrDefault("setAverage", false)),
358  average_(getAverage(dict, setAverage_)),
359  interpolationScheme_(interpolationCell<Type>::typeName)
360 {
361  if
362  (
363  mapper_.sampleDatabase()
364  && (
365  mapper_.mode() != mappedPatchBase::NEARESTPATCHFACE
366  && mapper_.mode() != mappedPatchBase::NEARESTPATCHFACEAMI
367  )
368  )
369  {
371  << "Mapping using the database only supported for "
372  << "sampleModes "
373  << mappedPatchBase::sampleModeNames_
374  [
375  mappedPatchBase::NEARESTPATCHFACE
376  ]
377  << " and "
378  << mappedPatchBase::sampleModeNames_
379  [
380  mappedPatchBase::NEARESTPATCHFACEAMI
381  ]
382  << exit(FatalError);
383  }
384 
385  if (mapper_.mode() == mappedPatchBase::NEARESTCELL)
386  {
387  dict.readEntry("interpolationScheme", interpolationScheme_);
388  }
389 
390  // Note: in database mode derived boundary conditions need to initialise
391  // fields
392 }
393 
394 
395 template<class Type>
397 (
398  const mappedPatchBase& mapper,
399  const fvPatchField<Type>& patchField,
400  const dictionary& dict,
401  const Field<Type>& fld
402 )
403 :
405 {
406  if (mapper_.sampleDatabase())
407  {
408  if (mapper_.mode() == mappedPatchBase::NEARESTPATCHFACE)
409  {
410  // Store my data on receive buffers so we have some initial data
411  initRetrieveField
412  (
413  patchField_.internalField().time(),
414  //patchField_.patch().boundaryMesh().mesh().name(),
415  mapper_.sampleRegion(),
416  //patchField_.patch().name(),
417  mapper_.samplePatch(),
418  mapper_.map().constructMap(),
419  patchField_.internalField().name(),
420  patchField_
421  );
422  }
423  else if (mapper_.mode() == mappedPatchBase::NEARESTPATCHFACEAMI)
424  {
425  // Depend on fall-back (sorting dummy field) in retrieveField
426  // since it would be too hard to determine the field that gives
427  // the wanted result after interpolation
428  }
429  }
430 }
431 
432 
433 template<class Type>
435 (
436  const mappedPatchBase& mapper,
437  const fvPatchField<Type>& patchField
438 )
439 :
440  mapper_(mapper),
441  patchField_(patchField),
442  fieldName_(patchField_.internalField().name()),
443  setAverage_(false),
444  average_(Zero),
445  interpolationScheme_(interpolationCell<Type>::typeName)
446 {}
447 
448 
449 template<class Type>
451 (
452  const mappedPatchFieldBase<Type>& mapper
453 )
454 :
455  mapper_(mapper.mapper_),
456  patchField_(mapper.patchField_),
457  fieldName_(mapper.fieldName_),
458  setAverage_(mapper.setAverage_),
459  average_(mapper.average_),
460  interpolationScheme_(mapper.interpolationScheme_)
461 {}
462 
463 
464 template<class Type>
466 (
467  const mappedPatchBase& mapper,
468  const fvPatchField<Type>& patchField,
469  const mappedPatchFieldBase<Type>& base
470 )
471 :
472  mapper_(mapper),
473  patchField_(patchField),
474  fieldName_(base.fieldName_),
475  setAverage_(base.setAverage_),
476  average_(base.average_),
477  interpolationScheme_(base.interpolationScheme_)
478 {}
479 
480 
481 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
482 
483 template<class Type>
484 template<class Type2>
487 {
489 
490  if (mapper_.sameRegion())
491  {
492  if (fieldName == patchField_.internalField().name())
493  {
494  // Optimisation: bypass field lookup
495  return
496  dynamic_cast<const fieldType&>
497  (
498  patchField_.internalField()
499  );
500  }
501  else
502  {
503  const fvMesh& thisMesh = patchField_.patch().boundaryMesh().mesh();
504  return thisMesh.template lookupObject<fieldType>(fieldName);
505  }
506  }
507 
508  const fvMesh& nbrMesh = refCast<const fvMesh>(mapper_.sampleMesh());
509  return nbrMesh.template lookupObject<fieldType>(fieldName);
510 }
511 
512 
513 template<class Type>
516 {
517  return sampleField<Type>(fieldName_);
518 }
519 
520 
521 template<class Type>
522 template<class T>
524 (
525  const word& fieldName,
526  Field<T>& fld
527 ) const
528 {
529  if (mapper_.sampleDatabase())
530  {
531  const label myComm = mapper_.getCommunicator(); // Get or create
532 
533  if (mapper_.mode() != mappedPatchBase::NEARESTPATCHFACEAMI)
534  {
535  // Store my data on send buffers
536  storeField
537  (
538  patchField_.internalField().time(),
539  patchField_.patch().boundaryMesh().mesh().name(),
540  patchField_.patch().name(),
541  myComm,
542  mapper_.map().subMap(),
543  fieldName,
544  fld
545  );
546  // Construct my data from receive buffers
547  fld.setSize(mapper_.map().constructSize());
548  retrieveField
549  (
550  true, // allow unset
551  patchField_.internalField().time(),
552  mapper_.sampleRegion(),
553  mapper_.samplePatch(),
554  myComm,
555  mapper_.map().constructMap(),
556  fieldName,
557  fld
558  );
559  }
560  else
561  {
562  const AMIPatchToPatchInterpolation& AMI = mapper_.AMI();
563 
564  // The AMI does an interpolateToSource/ToTarget. This is a
565  // mapDistribute (so using subMap/constructMap) and then a
566  // weighted sum. We'll store the sent data as before and
567  // do the weighted summation after the retrieveField
568 
569  if (mapper_.masterWorld())
570  {
571  // See AMIInterpolation::interpolateToSource. Use tgtMap,
572  // srcAddress, srcWeights
573  storeAndRetrieveField
574  (
575  fieldName,
576  myComm,
577  AMI.srcMap().subMap(),
578  AMI.tgtMap().constructSize(),
579  AMI.tgtMap().constructMap(),
580  AMI.srcAddress(),
581  AMI.srcWeights(),
582  fld
583  );
584  }
585  else
586  {
587  // See AMIInterpolation::interpolateToTarget.
588  // Use srcMap, tgtAddress, tgtWeights
589  storeAndRetrieveField
590  (
591  fieldName,
592  myComm,
593  AMI.tgtMap().subMap(),
594  AMI.srcMap().constructSize(),
595  AMI.srcMap().constructMap(),
596  AMI.tgtAddress(),
597  AMI.tgtWeights(),
598  fld
599  );
600  }
601  }
602  }
603  else
604  {
605  mapper_.distribute(fld);
606  }
607 }
608 
609 
610 template<class Type>
611 //template<class T>
614 (
615 // const GeometricField<T, fvPatchField, volMesh>& fld
616 ) const
617 {
619 
620  // Since we're inside initEvaluate/evaluate there might be processor
621  // comms underway. Change the tag we use.
622  int oldTag = UPstream::msgType();
623  UPstream::msgType() = oldTag + 1;
624 
625  const fvMesh& thisMesh = patchField_.patch().boundaryMesh().mesh();
626 
627  // Result of obtaining remote values
628  auto tnewValues = tmp<Field<Type>>::New();
629  auto& newValues = tnewValues.ref();
630 
631  switch (mapper_.mode())
632  {
633  case mappedPatchBase::NEARESTCELL:
634  {
635  const fieldType& fld = sampleField();
636  const mapDistribute& distMap = mapper_.map();
637 
638  if (interpolationScheme_ != interpolationCell<Type>::typeName)
639  {
640  if (!mapper_.sameWorld() || mapper_.sampleDatabase())
641  {
643  << "Interpolating cell values from different world"
644  << " or database currently not supported"
645  << exit(FatalError);
646  }
647 
648  const fvMesh& nbrMesh =
649  refCast<const fvMesh>(mapper_.sampleMesh());
650 
651  // Send back sample points to the processor that holds the cell
652  vectorField samples(mapper_.samplePoints());
653 
654  distMap.reverseDistribute
655  (
656  (
657  mapper_.sameRegion()
658  ? thisMesh.nCells()
659  : nbrMesh.nCells()
660  ),
661  point::max,
662  samples
663  );
664 
665  auto interpolator =
667  (
668  interpolationScheme_,
669  fld
670  );
671 
672  const auto& interp = *interpolator;
673 
674  newValues.setSize(samples.size(), pTraits<Type>::max);
675  forAll(samples, celli)
676  {
677  if (samples[celli] != point::max)
678  {
679  newValues[celli] = interp.interpolate
680  (
681  samples[celli],
682  celli
683  );
684  }
685  }
686  }
687  else
688  {
689  newValues = fld;
690  }
691 
692  distribute(fieldName_, newValues);
693 
694  break;
695  }
696  case mappedPatchBase::NEARESTPATCHFACE:
697  case mappedPatchBase::NEARESTPATCHFACEAMI:
698  {
699  if (mapper_.sameWorld())
700  {
701  const fvMesh& nbrMesh =
702  refCast<const fvMesh>(mapper_.sampleMesh());
703  const fieldType& fld = sampleField();
704 
705  const label nbrPatchID =
706  nbrMesh.boundaryMesh().findPatchID(mapper_.samplePatch());
707 
708  if (nbrPatchID < 0)
709  {
711  << "Unable to find sample patch " << mapper_.samplePatch()
712  << " in region " << mapper_.sampleRegion()
713  << " for patch " << patchField_.patch().name() << nl
714  << abort(FatalError);
715  }
716 
717  const auto& nbrField = fld;
718 
719  newValues = nbrField.boundaryField()[nbrPatchID];
720  }
721  else
722  {
723  // Start off from my patch values, let distribute function below
724  // do all the work
725  newValues = patchField_;
726  }
727  distribute(fieldName_, newValues);
728 
729  break;
730  }
731  case mappedPatchBase::NEARESTFACE:
732  {
733  Field<Type> allValues;
734  if (mapper_.sameWorld())
735  {
736  const fvMesh& nbrMesh =
737  refCast<const fvMesh>(mapper_.sampleMesh());
738  const fieldType& fld = sampleField();
739 
740  allValues.setSize(nbrMesh.nFaces(), Zero);
741 
742  const auto& nbrField = fld;
743 
744  for (const fvPatchField<Type>& pf : nbrField.boundaryField())
745  {
746  label faceStart = pf.patch().start();
747 
748  forAll(pf, facei)
749  {
750  allValues[faceStart++] = pf[facei];
751  }
752  }
753  }
754  else
755  {
756  // Start off from my patch values. Let distribute function below
757  // do all the work
758  allValues.setSize(thisMesh.nFaces(), Zero);
759 
760  const fieldType& thisFld = dynamic_cast<const fieldType&>
761  (
762  patchField_.internalField()
763  );
764 
765  for (const fvPatchField<Type>& pf : thisFld.boundaryField())
766  {
767  label faceStart = pf.patch().start();
768 
769  forAll(pf, facei)
770  {
771  allValues[faceStart++] = pf[facei];
772  }
773  }
774  }
775 
776  distribute(fieldName_, allValues);
777  newValues.transfer(allValues);
778 
779  break;
780  }
781  default:
782  {
784  << "Unknown sampling mode: " << mapper_.mode() << nl
785  << abort(FatalError);
786  }
787  }
788 
789  if (setAverage_)
790  {
791  Type averagePsi =
792  gSum(patchField_.patch().magSf()*newValues)
793  /gSum(patchField_.patch().magSf());
794 
795  if (mag(averagePsi) > 0.5*mag(average_))
796  {
797  newValues *= mag(average_)/mag(averagePsi);
798  }
799  else
800  {
801  newValues += (average_ - averagePsi);
802  }
803  }
804 
805  // Restore tag
806  UPstream::msgType() = oldTag;
807 
808  return tnewValues;
809 }
810 
811 
812 //template<class Type>
813 //Foam::tmp<Foam::Field<Type>>
814 //Foam::mappedPatchFieldBase<Type>::mappedField() const
815 //{
816 // const GeometricField<Type, fvPatchField, volMesh>& fld = sampleField();
817 // return mappedField<Type>(fld);
818 //}
819 
820 
821 template<class Type>
824 {
825  // Swap to obtain full local values of neighbour internal field
826  tmp<Field<Type>> tnbrIntFld(new Field<Type>());
827  Field<Type>& nbrIntFld = tnbrIntFld.ref();
828 
829  if (mapper_.sameWorld())
830  {
831  // Same world so lookup
832  const label nbrPatchID = mapper_.samplePolyPatch().index();
833  const auto& nbrField = this->sampleField();
834  nbrIntFld = nbrField.boundaryField()[nbrPatchID].patchInternalField();
835  }
836  else
837  {
838  // Different world so use my region,patch. Distribution below will
839  // do the reordering
840  nbrIntFld = patchField_.patchInternalField();
841  }
842 
843  // Since we're inside initEvaluate/evaluate there might be processor
844  // comms underway. Change the tag we use.
845  int oldTag = UPstream::msgType();
846  UPstream::msgType() = oldTag+1;
847 
848  distribute(fieldName_, nbrIntFld);
849 
850  // Restore tag
851  UPstream::msgType() = oldTag;
852 
853  return tnbrIntFld;
854 }
855 
856 
857 template<class Type>
860 {
861  // Swap to obtain full local values of neighbour internal field
862  tmp<scalarField> tnbrKDelta(new scalarField());
863  scalarField& nbrKDelta = tnbrKDelta.ref();
864 
865  if (mapper_.sameWorld())
866  {
867  // Same world so lookup
868  const auto& nbrMesh = refCast<const fvMesh>(this->mapper_.sampleMesh());
869  const label nbrPatchID = mapper_.samplePolyPatch().index();
870  const auto& nbrPatch = nbrMesh.boundary()[nbrPatchID];
871  nbrKDelta = nbrPatch.deltaCoeffs();
872  }
873  else
874  {
875  // Different world so use my region,patch. Distribution below will
876  // do the reordering
877  nbrKDelta = patchField_.patch().deltaCoeffs();
878  }
879 
880 
881  // Since we're inside initEvaluate/evaluate there might be processor
882  // comms underway. Change the tag we use.
883  const int oldTag = UPstream::msgType();
884  UPstream::msgType() = oldTag+1;
885 
886  distribute(fieldName_ + "_deltaCoeffs", nbrKDelta);
887 
888  // Restore tag
889  UPstream::msgType() = oldTag;
890 
891  return tnbrKDelta;
892 }
893 
894 
895 template<class Type>
897 (
898  const word& fieldName,
899  tmp<scalarField>& thisWeights,
900  tmp<scalarField>& nbrWeights
901 ) const
902 {
903  thisWeights = new scalarField(patchField_.patch().deltaCoeffs());
904  if (!fieldName.empty())
905  {
906  thisWeights.ref() *=
907  patchField_.patch().template lookupPatchField
908  <
910  scalar
911  >
912  (
913  fieldName
914  ).patchInternalField();
915  }
916 
917 
918  // Swap to obtain full local values of neighbour internal field
919 
920  if (mapper_.sameWorld())
921  {
922  // Same world so lookup
923  const auto& nbrMesh = refCast<const fvMesh>(mapper_.sampleMesh());
924  const label nbrPatchID = mapper_.samplePolyPatch().index();
925  const auto& nbrPatch = nbrMesh.boundary()[nbrPatchID];
926 
927  nbrWeights = new scalarField(nbrPatch.deltaCoeffs());
928 
929  if (!fieldName.empty())
930  {
931  // Weightfield is volScalarField
932  const auto& nbrWeightField =
933  nbrMesh.template lookupObject<volScalarField>(fieldName);
934  nbrWeights.ref() *=
935  nbrWeightField.boundaryField()[nbrPatchID].patchInternalField();
936  }
937  }
938  else
939  {
940  // Different world so use my region,patch. Distribution below will
941  // do the reordering
942  nbrWeights = new scalarField(thisWeights());
943  }
944 
945  // Since we're inside initEvaluate/evaluate there might be processor
946  // comms underway. Change the tag we use.
947  int oldTag = UPstream::msgType();
948  UPstream::msgType() = oldTag+1;
949 
950  distribute(fieldName_ + "_weights", nbrWeights.ref());
951 
952  // Restore tag
953  UPstream::msgType() = oldTag;
954 }
955 
956 
957 template<class Type>
959 (
960  const fvPatch& p,
962 )
963 {
964  if (!isA<mappedPatchBase>(p.patch()))
965  {
967  << "Incorrect patch type " << p.patch().type()
968  << " for patch " << p.patch().name()
969  << " of field " << iF.name()
970  << " in file " << iF.objectPath() << nl
971  << "Type should be a mappedPatch"
972  << exit(FatalError);
973  }
974  return refCast<const mappedPatchBase>(p.patch());
975 }
976 
977 
978 template<class Type>
979 template<class T>
981 (
982  const word& fieldName,
983  const Field<T>& fld
984 ) const
985 {
986  if (mapper_.sampleDatabase())
987  {
988  // Store my data on receive buffers (reverse of storeField;
989  // i.e. retrieveField will obtain patchField)
990  if (mapper_.mode() == mappedPatchBase::NEARESTPATCHFACE)
991  {
992  initRetrieveField
993  (
994  patchField_.internalField().time(),
995  mapper_.sampleRegion(),
996  mapper_.samplePatch(),
997  mapper_.map().constructMap(),
998  fieldName,
999  fld
1000  );
1001  }
1002  }
1003 }
1004 
1005 
1006 template<class Type>
1008 {
1009  os.writeEntry("field", fieldName_);
1010 
1011  if (setAverage_)
1012  {
1013  os.writeEntry("setAverage", "true");
1014  os.writeEntry("average", average_);
1015  }
1016 
1017  if (mapper_.mode() == mappedPatchBase::NEARESTCELL)
1018  {
1019  os.writeEntry("interpolationScheme", interpolationScheme_);
1020  }
1021 }
1022 
1023 
1024 // ************************************************************************* //
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::mapDistributeBase::subMap
const labelListList & subMap() const
From subsetted data back to original data.
Definition: mapDistributeBase.H:289
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::mappedPatchFieldBase::average_
const Type average_
Definition: mappedPatchFieldBase.H:131
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::mappedPatchFieldBase::mappedPatchFieldBase
mappedPatchFieldBase(const mappedPatchBase &mapper, const fvPatchField< Type > &patchField, const word &fieldName, const bool setAverage, const Type average, const word &interpolationScheme)
Construct from components.
Definition: mappedPatchFieldBase.C:321
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::mappedPatchFieldBase::patchField_
const fvPatchField< Type > & patchField_
Underlying patch field.
Definition: mappedPatchFieldBase.H:121
Foam::mappedPatchFieldBase::storeAndRetrieveField
bool storeAndRetrieveField(const word &fieldName, const label myComm, const labelListList &subMap, const label constructSize, const labelListList &constructMap, const labelListList &address, const scalarListList &weights, Field< T > &fld) const
Helper : storeField and retrieveField and interpolate. Leaves fld.
Definition: mappedPatchFieldBase.C:252
interpolationCell.H
Foam::mappedPatchFieldBase::mappedField
virtual tmp< Field< Type > > mappedField() const
Map sampleField onto *this patch.
Definition: mappedPatchFieldBase.C:614
mappedPatchFieldBase.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::mappedPatchFieldBase::storeField
void storeField(const objectRegistry &obr, const word &region, const word &patch, const label myComm, const labelListList &procToMap, const word &fieldName, const Field< T > &fld) const
Store elements of field onto (sub) registry.
Definition: mappedPatchFieldBase.C:54
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
Foam::AMIInterpolation::tgtAddress
const labelListList & tgtAddress() const
Return const access to target patch addressing.
Definition: AMIInterpolationI.H:195
Foam::AMIInterpolation::tgtWeights
const scalarListList & tgtWeights() const
Return const access to target patch weights.
Definition: AMIInterpolationI.H:207
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::interpolationCell
Uses the cell value for any point in the cell.
Definition: interpolationCell.H:50
Foam::mappedPatchFieldBase::fieldName_
word fieldName_
Name of field to sample.
Definition: mappedPatchFieldBase.H:124
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::mappedPatchFieldBase::retrieveField
bool retrieveField(const bool allowUnset, const objectRegistry &obr, const word &region, const word &patch, const label myComm, const labelListList &procToMap, const word &fieldName, Field< T > &fld) const
Construct field from registered elements.
Definition: mappedPatchFieldBase.C:105
Foam::AMIInterpolation::tgtMap
const mapDistribute & tgtMap() const
Definition: AMIInterpolationI.H:231
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::polyBoundaryMesh::mesh
const polyMesh & mesh() const noexcept
Return the mesh reference.
Definition: polyBoundaryMesh.H:152
Foam::Field< T >
Foam::fac::average
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:47
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::mappedPatchFieldBase::distribute
void distribute(const word &fieldName, Field< T > &newValues) const
Definition: mappedPatchFieldBase.C:524
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
samples
scalarField samples(nIntervals, Zero)
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:765
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
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
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mappedPatchFieldBase
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
Definition: mappedPatchFieldBase.H:105
Foam::mappedPatchFieldBase::mapper
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< Type, volMesh > &iF)
Check that patch is of correct type.
Definition: mappedPatchFieldBase.C:959
Foam::mappedPatchFieldBase::mappedWeightField
virtual tmp< scalarField > mappedWeightField() const
Map optional weightField onto *this patch. Default is deltaCoeffs.
Definition: mappedPatchFieldBase.C:859
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::mappedPatchFieldBase::initRetrieveField
void initRetrieveField(const objectRegistry &obr, const word &region, const word &patch, const labelListList &map, const word &fieldName, const Field< T > &fld) const
Construct field from registered elements.
Definition: mappedPatchFieldBase.C:201
os
OBJstream os(runTime.globalPath()/outputName)
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:216
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::mappedPatchFieldBase::interpolationScheme_
word interpolationScheme_
Interpolation scheme to use for nearestcell mode.
Definition: mappedPatchFieldBase.H:134
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::mapDistributeBase::constructSize
label constructSize() const
Constructed data size.
Definition: mapDistributeBase.H:277
Foam::multiplyWeightedOp
Definition: ops.H:250
Foam::mappedPatchFieldBase::write
virtual void write(Ostream &os) const
Write.
Definition: mappedPatchFieldBase.C:1007
Foam::AMIInterpolation::srcMap
const mapDistribute & srcMap() const
Definition: AMIInterpolationI.H:177
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::mappedPatchFieldBase::mapper_
const mappedPatchBase & mapper_
Mapping engine.
Definition: mappedPatchFieldBase.H:118
Foam::mappedPatchFieldBase::sampleField
const GeometricField< Type, fvPatchField, volMesh > & sampleField() const
Field to sample. Either on my or nbr mesh.
Definition: mappedPatchFieldBase.C:515
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< labelList >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:56
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::mappedPatchFieldBase::mappedInternalField
virtual tmp< Field< Type > > mappedInternalField() const
Map internal of sampleField onto *this patch.
Definition: mappedPatchFieldBase.C:823
Foam::AMIInterpolation
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
Definition: AMIInterpolation.H:79
Foam::mapDistribute::reverseDistribute
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
Definition: mapDistributeTemplates.C:182
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::plusEqOp
Definition: ops.H:72
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
Foam::mapDistributeBase::constructMap
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
Definition: mapDistributeBase.H:301
Foam::mappedPatchFieldBase::setAverage_
const bool setAverage_
If true adjust the mapped field to maintain average value average_.
Definition: mappedPatchFieldBase.H:127
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::AMIInterpolation::srcAddress
const labelListList & srcAddress() const
Return const access to source patch addressing.
Definition: AMIInterpolationI.H:129
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::AMIInterpolation::srcWeights
const scalarListList & srcWeights() const
Return const access to source patch weights.
Definition: AMIInterpolationI.H:141
mappedPatchBase.H