Loading...
Searching...
No Matches
cellCellStencilTemplates.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) 2018-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
29
30// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
31
32template<class T>
34(
36 const fvMesh& mesh,
37 const cellCellStencil& overlap,
38 const List<scalarList>& wghts
39)
40{
41 const labelListList& stencil = overlap.cellStencil();
42
43 if (stencil.size() != mesh.nCells())
44 {
45 return;
46 }
47
48 const mapDistribute& map = overlap.cellInterpolationMap();
49 const labelUList& cellIDs = overlap.interpolationCells();
50 const scalarList& factor = overlap.cellInterpolationWeight();
51
52 Field<T> work(psi);
53 map.mapDistributeBase::distribute(work, UPstream::msgType()+1);
54
55 forAll(cellIDs, i)
56 {
57 const label celli = cellIDs[i];
58
59 const scalarList& w = wghts[celli];
60 const labelList& nbrs = stencil[celli];
61 const scalar f = factor[celli];
62
63 if (nbrs.size() == 0 && f != 0.0)
64 {
65 FatalErrorInFunction << "problem: cell:" << celli
66 << " at:" << mesh.cellCentres()[celli]
67 << " type:" << overlap.cellTypes()[celli]
68 << " stencil:" << nbrs
69 << " factor:" << f << exit(FatalError);
70 }
71
73 forAll(nbrs, nbrI)
74 {
75 s += w[nbrI]*work[nbrs[nbrI]];
76 }
78 psi[celli] = (1.0-f)*psi[celli] + f*s;
79 }
80}
81
82
83template<class T>
84void Foam::cellCellStencil::interpolate(const fvMesh& mesh, Field<T>& psi) const
85{
86 const cellCellStencil& overlap = *this;
87
89 (
90 psi,
91 mesh,
93 overlap.cellInterpolationWeights()
94 );
95}
96
97
98template<class GeoField>
99void Foam::cellCellStencil::interpolate(GeoField& psi) const
100{
102 (
103 psi.mesh(),
104 psi.primitiveFieldRef()
105 );
106 psi.correctBoundaryConditions();
107}
108
109
110template<class GeoField>
112(
113 const fvMesh& mesh,
114 const wordHashSet& suppressed
115) const
116{
117 const cellCellStencil& overlap = *this;
118
119 for (const GeoField& field : mesh.thisDb().csorted<GeoField>())
120 {
121 const word& name = field.name();
122
123 if (!suppressed.found(baseName(name)))
124 {
125 if (debug)
126 {
127 Pout<< "cellCellStencil::interpolate: interpolating : "
128 << name << endl;
129 }
130
131 auto& fld = const_cast<GeoField&>(field);
132
134 (
135 fld.primitiveFieldRef(),
136 mesh,
137 overlap,
138 overlap.cellInterpolationWeights()
139 );
140 }
141 else
142 {
143 if (debug)
144 {
145 Pout<< "cellCellStencil::interpolate: skipping : "
146 << name << endl;
148 }
149 }
150}
151
152template<class Type>
153Foam::tmp<Foam::volScalarField>
155(
156 const fvMesh& mesh,
157 const word& name,
158 const UList<Type>& psi
159)
160{
161 auto tfld = volScalarField::New
162 (
163 name,
165 mesh,
168 );
169 auto& fld = tfld.ref();
170
171 forAll(psi, cellI)
172 {
173 fld[cellI] = psi[cellI];
174 }
175 return tfld;
176}
177
178
179template<class GeoField, class SuppressBC>
181(
182 GeoField& psi
183)
184{
185 // correctBoundaryConditions excluding oversetFvPatchFields
186
187 psi.boundaryFieldRef().evaluate_if
188 (
189 [](const auto& pfld)
190 {
191 return (!isA<SuppressBC>(pfld));
192 },
194 );
195}
196
197
198// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
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())
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
@ NO_REGISTER
Do not request registration (bool: false).
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
Calculation of interpolation stencils.
static word baseName(const word &name)
Helper: strip off trailing _0.
static void correctBoundaryConditions(GeoField &psi)
Version of correctBoundaryConditions that excludes 'overset' bcs.
static void interpolate(Field< T > &psi, const fvMesh &mesh, const cellCellStencil &overlap, const List< scalarList > &wghts)
Interpolation of acceptor cells from donor cells.
static tmp< volScalarField > createField(const fvMesh &mesh, const word &name, const UList< Type > &)
Helper: create volScalarField for postprocessing.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Class containing processor-to-processor mapping information.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & psi
rDeltaTY field()
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const cellCellStencilObject & overlap
Definition correctPhi.H:57
Namespace for handling debugging switches.
Definition debug.C:45
const dimensionSet dimless
Dimensionless.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299