Loading...
Searching...
No Matches
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-2025 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:
44 cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p)),
45 patchNeighbourFieldPtr_(nullptr)
46{}
47
48
49template<class Type>
51(
52 const fvPatch& p,
54 const dictionary& dict
55)
56:
58 coupledFvPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
59 cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p, dict)),
60 patchNeighbourFieldPtr_(nullptr)
61{
63 {
65 << " patch type '" << p.type()
66 << "' not constraint type '" << typeName << "'"
67 << "\n for patch " << p.name()
68 << " of field " << this->internalField().name()
69 << " in file " << this->internalField().objectPath()
71 }
72
73 if (cacheNeighbourField())
74 {
75 // Handle neighbour value first, before any evaluate()
76 const auto* hasNeighbValue =
77 dict.findEntry("neighbourValue", keyType::LITERAL);
78
79 if (hasNeighbValue)
80 {
81 patchNeighbourFieldPtr_.reset
82 (
83 new Field<Type>(*hasNeighbValue, p.size())
84 );
85 }
86 }
87
88 // Use 'value' supplied, or set to coupled field
89 if (!this->readValueEntry(dict) && this->coupled())
90 {
91 // Extra check: make sure that the non-overlap patch is before
92 // this so it has actually been read - evaluate will crash otherwise
93
94 const auto& fld =
95 static_cast<const GeometricField<Type, fvPatchField, volMesh>&>
96 (
97 this->primitiveField()
98 );
99
100 if (!fld.boundaryField().set(cyclicACMIPatch_.nonOverlapPatchID()))
101 {
103 << " patch " << p.name()
104 << " of field " << this->internalField().name()
105 << " refers to non-overlap patch "
106 << cyclicACMIPatch_.cyclicACMIPatch().nonOverlapPatchName()
107 << " which is not constructed yet." << nl
108 << " Either supply an initial value or change the ordering"
109 << " in the file"
110 << exit(FatalIOError);
111 }
112
113 // Tricky: avoid call to evaluate without call to initEvaluate.
114 // For now just disable the localConsistency to make it use the
115 // old logic (ultimately calls the fully self contained
116 // patchNeighbourField)
117
118 const auto oldConsistency = FieldBase::localBoundaryConsistency(0);
119
123 }
124}
125
126
127template<class Type>
129(
130 const cyclicACMIFvPatchField<Type>& ptf,
131 const fvPatch& p,
132 const DimensionedField<Type, volMesh>& iF,
133 const fvPatchFieldMapper& mapper
134)
135:
136 cyclicACMILduInterfaceField(),
137 coupledFvPatchField<Type>(ptf, p, iF, mapper),
138 cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p)),
139 patchNeighbourFieldPtr_(nullptr)
140{
141 //if (ptf.patchNeighbourFieldPtr_ && cacheNeighbourField())
142 //{
143 // patchNeighbourFieldPtr_.reset
144 // (
145 // new Field<Type>(ptf.patchNeighbourFieldPtr_(), mapper)
146 // );
147 //}
148
149 if (!isA<cyclicACMIFvPatch>(this->patch()))
150 {
152 << "\n patch type '" << p.type()
153 << "' not constraint type '" << typeName << "'"
154 << "\n for patch " << p.name()
155 << " of field " << this->internalField().name()
156 << " in file " << this->internalField().objectPath()
157 << exit(FatalError);
158 }
159 if (debug && !ptf.all_ready())
160 {
161 FatalErrorInFunction
162 << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
163 << abort(FatalError);
164 }
165}
166
167
168template<class Type>
170(
172)
173:
175 coupledFvPatchField<Type>(ptf),
176 cyclicACMIPatch_(ptf.cyclicACMIPatch_),
177 patchNeighbourFieldPtr_(nullptr)
178{
179 if (debug && !ptf.all_ready())
180 {
181 FatalErrorInFunction
182 << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
183 << abort(FatalError);
184 }
185}
186
187
188template<class Type>
190(
193)
194:
196 coupledFvPatchField<Type>(ptf, iF),
197 cyclicACMIPatch_(ptf.cyclicACMIPatch_),
198 patchNeighbourFieldPtr_(nullptr)
199{
200 if (debug && !ptf.all_ready())
201 {
202 FatalErrorInFunction
203 << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
204 << abort(FatalError);
206}
207
208
209// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210
211template<class Type>
213{
214 return cyclicACMIPatch_.coupled();
215}
216
217
218template<class Type>
219bool Foam::cyclicACMIFvPatchField<Type>::all_ready() const
220{
221 int done = 0;
222
223 if
224 (
226 (
227 recvRequests_.start(),
228 recvRequests_.size()
229 )
231 (
232 recvRequests1_.start(),
233 recvRequests1_.size()
234 )
235 )
236 {
237 recvRequests_.clear();
238 recvRequests1_.clear();
239 ++done;
240 }
241
242 if
243 (
245 (
246 sendRequests_.start(),
247 sendRequests_.size()
248 )
250 (
251 sendRequests1_.start(),
252 sendRequests1_.size()
253 )
254 )
255 {
256 sendRequests_.clear();
257 sendRequests1_.clear();
258 ++done;
260
261 return (done == 2);
262}
263
264
265template<class Type>
267{
268 if
269 (
271 (
272 recvRequests_.start(),
273 recvRequests_.size()
274 )
276 (
277 recvRequests1_.start(),
278 recvRequests1_.size()
279 )
280 )
281 {
282 recvRequests_.clear();
283 recvRequests1_.clear();
284
285 if
286 (
288 (
289 sendRequests_.start(),
290 sendRequests_.size()
291 )
293 (
294 sendRequests1_.start(),
295 sendRequests1_.size()
296 )
297 )
298 {
299 sendRequests_.clear();
300 sendRequests1_.clear();
301 }
302
303 return true;
304 }
305
306 return false;
307}
308
309
310template<class Type>
311Foam::tmp<Foam::Field<Type>>
312Foam::cyclicACMIFvPatchField<Type>::getNeighbourField
313(
314 const UList<Type>& internalData
315) const
316{
318 << "cyclicACMIFvPatchField::getNeighbourField(const UList<Type>&) :"
319 << " field:" << this->internalField().name()
320 << " patch:" << this->patch().name()
321 << endl;
322
323 // By pass polyPatch to get nbrId. Instead use cyclicACMIFvPatch virtual
324 // neighbPatch()
325 const auto& neighbPatch = cyclicACMIPatch_.neighbPatch();
326 const labelUList& nbrFaceCells = neighbPatch.faceCells();
327
328 tmp<Field<Type>> tpnf
329 (
330 cyclicACMIPatch_.interpolate
331 (
332 Field<Type>(internalData, nbrFaceCells)
333 )
334 );
335
336 if (doTransform())
337 {
338 transform(tpnf.ref(), forwardT(), tpnf());
339 }
340
341 return tpnf;
342}
343
344
345template<class Type>
346bool Foam::cyclicACMIFvPatchField<Type>::cacheNeighbourField()
347{
348 /*
349 return (FieldBase::localBoundaryConsistency() != 0);
350 */
351 return false;
352}
353
354
355template<class Type>
356Foam::tmp<Foam::Field<Type>>
357Foam::cyclicACMIFvPatchField<Type>::getPatchNeighbourField
358(
359 const bool checkCommunicator
360) const
361{
362 const auto& AMI = this->ownerAMI();
363
364 if
365 (
366 AMI.distributed() && cacheNeighbourField()
367 && (!checkCommunicator || AMI.comm() != -1)
368 )
369 {
370 if (!this->ready())
371 {
373 << "Outstanding recv request(s) on patch "
374 << cyclicACMIPatch_.name()
375 << " field " << this->internalField().name()
376 << abort(FatalError);
377 }
378
379 // Initialise if not done in construct-from-dictionary
380 if (!patchNeighbourFieldPtr_)
381 {
383 << "cyclicACMIFvPatchField::patchNeighbourField() :"
384 << " field:" << this->internalField().name()
385 << " patch:" << this->patch().name()
386 << " calculating&caching patchNeighbourField"
387 << endl;
388
389 // Do interpolation and store result
390 patchNeighbourFieldPtr_.reset
391 (
392 getNeighbourField(this->primitiveField()).ptr()
393 );
394 }
395 else
396 {
398 << "cyclicACMIFvPatchField::patchNeighbourField() :"
399 << " field:" << this->internalField().name()
400 << " patch:" << this->patch().name()
401 << " returning cached patchNeighbourField"
402 << endl;
403 }
404 return patchNeighbourFieldPtr_();
405 }
406 else
407 {
408 // Do interpolation
410 << "cyclicACMIFvPatchField::evaluate() :"
411 << " field:" << this->internalField().name()
412 << " patch:" << this->patch().name()
413 << " calculating up-to-date patchNeighbourField"
414 << endl;
415
416 return getNeighbourField(this->primitiveField());
417 }
418}
419
420
421template<class Type>
424{
425 return this->getPatchNeighbourField(true); // checkCommunicator = true
426}
427
428
429template<class Type>
431(
432 UList<Type>& pnf
433) const
434{
435 // checkCommunicator = false
436 auto tpnf = this->getPatchNeighbourField(false);
437 pnf.deepCopy(tpnf());
438}
439
440
441template<class Type>
444{
445 const auto& fld =
447 (
448 this->primitiveField()
449 );
450
452 (
453 fld.boundaryField()[cyclicACMIPatch_.neighbPatchID()]
454 );
455}
456
457
458template<class Type>
461{
462 const auto& fld =
464 (
465 this->primitiveField()
466 );
468 // WIP: Needs to re-direct nonOverlapPatchID to new patchId for assembly?
469 return fld.boundaryField()[cyclicACMIPatch_.nonOverlapPatchID()];
470}
471
472
473template<class Type>
475(
476 const Pstream::commsTypes commsType
477)
478{
479 if (!this->updated())
480 {
481 this->updateCoeffs();
482 }
483
484 const auto& AMI = this->ownerAMI();
485
486 if (AMI.distributed() && cacheNeighbourField() && AMI.comm() != -1)
487 {
488 if (commsType != UPstream::commsTypes::nonBlocking)
489 {
490 // Invalidate old field - or flag as fatal?
491 patchNeighbourFieldPtr_.reset(nullptr);
492 return;
493 }
494
495 // Start sending
497 << "cyclicACMIFvPatchField::initEvaluate() :"
498 << " field:" << this->internalField().name()
499 << " patch:" << this->patch().name()
500 << " starting send&receive"
501 << endl;
502
503 // Bypass polyPatch to get nbrId.
504 // - use cyclicACMIFvPatch::neighbPatch() virtual instead
505 const cyclicACMIFvPatch& neighbPatch = cyclicACMIPatch_.neighbPatch();
506 const labelUList& nbrFaceCells = neighbPatch.faceCells();
507 const Field<Type> pnf(this->primitiveField(), nbrFaceCells);
508
509 // Assert that all receives are known to have finished
510 if (!recvRequests_.empty() || !recvRequests1_.empty())
511 {
513 << "Outstanding recv request(s) on patch "
514 << cyclicACMIPatch_.name()
515 << " field " << this->internalField().name()
516 << abort(FatalError);
517 }
518
519 // Assume that sends are also OK
520 sendRequests_.clear();
521 sendRequests1_.clear();
522
523 cyclicACMIPatch_.initInterpolate
524 (
525 pnf,
526 sendRequests_,
527 recvRequests_,
528 sendBufs_,
529 recvBufs_,
530
531 sendRequests1_,
532 recvRequests1_,
533 sendBufs1_,
534 recvBufs1_
535 );
536 }
537}
538
539
540template<class Type>
542(
543 const Pstream::commsTypes commsType
544)
545{
546 if (!this->updated())
547 {
548 this->updateCoeffs();
549 }
550
551 const auto& AMI = this->ownerAMI();
552
553 if (AMI.distributed() && cacheNeighbourField() && AMI.comm() != -1)
554 {
555 // Calculate patchNeighbourField
556 if (commsType != UPstream::commsTypes::nonBlocking)
557 {
559 << "Can only evaluate distributed AMI with nonBlocking"
560 << exit(FatalError);
561 }
562
563 patchNeighbourFieldPtr_.reset(nullptr);
564
565 if (!this->ready())
566 {
568 << "Outstanding recv request(s) on patch "
569 << cyclicACMIPatch_.name()
570 << " field " << this->internalField().name()
571 << abort(FatalError);
572 }
573
574 patchNeighbourFieldPtr_.reset
575 (
576 cyclicACMIPatch_.interpolate
577 (
578 Field<Type>::null(), // Not used for distributed
579 recvRequests_,
580 recvBufs_,
581 recvRequests1_,
582 recvBufs1_
583 ).ptr()
584 );
585
586 // Receive requests all handled by last function call
587 recvRequests_.clear();
588 recvRequests1_.clear();
589
590 if (doTransform())
591 {
592 // In-place transform
593 auto& pnf = *patchNeighbourFieldPtr_;
594 transform(pnf, forwardT(), pnf);
595 }
596 }
598 // Use patchNeighbourField() and patchInternalField() to obtain face value
600}
601
602
603template<class Type>
605(
606 solveScalarField& result,
607 const bool add,
608 const lduAddressing& lduAddr,
609 const label patchId,
610 const solveScalarField& psiInternal,
611 const scalarField& coeffs,
612 const direction cmpt,
613 const Pstream::commsTypes commsType
614) const
615{
616 const auto& AMI = this->ownerAMI();
617
618 if (AMI.distributed() && AMI.comm() != -1)
619 {
620 // Start sending
621 if (commsType != UPstream::commsTypes::nonBlocking)
622 {
624 << "Can only evaluate distributed AMI with nonBlocking"
625 << exit(FatalError);
626 }
627
629 << "cyclicACMIFvPatchField::initInterfaceMatrixUpdate() :"
630 << " field:" << this->internalField().name()
631 << " patch:" << this->patch().name()
632 << " starting send&receive"
633 << endl;
634
635 const labelUList& nbrFaceCells =
636 lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
637
638 solveScalarField pnf(psiInternal, nbrFaceCells);
639
640 // Transform according to the transformation tensors
641 transformCoupleField(pnf, cmpt);
642
643 // Assert that all receives are known to have finished
644 if (!recvRequests_.empty() || !recvRequests1_.empty())
645 {
647 << "Outstanding recv request(s) on patch "
648 << cyclicACMIPatch_.name()
649 << " field " << this->internalField().name()
650 << abort(FatalError);
651 }
652
653 // Assume that sends are also OK
654 sendRequests_.clear();
655 sendRequests1_.clear();
656
657 cyclicACMIPatch_.initInterpolate
658 (
659 pnf,
660 sendRequests_,
661 recvRequests_,
662 scalarSendBufs_,
663 scalarRecvBufs_,
664
665 sendRequests1_,
666 recvRequests1_,
667 scalarSendBufs1_,
668 scalarRecvBufs1_
669 );
671
672 this->updatedMatrix(false);
673}
674
675
676template<class Type>
678(
679 solveScalarField& result,
680 const bool add,
681 const lduAddressing& lduAddr,
682 const label patchId,
683 const solveScalarField& psiInternal,
684 const scalarField& coeffs,
685 const direction cmpt,
686 const Pstream::commsTypes commsType
687) const
688{
689 DebugPout<< "cyclicACMIFvPatchField::updateInterfaceMatrix() :"
690 << " field:" << this->internalField().name()
691 << " patch:" << this->patch().name()
692 << endl;
693
694 // note: only applying coupled contribution
695
696 const labelUList& faceCells = lduAddr.patchAddr(patchId);
697
699
700 const auto& AMI = this->ownerAMI();
701
702 if (AMI.distributed() && AMI.comm() != -1)
703 {
704 if (commsType != UPstream::commsTypes::nonBlocking)
705 {
707 << "Can only evaluate distributed AMI with nonBlocking"
708 << exit(FatalError);
709 }
710
712 << "cyclicACMIFvPatchField::evaluate() :"
713 << " field:" << this->internalField().name()
714 << " patch:" << this->patch().name()
715 << " consuming received coupled neighbourfield"
716 << endl;
717
718 pnf =
719 cyclicACMIPatch_.interpolate
720 (
721 solveScalarField::null(), // Not used for distributed
722 recvRequests_,
723 scalarRecvBufs_,
724 recvRequests1_,
725 scalarRecvBufs1_
726 );
727
728 // Receive requests all handled by last function call
729 recvRequests_.clear();
730 recvRequests1_.clear();
731 }
732 else
733 {
734 const labelUList& nbrFaceCells =
735 lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
736
737 pnf = solveScalarField(psiInternal, nbrFaceCells);
738
739 // Transform according to the transformation tensors
740 transformCoupleField(pnf, cmpt);
741
742 pnf = cyclicACMIPatch_.interpolate(pnf);
743 }
744
745 this->addToInternalField(result, !add, faceCells, coeffs, pnf);
746
747 this->updatedMatrix(true);
748}
749
750
751template<class Type>
753(
754 Field<Type>& result,
755 const bool add,
756 const lduAddressing& lduAddr,
757 const label patchId,
758 const Field<Type>& psiInternal,
759 const scalarField& coeffs,
760 const Pstream::commsTypes commsType
761) const
762{
763 const auto& AMI = this->ownerAMI();
764
765 if (AMI.distributed() && AMI.comm() != -1)
766 {
767 if (commsType != UPstream::commsTypes::nonBlocking)
768 {
770 << "Can only evaluate distributed AMI with nonBlocking"
771 << exit(FatalError);
772 }
773
774 const labelUList& nbrFaceCells =
775 lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
776
777 Field<Type> pnf(psiInternal, nbrFaceCells);
778
779 // Transform according to the transformation tensors
780 transformCoupleField(pnf);
781
782 // Assert that all receives are known to have finished
783 if (!recvRequests_.empty() || !recvRequests1_.empty())
784 {
786 << "Outstanding recv request(s) on patch "
787 << cyclicACMIPatch_.name()
788 << " field " << this->internalField().name()
789 << abort(FatalError);
790 }
791
792 // Assume that sends are also OK
793 sendRequests_.clear();
794 sendRequests1_.clear();
795
796 cyclicACMIPatch_.initInterpolate
797 (
798 pnf,
799 sendRequests_,
800 recvRequests_,
801 sendBufs_,
802 recvBufs_,
803
804 sendRequests1_,
805 recvRequests1_,
806 sendBufs1_,
807 recvBufs1_
808 );
810
811 this->updatedMatrix(false);
812}
813
814
815template<class Type>
817(
818 Field<Type>& result,
819 const bool add,
820 const lduAddressing& lduAddr,
821 const label patchId,
822 const Field<Type>& psiInternal,
823 const scalarField& coeffs,
824 const Pstream::commsTypes commsType
825) const
826{
827 // note: only applying coupled contribution
828
829 const labelUList& faceCells = lduAddr.patchAddr(patchId);
830
831 const auto& AMI = this->ownerAMI();
832
833 Field<Type> pnf;
834
835 if (AMI.distributed() && AMI.comm() != -1)
836 {
837 if (commsType != UPstream::commsTypes::nonBlocking)
838 {
840 << "Can only evaluate distributed AMI with nonBlocking"
841 << exit(FatalError);
842 }
843
844 pnf =
845 cyclicACMIPatch_.interpolate
846 (
847 Field<Type>::null(), // Not used for distributed
848 recvRequests_,
849 recvBufs_,
850 recvRequests1_,
851 recvBufs1_
852 );
853
854 // Receive requests all handled by last function call
855 recvRequests_.clear();
856 recvRequests1_.clear();
857 }
858 else
859 {
860 const labelUList& nbrFaceCells =
861 lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
862
863 pnf = Field<Type>(psiInternal, nbrFaceCells);
864
865 // Transform according to the transformation tensors
866 transformCoupleField(pnf);
867
868 pnf = cyclicACMIPatch_.interpolate(pnf);
869 }
870
871 this->addToInternalField(result, !add, faceCells, coeffs, pnf);
872
873 this->updatedMatrix(true);
874}
875
876
877template<class Type>
879(
880 fvMatrix<Type>& matrix
881)
882{
883 const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
884
885 // Nothing to be done by the AMI, but re-direct to non-overlap patch
886 // with non-overlap patch weights
888
889 const_cast<fvPatchField<Type>&>(npf).manipulateMatrix(matrix, 1.0 - mask);
890}
891
892
893template<class Type>
895(
896 fvMatrix<Type>& matrix,
897 const label mat,
898 const direction cmpt
899)
900{
901 if (this->cyclicACMIPatch().owner())
902 {
903 label index = this->patch().index();
904
905 const label globalPatchID =
906 matrix.lduMeshAssembly().patchLocalToGlobalMap()[mat][index];
907
908 const Field<scalar> intCoeffsCmpt
909 (
910 matrix.internalCoeffs()[globalPatchID].component(cmpt)
911 );
912
913 const Field<scalar> boundCoeffsCmpt
914 (
915 matrix.boundaryCoeffs()[globalPatchID].component(cmpt)
916 );
917
918 tmp<Field<scalar>> tintCoeffs(coeffs(matrix, intCoeffsCmpt, mat));
919 tmp<Field<scalar>> tbndCoeffs(coeffs(matrix, boundCoeffsCmpt, mat));
920 const Field<scalar>& intCoeffs = tintCoeffs.ref();
921 const Field<scalar>& bndCoeffs = tbndCoeffs.ref();
922
923 const labelUList& u = matrix.lduAddr().upperAddr();
924 const labelUList& l = matrix.lduAddr().lowerAddr();
925
926 label subFaceI = 0;
927
928 const labelList& faceMap =
929 matrix.lduMeshAssembly().faceBoundMap()[mat][index];
930
931 forAll (faceMap, j)
932 {
933 label globalFaceI = faceMap[j];
934
935 const scalar boundCorr = -bndCoeffs[subFaceI];
936 const scalar intCorr = -intCoeffs[subFaceI];
937
938 matrix.upper()[globalFaceI] += boundCorr;
939 matrix.diag()[u[globalFaceI]] -= intCorr;
940 matrix.diag()[l[globalFaceI]] -= boundCorr;
941
942 if (matrix.asymmetric())
943 {
944 matrix.lower()[globalFaceI] += intCorr;
945 }
946 subFaceI++;
947 }
948
949 // Set internalCoeffs and boundaryCoeffs in the assembly matrix
950 // on cyclicACMI patches to be used in the individual matrix by
951 // matrix.flux()
952 if (matrix.psi(mat).mesh().fluxRequired(this->internalField().name()))
953 {
954 matrix.internalCoeffs().set
955 (
956 globalPatchID, intCoeffs*pTraits<Type>::one
957 );
958 matrix.boundaryCoeffs().set
959 (
960 globalPatchID, bndCoeffs*pTraits<Type>::one
961 );
962
963 const label nbrPathID =
964 cyclicACMIPatch_.cyclicACMIPatch().neighbPatchID();
965
966 const label nbrGlobalPatchID =
968 [mat][nbrPathID];
969
970 matrix.internalCoeffs().set
971 (
972 nbrGlobalPatchID, intCoeffs*pTraits<Type>::one
973 );
974 matrix.boundaryCoeffs().set
975 (
976 nbrGlobalPatchID, bndCoeffs*pTraits<Type>::one
977 );
978 }
979 }
980}
981
982
983template<class Type>
984Foam::tmp<Foam::Field<Foam::scalar>>
985Foam::cyclicACMIFvPatchField<Type>::coeffs
986(
987 fvMatrix<Type>& matrix,
988 const Field<scalar>& coeffs,
989 const label mat
990) const
991{
992 const label index(this->patch().index());
993
994 const label nSubFaces
995 (
996 matrix.lduMeshAssembly().cellBoundMap()[mat][index].size()
997 );
998
999 auto tmapCoeffs = tmp<Field<scalar>>::New(nSubFaces, Zero);
1000 auto& mapCoeffs = tmapCoeffs.ref();
1001
1002 const scalarListList& srcWeight =
1003 cyclicACMIPatch_.cyclicACMIPatch().AMI().srcWeights();
1004
1005 const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
1006
1007 const scalar tol = cyclicACMIPolyPatch::tolerance();
1008 label subFaceI = 0;
1009 forAll(mask, faceI)
1010 {
1011 const scalarList& w = srcWeight[faceI];
1012 for(label i=0; i<w.size(); i++)
1013 {
1014 if (mask[faceI] > tol)
1015 {
1016 const label localFaceId =
1017 matrix.lduMeshAssembly().facePatchFaceMap()
1018 [mat][index][subFaceI];
1019 mapCoeffs[subFaceI] = w[i]*coeffs[localFaceId];
1020 }
1021 subFaceI++;
1022 }
1024
1025 return tmapCoeffs;
1026}
1027
1028
1029template<class Type>
1031{
1032 // Update non-overlap patch - some will implement updateCoeffs, and
1033 // others will implement evaluate
1034
1035 // Pass in (1 - mask) to give non-overlap patch the chance to do
1036 // manipulation of non-face based data
1037
1038 const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
1040 const_cast<fvPatchField<Type>&>(npf).updateWeightedCoeffs(1.0 - mask);
1041}
1042
1043
1044template<class Type>
1046{
1049
1050 if (patchNeighbourFieldPtr_)
1051 {
1052 patchNeighbourFieldPtr_->writeEntry("neighbourValue", os);
1053 }
1054}
1055
1056
1057// ************************************************************************* //
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))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const noexcept
Return const reference to mesh.
static int localBoundaryConsistency() noexcept
Get flag for local boundary consistency checks.
Definition Field.H:133
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
static const Field< Type > & null() noexcept
Return a null Field (reference to a nullObject). Behaves like an empty Field.
Definition Field.H:192
constexpr Field() noexcept
Default construct.
Definition FieldI.H:24
Generic GeometricField class.
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
fileName objectPath() const
The complete path + object name.
Definition IOobjectI.H:313
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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 deepCopy(const UList< T > &list)
Copy elements of the given UList. Sizes must match!
Definition UList.C:95
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool finishedRequests(label pos, label len=-1)
Non-blocking comms: have all requests (from position onwards) finished? Corresponds to MPI_Testall().
commsTypes
Communications types.
Definition UPstream.H:81
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
Abstract base class for coupled patches.
coupledFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
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 bool doTransform() const
Does the patch field perform the transformation.
virtual void initInterfaceMatrixUpdate(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
Initialise neighbour matrix update.
const cyclicACMIFvPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic AMI patch.
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.
cyclicACMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
virtual void write(Ostream &os) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
virtual bool ready() const
Are all (receive) data available?
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
virtual const cyclicACMIFvPatch & neighbPatch() const
Return neighbour fvPatch.
Abstract base class for cyclic ACMI coupled interfaces.
void transformCoupleField(Field< Type > &f) const
Transform given patch field.
static scalar tolerance()
Overlap tolerance.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition fvMatrix.H:549
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition fvMatrix.H:567
const lduPrimitiveMeshAssembly & lduMeshAssembly()
Return optional lduAdressing.
Definition fvMatrix.H:478
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
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...
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
@ LITERAL
String literal.
Definition keyType.H:82
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.
bool updatedMatrix() const noexcept
Whether matrix has been updated.
void addToInternalField(Field< Type > &result, const bool add, const labelUList &faceCells, const scalarField &coeffs, const Field< Type > &vals) const
Add/subtract weighted contributions to internal field.
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition lduMatrix.H:769
bool asymmetric() const noexcept
Matrix is asymmetric (ie, full).
Definition lduMatrix.H:850
const scalarField & diag() const
Definition lduMatrix.C:195
const scalarField & upper() const
Definition lduMatrix.C:235
const scalarField & lower() const
Definition lduMatrix.C:306
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 traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
bool fluxRequired(const word &name) const
Get flux-required for given name, or default.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
label patchId(-1)
#define DebugPout
Report an information message using Foam::Pout.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< label > labelList
A List of labels.
Definition List.H:62
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
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
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dict add("bounds", meshBb)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Spatial transformation functions for primitive fields.