setAlphaField.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) 2016-2017 DHI
9  Copyright (C) 2017-2024 OpenCFD Ltd.
10  Copyright (C) 2017-2020 German Aerospace Center (DLR)
11  Copyright (C) 2020 Johan Roenby
12 -------------------------------------------------------------------------------
13 License
14  This file is part of OpenFOAM.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 Application
30  setAlphaField
31 
32 Description
33  Uses cutCellIso to create a volume fraction field from either a cylinder,
34  a sphere or a plane.
35 
36  Original code supplied by Johan Roenby, DHI (2016)
37  Modification Henning Scheufler, DLR (2019)
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
42 #include "timeSelector.H"
43 
44 #include "triSurface.H"
45 #include "triSurfaceTools.H"
46 
47 #include "implicitFunction.H"
48 
49 #include "cutCellIso.H"
50 #include "cutFaceIso.H"
51 #include "OBJstream.H"
52 
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 void isoFacesToFile
57 (
58  const UList<List<point>>& faces,
59  const word& fileName
60 )
61 {
62  // Writing isofaces to OBJ file for inspection in paraview
63  OBJstream os(fileName + ".obj");
64 
65  if (Pstream::parRun())
66  {
67  // Collect points from all the processors
68  List<List<List<point>>> allProcFaces(Pstream::nProcs());
69  allProcFaces[Pstream::myProcNo()] = faces;
70  Pstream::gatherList(allProcFaces);
71 
72  if (Pstream::master())
73  {
74  Info<< "Writing file: " << fileName << endl;
75 
76  for (const auto& procFaces : allProcFaces)
77  {
78  for (const auto& facePts : procFaces)
79  {
80  os.writeFace(facePts);
81  }
82  }
83  }
84  }
85  else
86  {
87  Info<< "Writing file: " << fileName << endl;
88 
89  for (const auto& facePts : faces)
90  {
91  os.writeFace(facePts);
92  }
93  }
94 }
95 
96 void setAlpha
97 (
99  DynamicList<List<point>>& facePts,
100  scalarField& f,
101  const bool writeOBJ
102 )
103 {
104  const fvMesh& mesh = alpha1.mesh();
105  cutCellIso cutCell(mesh, f);
106  cutFaceIso cutFace(mesh, f);
107 
108  forAll(alpha1, cellI)
109  {
110  cutCell.calcSubCell(cellI, 0.0);
111 
112  alpha1[cellI] = clamp(cutCell.VolumeOfFluid(), zero_one{});
113 
114  if (writeOBJ && (mag(cutCell.faceArea()) >= 1e-14))
115  {
116  facePts.append(cutCell.facePoints());
117  }
118  }
119 
120  // Setting boundary alpha1 values
121  forAll(mesh.boundary(), patchi)
122  {
123  if (mesh.boundary()[patchi].size() > 0)
124  {
125  const label start = mesh.boundary()[patchi].patch().start();
126  scalarField& alphap = alpha1.boundaryFieldRef()[patchi];
127  const scalarField& magSfp = mesh.magSf().boundaryField()[patchi];
128 
129  forAll(alphap, patchFacei)
130  {
131  const label facei = patchFacei + start;
132  cutFace.calcSubFace(facei, 0.0);
133  alphap[patchFacei] =
134  mag(cutFace.subFaceArea())/magSfp[patchFacei];
135  }
136  }
137  }
138 }
139 
140 
141 int main(int argc, char *argv[])
142 {
143  argList::noFunctionObjects(); // Never use function objects
144 
145  timeSelector::addOptions_singleTime(); // Single-time options
146 
147  argList::addNote
148  (
149  "Uses cutCellIso to create a volume fraction field from an "
150  "implicit function."
151  );
152 
153  argList::addOption
154  (
155  "dict",
156  "file",
157  "Alternative setAlphaFieldDict dictionary"
158  );
159 
160  #include "addRegionOption.H"
161  #include "setRootCase.H"
162  #include "createTime.H"
163 
164  // Set time from specified time options, or no-op
165  timeSelector::setTimeIfPresent(runTime, args);
166 
167  #include "createNamedMesh.H"
168 
169  const word dictName("setAlphaFieldDict");
170  #include "setSystemMeshDictionaryIO.H"
171 
172  IOdictionary setAlphaFieldDict(dictIO);
173 
174  Info<< "Reading " << setAlphaFieldDict.name() << endl;
175 
176  const word fieldName = setAlphaFieldDict.get<word>("field");
177  const bool invert = setAlphaFieldDict.getOrDefault("invert", false);
178  const bool writeOBJ = setAlphaFieldDict.getOrDefault("writeOBJ", false);
179 
180  Info<< "Reading field " << fieldName << nl << endl;
182  (
183  IOobject
184  (
185  fieldName,
186  runTime.timeName(),
187  mesh,
188  IOobject::MUST_READ,
189  IOobject::AUTO_WRITE
190  ),
191  mesh
192  );
193 
194  autoPtr<implicitFunction> func = implicitFunction::New
195  (
196  setAlphaFieldDict.get<word>("type"),
197  setAlphaFieldDict
198  );
199 
200  scalarField f(mesh.nPoints(), Zero);
201 
202  forAll(f, pi)
203  {
204  f[pi] = func->value(mesh.points()[pi]);
205  };
206 
207  DynamicList<List<point>> facePts;
208 
209  setAlpha(alpha1, facePts, f, writeOBJ);
210 
211  if (writeOBJ)
212  {
213  isoFacesToFile(facePts, fieldName + "0");
214  }
215 
216  ISstream::defaultPrecision(18);
217 
218  if (invert)
219  {
220  alpha1 = scalar(1) - alpha1;
221  }
222 
223  Info<< "Writing new alpha field " << alpha1.name() << endl;
224  alpha1.write();
225 
226  {
227  const auto& alpha = alpha1.primitiveField();
228 
229  auto limits = gMinMax(alpha);
230  auto total = gWeightedSum(mesh.V(), alpha);
231 
232  Info<< "sum(alpha*V):" << total
233  << ", 1-max(alpha1): " << 1 - limits.max()
234  << " min(alpha1): " << limits.min() << endl;
235  }
236 
237  Info<< "\nEnd\n" << endl;
238 
239  return 0;
240 }
241 
242 
243 // ************************************************************************* //
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const word dictName("faMeshDefinition")
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Required Classes.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
auto limits
Definition: setRDeltaT.H:186
Required Classes.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar pi(M_PI)
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
labelList f(nPoints)
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:28
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject dictIO
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Foam::argList args(argc, argv)
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1