Parallel Finite Element - General geometry Ewald-like Method
Part of Continuum-particle Simulation Suite under MICCOM
assemble_stokes.h
Go to the documentation of this file.
1 // Parallel Finite Element-General Geometry Ewald-like Method.
2 // Copyright (C) 2015-2016 Xujun Zhao, Jiyuan Li, Xikai Jiang
3 
4 // This code is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 
10 // This code is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 // Lesser General Public License for more details.
14 
15 
16 // You should have received a copy of the GNU General Public
17 // License along with this code; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 
20 
21 #pragma once
22 
23 // Local includes
24 #include "assemble_system.h"
25 
36 {
37 public:
42  AssembleStokes(EquationSystems& es);
43 
44 
49 
50 
58  void assemble_global_K(const std::string& system_name,
59  const std::string& option) override;
60 
61 
68  void assemble_global_F(const std::string& system_name,
69  const std::string& option) override;
70 
71 
78  void assemble_element_KIJ(const std::vector<Real>& JxW,
79  const std::vector<std::vector<RealGradient> >& dphi,
80  const unsigned int n_u_dofs,
81  const unsigned int I,
82  const unsigned int J,
83  DenseMatrix<Number>& Kij) override;
84 
85 
91  void assemble_element_QI(const std::vector<Real>& JxW,
92  const std::vector<std::vector<RealGradient> >& dphi,
93  const std::vector<std::vector<Real> >& psi,
94  const unsigned int n_v_dofs,
95  const unsigned int n_p_dofs,
96  const unsigned int I,
97  DenseMatrix<Number>& Qi);
98 
99 
104  void assemble_element_MIJ(const std::vector<Real>& JxW,
105  const std::vector<std::vector<Real> >& phi,
106  const unsigned int n_v_dofs,
107  DenseMatrix<Number>& Mij);
108 
109 
114  void compute_element_rhs(const Elem* elem,
115  const unsigned int n_u_dofs,
116  FEBase& fe_v,
117  const std::vector<std::size_t> n_list,
118  const bool& pf_flag,
119  const std::string& option,
120  const Real& alpha,
121  DenseVector<Number>& Fe) override;
122 
126  void select_boundary_side(const Elem* elem) override;
127 
131  void apply_bc_by_penalty(const Elem* elem,
132  const std::string& matrix_or_vector,
133  DenseMatrix<Number>& Ke,
134  DenseVector<Number>& Fe,
135  const std::string& option) override;
136 
137 
141  Real boundary_pressure_jump(const std::string& which_side) const;
142 
143 
147  Real boundary_traction(const std::string& which_side) const;
148 };
void assemble_element_QI(const std::vector< Real > &JxW, const std::vector< std::vector< RealGradient > > &dphi, const std::vector< std::vector< Real > > &psi, const unsigned int n_v_dofs, const unsigned int n_p_dofs, const unsigned int I, DenseMatrix< Number > &Qi)
Assemble the element matrices Q_I, i.e., kup, kvp, kwp, for pressure.
Definition: assemble_stokes.C:636
This base class provides the template for assembling the matrix and right-hand-side vector when solvi...
Definition: assemble_system.h:48
void select_boundary_side(const Elem *elem) override
select sides on the boundary for all elements
Definition: assemble_stokes.C:685
void assemble_element_MIJ(const std::vector< Real > &JxW, const std::vector< std::vector< Real > > &phi, const unsigned int n_v_dofs, DenseMatrix< Number > &Mij)
Assemble the element mass matrix M for preconditioning matrix.
Definition: assemble_stokes.C:662
Real boundary_traction(const std::string &which_side) const
Define the pressure jump at the inlet and outlet of the channel.
Definition: assemble_stokes.C:844
AssembleStokes(EquationSystems &es)
Constructor.
Definition: assemble_stokes.C:52
void compute_element_rhs(const Elem *elem, const unsigned int n_u_dofs, FEBase &fe_v, const std::vector< std::size_t > n_list, const bool &pf_flag, const std::string &option, const Real &alpha, DenseVector< Number > &Fe) override
Assemble function for the right-hand-side in Stokes equation.
Definition: assemble_stokes.C:493
void apply_bc_by_penalty(const Elem *elem, const std::string &matrix_or_vector, DenseMatrix< Number > &Ke, DenseVector< Number > &Fe, const std::string &option) override
Apply BCs by penalty method.
Definition: assemble_stokes.C:722
~AssembleStokes()
Destructor.
Definition: assemble_stokes.C:61
This class provides the basic components for assembling the matrix and vector when solving Stokes equ...
Definition: assemble_stokes.h:35
Real boundary_pressure_jump(const std::string &which_side) const
Define the pressure jump at the inlet and outlet of the channel.
Definition: assemble_stokes.C:818
void assemble_element_KIJ(const std::vector< Real > &JxW, const std::vector< std::vector< RealGradient > > &dphi, const unsigned int n_u_dofs, const unsigned int I, const unsigned int J, DenseMatrix< Number > &Kij) override
Assemble the element matrix K_IJ.
Definition: assemble_stokes.C:607
void assemble_global_K(const std::string &system_name, const std::string &option) override
Assemble the Global Matrix K.
Definition: assemble_stokes.C:69
void assemble_global_F(const std::string &system_name, const std::string &option) override
Assemble the Global force vector F.
Definition: assemble_stokes.C:308