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-2024 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 
38 static inline tmp<volScalarField> createField
39 (
40  const fvMesh& mesh,
41  const scalar val
42 )
43 {
44  return volScalarField::New
45  (
48  mesh,
49  val,
50  dimless
51  );
52 }
53 
54 } // End namespace Foam
55 
56 
57 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58 
59 template<class Type>
61 (
62  refPtr<fvMesh>& meshRef,
63  const autoPtr<Time>& runTimePtr
64 )
65 {
66  const Time& runTime = runTimePtr();
67  const word meshName(polyMesh::defaultRegion);
68 
69  // Fetch mesh from Time database
70  meshRef.cref
71  (
72  runTime.cfindObject<fvMesh>(meshName)
73  );
74 
75  if (!meshRef)
76  {
77  // Fallback: load mesh from disk and cache it
78  meshRef.reset
79  (
80  new fvMesh
81  (
82  IOobject
83  (
84  meshName,
85  runTime.timeName(),
86  runTime,
87  IOobject::MUST_READ,
88  IOobject::NO_WRITE,
89  IOobject::REGISTER
90  )
91  )
92  );
93  }
94 }
95 
96 
97 template<class Type>
99 (
100  const fvMesh& srcMesh,
101  const fvMesh& tgtMesh
102 )
103 {
104  if (consistent_)
105  {
106  interpPtr_.reset
107  (
108  new meshToMesh
109  (
110  srcMesh,
111  tgtMesh,
112  mapMethodName_,
113  patchMapMethodName_
114  )
115  );
116  }
117  else
118  {
119  interpPtr_.reset
120  (
121  new meshToMesh
122  (
123  srcMesh,
124  tgtMesh,
125  mapMethodName_,
126  patchMapMethodName_,
127  patchMap_,
128  cuttingPatches_
129  )
130  );
131  }
132 }
133 
134 
135 template<class Type>
136 template<class VolFieldType>
138 (
139  const fvMesh& thisMesh,
140  const word& fieldName
141 ) const
142 {
143  auto* ptr = thisMesh.getObjectPtr<VolFieldType>(fieldName);
144 
145  if (!ptr)
146  {
147  ptr = new VolFieldType
148  (
149  IOobject
150  (
151  fieldName,
152  thisMesh.time().timeName(),
153  thisMesh.thisDb(),
154  IOobject::MUST_READ,
155  IOobject::NO_WRITE,
156  IOobject::REGISTER
157  ),
158  thisMesh
159  );
160  regIOobject::store(ptr);
161  }
162 
163  return *ptr;
164 }
165 
166 
167 template<class Type>
169 {
170  const fvMesh& srcMesh = srcMeshPtr_();
171  const fvMesh& tgtMesh = mesh_;
172 
173  // Create mask fields
174  const volScalarField srcFld(createField(srcMesh, 1));
175  volScalarField tgtFld(createField(tgtMesh, 0));
176 
177  // Map the mask field of 1s onto the mask field of 0s
178  interpPtr_->mapSrcToTgt(srcFld, plusEqOp<scalar>(), tgtFld);
179 
180  // Identify and collect cell indices whose values were changed from 0 to 1
181  DynamicList<label> cells;
182  forAll(tgtFld, i)
183  {
184  if (tgtFld[i] != 0)
185  {
186  cells.append(i);
187  }
188  }
189 
190  return cells;
191 }
192 
193 
194 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
195 
196 template<class Type>
198 :
199  positionPtr_(),
200  directionPtr_(),
201  points_(),
202  origin_(),
203  normal_(),
204  active_(false)
205 {}
206 
207 
208 template<class Type>
210 (
211  const word& name,
212  const word& modelType,
213  const dictionary& dict,
214  const fvMesh& mesh
215 )
216 :
217  fv::option(name, modelType, dict, mesh),
218  transform_(),
219  srcTimePtr_(),
220  srcMeshPtr_(),
221  interpPtr_(),
222  patchMap_(),
223  cells_(),
224  cuttingPatches_(),
225  mapMethodName_(),
226  patchMapMethodName_(),
227  consistent_(false)
228 {
229  read(dict);
230 
231  setSourceMesh(srcMeshPtr_, srcTimePtr_);
232 
233  const fvMesh& srcMesh = srcMeshPtr_();
234  const fvMesh& tgtMesh = mesh_;
235  createInterpolation(srcMesh, tgtMesh);
236 
237  cells_ = tgtCellIDs();
238 
239  if (returnReduceAnd(cells_.empty()))
240  {
242  << "No cells selected!" << endl;
243  }
244 
245  transform_.initialize(srcMesh, dict);
246 }
247 
248 
249 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
250 
251 template<class Type>
253 (
254  const fvMesh& srcMesh,
255  const dictionary& dict
256 )
257 {
258  const dictionary* subDictPtr = dict.findDict("transform");
259 
260  if (!subDictPtr)
261  {
262  return false;
263  }
264 
265  positionPtr_.reset
266  (
268  (
269  "position",
270  *subDictPtr,
271  word::null,
272  &srcMesh
273  )
274  );
275 
276  directionPtr_.reset
277  (
279  (
280  "direction",
281  *subDictPtr,
282  word::null,
283  &srcMesh
284  )
285  );
286 
287  if (positionPtr_)
288  {
289  subDictPtr->readIfPresent("origin", origin_);
290  }
291 
292  if (directionPtr_)
293  {
294  subDictPtr->readIfPresent("normal", normal_);
295  normal_.normalise();
296  }
297 
298  points_ = srcMesh.points();
299 
300  active_ = true;
301 
302  return true;
303 }
304 
305 
306 template<class Type>
308 (
309  refPtr<fvMesh>& srcMeshPtr,
310  const scalar t
311 )
312 {
313  if (!positionPtr_)
314  {
315  return;
316  }
317 
318  const pointField translate
319  (
320  points_ + (positionPtr_->value(t) - origin_)
321  );
322 
323  fvMesh& srcMesh = srcMeshPtr.ref();
324  srcMesh.movePoints(translate);
325 }
326 
327 
328 template<class Type>
330 (
331  refPtr<fvMesh>& srcMeshPtr,
332  const scalar t
333 )
334 {
335  if (!directionPtr_)
336  {
337  return;
338  }
339 
340  const vector dir(normalised(directionPtr_->value(t)));
341 
342  const tensor rot(rotationTensor(normal_, dir));
343 
344  pointField rotate(points_);
345 
346  Foam::transform(rotate, rot, rotate);
348  fvMesh& srcMesh = srcMeshPtr.ref();
349  srcMesh.movePoints(rotate);
350 }
351 
352 
353 template<class Type>
355 {
356  if (!fv::option::read(dict))
357  {
358  return false;
359  }
360 
361  fieldNames_.resize(1, coeffs_.getWord("field"));
362 
364 
365  // Load the time database for the source mesh once per simulation
366  if (!srcTimePtr_)
367  {
368  fileName srcMesh(coeffs_.get<fileName>("srcMesh").expand());
369  srcMesh.clean();
370 
371  srcTimePtr_.reset(Time::New(srcMesh));
372 
373  // Set time-step of source database to an arbitrary yet safe value
374  srcTimePtr_().setDeltaT(1.0);
375  }
376 
377  coeffs_.readEntry("mapMethod", mapMethodName_);
378  if (!meshToMesh::interpolationMethodNames_.found(mapMethodName_))
379  {
380  FatalIOErrorInFunction(coeffs_)
381  << type() << " " << name() << ": unknown map method "
382  << mapMethodName_ << nl
383  << "Available methods include: "
385  << exit(FatalIOError);
386  }
387 
388  coeffs_.readIfPresent("consistent", consistent_);
389  coeffs_.readIfPresent("patchMap", patchMap_);
390  coeffs_.readIfPresent("cuttingPatches", cuttingPatches_);
391 
392  if (!coeffs_.readIfPresent("patchMapMethod", patchMapMethodName_))
393  {
395  (
397  );
398  patchMapMethodName_ = meshToMesh::interpolationMethodAMI(mapMethod);
399  }
401  return true;
402 }
403 
404 
405 template<class Type>
407 (
408  fvMatrix<Type>& eqn,
409  const label
410 )
411 {
412  DebugInfo
413  << "MapFieldConstraint<"
415  << ">::constrain for source " << name_ << endl;
416 
417  // Translate and/or rotate source mesh if requested
418  if (transform_.isActive())
419  {
420  // Use time from mesh_ since source mesh does not advance in time
421  const scalar t = mesh_.time().timeOutputValue();
422  transform_.translate(srcMeshPtr_, t);
423  transform_.rotate(srcMeshPtr_, t);
424  }
425 
426  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
427 
428  const word& fldName = fieldNames_[0];
429 
430  const fvMesh& srcMesh = srcMeshPtr_();
431  const fvMesh& tgtMesh = mesh_;
432 
433  // Fetch source and target fields
434  const VolFieldType& srcFld = getOrReadField<VolFieldType>(srcMesh, fldName);
435  VolFieldType& tgtFld = tgtMesh.lookupObjectRef<VolFieldType>(fldName);
436 
437  // When mesh/src changes, reinitilize mesh-to-mesh members
438  if (tgtMesh.changing() || transform_.isActive())
439  {
440  createInterpolation(srcMesh, tgtMesh);
441  cells_ = tgtCellIDs();
442  }
443 
444  // Map source-mesh field onto target-mesh field
445  interpPtr_->mapSrcToTgt(srcFld, plusEqOp<Type>(), tgtFld);
446 
447  // Constrain mapped field in target mesh to avoid overwrite by solver
448  eqn.setValues(cells_, UIndirectList<Type>(tgtFld, cells_));
449 }
450 
451 
452 // ************************************************************************* //
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:675
engineTime & runTime
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static tmp< volScalarField > createField(const fvMesh &mesh, const scalar val)
virtual void movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:929
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.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
const dimensionSet dimless
Dimensionless.
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from vector n1 to n2.
Definition: transform.H:47
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
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:406
static const word null
An empty word.
Definition: word.H:84
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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:768
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:637
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
A special matrix type and solver, designed for finite volume solutions of scalar equations.
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:972
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
bool found
virtual void constrain(fvMatrix< Type > &eqn, const label)
Set value on field.
Do not request registration (bool: false)
Namespace for OpenFOAM.
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 ...