Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
OpenFOAM
fields
GeometricFields
GeometricField
MapGeometricFields.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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
Class
25
Foam::MapInternalField
26
27
Description
28
Generic internal field mapper. For "real" mapping, add template
29
specialisations for mapping of internal fields depending on mesh
30
type.
31
32
\*---------------------------------------------------------------------------*/
33
34
#ifndef MapGeometricFields_H
35
#define MapGeometricFields_H
36
37
#include <
OpenFOAM/polyMesh.H
>
38
39
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40
41
namespace
Foam
42
{
43
44
template
<
class
Type,
class
MeshMapper,
class
GeoMesh>
45
class
MapInternalField
46
{
47
public
:
48
49
MapInternalField
()
50
{}
51
52
void
operator()
53
(
54
Field<Type>
& field,
55
const
MeshMapper& mapper
56
)
const
;
57
};
58
59
60
//- Generic Geometric field mapper.
61
// For "real" mapping, add template specialisations
62
// for mapping of internal fields depending on mesh type.
63
template
64
<
65
class
Type,
66
template
<
class
>
class
PatchField,
67
class
MeshMapper,
68
class
GeoMesh
69
>
70
void
MapGeometricFields
71
(
72
const
MeshMapper& mapper
73
)
74
{
75
HashTable<const GeometricField<Type, PatchField, GeoMesh>
*>
fields
76
(
77
mapper.thisDb().objectRegistry::lookupClass
78
<
GeometricField<Type, PatchField, GeoMesh>
>()
79
);
80
81
// It is necessary to enforce that all old-time fields are stored
82
// before the mapping is performed. Otherwise, if the
83
// old-time-level field is mapped before the field itself, sizes
84
// will not match.
85
86
for
87
(
88
typename
HashTable
<
const
GeometricField<Type, PatchField, GeoMesh>
*>::
89
iterator fieldIter = fields.
begin
();
90
fieldIter != fields.
end
();
91
++fieldIter
92
)
93
{
94
GeometricField<Type, PatchField, GeoMesh>
& field =
95
const_cast<
GeometricField<Type, PatchField, GeoMesh>
&
>
96
(*fieldIter());
97
98
//Note: check can be removed once pointFields are actually stored on
99
// the pointMesh instead of now on the polyMesh!
100
if
(&field.
mesh
() == &mapper.mesh())
101
{
102
field.
storeOldTimes
();
103
}
104
}
105
106
for
107
(
108
typename
HashTable
<
const
GeometricField<Type, PatchField, GeoMesh>
*>::
109
iterator fieldIter = fields.
begin
();
110
fieldIter != fields.
end
();
111
++fieldIter
112
)
113
{
114
GeometricField<Type, PatchField, GeoMesh>
& field =
115
const_cast<
GeometricField<Type, PatchField, GeoMesh>
&
>
116
(*fieldIter());
117
118
if
(&field.
mesh
() == &mapper.mesh())
119
{
120
if
(polyMesh::debug)
121
{
122
Info
<<
"Mapping "
<< field.
typeName
<<
' '
<< field.
name
()
123
<<
endl
;
124
}
125
126
// Map the internal field
127
MapInternalField<Type, MeshMapper, GeoMesh>
()
128
(
129
field.
internalField
(),
130
mapper
131
);
132
133
// Map the patch fields
134
forAll
(field.
boundaryField
(),
patchi
)
135
{
136
// Cannot check sizes for patch fields because of
137
// empty fields in FV and because point fields get their size
138
// from the patch which has already been resized
139
//
140
141
field.
boundaryField
()[
patchi
].autoMap
142
(
143
mapper.boundaryMap()[
patchi
]
144
);
145
}
146
147
field.
instance
() = field.
time
().
timeName
();
148
}
149
else
if
(polyMesh::debug)
150
{
151
Info
<<
"Not mapping "
<< field.
typeName
<<
' '
<< field.
name
()
152
<<
" since originating mesh differs from that of mapper."
153
<<
endl
;
154
}
155
}
156
}
157
158
159
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160
161
}
// End namespace Foam
162
163
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164
165
#endif
166
167
// ************************ vim: set sw=4 sts=4 et: ************************ //