Loading...
Searching...
No Matches
ITHACAassign.C
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/>.
24\*---------------------------------------------------------------------------*/
25#include "ITHACAassign.H"
26
28{
29
30template<class TypeField>
31PtrList<TypeField> averageSubtract(PtrList<TypeField>
32 fields, Eigen::MatrixXd ind, PtrList<TypeField>& ave)
33{
34 PtrList<TypeField> aveSubtracted;
35 Eigen::VectorXd newInd;
36 newInd.resize(ind.size() + 1);
37 newInd.head(ind.size()) = ind;
38 newInd(ind.size()) = fields.size();
39
40 for (label i = 0; i < ind.size(); i++)
41 {
42 TypeField aveTemp("nut", fields[0] * 0);
43
44 for (label j = newInd(i); j < newInd(i + 1); j++)
45 {
46 aveTemp += fields[j];
47 }
48
49 aveTemp /= newInd(i + 1) - newInd(i);
50 ave.append(aveTemp.clone());
51 }
52
53 for (label i = 0; i < ind.size(); i++)
54 {
55 for (label j = newInd(i); j < newInd(i + 1); j++)
56 {
57 TypeField newfield("nut", fields[0] * 0);
58 newfield = fields[j] - ave[i];
59 aveSubtracted.append(newfield.clone());
60 }
61 }
62
63 return aveSubtracted;
64}
65
66template PtrList<volScalarField> averageSubtract(
67 PtrList<volScalarField>
68 fields, Eigen::MatrixXd ind, PtrList<volScalarField>& ave);
69template PtrList<volVectorField> averageSubtract(
70 PtrList<volVectorField>
71 fields, Eigen::MatrixXd ind, PtrList<volVectorField>& ave);
72
73template<class TypeField>
74TypeField computeAverage(PtrList<TypeField>& fields)
75{
76 TypeField av(fields[0]);
77
78 for (label i = 1; i < fields.size(); i++)
79 {
80 av += fields[i];
81 }
82
83 av = av / fields.size();
84 return av;
85}
86
87template volVectorField computeAverage(
88 PtrList<volVectorField>& fields);
89template volScalarField computeAverage(
90 PtrList<volScalarField>& fields);
91
92
93template<typename Type>
94void assignIF(GeometricField<Type, fvPatchField, volMesh>& s,
95 Type value)
96{
97 for (label i = 0; i < s.internalField().size(); i++)
98 {
99 s.ref()[i] = value;
100 }
101}
102
103template void assignIF(
104 GeometricField<scalar, fvPatchField, volMesh>& field, scalar value);
105template void assignIF(
106 GeometricField<vector, fvPatchField, volMesh>& field, vector value);
107
108template<typename Type>
109void assignIF(GeometricField<Type, fvPatchField, volMesh>& s,
110 Type& value, List<label>& indices)
111{
112 for (label i = 0; i < indices.size(); i++)
113 {
114 s.ref()[indices[i]] = value;
115 }
116}
117
118template void assignIF(GeometricField<scalar, fvPatchField, volMesh>& s,
119 scalar& value, List<label>& indices);
120template void assignIF(GeometricField<vector, fvPatchField, volMesh>& s,
121 vector& value, List<label>& indices);
122
123template<typename Type>
124void assignIF(GeometricField<Type, fvPatchField, volMesh>& s,
125 Type& value, label index)
126{
127 s.ref()[index] = value;
128}
129
130template void assignIF(GeometricField<scalar, fvPatchField, volMesh>& field,
131 scalar& value, label index);
132template void assignIF(GeometricField<vector, fvPatchField, volMesh>& field,
133 vector& value, label index);
134
135void assignONE(volScalarField& s, List<label>& L)
136{
137 for (label i = 0; i < L.size(); i++)
138 {
139 s.ref()[L[i]] = 1;
140 }
141}
142
143void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
144 double value)
145{
146 label sizeBC = s.boundaryField()[BC_ind].size();
147 List<double> valueList(sizeBC);
148
149 for (label i = 0; i < sizeBC; i++)
150 {
151 valueList[i] = value;
152 }
153
154 assignBC(s, BC_ind, valueList);
155}
156
157void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
158 Eigen::MatrixXd valueVec)
159{
160 label sizeBC = s.boundaryField()[BC_ind].size();
161 M_Assert(sizeBC == valueVec.size(),
162 "The size of the given values matrix has to be equal to the dimension of the boundaryField");
163 List<double> valueList(sizeBC);
164
165 for (label i = 0; i < sizeBC; i++)
166 {
167 valueList[i] = valueVec(i);
168 }
169
170 assignBC(s, BC_ind, valueList);
171}
172
173// Assign a BC for a scalar field
174void assignBC(GeometricField<scalar, fvPatchField, volMesh>& s, label BC_ind,
175 List<double> valueList)
176{
177 word typeBC = s.boundaryField()[BC_ind].type();
178 label sizeBC = s.boundaryField()[BC_ind].size();
179 M_Assert(sizeBC == valueList.size(),
180 "The size of the given values list has to be equal to the dimension of the boundaryField");
182
183 if (typeBC == "fixedGradient")
184 {
185 fixedGradientFvPatchScalarField& Tpatch =
186 refCast<fixedGradientFvPatchScalarField>(s.boundaryFieldRef()[BC_ind]);
187 scalarField& gradTpatch = Tpatch.gradient();
188 forAll(gradTpatch, faceI)
189 {
190 double value = valueList[faceI];
191 gradTpatch[faceI] = value;
192 }
193 }
194 else if (typeBC == "freestream")
195 {
196 for (label i = 0; i < sizeBC; i++)
197 {
198 double value = valueList[i];
199 s.boundaryFieldRef()[BC_ind][i] = value;
200 }
201
202 freestreamFvPatchField<scalar>& Tpatch =
203 refCast<freestreamFvPatchField<scalar >> (s.boundaryFieldRef()[BC_ind]);
204 scalarField& gradTpatch = Tpatch.freestreamValue();
205 forAll(gradTpatch, faceI)
206 {
207 double value = valueList[faceI];
208 gradTpatch[faceI] = value;
209 }
210 }
211
212 else if (typeBC == "empty" || typeBC == "zeroGradient")
213 {}
214
215 else
216 {
217 try
218 {
219 if (typeBC != "fixedGradient" && typeBC != "freestream" && typeBC != "empty"
220 && typeBC != "zeroGradient" && typeBC != "fixedValue" && typeBC != "calculated"
221 && typeBC != "fixedFluxPressure" && typeBC != "processor"
222 && typeBC != "nutkWallFunction" && typeBC != "mixedEnergy")
223 {
224 word message = "Pay attention, your typeBC " + typeBC + " for " + s.name() +
225 " is not included into the developed ones. Your BC will be treated as a classical fixedValue.";
226 throw (message);
227 }
228 }
229 catch (const word message)
230 {
231 if (para->warnings)
232 {
233 WarningInFunction << message << endl;
234 }
235 }
236
237 for (label i = 0; i < sizeBC; i++)
238 {
239 double value = valueList[i];
240 s.boundaryFieldRef()[BC_ind][i] = value;
241 }
242 }
243}
244
245void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
246 vector value)
247{
248 M_Assert(value.size() == 3,
249 "The size of the given vector has to be equal to 3 for the 3 components");
250 label sizeBC = s.boundaryField()[BC_ind].size();
251 List<vector> valueList(sizeBC);
252
253 for (label i = 0; i < sizeBC; i++)
254 {
255 valueList[i] = value;
256 }
257
258 assignBC(s, BC_ind, valueList);
259}
260
261void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
262 tensor value)
263{
264 M_Assert(value.size() == 9,
265 "The size of the given vector has to be equal to 3x3 for the 3x3 components");
266 label sizeBC = s.boundaryField()[BC_ind].size();
267 List<tensor> valueList(sizeBC);
268
269 for (label i = 0; i < sizeBC; i++)
270 {
271 valueList[i] = value;
272 }
273
274 assignBC(s, BC_ind, valueList);
275}
276
277void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
278 Eigen::MatrixXd valueVec)
279{
280 label sizeBC = s.boundaryField()[BC_ind].size();
281 M_Assert(sizeBC * 3 == valueVec.size(),
282 "The size of the given values matrix has to be equal to 3 times the dimension of the boundaryField");
283 List<vector> valueList(sizeBC);
284
285 for (label i = 0; i < sizeBC; i++)
286 {
287 for (label j = 0; j < 3; j++)
288 {
289 valueList[i].component(j) = valueVec(i + sizeBC * j);
290 }
291 }
292
293 assignBC(s, BC_ind, valueList);
294}
295
296void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
297 Eigen::MatrixXd valueVec)
298{
299 label sizeBC = s.boundaryField()[BC_ind].size();
300 M_Assert(sizeBC * 9 == valueVec.size(),
301 "The size of the given values matrix has to be equal to 9 times the dimension of the boundaryField");
302 List<tensor> valueList(sizeBC);
303
304 for (label i = 0; i < sizeBC; i++)
305 {
306 for (label j = 0; j < 9; j++)
307 {
308 valueList[i].component(j) = valueVec(i + sizeBC * j);
309 }
310 }
311
312 assignBC(s, BC_ind, valueList);
313}
314
315void assignBC(GeometricField<scalar, fvsPatchField, surfaceMesh>& s,
316 label BC_ind, Eigen::MatrixXd valueVec)
317{
318 label sizeBC = s.boundaryField()[BC_ind].size();
319 M_Assert(sizeBC == valueVec.rows() && valueVec.cols() == 1,
320 "The given matrix must be a column one with the size equal to 3 times the dimension of the boundaryField");
321 List<scalar> valueList(sizeBC);
322
323 for (label i = 0; i < sizeBC; i++)
324 {
325 valueList[i] = valueVec(i);
326 }
327
328 assignBC(s, BC_ind, valueList);
329}
330
331void assignBC(GeometricField<vector, fvsPatchField, surfaceMesh>& s,
332 label BC_ind, Eigen::MatrixXd valueVec)
333{
334 label sizeBC = s.boundaryField()[BC_ind].size();
335 M_Assert(sizeBC * 3 == valueVec.rows() && valueVec.cols() == 1,
336 "The given matrix must be a column one with the size equal to the dimension of the boundaryField");
337 List<vector> valueList(sizeBC);
338
339 for (label i = 0; i < sizeBC; i++)
340 {
341 valueList[i].component(0) = valueVec(i);
342 valueList[i].component(1) = valueVec(i + sizeBC);
343 valueList[i].component(2) = valueVec(i + sizeBC * 2);
344 }
345
346 assignBC(s, BC_ind, valueList);
347}
348
349// Assign a BC for a vector field
350void assignBC(GeometricField<vector, fvPatchField, volMesh>& s, label BC_ind,
351 List<vector> valueList)
352{
353 word typeBC = s.boundaryField()[BC_ind].type();
354 label sizeBC = s.boundaryField()[BC_ind].size();
355 M_Assert(sizeBC == valueList.size(),
356 "The size of the given values list has to be equal to the dimension of the boundaryField");
358
359 if (s.boundaryField()[BC_ind].type() == "fixedGradient")
360 {
361 Info << "This Feature is not implemented for this boundary condition" << endl;
362 exit(0);
363 }
364 else if (s.boundaryField()[BC_ind].type() == "freestream")
365 {
366 for (label i = 0; i < s.boundaryField()[BC_ind].size(); i++)
367 {
368 s.boundaryFieldRef()[BC_ind][i] = valueList[i];
369 }
370
371 freestreamFvPatchField<vector>& Tpatch =
372 refCast<freestreamFvPatchField<vector >> (s.boundaryFieldRef()[BC_ind]);
373 vectorField& gradTpatch = Tpatch.freestreamValue();
374 forAll(gradTpatch, faceI)
375 {
376 gradTpatch[faceI] = valueList[faceI];
377 }
378 }
379
380 else if (s.boundaryField()[BC_ind].type() == "empty"
381 || s.boundaryField()[BC_ind].type() == "zeroGradient")
382 {}
383
384 else
385 {
386 try
387 {
388 if (typeBC != "fixedGradient" && typeBC != "freestream" && typeBC != "empty"
389 && typeBC != "zeroGradient" && typeBC != "fixedValue" && typeBC != "calculated"
390 && typeBC != "processor")
391 {
392 word message = "Pay attention, your typeBC " + typeBC + " for " + s.name() +
393 " is not included into the developed ones. Your BC will be treated as a classical fixedValue.";
394 throw (message);
395 }
396 }
397 catch (const word message)
398 {
399 if (para->warnings)
400 {
401 WarningInFunction << message << endl;
402 }
403 }
404
405 for (label i = 0; i < sizeBC; i++)
406 {
407 s.boundaryFieldRef()[BC_ind][i] = valueList[i];
408 }
409 }
410}
411
412
413// Assign a BC for a tensor field
414void assignBC(GeometricField<tensor, fvPatchField, volMesh>& s, label BC_ind,
415 List<tensor> valueList)
416{
417 word typeBC = s.boundaryField()[BC_ind].type();
418 label sizeBC = s.boundaryField()[BC_ind].size();
419 M_Assert(sizeBC == valueList.size(),
420 "The size of the given values list has to be equal to the dimension of the boundaryField");
422
423 if (s.boundaryField()[BC_ind].type() == "fixedGradient")
424 {
425 Info << "This Feature is not implemented for this boundary condition" << endl;
426 exit(0);
427 }
428 else if (s.boundaryField()[BC_ind].type() == "freestream")
429 {
430 for (label i = 0; i < s.boundaryField()[BC_ind].size(); i++)
431 {
432 s.boundaryFieldRef()[BC_ind][i] = valueList[i];
433 }
434
435 freestreamFvPatchField<tensor>& Tpatch =
436 refCast<freestreamFvPatchField<tensor >> (s.boundaryFieldRef()[BC_ind]);
437 tensorField& gradTpatch = Tpatch.freestreamValue();
438 forAll(gradTpatch, faceI)
439 {
440 gradTpatch[faceI] = valueList[faceI];
441 }
442 }
443
444 else if (s.boundaryField()[BC_ind].type() == "empty"
445 || s.boundaryField()[BC_ind].type() == "zeroGradient")
446 {}
447
448 else
449 {
450 try
451 {
452 if (typeBC != "fixedGradient" && typeBC != "freestream" && typeBC != "empty"
453 && typeBC != "zeroGradient" && typeBC != "fixedValue" && typeBC != "calculated"
454 && typeBC != "processor")
455 {
456 word message = "Pay attention, your typeBC " + typeBC + " for " + s.name() +
457 " is not included into the developed ones. Your BC will be treated as a classical fixedValue.";
458 throw (message);
459 }
460 }
461 catch (const word message)
462 {
463 if (para->warnings)
464 {
465 WarningInFunction << message << endl;
466 }
467 }
468
469 for (label i = 0; i < sizeBC; i++)
470 {
471 s.boundaryFieldRef()[BC_ind][i] = valueList[i];
472 }
473 }
474}
475
476
477template<typename Type>
478void assignBC(GeometricField<Type, fvsPatchField, surfaceMesh>& s,
479 label BC_ind,
480 List<Type>& valueList)
481{
482 word typeBC = s.boundaryField()[BC_ind].type();
483 label sizeBC = s.boundaryField()[BC_ind].size();
484 M_Assert(sizeBC == valueList.size(),
485 "The size of the given values list has to be equal to the dimension of the boundaryField");
487
488 if (s.boundaryField()[BC_ind].type() == "fixedGradient")
489 {
490 fixedGradientFvPatchField<Type>& Tpatch =
491 refCast<fixedGradientFvPatchField<Type >> (s.boundaryFieldRef()[BC_ind]);
492 Field<Type>& gradTpatch = Tpatch.gradient();
493 forAll(gradTpatch, faceI)
494 {
495 gradTpatch[faceI] = valueList[faceI];
496 }
497 }
498 else if (s.boundaryField()[BC_ind].type() == "freestream")
499 {
500 for (label i = 0; i < s.boundaryField()[BC_ind].size(); i++)
501 {
502 s.boundaryFieldRef()[BC_ind][i] = valueList[i];
503 }
504
505 freestreamFvPatchField<Type>& Tpatch =
506 refCast<freestreamFvPatchField<Type >> (s.boundaryFieldRef()[BC_ind]);
507 Field<Type>& gradTpatch = Tpatch.freestreamValue();
508 forAll(gradTpatch, faceI)
509 {
510 gradTpatch[faceI] = valueList[faceI];
511 }
512 }
513 else if (s.boundaryField()[BC_ind].type() == "empty"
514 || s.boundaryField()[BC_ind].type() == "zeroGradient")
515 {}
516 else
517 {
518 try
519 {
520 if (typeBC != "fixedGradient" && typeBC != "freestream" && typeBC != "empty"
521 || typeBC != "zeroGradient" && typeBC != "fixedValue" && typeBC != "calculated"
522 && typeBC != "fixedFluxPressure" && typeBC != "processor")
523 {
524 word message = "Pay attention, your typeBC " + typeBC + " for " + s.name() +
525 " is not included into the developed ones. Your BC will be treated as a classical fixedValue.";
526 throw (message);
527 }
528 }
529 catch (const word message)
530 {
531 if (para->warnings)
532 {
533 WarningInFunction << message << endl;
534 }
535 }
536
537 for (label i = 0; i < s.boundaryField()[BC_ind].size(); i++)
538 {
539 s.boundaryFieldRef()[BC_ind][i] = valueList[i];
540 }
541 }
542}
543
544template void assignBC(
545 GeometricField<scalar, fvsPatchField, surfaceMesh>& s, label BC_ind,
546 List<scalar>& valueList);
547template void assignBC(
548 GeometricField<vector, fvsPatchField, surfaceMesh>& s, label BC_ind,
549 List<vector>& valueList);
550
551template<typename Type>
552void assignBC(GeometricField<Type, fvsPatchField, surfaceMesh>& s,
553 label BC_ind,
554 Type& value)
555{
556 label sizeBC = s.boundaryField()[BC_ind].size();
557 List<Type> valueList(sizeBC);
558
559 for (label i = 0; i < sizeBC; i++)
560 {
561 valueList[i] = value;
562 }
563
564 assignBC(s, BC_ind, valueList);
565}
566
567template void assignBC(
568 GeometricField<scalar, fvsPatchField, surfaceMesh>& s, label BC_ind,
569 scalar& valueList);
570template void assignBC(
571 GeometricField<vector, fvsPatchField, surfaceMesh>& s, label BC_ind,
572 vector& valueList);
573
574template<typename Type>
575void changeNeumann2Dirichlet(GeometricField<Type, fvPatchField, volMesh>&
576 field,
577 Type& value)
578{
579 forAll(field.mesh().boundary(), i)
580 {
581 if (field.boundaryField()[i].type() == "zeroGradient" ||
582 field.boundaryField()[i].type() == "fixedGradient")
583 {
584 ITHACAutilities::changeBCtype(field, "fixedValue", i);
585 assignBC(field, i, value);
586 }
587 }
588}
589
591 GeometricField<scalar, fvPatchField, volMesh>& field, scalar& value);
593 GeometricField<vector, fvPatchField, volMesh>& field, vector& value);
594
595template<typename Type>
596void assignZeroDirichlet(GeometricField<Type, fvPatchField, volMesh>& field)
597{
598 Type v;
599 v = v * 0;
600 assignIF(field, v);
601 forAll(field.mesh().boundary(), i)
602 {
603 if (field.boundaryField()[i].type() == "fixedValue")
604 {
605 assignBC(field, i, v);
606 }
607 }
608
610}
611
613 GeometricField<vector, fvPatchField, volMesh>& field);
615 GeometricField<scalar, fvPatchField, volMesh>& field);
616
617template<typename Type>
618void setBoxToValue(GeometricField<Type, fvPatchField, volMesh>& field,
619 Eigen::MatrixXd Box, Type value)
620{
621 M_Assert(Box.rows() == 2
622 && Box.cols() == 3,
623 "The box must be a 2*3 matrix shaped in this way: \nBox = \t|x0, y0, z0|\n\t|x1, yi, z1|\n");
624
625 for (label i = 0; i < field.internalField().size(); i++)
626 {
627 auto cx = field.mesh().C()[i].component(vector::X);
628 auto cy = field.mesh().C()[i].component(vector::Y);
629 auto cz = field.mesh().C()[i].component(vector::Z);
630
631 if (cx >= Box(0, 0) && cy >= Box(0, 1) && cz >= Box(0, 2) && cx <= Box(1, 0)
632 && cy <= Box(1, 1) && cz <= Box(1, 2) )
633 {
634 field.ref()[i] = value;
635 }
636 }
637
638 for (label i = 0; i < field.boundaryField().size(); i++)
639 {
640 for (label j = 0; j < field.boundaryField()[i].size(); j++)
641 {
642 if (field.boundaryField()[i].type() == "fixedValue"
643 || field.boundaryField()[i].type() == "calculated")
644 {
645 auto cx = field.mesh().C().boundaryField()[i][j][0];
646 auto cy = field.mesh().C().boundaryField()[i][j][1];
647 auto cz = field.mesh().C().boundaryField()[i][j][2];
648
649 if (cx >= Box(0, 0) && cy >= Box(0, 1) && cz >= Box(0, 2) && cx <= Box(1, 0)
650 && cy <= Box(1, 1) && cz <= Box(1, 2) )
651 {
652 field.boundaryFieldRef()[i][j] = value;
653 }
654 }
655 }
656 }
657}
658
659template void setBoxToValue(GeometricField<scalar, fvPatchField, volMesh>&
660 field, Eigen::MatrixXd Box, scalar value);
661template void setBoxToValue(GeometricField<vector, fvPatchField, volMesh>&
662 field, Eigen::MatrixXd Box, vector value);
663
664template<typename Type>
665void setIndices2Value(labelList& ind2set, List<Type>& value2set,
666 labelList& movingIDS, List<Type>& originalList)
667{
668 M_Assert(ind2set.size() == value2set.size(),
669 "The size of the indices must be equal to the size of the values list");
670 M_Assert(originalList.size() >= value2set.size(),
671 "The size of the original list of values must be bigger than the size of the list of values you want to set");
672 labelList ind_ok(ind2set);
673
674 for (label i = 0; i < ind2set.size(); i++)
675 {
676 for (label k = 0; k < movingIDS.size(); k++)
677 {
678 if (ind2set[i] == movingIDS[k])
679 {
680 ind_ok[i] = k;
681 break;
682 }
683 }
684 }
685
686 for (label i = 0; i < ind2set.size(); i++)
687 {
688 originalList[ind_ok[i]] = value2set[i];
689 }
690}
691
692template void setIndices2Value(labelList& ind2set, List<scalar>& value2set,
693 labelList& movingIDS, List<scalar>& originalList);
694template void setIndices2Value(labelList& ind2set, List<vector>& value2set,
695 labelList& movingIDS, List<vector>& originalList);
696
697template<class Type>
699 GeometricField<Type, fvPatchField, volMesh>& field, word BCtype,
700 label BC_ind)
701{
702 field.boundaryFieldRef().set(BC_ind, fvPatchField<Type>::New(BCtype,
703 field.mesh().boundary()[BC_ind], field));
704}
705
707(GeometricField<scalar, fvPatchField, volMesh>& field, word BCtype,
708 label BC_ind);
710(GeometricField<vector, fvPatchField, volMesh>& field, word BCtype,
711 label BC_ind);
712
713template<typename Type>
715 GeometricField<Type, fvPatchField, volMesh>& field, label BC_ind,
716 List<Type>& value, List<Type>& grad, List<scalar>& valueFrac)
717{
718 std::string message = "Patch is NOT mixed. It is of type: " +
719 field.boundaryField()[BC_ind].type();
720 M_Assert(field.boundaryField()[BC_ind].type() == "mixed", message.c_str());
721
722 if (field.boundaryField()[BC_ind].type() == "mixed")
723 {
724 mixedFvPatchField<Type>& Tpatch =
725 refCast<mixedFvPatchField<Type >> (field.boundaryFieldRef()[BC_ind]);
726 Field<Type>& valueTpatch = Tpatch.refValue();
727 Field<Type>& gradTpatch = Tpatch.refGrad();
728 Field<scalar>& valueFracTpatch = Tpatch.valueFraction();
729 valueTpatch = value;
730 gradTpatch = grad;
731 valueFracTpatch = valueFrac;
732 }
733}
734
736 GeometricField<scalar, fvPatchField, volMesh>& field, label BC_ind,
737 List<scalar>& value, List<scalar>& grad, List<scalar>& valueFrac);
738
740 GeometricField<vector, fvPatchField, volMesh>& field, label BC_ind,
741 List<vector>& value, List<vector>& grad, List<scalar>& valueFrac);
742
743template<typename Type>
745 PtrList<GeometricField<Type, fvPatchField, volMesh >> & fields)
746{
748 word normType = para->ITHACAdict->lookupOrDefault<word>("normalizationNorm",
749 "L2");
750 M_Assert(normType == "L2" ||
751 normType == "Frobenius", "The normalizationNorm can be only L2 or Frobenius" );
752 Eigen::MatrixXd eigenFields = Foam2Eigen::PtrList2Eigen(fields);
753 List<Eigen::MatrixXd> eigenFieldsBC = Foam2Eigen::PtrList2EigenBC(fields);
754
755 for (label i = 0; i < fields.size(); i++)
756 {
757 double norm;
758
759 if (normType == "L2")
760 {
761 norm = L2Norm(fields[i]);
762 }
763 else if (normType == "Frobenius")
764 {
765 norm = frobNorm(fields[i]);
766 }
767
768 GeometricField<Type, fvPatchField, volMesh> tmp2(fields[0].name(),
769 fields[0] * 0);
770 Eigen::VectorXd vec = eigenFields.col(i) / norm;
771 tmp2 = Foam2Eigen::Eigen2field(tmp2, vec);
772
773 // Adjusting boundary conditions
774 for (label k = 0; k < tmp2.boundaryField().size(); k++)
775 {
776 Eigen::MatrixXd vec = eigenFieldsBC[k].col(i) / norm;
777 assignBC(tmp2, k, vec);
778 }
779
780 fields.set(i, tmp2.clone());
781 }
782}
783
784template void normalizeFields(
785 PtrList<GeometricField<scalar, fvPatchField, volMesh >> & fields);
786template void normalizeFields(
787 PtrList<GeometricField<vector, fvPatchField, volMesh >>& fields);
788
789
790template<typename Type>
791Eigen::MatrixXd getValues(GeometricField<Type, fvPatchField,
792 volMesh>& field, labelList& indices)
793{
794 List<Type> list(indices.size());
795 M_Assert(max(indices) < field.size(),
796 "The list indices are too large respect to field dimension");
797
798 for (label i = 0; i < indices.size(); i++)
799 {
800 list[i] = field[indices[i]];
801 }
802
803 return Foam2Eigen::field2Eigen(list);
804}
805
806template<>
807Eigen::MatrixXd getValues(GeometricField<vector, fvPatchField,
808 volMesh>& field, labelList& indices, labelList* xyz)
809{
810 M_Assert(max(indices) < field.size(),
811 "The list of indices is too large respect to field dimension. There is at least one value larger than the dimension of the list");
812
813 if (xyz != NULL)
814 {
815 List<scalar> list;
816 list.resize(indices.size());
817 M_Assert(max(* xyz) <= 2,
818 "The list of xyz positions contains at list one value larger than 2");
819 labelList l = * xyz;
820
821 for (label i = 0; i < indices.size(); i++)
822 {
823 list[i] = field[indices[i]][l[i]];
824 }
825
826 return Foam2Eigen::field2Eigen(list);
827 }
828 else
829 {
830 List<vector> list;
831 list.resize(indices.size());
832
833 for (label i = 0; i < indices.size(); i++)
834 {
835 list[i] = field[indices[i]];
836 }
837
838 return Foam2Eigen::field2Eigen(list);
839 }
840}
841
842template<>
843Eigen::MatrixXd getValues(GeometricField<scalar, fvPatchField,
844 volMesh>& field, labelList& indices, labelList* xyz)
845{
846 M_Assert(max(indices) < field.size(),
847 "The list of indices is too large respect to field dimension. There is at least one value larger than the dimension of the list");
848 List<scalar> list;
849 list.resize(indices.size());
850
851 for (label i = 0; i < indices.size(); i++)
852 {
853 list[i] = field[indices[i]];
854 }
855
856 return Foam2Eigen::field2Eigen(list);
857}
858
859template<typename T>
860Eigen::MatrixXd getValues(PtrList<GeometricField<T, fvPatchField,
861 volMesh >> & fields, labelList& indices, labelList* xyz)
862{
863 Eigen::MatrixXd out;
864 Eigen::MatrixXd a = getValues(fields[0], indices, xyz);
865 out.resize(a.rows(), fields.size());
866 out.col(0) = a;
867 for (label i = 1; i < fields.size(); i++)
868 {
869 out.col(i) = getValues(fields[i], indices, xyz);
870 }
871
872 return out;
873}
874
875
876template
877Eigen::MatrixXd getValues(PtrList<GeometricField<scalar, fvPatchField,
878 volMesh >>& fields, labelList& indices, labelList* xyz);
879template
880Eigen::MatrixXd getValues(PtrList<GeometricField<vector, fvPatchField,
881 volMesh >> & fields, labelList& indices, labelList* xyz);
882
883
884}
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
#define M_Assert(Expr, Msg)
Header file of the ITHACAassign file.
ITHACAparameters * para
Definition NLsolve.H:40
static GeometricField< scalar, PatchField, GeoMesh > Eigen2field(GeometricField< scalar, PatchField, GeoMesh > &field, Eigen::VectorXd &eigen_vector, bool correctBC=true)
Convert a vector in Eigen format into an OpenFOAM scalar GeometricField.
Definition Foam2Eigen.C:514
static Eigen::VectorXd field2Eigen(GeometricField< tensor, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
Definition Foam2Eigen.C:40
static List< Eigen::MatrixXd > PtrList2EigenBC(PtrList< GeometricField< scalar, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert an OpenFOAM scalar field to a List of Eigen Vectors, one for each boundary.
Definition Foam2Eigen.C:270
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field)
Definition Foam2Eigen.C:646
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance(fvMesh &mesh, Time &localTime)
Gets an instance of ITHACAparameters, to be used if the instance is not existing.
volScalarField & T
Definition createT.H:46
volScalarField & v
Namespace to implement some useful assign operation of OF fields.
template void changeBCtype< vector >(GeometricField< vector, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
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.
double L2Norm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:41
template void changeBCtype< scalar >(GeometricField< scalar, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
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.
double frobNorm(GeometricField< Type, PatchField, GeoMesh > &field)
Evaluate the Frobenius norm of a field.
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)
template void assignMixedBC< scalar >(GeometricField< scalar, fvPatchField, volMesh > &field, label BC_ind, List< scalar > &value, List< scalar > &grad, List< scalar > &valueFrac)
void changeBCtype(GeometricField< Type, fvPatchField, volMesh > &field, word BCtype, label BC_ind)
Change the boundary condition type for a GeometricField.
template void assignMixedBC< vector >(GeometricField< vector, fvPatchField, volMesh > &field, label BC_ind, List< vector > &value, List< vector > &grad, List< scalar > &valueFrac)
void setBoxToValue(GeometricField< Type, fvPatchField, volMesh > &field, Eigen::MatrixXd Box, Type value)
Set value of a volScalarField to a constant inside a given box.
label i
Definition pEqn.H:46