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-2023 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 
113  }
114 }
115 
116 
117 template<class Type>
119 (
120  const cyclicACMIFvPatchField<Type>& ptf,
121  const fvPatch& p,
122  const DimensionedField<Type, volMesh>& iF,
123  const fvPatchFieldMapper& mapper
124 )
125 :
126  cyclicACMILduInterfaceField(),
127  coupledFvPatchField<Type>(ptf, p, iF, mapper),
128  cyclicACMIPatch_(refCast<const cyclicACMIFvPatch>(p)),
129  patchNeighbourFieldPtr_(nullptr)
130 {
131  //if (ptf.patchNeighbourFieldPtr_ && cacheNeighbourField())
132  //{
133  // patchNeighbourFieldPtr_.reset
134  // (
135  // new Field<Type>(ptf.patchNeighbourFieldPtr_(), mapper)
136  // );
137  //}
138 
139  if (!isA<cyclicACMIFvPatch>(this->patch()))
140  {
142  << "\n patch type '" << p.type()
143  << "' not constraint type '" << typeName << "'"
144  << "\n for patch " << p.name()
145  << " of field " << this->internalField().name()
146  << " in file " << this->internalField().objectPath()
147  << exit(FatalError);
148  }
149  if (debug && !ptf.all_ready())
150  {
152  << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
154  }
155 }
156 
157 
158 template<class Type>
160 (
161  const cyclicACMIFvPatchField<Type>& ptf
162 )
163 :
164  cyclicACMILduInterfaceField(),
165  coupledFvPatchField<Type>(ptf),
166  cyclicACMIPatch_(ptf.cyclicACMIPatch_),
167  patchNeighbourFieldPtr_(nullptr)
168 {
169  if (debug && !ptf.all_ready())
170  {
172  << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
174  }
175 }
176 
177 
178 template<class Type>
180 (
181  const cyclicACMIFvPatchField<Type>& ptf,
182  const DimensionedField<Type, volMesh>& iF
183 )
184 :
185  cyclicACMILduInterfaceField(),
186  coupledFvPatchField<Type>(ptf, iF),
187  cyclicACMIPatch_(ptf.cyclicACMIPatch_),
188  patchNeighbourFieldPtr_(nullptr)
189 {
190  if (debug && !ptf.all_ready())
191  {
193  << "Outstanding request(s) on patch " << cyclicACMIPatch_.name()
194  << abort(FatalError);
195  }
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
201 template<class Type>
203 {
204  return cyclicACMIPatch_.coupled();
205 }
206 
207 
208 template<class Type>
210 {
211  int done = 0;
212 
213  if
214  (
216  (
217  recvRequests_.start(),
218  recvRequests_.size()
219  )
220  )
221  {
222  recvRequests_.clear();
223  ++done;
224  }
225 
226  if
227  (
229  (
230  sendRequests_.start(),
231  sendRequests_.size()
232  )
233  )
234  {
235  sendRequests_.clear();
236  ++done;
237  }
238 
239  return (done == 2);
240 }
241 
242 
243 template<class Type>
245 {
246  if
247  (
249  (
250  recvRequests_.start(),
251  recvRequests_.size()
252  )
253  )
254  {
255  recvRequests_.clear();
256 
257  if
258  (
260  (
261  sendRequests_.start(),
262  sendRequests_.size()
263  )
264  )
265  {
266  sendRequests_.clear();
267  }
268 
269  return true;
270  }
271 
272  return false;
273 }
274 
275 
276 template<class Type>
279 (
280  const Field<Type>& iField
281 ) const
282 {
283  DebugPout
284  << "cyclicACMIFvPatchField::patchNeighbourField(const Field<Type>&) :"
285  << " field:" << this->internalField().name()
286  << " patch:" << this->patch().name()
287  << endl;
288 
289  // By pass polyPatch to get nbrId. Instead use cyclicACMIFvPatch virtual
290  // neighbPatch()
291  const cyclicACMIFvPatch& neighbPatch = cyclicACMIPatch_.neighbPatch();
292  const labelUList& nbrFaceCells = neighbPatch.faceCells();
293 
294  tmp<Field<Type>> tpnf
295  (
296  cyclicACMIPatch_.interpolate
297  (
298  Field<Type>
299  (
300  iField,
301  nbrFaceCells
302  //cpp.neighbPatch().faceCells()
303  )
304  )
305  );
306 
307  if (doTransform())
308  {
309  transform(tpnf.ref(), forwardT(), tpnf());
310  }
311 
312  return tpnf;
313 }
314 
315 
316 template<class Type>
318 {
319  /*
320  return
321  (
322  GeometricField<Type, fvPatchField, volMesh>::Boundary::localConsistency
323  != 0
324  );
325  */
326  return false;
327 }
328 
329 
330 template<class Type>
333 {
334  if (this->ownerAMI().distributed() && cacheNeighbourField())
335  {
336  if (!this->ready())
337  {
339  << "Outstanding recv request(s) on patch "
340  << cyclicACMIPatch_.name()
341  << " field " << this->internalField().name()
342  << abort(FatalError);
343  }
344 
345  // Initialise if not done in construct-from-dictionary
346  if (!patchNeighbourFieldPtr_)
347  {
348  DebugPout
349  << "cyclicACMIFvPatchField::patchNeighbourField() :"
350  << " field:" << this->internalField().name()
351  << " patch:" << this->patch().name()
352  << " calculating&caching patchNeighbourField"
353  << endl;
354 
355  // Do interpolation and store result
356  patchNeighbourFieldPtr_.reset
357  (
358  patchNeighbourField(this->primitiveField()).ptr()
359  );
360  }
361  else
362  {
363  DebugPout
364  << "cyclicACMIFvPatchField::patchNeighbourField() :"
365  << " field:" << this->internalField().name()
366  << " patch:" << this->patch().name()
367  << " returning cached patchNeighbourField"
368  << endl;
369  }
370  return patchNeighbourFieldPtr_();
371  }
372  else
373  {
374  // Do interpolation
375  DebugPout
376  << "cyclicACMIFvPatchField::evaluate() :"
377  << " field:" << this->internalField().name()
378  << " patch:" << this->patch().name()
379  << " calculating up-to-date patchNeighbourField"
380  << endl;
381 
382  return patchNeighbourField(this->primitiveField());
383  }
384 }
385 
386 
387 template<class Type>
390 {
391  const GeometricField<Type, fvPatchField, volMesh>& fld =
392  static_cast<const GeometricField<Type, fvPatchField, volMesh>&>
393  (
394  this->primitiveField()
395  );
396 
397  return refCast<const cyclicACMIFvPatchField<Type>>
398  (
399  fld.boundaryField()[cyclicACMIPatch_.neighbPatchID()]
400  );
401 }
402 
403 
404 template<class Type>
407 {
410  (
411  this->primitiveField()
412  );
413 
414  // WIP: Needs to re-direct nonOverlapPatchID to new patchId for assembly?
415  return fld.boundaryField()[cyclicACMIPatch_.nonOverlapPatchID()];
416 }
417 
418 
419 template<class Type>
421 (
422  const Pstream::commsTypes commsType
423 )
424 {
425  if (!this->updated())
426  {
427  this->updateCoeffs();
428  }
429 
430  if (this->ownerAMI().distributed() && cacheNeighbourField())
431  {
432  if (commsType != UPstream::commsTypes::nonBlocking)
433  {
434  // Invalidate old field - or flag as fatal?
435  patchNeighbourFieldPtr_.reset(nullptr);
436  return;
437  }
438 
439  // Start sending
440  DebugPout
441  << "cyclicACMIFvPatchField::initEvaluate() :"
442  << " field:" << this->internalField().name()
443  << " patch:" << this->patch().name()
444  << " starting send&receive"
445  << endl;
446 
447  if (!this->ready())
448  {
450  << "Outstanding recv request(s) on patch "
451  << cyclicACMIPatch_.name()
452  << " field " << this->internalField().name()
453  << abort(FatalError);
454  }
455 
456  // By-pass polyPatch to get nbrId. Instead use cyclicACMIFvPatch virtual
457  // neighbPatch()
458  const cyclicACMIFvPatch& neighbPatch = cyclicACMIPatch_.neighbPatch();
459  const labelUList& nbrFaceCells = neighbPatch.faceCells();
460  const Field<Type> pnf(this->primitiveField(), nbrFaceCells);
461 
462  cyclicACMIPatch_.initInterpolate
463  (
464  pnf,
465  sendRequests_,
466  sendBufs_,
467  recvRequests_,
468  recvBufs_
469  );
470  }
471 }
472 
473 
474 template<class Type>
476 (
477  const Pstream::commsTypes commsType
478 )
479 {
480  if (!this->updated())
481  {
482  this->updateCoeffs();
483  }
484 
485  const auto& AMI = this->ownerAMI();
486 
487  if (AMI.distributed() && cacheNeighbourField())
488  {
489  // Calculate patchNeighbourField
490  if (commsType != UPstream::commsTypes::nonBlocking)
491  {
493  << "Can only evaluate distributed AMI with nonBlocking"
494  << exit(FatalError);
495  }
496 
497  patchNeighbourFieldPtr_.reset(nullptr);
498 
499  if (!this->ready())
500  {
502  << "Outstanding recv request(s) on patch "
503  << cyclicACMIPatch_.name()
504  << " field " << this->internalField().name()
505  << abort(FatalError);
506  }
507 
508  patchNeighbourFieldPtr_.reset
509  (
510  cyclicACMIPatch_.interpolate
511  (
512  Field<Type>::null(), // Not used for distributed
513  recvRequests_,
514  recvBufs_
515  ).ptr()
516  );
517 
518  auto& patchNeighbourField = patchNeighbourFieldPtr_.ref();
519 
520  if (doTransform())
521  {
522  // In-place transform
523  transform(patchNeighbourField, forwardT(), patchNeighbourField);
524  }
525  }
526 
527  // Use patchNeighbourField() and patchInternalField() to obtain face value
529 }
530 
531 
532 template<class Type>
534 (
535  solveScalarField& result,
536  const bool add,
537  const lduAddressing& lduAddr,
538  const label patchId,
539  const solveScalarField& psiInternal,
540  const scalarField& coeffs,
541  const direction cmpt,
542  const Pstream::commsTypes commsType
543 ) const
544 {
545  if (this->ownerAMI().distributed())
546  {
547  // Start sending
548  if (commsType != UPstream::commsTypes::nonBlocking)
549  {
551  << "Can only evaluate distributed AMI with nonBlocking"
552  << exit(FatalError);
553  }
554 
555  DebugPout
556  << "cyclicACMIFvPatchField::initInterfaceMatrixUpdate() :"
557  << " field:" << this->internalField().name()
558  << " patch:" << this->patch().name()
559  << " starting send&receive"
560  << endl;
561 
562  if (!this->ready())
563  {
565  << "Outstanding recv request(s) on patch "
566  << cyclicACMIPatch_.name()
567  << " field " << this->internalField().name()
568  << abort(FatalError);
569  }
570 
571  const labelUList& nbrFaceCells =
572  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
573 
574  solveScalarField pnf(psiInternal, nbrFaceCells);
575 
576  // Transform according to the transformation tensors
577  transformCoupleField(pnf, cmpt);
578 
579  cyclicACMIPatch_.initInterpolate
580  (
581  pnf,
582  sendRequests_,
583  scalarSendBufs_,
584  recvRequests_,
585  scalarRecvBufs_
586  );
587  }
588 }
589 
590 
591 template<class Type>
593 (
594  solveScalarField& result,
595  const bool add,
596  const lduAddressing& lduAddr,
597  const label patchId,
598  const solveScalarField& psiInternal,
599  const scalarField& coeffs,
600  const direction cmpt,
601  const Pstream::commsTypes commsType
602 ) const
603 {
604  DebugPout<< "cyclicACMIFvPatchField::updateInterfaceMatrix() :"
605  << " field:" << this->internalField().name()
606  << " patch:" << this->patch().name()
607  << endl;
608 
609  // note: only applying coupled contribution
610 
611  const labelUList& faceCells = lduAddr.patchAddr(patchId);
612 
613  solveScalarField pnf;
614 
615  if (this->ownerAMI().distributed())
616  {
617  if (commsType != UPstream::commsTypes::nonBlocking)
618  {
620  << "Can only evaluate distributed AMI with nonBlocking"
621  << exit(FatalError);
622  }
623 
624  DebugPout
625  << "cyclicACMIFvPatchField::evaluate() :"
626  << " field:" << this->internalField().name()
627  << " patch:" << this->patch().name()
628  << " consuming received coupled neighbourfield"
629  << endl;
630 
631  pnf =
632  cyclicACMIPatch_.interpolate
633  (
634  solveScalarField::null(), // Not used for distributed
635  recvRequests_,
636  scalarRecvBufs_
637  );
638  }
639  else
640  {
641  const labelUList& nbrFaceCells =
642  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
643 
644  pnf = solveScalarField(psiInternal, nbrFaceCells);
645 
646  // Transform according to the transformation tensors
647  transformCoupleField(pnf, cmpt);
648 
649  pnf = cyclicACMIPatch_.interpolate(pnf);
650  }
652  this->addToInternalField(result, !add, faceCells, coeffs, pnf);
653 }
654 
655 
656 template<class Type>
658 (
659  Field<Type>& result,
660  const bool add,
661  const lduAddressing& lduAddr,
662  const label patchId,
663  const Field<Type>& psiInternal,
664  const scalarField& coeffs,
665  const Pstream::commsTypes commsType
666 ) const
667 {
668  const auto& AMI = this->ownerAMI();
669 
670  if (AMI.distributed())
671  {
672  if (commsType != UPstream::commsTypes::nonBlocking)
673  {
675  << "Can only evaluate distributed AMI with nonBlocking"
676  << exit(FatalError);
677  }
678 
679  if (!this->ready())
680  {
682  << "Outstanding recv request(s) on patch "
683  << cyclicACMIPatch_.name()
684  << " field " << this->internalField().name()
685  << abort(FatalError);
686  }
687 
688  const labelUList& nbrFaceCells =
689  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
690 
691  Field<Type> pnf(psiInternal, nbrFaceCells);
692 
693  // Transform according to the transformation tensors
694  transformCoupleField(pnf);
695 
696  cyclicACMIPatch_.initInterpolate
697  (
698  pnf,
699  sendRequests_,
700  sendBufs_,
701  recvRequests_,
702  recvBufs_
703  );
704  }
705 }
706 
707 
708 template<class Type>
710 (
711  Field<Type>& result,
712  const bool add,
713  const lduAddressing& lduAddr,
714  const label patchId,
715  const Field<Type>& psiInternal,
716  const scalarField& coeffs,
717  const Pstream::commsTypes commsType
718 ) const
719 {
720  // note: only applying coupled contribution
721 
722  const labelUList& faceCells = lduAddr.patchAddr(patchId);
723 
724  const auto& AMI = this->ownerAMI();
725 
726  Field<Type> pnf;
727 
728  if (AMI.distributed())
729  {
730  if (commsType != UPstream::commsTypes::nonBlocking)
731  {
733  << "Can only evaluate distributed AMI with nonBlocking"
734  << exit(FatalError);
735  }
736 
737  pnf =
738  cyclicACMIPatch_.interpolate
739  (
740  Field<Type>::null(), // Not used for distributed
741  recvRequests_,
742  recvBufs_
743  );
744  }
745  else
746  {
747  const labelUList& nbrFaceCells =
748  lduAddr.patchAddr(cyclicACMIPatch_.neighbPatchID());
749 
750  pnf = Field<Type>(psiInternal, nbrFaceCells);
751 
752  // Transform according to the transformation tensors
753  transformCoupleField(pnf);
754 
755  pnf = cyclicACMIPatch_.interpolate(pnf);
756  }
758  this->addToInternalField(result, !add, faceCells, coeffs, pnf);
759 }
760 
761 
762 template<class Type>
764 (
765  fvMatrix<Type>& matrix
766 )
767 {
768  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
769 
770  // Nothing to be done by the AMI, but re-direct to non-overlap patch
771  // with non-overlap patch weights
772  const fvPatchField<Type>& npf = nonOverlapPatchField();
774  const_cast<fvPatchField<Type>&>(npf).manipulateMatrix(matrix, 1.0 - mask);
775 }
776 
777 
778 template<class Type>
780 (
781  fvMatrix<Type>& matrix,
782  const label mat,
783  const direction cmpt
784 )
785 {
786  if (this->cyclicACMIPatch().owner())
787  {
788  label index = this->patch().index();
789 
790  const label globalPatchID =
791  matrix.lduMeshAssembly().patchLocalToGlobalMap()[mat][index];
792 
793  const Field<scalar> intCoeffsCmpt
794  (
795  matrix.internalCoeffs()[globalPatchID].component(cmpt)
796  );
797 
798  const Field<scalar> boundCoeffsCmpt
799  (
800  matrix.boundaryCoeffs()[globalPatchID].component(cmpt)
801  );
802 
803  tmp<Field<scalar>> tintCoeffs(coeffs(matrix, intCoeffsCmpt, mat));
804  tmp<Field<scalar>> tbndCoeffs(coeffs(matrix, boundCoeffsCmpt, mat));
805  const Field<scalar>& intCoeffs = tintCoeffs.ref();
806  const Field<scalar>& bndCoeffs = tbndCoeffs.ref();
807 
808  const labelUList& u = matrix.lduAddr().upperAddr();
809  const labelUList& l = matrix.lduAddr().lowerAddr();
810 
811  label subFaceI = 0;
812 
813  const labelList& faceMap =
814  matrix.lduMeshAssembly().faceBoundMap()[mat][index];
815 
816  forAll (faceMap, j)
817  {
818  label globalFaceI = faceMap[j];
819 
820  const scalar boundCorr = -bndCoeffs[subFaceI];
821  const scalar intCorr = -intCoeffs[subFaceI];
822 
823  matrix.upper()[globalFaceI] += boundCorr;
824  matrix.diag()[u[globalFaceI]] -= intCorr;
825  matrix.diag()[l[globalFaceI]] -= boundCorr;
826 
827  if (matrix.asymmetric())
828  {
829  matrix.lower()[globalFaceI] += intCorr;
830  }
831  subFaceI++;
832  }
833 
834  // Set internalCoeffs and boundaryCoeffs in the assembly matrix
835  // on cyclicACMI patches to be used in the individual matrix by
836  // matrix.flux()
837  if (matrix.psi(mat).mesh().fluxRequired(this->internalField().name()))
838  {
839  matrix.internalCoeffs().set
840  (
841  globalPatchID, intCoeffs*pTraits<Type>::one
842  );
843  matrix.boundaryCoeffs().set
844  (
845  globalPatchID, bndCoeffs*pTraits<Type>::one
846  );
847 
848  const label nbrPathID =
849  cyclicACMIPatch_.cyclicACMIPatch().neighbPatchID();
850 
851  const label nbrGlobalPatchID =
853  [mat][nbrPathID];
854 
855  matrix.internalCoeffs().set
856  (
857  nbrGlobalPatchID, intCoeffs*pTraits<Type>::one
858  );
859  matrix.boundaryCoeffs().set
860  (
861  nbrGlobalPatchID, bndCoeffs*pTraits<Type>::one
862  );
863  }
864  }
865 }
866 
867 
868 template<class Type>
871 (
872  fvMatrix<Type>& matrix,
873  const Field<scalar>& coeffs,
874  const label mat
875 ) const
876 {
877  const label index(this->patch().index());
878 
879  const label nSubFaces
880  (
881  matrix.lduMeshAssembly().cellBoundMap()[mat][index].size()
882  );
883 
884  Field<scalar> mapCoeffs(nSubFaces, Zero);
885 
886  const scalarListList& srcWeight =
887  cyclicACMIPatch_.cyclicACMIPatch().AMI().srcWeights();
888 
889  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
890 
891  const scalar tol = cyclicACMIPolyPatch::tolerance();
892  label subFaceI = 0;
893  forAll(mask, faceI)
894  {
895  const scalarList& w = srcWeight[faceI];
896  for(label i=0; i<w.size(); i++)
897  {
898  if (mask[faceI] > tol)
899  {
900  const label localFaceId =
901  matrix.lduMeshAssembly().facePatchFaceMap()
902  [mat][index][subFaceI];
903  mapCoeffs[subFaceI] = w[i]*coeffs[localFaceId];
904  }
905  subFaceI++;
906  }
907  }
908 
909  return tmp<Field<scalar>>(new Field<scalar>(mapCoeffs));
910 }
911 
912 
913 template<class Type>
915 {
916  // Update non-overlap patch - some will implement updateCoeffs, and
917  // others will implement evaluate
918 
919  // Pass in (1 - mask) to give non-overlap patch the chance to do
920  // manipulation of non-face based data
921 
922  const scalarField& mask = cyclicACMIPatch_.cyclicACMIPatch().mask();
923  const fvPatchField<Type>& npf = nonOverlapPatchField();
924  const_cast<fvPatchField<Type>&>(npf).updateWeightedCoeffs(1.0 - mask);
925 }
926 
927 
928 template<class Type>
930 {
933 
934  if (patchNeighbourFieldPtr_)
935  {
936  patchNeighbourFieldPtr_->writeEntry("neighbourValue", os);
937  }
938 }
939 
940 
941 // ************************************************************************* //
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)
dictionary dict
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
"blocking" : (MPI_Bsend, MPI_Recv)
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:72
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:598
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). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
virtual label nonOverlapPatchID() const
Return neighbour.
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
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
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:375
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
scalarField & upper()
Definition: lduMatrix.C:208
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:137
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.
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:670
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
Definition: lduMatrix.H:791
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)
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 lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:734
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
scalarField & lower()
Definition: lduMatrix.C:179
static const Field< Type > & null()
Return nullObject reference Field.
Definition: FieldI.H:24
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...
Definition: areaFieldsFwd.H:42
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
Definition: fvPatchField.H:662
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" : (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
scalarField & diag()
Definition: lduMatrix.C:197
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
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