exactPatchDistMethod.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) 2018-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "exactPatchDistMethod.H"
31 #include "volFields.H"
32 #include "DynamicField.H"
33 #include "OBJstream.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace patchDistMethods
40 {
41  defineTypeNameAndDebug(exact, 0);
42  addToRunTimeSelectionTable(patchDistMethod, exact, dictionary);
43 }
44 }
45 
47 Foam::patchDistMethods::exact::patchSurface() const
48 {
49  if (!patchSurfPtr_)
50  {
51  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
52 
53  Random rndGen(0);
54 
55  treeBoundBox localBb(mesh_.points());
56 
57  // Determine mesh bounding boxes:
58  List<treeBoundBox> meshBb
59  (
60  1,
61  localBb.extend(rndGen, 1E-3)
62  );
63 
64  // Dummy bounds dictionary
66  dict.add("bounds", meshBb);
67  dict.add
68  (
69  "distributionType",
71  [
72  //distributedTriSurfaceMesh::FOLLOW // use mesh bb
73  //distributedTriSurfaceMesh::INDEPENDENT // master-only
75  ]
76  );
77  dict.add("mergeDistance", 1e-6*localBb.mag());
78 
79 
80  Info<< "Triangulating local patch faces" << nl << endl;
81 
82  labelList mapTriToGlobal;
83 
84  patchSurfPtr_.reset
85  (
86  new distributedTriSurfaceMesh
87  (
88  IOobject
89  (
90  "wallSurface.stl",
91  mesh_.time().constant(),// directory
92  "triSurface", // instance
93  mesh_.time(), // registry
96  ),
98  (
99  pbm,
100  patchIDs_,
101  mapTriToGlobal
102  ),
103  dict
104  )
105  );
106 
107  // Do redistribution
108  Info<< "Redistributing surface" << nl << endl;
109  autoPtr<mapDistribute> faceMap;
110  autoPtr<mapDistribute> pointMap;
111  patchSurfPtr_().distribute
112  (
113  meshBb,
114  false, //keepNonMapped,
115  faceMap,
116  pointMap
117  );
118  faceMap.clear();
119  pointMap.clear();
120  }
121  return patchSurfPtr_();
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
127 Foam::patchDistMethods::exact::exact
128 (
129  const dictionary& dict,
130  const fvMesh& mesh,
131  const labelHashSet& patchIDs
132 )
133 :
135 {}
136 
137 
138 Foam::patchDistMethods::exact::exact
139 (
140  const fvMesh& mesh,
141  const labelHashSet& patchIDs
142 )
143 :
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
151 {
152  return correct(y, const_cast<volVectorField&>(volVectorField::null()));
153 }
154 
155 
157 (
158  volScalarField& y,
160 )
161 {
162  const distributedTriSurfaceMesh& surf = patchSurface();
163 
164  List<pointIndexHit> info;
165  surf.findNearest
166  (
167  mesh_.cellCentres(),
168  scalarField(mesh_.nCells(), Foam::sqr(GREAT)),
169  info
170  );
171 
172  // Take over hits
173  // label nHits = 0;
174  forAll(info, celli)
175  {
176  if (info[celli].hit())
177  {
178  const point& cc = mesh_.cellCentres()[celli];
179  y[celli] = info[celli].point().dist(cc);
180  // ++nHits;
181  }
182  //else
183  //{
184  // // Miss. Do what? Not possible with GREAT hopefully ...
185  //}
186  }
187  y.correctBoundaryConditions();
188 
189  if (debug)
190  {
191  OBJstream str(mesh_.time().timePath()/"wallPoint.obj");
192  Info<< type() << ": dumping nearest wall point to "
193  << str.name() << endl;
194  forAll(mesh_.cellCentres(), celli)
195  {
196  const point& cc = mesh_.cellCentres()[celli];
197  str.writeLine(cc, info[celli].point());
198  }
199  }
200 
201 
202  // Only calculate n if the field is defined
203  if (notNull(n))
204  {
205  surf.getNormal(info, n.primitiveFieldRef());
206  n.correctBoundaryConditions();
207  }
208 
209  return true;
210 }
211 
212 
213 // ************************************************************************* //
static const Enum< distributionType > distributionTypeNames_
const labelList patchIDs(pbm.indices(polyPatchNames, true))
const polyBoundaryMesh & pbm
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, labelList &faceMap, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Macros for easy insertion into run-time selection tables.
const fvMesh & mesh_
Reference to the mesh.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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
scalar y
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
IOoject and searching on distributed triSurface. All processor hold (possibly overlapping) part of th...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
int debug
Static debugging option.
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
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/v2312/OpenFOAM-v2312/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
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
vector point
Point is a vector.
Definition: point.H:37
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
defineTypeNameAndDebug(advectionDiffusion, 0)
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Nothing to be read.
const labelHashSet patchIDs_
Set of patch IDs.
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
List< label > labelList
A List of labels.
Definition: List.H:62
Specialisation of patchDist for wall distance calculation.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:246
Namespace for OpenFOAM.