MapFieldConstraint.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) 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 \*---------------------------------------------------------------------------*/
27 
28 #include "MapFieldConstraint.H"
29 #include "fvMatrices.H"
30 #include "meshToMesh.H"
31 #include "Function1.H"
32 
33 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
39 static inline tmp<volScalarField> createField
40 (
41  const fvMesh& mesh,
42  const scalar val
43 )
44 {
46  (
47  IOobject
48  (
50  mesh.time().timeName(),
51  mesh,
55  ),
56  mesh,
58  );
59 }
60 } // End namespace fv
61 } // End namespace Foam
62 
63 
64 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
65 
66 template<class Type>
68 (
69  refPtr<fvMesh>& meshRef,
71 )
72 {
73  const Time& runTime = runTimePtr();
74  const word meshName(polyMesh::defaultRegion);
75 
76  // Fetch mesh from Time database
77  meshRef.cref
78  (
79  runTime.cfindObject<fvMesh>(meshName)
80  );
81 
82  if (!meshRef)
83  {
84  // Fallback: load mesh from disk and cache it
85  meshRef.reset
86  (
87  new fvMesh
88  (
89  IOobject
90  (
91  meshName,
92  runTime.timeName(),
93  runTime,
94  IOobject::MUST_READ,
95  IOobject::NO_WRITE,
96  IOobject::REGISTER
97  )
98  )
99  );
100  }
101 }
102 
103 
104 template<class Type>
106 (
107  const fvMesh& srcMesh,
108  const fvMesh& tgtMesh
109 )
110 {
111  if (consistent_)
112  {
113  interpPtr_.reset
114  (
115  new meshToMesh
116  (
117  srcMesh,
118  tgtMesh,
119  mapMethodName_,
120  patchMapMethodName_
121  )
122  );
123  }
124  else
125  {
126  interpPtr_.reset
127  (
128  new meshToMesh
129  (
130  srcMesh,
131  tgtMesh,
132  mapMethodName_,
133  patchMapMethodName_,
134  patchMap_,
135  cuttingPatches_
136  )
137  );
138  }
139 }
140 
141 
142 template<class Type>
143 template<class VolFieldType>
145 (
146  const fvMesh& thisMesh,
147  const word& fieldName
148 ) const
149 {
150  auto* ptr = thisMesh.getObjectPtr<VolFieldType>(fieldName);
151 
152  if (!ptr)
153  {
154  ptr = new VolFieldType
155  (
156  IOobject
157  (
158  fieldName,
159  thisMesh.time().timeName(),
160  thisMesh,
161  IOobject::MUST_READ,
162  IOobject::NO_WRITE,
163  IOobject::REGISTER
164  ),
165  thisMesh
166  );
167  thisMesh.objectRegistry::store(ptr);
168  }
169 
170  return *ptr;
171 }
172 
173 
174 template<class Type>
176 {
177  const fvMesh& srcMesh = srcMeshPtr_();
178  const fvMesh& tgtMesh = mesh_;
179 
180  // Create mask fields
181  const volScalarField srcFld(createField(srcMesh, 1));
182  volScalarField tgtFld(createField(tgtMesh, 0));
183 
184  // Map the mask field of 1s onto the mask field of 0s
185  interpPtr_->mapSrcToTgt(srcFld, plusEqOp<scalar>(), tgtFld);
186 
187  // Identify and collect cell indices whose values were changed from 0 to 1
188  DynamicList<label> cells;
189  forAll(tgtFld, i)
190  {
191  if (tgtFld[i] != 0)
192  {
193  cells.append(i);
194  }
195  }
196 
197  return cells;
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
202 
203 template<class Type>
205 :
206  positionPtr_(),
207  directionPtr_(),
208  points_(),
209  origin_(),
210  normal_(),
211  active_(false)
212 {}
213 
214 
215 template<class Type>
217 (
218  const word& name,
219  const word& modelType,
220  const dictionary& dict,
221  const fvMesh& mesh
222 )
223 :
224  fv::option(name, modelType, dict, mesh),
225  transform_(),
226  srcTimePtr_(),
227  srcMeshPtr_(),
228  interpPtr_(),
229  patchMap_(),
230  cells_(),
231  cuttingPatches_(),
232  mapMethodName_(),
233  patchMapMethodName_(),
234  consistent_(false)
235 {
236  read(dict);
237 
238  setSourceMesh(srcMeshPtr_, srcTimePtr_);
239 
240  const fvMesh& srcMesh = srcMeshPtr_();
241  const fvMesh& tgtMesh = mesh_;
242  createInterpolation(srcMesh, tgtMesh);
243 
244  cells_ = tgtCellIDs();
245 
246  if (returnReduceAnd(cells_.empty()))
247  {
249  << "No cells selected!" << endl;
250  }
251 
252  transform_.initialize(srcMesh, dict);
253 }
254 
255 
256 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
257 
258 template<class Type>
260 (
261  const fvMesh& srcMesh,
262  const dictionary& dict
263 )
264 {
265  const dictionary* subDictPtr = dict.findDict("transform");
266 
267  if (!subDictPtr)
268  {
269  return false;
270  }
271 
272  positionPtr_.reset
273  (
275  (
276  "position",
277  *subDictPtr,
278  word::null,
279  &srcMesh
280  )
281  );
282 
283  directionPtr_.reset
284  (
286  (
287  "direction",
288  *subDictPtr,
289  word::null,
290  &srcMesh
291  )
292  );
293 
294  if (positionPtr_)
295  {
296  subDictPtr->readIfPresent("origin", origin_);
297  }
298 
299  if (directionPtr_)
300  {
301  subDictPtr->readIfPresent("normal", normal_);
302  normal_.normalise();
303  }
304 
305  points_ = srcMesh.points();
306 
307  active_ = true;
308 
309  return true;
310 }
311 
312 
313 template<class Type>
315 (
316  refPtr<fvMesh>& srcMeshPtr,
317  const scalar t
318 )
319 {
320  if (!positionPtr_)
321  {
322  return;
323  }
324 
325  const pointField translate
326  (
327  points_ + (positionPtr_->value(t) - origin_)
328  );
329 
330  fvMesh& srcMesh = srcMeshPtr.ref();
331  srcMesh.movePoints(translate);
332 }
333 
334 
335 template<class Type>
337 (
338  refPtr<fvMesh>& srcMeshPtr,
339  const scalar t
340 )
341 {
342  if (!directionPtr_)
343  {
344  return;
345  }
346 
347  const vector dir(normalised(directionPtr_->value(t)));
348 
349  const tensor rot(rotationTensor(normal_, dir));
350 
351  pointField rotate(points_);
352 
353  Foam::transform(rotate, rot, rotate);
355  fvMesh& srcMesh = srcMeshPtr.ref();
356  srcMesh.movePoints(rotate);
357 }
358 
359 
360 template<class Type>
362 {
363  if (!fv::option::read(dict))
364  {
365  return false;
366  }
367 
368  fieldNames_.resize(1, coeffs_.getWord("field"));
369 
371 
372  // Load the time database for the source mesh once per simulation
373  if (!srcTimePtr_)
374  {
375  fileName srcMesh(coeffs_.get<fileName>("srcMesh").expand());
376  srcMesh.clean();
377 
378  srcTimePtr_.reset(Time::New(srcMesh));
379 
380  // Set time-step of source database to an arbitrary yet safe value
381  srcTimePtr_().setDeltaT(1.0);
382  }
383 
384  coeffs_.readEntry("mapMethod", mapMethodName_);
385  if (!meshToMesh::interpolationMethodNames_.found(mapMethodName_))
386  {
387  FatalIOErrorInFunction(coeffs_)
388  << type() << " " << name() << ": unknown map method "
389  << mapMethodName_ << nl
390  << "Available methods include: "
392  << exit(FatalIOError);
393  }
394 
395  coeffs_.readIfPresent("consistent", consistent_);
396  coeffs_.readIfPresent("patchMap", patchMap_);
397  coeffs_.readIfPresent("cuttingPatches", cuttingPatches_);
398 
399  if (!coeffs_.readIfPresent("patchMapMethod", patchMapMethodName_))
400  {
402  (
404  );
405  patchMapMethodName_ = meshToMesh::interpolationMethodAMI(mapMethod);
406  }
408  return true;
409 }
410 
411 
412 template<class Type>
414 (
415  fvMatrix<Type>& eqn,
416  const label
417 )
418 {
419  DebugInfo
420  << "MapFieldConstraint<"
422  << ">::constrain for source " << name_ << endl;
423 
424  // Translate and/or rotate source mesh if requested
425  if (transform_.isActive())
426  {
427  // Use time from mesh_ since source mesh does not advance in time
428  const scalar t = mesh_.time().timeOutputValue();
429  transform_.translate(srcMeshPtr_, t);
430  transform_.rotate(srcMeshPtr_, t);
431  }
432 
433  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
434 
435  const word& fldName = fieldNames_[0];
436 
437  const fvMesh& srcMesh = srcMeshPtr_();
438  const fvMesh& tgtMesh = mesh_;
439 
440  // Fetch source and target fields
441  const VolFieldType& srcFld = getOrReadField<VolFieldType>(srcMesh, fldName);
442  VolFieldType& tgtFld = tgtMesh.lookupObjectRef<VolFieldType>(fldName);
443 
444  // When mesh/src changes, reinitilize mesh-to-mesh members
445  if (tgtMesh.changing() || transform_.isActive())
446  {
447  createInterpolation(srcMesh, tgtMesh);
448  cells_ = tgtCellIDs();
449  }
450 
451  // Map source-mesh field onto target-mesh field
452  interpPtr_->mapSrcToTgt(srcFld, plusEqOp<Type>(), tgtFld);
453 
454  // Constrain mapped field in target mesh to avoid overwrite by solver
455  eqn.setValues(cells_, UIndirectList<Type>(tgtFld, cells_));
456 }
457 
458 
459 // ************************************************************************* //
dictionary dict
A class for handling file names.
Definition: fileName.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
static const Enum< interpolationMethod > interpolationMethodNames_
Definition: meshToMesh.H:77
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
engineTime & runTime
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:918
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: TimeNew.C:26
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const dimensionSet dimless
Dimensionless.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: refPtrI.H:216
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
A class for managing references or pointers (no reference counting)
Definition: HashPtrTable.H:49
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
Info<< "Create engine time\"<< endl;autoPtr< engineTime > runTimePtr(engineTime::New(Time::controlDictName, args.rootPath(), args.caseName()))
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
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:405
static const word null
An empty word.
Definition: word.H:84
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvOptionIO.C:48
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
Vector< scalar > vector
Definition: vector.H:57
#define DebugInfo
Report an information message using Foam::Info.
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
static word interpolationMethodAMI(const interpolationMethod method)
Conversion between mesh and patch interpolation methods.
Definition: meshToMesh.C:639
static autoPtr< Function1< Type > > NewIfPresent(const word &entryName, const dictionary &dict, const word &redirectType, const objectRegistry *obrPtr=nullptr)
An optional selector, with fallback redirection.
Definition: Function1New.C:209
bool changing() const noexcept
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:767
The MapFieldConstraint constrains values of given fields of Type with a source field from an external...
#define WarningInFunction
Report a warning using Foam::Warning.
string & expand(const bool allowEmpty=false)
Inplace expand initial tags, tildes, and all occurrences of environment variables as per stringOps::e...
Definition: string.C:166
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
Nothing to be read.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
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.
void setValues(const labelUList &cellLabels, const Type &value)
Set solution in given cells to the specified value and eliminate the corresponding equations from the...
Definition: fvMatrix.C:978
interpolationMethod
Enumeration specifying interpolation method.
Definition: meshToMesh.H:69
A class for managing temporary objects.
Definition: HashPtrTable.H:50
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
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
bool found
virtual void constrain(fvMatrix< Type > &eqn, const label)
Set value on field.
Do not request registration (bool: false)
Namespace for OpenFOAM.
static tmp< volScalarField > createField(const fvMesh &mesh, const scalar val)
void reset(T *p=nullptr) noexcept
Delete managed pointer and set to new given pointer.
Definition: refPtrI.H:314
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...