31 template<
class TransferType>
36 List<TransferType>& faceDist
41 label nChangedFaces = 0;
43 forAll(mesh.boundaryMesh(), patchI)
45 if (patchIDs.found(patchI))
47 const polyPatch& patch = mesh.boundaryMesh()[patchI];
49 const Field<Type>& patchField = initialPatchValuePtrs_[patchI];
51 forAll(patch.faceCentres(), patchFaceI)
53 label meshFaceI = patch.start() + patchFaceI;
55 changedFaces[nChangedFaces] = meshFaceI;
57 faceDist[nChangedFaces] =
60 patch.faceCentres()[patchFaceI],
61 patchField[patchFaceI],
73 template<
class TransferType>
76 const MeshWave<TransferType>& waveInfo
81 const List<TransferType>& cellInfo = waveInfo.allCellInfo();
82 const List<TransferType>& faceInfo = waveInfo.allFaceInfo();
87 distance_.setSize(cellInfo.size());
91 const TransferType & wpn = cellInfo[cellI];
93 scalar dist = wpn.distSqr();
95 if (cellInfo[cellI].valid(waveInfo.data()))
99 cellData_[cellI] = cellInfo[cellI].data();
106 distance_[cellI] =
mag(dist);
109 cellData_[cellI] = cellInfo[cellI].data();
116 forAll(patchDistance_, patchI)
118 const polyPatch& patch = mesh.boundaryMesh()[patchI];
123 patchDistance_.set(patchI, patchFieldPtr);
128 Field<Type>* patchDataFieldPtr =
new Field<Type>(patch.size());
130 patchData_.set(patchI, patchDataFieldPtr);
132 Field<Type>& patchDataField = *patchDataFieldPtr;
135 forAll(patchField, patchFaceI)
137 label meshFaceI = patch.start() + patchFaceI;
139 scalar dist = faceInfo[meshFaceI].distSqr();
141 if (faceInfo[meshFaceI].valid(waveInfo.data()))
145 patchField[patchFaceI] =
Foam::sqrt(dist) + SMALL;
147 patchDataField[patchFaceI] = faceInfo[meshFaceI].data();
153 patchField[patchFaceI] =
mag(dist);
156 patchDataField[patchFaceI] = faceInfo[meshFaceI].data();
170 template<
class TransferType>
176 const bool correctWalls
181 initialPatchValuePtrs_(initialPatchValuePtrs),
182 correctWalls_(correctWalls),
195 template<
class TransferType>
203 template<
class TransferType>
212 label nWalls = sumPatchSize(patchIDs_);
217 setChangedFaces(patchIDs_, changedFaces, faceDist);
228 mesh().globalData().nTotalCells()+1
236 nUnset_ = getValues(waveInfo);
247 correctBoundaryFaceCells
254 correctBoundaryPointCells
266 forAll(wallCells, wallCellI)
268 label cellI = wallCells[wallCellI];
270 label faceI = nearestFace[cellI];
272 cellData_[cellI] = faceInfo[faceI].
data();
dimensionedScalar sqrt(const dimensionedScalar &ds)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
label size() const
Return the number of elements in the PtrList.
dimensioned< scalar > mag(const dimensioned< Type > &)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil)-psi)*pSat, rhoMin);#1"/home/dm2/henry/OpenFOAM/OpenFOAM-3.0.x/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H"1{alphav=max(min((rho-rholSat)/(rhovSat-rholSat), scalar(1)), scalar(0));alphal=1.0-alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
T * data()
Return a pointer to the first data element,.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
List< Key > toc() const
Return the table of contents.
Pre-declare SubField and related Field type.
const List< Type > & allFaceInfo() const
Get allFaceInfo.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Mesh consisting of general polyhedral cells.
virtual void correct()
Correct for mesh geom/topo changes.
List< label > labelList
A List of labels.
Collection of functions used in wall distance calculation.
virtual ~patchDataWave()
Destructor.
patchDataWave(const polyMesh &mesh, const labelHashSet &patchIDs, const UPtrList< Field< Type > > &initialPatchValuePtrs, bool correctWalls=true)
Construct from mesh, information on patches to initialize and flag.
volScalarField scalarField(fieldObject, mesh)
Takes a set of patches to start MeshWave from.