G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorThermoBoundary.h
Go to the documentation of this file.
1
16#pragma once
17
19#include <gsCore/gsFuncData.h>
20
21namespace gismo
22{
23
24template <class T>
25class gsVisitorThermoBoundary
26{
27public:
28 gsVisitorThermoBoundary(boxSide s,
29 const gsFunctionSet<T> & temperatureField_)
30 : side(s),
31 temperatureField(temperatureField_) {}
32
33 void initialize(const gsBasisRefs<T>& basisRefs,
34 const index_t patchIndex,
35 const gsOptionList & options,
36 gsQuadRule<T> & rule)
37 {
38 GISMO_UNUSED(patchIndex);
39 // parametric dimension of the first displacement component
40 dim = basisRefs.front().dim();
41 // a quadrature rule is defined by the basis for the first displacement component.
42 rule = gsQuadrature::get(basisRefs.front(), options);
43 // saving necessary info
44 patch = patchIndex;
45 paramTemp = options.getSwitch("ParamTemp");
46 initTemp = options.getReal("InitTemp");
47 thermalExpCoef = options.getReal("ThExpCoef");
48 T E = options.getReal("YoungsModulus");
49 T pr = options.getReal("PoissonsRatio");
50 lambda = E * pr / ( ( 1. + pr ) * ( 1. - 2. * pr ) );
51 mu = E / ( 2. * ( 1. + pr ) );
52 // resize containers for global indices
53 globalIndices.resize(dim);
54 blockNumbers.resize(dim);
55 }
56
57 inline void evaluate(const gsBasisRefs<T> & basisRefs,
58 const gsGeometry<T> & geo,
59 const gsMatrix<T> & quNodes)
60 {
61 // store quadrature points of the element for geometry evaluation
62 md.points = quNodes;
63 // NEED_VALUE to get points in the physical domain for evaluation of the temperature gradients
64 // NEED_MEASURE to get the Jacobian determinant values for integration
65 if (paramTemp)
66 md.flags = NEED_MEASURE;
67 else
68 md.flags = NEED_VALUE | NEED_MEASURE;
69 // Compute image of the quadrature points plus gradient, jacobian and other necessary data
70 geo.computeMap(md);
71 // Evaluate temperature
72 if (paramTemp) // evaluate temperature in the parametric domain
73 temperatureField.piece(patch).eval_into(quNodes,tempValues);
74 else // evaluate temperature in the physical domain
75 temperatureField.eval_into(md.values[0],tempValues);
76 // find local indices of the displacement basis functions active on the element
77 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
78 N_D = localIndicesDisp.rows();
79 // Evaluate displacement basis functions on the element
80 basisRefs.front().eval_into(quNodes,basisValuesDisp);
81 }
82
83 inline void assemble(gsDomainIterator<T> & element,
84 const gsVector<T> & quWeights)
85 {
86 GISMO_UNUSED(element);
87 // Initialize local matrix/rhs
88 localRhs.setZero(dim*N_D,1);
89 // loop over the quadrature nodes
90 for (index_t q = 0; q < quWeights.rows(); ++q)
91 {
92 // Compute the outer normal vector on the side
93 // normal length equals to the local area measure
94 outerNormal(md, q, side, unormal);
95 // Collect the factors here: quadrature weight, geometry measure and time factor
96 const T weight = thermalExpCoef*(2*mu+dim*lambda)*quWeights[q] * (tempValues.at(q) - initTemp);
97
98 for (index_t d = 0; d < dim; ++d)
99 localRhs.middleRows(d*N_D,N_D).noalias() += weight * unormal(d,0) * basisValuesDisp.col(q) ;
100 }
101 }
102
103 inline void localToGlobal(const int patchIndex,
104 const std::vector<gsMatrix<T> > & eliminatedDofs,
105 gsSparseSystem<T> & system)
106 {
107 GISMO_UNUSED(eliminatedDofs);
108 // computes global indices for displacement components
109 for (short_t d = 0; d < dim; ++d)
110 {
111 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
112 blockNumbers.at(d) = d;
113 }
114 // push to global system
115 system.pushToRhs(localRhs,globalIndices,blockNumbers);
116 }
117
118protected:
119 // problem info
120 short_t dim;
121 // Lame and thermal expansion coefficients, initial temperature
122 T lambda, mu, thermalExpCoef, initTemp;
123 // geometry mapping
124 gsMapData<T> md;
125 // local components of the global linear system
126 gsMatrix<T> localRhs;
127 // local indices (at the current patch) of the displacement basis functions active at the current element
128 gsMatrix<index_t> localIndicesDisp;
129 // number of displacement basis functions active at the current element
130 index_t N_D;
131 // values of displacement basis functions at quadrature points at the current element stored as a N_D x numQuadPoints matrix;
132 gsMatrix<T> basisValuesDisp;
133 // Temperature info
134 index_t patch;
135 boxSide side;
136 const gsFunctionSet<T> & temperatureField;
137 // true if temperature field is defined in the parametric domain; false if in the physical
138 bool paramTemp;
139 // temperature gradient evaluated at the quadrature points or at their images in the physical domain;
140 // stored as a dim x numQuadPoints matrix
141 gsMatrix<T> tempGrads;
142 // temperature evaluated at the quadrature points or at their images in the physical domain;
143 // stored as a 1 x numQuadPoints matrix
144 gsMatrix<T> tempValues;
145
146 // all temporary matrices defined here for efficiency
147 gsVector<T> unormal;
148 // containers for global indices
149 std::vector< gsMatrix<index_t> > globalIndices;
150 gsVector<index_t> blockNumbers;
151
152}; //class definition ends
153
154} // namespace ends
155
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
This object is a cache for computed values from an evaluator.
Creates a variety of quadrature rules.
The G+Smo namespace, containing all definitions for the library.
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ NEED_MEASURE
The density of the measure pull back.
Definition gsForwardDeclarations.h:76
static gsQuadRule< T > get(const gsBasis< T > &basis, const gsOptionList &options, short_t fixDir=-1)
Constructs a quadrature rule based on input options.
Definition gsQuadrature.h:51