G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsVisitorThermoBoundary.h
Go to the documentation of this file.
1 
16 #pragma once
17 
19 #include <gsCore/gsFuncData.h>
20 
21 namespace gismo
22 {
23 
24 template <class T>
25 class gsVisitorThermoBoundary
26 {
27 public:
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  // parametric dimension of the first displacement component
39  dim = basisRefs.front().dim();
40  // a quadrature rule is defined by the basis for the first displacement component.
41  rule = gsQuadrature::get(basisRefs.front(), options);
42  // saving necessary info
43  patch = patchIndex;
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 ) );
51  // resize containers for global indices
52  globalIndices.resize(dim);
53  blockNumbers.resize(dim);
54  }
55 
56  inline void evaluate(const gsBasisRefs<T> & basisRefs,
57  const gsGeometry<T> & geo,
58  const gsMatrix<T> & quNodes)
59  {
60  // store quadrature points of the element for geometry evaluation
61  md.points = quNodes;
62  // NEED_VALUE to get points in the physical domain for evaluation of the temperature gradients
63  // NEED_MEASURE to get the Jacobian determinant values for integration
64  if (paramTemp)
65  md.flags = NEED_MEASURE;
66  else
67  md.flags = NEED_VALUE | NEED_MEASURE;
68  // Compute image of the quadrature points plus gradient, jacobian and other necessary data
69  geo.computeMap(md);
70  // Evaluate temperature
71  if (paramTemp) // evaluate temperature in the parametric domain
72  temperatureField.piece(patch).eval_into(quNodes,tempValues);
73  else // evaluate temperature in the physical domain
74  temperatureField.eval_into(md.values[0],tempValues);
75  // find local indices of the displacement basis functions active on the element
76  basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
77  N_D = localIndicesDisp.rows();
78  // Evaluate displacement basis functions on the element
79  basisRefs.front().eval_into(quNodes,basisValuesDisp);
80  }
81 
82  inline void assemble(gsDomainIterator<T> & element,
83  const gsVector<T> & quWeights)
84  {
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  // Compute the outer normal vector on the side
91  // normal length equals to the local area measure
92  outerNormal(md, q, side, unormal);
93  // Collect the factors here: quadrature weight, geometry measure and time factor
94  const T weight = thermalExpCoef*(2*mu+dim*lambda)*quWeights[q] * (tempValues.at(q) - initTemp);
95 
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) ;
98  }
99  }
100 
101  inline void localToGlobal(const int patchIndex,
102  const std::vector<gsMatrix<T> > & eliminatedDofs,
103  gsSparseSystem<T> & system)
104  {
105  // computes global indices for displacement components
106  for (short_t d = 0; d < dim; ++d)
107  {
108  system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
109  blockNumbers.at(d) = d;
110  }
111  // push to global system
112  system.pushToRhs(localRhs,globalIndices,blockNumbers);
113  }
114 
115 protected:
116  // problem info
117  short_t dim;
118  // Lame and thermal expansion coefficients, initial temperature
119  T lambda, mu, thermalExpCoef, initTemp;
120  // geometry mapping
121  gsMapData<T> md;
122  // local components of the global linear system
123  gsMatrix<T> localRhs;
124  // local indices (at the current patch) of the displacement basis functions active at the current element
125  gsMatrix<index_t> localIndicesDisp;
126  // number of displacement basis functions active at the current element
127  index_t N_D;
128  // values of displacement basis functions at quadrature points at the current element stored as a N_D x numQuadPoints matrix;
129  gsMatrix<T> basisValuesDisp;
130  // Temperature info
131  index_t patch;
132  boxSide side;
133  const gsFunctionSet<T> & temperatureField;
134  // true if temperature field is defined in the parametric domain; false if in the physical
135  bool paramTemp;
136  // temperature gradient evaluated at the quadrature points or at their images in the physical domain;
137  // stored as a dim x numQuadPoints matrix
138  gsMatrix<T> tempGrads;
139  // temperature evaluated at the quadrature points or at their images in the physical domain;
140  // stored as a 1 x numQuadPoints matrix
141  gsMatrix<T> tempValues;
142 
143  // all temporary matrices defined here for efficiency
144  gsVector<T> unormal;
145  // containers for global indices
146  std::vector< gsMatrix<index_t> > globalIndices;
147  gsVector<index_t> blockNumbers;
148 
149 }; //class definition ends
150 
151 } // namespace ends
152 
#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.