MapFieldConstraint.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) 2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::fv::MapFieldConstraint
28 
29 Group
30  grpFvOptionsConstraints
31 
32 Description
33  The \c MapFieldConstraint constrains values of given fields of \c Type
34  with a source field from an external mesh, where
35  \c <Type>={scalar,vector,sphericalTensor,symmTensor,tensor}.
36  Optionally, the source field can be translated and/or rotated as a function
37  of time.
38 
39 Usage
40  Minimal example by using \c constant/fvOptions:
41  \verbatim
42  <Type>MapFieldConstraint1
43  {
44  // Mandatory entries
45  type <Type>MapFieldConstraint;
46  field <word>;
47  srcMesh <fileName>;
48  mapMethod <word>;
49 
50  // Optional entries
51  consistent <bool>;
52  patchMapMethod <word>;
53  transform
54  {
55  // Optional entries
56  position <Function1<vector>>;
57  origin <vector>;
58 
59  direction <Function1<vector>>;
60  normal <vector>;
61  }
62 
63  // Conditional entries
64 
65  // when consistent=false
66  patchMap <HashTable<word>>; // (<patchSrc> <patchTgt>);
67  cuttingPatches <wordList>; // (<patchTgt1> ... <patchTgtN>);
68 
69  // Inherited entries
70  ...
71  }
72  \endverbatim
73 
74  where the entries mean:
75  \table
76  Property | Description | Type | Reqd | Deflt
77  type | Type name: <Type>MapFieldConstraint | word | yes | -
78  field | Name of operand field | word | yes | -
79  srcMesh | Directory path to mesh to map from | fileName | yes | -
80  mapMethod | Mapping method | word | yes | -
81  consistent | Flag to determine if meshes have consistent boundaries <!--
82  --> | bool | no | false
83  patchMapMethod | Name of patch-map method | word | no | -
84  patchMap | Coincident source/target patches in two cases <!--
85  --> | wordHashTable | no | -
86  cuttingPatches | Target patches cutting the source domain <!--
87  --> | wordList | no | -
88  transform | Transform settings for source mesh points <!--
89  --> | dict | no | -
90  position | Position of source mesh as a function of time <!--
91  --> | Function1<vector> | no | -
92  direction | Direction of source mesh as a function of time <!--
93  --> | Function1<vector> | no | -
94  origin | Origin of source mesh | vector | no | -
95  normal | Normal of reference plane representing source mesh <!--
96  --> | vector | no | -
97  \endtable
98 
99  The inherited entries are elaborated in:
100  - \link fvOption.H \endlink
101  - \link Function1.H \endlink
102 
103 SourceFiles
104  MapFieldConstraint.C
105 
106 \*---------------------------------------------------------------------------*/
107 
108 #ifndef Foam_fv_MapFieldConstraint_H
109 #define Foam_fv_MapFieldConstraint_H
110 
111 #include "fvOption.H"
112 #include "fvMesh.H"
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
115 
116 namespace Foam
117 {
118 
119 // Forward Declarations
120 class meshToMesh;
121 template<class Type> class Function1;
122 
123 namespace fv
124 {
125 
126 /*---------------------------------------------------------------------------*\
127  Class MapFieldConstraint Declaration
128 \*---------------------------------------------------------------------------*/
129 
130 template<class Type>
131 class MapFieldConstraint
132 :
133  public fv::option
134 {
135  // Private Classes
136 
137  class transform
138  {
139  // Private Data
140 
141  //- Position of source mesh as a function of time
142  autoPtr<Function1<point>> positionPtr_;
143 
144  //- Direction of source mesh as a function of time
145  autoPtr<Function1<point>> directionPtr_;
146 
147  //- Cached points of source mesh
148  pointField points_;
149 
150  //- Origin of source mesh
151  point origin_;
152 
153  //- Normal of reference plane representing source mesh
154  vector normal_;
155 
156  //- Flag to deduce if transformation is active
157  bool active_;
158 
159 
160  public:
161 
162  // Constructors
163 
164  //- Default construct
165  transform();
166 
167  //- No copy construct
168  transform(const transform&) = delete;
169 
170  //- No copy assignment
171  void operator=(const transform&) = delete;
172 
173 
174  // Member Functions
175 
176  // Access
177 
178  //- Return flag to deduce if transformation is active
179  bool isActive() const noexcept { return active_; }
180 
181 
182  // Evaluation
183 
184  //- Translate source mesh as a function of time
185  void translate(refPtr<fvMesh>& srcMeshPtr, const scalar time);
186 
187  //- Rotate source mesh as a function of time
188  void rotate(refPtr<fvMesh>& srcMeshPtr, const scalar time);
189 
190 
191  // I-O
192 
193  //- Initialize the class members
194  bool initialize(const fvMesh& srcMesh, const dictionary& dict);
195  };
196 
197 
198  // Private Data
199 
200  //- Transformation settings for source mesh
201  transform transform_;
202 
203  //- Time database for source mesh to map from
204  autoPtr<Time> srcTimePtr_;
205 
206  //- Source mesh to map from
207  refPtr<fvMesh> srcMeshPtr_;
208 
209  //- Mesh-to-mesh interpolation from source mesh to target mesh
210  autoPtr<meshToMesh> interpPtr_;
211 
212  //- List of coincident source/target patches in two cases
213  HashTable<word> patchMap_;
214 
215  //- Set of cells to apply source to
216  labelList cells_;
217 
218  //- List of names of target patches cutting the source domain
219  wordList cuttingPatches_;
220 
221  //- Name of map method
222  word mapMethodName_;
223 
224  //- Name of patch-map method
225  word patchMapMethodName_;
226 
227  //- Flag to determine if meshes have consistent boundaries
228  bool consistent_;
229 
230 
231  // Private Member Functions
232 
233  //- Helper function to set source mesh
234  // Fetch fvMesh from a given Time database
235  // Otherwise, load it from disk and cache it to the database
236  void setSourceMesh
237  (
238  refPtr<fvMesh>& meshRef,
240  );
241 
242  //- Helper function to create the mesh-to-mesh interpolation
243  void createInterpolation
244  (
245  const fvMesh& srcMesh,
246  const fvMesh& tgtMesh
247  );
248 
249  //- Return requested field from object registry
250  //- otherwise read it from disk and register it to the object registry
251  template<class VolFieldType>
252  VolFieldType& getOrReadField
253  (
254  const fvMesh& thisMesh,
255  const word& fieldName
256  ) const;
257 
258  //- Return the local cell indices of the target mesh
259  labelList tgtCellIDs() const;
260 
261 
262 public:
263 
264  //- Runtime type information
265  TypeName("MapFieldConstraint");
266 
267 
268  // Constructors
269 
270  //- Construct from components
272  (
273  const word& name,
274  const word& modelType,
275  const dictionary& dict,
276  const fvMesh& mesh
277  );
278 
279  //- No copy construct
280  MapFieldConstraint(const MapFieldConstraint&) = delete;
281 
282  //- No copy assignment
283  void operator=(const MapFieldConstraint&) = delete;
284 
285 
286  //- Destructor
287  virtual ~MapFieldConstraint() = default;
288 
289 
290  // Member Functions
291 
292  //- Read source dictionary
293  virtual bool read(const dictionary& dict);
294 
295  //- Set value on field
296  virtual void constrain(fvMatrix<Type>& eqn, const label);
297 };
298 
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 } // End namespace fv
303 } // End namespace Foam
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 #ifdef NoRepository
308  #include "MapFieldConstraint.C"
309 #endif
310 
311 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
312 
313 #endif
314 
315 // ************************************************************************* //
dictionary dict
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition: fvOptionI.H:30
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void operator=(const MapFieldConstraint &)=delete
No copy assignment.
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
virtual ~MapFieldConstraint()=default
Destructor.
Info<< "Create engine time\"<< endl;autoPtr< engineTime > runTimePtr(engineTime::New(Time::controlDictName, args.rootPath(), args.caseName()))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
const word & name() const noexcept
Return const access to the source name.
Definition: fvOptionI.H:24
Vector< scalar > vector
Definition: vector.H:57
const direction noexcept
Definition: Scalar.H:258
The MapFieldConstraint constrains values of given fields of Type with a source field from an external...
List< word > wordList
List of word.
Definition: fileName.H:59
vector point
Point is a vector.
Definition: point.H:37
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual bool isActive()
Is the source active?
Definition: fvOption.C:112
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
MapFieldConstraint(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
List< label > labelList
A List of labels.
Definition: List.H:62
virtual bool read(const dictionary &dict)
Read source dictionary.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
virtual void constrain(fvMatrix< Type > &eqn, const label)
Set value on field.
bool active_
Source active flag.
Definition: fvOption.H:167
TypeName("MapFieldConstraint")
Runtime type information.
Namespace for OpenFOAM.