G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsNsAssembler.hpp
Go to the documentation of this file.
1 
15 #pragma once
16 
18 
19 #include <gsElasticity/gsBasePde.h>
20 #include <gsUtils/gsPointGrid.h>
21 
22 // Element visitors
25 
26 namespace gismo
27 {
28 
29 template<class T>
31  gsMultiBasis<T> const & basisVel,
32  gsMultiBasis<T> const & basisPres,
33  gsBoundaryConditions<T> const & bconditions,
34  const gsFunction<T> & body_force)
35 {
36  // comment: same problems as in gsElasticityAssembler
37  gsPiecewiseFunction<T> rightHandSides;
38  rightHandSides.addPiece(body_force);
39  typename gsPde<T>::Ptr pde( new gsBasePde<T>(patches,bconditions,rightHandSides) );
40 
41  m_dim = body_force.targetDim();
42  for (short_t d = 0; d < m_dim; ++d)
43  m_bases.push_back(basisVel);
44  m_bases.push_back(basisPres);
45 
46  Base::initialize(pde, m_bases, defaultOptions());
47 }
48 
49 template <class T>
51 {
52  gsOptionList opt = Base::defaultOptions();
53  opt.addReal("Viscosity","Kinematic viscosity of the fluid",0.001);
54  opt.addReal("Density","Density of the fluid",1.);
55  opt.addReal("ForceScaling","Force scaling parameter",1.);
56  opt.addInt("Assembly","Type of the linear system to assemble",ns_assembly::newton_update);
57  return opt;
58 }
59 
60 template <class T>
62 {
63  // Pick up values from options
64  const T bdA = m_options.getReal("bdA");
65  const index_t bdB = m_options.getInt("bdB");
66  const T bdO = m_options.getReal("bdO");
67 
68  index_t deg = 0;
69  for (index_t d = 0; d < m_bases[0][0].dim(); ++d )
70  if (m_bases[0][0].degree(d) > deg)
71  deg = m_bases[0][0].degree(d);
72 
73  // m_dim velocity*velocity blocks + 1 pressure*velocity block (additioanal factor 2 for subgrid element)
74  index_t numElPerColumn = pow((bdA*deg+bdB),m_dim)*m_dim + pow((2*bdA*deg+bdB),m_dim);
75 
76  m_system.reserve(numElPerColumn*(1+bdO),1);
77 }
78 
79 template <class T>
81 {
82  GISMO_ENSURE(m_dim == m_pde_ptr->domain().parDim(), "The RHS dimension and the domain dimension don't match!");
83  GISMO_ENSURE(m_dim == 2 || m_dim == 3, "Only two- and three-dimenstion domains are supported!");
84 
85  std::vector<gsDofMapper> m_dofMappers(m_bases.size());
86  for (unsigned d = 0; d < m_bases.size(); d++)
87  m_bases[d].getMapper((dirichlet::strategy)m_options.getInt("DirichletStrategy"),
88  iFace::glue,m_pde_ptr->bc(),m_dofMappers[d],d,true);
89 
90  m_system = gsSparseSystem<T>(m_dofMappers, gsVector<index_t>::Ones(m_bases.size()));
91  reserve();
92 
93  for (unsigned d = 0; d < m_bases.size(); ++d)
94  Base::computeDirichletDofs(d);
95 }
96 
97 //--------------------- SYSTEM ASSEMBLY ----------------------------------//
98 
99 template<class T>
100 void gsNsAssembler<T>::assemble(bool saveEliminationMatrix)
101 {
102  m_system.matrix().setZero();
103  reserve();
104  m_system.rhs().setZero();
105 
106  gsVisitorStokes<T> visitor(*m_pde_ptr);
107  Base::template push<gsVisitorStokes<T> >(visitor);
108 
109  m_system.matrix().makeCompressed();
110 }
111 
112 template <class T>
113 bool gsNsAssembler<T>::assemble(const gsMatrix<T> & solutionVector,
114  const std::vector<gsMatrix<T> > & fixedDoFs)
115 {
116  gsMultiPatch<T> velocity, pressure;
117  constructSolution(solutionVector,fixedDoFs,velocity,pressure);
118  assemble(velocity,pressure);
119 
120  return true;
121 }
122 
123 template <class T>
125  const gsMultiPatch<T> & pressure)
126 {
127  m_system.matrix().setZero();
128  reserve();
129  m_system.rhs().setZero();
130 
131  gsVisitorNavierStokes<T> visitor(*m_pde_ptr,velocity,pressure);
132  Base::template push<gsVisitorNavierStokes<T> >(visitor);
133 
134  m_system.matrix().makeCompressed();
135 }
136 
137 //--------------------- SOLUTION CONSTRUCTION ----------------------------------//
138 
139 template <class T>
141  const std::vector<gsMatrix<T> > & fixedDoFs,
142  gsMultiPatch<T>& velocity) const
143 {
144  gsVector<index_t> unknowns(m_dim);
145  for (short_t d = 0; d < m_dim; ++d)
146  unknowns.at(d) = d;
147  Base::constructSolution(solVector,fixedDoFs,velocity,unknowns);
148 }
149 
150 template <class T>
152  const std::vector<gsMatrix<T> > & fixedDoFs,
153  gsMultiPatch<T> & velocity, gsMultiPatch<T> & pressure) const
154 {
155  // construct displacement
156  constructSolution(solVector,fixedDoFs,velocity);
157  // construct pressure
158  constructPressure(solVector,fixedDoFs,pressure);
159 }
160 
161 template <class T>
163  const std::vector<gsMatrix<T> > & fixedDoFs,
164  gsMultiPatch<T>& pressure) const
165 {
166  gsVector<index_t> unknowns(1);
167  unknowns.at(0) = m_dim;
168  Base::constructSolution(solVector,m_ddof,pressure,unknowns);
169 }
170 
171 //--------------------- SPECIALS ----------------------------------//
172 
173 template <class T>
175  const std::vector<std::pair<index_t,boxSide> > & bdrySides, bool split) const
176 {
177  // all temporary data structures
178  gsMatrix<T> quNodes, pressureValues, physGradJac, sigma;
179  gsVector<T> quWeights, normal;
180  // NEED_MEASURE for integration
181  // NEED_GRAD_TRANSFORM for velocity gradients transformation from parametric to physical domain
183  // NEED_DERIV for velocity gradients
184  gsMapData<T> mdVel(NEED_DERIV);
185 
186  gsMatrix<T> force;
187  force.setZero(m_dim, 2);
188  const T viscosity = m_options.getReal("Viscosity");
189  const T density = m_options.getReal("Density");
190 
191  // loop over bdry sides
192  for (std::vector<std::pair<index_t,boxSide> >::const_iterator it = bdrySides.begin();
193  it != bdrySides.end(); ++it)
194  {
195  // basis of the patch
196  const gsBasis<T> & basis = m_bases[0][it->first];
197  // setting quadrature rule for the boundary side
198  gsGaussRule<T> bdQuRule(basis,1.0,1,it->second.direction());
199  // loop over elements of the side
200  typename gsBasis<T>::domainIter elem = basis.makeDomainIterator(it->second);
201  for (; elem->good(); elem->next())
202  {
203  // mapping quadrature rule to the element
204  bdQuRule.mapTo(elem->lowerCorner(),elem->upperCorner(),quNodes,quWeights);
205  // evaluate geoemtry mapping at the quad points
206  mdGeo.points = quNodes;
207  m_pde_ptr->patches().patch(it->first).computeMap(mdGeo);
208  // evaluate velocity at the quad points
209  mdVel.points = quNodes;
210  velocity.patch(it->first).computeMap(mdVel);
211  // evaluate pressure at the quad points
212  pressure.patch(it->first).eval_into(quNodes,pressureValues);
213 
214  // loop over quad points
215  for (index_t q = 0; q < quWeights.rows(); ++q)
216  {
217  // transform gradients from parametric to physical
218  physGradJac = mdVel.jacobian(q)*(mdGeo.jacobian(q).cramerInverse());
219  // normal length is the local measure
220  outerNormal(mdGeo,q,it->second,normal);
221  // pressure contribution to the stress tensor
222  sigma = pressureValues.at(q)*gsMatrix<T>::Identity(m_dim,m_dim);
223  force.col(0) += quWeights[q] * sigma * normal;
224  sigma = -1*density*viscosity*(physGradJac + physGradJac.transpose());
225  force.col(1) += quWeights[q] * sigma * normal;
226  }
227  }
228  }
229  if (split)
230  return force;
231  else
232  return force.col(0) + force.col(1);
233 }
234 
235 } // namespace ends
void addPiece(const gsFunction< T > &func)
Add a piece.
Definition: gsPiecewiseFunction.h:113
virtual void reserve()
a custom reserve function to allocate memory for the sparse matrix
Definition: gsNsAssembler.hpp:61
virtual gsMatrix< T > computeForce(const gsMultiPatch< T > &velocity, const gsMultiPatch< T > &pressure, const std::vector< std::pair< index_t, boxSide > > &bdrySides, bool split=false) const
Definition: gsNsAssembler.hpp:174
Gradient of the object.
Definition: gsForwardDeclarations.h:73
virtual void constructSolution(const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDoFs, gsMultiPatch< T > &velocity) const
Construct velocity from computed solution vector and fixed degrees of freedom.
Definition: gsNsAssembler.hpp:140
#define short_t
Definition: gsConfig.h:35
Visitor class for volumetric integration of the Stokes system.
the gsMapData is a cache of pre-computed function (map) values.
Definition: gsFuncData.h:324
gsNsAssembler(const gsMultiPatch< T > &patches, const gsMultiBasis< T > &basisVel, const gsMultiBasis< T > &basisPres, const gsBoundaryConditions< T > &bconditions, const gsFunction< T > &body_force)
Constructor.
Definition: gsNsAssembler.hpp:30
#define index_t
Definition: gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition: gsDebug.h:102
A function from a n-dimensional domain to an m-dimensional image.
Definition: gsFunction.h:59
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition: gsNsAssembler.hpp:50
virtual short_t targetDim() const
Dimension of the target space.
Definition: gsFunctionSet.h:560
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition: gsOptionList.cpp:201
T at(index_t i) const
Returns the i-th element of the vectorization of the matrix.
Definition: gsMatrix.h:211
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition: gsVector.h:35
Holds a set of patch-wise bases and their topology information.
Definition: gsMultiBasis.h:36
virtual void refresh()
Refresh routine to set dof-mappers.
Definition: gsNsAssembler.hpp:80
gsGeometry< T > & patch(size_t i) const
Return the i-th patch.
Definition: gsMultiPatch.h:226
T at(index_t i) const
Returns the i-th element of the vector.
Definition: gsVector.h:177
Provides matrix and rhs assebmly for stationary and transient incompressible Stokes and Navier-Stokes...
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition: gsMultiPatch.h:33
IMHO, a useless class, but it is necessary to use the gsAssembler class. Contains proper information ...
The density of the measure pull back.
Definition: gsForwardDeclarations.h:76
Visitor class for volumetric integration of thetangential Navier-Stokes system.
Class containing a set of boundary conditions.
Definition: gsBoundaryConditions.h:341
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition: gsOptionList.cpp:211
Provides functions to generate structured point data.
virtual void constructPressure(const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDoFs, gsMultiPatch< T > &pressure) const
@ brief Construct pressure from computed solution vector
Definition: gsNsAssembler.hpp:162
gsMatrix< T > points
input (parametric) points
Definition: gsFuncData.h:348
Gradient transformation matrix.
Definition: gsForwardDeclarations.h:77
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
A function depending on an index i, typically referring to a patch/sub-domain. On each patch a differ...
Definition: gsPiecewiseFunction.h:28
virtual void assemble()
Main assemble routine, to be implemented in derived classes.
Definition: gsAssembler.hpp:51
Class that represents the (tensor) Gauss-Legendre quadrature rule.
Definition: gsGaussRule.h:27
virtual domainIter makeDomainIterator() const
Create a domain iterator for the computational mesh of this basis, that points to the first element o...
Definition: gsBasis.hpp:493
stationary point iteration, 1st order, yields a new solution at each iteration
Definition: gsBaseUtils.h:47
A basis represents a family of scalar basis functions defined over a common parameter domain...
Definition: gsBasis.h:78