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  patchIDs_ = pbm.patchSet(patchNames).sortedToc();
90 
91  Info<< " applying to " << patchIDs_.size() << " patches:" << nl;
92 
93  for (const label patchi : patchIDs_)
94  {
95  Info<< " " << pbm[patchi].name() << endl;
96  }
97  }
98  else
99  {
100  Info<< " applying to all patches" << endl;
101 
102  patchIDs_ = identity(nPatches);
103  }
104 
105  patchTransferredMasses_.resize(patchIDs_.size(), Zero);
106 
107  if (patchIDs_.empty())
108  {
110  << "No patches selected"
111  << exit(FatalError);
112  }
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
117 
118 VoFPatchTransfer::~VoFPatchTransfer()
119 {}
120 
121 
122 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
123 
125 (
126  scalarField& availableMass,
127  scalarField& massToTransfer
128 )
129 {
131 }
132 
133 
135 (
136  scalarField& availableMass,
137  scalarField& massToTransfer,
138  scalarField& energyToTransfer
139 )
140 {
141  // Do not correct if no patches selected
142  if (!patchIDs_.size()) return;
143 
144  const scalarField& delta = film().delta();
145  const scalarField& rho = film().rho();
146  const scalarField& magSf = film().magSf();
147 
148  const polyBoundaryMesh& pbm = film().regionMesh().boundaryMesh();
149 
150 
151  const twoPhaseMixtureThermo& thermo
152  (
153  film().primaryMesh().lookupObject<twoPhaseMixtureThermo>
154  (
156  )
157  );
158 
159  const volScalarField& heVoF = thermo.thermo1().he();
160  const volScalarField& TVoF = thermo.thermo1().T();
161  const volScalarField CpVoF(thermo.thermo1().Cp());
162  const volScalarField& rhoVoF = thermo.thermo1().rho()();
163  const volScalarField& alphaVoF = thermo.alpha1();
164 
165  forAll(patchIDs_, pidi)
166  {
167  const label patchi = patchIDs_[pidi];
168  label primaryPatchi = -1;
169 
170  forAll(film().intCoupledPatchIDs(), i)
171  {
172  const label filmPatchi = film().intCoupledPatchIDs()[i];
173 
174  if (filmPatchi == patchi)
175  {
176  primaryPatchi = film().primaryPatchIDs()[i];
177  }
178  }
179 
180  if (primaryPatchi != -1)
181  {
182  scalarField deltaCoeffs
183  (
184  film().primaryMesh().boundary()[primaryPatchi].deltaCoeffs()
185  );
186  film().toRegion(patchi, deltaCoeffs);
187 
188  scalarField hp(heVoF.boundaryField()[primaryPatchi]);
189  film().toRegion(patchi, hp);
190 
191  scalarField Tp(TVoF.boundaryField()[primaryPatchi]);
192  film().toRegion(patchi, Tp);
193 
194  scalarField Cpp(CpVoF.boundaryField()[primaryPatchi]);
195  film().toRegion(patchi, Cpp);
196 
197  scalarField rhop(rhoVoF.boundaryField()[primaryPatchi]);
198  film().toRegion(patchi, rhop);
199 
200  scalarField alphap(alphaVoF.boundaryField()[primaryPatchi]);
201  film().toRegion(patchi, alphap);
202 
203  scalarField Vp
204  (
205  film().primaryMesh().boundary()[primaryPatchi]
206  .patchInternalField(film().primaryMesh().V())
207  );
208  film().toRegion(patchi, Vp);
209 
210  const polyPatch& pp = pbm[patchi];
211  const labelList& faceCells = pp.faceCells();
212 
213  // Accumulate the total mass removed from patch
214  scalar dMassPatch = 0;
215 
216  forAll(faceCells, facei)
217  {
218  const label celli = faceCells[facei];
219 
220  scalar dMass = 0;
221 
222  if
223  (
224  delta[celli] > 2*deltaFactorToVoF_/deltaCoeffs[facei]
225  || alphap[facei] > alphaToVoF_
226  )
227  {
228  dMass =
229  transferRateCoeff_*delta[celli]*rho[celli]*magSf[celli];
230 
231  massToTransfer[celli] += dMass;
232  energyToTransfer[celli] += dMass*film().hs()[celli];
233  }
234 
235  if
236  (
237  alphap[facei] > 0
238  && delta[celli] > 0
239  && delta[celli] < 2*deltaFactorToFilm_/deltaCoeffs[facei]
240  && alphap[facei] < alphaToFilm_
241  )
242  {
243  dMass =
244  -transferRateCoeff_*alphap[facei]*rhop[facei]*Vp[facei];
245 
246  massToTransfer[celli] += dMass;
247  energyToTransfer[celli] += dMass*hp[facei];
248  }
249 
250  availableMass[celli] -= dMass;
251  dMassPatch += dMass;
252  }
253 
254  patchTransferredMasses_[pidi] += dMassPatch;
255  addToTransferredMass(dMassPatch);
256  }
257  }
258 
260 
261  if (writeTime())
262  {
263  scalarField patchTransferredMasses0
264  (
265  getModelProperty<scalarField>
266  (
267  "patchTransferredMasses",
268  scalarField(patchTransferredMasses_.size(), Zero)
269  )
270  );
271 
272  scalarField patchTransferredMassTotals(patchTransferredMasses_);
273  Pstream::listCombineGather
274  (
275  patchTransferredMassTotals,
276  plusEqOp<scalar>()
277  );
278  patchTransferredMasses0 += patchTransferredMassTotals;
279 
280  setModelProperty<scalarField>
281  (
282  "patchTransferredMasses",
283  patchTransferredMasses0
284  );
285 
286  patchTransferredMasses_ = 0;
287  }
288 }
289 
290 
291 void VoFPatchTransfer::patchTransferredMassTotals
292 (
293  scalarField& patchMasses
294 ) const
295 {
296  // Do not correct if no patches selected
297  if (!patchIDs_.size()) return;
298 
299  scalarField patchTransferredMasses
300  (
301  getModelProperty<scalarField>
302  (
303  "patchTransferredMasses",
304  scalarField(patchTransferredMasses_.size(), Zero)
305  )
306  );
307 
308  scalarField patchTransferredMassTotals(patchTransferredMasses_);
309  Pstream::listCombineGather(patchTransferredMassTotals, plusEqOp<scalar>());
310 
311  forAll(patchIDs_, pidi)
312  {
313  const label patchi = patchIDs_[pidi];
314  patchMasses[patchi] +=
315  patchTransferredMasses[pidi] + patchTransferredMassTotals[pidi];
316  }
317 }
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 } // End namespace surfaceFilmModels
323 } // End namespace regionModels
324 } // End namespace Foam
325 
326 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:63
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
dictName
const word dictName("faMeshDefinition")
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
twoPhaseMixtureThermo.H
NotImplemented
#define NotImplemented
Definition: error.H:553
nPatches
const label nPatches
Definition: printMeshSummary.H:24
Foam::Info
messageStream Info
rho
rho
Definition: readInitialConditions.H:88
correct
fvOptions correct(rho)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
patchNames
wordList patchNames(nPatches)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
VoFPatchTransfer.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Definition: POSIX.C:717
Foam::identity
labelList identity(const label len, label start=0)
Definition: labelList.C:31
boundary
faceListList boundary
Definition: createBlockMesh.H:4