PLaSK library
Loading...
Searching...
No Matches
therm3d.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_TFem_H
15
#define PLASK__MODULE_THERMAL_TFem_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
struct
PLASK_SOLVER_API
ThermalFem3DSolver
:
public
FemSolverWithMaskedMesh
<Geometry3D, RectangularMesh<3>> {
28
29
protected
:
30
31
int
loopno
;
32
double
maxT
;
33
double
toterr
;
34
35
DataVector<double>
temperatures
;
36
37
DataVector<double>
thickness
;
38
39
DataVector<Vec<3,double>
>
fluxes
;
40
50
void
setMatrix(
FemMatrix
& A,
DataVector<double>
& B,
51
const
BoundaryConditionsWithMesh
<
RectangularMesh<3>::Boundary
,
double
>& btemperature,
52
const
BoundaryConditionsWithMesh
<
RectangularMesh<3>::Boundary
,
double
>& bheatflux,
53
const
BoundaryConditionsWithMesh
<
RectangularMesh<3>::Boundary
,
Convection
>& bconvection,
54
const
BoundaryConditionsWithMesh
<
RectangularMesh<3>::Boundary
,
Radiation
>& bradiation
55
);
56
60
template
<
typename
MatrixT>
61
void
applyBC
(MatrixT& A,
DataVector<double>
& B,
const
BoundaryConditionsWithMesh
<
RectangularMesh<3>::Boundary
,
double
>& btemperature);
62
64
template
<
typename
MatrixT>
65
MatrixT
makeMatrix
();
66
68
void
saveHeatFluxes();
// [W/m^2]
69
71
void
onInitialize()
override
;
72
74
void
onInvalidate()
override
;
75
76
public
:
77
78
double
inittemp
;
79
double
maxerr
;
80
81
// Boundary conditions
82
BoundaryConditions<RectangularMesh<3>::Boundary
,
double
>
temperature_boundary
;
83
BoundaryConditions<RectangularMesh<3>::Boundary
,
double
>
heatflux_boundary
;
84
BoundaryConditions<RectangularMesh<3>::Boundary
,
Convection
>
convection_boundary
;
85
BoundaryConditions<RectangularMesh<3>::Boundary
,
Radiation
>
radiation_boundary
;
86
87
typename
ProviderFor<Temperature,Geometry3D>::Delegate
outTemperature
;
88
89
typename
ProviderFor<HeatFlux,Geometry3D>::Delegate
outHeatFlux
;
90
91
typename
ProviderFor<ThermalConductivity,Geometry3D>::Delegate
outThermalConductivity
;
92
93
ReceiverFor<Heat,Geometry3D>
inHeat
;
94
100
double
compute(
int
loops=1);
101
102
void
loadConfiguration(
XMLReader
& source,
Manager
& manager)
override
;
// for solver configuration (see: *.xpl file with structures)
103
104
ThermalFem3DSolver
(
const
std::string& name=
""
);
105
106
std::string
getClassName
()
const override
{
return
"thermal.Static3D"
; }
107
108
~ThermalFem3DSolver
();
109
111
double
getErr
()
const
{
return
toterr; }
112
114
FemMatrixAlgorithm
getAlgorithm
()
const
{
return
algorithm; }
115
120
void
setAlgorithm(
FemMatrixAlgorithm
alg);
121
122
protected
:
123
124
struct
ThermalConductivityData
:
public
LazyDataImpl
<Tensor2<double>> {
125
const
ThermalFem3DSolver
*
solver
;
126
shared_ptr<const MeshD<3>>
dest_mesh
;
127
InterpolationFlags
flags
;
128
LazyData<double>
temps
;
129
ThermalConductivityData
(
const
ThermalFem3DSolver
* solver,
const
shared_ptr<
const
MeshD<3>
>& dst_mesh);
130
Tensor2<double>
at(std::size_t i)
const override
;
131
std::size_t size()
const override
;
132
};
133
134
const
LazyData<double>
getTemperatures(
const
shared_ptr<
const
MeshD<3>
>& dst_mesh,
InterpolationMethod
method)
const
;
135
136
const
LazyData<Vec<3>
> getHeatFluxes(
const
shared_ptr<
const
MeshD<3>
>& dst_mesh,
InterpolationMethod
method);
137
138
const
LazyData<Tensor2<double>
> getThermalConductivity(
const
shared_ptr<
const
MeshD<3>
>& dst_mesh,
InterpolationMethod
method);
139
};
140
141
}}}
//namespaces
142
143
#endif
solvers
thermal
static
therm3d.hpp
Generated by
1.9.8