42#include "freestreamFvPatchField.H"
45#pragma GCC diagnostic push
46#pragma GCC diagnostic ignored "-Wold-style-cast"
48#pragma GCC diagnostic pop
51#include "polyMeshTools.H"
53#include "mixedFvPatchFields.H"
54#include "fvMeshSubset.H"
55using namespace std::placeholders;
76template<
class TypeField>
78 fields, Eigen::MatrixXd ind, PtrList<TypeField>& ave);
89template<
class TypeField>
100template<
typename Type>
101void assignIF(GeometricField<Type, fvPatchField, volMesh>& field,
113template<
typename Type>
114void assignIF(GeometricField<Type, fvPatchField, volMesh>& field,
115 Type& value, List<label>& indices);
126template<
typename Type>
127void assignIF(GeometricField<Type, fvPatchField, volMesh>& field,
128 Type& value, label index);
137void assignONE(volScalarField& s, List<label>& L);
147void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
157void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
158 Eigen::MatrixXd valueVec);
167void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
168 List<double> valueList);
177void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
187void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
197void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
198 Eigen::MatrixXd valueVec);
207void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
208 Eigen::MatrixXd valueVec);
217void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
218 List<vector> valueList);
227void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
228 List<tensor> valueList);
239void assignBC(GeometricField<scalar, fvsPatchField, surfaceMesh>& field,
240 label BC_ind, Eigen::MatrixXd value);
242void assignBC(GeometricField<vector, fvsPatchField, surfaceMesh>& field,
243 label BC_ind, Eigen::MatrixXd value);
254template<
typename Type>
255void assignBC(GeometricField<Type, fvsPatchField, surfaceMesh>& field,
256 label BC_ind, List<Type>& value);
269template<
typename Type>
270void assignBC(GeometricField<Type, fvsPatchField, surfaceMesh>& field,
271 label BC_ind, Type& value);
282template<
typename Type>
296template<
typename Type>
300#pragma GCC diagnostic push
301#pragma GCC diagnostic ignored "-Wcomment"
323template<
typename Type>
324void setBoxToValue(GeometricField<Type, fvPatchField, volMesh>& field,
325 Eigen::MatrixXd Box, Type value);
326#pragma GCC diagnostic pop
338void changeBCtype(GeometricField<Type, fvPatchField, volMesh>&
339 field, word BCtype, label BC_ind);
352template<
typename Type>
354 labelList& movingIDS, List<Type>& originalList);
367template<
typename Type>
369 label BC_ind, Eigen::MatrixXd& value,
370 Eigen::MatrixXd& grad, Eigen::MatrixXd& valueFrac);
384template<
typename Type>
385void assignMixedBC(GeometricField<Type, fvPatchField, volMesh>& field,
386 label BC_ind, List<Type>& value,
387 List<Type>& grad, List<scalar>& valueFrac);
397template<
typename Type>
399 PtrList<GeometricField<Type, fvPatchField, volMesh>>& fields);
412template<
typename Type>
413Eigen::MatrixXd
getValues(GeometricField<Type, fvPatchField,
414 volMesh>& field, labelList& indices, labelList* xyz = NULL);
427template<
typename Type>
428Eigen::MatrixXd
getValues(PtrList<GeometricField<Type, fvPatchField,
429 volMesh>>& field, labelList& indices, labelList* xyz = NULL);
Header file of the Foam2Eigen class.
Header file of the ITHACAerror file.
Simple header and source file of the Color::Modifier class to change color to the output stream.
Namespace to implement some useful assign operation of OF fields.
void assignMixedBC(GeometricField< Type, fvPatchField, volMesh > &field, label BC_ind, List< Type > &value, List< Type > &grad, List< scalar > &valueFrac)
Assign value of a boundary condition of type "mixed".
void assignIF(GeometricField< Type, fvPatchField, volMesh > &s, Type value)
Assign internal field.
void changeNeumann2Dirichlet(GeometricField< Type, fvPatchField, volMesh > &field, Type &value)
Change all Neumann boundary conditions to Dirichlet boundary conditions.
TypeField computeAverage(PtrList< TypeField > &fields)
Calculates the average of a list of fields.
void normalizeFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &fields)
Normalize list of Geometric fields.
void assignONE(volScalarField &s, List< label > &L)
Assign one to volScalarField.
void setIndices2Value(labelList &ind2set, List< Type > &value2set, labelList &movingIDS, List< Type > &originalList)
Sets some given Indices of a list of objects to given values.
PtrList< TypeField > averageSubtract(PtrList< TypeField > fields, Eigen::MatrixXd ind, PtrList< TypeField > &ave)
A function to compute time-averaged fields for a set of different parameter samples and also the fiel...
void assignZeroDirichlet(GeometricField< Type, fvPatchField, volMesh > &field)
Assign zero internal field.
void assignBC(GeometricField< scalar, fvPatchField, volMesh > &s, label BC_ind, double value)
Assign uniform Boundary Condition to a volScalarField.
Eigen::MatrixXd getValues(GeometricField< Type, fvPatchField, volMesh > &field, labelList &indices)
void changeBCtype(GeometricField< Type, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
Change the boundary condition type for a GeometricField.
void setBoxToValue(GeometricField< Type, fvPatchField, volMesh > &field, Eigen::MatrixXd Box, Type value)
Set value of a volScalarField to a constant inside a given box.