25 class gsVisitorThermoBoundary
28 gsVisitorThermoBoundary(boxSide s,
29 const gsFunctionSet<T> & temperatureField_)
31 temperatureField(temperatureField_) {}
33 void initialize(
const gsBasisRefs<T>& basisRefs,
35 const gsOptionList & options,
39 dim = basisRefs.front().dim();
44 paramTemp = options.getSwitch(
"ParamTemp");
45 initTemp = options.getReal(
"InitTemp");
46 thermalExpCoef = options.getReal(
"ThExpCoef");
47 T E = options.getReal(
"YoungsModulus");
48 T pr = options.getReal(
"PoissonsRatio");
49 lambda = E * pr / ( ( 1. + pr ) * ( 1. - 2. * pr ) );
50 mu = E / ( 2. * ( 1. + pr ) );
52 globalIndices.resize(dim);
53 blockNumbers.resize(dim);
56 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
57 const gsGeometry<T> & geo,
58 const gsMatrix<T> & quNodes)
72 temperatureField.piece(patch).eval_into(quNodes,tempValues);
74 temperatureField.eval_into(md.values[0],tempValues);
76 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
77 N_D = localIndicesDisp.rows();
79 basisRefs.front().eval_into(quNodes,basisValuesDisp);
82 inline void assemble(gsDomainIterator<T> & element,
83 const gsVector<T> & quWeights)
86 localRhs.setZero(dim*N_D,1);
88 for (
index_t q = 0; q < quWeights.rows(); ++q)
92 outerNormal(md, q, side, unormal);
94 const T weight = thermalExpCoef*(2*mu+dim*lambda)*quWeights[q] * (tempValues.at(q) - initTemp);
96 for (
index_t d = 0; d < dim; ++d)
97 localRhs.middleRows(d*N_D,N_D).noalias() += weight * unormal(d,0) * basisValuesDisp.col(q) ;
101 inline void localToGlobal(
const int patchIndex,
102 const std::vector<gsMatrix<T> > & eliminatedDofs,
103 gsSparseSystem<T> & system)
106 for (
short_t d = 0; d < dim; ++d)
108 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
109 blockNumbers.at(d) = d;
112 system.pushToRhs(localRhs,globalIndices,blockNumbers);
119 T lambda, mu, thermalExpCoef, initTemp;
123 gsMatrix<T> localRhs;
125 gsMatrix<index_t> localIndicesDisp;
129 gsMatrix<T> basisValuesDisp;
133 const gsFunctionSet<T> & temperatureField;
138 gsMatrix<T> tempGrads;
141 gsMatrix<T> tempValues;
146 std::vector< gsMatrix<index_t> > globalIndices;
147 gsVector<index_t> blockNumbers;
#define short_t
Definition: gsConfig.h:35
#define index_t
Definition: gsConfig.h:32
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:48
Creates a variety of quadrature rules.
The density of the measure pull back.
Definition: gsForwardDeclarations.h:76
Value of the object.
Definition: gsForwardDeclarations.h:72
This object is a cache for computed values from an evaluator.