createZeroField.H
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) 2007-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 
29 \*---------------------------------------------------------------------------*/
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 #ifndef createZeroField_H
34 #define createZeroField_H
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class Type>
44 (
45  const fvMesh& mesh,
46  const word& name,
47  const dimensionSet dims,
48  bool printAllocation = false
49 )
50 {
51  if (printAllocation)
52  {
53  Info<< "Allocating new volField " << name << nl << endl;
54  }
55 
56  return
57  (
59  (
60  IOobject
61  (
62  name,
63  mesh.time().timeName(),
64  mesh,
67  ),
68  mesh,
69  dimensioned<Type>(dims, Zero)
70  )
71  );
72 }
73 
74 
75 template<class Type>
76 autoPtr<typename GeometricField<Type, fvPatchField, volMesh>::Boundary>
78 (
79  const fvMesh& mesh,
80  bool printAllocation = false
81 )
82 {
83  if (printAllocation)
84  {
85  Info<< "Allocating new boundaryField " << nl << endl;
86  }
87 
89  Boundary;
90 
91  // Make sure that the patchFields to be generated will be of type
92  // calculated, even if they are of constraint type.
93  // Necessary to avoid unexpected behaviour when computing sensitivities
94  // on symmetry patches (not a good practice either way)
95  const fvBoundaryMesh& bm = mesh.boundary();
96  wordList actualPatchTypes(bm.size(), word::null);
97  forAll(actualPatchTypes, patchi)
98  {
99  if
100  (
101  fvPatchField<Type>::patchConstructorTablePtr_->contains
102  (
103  bm[patchi].type()
104  )
105  )
106  {
107  actualPatchTypes[patchi] = bm[patchi].type();
108  }
109  }
110 
111  autoPtr<Boundary> bPtr
112  (
113  new Boundary
114  (
115  mesh.boundary(),
118  actualPatchTypes
119  )
120  );
121 
122  // Values are not assigned! Assign manually
123  Boundary& bRef = bPtr();
124  forAll(bRef, pI)
125  {
126  bRef[pI] = pTraits<Type>::zero;
127  }
128 
129  return (bPtr);
130 }
131 
132 
133 template<class Type>
134 autoPtr<List<Field<Type>>>
136 (
137  const fvMesh& mesh,
138  bool printAllocation = false
139 )
140 {
141  if (printAllocation)
142  {
143  Info<< "Allocating new point boundaryField " << nl << endl;
144  }
145 
146  autoPtr<List<Field<Type>>> bPtr
147  (
148  new List<Field<Type>>
149  (
150  mesh.boundary().size()
151  )
152  );
153 
154  List<Field<Type>>& bRef = bPtr();
155  forAll(bRef, pI)
156  {
157  bRef[pI] =
158  Field<Type>
159  (
160  mesh.boundaryMesh()[pI].nPoints(),
161  pTraits<Type>::zero
162  );
163  }
164 
165  return (bPtr);
166 }
167 
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 } // End namespace Foam
172 
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
174 
175 #endif
176 
177 // ************************************************************************* //
autoPtr< List< Field< Type > > > createZeroBoundaryPointFieldPtr(const fvMesh &mesh, bool printAllocation=false)
static const DimensionedField< Type, GeoMesh > & null() noexcept
Return a null DimensionedField (reference to a nullObject).
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of boundary fields.
Generic dimensioned Type class.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
#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
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
static const word null
An empty word.
Definition: word.H:84
autoPtr< GeometricField< Type, fvPatchField, volMesh > > createZeroFieldPtr(const fvMesh &mesh, const word &name, const dimensionSet dims, bool printAllocation=false)
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
List< word > wordList
List of word.
Definition: fileName.H:59
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
messageStream Info
Information stream (stdout output on master, null elsewhere)
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127