G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorThermo.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 gsVisitorThermo
26{
27public:
28 gsVisitorThermo(const gsFunctionSet<T> & temperatureField_)
29 : dim(0), lambda(0), mu(0), thermalExpCoef(0), N_D(0), patch(0),
30 temperatureField(temperatureField_), paramTemp(0) {}
31
32 void initialize(const gsBasisRefs<T> & basisRefs,
33 const index_t patchIndex,
34 const gsOptionList & options,
35 gsQuadRule<T> & rule)
36 {
37 // parametric dimension of the first displacement component
38 dim = basisRefs.front().dim();
39 // a quadrature rule is defined by the basis for the first displacement component.
40 rule = gsQuadrature::get(basisRefs.front(), options);
41 // saving necessary info
42 patch = patchIndex;
43 paramTemp = options.getSwitch("ParamTemp");
44 thermalExpCoef = options.getReal("ThExpCoef");
45 T E = options.getReal("YoungsModulus");
46 T pr = options.getReal("PoissonsRatio");
47 lambda = E * pr / ( ( 1. + pr ) * ( 1. - 2. * pr ) );
48 mu = E / ( 2. * ( 1. + pr ) );
49 // resize containers for global indices
50 globalIndices.resize(dim);
51 blockNumbers.resize(dim);
52 }
53
54 inline void evaluate(const gsBasisRefs<T> & basisRefs,
55 const gsGeometry<T> & geo,
56 const gsMatrix<T> & quNodes)
57 {
58 // store quadrature points of the element for geometry evaluation
59 md.points = quNodes;
60 // NEED_VALUE to get points in the physical domain for evaluation of the temperature gradients
61 // NEED_MEASURE to get the Jacobian determinant values for integration
62 // NEED_GRAD_TRANSFORM to get the Jacobian matrix to transform temperature gradient from the parametric to physical domain
63 if (paramTemp)
65 else
66 md.flags = NEED_VALUE | NEED_MEASURE;
67 // Compute image of the quadrature points plus gradient, jacobian and other necessary data
68 geo.computeMap(md);
69 // Compute temperature gradients
70 if (paramTemp) // evaluate gradients in the parametric domain
71 temperatureField.piece(patch).deriv_into(quNodes,tempGrads);
72 else // evaluate gradients in the physical domain
73 temperatureField.piece(patch).deriv_into(md.values[0],tempGrads);
74 // find local indices of the displacement basis functions active on the element
75 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
76 N_D = localIndicesDisp.rows();
77 // Evaluate displacement basis functions on the element
78 basisRefs.front().eval_into(quNodes,basisValuesDisp);
79 }
80
81 inline void assemble(gsDomainIterator<T> & element,
82 const gsVector<T> & quWeights)
83 {
84 GISMO_UNUSED(element);
85 // Initialize local matrix/rhs
86 localRhs.setZero(dim*N_D, 1);
87 // Loop over the quadrature nodes
88 for (index_t q = 0; q < quWeights.rows(); ++q)
89 {
90 // Multiply quadrature weight by the geometry measure
91 const T weight = thermalExpCoef*(2*mu+dim*lambda)*quWeights[q]*md.measure(q);
92
93 if (paramTemp) // transform temperature gradients to the physical domain
94 {
95 // temperature gradient at one point in the physical domain, dim x 1
96 transformGradients(md,q,tempGrads,physGrad);
97 for (index_t d = 0; d < dim; ++d)
98 localRhs.middleRows(d*N_D,N_D).noalias() -= weight * physGrad(d,0) * basisValuesDisp.col(q);
99 }
100 else // use temperature gradients as they are
101 for (index_t d = 0; d < dim; ++d)
102 localRhs.middleRows(d*N_D,N_D).noalias() -= weight * tempGrads(d,q) * basisValuesDisp.col(q);
103 }
104 }
105
106 inline void localToGlobal(const int patchIndex,
107 const std::vector<gsMatrix<T> > & eliminatedDofs,
108 gsSparseSystem<T> & system)
109 {
110 GISMO_UNUSED(eliminatedDofs);
111 // computes global indices for displacement components
112 for (short_t d = 0; d < dim; ++d)
113 {
114 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
115 blockNumbers.at(d) = d;
116 }
117 // push to global system
118 system.pushToRhs(localRhs,globalIndices,blockNumbers);
119 }
120
121protected:
122 // problem info
123 short_t dim;
124 // Lame and thermal expansion coefficients
125 T lambda, mu, thermalExpCoef;
126 // geometry mapping
127 gsMapData<T> md;
128 // local components of the global linear system
129 gsMatrix<T> localRhs;
130 // local indices (at the current patch) of the displacement basis functions active at the current element
131 gsMatrix<index_t> localIndicesDisp;
132 // number of displacement basis functions active at the current element
133 index_t N_D;
134 // values of displacement basis functions at quadrature points at the current element stored as a N_D x numQuadPoints matrix;
135 gsMatrix<T> basisValuesDisp;
136 // Temperature info
137 index_t patch;
138 const gsFunctionSet<T> & temperatureField;
139 // true if temperature field is defined in the parametric domain; false if in the physical
140 bool paramTemp;
141 // temperature gradient evaluated at the quadrature points or at their images in the physical domain;
142 // stored as a dim x numQuadPoints matrix
143 gsMatrix<T> tempGrads;
144
145 // all temporary matrices defined here for efficiency
146 gsMatrix<T> physGrad;
147 // containers for global indices
148 std::vector< gsMatrix<index_t> > globalIndices;
149 gsVector<index_t> blockNumbers;
150
151}; //class definition ends
152
153} // namespace ends
#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_GRAD_TRANSFORM
Gradient transformation matrix.
Definition gsForwardDeclarations.h:77
@ 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