25class gsVisitorElasticityNeumann
29 gsVisitorElasticityNeumann(
const gsPde<T> & ,
30 const boundary_condition<T> & s)
31 : neumannFunction_ptr(s.function().get()),
32 patchSide(s.side()) {}
34 void initialize(
const gsBasisRefs<T> & basisRefs,
36 const gsOptionList & options,
41 dim = basisRefs.front().dim();
45 forceScaling = options.getReal(
"ForceScaling");
47 globalIndices.resize(dim);
48 blockNumbers.resize(dim);
51 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
52 const gsGeometry<T> & geo,
53 const gsMatrix<T> & quNodes)
63 neumannFunction_ptr->eval_into(md.values[0], neumannValues);
65 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
66 N_D = localIndicesDisp.rows();
68 basisRefs.front().eval_into(quNodes,basisValuesDisp);
71 inline void assemble(gsDomainIterator<T> & element,
72 const gsVector<T> & quWeights)
76 localRhs.setZero(dim*N_D,1);
78 for (
index_t q = 0; q < quWeights.rows(); ++q)
82 outerNormal(md, q, patchSide, unormal);
84 const T weight = quWeights[q] * unormal.norm() * forceScaling;
86 for (
short_t d = 0; d < dim; ++d)
87 localRhs.middleRows(d*N_D,N_D).noalias() += weight * neumannValues(d,q) * basisValuesDisp.col(q);
91 inline void localToGlobal(
const int patchIndex,
92 const std::vector<gsMatrix<T> > & eliminatedDofs,
93 gsSparseSystem<T> & system)
97 for (
short_t d = 0; d < dim; ++d)
99 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
100 blockNumbers.at(d) = d;
103 system.pushToRhs(localRhs,globalIndices,blockNumbers);
109 const gsFunctionSet<T> * neumannFunction_ptr;
115 gsMatrix<T> localRhs;
117 gsMatrix<index_t> localIndicesDisp;
121 gsMatrix<T> basisValuesDisp;
123 gsMatrix<T> neumannValues;
128 std::vector< gsMatrix<index_t> > globalIndices;
129 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