Loading...
Searching...
No Matches
ITHACAassign.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/>.
24Namespace
25 ITHACAutilites
26Description
27 Utilities to assign fields and BCs of OF fields
28SourceFiles
29 ITHACAassign.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef ITHACAassign_H
38#define ITHACAassign_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#include "ITHACAerror.H"
58#include "ITHACAparameters.H"
59
61namespace ITHACAutilities
62{
63//-----------------------------------------------------------------------------
76template<class TypeField>
77PtrList<TypeField> averageSubtract(PtrList<TypeField>
78 fields, Eigen::MatrixXd ind, PtrList<TypeField>& ave);
79
80//-----------------------------------------------------------------------------
89template<class TypeField>
90TypeField computeAverage(PtrList<TypeField>& fields);
91
92//--------------------------------------------------------------------------
100template<typename Type>
101void assignIF(GeometricField<Type, fvPatchField, volMesh>& field,
102 Type value);
103
104//--------------------------------------------------------------------------
113template<typename Type>
114void assignIF(GeometricField<Type, fvPatchField, volMesh>& field,
115 Type& value, List<label>& indices);
116
117//--------------------------------------------------------------------------
126template<typename Type>
127void assignIF(GeometricField<Type, fvPatchField, volMesh>& field,
128 Type& value, label index);
129
130//--------------------------------------------------------------------------
137void assignONE(volScalarField& s, List<label>& L);
138
139
140//--------------------------------------------------------------------------
147void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
148 double value);
149
150//--------------------------------------------------------------------------
157void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
158 Eigen::MatrixXd valueVec);
159
160//--------------------------------------------------------------------------
167void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
168 List<double> valueList);
169
170//--------------------------------------------------------------------------
177void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
178 vector value);
179
180//--------------------------------------------------------------------------
187void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
188 tensor value);
189
190//--------------------------------------------------------------------------
197void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
198 Eigen::MatrixXd valueVec);
199
200//--------------------------------------------------------------------------
207void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
208 Eigen::MatrixXd valueVec);
209
210//--------------------------------------------------------------------------
217void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
218 List<vector> valueList);
219
220//--------------------------------------------------------------------------
227void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
228 List<tensor> valueList);
229
230//------------------------------------------------------------------------------
239void assignBC(GeometricField<scalar, fvsPatchField, surfaceMesh>& field,
240 label BC_ind, Eigen::MatrixXd value);
241
242void assignBC(GeometricField<vector, fvsPatchField, surfaceMesh>& field,
243 label BC_ind, Eigen::MatrixXd value);
244
245//------------------------------------------------------------------------------
254template<typename Type>
255void assignBC(GeometricField<Type, fvsPatchField, surfaceMesh>& field,
256 label BC_ind, List<Type>& value);
257
258
259
260//------------------------------------------------------------------------------
269template<typename Type>
270void assignBC(GeometricField<Type, fvsPatchField, surfaceMesh>& field,
271 label BC_ind, Type& value);
272
273
274//--------------------------------------------------------------------------
282template<typename Type>
283void changeNeumann2Dirichlet(GeometricField<Type, fvPatchField, volMesh>& field,
284 Type& value);
285
286
287//--------------------------------------------------------------------------
296template<typename Type>
297void assignZeroDirichlet(GeometricField<Type, fvPatchField,
298 volMesh>& field);
299
300#pragma GCC diagnostic push
301#pragma GCC diagnostic ignored "-Wcomment"
302//------------------------------------------------------------------------------
323template<typename Type>
324void setBoxToValue(GeometricField<Type, fvPatchField, volMesh>& field,
325 Eigen::MatrixXd Box, Type value);
326#pragma GCC diagnostic pop
327
328//--------------------------------------------------------------------------
337template<class Type>
338void changeBCtype(GeometricField<Type, fvPatchField, volMesh>&
339 field, word BCtype, label BC_ind);
340
341
342//--------------------------------------------------------------------------
352template<typename Type>
353void setIndices2Value(labelList& ind2set, List<Type>& value2set,
354 labelList& movingIDS, List<Type>& originalList);
355
356//--------------------------------------------------------------------------
367template<typename Type>
368void assignMixedBC(GeometricField<Type, fvPatchField, volMesh>& field,
369 label BC_ind, Eigen::MatrixXd& value,
370 Eigen::MatrixXd& grad, Eigen::MatrixXd& valueFrac);
371
372
373//--------------------------------------------------------------------------
384template<typename Type>
385void assignMixedBC(GeometricField<Type, fvPatchField, volMesh>& field,
386 label BC_ind, List<Type>& value,
387 List<Type>& grad, List<scalar>& valueFrac);
388
389
390//--------------------------------------------------------------------------
397template<typename Type>
398void normalizeFields(
399 PtrList<GeometricField<Type, fvPatchField, volMesh>>& fields);
400
401//------------------------------------------------------------------------------
412template<typename Type>
413Eigen::MatrixXd getValues(GeometricField<Type, fvPatchField,
414 volMesh>& field, labelList& indices, labelList* xyz = NULL);
415
416//------------------------------------------------------------------------------
427template<typename Type>
428Eigen::MatrixXd getValues(PtrList<GeometricField<Type, fvPatchField,
429 volMesh>>& field, labelList& indices, labelList* xyz = NULL);
430
431}
432
433#endif
Header file of the Foam2Eigen class.
Header file of the ITHACAerror file.
Simple header and source file of the Color::Modifier class to change color to the output stream.
Namespace to implement some useful assign operation of OF fields.
void assignMixedBC(GeometricField< Type, fvPatchField, volMesh > &field, label BC_ind, List< Type > &value, List< Type > &grad, List< scalar > &valueFrac)
Assign value of a boundary condition of type "mixed".
void assignIF(GeometricField< Type, fvPatchField, volMesh > &s, Type value)
Assign internal field.
void changeNeumann2Dirichlet(GeometricField< Type, fvPatchField, volMesh > &field, Type &value)
Change all Neumann boundary conditions to Dirichlet boundary conditions.
TypeField computeAverage(PtrList< TypeField > &fields)
Calculates the average of a list of fields.
void normalizeFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &fields)
Normalize list of Geometric fields.
void assignONE(volScalarField &s, List< label > &L)
Assign one to volScalarField.
void setIndices2Value(labelList &ind2set, List< Type > &value2set, labelList &movingIDS, List< Type > &originalList)
Sets some given Indices of a list of objects to given values.
PtrList< TypeField > averageSubtract(PtrList< TypeField > fields, Eigen::MatrixXd ind, PtrList< TypeField > &ave)
A function to compute time-averaged fields for a set of different parameter samples and also the fiel...
void assignZeroDirichlet(GeometricField< Type, fvPatchField, volMesh > &field)
Assign zero internal field.
void assignBC(GeometricField< scalar, fvPatchField, volMesh > &s, label BC_ind, double value)
Assign uniform Boundary Condition to a volScalarField.
Eigen::MatrixXd getValues(GeometricField< Type, fvPatchField, volMesh > &field, labelList &indices)
void changeBCtype(GeometricField< Type, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
Change the boundary condition type for a GeometricField.
void setBoxToValue(GeometricField< Type, fvPatchField, volMesh > &field, Eigen::MatrixXd Box, Type value)
Set value of a volScalarField to a constant inside a given box.