fvMeshToolsTemplates.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) 2012-2022 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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvMeshTools.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 template<class GeoField>
33 void Foam::fvMeshTools::setPatchFields
34 (
35  typename GeoField::Mesh& mesh,
36  const label patchi,
37  const dictionary& patchFieldDict
38 )
39 {
40  objectRegistry& obr = const_cast<objectRegistry&>(mesh.thisDb());
41 
42  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
43 
44  forAllIter(typename HashTable<GeoField*>, flds, iter)
45  {
46  GeoField& fld = *iter();
47 
48  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
49 
50  if
51  (
52  patchFieldDict.found(fld.name())
53  || !fvPatch::constraintType(mesh.boundary()[patchi].type())
54  )
55  {
56  bfld.set
57  (
58  patchi,
60  (
61  mesh.boundary()[patchi],
62  fld(),
63  patchFieldDict.subDict(fld.name())
64  )
65  );
66  }
67  }
68 }
69 
70 
71 template<class GeoField>
72 void Foam::fvMeshTools::setPatchFields
73 (
74  typename GeoField::Mesh& mesh,
75  const label patchi,
76  const typename GeoField::value_type& value
77 )
78 {
79  objectRegistry& obr = const_cast<objectRegistry&>(mesh.thisDb());
80 
81  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
82 
83  forAllIter(typename HashTable<GeoField*>, flds, iter)
84  {
85  GeoField& fld = *iter();
86 
87  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
88 
89  bfld[patchi] == value;
90  }
91 }
92 
93 
94 // ************************************************************************* //
Foam::surfaceFields.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: fvPatch.C:61
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
label patchi