25class gsVisitorThermoBoundary
28 gsVisitorThermoBoundary(boxSide s,
29 const gsFunctionSet<T> & temperatureField_)
31 temperatureField(temperatureField_) {}
33 void initialize(
const gsBasisRefs<T>& basisRefs,
35 const gsOptionList & options,
40 dim = basisRefs.front().dim();
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 ) );
53 globalIndices.resize(dim);
54 blockNumbers.resize(dim);
57 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
58 const gsGeometry<T> & geo,
59 const gsMatrix<T> & quNodes)
73 temperatureField.piece(patch).eval_into(quNodes,tempValues);
75 temperatureField.eval_into(md.values[0],tempValues);
77 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
78 N_D = localIndicesDisp.rows();
80 basisRefs.front().eval_into(quNodes,basisValuesDisp);
83 inline void assemble(gsDomainIterator<T> & element,
84 const gsVector<T> & quWeights)
88 localRhs.setZero(dim*N_D,1);
90 for (
index_t q = 0; q < quWeights.rows(); ++q)
94 outerNormal(md, q, side, unormal);
96 const T weight = thermalExpCoef*(2*mu+dim*lambda)*quWeights[q] * (tempValues.at(q) - initTemp);
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) ;
103 inline void localToGlobal(
const int patchIndex,
104 const std::vector<gsMatrix<T> > & eliminatedDofs,
105 gsSparseSystem<T> & system)
109 for (
short_t d = 0; d < dim; ++d)
111 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
112 blockNumbers.at(d) = d;
115 system.pushToRhs(localRhs,globalIndices,blockNumbers);
122 T lambda, mu, thermalExpCoef, initTemp;
126 gsMatrix<T> localRhs;
128 gsMatrix<index_t> localIndicesDisp;
132 gsMatrix<T> basisValuesDisp;
136 const gsFunctionSet<T> & temperatureField;
141 gsMatrix<T> tempGrads;
144 gsMatrix<T> tempValues;
149 std::vector< gsMatrix<index_t> > globalIndices;
150 gsVector<index_t> blockNumbers;
#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