Loading...
Searching...
No Matches
cyclicAMIFvPatchField.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) 2011-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 "cyclicAMIPolyPatch.H"
30#include "mapDistributeBase.H"
32#include "fvMatrix.H"
33#include "volFields.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37template<class Type>
39(
40 const fvPatch& p,
42)
43:
46 cyclicAMIPatch_(refCast<const cyclicAMIFvPatch>(p)),
47 patchNeighbourFieldPtr_(nullptr)
48{}
49
50
51template<class Type>
53(
54 const fvPatch& p,
56 const dictionary& dict
57)
58:
60 coupledFvPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
61 cyclicAMIPatch_(refCast<const cyclicAMIFvPatch>(p, dict)),
62 patchNeighbourFieldPtr_(nullptr)
63{
65 {
67 << "\n patch type '" << p.type()
68 << "' not constraint type '" << typeName << "'"
69 << "\n for patch " << p.name()
70 << " of field " << this->internalField().name()
71 << " in file " << this->internalField().objectPath()
73 }
74
75 if (cacheNeighbourField())
76 {
77 // Handle neighbour value first, before any evaluate()
78 const auto* hasNeighbValue =
79 dict.findEntry("neighbourValue", keyType::LITERAL);
80
81 if (hasNeighbValue)
82 {
83 patchNeighbourFieldPtr_.reset
84 (
85 new Field<Type>(*hasNeighbValue, p.size())
86 );
87 }
88 }
89
90 // Use 'value' supplied, or evaluate (if coupled) or set to internal field
91 if (!this->readValueEntry(dict))
92 {
93 if (this->coupled())
94 {
95 // Tricky: avoid call to evaluate without call to initEvaluate.
96 // For now just disable the localConsistency to make it use the
97 // old logic (ultimately calls the fully self contained
98 // patchNeighbourField)
99
100 const auto oldConsistency = FieldBase::localBoundaryConsistency(0);
101
103
105 }
106 else
107 {
108 this->extrapolateInternal(); // Zero-gradient patch values
109 }
110 }
111}
112
113
114template<class Type>
116(
118 const fvPatch& p,
120 const fvPatchFieldMapper& mapper
121)
122:
124 coupledFvPatchField<Type>(ptf, p, iF, mapper),
125 cyclicAMIPatch_(refCast<const cyclicAMIFvPatch>(p)),
126 patchNeighbourFieldPtr_(nullptr)
127{
128 //if (ptf.patchNeighbourFieldPtr_ && cacheNeighbourField())
129 //{
130 // patchNeighbourFieldPtr_.reset
131 // (
132 // new Field<Type>(ptf.patchNeighbourFieldPtr_(), mapper)
133 // );
134 //}
135
136 if (!isA<cyclicAMIFvPatch>(this->patch()))
137 {
139 << "\n patch type '" << p.type()
140 << "' not constraint type '" << typeName << "'"
141 << "\n for patch " << p.name()
142 << " of field " << this->internalField().name()
143 << " in file " << this->internalField().objectPath()
144 << exit(FatalError);
145 }
146 if (debug && !ptf.all_ready())
147 {
148 FatalErrorInFunction
149 << "Outstanding request(s) on patch " << cyclicAMIPatch_.name()
150 << abort(FatalError);
151 }
152}
153
154
155template<class Type>
157(
159)
160:
162 coupledFvPatchField<Type>(ptf),
163 cyclicAMIPatch_(ptf.cyclicAMIPatch_),
164 patchNeighbourFieldPtr_(nullptr)
165{
166 if (debug && !ptf.all_ready())
167 {
168 FatalErrorInFunction
169 << "Outstanding request(s) on patch " << cyclicAMIPatch_.name()
170 << abort(FatalError);
171 }
172}
173
174
175template<class Type>
177(
180)
181:
183 coupledFvPatchField<Type>(ptf, iF),
184 cyclicAMIPatch_(ptf.cyclicAMIPatch_),
185 patchNeighbourFieldPtr_(nullptr)
186{
187 if (debug && !ptf.all_ready())
188 {
189 FatalErrorInFunction
190 << "Outstanding request(s) on patch " << cyclicAMIPatch_.name()
191 << abort(FatalError);
192 }
193}
194
195
196// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197
198template<class Type>
199bool Foam::cyclicAMIFvPatchField<Type>::all_ready() const
200{
201 int done = 0;
202
203 if
204 (
206 (
207 recvRequests_.start(),
208 recvRequests_.size()
209 )
211 (
212 recvRequests1_.start(),
213 recvRequests1_.size()
214 )
215 )
216 {
217 recvRequests_.clear();
218 recvRequests1_.clear();
219 ++done;
220 }
221
222 if
223 (
225 (
226 sendRequests_.start(),
227 sendRequests_.size()
228 )
230 (
231 sendRequests1_.start(),
232 sendRequests1_.size()
233 )
235 {
236 sendRequests_.clear();
237 sendRequests1_.clear();
238 ++done;
240
241 return (done == 2);
242}
243
244
245template<class Type>
247{
248 if
249 (
251 (
252 recvRequests_.start(),
253 recvRequests_.size()
254 )
256 (
257 recvRequests1_.start(),
258 recvRequests1_.size()
259 )
260 )
261 {
262 recvRequests_.clear();
263 recvRequests1_.clear();
264
265 if
266 (
268 (
269 sendRequests_.start(),
270 sendRequests_.size()
271 )
273 (
274 sendRequests1_.start(),
275 sendRequests1_.size()
276 )
277 )
278 {
279 sendRequests_.clear();
280 sendRequests1_.clear();
281 }
282
283 return true;
284 }
286 return false;
287}
288
289
290
291template<class Type>
293(
294 const fvPatchFieldMapper& mapper
295)
298 patchNeighbourFieldPtr_.reset(nullptr);
299}
300
301
302template<class Type>
304(
305 const fvPatchField<Type>& ptf,
306 const labelList& addr
307)
308{
310 patchNeighbourFieldPtr_.reset(nullptr);
311}
312
313
314template<class Type>
316Foam::cyclicAMIFvPatchField<Type>::getNeighbourField
317(
318 const UList<Type>& internalData
319) const
320{
321 // By pass polyPatch to get nbrId. Instead use cyclicAMIFvPatch virtual
322 // neighbPatch()
323 const auto& neighbPatch = cyclicAMIPatch_.neighbPatch();
324 const labelUList& nbrFaceCells = neighbPatch.faceCells();
325
326 Field<Type> pnf(internalData, nbrFaceCells);
327 Field<Type> defaultValues;
328
329 if (cyclicAMIPatch_.applyLowWeightCorrection())
330 {
331 defaultValues = Field<Type>(internalData, cyclicAMIPatch_.faceCells());
332 }
334 tmp<Field<Type>> tpnf = cyclicAMIPatch_.interpolate(pnf, defaultValues);
335
336 if (doTransform())
337 {
338 transform(tpnf.ref(), forwardT(), tpnf());
339 }
340
341 return tpnf;
342}
343
344
345template<class Type>
346bool Foam::cyclicAMIFvPatchField<Type>::cacheNeighbourField() const
347{
348 // const auto& AMI = this->ownerAMI();
350 // if (AMI.cacheActive())
351 // {
352 // return false;
353 // }
354 // else
355 {
357 }
358}
359
360
361template<class Type>
363Foam::cyclicAMIFvPatchField<Type>::getPatchNeighbourField
364(
365 bool checkCommunicator
366) const
367{
368 const auto& AMI = this->ownerAMI();
369
370 if
371 (
372 AMI.distributed() && cacheNeighbourField()
373 && (!checkCommunicator || AMI.comm() != -1)
375 {
376 if (!this->ready())
377 {
379 << "Outstanding recv request(s) on patch "
380 << cyclicAMIPatch_.name()
381 << " field " << this->internalField().name()
382 << abort(FatalError);
383 }
384
385 const auto& fvp = this->patch();
386 const auto& mesh = fvp.boundaryMesh().mesh();
387
388 if
390 patchNeighbourFieldPtr_
391 && !mesh.upToDatePoints(this->internalField())
392 )
393 {
394 //DebugPout
395 // << "cyclicAMIFvPatchField::patchNeighbourField() :"
396 // << " field:" << this->internalField().name()
397 // << " patch:" << fvp.name()
398 // << " CLEARING patchNeighbourField"
399 // << endl;
400 patchNeighbourFieldPtr_.reset(nullptr);
401 }
402
403 // Initialise if not done in construct-from-dictionary
404 if (!patchNeighbourFieldPtr_)
405 {
406 //DebugPout
407 // << "cyclicAMIFvPatchField::patchNeighbourField() :"
408 // << " field:" << this->internalField().name()
409 // << " patch:" << fvp.name()
410 // << " caching patchNeighbourField"
411 // << endl;
412
413 // Do interpolation and store result
414 patchNeighbourFieldPtr_.reset
415 (
416 getNeighbourField(this->primitiveField()).ptr()
417 );
418 }
419 else
420 {
421 // Have cached value. Check
422 //if (debug)
423 //{
424 // tmp<Field<Type>> tpnf
425 // (
426 // getNeighbourField(this->primitiveField())
427 // );
428 // if (tpnf() != patchNeighbourFieldPtr_())
429 // {
430 // FatalErrorInFunction
431 // << "On field " << this->internalField().name()
432 // << " patch " << fvp.name() << endl
433 // << "Cached patchNeighbourField :"
434 // << flatOutput(patchNeighbourFieldPtr_()) << endl
435 // << "Calculated patchNeighbourField:"
436 // << flatOutput(tpnf()) << exit(FatalError);
437 // }
438 //}
439 }
440
441 return patchNeighbourFieldPtr_();
442 }
443 else
444 {
445 // Do interpolation
446 return getNeighbourField(this->primitiveField());
447 }
448}
449
450
451template<class Type>
455 // checkCommunicator = true
456 return this->getPatchNeighbourField(true);
457}
458
459
460template<class Type>
462(
463 UList<Type>& pnf
464) const
465{
466 // checkCommunicator = false
467 auto tpnf = this->getPatchNeighbourField(false);
468 pnf.deepCopy(tpnf());
469}
470
471
472template<class Type>
475{
476 const auto& fld =
478 (
479 this->primitiveField()
480 );
481
484 fld.boundaryField()[cyclicAMIPatch_.neighbPatchID()]
485 );
486}
487
488
489template<class Type>
492 const Pstream::commsTypes commsType
493)
494{
495 if (!this->updated())
496 {
497 this->updateCoeffs();
498 }
499
500 const auto& AMI = this->ownerAMI();
501
502 if (AMI.distributed() && cacheNeighbourField() && AMI.comm() != -1)
503 {
504 //DebugPout
505 // << "*** cyclicAMIFvPatchField::initEvaluate() :"
506 // << " field:" << this->internalField().name()
507 // << " patch:" << this->patch().name()
508 // << " sending patchNeighbourField"
509 // << endl;
510
511 if (commsType != UPstream::commsTypes::nonBlocking)
512 {
513 // Invalidate old field - or flag as fatal?
514 patchNeighbourFieldPtr_.reset(nullptr);
515 return;
516 }
517
518 // Start sending
519
520 // Bypass polyPatch to get nbrId.
521 // - use cyclicACMIFvPatch::neighbPatch() virtual instead
522 const cyclicAMIFvPatch& neighbPatch = cyclicAMIPatch_.neighbPatch();
523 const labelUList& nbrFaceCells = neighbPatch.faceCells();
524 const Field<Type> pnf(this->primitiveField(), nbrFaceCells);
525
526 const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
527
528 // Assert that all receives are known to have finished
529 if (!recvRequests_.empty() || !recvRequests1_.empty())
530 {
532 << "Outstanding recv request(s) on patch "
533 << cyclicAMIPatch_.name()
534 << " field " << this->internalField().name()
535 << abort(FatalError);
536 }
537
538 // Assume that sends are also OK
539 sendRequests_.clear();
540 sendRequests1_.clear();
541
543 (
544 pnf,
545 sendRequests_,
546 recvRequests_,
547 sendBufs_,
548 recvBufs_,
549
550 sendRequests1_,
551 recvRequests1_,
552 sendBufs1_,
553 recvBufs1_
554 );
555 }
556}
557
558
559template<class Type>
561(
562 const Pstream::commsTypes commsType
563)
564{
565 if (!this->updated())
566 {
567 this->updateCoeffs();
568 }
569
570 const auto& AMI = this->ownerAMI();
571
572 if (AMI.distributed() && cacheNeighbourField() && AMI.comm() != -1)
573 {
574 // Calculate patchNeighbourField
575 if (commsType != UPstream::commsTypes::nonBlocking)
576 {
578 << "Can only evaluate distributed AMI with nonBlocking"
579 << exit(FatalError);
580 }
581
582 patchNeighbourFieldPtr_.reset(nullptr);
583
584 const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
585
586 Field<Type> defaultValues;
587 if (AMI.applyLowWeightCorrection())
588 {
589 defaultValues = this->patchInternalField();
590 }
591
592 //DebugPout
593 // << "*** cyclicAMIFvPatchField::evaluate() :"
594 // << " field:" << this->internalField().name()
595 // << " patch:" << this->patch().name()
596 // << " receiving&caching patchNeighbourField"
597 // << endl;
598
599 patchNeighbourFieldPtr_.reset
600 (
601 cpp.interpolate
602 (
603 Field<Type>::null(), // Not used for distributed
604 recvRequests_,
605 recvBufs_,
606 recvRequests1_,
607 recvBufs1_,
608 defaultValues
609 ).ptr()
610 );
611
612 // Receive requests all handled by last function call
613 recvRequests_.clear();
614 recvRequests1_.clear();
615
616 if (doTransform())
617 {
618 // In-place transform
619 auto& pnf = *patchNeighbourFieldPtr_;
620 transform(pnf, forwardT(), pnf);
621 }
622 }
624 // Use patchNeighbourField() and patchInternalField() to obtain face value
626}
627
628
629template<class Type>
631(
632 solveScalarField& result,
633 const bool add,
634 const lduAddressing& lduAddr,
635 const label patchId,
636 const solveScalarField& psiInternal,
637 const scalarField& coeffs,
638 const direction cmpt,
639 const Pstream::commsTypes commsType
640) const
641{
642 const auto& AMI = this->ownerAMI();
643
644 if (AMI.distributed() && AMI.comm() != -1)
645 {
646 // Start sending
647 if (commsType != UPstream::commsTypes::nonBlocking)
648 {
650 << "Can only evaluate distributed AMI with nonBlocking"
651 << exit(FatalError);
652 }
653
654 const labelUList& nbrFaceCells =
655 lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
656
657 solveScalarField pnf(psiInternal, nbrFaceCells);
658
659 // Transform according to the transformation tensors
660 transformCoupleField(pnf, cmpt);
661
662 const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
663
664 // Assert that all receives are known to have finished
665 if (!recvRequests_.empty() || !recvRequests1_.empty())
666 {
668 << "Outstanding recv request(s) on patch "
669 << cyclicAMIPatch_.name()
670 << " field " << this->internalField().name()
671 << abort(FatalError);
672 }
673
674 // Assume that sends are also OK
675 sendRequests_.clear();
676 sendRequests1_.clear();
677
678 cpp.initInterpolate
679 (
680 pnf,
681 sendRequests_,
682 recvRequests_,
683 scalarSendBufs_,
684 scalarRecvBufs_,
685
686 sendRequests1_,
687 recvRequests1_,
688 scalarSendBufs1_,
689 scalarRecvBufs1_
690 );
692
693 this->updatedMatrix(false);
694}
695
696
697template<class Type>
699(
700 solveScalarField& result,
701 const bool add,
702 const lduAddressing& lduAddr,
703 const label patchId,
704 const solveScalarField& psiInternal,
705 const scalarField& coeffs,
706 const direction cmpt,
707 const Pstream::commsTypes commsType
708) const
709{
710 //DebugPout<< "cyclicAMIFvPatchField::updateInterfaceMatrix() :"
711 // << " field:" << this->internalField().name()
712 // << " patch:" << this->patch().name()
713 // << endl;
714
715 const labelUList& faceCells = lduAddr.patchAddr(patchId);
716
717 const auto& AMI = this->ownerAMI();
718
720
721 if (AMI.distributed() && AMI.comm() != -1)
722 {
723 if (commsType != UPstream::commsTypes::nonBlocking)
724 {
726 << "Can only evaluate distributed AMI with nonBlocking"
727 << exit(FatalError);
728 }
729
730 solveScalarField defaultValues;
731 if (AMI.applyLowWeightCorrection())
732 {
733 defaultValues = solveScalarField(psiInternal, faceCells);
734 }
735
736 const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
737
738 pnf =
739 cpp.interpolate
740 (
741 solveScalarField::null(), // Not used for distributed
742 recvRequests_,
743 scalarRecvBufs_,
744 recvRequests1_,
745 scalarRecvBufs1_,
746 defaultValues
747 );
748
749 // Receive requests all handled by last function call
750 recvRequests_.clear();
751 recvRequests1_.clear();
752 }
753 else
754 {
755 solveScalarField defaultValues;
756 if (cyclicAMIPatch_.applyLowWeightCorrection())
757 {
758 defaultValues = solveScalarField(psiInternal, faceCells);
759 }
760
761 const labelUList& nbrFaceCells =
762 lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
763
764 pnf = solveScalarField(psiInternal, nbrFaceCells);
765
766 // Transform according to the transformation tensors
767 transformCoupleField(pnf, cmpt);
768
769 pnf = cyclicAMIPatch_.interpolate(pnf, defaultValues);
770 }
771
772 // Multiply the field by coefficients and add into the result
773 this->addToInternalField(result, !add, faceCells, coeffs, pnf);
774
775 this->updatedMatrix(true);
776}
777
778
779template<class Type>
781(
782 Field<Type>& result,
783 const bool add,
784 const lduAddressing& lduAddr,
785 const label patchId,
786 const Field<Type>& psiInternal,
787 const scalarField& coeffs,
788 const Pstream::commsTypes commsType
789) const
790{
791 const auto& AMI = this->ownerAMI();
792
793 if (AMI.distributed() && AMI.comm() != -1)
794 {
795 if (commsType != UPstream::commsTypes::nonBlocking)
796 {
798 << "Can only evaluate distributed AMI with nonBlocking"
799 << exit(FatalError);
800 }
801
802 const labelUList& nbrFaceCells =
803 lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
804
805 Field<Type> pnf(psiInternal, nbrFaceCells);
806
807 // Transform according to the transformation tensors
808 transformCoupleField(pnf);
809
810 const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
811
812 // Assert that all receives are known to have finished
813 if (!recvRequests_.empty() || !recvRequests1_.empty())
814 {
816 << "Outstanding recv request(s) on patch "
817 << cyclicAMIPatch_.name()
818 << " field " << this->internalField().name()
819 << abort(FatalError);
820 }
821
822 // Assume that sends are also OK
823 sendRequests_.clear();
824 sendRequests1_.clear();
825
826 cpp.initInterpolate
827 (
828 pnf,
829 sendRequests_,
830 recvRequests_,
831 sendBufs_,
832 recvBufs_,
833
834 sendRequests1_,
835 recvRequests1_,
836 sendBufs1_,
837 recvBufs1_
838 );
840
841 this->updatedMatrix(false);
842}
843
844
845template<class Type>
847(
848 Field<Type>& result,
849 const bool add,
850 const lduAddressing& lduAddr,
851 const label patchId,
852 const Field<Type>& psiInternal,
853 const scalarField& coeffs,
854 const Pstream::commsTypes commsType
855) const
856{
857 //DebugPout<< "cyclicAMIFvPatchField::updateInterfaceMatrix() :"
858 // << " field:" << this->internalField().name()
859 // << " patch:" << this->patch().name()
860 // << endl;
861
862 const labelUList& faceCells = lduAddr.patchAddr(patchId);
863
864 const auto& AMI = this->ownerAMI();
865
866 Field<Type> pnf;
867
868 if (AMI.distributed() && AMI.comm() != -1)
869 {
870 if (commsType != UPstream::commsTypes::nonBlocking)
871 {
873 << "Can only evaluate distributed AMI with nonBlocking"
874 << exit(FatalError);
875 }
876
877 const cyclicAMIPolyPatch& cpp = cyclicAMIPatch_.cyclicAMIPatch();
878
879 Field<Type> defaultValues;
880 if (AMI.applyLowWeightCorrection())
881 {
882 defaultValues = Field<Type>(psiInternal, faceCells);
883 }
884
885 pnf =
886 cpp.interpolate
887 (
888 Field<Type>::null(), // Not used for distributed
889 recvRequests_,
890 recvBufs_,
891 recvRequests1_,
892 recvBufs1_,
893 defaultValues
894 );
895
896 // Receive requests all handled by last function call
897 recvRequests_.clear();
898 recvRequests1_.clear();
899 }
900 else
901 {
902 const labelUList& nbrFaceCells =
903 lduAddr.patchAddr(cyclicAMIPatch_.neighbPatchID());
904
905 pnf = Field<Type>(psiInternal, nbrFaceCells);
906
907 // Transform according to the transformation tensors
908 transformCoupleField(pnf);
909
910 Field<Type> defaultValues;
911 if (cyclicAMIPatch_.applyLowWeightCorrection())
912 {
913 defaultValues = Field<Type>(psiInternal, faceCells);
914 }
915
916 pnf = cyclicAMIPatch_.interpolate(pnf, defaultValues);
917 }
918
919 // Multiply the field by coefficients and add into the result
920 this->addToInternalField(result, !add, faceCells, coeffs, pnf);
921
922 this->updatedMatrix(true);
923}
924
925
926template<class Type>
928(
929 fvMatrix<Type>& matrix,
930 const label mat,
931 const direction cmpt
932)
933{
934 if (this->cyclicAMIPatch().owner())
935 {
936 const label index = this->patch().index();
937
938 const label globalPatchID =
939 matrix.lduMeshAssembly().patchLocalToGlobalMap()[mat][index];
940
941 const Field<scalar> intCoeffsCmpt
942 (
943 matrix.internalCoeffs()[globalPatchID].component(cmpt)
944 );
945
946 const Field<scalar> boundCoeffsCmpt
947 (
948 matrix.boundaryCoeffs()[globalPatchID].component(cmpt)
949 );
950
951 tmp<Field<scalar>> tintCoeffs(coeffs(matrix, intCoeffsCmpt, mat));
952 tmp<Field<scalar>> tbndCoeffs(coeffs(matrix, boundCoeffsCmpt, mat));
953 const Field<scalar>& intCoeffs = tintCoeffs.ref();
954 const Field<scalar>& bndCoeffs = tbndCoeffs.ref();
955
956 const labelUList& u = matrix.lduAddr().upperAddr();
957 const labelUList& l = matrix.lduAddr().lowerAddr();
958
959 label subFaceI = 0;
960
961 const labelList& faceMap =
962 matrix.lduMeshAssembly().faceBoundMap()[mat][index];
963
964 forAll (faceMap, j)
965 {
966 label globalFaceI = faceMap[j];
967
968 const scalar boundCorr = -bndCoeffs[subFaceI];
969 const scalar intCorr = -intCoeffs[subFaceI];
970
971 matrix.upper()[globalFaceI] += boundCorr;
972 matrix.diag()[u[globalFaceI]] -= intCorr;
973 matrix.diag()[l[globalFaceI]] -= boundCorr;
974
975 if (matrix.asymmetric())
976 {
977 matrix.lower()[globalFaceI] += intCorr;
978 }
979 subFaceI++;
980 }
981
982 // Set internalCoeffs and boundaryCoeffs in the assembly matrix
983 // on cyclicAMI patches to be used in the individual matrix by
984 // matrix.flux()
985 if (matrix.psi(mat).mesh().fluxRequired(this->internalField().name()))
986 {
987 matrix.internalCoeffs().set
988 (
989 globalPatchID, intCoeffs*pTraits<Type>::one
990 );
991 matrix.boundaryCoeffs().set
992 (
993 globalPatchID, bndCoeffs*pTraits<Type>::one
994 );
995
996 const label nbrPathID =
997 cyclicAMIPatch_.cyclicAMIPatch().neighbPatchID();
998
999 const label nbrGlobalPatchID =
1001 [mat][nbrPathID];
1002
1003 matrix.internalCoeffs().set
1004 (
1005 nbrGlobalPatchID, intCoeffs*pTraits<Type>::one
1006 );
1007 matrix.boundaryCoeffs().set
1008 (
1009 nbrGlobalPatchID, bndCoeffs*pTraits<Type>::one
1010 );
1011 }
1012 }
1013}
1014
1015
1016template<class Type>
1017Foam::tmp<Foam::Field<Foam::scalar>>
1018Foam::cyclicAMIFvPatchField<Type>::coeffs
1019(
1020 fvMatrix<Type>& matrix,
1021 const Field<scalar>& coeffs,
1022 const label mat
1023) const
1024{
1025 const label index(this->patch().index());
1026
1027 const label nSubFaces
1028 (
1029 matrix.lduMeshAssembly().cellBoundMap()[mat][index].size()
1030 );
1031
1032 auto tmapCoeffs = tmp<Field<scalar>>::New(nSubFaces, Zero);
1033 auto& mapCoeffs = tmapCoeffs.ref();
1034
1035 const scalarListList& srcWeight =
1036 cyclicAMIPatch_.cyclicAMIPatch().AMI().srcWeights();
1037
1038 label subFaceI = 0;
1039 forAll(*this, faceI)
1040 {
1041 const scalarList& w = srcWeight[faceI];
1042 for(label i=0; i<w.size(); i++)
1043 {
1044 const label localFaceId =
1045 matrix.lduMeshAssembly().facePatchFaceMap()[mat][index][subFaceI];
1046 mapCoeffs[subFaceI] = w[i]*coeffs[localFaceId];
1047 subFaceI++;
1048 }
1049 }
1050
1051 return tmapCoeffs;
1052}
1053
1054
1055template<class Type>
1056template<class Type2>
1057void Foam::cyclicAMIFvPatchField<Type>::collectStencilData
1058(
1059 const refPtr<mapDistribute>& mapPtr,
1060 const labelListList& stencil,
1061 const Type2& data,
1062 List<Type2>& expandedData
1063)
1064{
1065 expandedData.resize_nocopy(stencil.size());
1066 if (mapPtr)
1067 {
1068 Type2 work(data);
1069 mapPtr().distribute(work);
1070
1071 forAll(stencil, facei)
1072 {
1073 const auto& slots = stencil[facei];
1074 expandedData[facei].push_back
1075 (
1077 );
1078 }
1079 }
1080 else
1081 {
1082 forAll(stencil, facei)
1083 {
1084 const auto& slots = stencil[facei];
1085 expandedData[facei].push_back
1086 (
1089 }
1090 }
1091}
1092
1093
1094template<class Type>
1096{
1099
1100 if (patchNeighbourFieldPtr_)
1101 {
1102 patchNeighbourFieldPtr_->writeEntry("neighbourValue", os);
1104}
1105
1106
1107// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1108
1109template<class Type>
1111(
1112 const fvPatchField<Type>& ptf
1113)
1114{
1116
1117 //Pout<< "cyclicAMIFvPatchField::operator= :"
1118 // << " field:" << this->internalField().name()
1119 // << " patch:" << this->patch().name()
1120 // << " copying from field:" << ptf.internalField().name()
1121 // << endl;
1122
1123 const auto* cycPtr = isA<cyclicAMIFvPatchField<Type>>(ptf);
1124
1125 if
1126 (
1127 cycPtr
1128 && cycPtr->patchNeighbourFieldPtr_
1129 && cycPtr->patchNeighbourFieldPtr_->size() == this->size()
1130 )
1131 {
1132 if (!patchNeighbourFieldPtr_)
1133 {
1134 patchNeighbourFieldPtr_ = autoPtr<Field<Type>>::New();
1135 }
1136
1137 // Copy values
1138 *patchNeighbourFieldPtr_ = *(cycPtr->patchNeighbourFieldPtr_);
1139 }
1140 else
1142 patchNeighbourFieldPtr_.reset(nullptr);
1143 }
1144}
1145
1146
1147template<class Type>
1149(
1150 const fvPatchField<Type>& ptf
1151)
1152{
1154
1155 //Pout<< "cyclicAMIFvPatchField::operator== :"
1156 // << " field:" << this->internalField().name()
1157 // << " patch:" << this->patch().name()
1158 // << " copying from field:" << ptf.internalField().name()
1159 // << endl;
1160
1161 const auto* cycPtr = isA<cyclicAMIFvPatchField<Type>>(ptf);
1162
1163 if
1164 (
1165 cycPtr
1166 && cycPtr->patchNeighbourFieldPtr_
1167 && cycPtr->patchNeighbourFieldPtr_->size() == this->size()
1168 )
1169 {
1170 if (!patchNeighbourFieldPtr_)
1171 {
1172 patchNeighbourFieldPtr_ = autoPtr<Field<Type>>::New();
1173 }
1174
1175 // Copy values
1176 *patchNeighbourFieldPtr_ = *(cycPtr->patchNeighbourFieldPtr_);
1177 }
1178 else
1179 {
1180 patchNeighbourFieldPtr_.reset(nullptr);
1181 }
1182}
1183
1184
1185// ************************************************************************* //
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
static const char *const typeName
Typename for Field.
Definition Field.H:93
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
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
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
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 List with indirect addressing. Like IndirectList but does not store addressing.
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
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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.
cyclicAMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const cyclicAMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
virtual void manipulateMatrix(fvMatrix< Type > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
virtual bool coupled() const
Return true if coupled. Note that the underlying 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.
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
virtual void write(Ostream &os) const
Write.
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const cyclicAMIFvPatch & cyclicAMIPatch() const noexcept
Return local reference cast into the cyclic AMI patch.
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 Arbitrary Mesh Interface (AMI).
virtual const cyclicAMIFvPatch & neighbPatch() const
Return a reference to the neighbour patch.
Abstract base class for cyclic AMI coupled interfaces.
void transformCoupleField(Field< Type > &f) const
Transform given patch field.
Cyclic patch for Arbitrary Mesh Interface (AMI).
tmp< Field< Type > > interpolate(const Field< Type > &fld, const UList< Type > &defaultValues=UList< Type >()) const
Interpolate field.
void initInterpolate(const Field< Type > &fld, labelRange &sendRequests, labelRange &recvRequests, PtrList< List< Type > > &sendBuffers, PtrList< List< Type > > &recvBuffers, labelRange &sendRequests1, labelRange &recvRequests1, PtrList< List< Type > > &sendBuffers1, PtrList< List< Type > > &recvBuffers1) const
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 autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
virtual void operator==(const fvPatchField< Type > &)
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
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 operator=(const UList< Type > &)
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
void extrapolateInternal()
Assign the patch field from the internal field.
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
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
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
dynamicFvMesh & mesh
#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)
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< labelList > labelListList
List of labelList.
Definition labelList.H:38
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")
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
dict add("bounds", meshBb)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299