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-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
35 template<class Type>
37 (
38  const fvPatch& p,
40 )
41 :
43  coupledFvPatchField<Type>(p, iF),
44  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p)),
45  patchNeighbourFieldPtr_(nullptr)
46 {}
47 
48 
49 template<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 {
62  if (!isA<cyclicACMIFvPatch>(p))
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()
70  << exit(FatalIOError);
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 GeometricField<Type, fvPatchField, volMesh>& fld =
95  static_cast<const GeometricField<Type, fvPatchField, volMesh>&>
96  (
97  this->primitiveField()
98  );
99  if (!fld.boundaryField().set(cyclicACMIPatch_.nonOverlapPatchID()))
100  {
102  << " patch " << p.name()
103  << " of field " << this->internalField().name()
104  << " refers to non-overlap patch "
105  << cyclicACMIPatch_.cyclicACMIPatch().nonOverlapPatchName()
106  << " which is not constructed yet." << nl
107  << " Either supply an initial value or change the ordering"
108  << " in the file"
109  << exit(FatalIOError);
110  }
111 
112  // Tricky: avoid call to evaluate without call to initEvaluate.
113  // For now just disable the localConsistency to make it use the
114  // old logic (ultimately calls the fully self contained
115  // patchNeighbourField)
116 
117  int& consistency =
118  GeometricField<Type, fvPatchField, volMesh>::
119  Boundary::localConsistency;
120 
121  const int oldConsistency = consistency;
122  consistency = 0;
123 
125 
126  consistency = oldConsistency;
127  }
128 }
129 
130 
131 template<class Type>
133 (
134  const cyclicACMIFvPatchField<Type>& ptf,
135  const fvPatch& p,
136  const DimensionedField<Type, volMesh>& iF,
137  const fvPatchFieldMapper& mapper
138 )
139 :
140  cyclicACMILduInterfaceField(),
141  coupledFvPatchField<Type>(ptf, p, iF, mapper),
142  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p)),
143  patchNeighbourFieldPtr_(nullptr)
144 {
145  //if (ptf.patchNeighbourFieldPtr_ && cacheNeighbourField())
146  //{
147  // patchNeighbourFieldPtr_.reset
148  // (
149  // new Field<Type>(ptf.patchNeighbourFieldPtr_(), mapper)
150  // );
151  //}
152 
153  if (!isA<cyclicACMIFvPatch>(this->patch()))
154  {
156  << "\n patch type '" << p.type()
157  << "' not constraint type '" << typeName << "'"
158  << "\n for patch " << p.name()
159  << " of field " << this->internalField().name()
160  << " in file " << this->internalField().objectPath()
161  << exit(FatalError);
162  }
163  if (debug && !ptf.all_ready())
164  {
166  << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
168  }
169 }
170 
171 
172 template<class Type>
174 (
175  const cyclicACMIFvPatchField<Type>& ptf
176 )
177 :
178  cyclicACMILduInterfaceField(),
179  coupledFvPatchField<Type>(ptf),
180  cyclicACMIPatch_(ptf.cyclicACMIPatch_),
181  patchNeighbourFieldPtr_(nullptr)
182 {
183  if (debug && !ptf.all_ready())
184  {
186  << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
188  }
189 }
190 
191 
192 template<class Type>
194 (
195  const cyclicACMIFvPatchField<Type>& ptf,
196  const DimensionedField<Type, volMesh>& iF
197 )
198 :
199  cyclicACMILduInterfaceField(),
200  coupledFvPatchField<Type>(ptf, iF),
201  cyclicACMIPatch_(ptf.cyclicACMIPatch_),
202  patchNeighbourFieldPtr_(nullptr)
203 {
204  if (debug && !ptf.all_ready())
205  {
207  << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
208  << abort(FatalError);
209  }
210 }
211 
212 
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214 
215 template<class Type>
217 {
218  return cyclicACMIPatch_.coupled();
219 }
220 
221 
222 template<class Type>
224 {
225  int done = 0;
226 
227  if
228  (
230  (
231  recvRequests_.start(),
232  recvRequests_.size()
233  )
234  )
235  {
236  recvRequests_.clear();
237  ++done;
238  }
239 
240  if
241  (
243  (
244  sendRequests_.start(),
245  sendRequests_.size()
246  )
247  )
248  {
249  sendRequests_.clear();
250  ++done;
251  }
252 
253  return (done == 2);
254 }
255 
256 
257 template<class Type>
259 {
260  if
261  (
263  (
264  recvRequests_.start(),
265  recvRequests_.size()
266  )
267  )
268  {
269  recvRequests_.clear();
270 
271  if
272  (
274  (
275  sendRequests_.start(),
276  sendRequests_.size()
277  )
278  )
279  {
280  sendRequests_.clear();
281  }
282 
283  return true;
284  }
285 
286  return false;
287 }
288 
289 
290 template<class Type>
293 (
294  const Field<Type>& iField
295 ) const
296 {
297  DebugPout
298  << "cyclicACMIFvPatchField::patchNeighbourField(const Field<Type>&) :"
299  << " field:" << this->internalField().name()
300  << " patch:" << this->patch().name()
301  << endl;
302 
303  // By pass polyPatch to get nbrId. Instead use cyclicACMIFvPatch virtual
304  // neighbPatch()
305  const cyclicACMIFvPatch& neighbPatch = cyclicACMIPatch_.neighbPatch();
306  const labelUList& nbrFaceCells = neighbPatch.faceCells();
307 
308  tmp<Field<Type>> tpnf
309  (
310  cyclicACMIPatch_.interpolate
311  (
312  Field<Type>
313  (
314  iField,
315  nbrFaceCells
316  //cpp.neighbPatch().faceCells()
317  )
318  )
319  );
320 
321  if (doTransform())
322  {
323  transform(tpnf.ref(), forwardT(), tpnf());
324  }
325 
326  return tpnf;
327 }
328 
329 
330 template<class Type>
332 {
333  /*
334  return
335  (
336  GeometricField<Type, fvPatchField, volMesh>::Boundary::localConsistency
337  != 0
338  );
339  */
340  return false;
341 }
342 
343 
344 template<class Type>
347 {
348  if (this->ownerAMI().distributed() && cacheNeighbourField())
349  {
350  if (!this->ready())
351  {
353  << "Outstanding recv request(s) on patch "
354  << cyclicACMIPatch_.name()
355  << " field " << this->internalField().name()
356  << abort(FatalError);
357  }
358 
359  // Initialise if not done in construct-from-dictionary
360  if (!patchNeighbourFieldPtr_)
361  {
362  DebugPout
363  << "cyclicACMIFvPatchField::patchNeighbourField() :"
364  << " field:" << this->internalField().name()
365  << " patch:" << this->patch().name()
366  << " calculating&caching patchNeighbourField"
367  << endl;
368 
369  // Do interpolation and store result
370  patchNeighbourFieldPtr_.reset
371  (
372  patchNeighbourField(this->primitiveField()).ptr()
373  );
374  }
375  else
376  {
377  DebugPout
378  << "cyclicACMIFvPatchField::patchNeighbourField() :"
379  << " field:" << this->internalField().name()
380  << " patch:" << this->patch().name()
381  << " returning cached patchNeighbourField"
382  << endl;
383  }
384  return patchNeighbourFieldPtr_();
385  }
386  else
387  {
388  // Do interpolation
389  DebugPout
390  << "cyclicACMIFvPatchField::evaluate() :"
391  << " field:" << this->internalField().name()
392  << " patch:" << this->patch().name()
393  << " calculating up-to-date patchNeighbourField"
394  << endl;
395 
396  return patchNeighbourField(this->primitiveField());
397  }
398 }
399 
400 
401 template<class Type>
404 {
405  const GeometricField<Type, fvPatchField, volMesh>& fld =
406  static_cast<const GeometricField<Type, fvPatchField, volMesh>&>
407  (
408  this->primitiveField()
409  );
410 
411  return refCast<const cyclicACMIFvPatchField<Type>>
412  (
413  fld.boundaryField()[cyclicACMIPatch_.neighbPatchID()]
414  );
415 }
416 
417 
418 template<class Type>
421 {
424  (
425  this->primitiveField()
426  );
427 
428  // WIP: Needs to re-direct nonOverlapPatchID to new patchId for assembly?
429  return fld.boundaryField()[cyclicACMIPatch_.nonOverlapPatchID()];
430 }
431 
432 
433 template<class Type>
435 (
436  const Pstream::commsTypes commsType
437 )
438 {
439  if (!this->updated())
440  {
441  this->updateCoeffs();
442  }
443 
444  if (this->ownerAMI().distributed() && cacheNeighbourField())
445  {
446  if (commsType != UPstream::commsTypes::nonBlocking)
447  {
448  // Invalidate old field - or flag as fatal?
449  patchNeighbourFieldPtr_.reset(nullptr);
450  return;
451  }
452 
453  // Start sending
454  DebugPout
455  << "cyclicACMIFvPatchField::initEvaluate() :"
456  << " field:" << this->internalField().name()
457  << " patch:" << this->patch().name()
458  << " starting send&receive"
459  << endl;
460 
461  // Bypass polyPatch to get nbrId.
462  // - use cyclicACMIFvPatch::neighbPatch() virtual instead
463  const cyclicACMIFvPatch& neighbPatch = cyclicACMIPatch_.neighbPatch();
464  const labelUList& nbrFaceCells = neighbPatch.faceCells();
465  const Field<Type> pnf(this->primitiveField(), nbrFaceCells);
466 
467  // Assert that all receives are known to have finished
468  if (!recvRequests_.empty())
469  {
471  << "Outstanding recv request(s) on patch "
472  << cyclicACMIPatch_.name()
473  << " field " << this->internalField().name()
474  << abort(FatalError);
475  }
476 
477  // Assume that sends are also OK
478  sendRequests_.clear();
479 
480  cyclicACMIPatch_.initInterpolate
481  (
482  pnf,
483  sendRequests_,
484  sendBufs_,
485  recvRequests_,
486  recvBufs_
487  );
488  }
489 }
490 
491 
492 template<class Type>
494 (
495  const Pstream::commsTypes commsType
496 )
497 {
498  if (!this->updated())
499  {
500  this->updateCoeffs();
501  }
502 
503  const auto& AMI = this->ownerAMI();
504 
505  if (AMI.distributed() && cacheNeighbourField())
506  {
507  // Calculate patchNeighbourField
508  if (commsType != UPstream::commsTypes::nonBlocking)
509  {
511  << "Can only evaluate distributed AMI with nonBlocking"
512  << exit(FatalError);
513  }
514 
515  patchNeighbourFieldPtr_.reset(nullptr);
516 
517  if (!this->ready())
518  {
520  << "Outstanding recv request(s) on patch "
521  << cyclicACMIPatch_.name()
522  << " field " << this->internalField().name()
523  << abort(FatalError);
524  }
525 
526  patchNeighbourFieldPtr_.reset
527  (
528  cyclicACMIPatch_.interpolate
529  (
530  Field<Type>::null(), // Not used for distributed
531  recvRequests_,
532  recvBufs_
533  ).ptr()
534  );
535 
536  // Receive requests all handled by last function call
537  recvRequests_.clear();
538 
539 
540  auto& patchNeighbourField = patchNeighbourFieldPtr_.ref();
541 
542  if (doTransform())
543  {
544  // In-place transform
545  transform(patchNeighbourField, forwardT(), patchNeighbourField);
546  }
547  }
548 
549  // Use patchNeighbourField() and patchInternalField() to obtain face value
551 }
552 
553 
554 template<class Type>
556 (
557  solveScalarField& result,
558  const bool add,
559  const lduAddressing& lduAddr,
560  const label patchId,
561  const solveScalarField& psiInternal,
562  const scalarField& coeffs,
563  const direction cmpt,
564  const Pstream::commsTypes commsType
565 ) const
566 {
567  if (this->ownerAMI().distributed())
568  {
569  // Start sending
570  if (commsType != UPstream::commsTypes::nonBlocking)
571  {
573  << "Can only evaluate distributed AMI with nonBlocking"
574  << exit(FatalError);
575  }
576 
577  DebugPout
578  << "cyclicACMIFvPatchField::initInterfaceMatrixUpdate() :"
579  << " field:" << this->internalField().name()
580  << " patch:" << this->patch().name()
581  << " starting send&receive"
582  << endl;
583 
584  const labelUList& nbrFaceCells =
585  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
586 
587  solveScalarField pnf(psiInternal, nbrFaceCells);
588 
589  // Transform according to the transformation tensors
590  transformCoupleField(pnf, cmpt);
591 
592  // Assert that all receives are known to have finished
593  if (!recvRequests_.empty())
594  {
596  << "Outstanding recv request(s) on patch "
597  << cyclicACMIPatch_.name()
598  << " field " << this->internalField().name()
599  << abort(FatalError);
600  }
601 
602  // Assume that sends are also OK
603  sendRequests_.clear();
604 
605  cyclicACMIPatch_.initInterpolate
606  (
607  pnf,
608  sendRequests_,
609  scalarSendBufs_,
610  recvRequests_,
611  scalarRecvBufs_
612  );
613  }
615  this->updatedMatrix(false);
616 }
617 
618 
619 template<class Type>
621 (
622  solveScalarField& result,
623  const bool add,
624  const lduAddressing& lduAddr,
625  const label patchId,
626  const solveScalarField& psiInternal,
627  const scalarField& coeffs,
628  const direction cmpt,
629  const Pstream::commsTypes commsType
630 ) const
631 {
632  DebugPout<< "cyclicACMIFvPatchField::updateInterfaceMatrix() :"
633  << " field:" << this->internalField().name()
634  << " patch:" << this->patch().name()
635  << endl;
636 
637  // note: only applying coupled contribution
638 
639  const labelUList& faceCells = lduAddr.patchAddr(patchId);
640 
641  solveScalarField pnf;
642 
643  if (this->ownerAMI().distributed())
644  {
645  if (commsType != UPstream::commsTypes::nonBlocking)
646  {
648  << "Can only evaluate distributed AMI with nonBlocking"
649  << exit(FatalError);
650  }
651 
652  DebugPout
653  << "cyclicACMIFvPatchField::evaluate() :"
654  << " field:" << this->internalField().name()
655  << " patch:" << this->patch().name()
656  << " consuming received coupled neighbourfield"
657  << endl;
658 
659  pnf =
660  cyclicACMIPatch_.interpolate
661  (
662  solveScalarField::null(), // Not used for distributed
663  recvRequests_,
664  scalarRecvBufs_
665  );
666 
667  // Receive requests all handled by last function call
668  recvRequests_.clear();
669  }
670  else
671  {
672  const labelUList& nbrFaceCells =
673  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
674 
675  pnf = solveScalarField(psiInternal, nbrFaceCells);
676 
677  // Transform according to the transformation tensors
678  transformCoupleField(pnf, cmpt);
679 
680  pnf = cyclicACMIPatch_.interpolate(pnf);
681  }
682 
683  this->addToInternalField(result, !add, faceCells, coeffs, pnf);
685  this->updatedMatrix(true);
686 }
687 
688 
689 template<class Type>
691 (
692  Field<Type>& result,
693  const bool add,
694  const lduAddressing& lduAddr,
695  const label patchId,
696  const Field<Type>& psiInternal,
697  const scalarField& coeffs,
698  const Pstream::commsTypes commsType
699 ) const
700 {
701  const auto& AMI = this->ownerAMI();
702 
703  if (AMI.distributed())
704  {
705  if (commsType != UPstream::commsTypes::nonBlocking)
706  {
708  << "Can only evaluate distributed AMI with nonBlocking"
709  << exit(FatalError);
710  }
711 
712  const labelUList& nbrFaceCells =
713  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
714 
715  Field<Type> pnf(psiInternal, nbrFaceCells);
716 
717  // Transform according to the transformation tensors
718  transformCoupleField(pnf);
719 
720  // Assert that all receives are known to have finished
721  if (!recvRequests_.empty())
722  {
724  << "Outstanding recv request(s) on patch "
725  << cyclicACMIPatch_.name()
726  << " field " << this->internalField().name()
727  << abort(FatalError);
728  }
729 
730  // Assume that sends are also OK
731  sendRequests_.clear();
732 
733  cyclicACMIPatch_.initInterpolate
734  (
735  pnf,
736  sendRequests_,
737  sendBufs_,
738  recvRequests_,
739  recvBufs_
740  );
741  }
743  this->updatedMatrix(false);
744 }
745 
746 
747 template<class Type>
749 (
750  Field<Type>& result,
751  const bool add,
752  const lduAddressing& lduAddr,
753  const label patchId,
754  const Field<Type>& psiInternal,
755  const scalarField& coeffs,
756  const Pstream::commsTypes commsType
757 ) const
758 {
759  // note: only applying coupled contribution
760 
761  const labelUList& faceCells = lduAddr.patchAddr(patchId);
762 
763  const auto& AMI = this->ownerAMI();
764 
765  Field<Type> pnf;
766 
767  if (AMI.distributed())
768  {
769  if (commsType != UPstream::commsTypes::nonBlocking)
770  {
772  << "Can only evaluate distributed AMI with nonBlocking"
773  << exit(FatalError);
774  }
775 
776  pnf =
777  cyclicACMIPatch_.interpolate
778  (
779  Field<Type>::null(), // Not used for distributed
780  recvRequests_,
781  recvBufs_
782  );
783 
784  // Receive requests all handled by last function call
785  recvRequests_.clear();
786  }
787  else
788  {
789  const labelUList& nbrFaceCells =
790  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
791 
792  pnf = Field<Type>(psiInternal, nbrFaceCells);
793 
794  // Transform according to the transformation tensors
795  transformCoupleField(pnf);
796 
797  pnf = cyclicACMIPatch_.interpolate(pnf);
798  }
799 
800  this->addToInternalField(result, !add, faceCells, coeffs, pnf);
802  this->updatedMatrix(true);
803 }
804 
805 
806 template<class Type>
808 (
809  fvMatrix<Type>& matrix
810 )
811 {
812  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
813 
814  // Nothing to be done by the AMI, but re-direct to non-overlap patch
815  // with non-overlap patch weights
816  const fvPatchField<Type>& npf = nonOverlapPatchField();
818  const_cast<fvPatchField<Type>&>(npf).manipulateMatrix(matrix, 1.0 - mask);
819 }
820 
821 
822 template<class Type>
824 (
825  fvMatrix<Type>& matrix,
826  const label mat,
827  const direction cmpt
828 )
829 {
830  if (this->cyclicACMIPatch().owner())
831  {
832  label index = this->patch().index();
833 
834  const label globalPatchID =
835  matrix.lduMeshAssembly().patchLocalToGlobalMap()[mat][index];
836 
837  const Field<scalar> intCoeffsCmpt
838  (
839  matrix.internalCoeffs()[globalPatchID].component(cmpt)
840  );
841 
842  const Field<scalar> boundCoeffsCmpt
843  (
844  matrix.boundaryCoeffs()[globalPatchID].component(cmpt)
845  );
846 
847  tmp<Field<scalar>> tintCoeffs(coeffs(matrix, intCoeffsCmpt, mat));
848  tmp<Field<scalar>> tbndCoeffs(coeffs(matrix, boundCoeffsCmpt, mat));
849  const Field<scalar>& intCoeffs = tintCoeffs.ref();
850  const Field<scalar>& bndCoeffs = tbndCoeffs.ref();
851 
852  const labelUList& u = matrix.lduAddr().upperAddr();
853  const labelUList& l = matrix.lduAddr().lowerAddr();
854 
855  label subFaceI = 0;
856 
857  const labelList& faceMap =
858  matrix.lduMeshAssembly().faceBoundMap()[mat][index];
859 
860  forAll (faceMap, j)
861  {
862  label globalFaceI = faceMap[j];
863 
864  const scalar boundCorr = -bndCoeffs[subFaceI];
865  const scalar intCorr = -intCoeffs[subFaceI];
866 
867  matrix.upper()[globalFaceI] += boundCorr;
868  matrix.diag()[u[globalFaceI]] -= intCorr;
869  matrix.diag()[l[globalFaceI]] -= boundCorr;
870 
871  if (matrix.asymmetric())
872  {
873  matrix.lower()[globalFaceI] += intCorr;
874  }
875  subFaceI++;
876  }
877 
878  // Set internalCoeffs and boundaryCoeffs in the assembly matrix
879  // on cyclicACMI patches to be used in the individual matrix by
880  // matrix.flux()
881  if (matrix.psi(mat).mesh().fluxRequired(this->internalField().name()))
882  {
883  matrix.internalCoeffs().set
884  (
885  globalPatchID, intCoeffs*pTraits<Type>::one
886  );
887  matrix.boundaryCoeffs().set
888  (
889  globalPatchID, bndCoeffs*pTraits<Type>::one
890  );
891 
892  const label nbrPathID =
893  cyclicACMIPatch_.cyclicACMIPatch().neighbPatchID();
894 
895  const label nbrGlobalPatchID =
897  [mat][nbrPathID];
898 
899  matrix.internalCoeffs().set
900  (
901  nbrGlobalPatchID, intCoeffs*pTraits<Type>::one
902  );
903  matrix.boundaryCoeffs().set
904  (
905  nbrGlobalPatchID, bndCoeffs*pTraits<Type>::one
906  );
907  }
908  }
909 }
910 
911 
912 template<class Type>
915 (
916  fvMatrix<Type>& matrix,
917  const Field<scalar>& coeffs,
918  const label mat
919 ) const
920 {
921  const label index(this->patch().index());
922 
923  const label nSubFaces
924  (
925  matrix.lduMeshAssembly().cellBoundMap()[mat][index].size()
926  );
927 
928  auto tmapCoeffs = tmp<Field<scalar>>::New(nSubFaces, Zero);
929  auto& mapCoeffs = tmapCoeffs.ref();
930 
931  const scalarListList& srcWeight =
932  cyclicACMIPatch_.cyclicACMIPatch().AMI().srcWeights();
933 
934  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
935 
936  const scalar tol = cyclicACMIPolyPatch::tolerance();
937  label subFaceI = 0;
938  forAll(mask, faceI)
939  {
940  const scalarList& w = srcWeight[faceI];
941  for(label i=0; i<w.size(); i++)
942  {
943  if (mask[faceI] > tol)
944  {
945  const label localFaceId =
946  matrix.lduMeshAssembly().facePatchFaceMap()
947  [mat][index][subFaceI];
948  mapCoeffs[subFaceI] = w[i]*coeffs[localFaceId];
949  }
950  subFaceI++;
951  }
952  }
953 
954  return tmapCoeffs;
955 }
956 
957 
958 template<class Type>
960 {
961  // Update non-overlap patch - some will implement updateCoeffs, and
962  // others will implement evaluate
963 
964  // Pass in (1 - mask) to give non-overlap patch the chance to do
965  // manipulation of non-face based data
966 
967  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
968  const fvPatchField<Type>& npf = nonOverlapPatchField();
969  const_cast<fvPatchField<Type>&>(npf).updateWeightedCoeffs(1.0 - mask);
970 }
971 
972 
973 template<class Type>
975 {
978 
979  if (patchNeighbourFieldPtr_)
980  {
981  patchNeighbourFieldPtr_->writeEntry("neighbourValue", os);
982  }
983 }
984 
985 
986 // ************************************************************************* //
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
label patchId(-1)
const scalarField & diag() const
Definition: lduMatrix.C:195
dictionary dict
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
uint8_t direction
Definition: direction.H:46
const labelListList & patchLocalToGlobalMap() const
Return patchLocalToGlobalMap.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
Field< solveScalar > solveScalarField
commsTypes
Communications types.
Definition: UPstream.H:77
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
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 tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
const lduPrimitiveMeshAssembly & lduMeshAssembly()
Return optional lduAdressing.
Definition: fvMatrix.H:476
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition: typeInfo.H:172
virtual label nonOverlapPatchID() const
Return neighbour.
static const Field< Type > & null() noexcept
Return a null Field (reference to a nullObject). Behaves like an empty Field.
Definition: Field.H:137
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Generic GeometricField class.
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:403
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static const char *const typeName
Typename for Field.
Definition: Field.H:86
Spatial transformation functions for primitive fields.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual bool ready() const
Are all (receive) data available?
virtual void write(Ostream &os) const
Write.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:130
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const cyclicACMIFvPatch & neighbPatch() const
Return neighbour fvPatch.
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
const scalarField & lower() const
Definition: lduMatrix.C:306
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
Definition: fvPatchField.H:704
String literal.
Definition: keyType.H:82
const cyclicACMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
cyclicACMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Abstract base class for coupled patches.
Abstract base class for cyclic ACMI coupled interfaces.
const word & nonOverlapPatchName() const
Non-overlapping patch name.
bool asymmetric() const noexcept
Matrix is asymmetric (ie, full)
Definition: lduMatrix.H:846
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
int debug
Static debugging option.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
OBJstream os(runTime.globalPath()/outputName)
const labelListListList & faceBoundMap() const
Return boundary face map.
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.
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
const cyclicACMIPolyPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic patch.
const scalarField & upper() const
Definition: lduMatrix.C:235
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:765
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:637
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:696
const std::string patch
OpenFOAM patch number as a std::string.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
#define DebugPout
Report an information message using Foam::Pout.
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition: fvMatrix.H:547
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition: fvMatrix.H:565
virtual const word & name() const
Return name.
Definition: fvPatch.H:210
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
"buffered" : (MPI_Bsend, MPI_Recv)
static scalar tolerance()
Overlap tolerance.
static bool finishedRequests(const label pos, label len=-1)
Non-blocking comms: have all requests (from position onwards) finished? Corresponds to MPI_Testall() ...
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127