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 tmp<volScalarField> trhoVoF(thermo.thermo1().rho());
164  const volScalarField& rhoVoF = trhoVoF();
165  const volScalarField& alphaVoF = thermo.alpha1();
166 
167  forAll(patchIDs_, pidi)
168  {
169  const label patchi = patchIDs_[pidi];
170  label primaryPatchi = -1;
171 
172  forAll(film().intCoupledPatchIDs(), i)
173  {
174  const label filmPatchi = film().intCoupledPatchIDs()[i];
175 
176  if (filmPatchi == patchi)
177  {
178  primaryPatchi = film().primaryPatchIDs()[i];
179  }
180  }
181 
182  if (primaryPatchi != -1)
183  {
184  scalarField deltaCoeffs
185  (
186  film().primaryMesh().boundary()[primaryPatchi].deltaCoeffs()
187  );
188  film().toRegion(patchi, deltaCoeffs);
189 
190  scalarField hp(heVoF.boundaryField()[primaryPatchi]);
191  film().toRegion(patchi, hp);
192 
193  scalarField Tp(TVoF.boundaryField()[primaryPatchi]);
194  film().toRegion(patchi, Tp);
195 
196  scalarField Cpp(CpVoF.boundaryField()[primaryPatchi]);
197  film().toRegion(patchi, Cpp);
198 
199  scalarField rhop(rhoVoF.boundaryField()[primaryPatchi]);
200  film().toRegion(patchi, rhop);
201 
202  scalarField alphap(alphaVoF.boundaryField()[primaryPatchi]);
203  film().toRegion(patchi, alphap);
204 
205  scalarField Vp
206  (
207  film().primaryMesh().boundary()[primaryPatchi]
208  .patchInternalField(film().primaryMesh().V())
209  );
210  film().toRegion(patchi, Vp);
211 
212  const polyPatch& pp = pbm[patchi];
213  const labelList& faceCells = pp.faceCells();
214 
215  // Accumulate the total mass removed from patch
216  scalar dMassPatch = 0;
217 
218  forAll(faceCells, facei)
219  {
220  const label celli = faceCells[facei];
221 
222  scalar dMass = 0;
223 
224  if
225  (
226  delta[celli] > 2*deltaFactorToVoF_/deltaCoeffs[facei]
227  || alphap[facei] > alphaToVoF_
228  )
229  {
230  dMass =
231  transferRateCoeff_*delta[celli]*rho[celli]*magSf[celli];
232 
233  massToTransfer[celli] += dMass;
234  energyToTransfer[celli] += dMass*film().hs()[celli];
235  }
236 
237  if
238  (
239  alphap[facei] > 0
240  && delta[celli] > 0
241  && delta[celli] < 2*deltaFactorToFilm_/deltaCoeffs[facei]
242  && alphap[facei] < alphaToFilm_
243  )
244  {
245  dMass =
246  -transferRateCoeff_*alphap[facei]*rhop[facei]*Vp[facei];
247 
248  massToTransfer[celli] += dMass;
249  energyToTransfer[celli] += dMass*hp[facei];
250  }
251 
252  availableMass[celli] -= dMass;
253  dMassPatch += dMass;
254  }
255 
256  patchTransferredMasses_[pidi] += dMassPatch;
257  addToTransferredMass(dMassPatch);
258  }
259  }
260 
262 
263  if (writeTime())
264  {
265  scalarField patchTransferredMasses0
266  (
267  getModelProperty<scalarField>
268  (
269  "patchTransferredMasses",
270  scalarField(patchTransferredMasses_.size(), Zero)
271  )
272  );
273 
274  scalarField patchTransferredMassTotals(patchTransferredMasses_);
275  Pstream::listCombineGather
276  (
277  patchTransferredMassTotals,
278  plusEqOp<scalar>()
279  );
280  patchTransferredMasses0 += patchTransferredMassTotals;
281 
282  setModelProperty<scalarField>
283  (
284  "patchTransferredMasses",
285  patchTransferredMasses0
286  );
287 
288  patchTransferredMasses_ = 0;
289  }
290 }
291 
292 
293 void VoFPatchTransfer::patchTransferredMassTotals
294 (
295  scalarField& patchMasses
296 ) const
297 {
298  // Do not correct if no patches selected
299  if (!patchIDs_.size()) return;
300 
301  scalarField patchTransferredMasses
302  (
303  getModelProperty<scalarField>
304  (
305  "patchTransferredMasses",
306  scalarField(patchTransferredMasses_.size(), Zero)
307  )
308  );
309 
310  scalarField patchTransferredMassTotals(patchTransferredMasses_);
311  Pstream::listCombineGather(patchTransferredMassTotals, plusEqOp<scalar>());
312 
313  forAll(patchIDs_, pidi)
314  {
315  const label patchi = patchIDs_[pidi];
316  patchMasses[patchi] +=
317  patchTransferredMasses[pidi] + patchTransferredMassTotals[pidi];
318  }
319 }
320 
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 } // End namespace surfaceFilmModels
325 } // End namespace regionModels
326 } // End namespace Foam
327 
328 // ************************************************************************* //
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:608
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:72
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/v2406/OpenFOAM-v2406/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:696
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