Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
bf8df5f5
Commit
bf8df5f5
authored
Apr 05, 2013
by
andy
Browse files
BUG: Corrected sampling - handling of field regExps and field surface output
parent
42f49505
Changes
4
Hide whitespace changes
Inline
Side-by-side
src/sampling/sampledSurface/sampledSurface/sampledSurface.C
View file @
bf8df5f5
...
...
@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-201
2
OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-201
3
OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -259,6 +259,7 @@ Foam::tmp<Foam::vectorField> Foam::sampledSurface::sample
return
tmp
<
vectorField
>
(
NULL
);
}
Foam
::
tmp
<
Foam
::
sphericalTensorField
>
Foam
::
sampledSurface
::
sample
(
const
surfaceSphericalTensorField
&
sField
...
...
src/sampling/sampledSurface/sampledSurfaces/sampledSurfaces.C
View file @
bf8df5f5
...
...
@@ -171,17 +171,19 @@ void Foam::sampledSurfaces::write()
writeGeometry
();
}
sampleAndWrite
<
volScalarField
>
();
sampleAndWrite
<
volVectorField
>
();
sampleAndWrite
<
volSphericalTensorField
>
();
sampleAndWrite
<
volSymmTensorField
>
();
sampleAndWrite
<
volTensorField
>
();
sampleAndWrite
<
surfaceScalarField
>
();
sampleAndWrite
<
surfaceVectorField
>
();
sampleAndWrite
<
surfaceSphericalTensorField
>
();
sampleAndWrite
<
surfaceSymmTensorField
>
();
sampleAndWrite
<
surfaceTensorField
>
();
const
IOobjectList
objects
(
mesh_
,
mesh_
.
time
().
timeName
());
sampleAndWrite
<
volScalarField
>
(
objects
);
sampleAndWrite
<
volVectorField
>
(
objects
);
sampleAndWrite
<
volSphericalTensorField
>
(
objects
);
sampleAndWrite
<
volSymmTensorField
>
(
objects
);
sampleAndWrite
<
volTensorField
>
(
objects
);
sampleAndWrite
<
surfaceScalarField
>
(
objects
);
sampleAndWrite
<
surfaceVectorField
>
(
objects
);
sampleAndWrite
<
surfaceSphericalTensorField
>
(
objects
);
sampleAndWrite
<
surfaceSymmTensorField
>
(
objects
);
sampleAndWrite
<
surfaceTensorField
>
(
objects
);
}
}
...
...
src/sampling/sampledSurface/sampledSurfaces/sampledSurfaces.H
View file @
bf8df5f5
...
...
@@ -89,6 +89,7 @@ class sampledSurfaces
//- Tolerance for merging points (fraction of mesh bounding box)
static
scalar
mergeTol_
;
// Private data
//- Name of this set of surfaces,
...
...
@@ -113,6 +114,7 @@ class sampledSurfaces
//- Interpolation scheme to use
word
interpolationScheme_
;
// surfaces
//- Information for merging surfaces
...
...
@@ -159,7 +161,7 @@ class sampledSurfaces
);
//- Sample and write all sampled fields
template
<
class
Type
>
void
sampleAndWrite
();
template
<
class
Type
>
void
sampleAndWrite
(
const
IOobjectList
&
objects
);
//- Disallow default bitwise copy construct and assignment
sampledSurfaces
(
const
sampledSurfaces
&
);
...
...
@@ -233,7 +235,6 @@ public:
//- Update for changes of mesh due to readUpdate - expires the surfaces
virtual
void
readUpdate
(
const
polyMesh
::
readUpdateState
state
);
};
...
...
src/sampling/sampledSurface/sampledSurfaces/sampledSurfacesTemplates.C
View file @
bf8df5f5
...
...
@@ -166,36 +166,26 @@ void Foam::sampledSurfaces::sampleAndWrite
template
<
class
GeoField
>
void
Foam
::
sampledSurfaces
::
sampleAndWrite
()
void
Foam
::
sampledSurfaces
::
sampleAndWrite
(
const
IOobjectList
&
objects
)
{
DynamicList
<
wordRe
>
fields
;
bool
found
=
false
;
forAll
(
fieldSelection_
,
fieldI
)
if
(
loadFromFiles_
)
{
const
wordRe
field
=
fieldSelection_
[
fieldI
];
if
(
mesh_
.
thisDb
().
foundObject
<
GeoField
>
(
field
))
{
found
=
true
;
fields
.
append
(
field
);
}
}
if
(
fields
.
size
()
&&
found
)
{
forAll
(
fields
,
fieldI
)
IOobjectList
fieldObjects
(
objects
.
lookupClass
(
GeoField
::
typeName
));
forAll
(
fieldSelection_
,
fieldI
)
{
const
wordRe
field
=
fields
[
fieldI
];
if
((
Pstream
::
master
())
&&
verbose_
)
{
Pout
<<
"sampleAndWrite: "
<<
field
<<
endl
;
}
const
wordRe
fieldNameRe
=
fieldSelection_
[
fieldI
];
IOobjectList
fieldIO
=
fieldObjects
.
lookupRe
(
fieldNameRe
);
if
(
loadFromFiles_
)
forAllConstIter
(
IOobjectList
,
fieldIO
,
iter
)
{
const
GeoField
geoField
const
word
&
fieldName
=
iter
()
->
name
();
const
GeoField
fld
(
IOobject
(
field
,
field
Name
,
mesh_
.
time
().
timeName
(),
mesh_
,
IOobject
::
MUST_READ
...
...
@@ -203,13 +193,38 @@ void Foam::sampledSurfaces::sampleAndWrite()
mesh_
);
sampleAndWrite
(
geoField
);
if
((
Pstream
::
master
())
&&
verbose_
)
{
Pout
<<
"sampleAndWrite: "
<<
fieldName
<<
endl
;
}
sampleAndWrite
(
fld
);
}
else
}
}
else
{
forAll
(
fieldSelection_
,
fieldI
)
{
const
wordRe
&
fieldNameRe
=
fieldSelection_
[
fieldI
];
const
wordList
dbFields
(
mesh_
.
thisDb
().
foundObjectRe
<
GeoField
>
(
fieldNameRe
)
);
forAll
(
dbFields
,
i
)
{
const
word
&
fieldName
=
dbFields
[
i
];
if
((
Pstream
::
master
())
&&
verbose_
)
{
Pout
<<
"sampleAndWrite: "
<<
fieldName
<<
endl
;
}
sampleAndWrite
(
mesh_
.
thisDb
().
lookupObject
<
GeoField
>
(
field
)
mesh_
.
thisDb
().
lookupObject
<
GeoField
>
(
field
Name
)
);
}
}
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment