Go to the documentation of this file.
54 intFld.setSize(
mesh.nCells());
65 GeometricBoundaryField& bfld =
fld.boundaryField();
74 label unusedPatchI = 0;
76 forAll(oldPatchMap, patchI)
78 label newPatchI = oldPatchMap[patchI];
86 label nUsedPatches = unusedPatchI;
91 forAll(oldPatchMap, patchI)
93 label newPatchI = oldPatchMap[patchI];
97 oldToNew[patchI] = newPatchI;
101 oldToNew[patchI] = unusedPatchI++;
107 bfld.reorder(oldToNew);
109 bfld.setSize(
mesh.boundaryMesh().size());
113 label newPatchI = nUsedPatches;
114 newPatchI < bfld.size();
118 bfld.set(newPatchI, NULL);
125 forAll(oldPatchMap, patchI)
127 label newPatchI = oldPatchMap[patchI];
135 oldPatchStarts[patchI],
136 oldPatchSizes[patchI],
138 mesh.boundaryMesh()[newPatchI],
159 mesh.boundary()[newPatchI],
160 fld.dimensionedInternalField(),
177 forAll(addedPatchMap, patchI)
179 label newPatchI = addedPatchMap[patchI];
185 fldToAdd.mesh().boundaryMesh()[patchI];
187 if (!bfld(newPatchI))
213 mesh.boundary()[newPatchI],
214 fld.dimensionedInternalField(),
224 labelList addedToNew(oldPatch.size(), -1);
230 if (patchFaceI >= 0 && patchFaceI < newPatch.size())
232 addedToNew[i] = patchFaceI;
258 mesh.objectRegistry::lookupClass
265 meshToAdd.objectRegistry::lookupClass
278 iterator fieldIter =
fields.begin();
279 fieldIter !=
fields.end();
285 Pout<<
"MapVolFields : Storing old time for " << fieldIter()->
name()
297 iterator fieldIter =
fields.begin();
298 fieldIter !=
fields.end();
311 *fieldsToAdd[
fld.name()];
315 Pout<<
"MapVolFields : mapping " <<
fld.name()
316 <<
" and " << fldToAdd.name() <<
endl;
319 MapVolField<Type>(meshMap,
fld, fldToAdd);
324 <<
"Not mapping field " <<
fld.name()
325 <<
" since not present on mesh to add"
345 GeometricBoundaryField& bfld =
fld.boundaryField();
357 intFld.setSize(
mesh.nInternalFaces());
370 label start = oldPatchStarts[patchI];
376 if (newFaceI >= 0 && newFaceI <
mesh.nInternalFaces())
378 intFld[newFaceI] = pf[i];
394 label unusedPatchI = 0;
396 forAll(oldPatchMap, patchI)
398 label newPatchI = oldPatchMap[patchI];
406 label nUsedPatches = unusedPatchI;
411 forAll(oldPatchMap, patchI)
413 label newPatchI = oldPatchMap[patchI];
417 oldToNew[patchI] = newPatchI;
421 oldToNew[patchI] = unusedPatchI++;
427 bfld.reorder(oldToNew);
429 bfld.setSize(
mesh.boundaryMesh().size());
433 label newPatchI = nUsedPatches;
434 newPatchI < bfld.size();
438 bfld.set(newPatchI, NULL);
445 forAll(oldPatchMap, patchI)
447 label newPatchI = oldPatchMap[patchI];
455 oldPatchStarts[patchI],
456 oldPatchSizes[patchI],
458 mesh.boundaryMesh()[newPatchI],
478 mesh.boundary()[newPatchI],
479 fld.dimensionedInternalField(),
496 forAll(addedPatchMap, patchI)
498 label newPatchI = addedPatchMap[patchI];
504 fldToAdd.mesh().boundaryMesh()[patchI];
506 if (!bfld(newPatchI))
532 mesh.boundary()[newPatchI],
533 fld.dimensionedInternalField(),
543 labelList addedToNew(oldPatch.size(), -1);
549 if (patchFaceI >= 0 && patchFaceI < newPatch.size())
551 addedToNew[i] = patchFaceI;
579 mesh.objectRegistry::lookupClass<fldType>()
584 meshToAdd.objectRegistry::lookupClass<fldType>()
595 iterator fieldIter =
fields.begin();
596 fieldIter !=
fields.end();
602 Pout<<
"MapSurfaceFields : Storing old time for "
606 const_cast<fldType*
>(fieldIter())->storeOldTimes();
613 iterator fieldIter =
fields.begin();
614 fieldIter !=
fields.end();
618 fldType&
fld =
const_cast<fldType&
>(*fieldIter());
622 const fldType& fldToAdd = *fieldsToAdd[
fld.name()];
626 Pout<<
"MapSurfaceFields : mapping " <<
fld.name()
627 <<
" and " << fldToAdd.name() <<
endl;
630 MapSurfaceField<Type>(meshMap,
fld, fldToAdd);
635 <<
"Not mapping field " <<
fld.name()
636 <<
" since not present on mesh to add"
678 mesh.objectRegistry::lookupClass<fldType>(
true)
683 meshToAdd.objectRegistry::lookupClass<fldType>(
true)
689 iterator fieldIter =
fields.begin();
690 fieldIter !=
fields.end();
694 fldType&
fld =
const_cast<fldType&
>(*fieldIter());
698 const fldType& fldToAdd = *fieldsToAdd[
fld.name()];
702 Pout<<
"MapDimFields : mapping " <<
fld.name()
703 <<
" and " << fldToAdd.name() <<
endl;
706 MapDimField<Type>(meshMap,
fld, fldToAdd);
711 <<
"Not mapping field " <<
fld.name()
712 <<
" since not present on mesh to add"
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
const labelList & addedCellMap() const
#define forAll(list, i)
Loop across all elements in list.
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const labelList & addedFaceMap() const
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A patch is a list of labels that address the faces in the global face list.
virtual const fileName & name() const
Return the name of the stream.
static void MapVolFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all volFields of Type.
InternalField & internalField()
Return internal field.
static void MapSurfaceField(const mapAddedPolyMesh &meshMap, GeometricField< Type, fvsPatchField, surfaceMesh > &fld, const GeometricField< Type, fvsPatchField, surfaceMesh > &fldToAdd)
Update single surfaceField.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
const labelList & oldCellMap() const
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool found(const Key &) const
Return true if hashedEntry is found in table.
Mesh data needed to do the Finite Volume discretisation.
label start() const
Return start label of this patch in the polyMesh face list.
const labelList & oldFaceMap() const
const labelList & addedPatchMap() const
From added mesh patch index to new patch index or -1 if.
An STL-conforming hash table.
static void MapDimField(const mapAddedPolyMesh &meshMap, DimensionedField< Type, volMesh > &fld, const DimensionedField< Type, volMesh > &fldToAdd)
Update single dimensionedField.
prefixOSstream Pout(cout, "Pout")
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const labelList & oldPatchSizes() const
Return list of the old patch sizes.
static void MapDimFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all DimensionedFields of Type.
static void MapVolField(const mapAddedPolyMesh &meshMap, GeometricField< Type, fvPatchField, volMesh > &fld, const GeometricField< Type, fvPatchField, volMesh > &fldToAdd)
Update single volField.
direct fvPatchFieldMapper
void size(const label)
Override size to be inconsistent with allocated storage.
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh ('added m...
Generic GeometricField class.
#define WarningInFunction
Report a warning using Foam::Warning.
const labelList & oldPatchMap() const
From old patch index to new patch index or -1 if patch.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static void MapSurfaceFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all surfaceFields of Type.