VoFPatchTransfer.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) 2017-2021 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 
26 #include "VoFPatchTransfer.H"
28 #include "thermoSingleLayer.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace regionModels
36 {
37 namespace surfaceFilmModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(VoFPatchTransfer, 0);
43 addToRunTimeSelectionTable(transferModel, VoFPatchTransfer, dictionary);
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  surfaceFilmRegionModel& film,
50  const dictionary& dict
51 )
52 :
53  transferModel(type(), film, dict),
54  deltaFactorToVoF_
55  (
56  coeffDict_.lookupOrDefault<scalar>("deltaFactorToVoF", 1.0)
57  ),
58  deltaFactorToFilm_
59  (
60  coeffDict_.lookupOrDefault<scalar>("deltaFactorToFilm", 0.5)
61  ),
62  alphaToVoF_
63  (
64  coeffDict_.lookupOrDefault<scalar>("alphaToVoF", 0.5)
65  ),
66  alphaToFilm_
67  (
68  coeffDict_.lookupOrDefault<scalar>("alphaToFilm", 0.1)
69  ),
70  transferRateCoeff_
71  (
72  coeffDict_.lookupOrDefault<scalar>("transferRateCoeff", 0.1)
73  )
74 {
75  const polyBoundaryMesh& pbm = film.regionMesh().boundaryMesh();
76  patchIDs_.setSize
77  (
78  pbm.size() - film.regionMesh().globalData().processorPatches().size()
79  );
80 
81  if (coeffDict_.found("patches"))
82  {
83  const wordReList patchNames(coeffDict_.lookup("patches"));
84  const labelHashSet patchSet = pbm.patchSet(patchNames);
85 
86  Info<< " applying to patches:" << nl;
87 
88  label pidi = 0;
89  forAllConstIter(labelHashSet, patchSet, iter)
90  {
91  const label patchi = iter.key();
92  patchIDs_[pidi++] = patchi;
93  Info<< " " << pbm[patchi].name() << endl;
94  }
95  patchIDs_.setSize(pidi);
96  patchTransferredMasses_.setSize(pidi, 0);
97  }
98  else
99  {
100  Info<< " applying to all patches" << endl;
101 
102  forAll(patchIDs_, patchi)
103  {
104  patchIDs_[patchi] = patchi;
105  }
106 
107  patchTransferredMasses_.setSize(patchIDs_.size(), 0);
108  }
109 
110  if (!patchIDs_.size())
111  {
113  << "No patches selected"
114  << exit(FatalError);
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
120 
121 VoFPatchTransfer::~VoFPatchTransfer()
122 {}
123 
124 
125 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
126 
128 (
129  scalarField& availableMass,
130  scalarField& massToTransfer,
131  vectorField& momentumToTransfer
132 )
133 {
135 }
136 
137 
139 (
140  scalarField& availableMass,
141  scalarField& massToTransfer,
142  vectorField& momentumToTransfer,
143  scalarField& energyToTransfer
144 )
145 {
146  // Do not correct if no patches selected
147  if (!patchIDs_.size()) return;
148 
149  const thermoSingleLayer& film = filmType<thermoSingleLayer>();
150 
151  const scalarField& delta = film.delta();
152  const scalarField& rho = film.rho();
153  const scalarField& magSf = film.magSf();
154 
155  const polyBoundaryMesh& pbm = film.regionMesh().boundaryMesh();
156 
157 
158  const compressibleTwoPhaseMixture& thermo
159  (
160  film.primaryMesh().lookupObject<compressibleTwoPhaseMixture>
161  (
162  "phaseProperties"
163  )
164  );
165 
166  const volVectorField& UVoF
167  (
168  film.primaryMesh().lookupObject<volVectorField>("U")
169  );
170 
171  const volScalarField& alphaVoF = thermo.alpha1();
172  const volScalarField& rhoVoF = thermo.thermo1().rho()();
173  const volScalarField& heVoF = thermo.thermo1().he();
174  const volScalarField& TVoF = thermo.thermo1().T();
175  const volScalarField CpVoF(thermo.thermo1().Cp());
176 
177  forAll(patchIDs_, pidi)
178  {
179  const label patchi = patchIDs_[pidi];
180  label primaryPatchi = -1;
181 
182  forAll(film.intCoupledPatchIDs(), i)
183  {
184  const label filmPatchi = film.intCoupledPatchIDs()[i];
185 
186  if (filmPatchi == patchi)
187  {
188  primaryPatchi = film.primaryPatchIDs()[i];
189  }
190  }
191 
192  if (primaryPatchi != -1)
193  {
194  scalarField deltaCoeffs
195  (
196  film.primaryMesh().boundary()[primaryPatchi].deltaCoeffs()
197  );
198  film.toRegion(patchi, deltaCoeffs);
199 
200  scalarField alphap(alphaVoF.boundaryField()[primaryPatchi]);
201  film.toRegion(patchi, alphap);
202 
203  scalarField rhop(rhoVoF.boundaryField()[primaryPatchi]);
204  film.toRegion(patchi, rhop);
205 
206  vectorField Up(UVoF.boundaryField()[primaryPatchi]);
207  film.toRegion(patchi, Up);
208 
209  scalarField hp(heVoF.boundaryField()[primaryPatchi]);
210  film.toRegion(patchi, hp);
211 
212  scalarField Tp(TVoF.boundaryField()[primaryPatchi]);
213  film.toRegion(patchi, Tp);
214 
215  scalarField Cpp(CpVoF.boundaryField()[primaryPatchi]);
216  film.toRegion(patchi, Cpp);
217 
218  scalarField Vp
219  (
220  film.primaryMesh().boundary()[primaryPatchi]
221  .patchInternalField(film.primaryMesh().V())
222  );
223  film.toRegion(patchi, Vp);
224 
225  const polyPatch& pp = pbm[patchi];
226  const labelList& faceCells = pp.faceCells();
227 
228  const vectorField& U = film.U();
229  const scalarField& he = film.thermo().he();
230 
231  // Accumulate the total mass removed from patch
232  scalar dMassPatch = 0;
233 
234  forAll(faceCells, facei)
235  {
236  const label celli = faceCells[facei];
237 
238  scalar dMass = 0;
239 
240  if
241  (
242  delta[celli] > 2*deltaFactorToVoF_/deltaCoeffs[facei]
243  || alphap[facei] > alphaToVoF_
244  )
245  {
246  dMass =
247  transferRateCoeff_*delta[celli]*rho[celli]*magSf[celli];
248 
249  massToTransfer[celli] += dMass;
250  momentumToTransfer[celli] += dMass*U[celli];
251  energyToTransfer[celli] += dMass*he[celli];
252  }
253 
254  if
255  (
256  alphap[facei] > 0
257  && delta[celli] < 2*deltaFactorToFilm_/deltaCoeffs[facei]
258  && alphap[facei] < alphaToFilm_
259  )
260  {
261  dMass =
262  -transferRateCoeff_*alphap[facei]*rhop[facei]*Vp[facei];
263 
264  massToTransfer[celli] += dMass;
265  momentumToTransfer[celli] += dMass*Up[facei];
266  energyToTransfer[celli] += dMass*hp[facei];
267  }
268 
269  availableMass[celli] -= dMass;
270  dMassPatch += dMass;
271  }
272 
273  patchTransferredMasses_[pidi] += dMassPatch;
274  addToTransferredMass(dMassPatch);
275  }
276  }
277 
279 
280  if (writeTime())
281  {
282  scalarField patchTransferredMasses0
283  (
284  getModelProperty<scalarField>
285  (
286  "patchTransferredMasses",
287  scalarField(patchTransferredMasses_.size(), 0)
288  )
289  );
290 
291  scalarField patchTransferredMassTotals(patchTransferredMasses_);
293  (
294  patchTransferredMassTotals,
295  plusEqOp<scalar>()
296  );
297  patchTransferredMasses0 += patchTransferredMassTotals;
298 
299  setModelProperty<scalarField>
300  (
301  "patchTransferredMasses",
302  patchTransferredMasses0
303  );
304 
305  patchTransferredMasses_ = 0;
306  }
307 }
308 
309 
310 void VoFPatchTransfer::patchTransferredMassTotals
311 (
312  scalarField& patchMasses
313 ) const
314 {
315  // Do not correct if no patches selected
316  if (!patchIDs_.size()) return;
317 
318  scalarField patchTransferredMasses
319  (
320  getModelProperty<scalarField>
321  (
322  "patchTransferredMasses",
323  scalarField(patchTransferredMasses_.size(), 0)
324  )
325  );
326 
327  scalarField patchTransferredMassTotals(patchTransferredMasses_);
328  Pstream::listCombineGather(patchTransferredMassTotals, plusEqOp<scalar>());
329 
330  forAll(patchIDs_, pidi)
331  {
332  const label patchi = patchIDs_[pidi];
333  patchMasses[patchi] +=
334  patchTransferredMasses[pidi] + patchTransferredMassTotals[pidi];
335  }
336 }
337 
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 } // End namespace surfaceFilmModels
342 } // End namespace regionModels
343 } // End namespace Foam
344 
345 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
Definition: createFields.H:28
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
VoFPatchTransfer(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model.
Macros for easy insertion into run-time selection tables.
addToRunTimeSelectionTable(ejectionModel, BrunDrippingEjection, dictionary)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList patchNames(nPatches)
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const char nl
Definition: Ostream.H:260
label patchi
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
messageStream Info
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.