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-2018 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"
27 #include "twoPhaseMixtureThermo.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace regionModels
35 {
36 namespace surfaceFilmModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(VoFPatchTransfer, 0);
42 addToRunTimeSelectionTable(transferModel, VoFPatchTransfer, dictionary);
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 VoFPatchTransfer::VoFPatchTransfer
47 (
48  surfaceFilmRegionModel& film,
49  const dictionary& dict
50 )
51 :
52  transferModel(type(), film, dict),
53  deltaFactorToVoF_
54  (
55  coeffDict_.lookupOrDefault<scalar>("deltaFactorToVoF", 1.0)
56  ),
57  deltaFactorToFilm_
58  (
59  coeffDict_.lookupOrDefault<scalar>("deltaFactorToFilm", 0.5)
60  ),
61  alphaToVoF_
62  (
63  coeffDict_.lookupOrDefault<scalar>("alphaToVoF", 0.5)
64  ),
65  alphaToFilm_
66  (
67  coeffDict_.lookupOrDefault<scalar>("alphaToFilm", 0.1)
68  ),
69  transferRateCoeff_
70  (
71  coeffDict_.lookupOrDefault<scalar>("transferRateCoeff", 0.1)
72  )
73 {
74  const polyBoundaryMesh& pbm = film.regionMesh().boundaryMesh();
75  patchIDs_.setSize
76  (
77  pbm.size() - film.regionMesh().globalData().processorPatches().size()
78  );
79 
80  if (coeffDict_.found("patches"))
81  {
82  const wordReList patchNames(coeffDict_.lookup("patches"));
83  const labelHashSet patchSet = pbm.patchSet(patchNames);
84 
85  Info<< " applying to patches:" << nl;
86 
87  label pidi = 0;
88  forAllConstIter(labelHashSet, patchSet, iter)
89  {
90  const label patchi = iter.key();
91  patchIDs_[pidi++] = patchi;
92  Info<< " " << pbm[patchi].name() << endl;
93  }
94  patchIDs_.setSize(pidi);
95  patchTransferredMasses_.setSize(pidi, 0);
96  }
97  else
98  {
99  Info<< " applying to all patches" << endl;
100 
101  forAll(patchIDs_, patchi)
102  {
103  patchIDs_[patchi] = patchi;
104  }
105 
106  patchTransferredMasses_.setSize(patchIDs_.size(), 0);
107  }
108 
109  if (!patchIDs_.size())
110  {
112  << "No patches selected"
113  << exit(FatalError);
114  }
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
119 
120 VoFPatchTransfer::~VoFPatchTransfer()
121 {}
122 
123 
124 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125 
127 (
128  scalarField& availableMass,
129  scalarField& massToTransfer
130 )
131 {
133 }
134 
135 
137 (
138  scalarField& availableMass,
139  scalarField& massToTransfer,
140  scalarField& energyToTransfer
141 )
142 {
143  // Do not correct if no patches selected
144  if (!patchIDs_.size()) return;
145 
146  const scalarField& delta = film().delta();
147  const scalarField& rho = film().rho();
148  const scalarField& magSf = film().magSf();
149 
150  const polyBoundaryMesh& pbm = film().regionMesh().boundaryMesh();
151 
152 
153  const twoPhaseMixtureThermo& thermo
154  (
155  film().primaryMesh().lookupObject<twoPhaseMixtureThermo>
156  (
158  )
159  );
160 
161  const volScalarField& heVoF = thermo.thermo1().he();
162  const volScalarField& TVoF = thermo.thermo1().T();
163  const volScalarField CpVoF(thermo.thermo1().Cp());
164  const volScalarField& rhoVoF = thermo.thermo1().rho()();
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] < 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(), 0)
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(), 0)
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 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
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:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
rhoReactionThermo & thermo
Definition: createFields.H:28
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
wordList patchNames(nPatches)
List< label > labelList
A List of labels.
Definition: labelList.H:56
const word dictName("particleTrackDict")
faceListList boundary(nPatches)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
volScalarField scalarField(fieldObject, mesh)
static const char nl
Definition: Ostream.H:265
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/ubuntu/OpenFOAM-6/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:72
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
label patchi
messageStream Info
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:366
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.