Loading...
Searching...
No Matches
inverseLaplacianProblem_CG.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
13 License
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
33
34
36
37// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
38
39// Constructors
41
43 :
45{
46 Foam::Time& runTime = _runTime();
47 Foam::fvMesh& mesh = _mesh();
48#include "createFields.H"
49}
50
51// * * * * * * * * * * * * * * Full Order Methods * * * * * * * * * * * * * * //
52
54{
55 fvMesh& mesh = _mesh();
56 valueFraction.resize(mesh.boundaryMesh()["coldSide"].size());
57 homogeneousBCcoldSide.resize(mesh.boundaryMesh()["coldSide"].size());
58 valueFractionAdj.resize(mesh.boundaryMesh()["coldSide"].size());
59 Eigen::VectorXd faceCellDist =
61 forAll (valueFraction, faceI)
62 {
63 scalar faceDist = faceCellDist(faceI);
64 valueFraction[faceI] = 1.0 / (1.0 + (k / H / faceDist));
65 valueFractionAdj[faceI] = 1 / (1 + (1 / k / H / faceDist));
66 homogeneousBCcoldSide[faceI] = 0;
67 }
69}
70
72{
73 fvMesh& mesh = _mesh();
74 volScalarField& lambda = _lambda();
76 forAll(mesh.boundaryMesh(), patchI)
77 {
78 if (patchI == mesh.boundaryMesh().findPatchID("coldSide"))
79 {
82 }
83 else if (patchI == mesh.boundaryMesh().findPatchID("hotSide"))
84 {
86 }
87 else
88 {
90 }
91 }
93}
94
96{
97 volScalarField& lambda = _lambda();
98 fvMesh& mesh = _mesh();
100 dimensionedScalar sourceDim("sourceDim", dimensionSet(1, -1, -3, -1, 0, 0, 0),
101 1);
102 autoPtr<volScalarField> f_
103 (
104 new volScalarField("f", lambda)
105 );
106 volScalarField& f = f_();
107
108 if (interpolation)
109 {
111 {
112 f.ref()[interpolationPlane.cellID [cellI]] =
113 interpolationPlane.Tdiff[cellI] * k;
114 }
115 }
116 else
117 {
118 for (int i = 0; i < thermocouplesCellID.size(); i++)
119 {
120 f.ref()[thermocouplesCellID [i]] = Tdiff(i) * k /
122 }
123 }
124
125 return (f * sourceDim).ref();
126}
127
129{
130 fvMesh& mesh = _mesh();
131 volScalarField& deltaT = _deltaT();
133 forAll(mesh.boundaryMesh(), patchI)
134 {
135 if (patchI == mesh.boundaryMesh().findPatchID("coldSide"))
136 {
139 }
140 else if (patchI == mesh.boundaryMesh().findPatchID("hotSide"))
141 {
142 ITHACAutilities::assignBC(deltaT, patchI, - P / k);
143 }
144 else
145 {
147 }
148 }
149}
150
152{
153 restart();
154 volScalarField& lambda = _lambda();
155 Foam::Time& runTime = _runTime();
156 volScalarField f = assignAdjointBCandSource();
157 simpleControl& simple = _simple();
158#if defined(OFVER) && (OFVER == 6)
159
160 while (simple.loop(runTime))
161#else
162 while (simple.loop())
163#endif
164 {
165 //Info << "Time = " << runTime.timeName() << nl << endl;
166 while (simple.correctNonOrthogonal())
167 {
168 fvScalarMatrix TEqn
169 (
170 fvm::laplacian(DT, lambda) == - f
171 );
172 TEqn.solve();
173 }
174 }
175}
176
177
179{
180 restart();
182 solve("sensitivity");
183}
184
185void inverseLaplacianProblem_CG::solve(const char* problemID)
186{
187 volScalarField& T = _T();
188 volScalarField& deltaT = _deltaT();
189 simpleControl& simple = _simple();
190 Foam::Time& runTime = _runTime();
191
192 if (strcmp( problemID, "direct") == 0)
193 {
195 }
196 else if (strcmp( problemID, "sensitivity") == 0)
197 {
199 }
200 else
201 {
202 Info << "Problem name should be direct or sensitivity" << endl;
203 exit(10);
204 }
205
206#if defined(OFVER) && (OFVER == 6)
207
208 while (simple.loop(runTime))
209#else
210 while (simple.loop())
211#endif
212 {
213 while (simple.correctNonOrthogonal())
214 {
215 if (strcmp( problemID, "direct") == 0)
216 {
217 fvScalarMatrix TEqn
218 (
219 fvm::laplacian(DT, T)
220 );
221 TEqn.solve();
222 }
223 else if (strcmp( problemID, "sensitivity") == 0)
224 {
225 fvScalarMatrix TEqn
226 (
227 fvm::laplacian(DT, deltaT)
228 );
229 TEqn.solve();
230 }
231 }
232 }
233}
234
236{
237 Info << "Defining the plane for measurements interpolation" << endl;
238 fvMesh& mesh = _mesh();
239 //Define thermocouples plane
240 bool firstCell = 1;
242 {
243 if (thermocouplesCellID[cellI] != -1)
244 {
245 thermocouplesCellProc[cellI] = Pstream::myProcNo();
246
247 if (firstCell)
248 {
249 interpolationPlane.minX = mesh.C()[thermocouplesCellID[cellI]].component(0);
250 interpolationPlane.maxX = mesh.C()[thermocouplesCellID[cellI]].component(0);
251 interpolationPlane.Y = mesh.C()[thermocouplesCellID[cellI]].component(1);
252 interpolationPlane.minZ = mesh.C()[thermocouplesCellID[cellI]].component(2);
253 interpolationPlane.maxZ = mesh.C()[thermocouplesCellID[cellI]].component(2);
254 firstCell = 0;
255 }
256
257 if (mesh.C()[thermocouplesCellID[cellI]].component(0) < interpolationPlane.minX)
258 {
259 interpolationPlane.minX = mesh.C()[thermocouplesCellID[cellI]].component(0);
260 }
261
262 if (mesh.C()[thermocouplesCellID[cellI]].component(0) > interpolationPlane.maxX)
263 {
264 interpolationPlane.maxX = mesh.C()[thermocouplesCellID[cellI]].component(0) ;
265 }
266
267 if (mesh.C()[thermocouplesCellID[cellI]].component(2) < interpolationPlane.minZ)
268 {
269 interpolationPlane.minZ = mesh.C()[thermocouplesCellID[cellI]].component(2);
270 }
271
272 if (mesh.C()[thermocouplesCellID[cellI]].component(2) > interpolationPlane.maxZ)
273 {
274 interpolationPlane.maxZ = mesh.C()[thermocouplesCellID[cellI]].component(2) ;
275 }
276 }
277 else
278 {
279 Tmeas (cellI) = 0;
280 thermocouplesCellProc[cellI] = -1;
281 }
282 }
283}
284
286{
287 volScalarField& deltaT = _deltaT();
288
289 if (interpolation)
290 {
292 {
294 deltaT.internalField()[interpolationPlane.cellID [cellI]];
295 }
296 }
297 else
298 {
300 }
301}
302
304{
305 set_g();
307 cgIter = 0;
308 J = 0;
309 P = g;
310 gradJ = g; //Gradient of the cost function [W/m2]
311 gamma = 0.0;
312 gamma_den = 0.0;
313 label sampleI = 1;
314 gList.resize(0);
315 Tfield.resize(0);
316 lambdaField.resize(0);
317 deltaTfield.resize(0);
318
319 while (cgIter < cgIterMax)
320 {
321 Info << "Iteration " << cgIter + 1 << endl;
322 restart();
323 solveDirect();
324
325 if (saveSolInLists && cgIter == 0)
326 {
327 gList.append(g.clone());
328 }
329
330 volScalarField& T = _T();
331 //ITHACAstream::exportSolution(T, std::to_string(sampleI),
332 // "./ITHACAoutput/CGtest/", T.name());
334
336 {
337 Jlist.conservativeResize(cgIter + 1, 1);
338 Jlist(cgIter) = J;
339 ITHACAstream::exportMatrix(Jlist, "costFunctionFull", "eigen", "./");
340 return (1);
341 }
342
343 Jlist.conservativeResize(cgIter + 1, 1);
344 Jlist(cgIter) = J;
345 solveAdjoint();
346 volScalarField& lambda = _lambda();
347 //ITHACAstream::exportSolution(lambda, std::to_string(sampleI),
348 // "./ITHACAoutput/CGtest/", lambda.name());
349 computeGradJ();
352 volScalarField& deltaT = _deltaT();
353 //ITHACAstream::exportSolution(deltaT, std::to_string(sampleI),
354 // "./ITHACAoutput/CGtest/", deltaT.name());
358
359 if (saveSolInLists)
360 {
361 volScalarField& T = _T();
362 volScalarField& lambda = _lambda();
363 volScalarField& deltaT = _deltaT();
364 gList.append(g.clone());
365 Tfield.append(T.clone());
366 lambdaField.append(lambda.clone());
367 deltaTfield.append(deltaT.clone());
368 sampleI++;
369 }
370
371 cgIter++;
372 }
373
374 return (0);
375}
376
378{
379 fvMesh& mesh = _mesh();
380 volScalarField& lambda = _lambda();
381 gradJ_L2norm = 0;
382 forAll (lambda.boundaryField()[hotSide_ind], faceI)
383 {
384 gradJ [faceI] = - lambda.boundaryField()[hotSide_ind][faceI];
385 gradJ_L2norm += gradJ[faceI] * gradJ[faceI] *
386 mesh.magSf().boundaryField()[hotSide_ind][faceI];
387 }
388 gradJ_L2norm = Foam::sqrt(gradJ_L2norm);
389 Info << "gradJ L2norm = " << gradJ_L2norm << endl;
390}
391
393{
394 fvMesh& mesh = _mesh();
395 gamma = 0.0;
396 scalar gammaNum = 0;
397 gammaNum = gradJ_L2norm;
398
399 if (cgIter > 0)
400 {
401 //reduce(gamma, sumOp<double>());
402 gamma = gammaNum / gamma_den;
403 Info << "gamma = " << gamma << endl;
404 }
405
406 P = gradJ + gamma * P; //Updating P
407 gamma_den = gammaNum;
408 //reduce(gamma_den, sumOp<double>());
409}
410
412{
413 if (interpolation)
414 {
415 List<scalar> temp = interpolationPlane.Tdiff *
417 beta = 0.0;
419 {
420 beta += interpolationPlane.cellVol [cellI] * temp [cellI];
421 }
422 //reduce(beta, sumOp<double>());
424 scalar betaDiv = 0.0;
426 {
427 betaDiv += interpolationPlane.cellVol [cellI] * temp [cellI];
428 }
429 //reduce(betaDiv, sumOp<double>());
430 beta = beta / betaDiv;
431 temp.clear();
432 }
433 else
434 {
435 beta = Tdiff.dot(Tsens);
436 //reduce(beta, sumOp<double>());
437 double betaDiv = Tsens.dot(Tsens);
438 //reduce(betaDiv, sumOp<double>());
439 beta = beta / betaDiv;
440 }
441
442 Info << "beta = " << beta << endl;
443}
444
449
450
452{
453 double Jold = J;
454
455 if (interpolation)
456 {
457 List<scalar> sqTdiff;
459 J = 0.0;
460 forAll(sqTdiff, cellI)
461 {
462 J += 0.5 * sqTdiff[cellI] * interpolationPlane.cellVol[cellI];
463 }
464 sqTdiff.clear();
465 }
466 else
467 {
468 J = 0.5 * Tdiff.dot(Tdiff);
469 }
470
471 //reduce(J, sumOp<double>());
472 Info << "J = " << J << endl;
473
474 if (J <= Jtol)
475 {
476 Info << "Convergence reached in " << cgIter << " iterations" << endl;
477 return (1);
478 }
479 else if (Foam::mag((Jold - J) / J) <= JtolRel)
480 {
481 Info << "Relative tolerance criteria meet in " << cgIter << " iterations" <<
482 endl;
483 Info << "|Jold - J| / |J| = " << Foam::mag((Jold - J) / J) << endl;
484 return (1);
485 }
486 else
487 {
488 return (0);
489 }
490}
491
492
493
494int inverseLaplacianProblem_CG::isInPlane(double cx, double cy, double cz,
495 Foam::vector thermocoupleCellDim)
496{
497 return (cx >= interpolationPlane.minX - thermocoupleCellDim[0] / 4 &&
498 cy >= interpolationPlane.Y - thermocoupleCellDim[1] / 4 &&
499 cy <= interpolationPlane.Y + thermocoupleCellDim[1] / 4 &&
500 cz >= interpolationPlane.minZ - thermocoupleCellDim[2] / 4 &&
501 cx <= interpolationPlane.maxX + thermocoupleCellDim[0] / 4 &&
502 cz <= interpolationPlane.maxZ + thermocoupleCellDim[2] / 4
503 );
504}
505
507 const char* folder)
508{
509 fvMesh& mesh = _mesh();
510 volScalarField& T = _T();
511 volScalarField& lambda = _lambda();
512 volScalarField& deltaT = _deltaT();
513 autoPtr<volScalarField> gVolField_
514 (
515 new volScalarField("g", T)
516 );
517 volScalarField& gVolField = gVolField_();
519 //Access the mesh information for the boundary
520 const polyPatch& cPatch = mesh.boundaryMesh()[hotSide_ind];
521 //List of cells close to a boundary
522 const labelUList& faceCells = cPatch.faceCells();
523 forAll(cPatch, faceI)
524 {
525 //id of the owner cell having the face
526 label faceOwner = faceCells[faceI] ;
527 gVolField[faceOwner] = g[faceI];
528 }
529 ITHACAstream::exportSolution(T, std::to_string(folderNumber + 1), folder,
530 "T");
531 ITHACAstream::exportSolution(lambda, std::to_string(folderNumber + 1), folder,
532 "lambda");
533 ITHACAstream::exportSolution(deltaT, std::to_string(folderNumber + 1), folder,
534 "deltaT");
535 ITHACAstream::exportSolution(gVolField, std::to_string(folderNumber + 1),
536 folder, "g");
537}
538
539//Interpolates the values in Tmeas on the interpolation plane defined in readThermocouples()
541{
542 fvMesh& mesh = _mesh();
543 volScalarField& T = _T();
544 DataTable thermocouplesSamples;
545 DenseVector x(2);
546 unsigned i = 0;
547
549 {
552
553 //find the dimensions of the cells containing the thermocouples
554 //I am assuming all thermocouples are a the same y coordinate
555 //I am assuming all cells have the same dimensions
556 if ( Pstream::master() == true )
557 {
558 unsigned flag = 1;
559
560 while (flag == 1)
561 {
562 if (thermocouplesCellID [i] != -1)
563 {
564 labelList pLabels(mesh.cells()[thermocouplesCellID [i]].labels(mesh.faces()));
565 pointField pLocal(pLabels.size(), Foam::vector::zero);
567 mesh.points(),
568 mesh.cells()[thermocouplesCellID [i]],
569 pLabels,
570 pLocal);
571 flag = 0;
572 }
573
574 i++;
575 }
576 }
577
578 //reduce(interpolationPlane.thermocoupleCellDim, sumOp<vector>());
579 forAll(thermocouplesCellID, thermocoupleI)
580 {
581 if (thermocouplesCellID[thermocoupleI] != -1)
582 {
583 interpolationPlane.thermocoupleX [thermocoupleI] =
584 mesh.C()[thermocouplesCellID [thermocoupleI]].component(0);
585 interpolationPlane.thermocoupleZ [thermocoupleI] =
586 mesh.C()[thermocouplesCellID [thermocoupleI]].component(2);
587 }
588 else
589 {
590 interpolationPlane.thermocoupleX [thermocoupleI] = 0;
591 interpolationPlane.thermocoupleZ [thermocoupleI] = 0;
592 }
593 }
594 //reduce(interpolationPlane.thermocoupleX, sumOp<List<scalar>>());
595 //reduce(interpolationPlane.thermocoupleZ, sumOp<List<scalar>>());
596 }
597
598 forAll(thermocouplesCellID, thermocoupleI)
599 {
600 x(0) = interpolationPlane.thermocoupleX[thermocoupleI];
601 x(1) = interpolationPlane.thermocoupleZ[thermocoupleI];
602 thermocouplesSamples.addSample(x, Tmeas(thermocoupleI));
603 }
604 std::cout << Tmeas << std::endl;
605 RBFSpline rbfspline(thermocouplesSamples, RadialBasisFunctionType::GAUSSIAN);
606 auto inPlaneCellID = 0;
607 forAll(T.internalField(), cellI)
608 {
609 auto cx = mesh.C()[cellI].component(Foam::vector::X);
610 auto cy = mesh.C()[cellI].component(Foam::vector::Y);
611 auto cz = mesh.C()[cellI].component(Foam::vector::Z);
612
614 {
616 {
617 auto planeSize = interpolationPlane.cellID.size() + 1;
618 interpolationPlane.cellID.resize (planeSize);
619 interpolationPlane.Tmeas.resize (planeSize);
620 interpolationPlane.Tdirect.resize(planeSize);
621 interpolationPlane.Tdiff.resize (planeSize);
622 interpolationPlane.Tsens.resize (planeSize);
623 interpolationPlane.cellVol.resize(planeSize);
624 x(0) = cx;
625 x(1) = cz;
626 interpolationPlane.cellID [planeSize - 1] = cellI;
627 interpolationPlane.Tmeas [planeSize - 1] =
628 rbfspline.eval(x);
629 interpolationPlane.cellVol[planeSize - 1] =
630 mesh.V()[cellI];
631 }
632 }
633 else
634 {
636 {
637 x(0) = cx;
638 x(1) = cz;
639 interpolationPlane.Tmeas [inPlaneCellID] =
640 rbfspline.eval(x);
641 inPlaneCellID++;
642 }
643 }
644 }
646}
647
648//Interpolates the values in Tmeas on the interpolation plane defined in readThermocouples()
650 DenseMatrix& RBFweights, DenseMatrix& RBFbasis)
651{
652 fvMesh& mesh = _mesh();
653 volScalarField& T = _T();
654 DataTable thermocouplesSamples;
655 DenseVector x(2);
656 unsigned i = 0;
657
659 {
662
663 //find the dimensions of the cells containing the thermocouples
664 //I am assuming all thermocouples are a the same y coordinate
665 //I am assuming all cells have the same dimensions
666 if ( Pstream::master() == true )
667 {
668 unsigned flag = 1;
669
670 while (flag == 1)
671 {
672 if (thermocouplesCellID [i] != -1)
673 {
674 labelList pLabels(mesh.cells()[thermocouplesCellID [i]].labels(mesh.faces()));
675 pointField pLocal(pLabels.size(), Foam::vector::zero);
677 mesh.points(),
678 mesh.cells()[thermocouplesCellID [i]],
679 pLabels,
680 pLocal);
681 flag = 0;
682 }
683
684 i++;
685 }
686 }
687
688 //reduce(interpolationPlane.thermocoupleCellDim, sumOp<vector>());
689 forAll(thermocouplesCellID, thermocoupleI)
690 {
691 if (thermocouplesCellID[thermocoupleI] != -1)
692 {
693 interpolationPlane.thermocoupleX [thermocoupleI] =
694 mesh.C()[thermocouplesCellID [thermocoupleI]].component(0);
695 interpolationPlane.thermocoupleZ [thermocoupleI] =
696 mesh.C()[thermocouplesCellID [thermocoupleI]].component(2);
697 }
698 else
699 {
700 interpolationPlane.thermocoupleX [thermocoupleI] = 0;
701 interpolationPlane.thermocoupleZ [thermocoupleI] = 0;
702 }
703 }
704 //reduce(interpolationPlane.thermocoupleX, sumOp<List<scalar>>());
705 //reduce(interpolationPlane.thermocoupleZ, sumOp<List<scalar>>());
706 }
707
708 forAll(thermocouplesCellID, thermocoupleI)
709 {
710 x(0) = interpolationPlane.thermocoupleX[thermocoupleI];
711 x(1) = interpolationPlane.thermocoupleZ[thermocoupleI];
712 thermocouplesSamples.addSample(x, Tmeas(thermocoupleI));
713 }
714 std::cout << Tmeas << std::endl;
715 RBFSpline rbfspline(thermocouplesSamples, RadialBasisFunctionType::GAUSSIAN);
716 auto inPlaneCellID = 0;
717 forAll(T.internalField(), cellI)
718 {
719 auto cx = mesh.C()[cellI].component(Foam::vector::X);
720 auto cy = mesh.C()[cellI].component(Foam::vector::Y);
721 auto cz = mesh.C()[cellI].component(Foam::vector::Z);
722
724 {
726 {
727 auto planeSize = interpolationPlane.cellID.size() + 1;
728 interpolationPlane.cellID.resize (planeSize);
729 interpolationPlane.Tmeas.resize (planeSize);
730 interpolationPlane.Tdirect.resize(planeSize);
731 interpolationPlane.Tdiff.resize (planeSize);
732 interpolationPlane.Tsens.resize (planeSize);
733 interpolationPlane.cellVol.resize(planeSize);
734 x(0) = cx;
735 x(1) = cz;
736 interpolationPlane.cellID [planeSize - 1] = cellI;
737 interpolationPlane.Tmeas [planeSize - 1] =
738 rbfspline.eval(x);
739 interpolationPlane.cellVol[planeSize - 1] =
740 mesh.V()[cellI];
741 }
742 }
743 else
744 {
746 {
747 x(0) = cx;
748 x(1) = cz;
749 interpolationPlane.Tmeas [inPlaneCellID] =
750 rbfspline.eval(x);
751 inPlaneCellID++;
752 }
753 }
754 }
756}
757
759{
760 volScalarField& T = _T();
761
762 if (interpolation)
763 {
765 {
767 T.internalField()[interpolationPlane.cellID [cellI]];
768 }
771 }
772 else
773 {
775 Tdiff = Tdirect - Tmeas;
776 }
777}
778
780{
781 Info << "\nResetting time, mesh and fields: " << fieldName << "\n" << endl;
782 _simple.clear();
783
784 if (fieldName == "T" || fieldName == "all")
785 {
786 _T.clear();
787 }
788
789 if (fieldName == "lambda" || fieldName == "all")
790 {
791 _lambda.clear();
792 }
793
794 if (fieldName == "deltaT" || fieldName == "all")
795 {
796 _deltaT.clear();
797 }
798
799 argList& args = _args();
800 Time& runTime = _runTime();
801 //Reinitializing runTime
802 instantList Times = runTime.times();
803 runTime.setTime(Times[1], 1);
804 Foam::fvMesh& mesh = _mesh();
805 _simple = autoPtr<simpleControl>
806 (
807 new simpleControl
808 (
809 mesh
810 )
811 );
812
813 if (fieldName == "T" || fieldName == "all")
814 {
815 //Info << "ReReading field T\n" << endl;
816 _T = autoPtr<volScalarField>
817 (
818 new volScalarField
819 (
820 IOobject
821 (
822 "T",
823 runTime.timeName(),
824 mesh,
825 IOobject::MUST_READ,
826 IOobject::AUTO_WRITE
827 ),
828 mesh
829 )
830 );
831 }
832
833 if (fieldName == "deltaT" || fieldName == "all")
834 {
835 //Info << "ReReading field deltaT\n" << endl;
836 _deltaT = autoPtr<volScalarField>
837 (
838 new volScalarField
839 (
840 IOobject
841 (
842 "deltaT",
843 runTime.timeName(),
844 mesh,
845 IOobject::MUST_READ,
846 IOobject::AUTO_WRITE
847 ),
848 mesh
849 )
850 );
851 }
852
853 if (fieldName == "lambda" || fieldName == "all")
854 {
855 //Info << "ReReading field lambda\n" << endl;
856 _lambda = autoPtr<volScalarField>
857 (
858 new volScalarField
859 (
860 IOobject
861 (
862 "lambda",
863 runTime.timeName(),
864 mesh,
865 IOobject::MUST_READ,
866 IOobject::AUTO_WRITE
867 ),
868 mesh
869 )
870 );
871 }
872
873 Info << "Ready for new computation" << endl;
874}
875
876
877//void inverseLaplacianProblem_CG::offlineSolve()
878//{
879// volScalarField& T = _T();
880// volScalarField& lambda = _lambda();
881// volScalarField& deltaT = _deltaT();
882//
883// if (offline)
884// {
885// ITHACAstream::read_fields(Tfield, T, "./ITHACAoutput/Offline/");
886// ITHACAstream::read_fields(lambdaField, lambda, "./ITHACAoutput/Offline/");
887// ITHACAstream::read_fields(deltaTfield, deltaT, "./ITHACAoutput/Offline/");
888// }
889// else
890// {
891// offlinePhase = 1;
892// set_g();
893// readThermocouples();
894// set_valueFraction();
895// saveSolInLists = 0;
896// offlineSolutionI = 1;
897//
898// for (label i = 0; i < muSamples.cols(); i++)
899// {
900// Info << "Sample number " << i + 1 << endl;
901// Tmeas = muSamples.col(i);
902//
903// if (interpolation)
904// {
905// Info << "Interpolating thermocouples measurements in the " <<
906// "plane defined by them" << endl;
907// thermocouplesInterpolation();
908// }
909// else
910// {
911// Info << "NOT interpolating thermocouples measurements" << endl;
912// }
913//
914// if (!conjugateGradient())
915// {
916// Info << "Conjugate gradient method did not converge" <<
917// "Exiting" << endl;
918// exit(10);
919// }
920// }
921//
922// saveSolInLists = 0;
923// offlinePhase = 0;
924// }
925//}
926
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
TEqn solve()
int cgIterMax
Maximum CG iterations.
Eigen::VectorXd Tsens
Vector of solutions of the sensitivity problem at the thermocouples points.
Eigen::MatrixXd Jlist
Vector to store the const function J.
void set_valueFraction()
Set valueFraction list values for Robin condition.
PtrList< volScalarField > lambdaField
List of adjoint solutions.
double gamma
Conjugate coefficient.
bool interpolationPlaneDefined
Flag for definition of interpolation plane.
bool saveSolInLists
Flag to save solutions in lists.
double Jtol
Absolute stopping criterion for the CG.
void defineThermocouplesPlane()
Identifies the plane defined by the thermocouples.
void assignSensitivityBC()
Set BC and IF of the sensitivity problem.
thermocouplesPlane interpolationPlane
Interpolation plane.
void thermocouplesInterpolation()
Interpolates the thermocouples measurements in the plane defined in readThermocouples() using radial ...
int isInPlane(double cx, double cy, double cz, Foam::vector thermocoupleCellDim)
Checks if a cell crosses the interpolation plane.
int cgIter
Conjugate Gradient (CG) interations counter.
double gamma_den
Denoinator of the conjugate coefficient.
void sensibilitySolAtThermocouplesLocations()
Fill the Foam::vector containing the values of the sensitivity solution at the thermocouples location...
void assignAdjointBC()
Set BC of the adjoint problem.
void updateHeatFlux()
Updates the heat flux in the conjugate gradient iterations.
int conjugateGradient()
Conjugate gradient method.
void writeFields(label folderNumber, const char *folder)
Writes fields to file.
List< scalar > gradJ
Gradient of the cost function.
void solve(const char *problemID)
Solve for direct and sensitivity.
bool interpolation
Flag for interpolation of the temperature measurements.
List< scalar > P
Search direction.
void solveSensitivity()
Solve sensibility problem.
double gradJ_L2norm
L2 norm of the gradient of J.
PtrList< volScalarField > Tfield
List of temperature solutions.
double JtolRel
Relative stopping criterion for the CG.
autoPtr< volScalarField > _deltaT
Sensibility temperature field.
void computeGradJ()
Computes the gradient of cost function J and its L2 norm.
void solveAdjoint()
Solve adjoint problem.
void computeSearchStep()
Compute the search step beta.
PtrList< volScalarField > deltaTfield
List of sensitivity solutions.
volScalarField assignAdjointBCandSource()
Assign the BC for the adjoint problem and returs the source field.
List< scalar > valueFractionAdj
Value fraction for the adjoint Robin BC.
autoPtr< volScalarField > _lambda
Adjoint field.
void searchDirection()
Computes the search direction P.
void differenceBetweenDirectAndMeasure()
Computes the difference between direct problem solution and measure.
int conjugateGradientConvergenceCheck()
Convergence cher for the conjugate gradient method.
Implementation of a inverse Laplacian problem .
Eigen::VectorXd Tdiff
Difference between computed and measured temperatures at the thermocouples.
label hotSide_ind
Index of the hotSide patch.
dimensionedScalar DT
Dummy thermal conductivity with unitary value.
double H
Heat transfer coefficient [W/(m2 K)].
autoPtr< simpleControl > _simple
simpleControl
List< scalar > g
Heat flux at hotSide [W/m2].
double J
Cost funtion [K^2].
Eigen::VectorXd Tdirect
Vector of computed temperatures at the thermocouples locations [K].
List< scalar > homogeneousBCcoldSide
List of zeros of the size of coldSide patch.
List< scalar > valueFraction
Value fraction for the Robin BC.
scalar homogeneousBC
Homogenerous BC.
Eigen::VectorXd Tmeas
Vector of measured temperatures at the thermocouples locations [K].
List< List< scalar > > gList
List of boundary heat fluxes.
Foam::vector cellDim(const faceList &ff, const pointField &pp, const cell &cc, labelList pLabels, pointField pLocal)
Compute maximum cell dimension in x, y and z.
void set_g()
Set the right g size and fills it with zeros.
label coldSide_ind
Index of the coldSide patch.
autoPtr< fvMesh > _mesh
Mesh.
Eigen::VectorXd fieldValueAtThermocouples(volScalarField &field)
Interpolates the field value at the thermocouples points.
List< scalar > refGrad
Reference gradient for the Robin BC.
List< int > thermocouplesCellProc
List of incedes of the processors contining each thermocouple.
int thermocouplesNum
Number of thermocouples.
double k
Thermal diffusivity [W/(m K)].
void solveDirect()
Solve direct problem.
autoPtr< volScalarField > _T
Temperature field.
List< int > thermocouplesCellID
List of cells indices containing a thermocouple.
autoPtr< argList > _args
argList
volScalarField & T
Definition createT.H:46
Header file of the inverseLaplacianProblem_CG class.
fvScalarMatrix & TEqn
Definition TEqn.H:15
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
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 assignBC(GeometricField< scalar, fvPatchField, volMesh > &s, label BC_ind, double value)
Assign uniform Boundary Condition to a volScalarField.
Eigen::VectorXd boudaryFaceToCellDistance(fvMesh &mesh, label BC_ind)
Compute the distance between the boundary face center and the boundary cell center.
simpleControl simple(mesh)
label i
Definition pEqn.H:46
List< scalar > thermocoupleZ
List< scalar > thermocoupleX
Foam::vector thermocoupleCellDim