64 const scalar
f = factor[celli];
69 s += w[nbrI]*work[nbrs[nbrI]];
79template<
class GeoField>
87template<
class GeoField>
90 auto flds(this->lookupClass<GeoField>());
91 for (
auto fldPtr : flds)
98 Pout<<
"dynamicOversetFvMesh::interpolate: interpolating : "
107 Pout<<
"dynamicOversetFvMesh::interpolate: skipping : " <<
name
115template<
class GeoField,
class PatchType>
118 typename GeoField::Boundary& bfld,
126 if (typeOnly == (isA<PatchType>(bfld[patchi]) !=
nullptr))
144 if (typeOnly == (isA<PatchType>(bfld[patchi]) !=
nullptr))
167 for (
direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
169 forAll(internalCoeffs, patchi)
171 const labelUList& fc = lduAddr().patchAddr(patchi);
172 const Field<Type>& intCoeffs = internalCoeffs[patchi];
176 norm[fc[i]] += cmptCoeffs[i];
185 const scalar&
n = norm[celli];
199 Pout<<
"For field " << m.
psi().name() <<
" have zero diagonals for "
200 << nZeroDiag <<
" cells" <<
endl;
218 extrapolatedNorm[celli] = -GREAT;
224 for (label facei = 0; facei < nInternalFaces(); facei++)
226 label ownType = types[own[facei]];
227 label neiType = types[nei[facei]];
239 for (label facei = nInternalFaces(); facei <
nFaces(); facei++)
241 label ownType = types[own[facei]];
242 label neiType = nbrTypes[facei-nInternalFaces()];
263 for (
const label facei : isFront)
265 if (extrapolatedNorm[own[facei]] == -GREAT)
268 newNorm[own[facei]] = cellAverage
281 isInternalFace(facei)
282 && extrapolatedNorm[nei[facei]] == -GREAT
286 newNorm[nei[facei]] = cellAverage
307 isFront.transfer(newIsFront);
314 scalar&
n = norm[celli];
324 n = extrapolatedNorm[celli];
363 if (!isA<fvMeshPrimitiveLduAddressing>(addr))
366 <<
"Problem : addressing is not fvMeshPrimitiveLduAddressing"
373 upper.setSize(upperAddr.
size(), 0.0);
375 lower.setSize(lowerAddr.
size(), 0.0);
383 if (interfaces.
size() > nOldInterfaces)
392 label patchi = nOldInterfaces;
393 patchi < interfaces.
size();
397 const labelUList& fc = interfaces[patchi].faceCells();
406 typename GeoField::Boundary& bfld =
407 const_cast<GeoField&
>(m.
psi()).boundaryFieldRef();
416 for (label patchi = 0; patchi < nOldInterfaces; patchi++)
418 if (isA<processorFvPatch>(bfld[patchi].patch()))
427 label patchi = nOldInterfaces;
428 patchi < interfaces.
size();
438 bfld[addPatchi].patch(),
463 label celli = upperAddr[facei];
464 lower[facei] *= 1.0-factor[celli];
469 label celli = lowerAddr[facei];
470 upper[facei] *= 1.0-factor[celli];
474 for (label patchi = 0; patchi < nOldInterfaces; ++patchi)
485 scalar
f = factor[celli];
486 intCoeffs[i] *= 1.0-
f;
487 bouCoeffs[i] *= 1.0-
f;
508 label startLabel = ownerStartAddr[celli];
509 label endLabel = ownerStartAddr[celli + 1];
511 for (label facei = startLabel; facei < endLabel; facei++)
519 for (label i = startLabel; i < endLabel; i++)
521 label facei = losortAddr[i];
525 diag[celli] = normalisation[celli];
526 source[celli] = normalisation[celli]*m.
psi()[celli];
544 const scalar
f = factor[celli];
547 const labelList& nbrFaces = stencilFaces_[celli];
548 const labelList& nbrPatches = stencilPatches_[celli];
553 <<
" at:" <<
C()[celli]
554 <<
" . Should this be in interpolationCells()????"
561 diag[celli] *= (1.0-
f);
562 source[celli] *= (1.0-
f);
566 label patchi = nbrPatches[nbri];
567 label facei = nbrFaces[nbri];
571 label nbrCelli = nbrs[nbri];
574 const scalar
s = normalisation[celli]*
f*w[nbri];
576 scalar& u = upper[facei];
577 scalar& l = lower[facei];
578 if (celli < nbrCelli)
606 const scalar
s = normalisation[celli]*
f*w[nbri];
657 typename GeoField::Boundary& bpsi =
658 const_cast<GeoField&
>(m.
psi()).boundaryFieldRef();
660 bool hasOverset =
false;
674 Pout<<
"dynamicOversetFvMesh::solve() :"
675 <<
" bypassing matrix adjustment for field " << m.
psi().name()
684 Pout<<
"dynamicOversetFvMesh::solve() :"
685 <<
" adjusting matrix for interpolation for field "
698 m.
psi().name() +
"_normalisation",
699 this->time().timeName(),
718 Pout<<
"dynamicOversetFvMesh::solve() :"
719 <<
" writing matrix normalisation for field " << m.
psi().name()
734 const label oldSize = bpsi.size();
736 addInterpolation(m, norm);
739 correctBoundaryConditions<GeoField, calculatedProcessorFvPatchField<Type>>
759 bpsi.setSize(oldSize);
783 os <<
"** Matrix **" <<
endl;
799 forAll(interfaces, patchi)
801 if (interfaces.
set(patchi))
803 const labelUList& fc = interfaces[patchi].faceCells();
807 cellToPatch[fc[i]].
append(patchi);
808 cellToPatchFace[fc[i]].
append(i);
816 os <<
"cell:" << celli <<
" diag:" <<
diag[celli]
817 <<
" source:" << source[celli] <<
endl;
819 label startLabel = ownerStartAddr[celli];
820 label endLabel = ownerStartAddr[celli + 1];
822 for (label facei = startLabel; facei < endLabel; facei++)
824 os <<
" face:" << facei
825 <<
" nbr:" << upperAddr[facei] <<
" coeff:" << upper[facei]
832 for (label i = startLabel; i < endLabel; i++)
834 label facei = losortAddr[i];
835 os <<
" face:" << facei
836 <<
" nbr:" << lowerAddr[facei] <<
" coeff:" << lower[facei]
840 forAll(cellToPatch[celli], i)
842 label patchi = cellToPatch[celli][i];
843 label patchFacei = cellToPatchFace[celli][i];
845 os <<
" patch:" << patchi
846 <<
" patchface:" << patchFacei
858 os <<
"patch:" << patchi
863 interfaces.
set(patchi)
864 && isA<processorLduInterface>(interfaces[patchi])
868 refCast<const processorLduInterface>(interfaces[patchi]);
876 os <<
" cell:" << fc[i]
886 m.
psi().boundaryField().scalarInterfaces();
887 forAll(interfaceFields, inti)
889 if (interfaceFields.
set(inti))
891 os <<
"interface:" << inti
892 <<
" if type:" << interfaceFields[inti].interface().type()
893 <<
" fld type:" << interfaceFields[inti].type() <<
endl;
897 os <<
"** End of Matrix **" <<
endl;
901template<
class GeoField>
904 typename GeoField::Boundary& bfld =
fld.boundaryFieldRef();
938template<
class GeoField>
941 Pout<<
"** starting checking of " <<
fld.name() <<
endl;
943 GeoField fldCorr(
fld.name()+
"_correct",
fld);
944 correctCoupledBoundaryConditions(fldCorr);
946 const typename GeoField::Boundary& bfld =
fld.boundaryField();
947 const typename GeoField::Boundary& bfldCorr = fldCorr.boundaryField();
951 const typename GeoField::Patch& pfld = bfld[patchi];
952 const typename GeoField::Patch& pfldCorr = bfldCorr[patchi];
958 if (pfld[i] != pfldCorr[i])
960 Pout<<
" " << i <<
" orig:" << pfld[i]
961 <<
" corrected:" << pfldCorr[i] <<
endl;
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Graphite solid properties.
const Field< Type > & field() const
Return field.
A field of fields is a PtrList of fields with reference counting.
Generic templated field type.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
void correctBoundaryConditions()
Correct boundary field.
bool found(const Key &key) const
Return true if hashed entry is found in table.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
void transfer(List< T > &list)
void setSize(const label n)
Alias for resize()
void append(const T &val)
Append an element at the end of the list.
virtual const fileName & name() const
Get the name of the stream.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
void size(const label n)
Older name for setAddressableSize.
@ nonBlocking
"nonBlocking"
static int & msgType() noexcept
Message tag of standard messages.
static label nRequests()
Get number of outstanding requests.
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
static commsTypes defaultCommsType
Default commsType.
static bool & parRun() noexcept
Test if this a parallel run.
const T * set(const label i) const
label size() const noexcept
The number of elements in the list.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
void set(const bitSet &bitset)
Set specified bits from another bitset.
A processorFvPatchField type bypassing fvPatch.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
virtual const labelUList & cellTypes() const
Return the cell type list.
virtual const List< scalarList > & cellInterpolationWeights() const
Weights for cellStencil.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
Calculation of interpolation stencils.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
static void checkCoupledBC(const GeoField &fld)
Debug: check halo swap is ok.
void addInterpolation(fvMatrix< Type > &, const scalarField &norm) const
Add interpolation to matrix (coefficients)
tmp< scalarField > normalisation(const fvMatrix< Type > &m) const
Freeze values at holes.
static void correctCoupledBoundaryConditions(GeoField &fld)
Debug: correct coupled bc.
lduInterfacePtrsList interfaces() const
virtual bool write()
Write the output fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
const FieldField< Field, Type > & internalCoeffs() const noexcept
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Field< Type > & source() noexcept
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const labelUList & ownerStartAddr() const
Return owner start addressing.
const labelUList & losortStartAddr() const
Return losort start addressing.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
label size() const
Return number of equations.
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.
const labelUList & losortAddr() const
Return losort addressing.
Class containing processor-to-processor mapping information.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
label nCells() const noexcept
Number of mesh cells.
An abstract base class for processor coupled interfaces.
virtual int neighbProcNo() const =0
Return neighbour processor number (rank in communicator)
virtual int myProcNo() const =0
Return processor number (rank in communicator)
virtual bool write(const bool valid=true) const
Write using setting from DB.
bool interpolate() const noexcept
Same as isPointData()
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
const volScalarField & psi
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const cellCellStencilObject & overlap
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
static constexpr const zero Zero
Global zero (0)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
cellMask correctBoundaryConditions()
#define forAll(list, i)
Loop across all elements in list.
A non-counting (dummy) refCount.