NURBS3DVolumeCylindrical.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2007-2023 PCOpt/NTUA
9  Copyright (C) 2013-2023 FOSS GP
10  Copyright (C) 2019 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 \*---------------------------------------------------------------------------*/
29 
32 #include "pointMesh.H"
33 #include "pointPatchField.H"
34 #include "pointPatchFieldsFwd.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41  defineTypeNameAndDebug(NURBS3DVolumeCylindrical, 0);
43  (
44  NURBS3DVolume,
45  NURBS3DVolumeCylindrical,
46  dictionary
47  );
48 }
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 (
55  const vector& localSystemCoordinates
56 ) const
57 {
58  vector cartesianCoors
59  (
60  localSystemCoordinates.x()*cos(localSystemCoordinates.y()),
61  localSystemCoordinates.x()*sin(localSystemCoordinates.y()),
62  localSystemCoordinates.z()
63  );
64  cartesianCoors += origin_;
65 
66  return cartesianCoors;
67 }
68 
69 
71 (
72  label globalPointIndex
73 )
74 {
75  const vector& localCoors = localSystemCoordinates_[globalPointIndex];
76  tensor transformTensor
77  (
78  cos(localCoors.y()), -localCoors.x()*sin(localCoors.y()), 0,
79  sin(localCoors.y()), localCoors.x()*cos(localCoors.y()), 0,
80  0, 0, 1
81  );
82 
83  return transformTensor;
84 }
85 
86 
88 (
89  const vectorField& cartesianPoints
90 )
91 {
92  forAll(cartesianPoints, pI)
93  {
94  const vector point(cartesianPoints[pI] - origin_);
95  vector cylindricalCoors(Zero);
96 
97  const scalar R(Foam::sqrt(sqr(point.x()) + sqr(point.y())));
98  const scalar theta(atan2(point.y(), point.x()));
99  cylindricalCoors.x() = R;
100  cylindricalCoors.y() = theta;
101  cylindricalCoors.z() = cartesianPoints[pI].z();
102  localSystemCoordinates_[pI] = cylindricalCoors;
103  }
104 
105  pointVectorField cylindricalCoors
106  (
107  IOobject
108  (
109  "cylindricalCoors" + name_,
110  mesh_.time().timeName(),
111  mesh_,
114  ),
115  pointMesh::New(mesh_),
116  vector::zero
117  );
118  cylindricalCoors.primitiveFieldRef() = localSystemCoordinates_;
119  cylindricalCoors.write();
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124 
125 Foam::NURBS3DVolumeCylindrical::NURBS3DVolumeCylindrical
126 (
127  const dictionary& dict,
128  const fvMesh& mesh,
129  bool computeParamCoors
130 )
131 :
132  NURBS3DVolume(dict, mesh, computeParamCoors),
133  origin_(dict.get<vector>("origin"))
134 {
136  writeCps("cpsBsplines" + mesh_.time().timeName());
137  if (computeParamCoors)
138  {
141  }
142 }
143 
144 
145 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
void updateLocalCoordinateSystem(const vectorField &cartesianPoints)
Update coordinates in the local system based on the cartesian points.
vector transformPointToCartesian(const vector &localCoordinates) const
Transform a point from its coordinate system to a cartesian system.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tensor transformationTensorDxDb(label globalPointIndex)
Transformation tensor for dxdb, from local coordinate system to cartesian.
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const fvMesh & mesh() const
Get mesh.
dimensionedScalar sqrt(const dimensionedScalar &ds)
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
vector origin_
Translation vector.
const pointVectorField & getParametricCoordinates()
Get parametric coordinates.
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:69
Macros for easy insertion into run-time selection tables.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void determineActiveDesignVariablesAndPoints()
Create lists with active design variables and control points.
const fvMesh & mesh_
Definition: NURBS3DVolume.H:79
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar sin(const dimensionedScalar &ds)
void writeCps(const fileName &="cpsFile", const bool transform=true) const
Write control points on a cartesian coordinates system for visualization.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
defineTypeNameAndDebug(combustionModel, 0)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
vector point
Point is a vector.
Definition: point.H:37
#define R(A, B, C, D, E, F, K, M)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Nothing to be read.
Tensor of scalars, i.e. Tensor<scalar>.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127