cyclicACMIFvPatchField.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) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "fvMatrix.H"
31#include "transformField.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
35template<class Type>
37(
38 const fvPatch& p,
40)
41:
43 coupledFvPatchField<Type>(p, iF),
44 cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p))
45{}
46
47
48template<class Type>
50(
51 const fvPatch& p,
53 const dictionary& dict
54)
55:
57 coupledFvPatchField<Type>(p, iF, dict, dict.found("value")),
58 cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p, dict))
59{
60 if (!isA<cyclicACMIFvPatch>(p))
61 {
63 << " patch type '" << p.type()
64 << "' not constraint type '" << typeName << "'"
65 << "\n for patch " << p.name()
66 << " of field " << this->internalField().name()
67 << " in file " << this->internalField().objectPath()
69 }
70
71 if (!dict.found("value") && this->coupled())
72 {
73 // Extra check: make sure that the non-overlap patch is before
74 // this so it has actually been read - evaluate will crash otherwise
75
78 (
79 this->primitiveField()
80 );
81 if (!fld.boundaryField().set(cyclicACMIPatch_.nonOverlapPatchID()))
82 {
84 << " patch " << p.name()
85 << " of field " << this->internalField().name()
86 << " refers to non-overlap patch "
87 << cyclicACMIPatch_.cyclicACMIPatch().nonOverlapPatchName()
88 << " which is not constructed yet." << nl
89 << " Either supply an initial value or change the ordering"
90 << " in the file"
92 }
93
95 }
96}
97
98
99template<class Type>
101(
103 const fvPatch& p,
105 const fvPatchFieldMapper& mapper
106)
107:
109 coupledFvPatchField<Type>(ptf, p, iF, mapper),
110 cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p))
111{
112 if (!isA<cyclicACMIFvPatch>(this->patch()))
113 {
115 << "' not constraint type '" << typeName << "'"
116 << "\n for patch " << p.name()
117 << " of field " << this->internalField().name()
118 << " in file " << this->internalField().objectPath()
119 << exit(FatalError);
120 }
121}
122
123
124
125template<class Type>
127(
129)
130:
132 coupledFvPatchField<Type>(ptf),
133 cyclicACMIPatch_(ptf.cyclicACMIPatch_)
134{}
135
136
137template<class Type>
139(
142)
143:
145 coupledFvPatchField<Type>(ptf, iF),
146 cyclicACMIPatch_(ptf.cyclicACMIPatch_)
147{}
148
149
150// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151
152template<class Type>
154{
155 return cyclicACMIPatch_.coupled();
156}
157
158
159template<class Type>
162{
163 const Field<Type>& iField = this->primitiveField();
164 //const cyclicACMIPolyPatch& cpp = cyclicACMIPatch_.cyclicACMIPatch();
165
166 // By pass polyPatch to get nbrId. Instead use cyclicAMIFvPatch virtual
167 // neighbPatch()
168 const cyclicACMIFvPatch& neighbPatch = cyclicACMIPatch_.neighbPatch();
169 const labelUList& nbrFaceCells = neighbPatch.faceCells();
170
171 tmp<Field<Type>> tpnf
172 (
173 cyclicACMIPatch_.interpolate
174 (
176 (
177 iField,
178 nbrFaceCells
179 //cpp.neighbPatch().faceCells()
180 )
181 )
182 );
183
184 if (doTransform())
185 {
186 tpnf.ref() = transform(forwardT(), tpnf());
187 }
188
189 return tpnf;
190}
191
192
193template<class Type>
196{
199 (
200 this->primitiveField()
201 );
202
203 return refCast<const cyclicACMIFvPatchField<Type>>
204 (
205 fld.boundaryField()[cyclicACMIPatch_.neighbPatchID()]
206 );
207}
208
209
210template<class Type>
213{
216 (
217 this->primitiveField()
218 );
219
220 // WIP: Needs to re-direct nonOverlapPatchID to new patchId for assembly?
221 return fld.boundaryField()[cyclicACMIPatch_.nonOverlapPatchID()];
222}
223
224
225template<class Type>
227(
228 solveScalarField& result,
229 const bool add,
230 const lduAddressing& lduAddr,
231 const label patchId,
232 const solveScalarField& psiInternal,
233 const scalarField& coeffs,
234 const direction cmpt,
236) const
237{
238 // note: only applying coupled contribution
239
240// const labelUList& nbrFaceCellsCoupled =
241// lduAddr.patchAddr
242// (
243// cyclicACMIPatch_.cyclicACMIPatch().neighbPatchID()
244// );
245
246 const labelUList& nbrFaceCellsCoupled =
247 lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
248
249 solveScalarField pnf(psiInternal, nbrFaceCellsCoupled);
250
251 // Transform according to the transformation tensors
252 transformCoupleField(pnf, cmpt);
253
254 pnf = cyclicACMIPatch_.interpolate(pnf);
255
256 const labelUList& faceCells = lduAddr.patchAddr(patchId);
257
258 this->addToInternalField(result, !add, faceCells, coeffs, pnf);
259}
260
261
262template<class Type>
264(
265 Field<Type>& result,
266 const bool add,
267 const lduAddressing& lduAddr,
268 const label patchId,
269 const Field<Type>& psiInternal,
270 const scalarField& coeffs,
272) const
273{
274 // note: only applying coupled contribution
275
276 const labelUList& nbrFaceCellsCoupled =
277 lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
278
279 Field<Type> pnf(psiInternal, nbrFaceCellsCoupled);
280
281 // Transform according to the transformation tensors
282 transformCoupleField(pnf);
283
284 pnf = cyclicACMIPatch_.interpolate(pnf);
285
286 const labelUList& faceCells = lduAddr.patchAddr(patchId);
287
288 this->addToInternalField(result, !add, faceCells, coeffs, pnf);
289}
290
291
292template<class Type>
294(
295 fvMatrix<Type>& matrix
296)
297{
298 const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
299
300 // Nothing to be done by the AMI, but re-direct to non-overlap patch
301 // with non-overlap patch weights
302 const fvPatchField<Type>& npf = nonOverlapPatchField();
303
304 const_cast<fvPatchField<Type>&>(npf).manipulateMatrix(matrix, 1.0 - mask);
305}
306
307
308template<class Type>
310(
311 fvMatrix<Type>& matrix,
312 const label mat,
313 const direction cmpt
314)
315{
316 if (this->cyclicACMIPatch().owner())
317 {
318 label index = this->patch().index();
319
320 const label globalPatchID =
321 matrix.lduMeshAssembly().patchLocalToGlobalMap()[mat][index];
322
323 const Field<scalar> intCoeffsCmpt
324 (
325 matrix.internalCoeffs()[globalPatchID].component(cmpt)
326 );
327
328 const Field<scalar> boundCoeffsCmpt
329 (
330 matrix.boundaryCoeffs()[globalPatchID].component(cmpt)
331 );
332
333 tmp<Field<scalar>> tintCoeffs(coeffs(matrix, intCoeffsCmpt, mat));
334 tmp<Field<scalar>> tbndCoeffs(coeffs(matrix, boundCoeffsCmpt, mat));
335 const Field<scalar>& intCoeffs = tintCoeffs.ref();
336 const Field<scalar>& bndCoeffs = tbndCoeffs.ref();
337
338 const labelUList& u = matrix.lduAddr().upperAddr();
339 const labelUList& l = matrix.lduAddr().lowerAddr();
340
341 label subFaceI = 0;
342
343 const labelList& faceMap =
344 matrix.lduMeshAssembly().faceBoundMap()[mat][index];
345
346 forAll (faceMap, j)
347 {
348 label globalFaceI = faceMap[j];
349
350 const scalar boundCorr = -bndCoeffs[subFaceI];
351 const scalar intCorr = -intCoeffs[subFaceI];
352
353 matrix.upper()[globalFaceI] += boundCorr;
354 matrix.diag()[u[globalFaceI]] -= intCorr;
355 matrix.diag()[l[globalFaceI]] -= boundCorr;
356
357 if (matrix.asymmetric())
358 {
359 matrix.lower()[globalFaceI] += intCorr;
360 }
361 subFaceI++;
362 }
363
364 // Set internalCoeffs and boundaryCoeffs in the assembly matrix
365 // on clyclicAMI patches to be used in the individual matrix by
366 // matrix.flux()
367 if (matrix.psi(mat).mesh().fluxRequired(this->internalField().name()))
368 {
369 matrix.internalCoeffs().set
370 (
371 globalPatchID, intCoeffs*pTraits<Type>::one
372 );
373 matrix.boundaryCoeffs().set
374 (
375 globalPatchID, bndCoeffs*pTraits<Type>::one
376 );
377
378 const label nbrPathID =
379 cyclicACMIPatch_.cyclicACMIPatch().neighbPatchID();
380
381 const label nbrGlobalPatchID =
383 [mat][nbrPathID];
384
385 matrix.internalCoeffs().set
386 (
387 nbrGlobalPatchID, intCoeffs*pTraits<Type>::one
388 );
389 matrix.boundaryCoeffs().set
390 (
391 nbrGlobalPatchID, bndCoeffs*pTraits<Type>::one
392 );
393 }
394 }
395}
396
397
398template<class Type>
401(
402 fvMatrix<Type>& matrix,
403 const Field<scalar>& coeffs,
404 const label mat
405) const
406{
407 const label index(this->patch().index());
408
409 const label nSubFaces
410 (
411 matrix.lduMeshAssembly().cellBoundMap()[mat][index].size()
412 );
413
414 Field<scalar> mapCoeffs(nSubFaces, Zero);
415
416 const scalarListList& srcWeight =
417 cyclicACMIPatch_.cyclicACMIPatch().AMI().srcWeights();
418
419 const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
420
421 const scalar tol = cyclicACMIPolyPatch::tolerance();
422 label subFaceI = 0;
423 forAll(mask, faceI)
424 {
425 const scalarList& w = srcWeight[faceI];
426 for(label i=0; i<w.size(); i++)
427 {
428 if (mask[faceI] > tol)
429 {
430 const label localFaceId =
432 [mat][index][subFaceI];
433 mapCoeffs[subFaceI] = w[i]*coeffs[localFaceId];
434 }
435 subFaceI++;
436 }
437 }
438
439 return tmp<Field<scalar>>(new Field<scalar>(mapCoeffs));
440}
441
442
443template<class Type>
445{
446 // Update non-overlap patch - some will implement updateCoeffs, and
447 // others will implement evaluate
448
449 // Pass in (1 - mask) to give non-overlap patch the chance to do
450 // manipulation of non-face based data
451
452 const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
453 const fvPatchField<Type>& npf = nonOverlapPatchField();
454 const_cast<fvPatchField<Type>&>(npf).updateWeightedCoeffs(1.0 - mask);
455}
456
457
458template<class Type>
460{
462 this->writeEntry("value", os);
463}
464
465
466// ************************************************************************* //
bool found
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static const char *const typeName
Typename for Field.
Definition: FieldBase.H:59
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
const T * set(const label i) const
Definition: PtrList.H:138
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
commsTypes
Types of communications.
Definition: UPstream.H:67
const dictionary & coeffs() const
Return const dictionary of the model.
Abstract base class for coupled patches.
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
virtual const labelUList & faceCells() const
Return faceCell addressing.
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
const cyclicACMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
virtual void updateInterfaceMatrix(solveScalarField &result, const bool add, const lduAddressing &lduAddr, const label patchId, const solveScalarField &psiInternal, const scalarField &coeffs, const direction cmpt, const Pstream::commsTypes commsType) const
Update result field based on interface functionality.
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
const cyclicACMIPolyPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic patch.
virtual const cyclicACMIFvPatch & neighbPatch() const
Return neighbour fvPatch.
virtual label nonOverlapPatchID() const
Return neighbour.
Abstract base class for cyclic ACMI coupled interfaces.
static scalar tolerance()
Overlap tolerance.
const word & nonOverlapPatchName() const
Non-overlapping patch name.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
virtual bool write()
Write the output fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
const FieldField< Field, Type > & internalCoeffs() const noexcept
Definition: fvMatrix.H:470
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
Definition: fvMatrix.H:484
const lduPrimitiveMeshAssembly & lduMeshAssembly()
Return optional lduAdressing.
Definition: fvMatrix.H:405
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:412
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
const Field< Type > & primitiveField() const
Return internal field reference.
Definition: fvPatchField.H:374
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:368
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
scalarField & upper()
Definition: lduMatrix.C:203
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:578
scalarField & lower()
Definition: lduMatrix.C:174
scalarField & diag()
Definition: lduMatrix.C:192
bool asymmetric() const
Definition: lduMatrix.H:633
const labelListListList & cellBoundMap() const
Return patch local sub-face to nbrCellId map.
const labelListListList & faceBoundMap() const
Return boundary face map.
const labelListListList & facePatchFaceMap() const
Return patch local sub-face to local patch face map.
const labelListList & patchLocalToGlobalMap() const
Return patchLocalToGlobalMap.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
label patchId(-1)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
uint8_t direction
Definition: direction.H:56
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Spatial transformation functions for primitive fields.