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
112 errField = autoPtr<GeometricField<Type, PatchField, GeoMesh>>
113 (new GeometricField<Type, PatchField, GeoMesh>(field1_S() - field2_S()));
114
115 if (frobNorm(field1) <= 1e-6)
116 {
117 err = 0;
118 }
119 else
120 {
121 err = frobNorm(errField()) / frobNorm(field1_S());
122 }
123
124 return err;
125}
126
127
128template double errorFrobRel(GeometricField<scalar, fvPatchField, volMesh>&
129 field1,
130 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
131template double errorFrobRel(GeometricField<vector, fvPatchField, volMesh>&
132 field1,
133 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
134
135template double errorFrobRel(GeometricField<scalar, fvsPatchField, surfaceMesh>&
136 field1,
137 GeometricField<scalar, fvsPatchField, surfaceMesh>& field2,
138 List<label>* labels);
139
140
141template<typename T>
142double errorLinfRel(GeometricField<T, fvPatchField, volMesh>& field1,
143 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels)
144{
145 double err;
146 autoPtr<GeometricField<T, fvPatchField, volMesh>> errField;
147 autoPtr<GeometricField<T, fvPatchField, volMesh>> field1_S;
148 autoPtr<GeometricField<T, fvPatchField, volMesh>> field2_S;
149 autoPtr<fvMeshSubset> submesh;
150
151 if (labels != NULL)
152 {
153 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(field1.mesh()));
154#if OPENFOAM >= 1812
155 submesh->setCellSubset(*labels);
156#else
157 submesh->setLargeCellSubset(*labels);
158#endif
159 GeometricField<T, fvPatchField, volMesh> field1tmp(submesh->interpolate(
160 field1));
161 GeometricField<T, fvPatchField, volMesh> field2tmp(submesh->interpolate(
162 field2));
163 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
164 (new GeometricField<T, fvPatchField, volMesh>(field1tmp.clone()));
165 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
166 (new GeometricField<T, fvPatchField, volMesh>(field2tmp.clone()));
167 }
168 else
169 {
170 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
171 (new GeometricField<T, fvPatchField, volMesh>(field1));
172 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
173 (new GeometricField<T, fvPatchField, volMesh>(field2));
174 }
175
176 errField = autoPtr<GeometricField<T, fvPatchField, volMesh>>
177 (new GeometricField<T, fvPatchField, volMesh>(field1_S() - field2_S()));
178
179 if (LinfNorm(field1_S()) <= 1e-6)
180 {
181 err = 0;
182 }
183 else
184 {
185 err = LinfNorm(errField()) / LinfNorm(field1_S());
186 }
187
188 return err;
189}
190
191
192template double errorLinfRel(GeometricField<scalar, fvPatchField, volMesh>&
193 field1,
194 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
195template double errorLinfRel(GeometricField<vector, fvPatchField, volMesh>&
196 field1,
197 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
198
199template<>
200double errorL2Abs(GeometricField<vector, fvPatchField, volMesh>& field1,
201 GeometricField<vector, fvPatchField, volMesh>& field2, volScalarField& Volumes)
202
203{
204 volScalarField diffFields2 = (((field1 - field2) & (field1 - field2)) *
205 Volumes).ref();
206 double err = Foam::sqrt(gSum(diffFields2));
207 return err;
208}
209
210template<>
211double errorL2Abs(GeometricField<scalar, fvPatchField, volMesh>& field1,
212 GeometricField<scalar, fvPatchField, volMesh>& field2, volScalarField& Volumes)
213
214{
215 volScalarField diffFields2 = (((field1 - field2) * (field1 - field2)) *
216 Volumes).ref();
217 double err = Foam::sqrt(gSum(diffFields2));
218 return err;
219}
220
221template<typename T>
222double errorL2Abs(GeometricField<T, fvPatchField, volMesh>& field1,
223 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels)
224{
225 autoPtr<GeometricField<T, fvPatchField, volMesh>> errField;
226 autoPtr<GeometricField<T, fvPatchField, volMesh>> field1_S;
227 autoPtr<GeometricField<T, fvPatchField, volMesh>> field2_S;
228 autoPtr<fvMeshSubset> submesh;
229
230 if (labels != NULL)
231 {
232 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(field1.mesh()));
233#if OPENFOAM >= 1812
234 submesh->setCellSubset(*labels);
235#else
236 submesh->setLargeCellSubset(*labels);
237#endif
238 GeometricField<T, fvPatchField, volMesh> field1tmp(submesh->interpolate(
239 field1));
240 GeometricField<T, fvPatchField, volMesh> field2tmp(submesh->interpolate(
241 field2));
242 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
243 (new GeometricField<T, fvPatchField, volMesh>(field1tmp.clone()));
244 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
245 (new GeometricField<T, fvPatchField, volMesh>(field2tmp.clone()));
246 }
247 else
248 {
249 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
250 (new GeometricField<T, fvPatchField, volMesh>(field1));
251 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
252 (new GeometricField<T, fvPatchField, volMesh>(field2));
253 }
254
255 errField = autoPtr<GeometricField<T, fvPatchField, volMesh>>
256 (new GeometricField<T, fvPatchField, volMesh>(field1_S() - field2_S()));
257 double err = L2Norm(errField());
258 return err;
259}
260
261template double errorL2Abs(GeometricField<scalar, fvPatchField, volMesh>&
262 field1,
263 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
264template double errorL2Abs(GeometricField<vector, fvPatchField, volMesh>&
265 field1,
266 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
267
268
269template<class T, template<class> class PatchField, class GeoMesh>
270Eigen::MatrixXd errorFrobRel(PtrList<GeometricField<T, PatchField, GeoMesh>>&
271 fields1,
272 PtrList<GeometricField<T, PatchField, GeoMesh>>& fields2, List<label>* labels)
273{
274 Eigen::VectorXd err;
275
276 if (fields1.size() != fields2.size())
277 {
278 Info << "The two fields do not have the same size, code will abort" << endl;
279 exit(0);
280 }
281
282 err.resize(fields1.size(), 1);
283
284 for (label k = 0; k < fields1.size(); k++)
285 {
286 err(k, 0) = errorFrobRel(fields1[k], fields2[k], labels);
287 Info << " Error is " << err[k] << endl;
288 }
289
290 return err;
291}
292
293template Eigen::MatrixXd errorFrobRel(
294 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields1,
295 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields2,
296 List<label>* labels);
297template Eigen::MatrixXd errorFrobRel(
298 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields1,
299 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields2,
300 List<label>* labels);
301template Eigen::MatrixXd errorFrobRel(
302 PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh>>& fields1,
303 PtrList<GeometricField<scalar, fvsPatchField, surfaceMesh>>& fields2,
304 List<label>* labels);
305
306
307template<typename T>
308Eigen::MatrixXd errorL2Abs(
309 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields1,
310 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
311 PtrList<volScalarField>& Volumes)
312{
313 M_Assert(fields1.size() == fields2.size(),
314 "The two fields do not have the same size, code will abort");
315 M_Assert(fields1.size() == Volumes.size(),
316 "The volumes field and the two solution fields do not have the same size, code will abort");
317 Eigen::VectorXd err;
318 err.resize(fields1.size(), 1);
319
320 for (label k = 0; k < fields1.size(); k++)
321 {
322 err(k, 0) = errorL2Abs(fields1[k], fields2[k], Volumes[k]);
323 Info << " Error is " << err[k] << endl;
324 }
325
326 return err;
327}
328
329template Eigen::MatrixXd errorL2Abs(
330 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields1,
331 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields2,
332 PtrList<volScalarField>& Volumes);
333template Eigen::MatrixXd errorL2Abs(
334 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields1,
335 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields2,
336 PtrList<volScalarField>& Volumes);
337
338template<typename T>
339Eigen::MatrixXd errorL2Abs(PtrList<GeometricField<T, fvPatchField, volMesh>>&
340 fields1,
341 PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2, List<label>* labels)
342{
343 Eigen::VectorXd err;
344
345 if (fields1.size() != fields2.size())
346 {
347 Info << "The two fields do not have the same size, code will abort" << endl;
348 exit(0);
349 }
350
351 err.resize(fields1.size(), 1);
352
353 for (label k = 0; k < fields1.size(); k++)
354 {
355 err(k, 0) = errorL2Abs(fields1[k], fields2[k], labels);
356 Info << " Error is " << err[k] << endl;
357 }
358
359 return err;
360}
361
362template Eigen::MatrixXd errorL2Abs(
363 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields1,
364 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields2,
365 List<label>* labels);
366template Eigen::MatrixXd errorL2Abs(
367 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields1,
368 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields2,
369 List<label>* labels);
370
371template<typename T>
372double errorL2Rel(GeometricField<T, fvPatchField, volMesh>& field1,
373 GeometricField<T, fvPatchField, volMesh>& field2, List<label>* labels)
374{
375 double err;
376 autoPtr<GeometricField<T, fvPatchField, volMesh>> errField;
377 autoPtr<GeometricField<T, fvPatchField, volMesh>> field1_S;
378 autoPtr<GeometricField<T, fvPatchField, volMesh>> field2_S;
379 autoPtr<fvMeshSubset> submesh;
380
381 if (labels != NULL)
382 {
383 submesh = autoPtr<fvMeshSubset>(new fvMeshSubset(field1.mesh()));
384#if OPENFOAM >= 1812
385 submesh->setCellSubset(*labels);
386#else
387 submesh->setLargeCellSubset(*labels);
388#endif
389 GeometricField<T, fvPatchField, volMesh> field1tmp(submesh->interpolate(
390 field1));
391 GeometricField<T, fvPatchField, volMesh> field2tmp(submesh->interpolate(
392 field2));
393 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
394 (new GeometricField<T, fvPatchField, volMesh>(field1tmp.clone()));
395 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
396 (new GeometricField<T, fvPatchField, volMesh>(field2tmp.clone()));
397 }
398 else
399 {
400 field1_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
401 (new GeometricField<T, fvPatchField, volMesh>(field1));
402 field2_S = autoPtr<GeometricField<T, fvPatchField, volMesh>>
403 (new GeometricField<T, fvPatchField, volMesh>(field2));
404 }
405
406 errField = autoPtr<GeometricField<T, fvPatchField, volMesh>>
407 (new GeometricField<T, fvPatchField, volMesh>(field1_S() - field2_S()));
408
409 if (L2Norm(field1) <= 1e-6)
410 {
411 err = 0;
412 }
413 else
414 {
415 err = L2Norm(errField()) / L2Norm(
416 field1_S());
417 }
418
419 return err;
420}
421
422template double errorL2Rel(GeometricField<scalar, fvPatchField, volMesh>&
423 field1,
424 GeometricField<scalar, fvPatchField, volMesh>& field2, List<label>* labels);
425template double errorL2Rel(GeometricField<vector, fvPatchField, volMesh>&
426 field1,
427 GeometricField<vector, fvPatchField, volMesh>& field2, List<label>* labels);
428
429template<typename T>
430Eigen::MatrixXd errorL2Rel(PtrList<GeometricField<T, fvPatchField, volMesh>>&
431 fields1, PtrList<GeometricField<T, fvPatchField, volMesh>>& fields2,
432 List<label>* labels)
433{
434 Eigen::VectorXd err;
435
436 if (fields1.size() != fields2.size())
437 {
438 Info << "The two fields do not have the same size, code will abort" << endl;
439 exit(0);
440 }
441
442 err.resize(fields1.size(), 1);
443
444 for (label k = 0; k < fields1.size(); k++)
445 {
446 err(k, 0) = errorL2Rel(fields1[k], fields2[k], labels);
447 Info << " Error is " << err[k] << endl;
448 }
449
450 return err;
451}
452
453template Eigen::MatrixXd errorL2Rel(
454 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields1,
455 PtrList<GeometricField<scalar, fvPatchField, volMesh>>& fields2,
456 List<label>* labels);
457template Eigen::MatrixXd errorL2Rel(
458 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields1,
459 PtrList<GeometricField<vector, fvPatchField, volMesh>>& fields2,
460 List<label>* labels);
461
462template<>
463double H1Seminorm(GeometricField<scalar, fvPatchField, volMesh>& field)
464{
465 double a;
466 a = Foam::sqrt(fvc::domainIntegrate(fvc::grad(field) & fvc::grad(
467 field)).value());
468 return a;
469}
470
471template<>
472double H1Seminorm(GeometricField<vector, fvPatchField, volMesh>& field)
473{
474 double a;
475 a = Foam::sqrt(fvc::domainIntegrate(fvc::grad(field)
476 && fvc::grad(field)).value());
477 return a;
478}
479
480
481template<class Type, template<class> class PatchField, class GeoMesh>
482double frobNorm(GeometricField<Type, PatchField, GeoMesh>& field)
483{
484 double norm(0);
485 Eigen::VectorXd vF = Foam2Eigen::field2Eigen(field);
486 norm = vF.norm();
487 return norm;
488}
489
490template double frobNorm(GeometricField<scalar, fvPatchField, volMesh>& field);
491template double frobNorm(GeometricField<vector, fvPatchField, volMesh>& field);
492
493double L2normOnPatch(fvMesh& mesh, volScalarField& field,
494 word patch)
495{
496 double L2 = 0;
497 //Access the mesh information for the boundary
498 label patchID = mesh.boundaryMesh().findPatchID(patch);
499 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
500 //List of cells close to a boundary
501 const labelUList& faceCells = cPatch.faceCells();
502 forAll(cPatch, faceI)
503 {
504 //id of the owner cell having the face
505 label faceOwner = faceCells[faceI] ;
506 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
507 L2 += faceArea * field[faceOwner] * field[faceOwner];
508 }
509 return Foam::sqrt(L2);
510}
511
512double L2productOnPatch(fvMesh& mesh, List<scalar>& field1,
513 List<scalar>& field2, word patch)
514{
515 M_Assert(field1.size() == field2.size(),
516 "The two fields do not have the same size, code will abort");
517 double L2 = 0;
518 //Access the mesh information for the boundary
519 label patchID = mesh.boundaryMesh().findPatchID(patch);
520 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
521 M_Assert(field1.size() == cPatch.size(),
522 "The filds must have the same size of the patch, code will abort");
523 forAll(cPatch, faceI)
524 {
525 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
526 L2 += faceArea * field1[faceI] * field2[faceI];
527 }
528 return L2;
529}
530
531double LinfNormOnPatch(fvMesh& mesh, volScalarField& field,
532 word patch)
533{
534 double Linf = 0;
535 //Access the mesh information for the boundary
536 label patchID = mesh.boundaryMesh().findPatchID(patch);
537 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
538 //List of cells close to a boundary
539 const labelUList& faceCells = cPatch.faceCells();
540 forAll(cPatch, faceI)
541 {
542 label faceOwner = faceCells[faceI] ;
543
544 if (faceI == 0)
545 {
546 Linf = std::abs(field[faceOwner]);
547 }
548 else if (std::abs(field[faceOwner]) > Linf)
549 {
550 Linf = std::abs(field[faceOwner]);
551 }
552 }
553 return Linf;
554}
555
556double integralOnPatch(fvMesh& mesh, volScalarField& field,
557 word patch)
558{
559 double integral = 0;
560 //Access the mesh information for the boundary
561 label patchID = mesh.boundaryMesh().findPatchID(patch);
562 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
563 //List of cells close to a boundary
564 const labelUList& faceCells = cPatch.faceCells();
565 forAll(cPatch, faceI)
566 {
567 //id of the owner cell having the face
568 label faceOwner = faceCells[faceI] ;
569 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
570 integral += faceArea * field[faceOwner];
571 }
572 return integral;
573}
574
575double integralOnPatch(fvMesh& mesh, List<scalar> field,
576 word patch)
577{
578 double integral = 0;
579 //Access the mesh information for the boundary
580 label patchID = mesh.boundaryMesh().findPatchID(patch);
581 const polyPatch& cPatch = mesh.boundaryMesh()[patchID];
582 //List of cells close to a boundary
583 const labelUList& faceCells = cPatch.faceCells();
584 int patchSize = mesh.magSf().boundaryField()[patchID].size();
585 std::string message = "The input list (size = " + std::to_string(field.size())
586 + ") must have the same size of the patch mesh (size = "
587 + std::to_string(patchSize) + ")";
588 M_Assert( patchSize == field.size(), message.c_str());
589 forAll(cPatch, faceI)
590 {
591 //id of the owner cell having the face
592 label faceOwner = faceCells[faceI] ;
593 scalar faceArea = mesh.magSf().boundaryField()[patchID][faceI];
594 integral += faceArea * field[faceI];
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)