PLaSK library
Loading...
Searching...
No Matches
therm2d.hpp
Go to the documentation of this file.
1
/*
2
* This file is part of PLaSK (https://plask.app) by Photonics Group at TUL
3
* Copyright (c) 2022 Lodz University of Technology
4
*
5
* This program is free software: you can redistribute it and/or modify
6
* it under the terms of the GNU General Public License as published by
7
* the Free Software Foundation, version 3.
8
*
9
* This program is distributed in the hope that it will be useful,
10
* but WITHOUT ANY WARRANTY; without even the implied warranty of
11
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12
* GNU General Public License for more details.
13
*/
14
#ifndef PLASK__MODULE_THERMAL_THERM2D_H
15
#define PLASK__MODULE_THERMAL_THERM2D_H
16
17
#include <
plask/plask.hpp
>
18
#include <
plask/common/fem.hpp
>
19
20
#include "
common.hpp
"
21
22
namespace
plask
{
namespace
thermal {
namespace
tstatic {
23
27
template
<
typename
Geometry2DType>
28
struct
PLASK_SOLVER_API
ThermalFem2DSolver
:
public
FemSolverWithMaskedMesh
<Geometry2DType, RectangularMesh<2>> {
29
protected
:
30
int
loopno
;
31
double
maxT
;
32
double
toterr
;
34
35
DataVector<double>
temperatures
;
36
37
DataVector<double>
thickness
;
38
39
DataVector<Vec<2, double>
>
fluxes
;
40
42
void
setMatrix
(
FemMatrix
& A,
43
DataVector<double>
& B,
44
const
BoundaryConditionsWithMesh
<
RectangularMesh<2>::Boundary
,
double
>& btemperature,
45
const
BoundaryConditionsWithMesh
<
RectangularMesh<2>::Boundary
,
double
>& bheatflux,
46
const
BoundaryConditionsWithMesh
<
RectangularMesh<2>::Boundary
,
Convection
>& bconvection,
47
const
BoundaryConditionsWithMesh
<
RectangularMesh<2>::Boundary
,
Radiation
>& bradiation);
48
50
template
<
typename
MatrixT> MatrixT
makeMatrix
();
51
53
void
saveHeatFluxes();
// [W/m^2]
54
56
void
onInitialize()
override
;
57
59
void
onInvalidate()
override
;
60
61
public
:
62
// Boundary conditions
63
BoundaryConditions<RectangularMesh<2>::Boundary
,
double
>
64
temperature_boundary
;
65
BoundaryConditions<RectangularMesh<2>::Boundary
,
double
>
66
heatflux_boundary
;
67
BoundaryConditions<RectangularMesh<2>::Boundary
,
Convection
>
convection_boundary
;
68
BoundaryConditions<RectangularMesh<2>::Boundary
,
Radiation
>
radiation_boundary
;
69
70
typename
ProviderFor<Temperature, Geometry2DType>::Delegate
outTemperature
;
71
72
typename
ProviderFor<HeatFlux, Geometry2DType>::Delegate
outHeatFlux
;
73
74
typename
ProviderFor<ThermalConductivity, Geometry2DType>::Delegate
outThermalConductivity
;
75
76
ReceiverFor<Heat, Geometry2DType>
inHeat
;
77
78
double
maxerr
;
79
double
inittemp
;
80
86
double
compute(
int
loops = 1);
87
89
double
getErr
()
const
{
return
toterr; }
90
91
void
loadConfiguration(
XMLReader
& source,
92
Manager
& manager)
override
;
// for solver configuration (see: *.xpl file with structures)
93
94
ThermalFem2DSolver
(
const
std::string& name =
""
);
95
96
std::string
getClassName
()
const override
;
97
98
~ThermalFem2DSolver
();
99
100
public
:
101
struct
ThermalConductivityData
:
public
LazyDataImpl
<Tensor2<double>> {
102
const
ThermalFem2DSolver
*
solver
;
103
shared_ptr<const MeshD<2>>
dest_mesh
;
104
InterpolationFlags
flags
;
105
LazyData<double>
temps
;
106
ThermalConductivityData
(
const
ThermalFem2DSolver
* solver,
const
shared_ptr<
const
MeshD<2>
>& dst_mesh);
107
Tensor2<double>
at(std::size_t i)
const override
;
108
std::size_t size()
const override
;
109
};
110
111
const
LazyData<double>
getTemperatures(
const
shared_ptr<
const
MeshD<2>
>& dest_mesh,
InterpolationMethod
method)
const
;
112
113
const
LazyData<Vec<2>
> getHeatFluxes(
const
shared_ptr<
const
MeshD<2>
>& dest_mesh,
InterpolationMethod
method);
114
115
const
LazyData<Tensor2<double>
> getThermalConductivity(
const
shared_ptr<
const
MeshD<2>
>& dest_mesh,
InterpolationMethod
method);
116
};
117
118
}}}
// namespace plask::thermal::tstatic
119
120
#endif
solvers
thermal
static
therm2d.hpp
Generated by
1.9.8