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-2019 PCOpt/NTUA
9  Copyright (C) 2013-2019 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_),
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  {
140  }
141 }
142 
143 
144 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::pointVectorField
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
Definition: pointFieldsFwd.H:57
Foam::Tensor
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: complexI.H:268
Foam::polyMesh::points
virtual const pointField & points() const
Definition: polyMesh.C:1062
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
pointPatchField.H
Foam::NURBS3DVolumeCylindrical::transformPointToCartesian
vector transformPointToCartesian(const vector &localCoordinates) const
Definition: NURBS3DVolumeCylindrical.C:47
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:257
Foam::NURBS3DVolume
NURBS3DVolume morpher. Includes support functions for gradient computations Base class providing supp...
Definition: NURBS3DVolume.H:68
Foam::atan2
dimensionedScalar atan2(const dimensionedScalar &x, const dimensionedScalar &y)
Definition: dimensionedScalar.C:305
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Definition: MeshObject.C:41
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
pointPatchFieldsFwd.H
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::NURBS3DVolume::getParametricCoordinates
const pointVectorField & getParametricCoordinates()
Definition: NURBS3DVolume.C:1597
R
#define R(A, B, C, D, E, F, K, M)
Foam::NURBS3DVolume::mesh_
const fvMesh & mesh_
Definition: NURBS3DVolume.H:76
Foam::Field
Generic templated field type.
Definition: Field.H:59
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:66
Foam::NURBS3DVolumeCylindrical::updateLocalCoordinateSystem
void updateLocalCoordinateSystem(const vectorField &cartesianPoints)
Definition: NURBS3DVolumeCylindrical.C:81
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::NURBS3DVolumeCylindrical::origin_
vector origin_
Definition: NURBS3DVolumeCylindrical.H:58
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:81
Foam
Definition: atmBoundaryLayer.C:26
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:78
Foam::NURBS3DVolumeCylindrical::transformationTensorDxDb
tensor transformationTensorDxDb(label globalPointIndex)
Definition: NURBS3DVolumeCylindrical.C:64
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Definition: GeometricField.C:759
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Definition: DimensionedFieldReuseFunctions.H:100
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:137
NURBS3DVolumeCylindrical.H
Foam::NURBS3DVolume::writeCps
void writeCps(const fileName &="cpsFile", const bool transform=true) const
Definition: NURBS3DVolume.C:1773
Foam::PtrListOps::get
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Form zero
Definition: VectorSpace.H:111
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:72
Foam::NURBS3DVolume::mesh
const fvMesh & mesh() const
Definition: NURBS3DVolumeI.H:93
Foam::point
vector point
Point is a vector.
Definition: point.H:37
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:49
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pointMesh.H
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:258