layerParameters.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) 2011-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 "layerParameters.H"
27 #include "polyBoundaryMesh.H"
28 #include "unitConversion.H"
29 #include "refinementSurfaces.H"
30 #include "searchableSurfaces.H"
31 #include "medialAxisMeshMover.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 const Foam::scalar Foam::layerParameters::defaultConcaveAngle = 90;
36 
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 Foam::scalar Foam::layerParameters::layerExpansionRatio
41 (
42  const label n,
43  const scalar totalOverFirst
44 ) const
45 {
46  if (n <= 1)
47  {
48  return 1;
49  }
50 
51  const scalar tol = 1e-8;
52 
53  if (mag(n - totalOverFirst) < tol)
54  {
55  return 1;
56  }
57 
58  const label maxIters = 100;
59 
60  // Calculate the bounds of the solution
61  scalar minR;
62  scalar maxR;
63 
64  if (totalOverFirst < n)
65  {
66  minR = 0;
67  maxR = pow(totalOverFirst/n, 1/(n-1));
68  }
69  else
70  {
71  minR = pow(totalOverFirst/n, 1/(n-1));
72  maxR = totalOverFirst/(n - 1);
73  }
74 
75  // Starting guess
76  scalar r = 0.5*(minR + maxR);
77 
78  for (label i = 0; i < maxIters; ++i)
79  {
80  const scalar prevr = r;
81  const scalar fx = pow(r, n) - totalOverFirst*r - (1 - totalOverFirst);
82  const scalar dfx = n*pow(r, n - 1) - totalOverFirst;
83 
84  r -= fx/dfx;
85 
86  if (mag(r - prevr) < tol)
87  {
88  break;
89  }
90  }
91 
92  return r;
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
98 Foam::layerParameters::layerParameters
99 (
100  const dictionary& dict,
101  const polyBoundaryMesh& boundaryMesh
102 )
103 :
104  dict_(dict),
105  numLayers_(boundaryMesh.size(), -1),
106  relativeSizes_(dict.lookup("relativeSizes")),
107  layerSpec_(ILLEGAL),
108  firstLayerThickness_(boundaryMesh.size(), -123),
109  finalLayerThickness_(boundaryMesh.size(), -123),
110  thickness_(boundaryMesh.size(), -123),
111  expansionRatio_(boundaryMesh.size(), -123),
112  minThickness_
113  (
114  boundaryMesh.size(),
115  readScalar(dict.lookup("minThickness"))
116  ),
117  featureAngle_(readScalar(dict.lookup("featureAngle"))),
118  concaveAngle_
119  (
120  dict.lookupOrDefault("concaveAngle", defaultConcaveAngle)
121  ),
122  nGrow_(readLabel(dict.lookup("nGrow"))),
123  maxFaceThicknessRatio_
124  (
125  readScalar(dict.lookup("maxFaceThicknessRatio"))
126  ),
127  nBufferCellsNoExtrude_
128  (
129  readLabel(dict.lookup("nBufferCellsNoExtrude"))
130  ),
131  nLayerIter_(readLabel(dict.lookup("nLayerIter"))),
132  nRelaxedIter_(labelMax),
133  additionalReporting_(dict.lookupOrDefault("additionalReporting", false)),
134  meshShrinker_
135  (
136  dict.lookupOrDefault
137  (
138  "meshShrinker",
139  medialAxisMeshMover::typeName
140  )
141  )
142 {
143  // Detect layer specification mode
144 
145  label nSpec = 0;
146 
147  bool haveFirst = dict.found("firstLayerThickness");
148  if (haveFirst)
149  {
150  firstLayerThickness_ = scalarField
151  (
152  boundaryMesh.size(),
153  readScalar(dict.lookup("firstLayerThickness"))
154  );
155  nSpec++;
156  }
157  bool haveFinal = dict.found("finalLayerThickness");
158  if (haveFinal)
159  {
160  finalLayerThickness_ = scalarField
161  (
162  boundaryMesh.size(),
163  readScalar(dict.lookup("finalLayerThickness"))
164  );
165  nSpec++;
166  }
167  bool haveTotal = dict.found("thickness");
168  if (haveTotal)
169  {
170  thickness_ = scalarField
171  (
172  boundaryMesh.size(),
173  readScalar(dict.lookup("thickness"))
174  );
175  nSpec++;
176  }
177  bool haveExp = dict.found("expansionRatio");
178  if (haveExp)
179  {
180  expansionRatio_ = scalarField
181  (
182  boundaryMesh.size(),
183  readScalar(dict.lookup("expansionRatio"))
184  );
185  nSpec++;
186  }
187 
188 
189  if (haveFirst && haveTotal)
190  {
191  layerSpec_ = FIRST_AND_TOTAL;
192  Info<< "Layer thickness specified as first layer and overall thickness."
193  << endl;
194  }
195  else if (haveFirst && haveExp)
196  {
197  layerSpec_ = FIRST_AND_EXPANSION;
198  Info<< "Layer thickness specified as first layer and expansion ratio."
199  << endl;
200  }
201  else if (haveFinal && haveTotal)
202  {
203  layerSpec_ = FINAL_AND_TOTAL;
204  Info<< "Layer thickness specified as final layer and overall thickness."
205  << endl;
206  }
207  else if (haveFinal && haveExp)
208  {
209  layerSpec_ = FINAL_AND_EXPANSION;
210  Info<< "Layer thickness specified as final layer and expansion ratio."
211  << endl;
212  }
213  else if (haveTotal && haveExp)
214  {
215  layerSpec_ = TOTAL_AND_EXPANSION;
216  Info<< "Layer thickness specified as overall thickness"
217  << " and expansion ratio." << endl;
218  }
219 
220 
221  if (layerSpec_ == ILLEGAL || nSpec != 2)
222  {
224  (
225  dict
226  ) << "Over- or underspecified layer thickness."
227  << " Please specify" << nl
228  << " first layer thickness ('firstLayerThickness')"
229  << " and overall thickness ('thickness') or" << nl
230  << " first layer thickness ('firstLayerThickness')"
231  << " and expansion ratio ('expansionRatio') or" << nl
232  << " final layer thickness ('finalLayerThickness')"
233  << " and expansion ratio ('expansionRatio') or" << nl
234  << " final layer thickness ('finalLayerThickness')"
235  << " and overall thickness ('thickness') or" << nl
236  << " overall thickness ('thickness')"
237  << " and expansion ratio ('expansionRatio'"
238  << exit(FatalIOError);
239  }
240 
241 
242  dict.readIfPresent("nRelaxedIter", nRelaxedIter_);
243 
244  if (nLayerIter_ < 0 || nRelaxedIter_ < 0)
245  {
247  << "Layer iterations should be >= 0." << endl
248  << "nLayerIter:" << nLayerIter_
249  << " nRelaxedIter:" << nRelaxedIter_
250  << exit(FatalIOError);
251  }
252 
253 
254  const dictionary& layersDict = dict.subDict("layers");
255 
256  forAllConstIter(dictionary, layersDict, iter)
257  {
258  if (iter().isDict())
259  {
260  const keyType& key = iter().keyword();
261  const labelHashSet patchIDs
262  (
263  boundaryMesh.patchSet(List<wordRe>(1, wordRe(key)))
264  );
265 
266  if (patchIDs.size() == 0)
267  {
268  IOWarningInFunction(layersDict)
269  << "Layer specification for " << key
270  << " does not match any patch." << endl
271  << "Valid patches are " << boundaryMesh.names() << endl;
272  }
273  else
274  {
275  const dictionary& layerDict = iter().dict();
276 
277  forAllConstIter(labelHashSet, patchIDs, patchiter)
278  {
279  const label patchi = patchiter.key();
280 
281  numLayers_[patchi] =
282  readLabel(layerDict.lookup("nSurfaceLayers"));
283 
284  switch (layerSpec_)
285  {
286  case FIRST_AND_TOTAL:
287  layerDict.readIfPresent
288  (
289  "firstLayerThickness",
290  firstLayerThickness_[patchi]
291  );
292  layerDict.readIfPresent
293  (
294  "thickness",
295  thickness_[patchi]
296  );
297  break;
298 
299  case FIRST_AND_EXPANSION:
300  layerDict.readIfPresent
301  (
302  "firstLayerThickness",
303  firstLayerThickness_[patchi]
304  );
305  layerDict.readIfPresent
306  (
307  "expansionRatio",
308  expansionRatio_[patchi]
309  );
310  break;
311 
312  case FINAL_AND_TOTAL:
313  layerDict.readIfPresent
314  (
315  "finalLayerThickness",
316  finalLayerThickness_[patchi]
317  );
318  layerDict.readIfPresent
319  (
320  "thickness",
321  thickness_[patchi]
322  );
323  break;
324 
325  case FINAL_AND_EXPANSION:
326  layerDict.readIfPresent
327  (
328  "finalLayerThickness",
329  finalLayerThickness_[patchi]
330  );
331  layerDict.readIfPresent
332  (
333  "expansionRatio",
334  expansionRatio_[patchi]
335  );
336  break;
337 
338  case TOTAL_AND_EXPANSION:
339  layerDict.readIfPresent
340  (
341  "thickness",
342  thickness_[patchi]
343  );
344  layerDict.readIfPresent
345  (
346  "expansionRatio",
347  expansionRatio_[patchi]
348  );
349  break;
350 
351  default:
353  (
354  dict
355  ) << "problem." << exit(FatalIOError);
356  break;
357  }
358 
359  layerDict.readIfPresent
360  (
361  "minThickness",
362  minThickness_[patchi]
363  );
364  }
365  }
366  }
367  }
368 }
369 
370 
371 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
372 
374 (
375  const label nLayers,
376  const scalar firstLayerThickess,
377  const scalar finalLayerThickess,
378  const scalar totalThickness,
379  const scalar expansionRatio
380 ) const
381 {
382  switch (layerSpec_)
383  {
384  case FIRST_AND_TOTAL:
385  case FINAL_AND_TOTAL:
386  case TOTAL_AND_EXPANSION:
387  {
388  return totalThickness;
389  }
390  break;
391 
392  case FIRST_AND_EXPANSION:
393  {
394  if (mag(expansionRatio-1) < small)
395  {
396  return firstLayerThickess * nLayers;
397  }
398  else
399  {
400  return firstLayerThickess
401  *(1 - pow(expansionRatio, nLayers))
402  /(1 - expansionRatio);
403  }
404  }
405  break;
406 
407  case FINAL_AND_EXPANSION:
408  {
409  if (mag(expansionRatio-1) < small)
410  {
411  return finalLayerThickess * nLayers;
412  }
413  else
414  {
415  const scalar invExpansion = 1.0/expansionRatio;
416 
417  return finalLayerThickess
418  *(1 - pow(invExpansion, nLayers))
419  /(1 - invExpansion);
420  }
421  }
422  break;
423 
424  default:
425  {
427  << exit(FatalError);
428  return -vGreat;
429  }
430  }
431 }
432 
433 
434 Foam::scalar Foam::layerParameters::layerExpansionRatio
435 (
436  const label nLayers,
437  const scalar firstLayerThickess,
438  const scalar finalLayerThickess,
439  const scalar totalThickness,
440  const scalar expansionRatio
441 ) const
442 {
443  switch (layerSpec_)
444  {
445  case FIRST_AND_EXPANSION:
446  case FINAL_AND_EXPANSION:
447  case TOTAL_AND_EXPANSION:
448  {
449  return expansionRatio;
450  }
451  break;
452 
453  case FIRST_AND_TOTAL:
454  {
455  return layerExpansionRatio
456  (
457  nLayers,
458  totalThickness/firstLayerThickess
459  );
460  }
461  break;
462 
463  case FINAL_AND_TOTAL:
464  {
465  return
466  1.0
467  /layerExpansionRatio
468  (
469  nLayers,
470  totalThickness/finalLayerThickess
471  );
472  }
473  break;
474 
475  default:
476  {
478  << "Illegal thickness specification" << exit(FatalError);
479  return -vGreat;
480  }
481  }
482 }
483 
484 
486 (
487  const label nLayers,
488  const scalar firstLayerThickess,
489  const scalar finalLayerThickess,
490  const scalar totalThickness,
491  const scalar expansionRatio
492 ) const
493 {
494  switch (layerSpec_)
495  {
496  case FIRST_AND_EXPANSION:
497  case FIRST_AND_TOTAL:
498  {
499  return firstLayerThickess;
500  }
501 
502  case FINAL_AND_EXPANSION:
503  {
504  return finalLayerThickess*pow(1.0/expansionRatio, nLayers-1);
505  }
506  break;
507 
508  case FINAL_AND_TOTAL:
509  {
510  const scalar r = layerExpansionRatio
511  (
512  nLayers,
513  firstLayerThickess,
514  finalLayerThickess,
515  totalThickness,
516  expansionRatio
517  );
518 
519  return finalLayerThickess/pow(r, nLayers-1);
520  }
521  break;
522 
523  case TOTAL_AND_EXPANSION:
524  {
525  const scalar r = finalLayerThicknessRatio
526  (
527  nLayers,
528  expansionRatio
529  );
530 
531  const scalar finalThickness = r*totalThickness;
532 
533  return finalThickness/pow(expansionRatio, nLayers-1);
534  }
535  break;
536 
537  default:
538  {
540  << "Illegal thickness specification" << exit(FatalError);
541  return -vGreat;
542  }
543  }
544 }
545 
546 
548 (
549  const label nLayers,
550  const scalar expansionRatio
551 ) const
552 {
553  if (nLayers > 0)
554  {
555  if (mag(expansionRatio-1) < small)
556  {
557  return 1.0/nLayers;
558  }
559  else
560  {
561  return
562  pow(expansionRatio, nLayers - 1)
563  *(1 - expansionRatio)
564  /(1 - pow(expansionRatio, nLayers));
565  }
566  }
567  else
568  {
569  return 0;
570  }
571 }
572 
573 
574 // ************************************************************************* //
A class for handling keywords in dictionaries.
Definition: keyType.H:64
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Unit conversion functions.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
scalar finalLayerThicknessRatio(const label nLayers, const scalar expansionRatio) const
Determine ratio of final layer thickness to.
scalar layerThickness(const label nLayers, const scalar firstLayerThickess, const scalar finalLayerThickess, const scalar totalThickness, const scalar expansionRatio) const
Determine overall thickness. Uses two of the four parameters.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
wordList names() const
Return a list of patch names.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
const scalarField & firstLayerThickness() const
Wanted thickness of the layer nearest to the wall.
static const label labelMax
Definition: label.H:62
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:74
Foam::polyBoundaryMesh.
static const char nl
Definition: Ostream.H:265
const scalarField & expansionRatio() const
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
IOerror FatalIOError