Loading...
Searching...
No Matches
ITHACAerror.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-------------------------------------------------------------------------------
12
13License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28\*---------------------------------------------------------------------------*/
29
30#include "ITHACAerror.H"
31
34
35// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36
37namespace ITHACAutilities
38{
39
40template<>
41double L2Norm(GeometricField<scalar, fvPatchField, volMesh>& field)
42{
43 double a;
44 a = Foam::sqrt(fvc::domainIntegrate(field * field).value());
45 return a;
46}
47
48template<>
49double L2Norm(GeometricField<vector, fvPatchField, volMesh>& field)
50{
51 double a;
52 a = Foam::sqrt(fvc::domainIntegrate(field & field).value());
53 return a;
54}
55
56template<>
57double LinfNorm(GeometricField<scalar, fvPatchField, volMesh>& field)
58{
59 double a;
60 a = Foam::max(Foam::sqrt(field.internalField() *
61 field.internalField())).value();
62 return a;
63}
64
65template<>
66double LinfNorm(GeometricField<vector, fvPatchField, volMesh>& field)
67{
68 double a;
69 Info << "LinfNorm(GeometricField<vector, fvPatchField, volMesh>& field) is still to be implemented"
70 << endl;
71 exit(12);
72 return a;
73}
74
75
76
77template<class Type, template<class> class PatchField, class GeoMesh>
78double errorFrobRel(GeometricField<Type, PatchField, GeoMesh>& field1,
79 GeometricField<Type, PatchField, GeoMesh>& field2, List<label>* labels)
80{
81 double err;
82 autoPtr<GeometricField<Type, PatchField, GeoMesh >> errField;
83 autoPtr<GeometricField<Type, PatchField, GeoMesh >> field1_S;
84 autoPtr<GeometricField<Type, PatchField, GeoMesh >> field2_S;
85 autoPtr<fvMeshSubset> submesh;
86
87 if (labels != NULL)
88 {
89 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(field1.mesh()));
90#if OPENFOAM >= 1812
91 submesh->setCellSubset(* labels);
92#else
93 submesh->setLargeCellSubset(* labels);
94#endif
95 GeometricField<Type, PatchField, GeoMesh> field1tmp(submesh->interpolate(
96 field1));
97 GeometricField<Type, PatchField, GeoMesh> field2tmp(submesh->interpolate(
98 field2));
99 field1_S = autoPtr<GeometricField<Type, PatchField, GeoMesh >>
100 (new GeometricField<Type, PatchField, GeoMesh>(field1tmp.clone()));
101 field2_S = autoPtr<GeometricField<Type, PatchField, GeoMesh >>
102 (new GeometricField<Type, PatchField, GeoMesh>(field2tmp.clone()));
103 }
104 else
105 {
106 field1_S = autoPtr<GeometricField<Type, PatchField, GeoMesh >>
107 (new GeometricField<Type, PatchField, GeoMesh>(field1));
108 field2_S = autoPtr<GeometricField<Type, PatchField, GeoMesh >>
109 (new GeometricField<Type, PatchField, GeoMesh>(field2));
110 }
111 errField = autoPtr<GeometricField<Type, PatchField, GeoMesh >>
112 (new GeometricField<Type, PatchField, GeoMesh>(field1_S() - field2_S()));
113
114 if (frobNorm(field1) <= 1e-6)
115 {
116 err = 0;
117 }
118 else
119 {
120 err = frobNorm(errField()) / frobNorm(field1_S());
121 }
122 return err;
123}
124
125
126template double errorFrobRel(GeometricField<scalar, fvPatchField, volMesh>&
127 field1,
128 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
129template double errorFrobRel(GeometricField<vector, fvPatchField, volMesh>&
130 field1,
131 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
132
133template double errorFrobRel(GeometricField<scalar, fvsPatchField, surfaceMesh>
134 &
135 field1,
136 GeometricField<scalar, fvsPatchField, surfaceMesh>& field2,
137 List<label>* labels);
138
139
140template<typename T>
141double errorLinfRel(GeometricField<T, fvPatchField, volMesh>& field1,
142 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels)
143{
144 double err;
145 autoPtr<GeometricField<T, fvPatchField, volMesh >> errField;
146 autoPtr<GeometricField<T, fvPatchField, volMesh >> field1_S;
147 autoPtr<GeometricField<T, fvPatchField, volMesh >> field2_S;
148 autoPtr<fvMeshSubset> submesh;
149
150 if (labels != NULL)
151 {
152 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(field1.mesh()));
153#if OPENFOAM >= 1812
154 submesh->setCellSubset(* labels);
155#else
156 submesh->setLargeCellSubset(* labels);
157#endif
158 GeometricField<T, fvPatchField, volMesh> field1tmp(submesh->interpolate(
159 field1));
160 GeometricField<T, fvPatchField, volMesh> field2tmp(submesh->interpolate(
161 field2));
162 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
163 (new GeometricField<T, fvPatchField, volMesh>(field1tmp.clone()));
164 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
165 (new GeometricField<T, fvPatchField, volMesh>(field2tmp.clone()));
166 }
167 else
168 {
169 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
170 (new GeometricField<T, fvPatchField, volMesh>(field1));
171 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
172 (new GeometricField<T, fvPatchField, volMesh>(field2));
173 }
174 errField = autoPtr<GeometricField<T, fvPatchField, volMesh >>
175 (new GeometricField<T, fvPatchField, volMesh>(field1_S() - field2_S()));
176
177 if (LinfNorm(field1_S()) <= 1e-6)
178 {
179 err = 0;
180 }
181 else
182 {
183 err = LinfNorm(errField()) / LinfNorm(field1_S());
184 }
185 return err;
186}
187
188
189template double errorLinfRel(GeometricField<scalar, fvPatchField, volMesh>&
190 field1,
191 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
192template double errorLinfRel(GeometricField<vector, fvPatchField, volMesh>&
193 field1,
194 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
195
196template<>
197double errorL2Abs(GeometricField<vector, fvPatchField, volMesh>& field1,
198 GeometricField<vector, fvPatchField, volMesh>& field2, volScalarField& Volumes)
199
200{
201 volScalarField diffFields2 = (((field1 - field2) & (field1 - field2)) *
202 Volumes).ref();
203 double err = Foam::sqrt(gSum(diffFields2));
204 return err;
205}
206
207template<>
208double errorL2Abs(GeometricField<scalar, fvPatchField, volMesh>& field1,
209 GeometricField<scalar, fvPatchField, volMesh>& field2, volScalarField& Volumes)
210
211{
212 volScalarField diffFields2 = (((field1 - field2) * (field1 - field2)) *
213 Volumes).ref();
214 double err = Foam::sqrt(gSum(diffFields2));
215 return err;
216}
217
218template<typename T>
219double errorL2Abs(GeometricField<T, fvPatchField, volMesh>& field1,
220 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels)
221{
222 autoPtr<GeometricField<T, fvPatchField, volMesh >> errField;
223 autoPtr<GeometricField<T, fvPatchField, volMesh >> field1_S;
224 autoPtr<GeometricField<T, fvPatchField, volMesh >> field2_S;
225 autoPtr<fvMeshSubset> submesh;
226
227 if (labels != NULL)
228 {
229 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(field1.mesh()));
230#if OPENFOAM >= 1812
231 submesh->setCellSubset(* labels);
232#else
233 submesh->setLargeCellSubset(* labels);
234#endif
235 GeometricField<T, fvPatchField, volMesh> field1tmp(submesh->interpolate(
236 field1));
237 GeometricField<T, fvPatchField, volMesh> field2tmp(submesh->interpolate(
238 field2));
239 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
240 (new GeometricField<T, fvPatchField, volMesh>(field1tmp.clone()));
241 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
242 (new GeometricField<T, fvPatchField, volMesh>(field2tmp.clone()));
243 }
244 else
245 {
246 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
247 (new GeometricField<T, fvPatchField, volMesh>(field1));
248 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
249 (new GeometricField<T, fvPatchField, volMesh>(field2));
250 }
251 errField = autoPtr<GeometricField<T, fvPatchField, volMesh >>
252 (new GeometricField<T, fvPatchField, volMesh>(field1_S() - field2_S()));
253 double err = L2Norm(errField());
254 return err;
255}
256
257template double errorL2Abs(GeometricField<scalar, fvPatchField, volMesh>&
258 field1,
259 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
260template double errorL2Abs(GeometricField<vector, fvPatchField, volMesh>&
261 field1,
262 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
263
264
265template<class T, template<class> class PatchField, class GeoMesh>
266Eigen::MatrixXd errorFrobRel(PtrList<GeometricField<T, PatchField, GeoMesh >> &
267 fields1,
268 PtrList<GeometricField<T, PatchField, GeoMesh >>& fields2, List<label>* labels)
269{
270 Eigen::VectorXd err;
271
272 if (fields1.size() != fields2.size())
273 {
274 Info << "The two fields do not have the same size, code will abort" << endl;
275 exit(0);
276 }
277
278 err.resize(fields1.size(), 1);
279
280 for (label k = 0; k < fields1.size(); k++)
281 {
282 err(k, 0) = errorFrobRel(fields1[k], fields2[k], labels);
283 Info << " Error is " << err[k] << endl;
284 }
285
286 return err;
287}
288
289template Eigen::MatrixXd errorFrobRel(
290 PtrList<GeometricField<scalar, fvPatchField, volMesh >> & fields1,
291 PtrList<GeometricField<scalar, fvPatchField, volMesh >>& fields2,
292 List<label>* labels);
293template Eigen::MatrixXd errorFrobRel(
294 PtrList<GeometricField<vector, fvPatchField, volMesh >> & fields1,
295 PtrList<GeometricField<vector, fvPatchField, volMesh >>& fields2,
296 List<label>* labels);
297template Eigen::MatrixXd errorFrobRel(
298 PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh >> & fields1,
299 PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh >>& fields2,
300 List<label>* labels);
301
302
303template<typename T>
304Eigen::MatrixXd errorL2Abs(
305 PtrList<GeometricField<T, fvPatchField, volMesh >> & fields1,
306 PtrList<GeometricField<T, fvPatchField, volMesh >>& fields2,
307 PtrList<volScalarField>& Volumes)
308{
309 M_Assert(fields1.size() == fields2.size(),
310 "The two fields do not have the same size, code will abort");
311 M_Assert(fields1.size() == Volumes.size(),
312 "The volumes field and the two solution fields do not have the same size, code will abort");
313 Eigen::VectorXd err;
314 err.resize(fields1.size(), 1);
315
316 for (label k = 0; k < fields1.size(); k++)
317 {
318 err(k, 0) = errorL2Abs(fields1[k], fields2[k], Volumes[k]);
319 Info << " Error is " << err[k] << endl;
320 }
321
322 return err;
323}
324
325template Eigen::MatrixXd errorL2Abs(
326 PtrList<GeometricField<scalar, fvPatchField, volMesh >> & fields1,
327 PtrList<GeometricField<scalar, fvPatchField, volMesh >>& fields2,
328 PtrList<volScalarField>& Volumes);
329template Eigen::MatrixXd errorL2Abs(
330 PtrList<GeometricField<vector, fvPatchField, volMesh >> & fields1,
331 PtrList<GeometricField<vector, fvPatchField, volMesh >>& fields2,
332 PtrList<volScalarField>& Volumes);
333
334template<typename T>
335Eigen::MatrixXd errorL2Abs(PtrList<GeometricField<T, fvPatchField, volMesh >> &
336 fields1,
337 PtrList<GeometricField<T, fvPatchField, volMesh >>& fields2,
338 List<label>* labels)
339{
340 Eigen::VectorXd err;
341
342 if (fields1.size() != fields2.size())
343 {
344 Info << "The two fields do not have the same size, code will abort" << endl;
345 exit(0);
346 }
347
348 err.resize(fields1.size(), 1);
349
350 for (label k = 0; k < fields1.size(); k++)
351 {
352 err(k, 0) = errorL2Abs(fields1[k], fields2[k], labels);
353 Info << " Error is " << err[k] << endl;
354 }
355
356 return err;
357}
358
359template Eigen::MatrixXd errorL2Abs(
360 PtrList<GeometricField<scalar, fvPatchField, volMesh >> & fields1,
361 PtrList<GeometricField<scalar, fvPatchField, volMesh >>& fields2,
362 List<label>* labels);
363template Eigen::MatrixXd errorL2Abs(
364 PtrList<GeometricField<vector, fvPatchField, volMesh >> & fields1,
365 PtrList<GeometricField<vector, fvPatchField, volMesh >>& fields2,
366 List<label>* labels);
367
368template<typename T>
369double errorL2Rel(GeometricField<T, fvPatchField, volMesh>& field1,
370 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels)
371{
372 double err;
373 autoPtr<GeometricField<T, fvPatchField, volMesh >> errField;
374 autoPtr<GeometricField<T, fvPatchField, volMesh >> field1_S;
375 autoPtr<GeometricField<T, fvPatchField, volMesh >> field2_S;
376 autoPtr<fvMeshSubset> submesh;
377
378 if (labels != NULL)
379 {
380 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(field1.mesh()));
381#if OPENFOAM >= 1812
382 submesh->setCellSubset(* labels);
383#else
384 submesh->setLargeCellSubset(* labels);
385#endif
386 GeometricField<T, fvPatchField, volMesh> field1tmp(submesh->interpolate(
387 field1));
388 GeometricField<T, fvPatchField, volMesh> field2tmp(submesh->interpolate(
389 field2));
390 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
391 (new GeometricField<T, fvPatchField, volMesh>(field1tmp.clone()));
392 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
393 (new GeometricField<T, fvPatchField, volMesh>(field2tmp.clone()));
394 }
395 else
396 {
397 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
398 (new GeometricField<T, fvPatchField, volMesh>(field1));
399 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh >>
400 (new GeometricField<T, fvPatchField, volMesh>(field2));
401 }
402 errField = autoPtr<GeometricField<T, fvPatchField, volMesh >>
403 (new GeometricField<T, fvPatchField, volMesh>(field1_S() - field2_S()));
404
405 if (L2Norm(field1) <= 1e-6)
406 {
407 err = 0;
408 }
409 else
410 {
411 err = L2Norm(errField()) / L2Norm(
412 field1_S());
413 }
414 return err;
415}
416
417template double errorL2Rel(GeometricField<scalar, fvPatchField, volMesh>&
418 field1,
419 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
420template double errorL2Rel(GeometricField<vector, fvPatchField, volMesh>&
421 field1,
422 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
423
424template<typename T>
425Eigen::MatrixXd errorL2Rel(PtrList<GeometricField<T, fvPatchField, volMesh >> &
426 fields1, PtrList<GeometricField<T, fvPatchField, volMesh >>& fields2,
427 List<label>* labels)
428{
429 Eigen::VectorXd err;
430
431 if (fields1.size() != fields2.size())
432 {
433 Info << "The two fields do not have the same size, code will abort" << endl;
434 exit(0);
435 }
436
437 err.resize(fields1.size(), 1);
438
439 for (label k = 0; k < fields1.size(); k++)
440 {
441 err(k, 0) = errorL2Rel(fields1[k], fields2[k], labels);
442 Info << " Error is " << err[k] << endl;
443 }
444
445 return err;
446}
447
448template Eigen::MatrixXd errorL2Rel(
449 PtrList<GeometricField<scalar, fvPatchField, volMesh >> & fields1,
450 PtrList<GeometricField<scalar, fvPatchField, volMesh >>& fields2,
451 List<label>* labels);
452template Eigen::MatrixXd errorL2Rel(
453 PtrList<GeometricField<vector, fvPatchField, volMesh >> & fields1,
454 PtrList<GeometricField<vector, fvPatchField, volMesh >>& fields2,
455 List<label>* labels);
456
457template<>
458double H1Seminorm(GeometricField<scalar, fvPatchField, volMesh>& field)
459{
460 double a;
461 a = Foam::sqrt(fvc::domainIntegrate(fvc::grad(field) & fvc::grad(
462 field)).value());
463 return a;
464}
465
466template<>
467double H1Seminorm(GeometricField<vector, fvPatchField, volMesh>& field)
468{
469 double a;
470 a = Foam::sqrt(fvc::domainIntegrate(fvc::grad(field)
471 && fvc::grad(field)).value());
472 return a;
473}
474
475
476template<class Type, template<class> class PatchField, class GeoMesh>
477double frobNorm(GeometricField<Type, PatchField, GeoMesh>& field)
478{
479 double norm(0);
480 Eigen::VectorXd vF = Foam2Eigen::field2Eigen(field);
481 norm = vF.norm();
482 return norm;
483}
484
485template double frobNorm(GeometricField<scalar, fvPatchField, volMesh>& field);
486template double frobNorm(GeometricField<vector, fvPatchField, volMesh>& field);
487
488double L2normOnPatch(fvMesh& mesh, volScalarField& field,
489 word patch)
490{
491 double L2 = 0;
492 //Access the mesh information for the boundary
493 label patchID = mesh.boundaryMesh().findPatchID(patch);
494 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
495 //List of cells close to a boundary
496 const labelUList& faceCells = cPatch.faceCells();
497 forAll(cPatch, faceI)
498 {
499 //id of the owner cell having the face
500 label faceOwner = faceCells[faceI] ;
501 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
502 L2 += faceArea * field[faceOwner] * field[faceOwner];
503 }
504
505 return Foam::sqrt(L2);
506}
507
508double L2productOnPatch(fvMesh& mesh, List<scalar>& field1,
509 List<scalar>& field2, word patch)
510{
511 M_Assert(field1.size() == field2.size(),
512 "The two fields do not have the same size, code will abort");
513 double L2 = 0;
514 //Access the mesh information for the boundary
515 label patchID = mesh.boundaryMesh().findPatchID(patch);
516 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
517 M_Assert(field1.size() == cPatch.size(),
518 "The filds must have the same size of the patch, code will abort");
519 forAll(cPatch, faceI)
520 {
521 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
522 L2 += faceArea * field1[faceI] * field2[faceI];
523 }
524
525 return L2;
526}
527
528double LinfNormOnPatch(fvMesh& mesh, volScalarField& field,
529 word patch)
530{
531 double Linf = 0;
532 //Access the mesh information for the boundary
533 label patchID = mesh.boundaryMesh().findPatchID(patch);
534 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
535 //List of cells close to a boundary
536 const labelUList& faceCells = cPatch.faceCells();
537 forAll(cPatch, faceI)
538 {
539 label faceOwner = faceCells[faceI] ;
540
541 if (faceI == 0)
542 {
543 Linf = std::abs(field[faceOwner]);
544 }
545 else if (std::abs(field[faceOwner]) > Linf)
546 {
547 Linf = std::abs(field[faceOwner]);
548 }
549 }
550
551 return Linf;
552}
553
554double integralOnPatch(fvMesh& mesh, volScalarField& field,
555 word patch)
556{
557 double integral = 0;
558 //Access the mesh information for the boundary
559 label patchID = mesh.boundaryMesh().findPatchID(patch);
560 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
561 //List of cells close to a boundary
562 const labelUList& faceCells = cPatch.faceCells();
563 forAll(cPatch, faceI)
564 {
565 //id of the owner cell having the face
566 label faceOwner = faceCells[faceI] ;
567 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
568 integral += faceArea * field[faceOwner];
569 }
570
571 return integral;
572}
573
574double integralOnPatch(fvMesh& mesh, List<scalar> field,
575 word patch)
576{
577 double integral = 0;
578 //Access the mesh information for the boundary
579 label patchID = mesh.boundaryMesh().findPatchID(patch);
580 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
581 //List of cells close to a boundary
582 const labelUList& faceCells = cPatch.faceCells();
583 int patchSize = mesh.magSf().boundaryField()[patchID].size();
584 std::string message = "The input list (size = " + std::to_string(field.size())
585 + ") must have the same size of the patch mesh (size = "
586 + std::to_string(patchSize) + ")";
587 M_Assert( patchSize == field.size(), message.c_str());
588 forAll(cPatch, faceI)
589 {
590 //id of the owner cell having the face
591 label faceOwner = faceCells[faceI] ;
592 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
593 integral += faceArea * field[faceI];
594 }
595
596 return integral;
597}
598
599}
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
#define M_Assert(Expr, Msg)
Header file of the ITHACAerror file.
static Eigen::VectorXd field2Eigen(GeometricField< tensor, PatchField, GeoMesh > &field)
Convert a vector OpenFOAM field into an Eigen Vector.
Definition Foam2Eigen.C:40
Namespace to implement some useful assign operation of OF fields.
double LinfNorm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:57
double errorLinfRel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in Linf norm.
double LinfNormOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the Linf norm of a field on a boundary patch.
double H1Seminorm(GeometricField< scalar, fvPatchField, volMesh > &field)
double L2Norm(GeometricField< scalar, fvPatchField, volMesh > &field)
Definition ITHACAerror.C:41
double integralOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the integral on a patch.
double L2normOnPatch(fvMesh &mesh, volScalarField &field, word patch)
Evaluate the L2 norm of a field on a boundary patch.
double errorL2Rel(GeometricField< T, fvPatchField, volMesh > &field1, GeometricField< T, fvPatchField, volMesh > &field2, List< label > *labels)
Computes the relative error between two geometric Fields in L2 norm.
double frobNorm(GeometricField< Type, PatchField, GeoMesh > &field)
Evaluate the Frobenius norm of a field.
double errorFrobRel(GeometricField< Type, PatchField, GeoMesh > &field1, GeometricField< Type, PatchField, GeoMesh > &field2, List< label > *labels)
Computes the relative error between two Fields in the Frobenius norm.
Definition ITHACAerror.C:78
double L2productOnPatch(fvMesh &mesh, List< scalar > &field1, List< scalar > &field2, word patch)
Evaluate the L2 inner product between two scalarLists.
double errorL2Abs(GeometricField< vector, fvPatchField, volMesh > &field1, GeometricField< vector, fvPatchField, volMesh > &field2, volScalarField &Volumes)