Test-mesh.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 | Copyright (C) 2011-2015 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 "argList.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 using namespace Foam;
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 // Main program:
35 
36 int main(int argc, char *argv[])
37 {
38 
39  #include "setRootCase.H"
40  #include "createTime.H"
41 
42  Info<< "Create mesh, no clear-out\n" << endl;
43  fvMesh mesh
44  (
45  IOobject
46  (
48  runTime.timeName(),
49  runTime,
51  )
52  );
53 
54  Info<< mesh.C() << endl;
55  Info<< mesh.V() << endl;
56 
57  surfaceVectorField Cf = mesh.Cf();
58 
59  Info<< Cf << endl;
60 
61  // Test construct from cellShapes
62  {
64  cellShapeList shapes(mesh.cellShapes());
65 
66  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
67 
68  faceListList boundaryFaces(pbm.size());
69  forAll(pbm, patchI)
70  {
71  boundaryFaces[patchI] = pbm[patchI];
72  }
73  wordList boundaryPatchNames(pbm.names());
74  PtrList<dictionary> boundaryDicts(pbm.size());
75  forAll(pbm, patchI)
76  {
77  OStringStream os;
78  os << pbm[patchI];
79  IStringStream is(os.str());
80  boundaryDicts.set(patchI, new dictionary(is));
81  }
82 
83  word defaultBoundaryPatchName = "defaultFaces";
84  word defaultBoundaryPatchType = emptyPolyPatch::typeName;
85 
86  fvMesh newMesh
87  (
88  IOobject
89  (
90  "newMesh",
91  runTime.timeName(),
92  runTime,
94  ),
96  shapes,
97  boundaryFaces,
98  boundaryPatchNames,
99  boundaryDicts,
100  defaultBoundaryPatchName,
101  defaultBoundaryPatchType
102  );
103 
104  Info<< newMesh.C() << endl;
105  Info<< newMesh.V() << endl;
106 
107  surfaceVectorField Cf = newMesh.Cf();
108 
109  Info<< Cf << endl;
110  }
111 
112 
113  Info<< "End\n" << endl;
114 
115  return 0;
116 }
117 
118 
119 // ************************************************************************* //
volFields.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fvMesh::Cf
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Definition: fvMeshGeometry.C:380
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::polyMesh::defaultRegion
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:527
Foam::OStringStream::str
string str() const
Return the string.
Definition: OStringStream.H:107
Foam::PtrList::set
bool set(const label) const
Is element set.
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:199
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
argList.H
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::primitiveMesh::cellShapes
const cellShapeList & cellShapes() const
Return cell shapes.
Definition: primitiveMesh.C:230
setRootCase.H
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::List< cellShape >
Foam::OStringStream
Output to memory buffer stream.
Definition: OStringStream.H:49
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
createTime.H
main
int main(int argc, char *argv[])
Definition: Test-mesh.C:36
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52