Go to the documentation of this file.
48 template<
class ReadGeoField,
class MappedGeoField>
55 const typename MappedGeoField::value_type& nullValue,
59 typedef typename MappedGeoField::value_type Type;
64 tetFields.
setSize(fieldObjects.size());
69 Info<<
"Converting " << ReadGeoField::typeName <<
' ' << iter.key()
72 ReadGeoField readField(*iter(),
mesh);
87 readField.registerObject()
97 fld.setSize(map.size(), nullValue);
100 label index = map[pointi];
104 label celli = index-1;
105 fld[pointi] = readField[celli];
109 label facei = -index-1;
118 fld[pointi] = readField.boundaryField()[patchi][localFacei];
136 tetFields[i].correctBoundaryConditions();
143 int main(
int argc,
char *argv[])
147 "Convert polyMesh results to tetDualMesh"
164 Info<<
"Create tetDualMesh for time = "
182 "pointDualAddressing",
190 if (pointDualAddressing.size() != tetDualMesh.
nPoints())
193 <<
"Size " << pointDualAddressing.size()
194 <<
" of addressing map " << pointDualAddressing.objectPath()
195 <<
" differs from number of points in mesh "
203 label nPatchFaces = 0;
205 forAll(pointDualAddressing, pointi)
207 label index = pointDualAddressing[pointi];
219 label facei = -index-1;
223 <<
"Face " << facei <<
" from index " << index
224 <<
" is not a boundary face."
238 Info<<
"tetDualMesh points : " << tetDualMesh.
nPoints()
239 <<
" of which mapped to" <<
nl
240 <<
" cells : " << nCells <<
nl
241 <<
" patch faces : " << nPatchFaces <<
nl
242 <<
" not mapped : " << nUnmapped <<
nl
251 ReadAndMapFields<volScalarField, pointScalarField>
262 ReadAndMapFields<volVectorField, pointVectorField>
273 ReadAndMapFields<volSphericalTensorField, pointSphericalTensorField>
284 ReadAndMapFields<volSymmTensorField, pointSymmTensorField>
295 ReadAndMapFields<volTensorField, pointTensorField>
305 tetDualMesh.objectRegistry::write();
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static constexpr const zero Zero
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
static word timeName(const scalar t, const int precision=precision_)
static void addNote(const string ¬e)
const fileName & facesInstance() const
const polyBoundaryMesh & boundaryMesh() const
Ostream & endl(Ostream &os)
const T * set(const label i) const
label nPoints() const noexcept
const labelList & patchID() const
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
instantList times() const
void setSize(const label newLen)
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
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))
Mesh data needed to do the Finite Volume discretisation.
List of IOobjects with searching and retrieving facilities.
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
forAllConstIters(mixture.phases(), phase)
label nInternalFaces() const noexcept
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
virtual void setTime(const Time &t)
IOobjectList lookupClass(const char *clsName) const