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-2020 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 
43 #include "triSurface.H"
44 #include "triSurfaceTools.H"
45 
46 #include "implicitFunction.H"
47 
48 #include "cutCellIso.H"
49 #include "cutFaceIso.H"
50 #include "OBJstream.H"
51 
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 void isoFacesToFile
56 (
57  const DynamicList<List<point>>& faces,
58  const word& fileName
59 )
60 {
61  // Writing isofaces to OBJ file for inspection in paraview
62  OBJstream os(fileName + ".obj");
63 
64  if (Pstream::parRun())
65  {
66  // Collect points from all the processors
67  List<DynamicList<List<point>>> allProcFaces(Pstream::nProcs());
68  allProcFaces[Pstream::myProcNo()] = faces;
69  Pstream::gatherList(allProcFaces);
70 
71  if (Pstream::master())
72  {
73  Info<< "Writing file: " << fileName << endl;
74 
75  for (const DynamicList<List<point>>& procFaces : allProcFaces)
76  {
77  for (const List<point>& facePts : procFaces)
78  {
79  os.write(face(identity(facePts.size())), facePts);
80  }
81  }
82  }
83  }
84  else
85  {
86  Info<< "Writing file: " << fileName << endl;
87 
88  for (const List<point>& facePts : faces)
89  {
90  os.write(face(identity(facePts.size())), facePts);
91  }
92  }
93 }
94 
95 void setAlpha
96 (
98  DynamicList<List<point>>& facePts,
99  scalarField& f,
100  const bool writeOBJ
101 )
102 {
103  const fvMesh& mesh = alpha1.mesh();
104  cutCellIso cutCell(mesh, f);
105  cutFaceIso cutFace(mesh, f);
106 
107  forAll(alpha1, cellI)
108  {
109  cutCell.calcSubCell(cellI, 0.0);
110 
111  alpha1[cellI] = max(min(cutCell.VolumeOfFluid(), 1), 0);
112 
113  if (writeOBJ && (mag(cutCell.faceArea()) >= 1e-14))
114  {
115  facePts.append(cutCell.facePoints());
116  }
117  }
118 
119  // Setting boundary alpha1 values
120  forAll(mesh.boundary(), patchi)
121  {
122  if (mesh.boundary()[patchi].size() > 0)
123  {
124  const label start = mesh.boundary()[patchi].patch().start();
125  scalarField& alphap = alpha1.boundaryFieldRef()[patchi];
126  const scalarField& magSfp = mesh.magSf().boundaryField()[patchi];
127 
128  forAll(alphap, patchFacei)
129  {
130  const label facei = patchFacei + start;
131  cutFace.calcSubFace(facei, 0.0);
132  alphap[patchFacei] =
133  mag(cutFace.subFaceArea())/magSfp[patchFacei];
134  }
135  }
136  }
137 }
138 
139 
140 int main(int argc, char *argv[])
141 {
142  argList::addNote
143  (
144  "Uses cutCellIso to create a volume fraction field from an "
145  "implicit function."
146  );
147 
148  #include "addDictOption.H"
149  #include "addRegionOption.H"
150  #include "setRootCase.H"
151  #include "createTime.H"
152  #include "createNamedMesh.H"
153 
154  const word dictName("setAlphaFieldDict");
155  #include "setSystemMeshDictionaryIO.H"
156 
157  IOdictionary setAlphaFieldDict(dictIO);
158 
159  Info<< "Reading " << setAlphaFieldDict.name() << endl;
160 
161  const word fieldName = setAlphaFieldDict.get<word>("field");
162  const bool invert = setAlphaFieldDict.getOrDefault("invert", false);
163  const bool writeOBJ = setAlphaFieldDict.getOrDefault("writeOBJ", false);
164 
165  Info<< "Reading field " << fieldName << nl << endl;
167  (
168  IOobject
169  (
170  fieldName,
171  runTime.timeName(),
172  mesh,
173  IOobject::MUST_READ,
174  IOobject::AUTO_WRITE
175  ),
176  mesh
177  );
178 
179  autoPtr<implicitFunction> func = implicitFunction::New
180  (
181  setAlphaFieldDict.get<word>("type"),
182  setAlphaFieldDict
183  );
184 
185  scalarField f(mesh.nPoints(), Zero);
186 
187  forAll(f, pi)
188  {
189  f[pi] = func->value(mesh.points()[pi]);
190  };
191 
192  DynamicList<List<point>> facePts;
193 
194  setAlpha(alpha1, facePts, f, writeOBJ);
195 
196  if (writeOBJ)
197  {
198  isoFacesToFile(facePts, fieldName + "0");
199  }
200 
201  ISstream::defaultPrecision(18);
202 
203  if (invert)
204  {
205  alpha1 = scalar(1) - alpha1;
206  }
207 
208  Info<< "Writing new alpha field " << alpha1.name() << endl;
209  alpha1.write();
210 
211  const scalarField& alpha = alpha1.internalField();
212 
213  Info<< "sum(alpha*V):" << gSum(mesh.V()*alpha)
214  << ", 1-max(alpha1): " << 1 - gMax(alpha)
215  << " min(alpha1): " << gMin(alpha) << endl;
216 
217  Info<< "\nEnd\n" << endl;
218 
219  return 0;
220 }
221 
222 
223 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
dictName
const word dictName("faMeshDefinition")
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
alpha1
const volScalarField & alpha1
Definition: setRegionFluidFields.H:8
triSurface.H
cutFaceIso.H
setSystemMeshDictionaryIO.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::invert
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
addDictOption.H
cutCellIso.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
addRegionOption.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::func
void func(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
createNamedMesh.H
Required Variables.
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
dictIO
IOobject dictIO
Definition: setConstantMeshDictionaryIO.H:1
implicitFunction.H
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
setRootCase.H
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
f
labelList f(nPoints)
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
createTime.H
fvCFD.H
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
triSurfaceTools.H
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
OBJstream.H