displacementMethodvelocityLaplacian.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-2020 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 
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 defineTypeNameAndDebug(displacementMethodvelocityLaplacian, 1);
43 (
44  displacementMethod,
45  displacementMethodvelocityLaplacian,
47 );
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 displacementMethodvelocityLaplacian::displacementMethodvelocityLaplacian
53 (
54  fvMesh& mesh,
55  const labelList& patchIDs
56 )
57 :
59  pointMotionU_
60  (
61  refCast<velocityLaplacianFvMotionSolver>(motionPtr_()).pointMotionU()
62  ),
63  cellMotionU_
64  (
65  refCast<velocityLaplacianFvMotionSolver>(motionPtr_()).cellMotionU()
66  ),
67  /*
68  // Getting a reference from the database is dangerous since multiple
69  // fields with the same name might be registered
70  pointMotionU_
71  (
72  const_cast<pointVectorField&>
73  (
74  mesh_.lookupObject<pointVectorField>("pointMotionU")
75  )
76  ),
77  cellMotionU_
78  (
79  const_cast<volVectorField&>
80  (
81  mesh_.lookupObject<volVectorField>("cellMotionU")
82  )
83  ),
84  */
85  resetFields_
86  (
87  IOdictionary::readContents
88  (
89  IOobject
90  (
91  "dynamicMeshDict",
92  mesh.time().constant(),
93  mesh,
94  IOobject::MUST_READ
95  )
96  ).subDict("velocityLaplacianCoeffs").getOrDefault<bool>
97  (
98  "resetFields",
99  true
100  )
101  )
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 (
109  const pointVectorField& pointMovement
110 )
111 {
112  if (resetFields_)
113  {
114  pointMotionU_.primitiveFieldRef() = vector::zero;
115  cellMotionU_.primitiveFieldRef() = vector::zero;
117  }
118 
119  // Set boundary mesh movement and calculate
120  // max current boundary displacement
121  for (label patchI : patchIDs_)
122  {
123  // Set boundary field. Needed for the motionSolver class
124  pointMotionU_.boundaryFieldRef()[patchI] ==
125  pointMovement.boundaryField()[patchI].patchInternalField()();
126 
127  // Set boundary field values as seen from the internalField!
128  // Needed for determining the max displacement
129  pointMotionU_.boundaryFieldRef()[patchI].setInInternalField
130  (
132  pointMovement.boundaryField()[patchI].patchInternalField()()
133  );
134 
135  // Find max value
137  max
138  (
140  gMax
141  (
142  mag
143  (
144  pointMotionU_.boundaryField()[patchI].
145  patchInternalField()
146  )
147  )
148  );
149  }
150 }
151 
152 
154 (
155  const volVectorField& cellMovement
156 )
157 {
159 
160  /*
161  // Set boundary mesh movement and calculate max current boundary
162  // displacement
163  forAll(patchIDs_, pI)
164  {
165  label patchI = patchIDs_[pI];
166  cellMotionU_.boundaryField()[patchI] ==
167  cellMovement.boundaryField()[patchI];
168 
169  // Find max value
170  maxDisplacement_ =
171  max
172  (
173  maxDisplacement_,
174  gMax
175  (
176  mag(cellMovement.boundaryField()[patchI])
177  )
178  );
179  }
180  */
181 }
182 
183 
185 (
186  const vectorField& controlField
187 )
188 {
190 }
191 
192 
194 (
195  const scalarField& controlField
196 )
197 {
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 } // End namespace Foam
205 
206 // ************************************************************************* //
const labelList patchIDs(pbm.indices(polyPatchNames, true))
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Type & refCast(U &obj)
A dynamic_cast (for references). Generates a FatalError on failed casts and uses the virtual type() m...
Definition: typeInfo.H:159
void setControlField(const vectorField &controlField)
Set control field as a vectorField. For methods working with parameters (RBF etc) ...
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
dynamicFvMesh & mesh
Abstract base class for displacement methods, which are a set or wrapper classes allowing to change t...
Type gMax(const FieldField< Field, Type > &f)
defineTypeNameAndDebug(combustionModel, 0)
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Mesh motion solver for an fvMesh. Based on solving the cell-centre Laplacian for the motion velocity...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
labelList patchIDs_
IDs of the patches to be moved.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
void setMotionField(const pointVectorField &pointMovement)
Set motion filed related to model based on given motion.