Loading...
Searching...
No Matches
ITHACAgeometry.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24File
25 ITHACAgeometry
26Description
27 set of methods for geometrical parametrization
28SourceFiles
29 ITHACAgeometry.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ITHACAgeometry_H
38#define ITHACAgeometry_H
39
40#include "fvCFD.H"
41#include "IOmanip.H"
42#include "freestreamFvPatchField.H"
43#include <sys/stat.h>
44#include <unistd.h>
45#pragma GCC diagnostic push
46#pragma GCC diagnostic ignored "-Wold-style-cast"
47#include <Eigen/Eigen>
48#pragma GCC diagnostic pop
49#include <functional>
50#include "./colormod.H"
51#include "polyMeshTools.H"
52#include <chrono>
53#include "mixedFvPatchFields.H"
54#include "fvMeshSubset.H"
55using namespace std::placeholders;
56#include "Foam2Eigen.H"
57
58namespace ITHACAutilities
59{
73labelList getIndicesFromBox(const fvMesh& mesh, List<label> indices,
74 Eigen::MatrixXd Box, List<vector>& points2Move);
75
76//--------------------------------------------------------------------------
85List<label> getIndices(const fvMesh& mesh, int index, int layers);
86
87//--------------------------------------------------------------------------
97List<label> getIndices(const fvMesh& mesh, int index_row, int index_col,
98 int layers);
99
100//--------------------------------------------------------------------------
108void getPointsFromPatch(fvMesh& mesh, label ind, List<vector>& points,
109 labelList& indices);
110
111//--------------------------------------------------------------------------
122List<vector> displacedSegment(List<vector> x0, double mux1, double mux2,
123 double muy1, double muy2);
124//----------------------------------------------------------------------
139vector displacePolabel(vector x0, vector x_low, vector x_up,
140 double mux_low, double mux_up,
141 double muy_low, double muy_up, double muz_low = 0.0, double muz_up = 0.0);
142
143//--------------------------------------------------------------------------
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);
164
165//--------------------------------------------------------------------------
173Eigen::MatrixXd rotationMatrix(vector AxisOfRotation,
174 double AngleOfRotation);
175//--------------------------------------------------------------------------
183Eigen::VectorXd boudaryFaceToCellDistance(fvMesh& mesh, label BC_ind);
184
185//--------------------------------------------------------------------------
197List<label> getIndicesFromDisc(const fvMesh& mesh, double radius,
198 vector origin, vector axis, List<double>& radii);
199
200//--------------------------------------------------------------------------
210template<typename type_f>
211List<label> getIndicesFromBox(
212 GeometricField<type_f, fvPatchField, volMesh>& field, Eigen::MatrixXd Box);
213
214//--------------------------------------------------------------------------
224template<typename type_f>
225fvMeshSubset* getSubMeshFromBox(
226 GeometricField<type_f, fvPatchField, volMesh>& field, Eigen::MatrixXd Box);
227
228//--------------------------------------------------------------------------
236volScalarField meshNonOrtho(fvMesh& mesh, volScalarField& NonOrtho);
237
238
239//--------------------------------------------------------------------------
248List<vector> rotatePoints(const List<vector>& originalPoints,
249 vector AxisOfRotation, double AngleOfRotation);
250
251}
252
253#endif
Header file of the Foam2Eigen class.
Foam::fvMesh & mesh
Definition createMesh.H:47
Simple header and source file of the Color::Modifier class to change color to the output stream.
volScalarField & v
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.