shapeSensitivitiesBaseTemplates.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-------------------------------------------------------------------------------
12License
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
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36
37// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38
39template<class Type>
41(
42 const autoPtr
43 <
45 >& sensFieldPtr,
46 const word& name
47) const
48{
50 (
52 (
53 name,
58 ),
61 );
62
63 for (const label patchI : sensitivityPatchIDs_)
64 {
65 volSensField.boundaryFieldRef()[patchI] = sensFieldPtr()[patchI];
66 }
67
68 volSensField.write();
69}
70
71
72template<class Type>
74(
75 const autoPtr<List<Field<Type>>>& sensFieldPtr,
76 const word& name
77) const
78{
80 (
82 (
83 name,
88 ),
91 //fixedValuePointPatchField<Type>::typeName
92 );
93
94 for (const label patchI : sensitivityPatchIDs_)
95 {
96 // Paraview does not visualise values on fixedValuePointPatchFields
97 // Only option is to store in internalField
98 //pointSensField.boundaryFieldRef()[patchI] == sensFieldPtr()[patchI];
99 pointSensField.boundaryField()[patchI].setInInternalField
100 (
101 pointSensField.primitiveFieldRef(),
102 sensFieldPtr()[patchI]
103 );
104 }
105
106 pointSensField.write();
107}
108
109
110template<class Type>
113(
114 const autoPtr
115 <
117 >& sensFieldPtr,
118 const word& name
119) const
120{
122 (
124 (
126 (
127 name,
132 ),
135 )
136 );
138 tVolSensField.ref();
139
141 volSensFieldbf = volSensField.boundaryFieldRef();
142
144 {
145 const label patchI = sensitivityPatchIDs_[pI];
146 volSensFieldbf[patchI] = sensFieldPtr()[patchI];
147 }
148
149 return tVolSensField;
150}
151
152
153// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154
155} // End namespace Foam
156
157// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
Generic GeometricBoundaryField class.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Generic dimensioned Type class.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
virtual bool write(const bool valid=true) const
Write using setting from DB.
tmp< GeometricField< Type, fvPatchField, volMesh > > constructVolSensitivtyField(const autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > &sensFieldPtr, const word &name) const
Constructs volField based on boundaryField and writes it.
void constructAndWriteSensitivityField(const autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > &sensFieldPtr, const word &name) const
Constructs volField based on boundaryField and writes it.
void constructAndWriteSensitivtyPointField(const autoPtr< List< Field< Type > > > &sensFieldPtr, const word &name) const
Constructs pointField based on boundaryField and writes it.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59