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