createNonConformalCouples.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
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 Application
25  createNonConformalCouples
26 
27 Description
28  Utility to create non-conformal couples between non-coupled patches.
29 
30 Usage
31  \b createNonConformalCouples <patch1> <patch2>
32 
33  Options:
34  - \par -overwrite \n
35  Replace the old mesh with the new one, rather than writing the new one
36  into a separate time directory
37 
38  - \par -fields \n
39  Add non-conformal boundary conditions to the fields.
40 
41 Note
42  If run with two arguments, these arguments specify the patches between
43  which a single couple is to be created. The resulting couple will not have
44  a transformation.
45 
46 Usage
47  \b createNonConformalCouples
48 
49 Note
50  If run without arguments then settings are read from a \b
51  system/createNonConformalCouplesDict dictionary (or from a different
52  dictionary specified by the \b -dict option). This dictionary can specify
53  the creation of multiple couples and/or couples with transformations.
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #include "argList.H"
59 #include "fvMeshTools.H"
60 #include "IOobjectList.H"
64 #include "polyMesh.H"
65 #include "processorPolyPatch.H"
66 #include "systemDict.H"
67 #include "Time.H"
68 
69 #include "ReadFields.H"
70 #include "volFields.H"
71 #include "surfaceFields.H"
72 #include "pointFields.H"
73 
74 using namespace Foam;
75 
76 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78 int main(int argc, char *argv[])
79 {
80  #include "addOverwriteOption.H"
81  #include "addRegionOption.H"
82  #include "addDictOption.H"
83 
84  const bool haveArgs = argList::hasArgs(argc, argv);
85  if (haveArgs)
86  {
87  argList::validArgs.append("patch1");
88  argList::validArgs.append("patch2");
90  (
91  "fields",
92  "add non-conformal boundary conditions to the fields"
93  );
94  }
95 
96  #include "setRootCase.H"
98 
99  // Flag to determine whether or not patches are added to fields
100  bool fields;
101 
102  // Patch names between which to create couples, field dictionaries, the
103  // associated cyclic name prefix and transformation (if any)
105  List<Pair<dictionary>> patchFieldDicts;
106  wordList cyclicNames;
108 
109  // If there are patch name arguments, then we assume fields are not being
110  // changed, the cyclic name is just the cyclic typename, and that there is
111  // no transformation. If there are no arguments then get all this
112  // information from the system dictionary.
113  if (haveArgs)
114  {
115  fields = args.optionFound("fields");
116 
118  patchFieldDicts.append(Pair<dictionary>());
119  cyclicNames.append(nonConformalCyclicPolyPatch::typeName);
120  transforms.append(cyclicTransform(true));
121  }
122  else
123  {
124  static const word dictName("createNonConformalCouplesDict");
125 
126  IOdictionary dict(systemDict(dictName, args, runTime));
127 
128  fields = dict.lookupOrDefault<bool>("fields", false);
129 
130  const dictionary& couplesDict =
131  dict.optionalSubDict("nonConformalCouples");
132 
133  forAllConstIter(dictionary, couplesDict, iter)
134  {
135  if (!iter().isDict()) continue;
136 
137  const dictionary& subDict = iter().dict();
138 
139  const bool havePatches = subDict.found("patches");
140  const bool haveOwnerNeighbour =
141  subDict.found("owner") || subDict.found("neighbour");
142 
143  if (havePatches == haveOwnerNeighbour)
144  {
145  FatalIOErrorInFunction(subDict)
146  << "Patches should be specified with either a single "
147  << "\"patches\" entry with a pair of patch names, or with "
148  << "two sub-dictionaries named \"owner\" and "
149  << "\"neighbour\"." << exit(FatalIOError);
150  }
151 
152  if (havePatches)
153  {
154  patchNames.append(subDict.lookup<Pair<word>>("patches"));
155  patchFieldDicts.append(Pair<dictionary>());
156  }
157 
158  if (haveOwnerNeighbour)
159  {
160  const dictionary& ownerDict = subDict.subDict("owner");
161  const dictionary& neighbourDict = subDict.subDict("neighbour");
162 
164  (
165  Pair<word>
166  (
167  ownerDict["patch"],
168  neighbourDict["patch"]
169  )
170  );
171  patchFieldDicts.append
172  (
174  (
175  ownerDict.subOrEmptyDict("patchFields"),
176  neighbourDict.subOrEmptyDict("patchFields")
177  )
178  );
179  }
180 
181  cyclicNames.append(subDict.dictName());
182  transforms.append(cyclicTransform(subDict, true));
183  }
184  }
185 
186  Foam::word meshRegionName = polyMesh::defaultRegion;
187  args.optionReadIfPresent("region", meshRegionName);
188 
189  #include "createNamedMesh.H"
190 
191  const polyBoundaryMesh& patches = mesh.boundaryMesh();
192 
193  const bool overwrite = args.optionFound("overwrite");
194 
195  const word oldInstance = mesh.pointsInstance();
196 
197  // Read the fields
198  IOobjectList objects(mesh, runTime.name());
199  if (fields) Info<< "Reading geometric fields" << nl << endl;
200  #include "readVolFields.H"
201  #include "readSurfaceFields.H"
202  #include "readPointFields.H"
203  if (fields) Info<< endl;
204 
205  // Make sure the mesh is not connected before couples are added
206  mesh.conform();
207 
208  // Start building lists of patches and patch-fields to add
209  List<polyPatch*> newPatches;
210  List<dictionary> newPatchFieldDicts;
211 
212  // Find the first processor patch and face
213  label firstProcPatchi = patches.size(), firstProcFacei = mesh.nFaces();
215  {
216  const polyPatch& pp = patches[patchi];
217 
218  if (isA<processorPolyPatch>(pp) && firstProcPatchi == patches.size())
219  {
220  firstProcPatchi = patchi;
221  firstProcFacei = pp.start();
222  }
223 
224  if (!isA<processorPolyPatch>(pp) && firstProcPatchi != patches.size())
225  {
227  << "Processor patches do not follow boundary patches"
228  << exit(FatalError);
229  }
230  }
231 
232  // Clone the non-processor patches
233  for (label patchi = 0; patchi < firstProcPatchi; ++ patchi)
234  {
235  const polyPatch& pp = patches[patchi];
236 
237  newPatches.append
238  (
239  pp.clone(patches, patchi, pp.size(), pp.start()).ptr()
240  );
241  newPatchFieldDicts.append
242  (
243  dictionary()
244  );
245  }
246 
247  // Convenience function to generate patch names for the owner or neighbour
248  auto nccPatchNames = [&](const label i)
249  {
250  return
251  Pair<word>
252  (
253  cyclicNames[i] + "_on_" + patchNames[i][0],
254  cyclicNames[i] + "_on_" + patchNames[i][1]
255  );
256  };
257 
258  // Add the cyclic patches
259  forAll(patchNames, i)
260  {
261  Info<< indent << "Adding "
262  << nonConformalCyclicPolyPatch::typeName
263  << " interfaces between patches: " << incrIndent << nl
264  << indent << patchNames[i] << decrIndent << nl
265  << indent << "Named:" << incrIndent << nl
266  << indent << Pair<word>(nccPatchNames(i)) << decrIndent << nl
267  << indent << "With transform: " << incrIndent << nl;
268  transforms[i].write(Info);
269  Info<< decrIndent << nl;
270 
271  newPatches.append
272  (
274  (
275  nccPatchNames(i)[0],
276  0,
277  firstProcFacei,
278  newPatches.size(),
279  patches,
280  nonConformalCyclicPolyPatch::typeName,
281  nccPatchNames(i)[1],
282  patchNames[i][0],
283  transforms[i]
284  )
285  );
286  newPatchFieldDicts.append
287  (
288  patchFieldDicts[i][0]
289  );
290  newPatches.append
291  (
293  (
294  nccPatchNames(i)[1],
295  0,
296  firstProcFacei,
297  newPatches.size(),
298  patches,
299  nonConformalCyclicPolyPatch::typeName,
300  nccPatchNames(i)[0],
301  patchNames[i][1],
302  inv(transforms[i])
303  )
304  );
305  newPatchFieldDicts.append
306  (
307  patchFieldDicts[i][1]
308  );
309  }
310 
311  // Add the error patches. Note there is only one for each source patch,
312  // regardless of how many interfaces are attached to that patch.
313  auto appendErrorPatches = [&](const bool owner)
314  {
315  wordHashSet patchANames;
316  forAll(patchNames, i)
317  {
318  patchANames.insert(patchNames[i][!owner]);
319  }
320  forAllConstIter(wordHashSet, patchANames, iter)
321  {
322  newPatches.append
323  (
325  (
326  nonConformalErrorPolyPatch::typeName + "_on_" + iter.key(),
327  0,
328  firstProcFacei,
329  newPatches.size(),
330  patches,
331  nonConformalErrorPolyPatch::typeName,
332  iter.key()
333  )
334  );
335  newPatchFieldDicts.append
336  (
337  dictionary()
338  );
339  }
340  };
341  appendErrorPatches(true);
342  appendErrorPatches(false);
343 
344  // Clone the processor patches
345  for (label patchi = firstProcPatchi; patchi < patches.size(); ++ patchi)
346  {
347  const polyPatch& pp = patches[patchi];
348 
349  newPatches.append
350  (
351  pp.clone(patches, newPatches.size(), pp.size(), pp.start()).ptr()
352  );
353  newPatchFieldDicts.append
354  (
355  dictionary()
356  );
357  }
358 
359  // Add the processor cyclic patches
360  if (Pstream::parRun())
361  {
362  forAll(patchNames, i)
363  {
364  const polyPatch& patch1 = patches[patchNames[i][0]];
365  const polyPatch& patch2 = patches[patchNames[i][1]];
366 
367  boolList procHasPatch1(Pstream::nProcs(), false);
368  procHasPatch1[Pstream::myProcNo()] = !patch1.empty();
369  Pstream::gatherList(procHasPatch1);
370  Pstream::scatterList(procHasPatch1);
371 
372  boolList procHasPatch2(Pstream::nProcs(), false);
373  procHasPatch2[Pstream::myProcNo()] = !patch2.empty();
374  Pstream::gatherList(procHasPatch2);
375  Pstream::scatterList(procHasPatch2);
376 
377  // Multiple cyclic interfaces must be ordered in a specific way for
378  // processor communication to function correctly.
379  //
380  // A communication that is sent from the cyclic owner is received
381  // on the cyclic neighbour and vice versa. Therefore, in a coupled
382  // pair of processors if one sends the owner first the other must
383  // receive the neighbour first.
384  //
385  // We ensure the above by ordering the patches so that for the
386  // lower indexed processor the owner interface comes first, and for
387  // the higher indexed processor the neighbour comes first.
388 
389  auto appendProcPatches = [&](const bool owner, const bool first)
390  {
391  const boolList& procHasPatchA =
392  owner ? procHasPatch1 : procHasPatch2;
393  const boolList& procHasPatchB =
394  owner ? procHasPatch2 : procHasPatch1;
395 
396  if (procHasPatchA[Pstream::myProcNo()])
397  {
398  forAll(procHasPatchB, proci)
399  {
400  if
401  (
402  (
403  (first && proci < Pstream::myProcNo())
404  || (!first && proci > Pstream::myProcNo())
405  )
406  && procHasPatchB[proci]
407  )
408  {
409  newPatches.append
410  (
412  (
413  0,
414  mesh.nFaces(),
415  newPatches.size(),
416  patches,
418  proci,
419  nccPatchNames(i)[!owner],
420  patchNames[i][!owner]
421  )
422  );
423  newPatchFieldDicts.append
424  (
425  patchFieldDicts[i][!owner]
426  );
427  }
428  }
429  }
430  };
431 
432  appendProcPatches(true, true);
433  appendProcPatches(false, true);
434  appendProcPatches(false, false);
435  appendProcPatches(true, false);
436  }
437  }
438 
439  // Re-patch the mesh. Note that new patches are all constraints, so the
440  // dictionary and patch type do not get used. Overrides will be handled
441  // later, once all patches have been added and the mesh has been stitched.
442  forAll(newPatches, newPatchi)
443  {
445  (
446  mesh,
447  *newPatches[newPatchi],
448  dictionary(),
450  false
451  );
452  }
453 
454  // Connect the mesh so that the new stitching topology gets written out
455  fvMeshStitchers::stationary(mesh).connect(false, false, false);
456 
457  // Set the fields on the new patches. All new patches are constraints, so
458  // this should only be creating overrides; e.g., jump cyclics.
459  forAll(newPatches, newPatchi)
460  {
461  if (!newPatchFieldDicts[newPatchi].empty())
462  {
463  fvMeshTools::setPatchFields
464  (
465  mesh,
466  newPatchi,
467  newPatchFieldDicts[newPatchi]
468  );
469  }
470  }
471 
472  mesh.setInstance(runTime.name());
473 
474  // Set the precision of the points data to 10
476 
477  if (!overwrite)
478  {
479  runTime++;
480  }
481  else
482  {
483  mesh.setInstance(oldInstance);
484  }
485 
486  // Write resulting mesh
487  Info<< "Writing mesh to " << runTime.name() << nl << endl;
488  mesh.write();
489 
490  Info<< "End\n" << endl;
491 
492  return 0;
493 }
494 
495 
496 // ************************************************************************* //
Field reading functions for post-processing utilities.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
A HashTable with keys but without contents.
Definition: HashSet.H:62
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
static unsigned int defaultPrecision()
Return the default precision.
Definition: IOstream.H:458
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
virtual Ostream & write(const char)=0
Write character.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static bool hasArgs(int argc, char *argv[])
Return true if there are arguments.
Definition: argList.C:298
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
Cyclic plane transformation.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:121
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:1049
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1076
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
virtual bool connect(const bool changing, const bool geometric, const bool load)
Connect the mesh by adding faces into the nonConformalCyclics.
Mesh stitcher for stationary meshes.
static label addPatch(fvMesh &mesh, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add patch. Inserts patch before all processor patches.
Definition: fvMeshTools.C:34
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
Non-conformal error poly patch. As nonConformalPolyPatch. This patch is a non-coupled non-conformal p...
Non-conformal processor cyclic poly patch. As nonConformalCyclicPolyPatch, but the neighbouring patch...
Foam::polyBoundaryMesh.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:215
A class for handling words, derived from string.
Definition: word.H:62
transformer transforms
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
const fvPatchList & patches
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
IOdictionary systemDict(const word &dictName, const argList &args, const objectRegistry &ob, const word &regionName=polyMesh::defaultRegion)
Definition: systemDict.C:92
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
static const char nl
Definition: Ostream.H:260
wordList patchNames(nPatches)
objects
dictionary dict
Foam::argList args(argc, argv)
Foam::surfaceFields.