VoFPatchTransfer.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) 2017 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 "VoFPatchTransfer.H"
30 #include "twoPhaseMixtureThermo.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace surfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(VoFPatchTransfer, 0);
45 addToRunTimeSelectionTable(transferModel, VoFPatchTransfer, dictionary);
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 VoFPatchTransfer::VoFPatchTransfer
50 (
51  surfaceFilmRegionModel& film,
52  const dictionary& dict
53 )
54 :
55  transferModel(type(), film, dict),
56  deltaFactorToVoF_
57  (
58  coeffDict_.getOrDefault<scalar>("deltaFactorToVoF", 1.0)
59  ),
60  deltaFactorToFilm_
61  (
62  coeffDict_.getOrDefault<scalar>("deltaFactorToFilm", 0.5)
63  ),
64  alphaToVoF_
65  (
66  coeffDict_.getOrDefault<scalar>("alphaToVoF", 0.5)
67  ),
68  alphaToFilm_
69  (
70  coeffDict_.getOrDefault<scalar>("alphaToFilm", 0.1)
71  ),
72  transferRateCoeff_
73  (
74  coeffDict_.getOrDefault<scalar>("transferRateCoeff", 0.1)
75  ),
76  patchIDs_(),
77  patchTransferredMasses_()
78 {
79  const polyBoundaryMesh& pbm = film.regionMesh().boundaryMesh();
80 
81  const label nPatches
82  (
83  pbm.size() - film.regionMesh().globalData().processorPatches().size()
84  );
85 
86  wordRes patchNames;
87  if (coeffDict_.readIfPresent("patches", patchNames))
88  {
89  // Can also use pbm.indices(), but no warnings...
90  patchIDs_ = pbm.patchSet(patchNames).sortedToc();
91 
92  Info<< " applying to " << patchIDs_.size() << " patches:" << nl;
93 
94  for (const label patchi : patchIDs_)
95  {
96  Info<< " " << pbm[patchi].name() << endl;
97  }
98  }
99  else
100  {
101  Info<< " applying to all patches" << endl;
102 
103  patchIDs_ = identity(nPatches);
104  }
105 
106  patchTransferredMasses_.resize(patchIDs_.size(), Zero);
107 
108  if (patchIDs_.empty())
109  {
111  << "No patches selected"
112  << exit(FatalError);
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
118 
119 VoFPatchTransfer::~VoFPatchTransfer()
120 {}
121 
122 
123 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
124 
126 (
127  scalarField& availableMass,
128  scalarField& massToTransfer
129 )
130 {
132 }
133 
134 
136 (
137  scalarField& availableMass,
138  scalarField& massToTransfer,
139  scalarField& energyToTransfer
140 )
141 {
142  // Do not correct if no patches selected
143  if (!patchIDs_.size()) return;
144 
145  const scalarField& delta = film().delta();
146  const scalarField& rho = film().rho();
147  const scalarField& magSf = film().magSf();
148 
149  const polyBoundaryMesh& pbm = film().regionMesh().boundaryMesh();
150 
151 
152  const twoPhaseMixtureThermo& thermo
153  (
154  film().primaryMesh().lookupObject<twoPhaseMixtureThermo>
155  (
157  )
158  );
159 
160  const volScalarField& heVoF = thermo.thermo1().he();
161  const volScalarField& TVoF = thermo.thermo1().T();
162  const volScalarField CpVoF(thermo.thermo1().Cp());
163  const volScalarField& rhoVoF = thermo.thermo1().rho()();
164  const volScalarField& alphaVoF = thermo.alpha1();
165 
166  forAll(patchIDs_, pidi)
167  {
168  const label patchi = patchIDs_[pidi];
169  label primaryPatchi = -1;
170 
171  forAll(film().intCoupledPatchIDs(), i)
172  {
173  const label filmPatchi = film().intCoupledPatchIDs()[i];
174 
175  if (filmPatchi == patchi)
176  {
177  primaryPatchi = film().primaryPatchIDs()[i];
178  }
179  }
180 
181  if (primaryPatchi != -1)
182  {
183  scalarField deltaCoeffs
184  (
185  film().primaryMesh().boundary()[primaryPatchi].deltaCoeffs()
186  );
187  film().toRegion(patchi, deltaCoeffs);
188 
189  scalarField hp(heVoF.boundaryField()[primaryPatchi]);
190  film().toRegion(patchi, hp);
191 
192  scalarField Tp(TVoF.boundaryField()[primaryPatchi]);
193  film().toRegion(patchi, Tp);
194 
195  scalarField Cpp(CpVoF.boundaryField()[primaryPatchi]);
196  film().toRegion(patchi, Cpp);
197 
198  scalarField rhop(rhoVoF.boundaryField()[primaryPatchi]);
199  film().toRegion(patchi, rhop);
200 
201  scalarField alphap(alphaVoF.boundaryField()[primaryPatchi]);
202  film().toRegion(patchi, alphap);
203 
204  scalarField Vp
205  (
206  film().primaryMesh().boundary()[primaryPatchi]
207  .patchInternalField(film().primaryMesh().V())
208  );
209  film().toRegion(patchi, Vp);
210 
211  const polyPatch& pp = pbm[patchi];
212  const labelList& faceCells = pp.faceCells();
213 
214  // Accumulate the total mass removed from patch
215  scalar dMassPatch = 0;
216 
217  forAll(faceCells, facei)
218  {
219  const label celli = faceCells[facei];
220 
221  scalar dMass = 0;
222 
223  if
224  (
225  delta[celli] > 2*deltaFactorToVoF_/deltaCoeffs[facei]
226  || alphap[facei] > alphaToVoF_
227  )
228  {
229  dMass =
230  transferRateCoeff_*delta[celli]*rho[celli]*magSf[celli];
231 
232  massToTransfer[celli] += dMass;
233  energyToTransfer[celli] += dMass*film().hs()[celli];
234  }
235 
236  if
237  (
238  alphap[facei] > 0
239  && delta[celli] > 0
240  && delta[celli] < 2*deltaFactorToFilm_/deltaCoeffs[facei]
241  && alphap[facei] < alphaToFilm_
242  )
243  {
244  dMass =
245  -transferRateCoeff_*alphap[facei]*rhop[facei]*Vp[facei];
246 
247  massToTransfer[celli] += dMass;
248  energyToTransfer[celli] += dMass*hp[facei];
249  }
250 
251  availableMass[celli] -= dMass;
252  dMassPatch += dMass;
253  }
254 
255  patchTransferredMasses_[pidi] += dMassPatch;
256  addToTransferredMass(dMassPatch);
257  }
258  }
259 
261 
262  if (writeTime())
263  {
264  scalarField patchTransferredMasses0
265  (
266  getModelProperty<scalarField>
267  (
268  "patchTransferredMasses",
269  scalarField(patchTransferredMasses_.size(), Zero)
270  )
271  );
272 
273  scalarField patchTransferredMassTotals(patchTransferredMasses_);
274  Pstream::listCombineGather
275  (
276  patchTransferredMassTotals,
277  plusEqOp<scalar>()
278  );
279  patchTransferredMasses0 += patchTransferredMassTotals;
280 
281  setModelProperty<scalarField>
282  (
283  "patchTransferredMasses",
284  patchTransferredMasses0
285  );
286 
287  patchTransferredMasses_ = 0;
288  }
289 }
290 
291 
292 void VoFPatchTransfer::patchTransferredMassTotals
293 (
294  scalarField& patchMasses
295 ) const
296 {
297  // Do not correct if no patches selected
298  if (!patchIDs_.size()) return;
299 
300  scalarField patchTransferredMasses
301  (
302  getModelProperty<scalarField>
303  (
304  "patchTransferredMasses",
305  scalarField(patchTransferredMasses_.size(), Zero)
306  )
307  );
308 
309  scalarField patchTransferredMassTotals(patchTransferredMasses_);
310  Pstream::listCombineGather(patchTransferredMassTotals, plusEqOp<scalar>());
311 
312  forAll(patchIDs_, pidi)
313  {
314  const label patchi = patchIDs_[pidi];
315  patchMasses[patchi] +=
316  patchTransferredMasses[pidi] + patchTransferredMassTotals[pidi];
317  }
318 }
319 
320 
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 
323 } // End namespace surfaceFilmModels
324 } // End namespace regionModels
325 } // End namespace Foam
326 
327 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
faceListList boundary
scalar delta
const polyBoundaryMesh & pbm
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const word dictName("faMeshDefinition")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList patchNames(nPatches)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
defineTypeNameAndDebug(kinematicSingleLayer, 0)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127