37#ifndef ITHACAgeometry_H
38#define ITHACAgeometry_H
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;
74 Eigen::MatrixXd Box, List<vector>& points2Move);
97List<label>
getIndices(
const fvMesh&
mesh,
int index_row,
int index_col,
123 double muy1,
double muy2);
140 double mux_low,
double mux_up,
141 double muy_low,
double muy_up,
double muz_low = 0.0,
double muz_up = 0.0);
161Field<vector>
rotateMesh(fvMesh&
mesh,
double r1,
double r2, vector axis,
162 double alpha, labelList movingPointsIDs,
163 List<double> radii, word angleVariationMethod =
"Linear",
double v = 3);
174 double AngleOfRotation);
198 vector origin, vector axis, List<double>& radii);
210template<
typename type_f>
212 GeometricField<type_f, fvPatchField, volMesh>& field, Eigen::MatrixXd Box);
224template<
typename type_f>
226 GeometricField<type_f, fvPatchField, volMesh>& field, Eigen::MatrixXd Box);
248List<vector>
rotatePoints(
const List<vector>& originalPoints,
249 vector AxisOfRotation,
double AngleOfRotation);
Header file of the Foam2Eigen class.
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.
Field< vector > rotateMesh(fvMesh &mesh, double r1, double r2, vector axis, double alpha, labelList movingPointsIDs, List< double > radii, word angleVariationMethod, double v)
A function that rotates the mesh rigidly, by a fixed specified angle for the points lying within a di...
List< label > getIndices(const fvMesh &mesh, int index, int layers)
Gets the indices of the cells around a certain cell.
vector displacePolabel(vector x0, vector x_low, vector x_up, double mux_low, double mux_up, double muy_low, double muy_up, double muz_low, double muz_up)
Displace a Polabel belonging to a given segment.
fvMeshSubset * getSubMeshFromBox(GeometricField< type_f, fvPatchField, volMesh > &field, Eigen::MatrixXd Box)
Gets a subMesh from a box of coordinates and a given field (used only for the mesh).
Eigen::MatrixXd rotationMatrix(vector AxisOfRotation, double AngleOfRotation)
Functions that return a Rotation Matrix given an axis of rotation and an angle in degrees.
List< vector > displacedSegment(List< vector > x0, double mux1, double mux2, double muy1, double muy2)
Get position of displaced segment of points given the displacements of the end points.
List< label > getIndicesFromDisc(const fvMesh &mesh, double radius, vector origin, vector axis, List< double > &radii)
Gets the indices of the points which are inside a disc of a specified radius before carrying out a ro...
volScalarField meshNonOrtho(fvMesh &mesh, volScalarField &NonOrtho)
Returns a scalarField that containes the non-orthogonality value of a given mesh.
void getPointsFromPatch(fvMesh &mesh, label ind, List< vector > &points, labelList &indices)
Get the polabel coordinates and indices from patch.
labelList getIndicesFromBox(const fvMesh &mesh, List< label > indices, Eigen::MatrixXd Box, List< vector > &points2Move)
Gives the indices conteined into a defined box.
Eigen::VectorXd boudaryFaceToCellDistance(fvMesh &mesh, label BC_ind)
Compute the distance between the boundary face center and the boundary cell center.
List< vector > rotatePoints(const List< vector > &originalPoints, vector AxisOfRotation, double AngleOfRotation)
Rotate a list of points in clockwise direction given an axis of rotation and an angle in degrees.