BlendedInterfacialModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2014-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "phaseSystem.H"
30 #include "surfaceInterpolate.H"
31 #include "zeroDimensionalFvMesh.H"
32 #include "mathematicalConstants.H"
33 #include "writeFile.H"
34 #include "triFace.H"
35 #include "noSetWriter.H"
36 #include "noSurfaceWriter.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace blendedInterfacialModel
43 {
44 
45 template<class GeoField>
47 
48 template<>
50 {
51  return f;
52 }
53 
54 template<>
56 {
57  return fvc::interpolate(f);
58 }
59 
60 template<class ModelType>
61 inline bool valid(const PtrList<ModelType>& l)
62 {
63  forAll(l, i)
64  {
65  if (l.set(i))
66  {
67  return true;
68  }
69  }
70  return false;
71 }
72 
73 } // End namespace blendedInterfacialModel
74 } // End namespace Foam
75 
76 
77 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
78 
79 template<class ModelType>
81 {
82  // Only generate warnings once per timestep
83  const label timeIndex = interface_.mesh().time().timeIndex();
84  if (checkTimeIndex_ == timeIndex) return;
85  checkTimeIndex_ = timeIndex;
86 
87  const phaseModel& phase1 = interface_.phase1();
88  const phaseModel& phase2 = interface_.phase2();
89 
90  const bool can1In2 = blending_->canBeContinuous(1);
91  const bool can2In1 = blending_->canBeContinuous(0);
92  const bool canS = blending_->canSegregate();
93 
94  // Warnings associated with redundant model specification
95 
96  if
97  (
98  !can1In2
99  && (
100  model1DispersedIn2_.valid()
101  || blendedInterfacialModel::valid(models1DispersedIn2Displaced_)
102  )
103  )
104  {
106  << "A " << ModelType::typeName << " was provided for "
107  << dispersedPhaseInterface(phase1, phase2).name()
108  << " but the associated blending does not permit " << phase2.name()
109  << " to be continuous so this model will not be used" << endl;
110  }
111 
112  if
113  (
114  !can2In1
115  && (
116  model2DispersedIn1_.valid()
117  || blendedInterfacialModel::valid(models2DispersedIn1Displaced_)
118  )
119  )
120  {
122  << "A " << ModelType::typeName << " was provided for "
123  << dispersedPhaseInterface(phase2, phase1).name()
124  << " but the associated blending does not permit " << phase1.name()
125  << " to be continuous so this model will not be used" << endl;
126  }
127 
128  if
129  (
130  !canS
131  && (
132  model1SegregatedWith2_.valid()
133  || blendedInterfacialModel::valid(models1SegregatedWith2Displaced_)
134  )
135  )
136  {
138  << "A " << ModelType::typeName << " was provided for "
139  << segregatedPhaseInterface(phase1, phase2).name()
140  << " but the associated blending does not permit segregation"
141  << " so this model will not be used" << endl;
142  }
143 
144  if
145  (
146  modelGeneral_.valid()
147  && (can1In2 || can2In1 || canS)
148  && (!can1In2 || model1DispersedIn2_.valid())
149  && (!can2In1 || model2DispersedIn1_.valid())
150  && (!canS || model1SegregatedWith2_.valid())
151  )
152  {
154  << "A " << ModelType::typeName << " was provided for "
155  << phaseInterface(phase1, phase2).name()
156  << " but other displaced and/or segregated models apply"
157  << " across the entire phase fraction space so this model"
158  << " will not be used" << endl;
159  }
160 
161  forAll(interface_.fluid().phases(), phasei)
162  {
163  const phaseModel& phaseD = interface_.fluid().phases()[phasei];
164 
165  if
166  (
167  modelsGeneralDisplaced_.set(phasei)
168  && (can1In2 || can2In1 || canS)
169  && (!can1In2 || models1DispersedIn2Displaced_.set(phasei))
170  && (!can2In1 || models2DispersedIn1Displaced_.set(phasei))
171  && (!canS || models1SegregatedWith2Displaced_.set(phasei))
172  )
173  {
175  << "A " << ModelType::typeName << " was provided for "
176  << displacedPhaseInterface(phase1, phase2, phaseD).name()
177  << " but other displaced and/or segregated models apply"
178  << " across the entire phase fraction space so this model"
179  << " will not be used" << endl;
180  }
181  }
182 
183  // Warnings associated with gaps in the blending space
184 
185  if (can1In2 && !modelGeneral_.valid() && !model1DispersedIn2_.valid())
186  {
188  << "Blending for " << ModelType::typeName << "s permits "
189  << phase2.name() << " to become continuous, but no model was "
190  << "provided for " << dispersedPhaseInterface(phase1, phase2).name()
191  << ". Consider adding a model for this configuration (or for "
192  << phaseInterface(phase1, phase2).name() << "), or if no model is "
193  << "needed then add a \"none\" model to suppress this warning."
194  << endl;
195  }
196 
197  if (can2In1 && !modelGeneral_.valid() && !model2DispersedIn1_.valid())
198  {
200  << "Blending for " << ModelType::typeName << "s permits "
201  << phase1.name() << " to become continuous, but no model was "
202  << "provided for " << dispersedPhaseInterface(phase2, phase1).name()
203  << ". Consider adding a model for this configuration (or for "
204  << phaseInterface(phase1, phase2).name() << "), or if no model is "
205  << "needed then add a \"none\" model to suppress this warning."
206  << endl;
207  }
208 
209  if (canS && !modelGeneral_.valid() && !model1SegregatedWith2_.valid())
210  {
212  << "Blending for " << ModelType::typeName << "s permits "
213  << "segregation but no model was provided for "
214  << segregatedPhaseInterface(phase2, phase1).name()
215  << ". Consider adding a model for this configuration (or for "
216  << phaseInterface(phase1, phase2).name() << "), or if no model is "
217  << "needed then add a \"none\" model to suppress this warning."
218  << endl;
219  }
220 }
221 
222 
223 template<class ModelType>
224 template<class GeoMesh>
226 (
227  const UPtrList<const volScalarField>& alphas,
228  tmp<GeometricField<scalar, GeoMesh>>& fG,
229  tmp<GeometricField<scalar, GeoMesh>>& f1D2,
230  tmp<GeometricField<scalar, GeoMesh>>& f2D1,
231  tmp<GeometricField<scalar, GeoMesh>>& fS,
232  PtrList<GeometricField<scalar, GeoMesh>>& fGD,
233  PtrList<GeometricField<scalar, GeoMesh>>& f1D2D,
234  PtrList<GeometricField<scalar, GeoMesh>>& f2D1D,
235  PtrList<GeometricField<scalar, GeoMesh>>& fSD
236 ) const
237 {
238  typedef GeometricField<scalar, GeoMesh> scalarGeoField;
239 
240  // Create a constant field
241  auto constant = [&](const scalar k)
242  {
243  return
245  (
246  Foam::name(k),
247  alphas.first().mesh(),
249  );
250  };
251 
252  // Get the dispersed blending functions
253  tmp<scalarGeoField> F1D2, F2D1;
254  if
255  (
256  model1DispersedIn2_.valid()
257  || model1SegregatedWith2_.valid()
258  || blendedInterfacialModel::valid(models1DispersedIn2Displaced_)
259  || blendedInterfacialModel::valid(models1SegregatedWith2Displaced_)
260  )
261  {
262  F1D2 =
263  blendedInterfacialModel::interpolate<scalarGeoField>
264  (
265  blending_->f1DispersedIn2(alphas)
266  );
267  }
268  if
269  (
270  model2DispersedIn1_.valid()
271  || model1SegregatedWith2_.valid()
272  || blendedInterfacialModel::valid(models2DispersedIn1Displaced_)
273  || blendedInterfacialModel::valid(models1SegregatedWith2Displaced_)
274  )
275  {
276  F2D1 =
277  blendedInterfacialModel::interpolate<scalarGeoField>
278  (
279  blending_->f2DispersedIn1(alphas)
280  );
281  }
282 
283  // Construct non-displaced coefficients
284  {
285  if (model1DispersedIn2_.valid())
286  {
287  f1D2 = F1D2().clone();
288  }
289 
290  if (model2DispersedIn1_.valid())
291  {
292  f2D1 = F2D1().clone();
293  }
294 
295  if (model1SegregatedWith2_.valid())
296  {
297  fS = constant(1);
298  if (f1D2.valid()) { fS.ref() -= f1D2(); }
299  if (f2D1.valid()) { fS.ref() -= f2D1(); }
300  }
301 
302  if (modelGeneral_.valid())
303  {
304  fG = constant(1);
305  if (f1D2.valid()) { fG.ref() -= f1D2(); }
306  if (f2D1.valid()) { fG.ref() -= f2D1(); }
307  if (fS.valid()) { fG.ref() -= fS(); }
308  }
309  }
310 
311  // Get the displaced blending function
312  tmp<scalarGeoField> FD;
313  if
314  (
315  blendedInterfacialModel::valid(modelsGeneralDisplaced_)
316  || blendedInterfacialModel::valid(models1DispersedIn2Displaced_)
317  || blendedInterfacialModel::valid(models2DispersedIn1Displaced_)
318  || blendedInterfacialModel::valid(models1SegregatedWith2Displaced_)
319  )
320  {
321  FD =
322  blendedInterfacialModel::interpolate<scalarGeoField>
323  (
324  blending_->fDisplaced(alphas)
325  );
326  }
327 
328  // Construct displaced coefficients
329  tmp<scalarGeoField> fDSum;
330  if
331  (
332  blendedInterfacialModel::valid(modelsGeneralDisplaced_)
333  || blendedInterfacialModel::valid(models1DispersedIn2Displaced_)
334  || blendedInterfacialModel::valid(models2DispersedIn1Displaced_)
335  || blendedInterfacialModel::valid(models1SegregatedWith2Displaced_)
336  )
337  {
338  fDSum = constant(0);
339 
340  forAll(alphas, phasei)
341  {
342  const phaseModel& phaseD = interface_.fluid().phases()[phasei];
343 
344  if (interface_.contains(phaseD)) continue;
345 
346  // Weight the contribution of a dispersed model for this phase
347  // according to the phase fraction in the subset of phases that are
348  // not part of this interface. This seems like a reasonable
349  // assumption if a stochastic viewpoint is taken. However, it is
350  // possible that other forms may be desired in which case this
351  // could be abstracted and controlled by the blending method.
352  tmp<scalarGeoField> alpha;
353  if
354  (
355  models1DispersedIn2Displaced_.set(phasei)
356  || models2DispersedIn1Displaced_.set(phasei)
357  || models1SegregatedWith2Displaced_.set(phasei)
358  || modelsGeneralDisplaced_.set(phasei)
359  )
360  {
361  alpha =
362  blendedInterfacialModel::interpolate<scalarGeoField>
363  (
364  alphas[phasei]
365  /max
366  (
367  1
368  - alphas[interface_.phase1().index()]
369  - alphas[interface_.phase2().index()],
370  phaseD.residualAlpha()
371  )
372  );
373  }
374 
375  if (models1DispersedIn2Displaced_.set(phasei))
376  {
377  f1D2D.set(phasei, alpha()*FD()*F1D2());
378  fDSum.ref() += f1D2D[phasei];
379  }
380 
381  if (models2DispersedIn1Displaced_.set(phasei))
382  {
383  f2D1D.set(phasei, alpha()*FD()*F2D1());
384  fDSum.ref() += f2D1D[phasei];
385  }
386 
387  if (models1SegregatedWith2Displaced_.set(phasei))
388  {
389  fSD.set(phasei, alpha()*FD());
390  if (f1D2D.set(phasei)) fSD[phasei] -= f1D2D[phasei];
391  if (f2D1D.set(phasei)) fSD[phasei] -= f2D1D[phasei];
392  fDSum.ref() += fSD[phasei];
393  }
394 
395  if (modelsGeneralDisplaced_.set(phasei))
396  {
397  fGD.set(phasei, alpha()*FD());
398  if (f1D2D.set(phasei)) fGD[phasei] -= f1D2D[phasei];
399  if (f2D1D.set(phasei)) fGD[phasei] -= f2D1D[phasei];
400  if (fSD.set(phasei)) fGD[phasei] -= fSD[phasei];
401  fDSum.ref() += fGD[phasei];
402  }
403  }
404  }
405 
406  // Remove the displaced part from the non-displaced coefficients
407  if (fDSum.valid())
408  {
409  if (f1D2.valid()) f1D2.ref() *= 1 - fDSum();
410  if (f2D1.valid()) f2D1.ref() *= 1 - fDSum();
411  if (fS.valid()) fS.ref() *= 1 - fDSum();
412  if (fG.valid()) fG.ref() *= 1 - fDSum();
413  }
414 }
415 
416 
417 template<class ModelType>
419 (
420  const word& format
421 ) const
422 {
423  using namespace constant::mathematical;
424 
425  const phaseSystem& fluid = interface_.fluid();
426  const label nPhases = fluid.phases().size();
427 
428  // Construct geometry and phase fraction values
430  faceList faces;
431  wordList fieldNames(nPhases);
432  PtrList<scalarField> fields(nPhases);
433  if (nPhases == 1)
434  {
435  // Single value
436  fieldNames[0] = fluid.phases()[0].volScalarField::name();
437  fields.set(0, new scalarField(1, 1));
438  }
439  else if (nPhases == 2)
440  {
441  // Single axis, using the first phase as the x-coordinate
442  static const label nDivisions = 128;
443  const scalarField divisionis(scalarList(identityMap(nDivisions)));
444  fieldNames[0] = fluid.phases()[0].volScalarField::name();
445  fieldNames[1] = fluid.phases()[1].volScalarField::name();
446  fields.set(0, new scalarField(divisionis/(nDivisions - 1)));
447  fields.set(1, new scalarField(1 - fields[0]));
448  }
449  else
450  {
451  // Polygon with as many vertices as there are phases. Each phase
452  // fraction equals one at a unique vertex, they vary linearly between
453  // vertices, and vary smoothly in the interior of the polygon. The sum
454  // of phase fractions is always one throughout the polygon (partition of
455  // unity). This is done using Wachspress coordinates. This does not
456  // cover the entire N-dimensional phase fraction space for N >= 4 (it is
457  // hard to imagine how that could be visualised) but it provides enough
458  // to give a good indication of what is going on.
459 
460  // Create the nodes of the blending space polygon
461  List<point> phaseNodes(nPhases);
462  forAll(phaseNodes, phasei)
463  {
464  const scalar theta = 2*pi*phasei/nPhases;
465  phaseNodes[phasei] = point(cos(theta), sin(theta), 0);
466  }
467 
468  // Create points within the polygon
469  static const label nDivisions = 32;
470  points.append(point::zero);
471  for (label divi = 0; divi < nDivisions; ++ divi)
472  {
473  const scalar s = scalar(divi + 1)/nDivisions;
474 
475  forAll(phaseNodes, phasei)
476  {
477  for (label i = 0; i < divi + 1; ++ i)
478  {
479  const scalar t = scalar(i)/(divi + 1);
480 
481  points.append
482  (
483  s*(1 - t)*phaseNodes[phasei]
484  + s*t*phaseNodes[(phasei + 1) % nPhases]
485  );
486  }
487  }
488  }
489 
490  // Create triangles within the polygon
491  forAll(phaseNodes, phasei)
492  {
493  faces.append(triFace(0, phasei + 1, (phasei + 1) % nPhases + 1));
494  }
495  for (label divi = 1; divi < nDivisions; ++ divi)
496  {
497  const label pointi0 = nPhases*(divi - 1)*divi/2 + 1;
498  const label pointi1 = nPhases*divi*(divi + 1)/2 + 1;
499 
500  forAll(phaseNodes, phasei)
501  {
502  for (label i = 0; i < divi + 1; ++ i)
503  {
504  const label pi00 =
505  pointi0
506  + ((phasei*divi + i) % (nPhases*divi));
507  const label pi01 =
508  pointi0
509  + ((phasei*divi + i + 1) % (nPhases*divi));
510  const label pi10 =
511  pointi1
512  + ((phasei*(divi + 1) + i) % (nPhases*(divi + 1)));
513  const label pi11 =
514  pointi1
515  + ((phasei*(divi + 1) + i + 1) % (nPhases*(divi + 1)));
516 
517  faces.append(triFace({pi00, pi10, pi11}));
518  if (i < divi) faces.append(triFace({pi00, pi11, pi01}));
519  }
520  }
521  }
522 
523  // Create phase fraction fields
524  forAll(fluid.phases(), phasei)
525  {
526  fieldNames[phasei] = fluid.phases()[phasei].volScalarField::name();
527  fields.set(phasei, new scalarField(points.size(), 0));
528 
529  const label phasei0 = (phasei + nPhases - 1) % nPhases;
530  const label phasei1 = (phasei + 1) % nPhases;
531 
532  const point& node0 = phaseNodes[phasei0];
533  const point& node = phaseNodes[phasei];
534  const point& node1 = phaseNodes[phasei1];
535 
536  forAll(points, i)
537  {
538  const scalar A = mag((node - node0) ^ (node1 - node0));
539  const scalar A0 = mag((node - node0) ^ (points[i] - node0));
540  const scalar A1 = mag((node1 - node) ^ (points[i] - node));
541 
542  if (A0 < rootSmall)
543  {
544  fields[phasei][i] =
545  great*mag(points[i] - node0)/mag(node - node0);
546  }
547  else if (A1 < rootSmall)
548  {
549  fields[phasei][i] =
550  great*mag(points[i] - node1)/mag(node - node1);
551  }
552  else
553  {
554  fields[phasei][i] = A/A0/A1;
555  }
556  }
557  }
558  forAll(points, i)
559  {
560  scalar s = 0;
561  forAll(fluid.phases(), phasei)
562  {
563  s += fields[phasei][i];
564  }
565  forAll(fluid.phases(), phasei)
566  {
567  fields[phasei][i] /= s;
568  }
569  }
570  }
571 
572  // Add the model coefficient fields
573  {
574  // Create alpha fields on a zero-dimensional one-cell mesh
575  const fvMesh mesh(zeroDimensionalFvMesh(fluid.mesh()));
576  PtrList<volScalarField> alphas(nPhases);
577  forAll(fluid.phases(), phasei)
578  {
579  alphas.set
580  (
581  phasei,
582  new volScalarField
583  (
584  IOobject
585  (
586  fluid.phases()[phasei].volScalarField::name(),
587  mesh.time().name(),
588  mesh
589  ),
590  mesh,
591  dimensionedScalar(dimless, fields[phasei][0])
592  )
593  );
594  }
595 
596  // Construct blending coefficient fields
597  {
598  tmp<volScalarField> fG, f1D2, f2D1, fS;
599  PtrList<volScalarField> fGD(nPhases);
600  PtrList<volScalarField> f1D2D(nPhases);
601  PtrList<volScalarField> f2D1D(nPhases);
602  PtrList<volScalarField> fSD(nPhases);
603  calculateBlendingCoeffs
604  (
605  alphas.convert<const volScalarField>(),
606  fG, f1D2, f2D1, fS,
607  fGD, f1D2D, f2D1D, fSD
608  );
609 
610  const phaseModel& phase1 = interface_.phase1();
611  const phaseModel& phase2 = interface_.phase2();
612 
613  auto addField = [&](const phaseInterface& interface)
614  {
615  fieldNames.append(interface.name());
616  fields.append(new scalarField(fields.first().size()));
617  };
618 
619  if (fG.valid())
620  {
621  addField(phaseInterface(phase1, phase2));
622  }
623  if (f1D2.valid())
624  {
625  addField(dispersedPhaseInterface(phase1, phase2));
626  }
627  if (f2D1.valid())
628  {
629  addField(dispersedPhaseInterface(phase2, phase1));
630  }
631  if (fS.valid())
632  {
633  addField(segregatedPhaseInterface(phase2, phase1));
634  }
635 
636  forAll(fluid.phases(), phasei)
637  {
638  const phaseModel& phaseD = fluid.phases()[phasei];
639 
640  if (fGD.set(phasei))
641  {
642  addField(displacedPhaseInterface(phase1, phase2, phaseD));
643  }
644  if (f1D2D.set(phasei))
645  {
646  addField
647  (
648  dispersedDisplacedPhaseInterface
649  (
650  phase1,
651  phase2,
652  phaseD
653  )
654  );
655  }
656  if (f2D1D.set(phasei))
657  {
658  addField
659  (
660  dispersedDisplacedPhaseInterface
661  (
662  phase2,
663  phase1,
664  phaseD
665  )
666  );
667  }
668  if (fSD.set(phasei))
669  {
670  addField
671  (
672  segregatedDisplacedPhaseInterface
673  (
674  phase1,
675  phase2,
676  phaseD
677  )
678  );
679  }
680  }
681  }
682 
683  // Populate blending coefficient fields
684  forAll(fields.first(), i)
685  {
686  forAll(fluid.phases(), phasei)
687  {
688  alphas[phasei] = fields[phasei][i];
689  }
690 
691  tmp<volScalarField> fG, f1D2, f2D1, fS;
692  PtrList<volScalarField> fGD(nPhases);
693  PtrList<volScalarField> f1D2D(nPhases);
694  PtrList<volScalarField> f2D1D(nPhases);
695  PtrList<volScalarField> fSD(nPhases);
696  calculateBlendingCoeffs
697  (
698  alphas.convert<const volScalarField>(),
699  fG, f1D2, f2D1, fS,
700  fGD, f1D2D, f2D1D, fSD
701  );
702 
703  label fieldi = nPhases;
704 
705  if (fG.valid())
706  {
707  fields[fieldi ++][i] = fG()[0];
708  }
709  if (f1D2.valid())
710  {
711  fields[fieldi ++][i] = f1D2()[0];
712  }
713  if (f2D1.valid())
714  {
715  fields[fieldi ++][i] = f2D1()[0];
716  }
717  if (fS.valid())
718  {
719  fields[fieldi ++][i] = fS()[0];
720  }
721 
722  forAll(fluid.phases(), phasei)
723  {
724  if (fGD.set(phasei))
725  {
726  fields[fieldi ++][i] = fGD[phasei][0];
727  }
728  if (f1D2D.set(phasei))
729  {
730  fields[fieldi ++][i] = f1D2D[phasei][0];
731  }
732  if (f2D1D.set(phasei))
733  {
734  fields[fieldi ++][i] = f2D1D[phasei][0];
735  }
736  if (fSD.set(phasei))
737  {
738  fields[fieldi ++][i] = fSD[phasei][0];
739  }
740  }
741  }
742  }
743 
744  // Write
745  const fileName path =
746  fluid.mesh().time().globalPath()
747  /functionObjects::writeFile::outputPrefix
748  /ModelType::typeName;
749  Info<< "Writing blending coefficients to " << path/interface_.name()
750  << endl;
751  if (nPhases <= 2)
752  {
753  // Strip out the first field and shuffle everything else up
754  const autoPtr<scalarField> field0 = fields.set(0, nullptr);
755  const word field0Name = fieldNames[0];
756  for (label fieldi = 1; fieldi < fields.size(); ++ fieldi)
757  {
758  fields.set(fieldi - 1, fields.set(fieldi, nullptr).ptr());
759  fieldNames[fieldi - 1] = fieldNames[fieldi];
760  }
761  fields.resize(fields.size() - 1);
763 
765  (
766  format,
767  IOstream::ASCII,
768  IOstream::UNCOMPRESSED
769  )->write
770  (
771  path,
772  interface_.name(),
773  coordSet(true, field0Name, field0),
774  fieldNames,
775  fields
776  );
777  }
778  else
779  {
781  (
782  format,
783  IOstream::ASCII,
784  IOstream::UNCOMPRESSED
785  )->write
786  (
787  path,
788  interface_.name(),
789  points,
790  faces,
791  true,
792  fieldNames,
793  fields
794  );
795  }
796 }
797 
798 
799 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
800 
801 template<class ModelType>
802 template<class Type, class GeoMesh, class ... Args>
805 (
806  tmp<GeometricField<Type, GeoMesh>>(ModelType::*method)(Args ...) const,
807  const word& name,
808  const dimensionSet& dims,
809  Args ... args
810 ) const
811 {
812  check();
813 
814  typedef GeometricField<scalar, GeoMesh> scalarGeoField;
815  typedef GeometricField<Type, GeoMesh> typeGeoField;
816 
817  // Get the blending coefficients
818  const label nPhases = interface_.fluid().phases().size();
819  tmp<scalarGeoField> fG, f1D2, f2D1, fS;
820  PtrList<scalarGeoField> fGD(nPhases);
821  PtrList<scalarGeoField> f1D2D(nPhases);
822  PtrList<scalarGeoField> f2D1D(nPhases);
823  PtrList<scalarGeoField> fSD(nPhases);
824  calculateBlendingCoeffs
825  (
826  interface_.fluid().phases()
827  .PtrList<phaseModel>::convert<const volScalarField>(),
828  fG, f1D2, f2D1, fS,
829  fGD, f1D2D, f2D1D, fSD
830  );
831 
832  // Construct the result
835  (
836  ModelType::typeName + ":"
837  + IOobject::groupName(name, interface_.name()),
838  interface_.mesh(),
839  dimensioned<Type>(dims, Zero)
840  );
841 
842  // Add the model contributions to the result
843  if (modelGeneral_.valid())
844  {
845  x.ref() += fG*(modelGeneral_().*method)(args ...);
846  }
847  if (model1DispersedIn2_.valid())
848  {
849  x.ref() += f1D2*(model1DispersedIn2_().*method)(args ...);
850  }
851  if (model2DispersedIn1_.valid())
852  {
853  x.ref() += f2D1*(model2DispersedIn1_().*method)(args ...);
854  }
855  if (model1SegregatedWith2_.valid())
856  {
857  x.ref() += fS*(model1SegregatedWith2_().*method)(args ...);
858  }
859  forAll(interface_.fluid().phases(), phasei)
860  {
861  if (modelsGeneralDisplaced_.set(phasei))
862  {
863  x.ref() +=
864  fGD[phasei]
865  *(modelsGeneralDisplaced_[phasei].*method)(args ...);
866  }
867  if (models1DispersedIn2Displaced_.set(phasei))
868  {
869  x.ref() +=
870  f1D2D[phasei]
871  *(models1DispersedIn2Displaced_[phasei].*method)(args ...);
872  }
873  if (models2DispersedIn1Displaced_.set(phasei))
874  {
875  x.ref() +=
876  f2D1D[phasei]
877  *(models2DispersedIn1Displaced_[phasei].*method)(args ...);
878  }
879  if (models1SegregatedWith2Displaced_.set(phasei))
880  {
881  x.ref() +=
882  fSD[phasei]
883  *(models1SegregatedWith2Displaced_[phasei].*method)(args ...);
884  }
885  }
886 
887  return x;
888 }
889 
890 
891 template<class ModelType>
892 template<class Type, class GeoMesh, class ... Args>
895 (
897  (ModelType::*method)(Args ...) const,
898  const word& name,
899  const dimensionSet& dims,
900  Args ... args
901 ) const
902 {
903  check();
904 
905  typedef GeometricField<scalar, GeoMesh> scalarGeoField;
906  typedef GeometricField<Type, GeoMesh> typeGeoField;
907 
908  // Get the blending coefficients
909  const label nPhases = interface_.fluid().phases().size();
910  tmp<scalarGeoField> fG, f1D2, f2D1, fS;
911  PtrList<scalarGeoField> fGD(nPhases);
912  PtrList<scalarGeoField> f1D2D(nPhases);
913  PtrList<scalarGeoField> f2D1D(nPhases);
914  PtrList<scalarGeoField> fSD(nPhases);
915  calculateBlendingCoeffs
916  (
917  interface_.fluid().phases()
918  .PtrList<phaseModel>::convert<const volScalarField>(),
919  fG, f1D2, f2D1, fS,
920  fGD, f1D2D, f2D1D, fSD
921  );
922 
923  // Construct the result
925 
926  // Add the model contributions to the result
927  auto addToXs = [&]
928  (
929  const scalarGeoField& f,
930  const HashPtrTable<typeGeoField>& dxs
931  )
932  {
933  forAllConstIter(typename HashPtrTable<typeGeoField>, dxs, dxIter)
934  {
935  if (xs.found(dxIter.key()))
936  {
937  *xs[dxIter.key()] += f**dxIter();
938  }
939  else
940  {
941  xs.insert
942  (
943  dxIter.key(),
945  (
946  ModelType::typeName + ':'
947  + IOobject::groupName
948  (
949  IOobject::groupName(name, dxIter.key()),
950  interface_.name()
951  ),
952  f**dxIter()
953  ).ptr()
954  );
955  }
956  }
957  };
958  if (modelGeneral_.valid())
959  {
960  addToXs(fG, (modelGeneral_().*method)(args ...));
961  }
962  if (model1DispersedIn2_.valid())
963  {
964  addToXs(f1D2, (model1DispersedIn2_().*method)(args ...));
965  }
966  if (model2DispersedIn1_.valid())
967  {
968  addToXs(f2D1, (model2DispersedIn1_().*method)(args ...));
969  }
970  if (model1SegregatedWith2_.valid())
971  {
972  addToXs(fS, (model1SegregatedWith2_().*method)(args ...));
973  }
974  forAll(interface_.fluid().phases(), phasei)
975  {
976  if (modelsGeneralDisplaced_.set(phasei))
977  {
978  addToXs
979  (
980  fGD[phasei],
981  (modelsGeneralDisplaced_[phasei].*method)(args ...)
982  );
983  }
984  if (models1DispersedIn2Displaced_.set(phasei))
985  {
986  addToXs
987  (
988  f1D2D[phasei],
989  (models1DispersedIn2Displaced_[phasei].*method)(args ...)
990  );
991  }
992  if (models2DispersedIn1Displaced_.set(phasei))
993  {
994  addToXs
995  (
996  f2D1D[phasei],
997  (models2DispersedIn1Displaced_[phasei].*method)(args ...)
998  );
999  }
1000  if (models1SegregatedWith2Displaced_.set(phasei))
1001  {
1002  addToXs
1003  (
1004  fSD[phasei],
1005  (models1SegregatedWith2Displaced_[phasei].*method)(args ...)
1006  );
1007  }
1008  }
1009 
1010  return xs;
1011 }
1012 
1013 
1014 template<class ModelType>
1015 template<class ... Args>
1017 (
1018  bool (ModelType::*method)(Args ...) const,
1019  Args ... args
1020 ) const
1021 {
1022  check();
1023 
1024  bool result = false;
1025 
1026  if (modelGeneral_.valid())
1027  {
1028  result = result || (modelGeneral_().*method)(args ...);
1029  }
1030  if (model1DispersedIn2_.valid())
1031  {
1032  result = result || (model1DispersedIn2_().*method)(args ...);
1033  }
1034  if (model2DispersedIn1_.valid())
1035  {
1036  result = result || (model2DispersedIn1_().*method)(args ...);
1037  }
1038  if (model1SegregatedWith2_.valid())
1039  {
1040  result = result || (model1SegregatedWith2_().*method)(args ...);
1041  }
1042 
1043  forAll(interface_.fluid().phases(), phasei)
1044  {
1045  if (modelsGeneralDisplaced_.set(phasei))
1046  {
1047  result =
1048  result
1049  || (modelsGeneralDisplaced_[phasei].*method)(args ...);
1050  }
1051  if (models1DispersedIn2Displaced_.set(phasei))
1052  {
1053  result =
1054  result
1055  || (models1DispersedIn2Displaced_[phasei].*method)(args ...);
1056  }
1057  if (models2DispersedIn1Displaced_.set(phasei))
1058  {
1059  result =
1060  result
1061  || (models2DispersedIn1Displaced_[phasei].*method)(args ...);
1062  }
1063  if (models1SegregatedWith2Displaced_.set(phasei))
1064  {
1065  result =
1066  result
1067  || (models1SegregatedWith2Displaced_[phasei].*method)(args ...);
1068  }
1069  }
1070 
1071  return result;
1072 }
1073 
1074 
1075 template<class ModelType>
1076 template<class ... Args>
1078 (
1079  const hashedWordList& (ModelType::*method)(Args ...) const,
1080  Args ... args
1081 ) const
1082 {
1083  check();
1084 
1085  wordList result;
1086 
1087  if (modelGeneral_.valid())
1088  {
1089  result.append((modelGeneral_().*method)(args ...));
1090  }
1091  if (model1DispersedIn2_.valid())
1092  {
1093  result.append((model1DispersedIn2_().*method)(args ...));
1094  }
1095  if (model2DispersedIn1_.valid())
1096  {
1097  result.append((model2DispersedIn1_().*method)(args ...));
1098  }
1099  if (model1SegregatedWith2_.valid())
1100  {
1101  result.append((model1SegregatedWith2_().*method)(args ...));
1102  }
1103 
1104  forAll(interface_.fluid().phases(), phasei)
1105  {
1106  if (modelsGeneralDisplaced_.set(phasei))
1107  {
1108  result.append
1109  (
1110  (modelsGeneralDisplaced_[phasei].*method)(args ...)
1111  );
1112  }
1113  if (models1DispersedIn2Displaced_.set(phasei))
1114  {
1115  result.append
1116  (
1117  (models1DispersedIn2Displaced_[phasei].*method)(args ...)
1118  );
1119  }
1120  if (models2DispersedIn1Displaced_.set(phasei))
1121  {
1122  result.append
1123  (
1124  (models2DispersedIn1Displaced_[phasei].*method)(args ...)
1125  );
1126  }
1127  if (models1SegregatedWith2Displaced_.set(phasei))
1128  {
1129  result.append
1130  (
1131  (models1SegregatedWith2Displaced_[phasei].*method)(args ...)
1132  );
1133  }
1134  }
1135 
1136  return hashedWordList(move(result));
1137 }
1138 
1139 
1140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1141 
1142 template<class ModelType>
1143 template<class ... Args>
1145 (
1146  const dictionary& dict,
1147  const phaseInterface& interface,
1148  const dictionary& blendingDict,
1149  const Args& ... args
1150 )
1151 :
1152  regIOobject
1153  (
1154  IOobject
1155  (
1156  IOobject::groupName(typeName, interface.name()),
1157  interface.fluid().mesh().time().name(),
1158  interface.fluid().mesh()
1159  )
1160  ),
1161  interface_(interface),
1162  blending_
1163  (
1164  blendingMethod::New(ModelType::typeName, blendingDict, interface)
1165  ),
1166  modelGeneral_(),
1167  model1DispersedIn2_(),
1168  model2DispersedIn1_(),
1169  model1SegregatedWith2_(),
1170  modelsGeneralDisplaced_(interface.fluid().phases().size()),
1171  models1DispersedIn2Displaced_(interface.fluid().phases().size()),
1172  models2DispersedIn1Displaced_(interface.fluid().phases().size()),
1173  models1SegregatedWith2Displaced_(interface.fluid().phases().size()),
1174  checkTimeIndex_(-1)
1175 {
1176  // Construct the models
1177  PtrList<phaseInterface> interfaces;
1178  PtrList<ModelType> models;
1180  <
1181  ModelType,
1188  >
1189  (
1190  interfaces,
1191  models,
1192  interface.fluid(),
1193  dict,
1194  wordHashSet({"blending"}),
1195  interface,
1196  args ...
1197  );
1198 
1199  // Unpack the interface and model lists to populate the models used for the
1200  // different parts of the blending space
1201  forAll(interfaces, i)
1202  {
1203  const phaseInterface& interface = interfaces[i];
1204 
1205  autoPtr<ModelType>* modelPtrPtr;
1206  PtrList<ModelType>* modelPtrsPtr;
1207 
1208  if (isA<dispersedPhaseInterface>(interface))
1209  {
1210  const phaseModel& dispersed =
1211  refCast<const dispersedPhaseInterface>(interface).dispersed();
1212 
1213  modelPtrPtr =
1214  interface_.index(dispersed) == 0
1215  ? &model1DispersedIn2_
1216  : &model2DispersedIn1_;
1217  modelPtrsPtr =
1218  interface_.index(dispersed) == 0
1219  ? &models1DispersedIn2Displaced_
1220  : &models2DispersedIn1Displaced_;
1221  }
1222  else if (isA<segregatedPhaseInterface>(interface))
1223  {
1224  modelPtrPtr = &model1SegregatedWith2_;
1225  modelPtrsPtr = &models1SegregatedWith2Displaced_;
1226  }
1227  else
1228  {
1229  modelPtrPtr = &modelGeneral_;
1230  modelPtrsPtr = &modelsGeneralDisplaced_;
1231  }
1232 
1233  if (!isA<displacedPhaseInterface>(interface))
1234  {
1235  *modelPtrPtr = models.set(i, nullptr);
1236  }
1237  else
1238  {
1239  const phaseModel& displacing =
1240  refCast<const displacedPhaseInterface>(interface).displacing();
1241 
1242  modelPtrsPtr->set(displacing.index(), models.set(i, nullptr));
1243  }
1244  }
1245 
1246  // Write out the blending space if needed
1247  if (blendingDict.found("format"))
1248  {
1249  const word format = blendingDict.lookup<word>("format");
1250 
1251  const label nPhases = interface_.fluid().phases().size();
1252 
1253  if
1254  (
1255  (nPhases <= 2 && format != noSetWriter::typeName)
1256  || (nPhases > 2 && format != noSurfaceWriter::typeName)
1257  )
1258  {
1259  postProcessBlendingCoefficients(format);
1260  }
1261  }
1262 }
1263 
1264 
1265 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1266 
1267 template<class ModelType>
1269 {}
1270 
1271 
1272 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
1273 
1274 template<class ModelType>
1275 const Foam::phaseInterface&
1277 {
1278  return interface_;
1279 }
1280 
1281 
1282 template<class ModelType>
1284 {
1285  return os.good();
1286 }
1287 
1288 
1289 // ************************************************************************* //
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Wrapper class for interfacial models for which multiple instances of the model are used for different...
BlendedInterfacialModel(const dictionary &dict, const phaseInterface &interface, const dictionary &blendingDict, const Args &... args)
Construct from a dictionary, an interface and a blending dictionary.
const phaseInterface & interface() const
Access the interface.
bool writeData(Ostream &os) const
Dummy write for regIOobject.
tmp< GeometricField< Type, GeoMesh > > evaluate(tmp< GeometricField< Type, GeoMesh >>(ModelType::*method)(Args ...) const, const word &name, const dimensionSet &dims, Args ... args) const
Return a blended field.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
Abstract base class for functions that are used to combine interfacial sub-models according to the vo...
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
Dimension set for the base types.
Definition: dimensionSet.H:125
Generic dimensioned Type class.
const word & name() const
Return const reference to name.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Class to represent an interface between phases which has been displaced to some extent by a third pha...
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
A wordList with hashed indices for faster lookup by name.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseSystem & fluid() const
Return the phase system.
virtual word name() const
Name.
label index(const phaseModel &phase) const
Return the index of the given phase. Generates a FatalError if.
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:163
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
Class to represent a interface between phases where the two phases are considered to be segregated,...
Class to represent a interface between phases where the two phases are considered to be segregated; t...
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
volScalarField scalarField(fieldObject, mesh)
static List< word > fieldNames
Definition: globalFoam.H:46
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:234
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
bool valid(const PtrList< ModelType > &l)
tmp< GeoField > interpolate(tmp< volScalarField > f)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
static const coefficient A("A", dimPressure, 611.21)
Namespace for OpenFOAM.
void generateInterfacialModels(PtrList< phaseInterface > &interfaces, PtrList< ModelType > &models, const phaseSystem &fluid, const dictionary &dict, const wordHashSet &ignoreKeys, const phaseInterface &interface, const Args &... args)
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
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:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:265
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
vector point
Point is a vector.
Definition: point.H:41
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dictionary & blendingDict(const phaseSystem &fluid, const dictionary &dict)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
HashSet wordHashSet
A HashSet with word keys.
Definition: HashSet.H:210
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
fvMesh zeroDimensionalFvMesh(const objectRegistry &db)
Construct a zero-dimensional unit-cube FV mesh.
List< face > faceList
Definition: faceListFwd.H:41
dimensionedScalar cos(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label timeIndex
Definition: getTimeIndex.H:4
word format(conversionProperties.lookup("format"))
face triFace(3)
labelList f(nPoints)
dictionary dict
Foam::argList args(argc, argv)