Electrowetting on dielectric
Classes | Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes
Velocity< dim > Class Template Reference
Inheritance diagram for Velocity< dim >:
Inheritance graph
[legend]
Collaboration diagram for Velocity< dim >:
Collaboration graph
[legend]

List of all members.

Classes

struct  ScratchData

Public Member Functions

 Velocity (const Triangulation< dim > &tria, Material_Parameters &params, const unsigned deg, const unsigned u_prec, const double ee, const double thres, const unsigned v_sweeps, const unsigned K_size)
virtual ~Velocity ()
void advance ()
double it_error ()

Public Attributes

AsFunction< dim > get_velocity
 velocity

Protected Types

typedef IteratorGroup
< 7, typename hp::DoFHandler
< dim >::active_cell_iterator > 
Iterator

Protected Member Functions

virtual void InitLADataSuffix ()
virtual void SetupDoFsSuffix ()
virtual void AssembleSystem (std::vector< AsFunction< dim > * > &data, const double time, const bool m_threaded=true)
void AssembleCell (const Iterator Its, ScratchData &scratch, typename Problem< dim >::PerTaskData &data) const
virtual void SolveSystemPreffix ()
virtual void ReinitPrec ()
virtual void DoSolve (SolverControl &control)
virtual void SetInitialData ()
virtual void PreRefinementPreffix ()
virtual void PostRefinementSuffix (const std::vector< TrilinosWrappers::Vector > &sol_tmp)

Protected Attributes

const material_functiondensity
 density
const material_functionviscosity
 viscosity
const material_functionslip_coeff
 slip coefficient
InterfaceEnergygamma_fs
 interface free energy
const double alpha
 surface transport parameter
const double lambda
 charge coupling parameter
const double gamma
 surface tension
const double delta
 interface thickness
const unsigned Krylov_size
 dimension of the Krylov subspace
TrilinosWrappers::PreconditionAMG prec
TrilinosWrappers::PreconditionAMG::AdditionalData prec_data
 Preconditioner data.
TrilinosWrappers::Vector oldsol
 the old solution
TrilinosWrappers::Vector olditsol
 the old iterative solution

Detailed Description

template<int dim>
class Velocity< dim >

The vector valued convection diffusion problem that defines the evolution of the velocity

\[ \frac{ D(\rho(\phi){\mathbf{u}}) }{Dt} - \nabla\cdot\left( \eta(\phi) \mathbf{S}({\mathbf{u}}) \right) + \nabla p = \mu \nabla \phi - q \nabla \left( V + \lambda q \right) + \frac12 \rho'(\phi)\phi_t {\mathbf{u}}, \quad \text{on }\Omega, \]

where

\[ \frac{ D(\rho(\phi){\mathbf{u}}) }{Dt} = \sigma(\sigma \mathbf{u})_t + \frac12 \nabla\cdot(\rho \mathbf{u})\mathbf{u}, \quad \mathbf{S}(\mathbf{u}) = \frac12 \left( \nabla \mathbf{u} + {\nabla \mathbf{u}}^T \right). \]

This system is supplemented with the boundary conditions

\[ {\mathbf{u}}\cdot{\mathbf{n}} = 0, \qquad \beta(\phi) {\mathbf{u}}_{\boldsymbol{\tau}} + \eta(\phi) \mathbf{S}({\mathbf{u}})_{{\mathbf{n}}{\boldsymbol{\tau}}} = \gamma\left( \Theta_{fs}'(\phi) + \delta \partial_{\mathbf{n}} \phi \right) \partial_{\boldsymbol{\tau}} \phi, \quad \text{on }\Gamma, \]

Definition at line 164 of file NSE.h.


Member Typedef Documentation

template<int dim>
typedef IteratorGroup<7,typename hp::DoFHandler<dim>::active_cell_iterator> Velocity< dim >::Iterator [protected]

We group the iterators of this problem and the problems this one takes information from as follows:

Definition at line 234 of file NSE.h.


Constructor & Destructor Documentation

template<int dim>
Velocity< dim >::Velocity ( const Triangulation< dim > &  tria,
Material_Parameters params,
const unsigned  deg,
const unsigned  u_prec,
const double  ee,
const double  thres,
const unsigned  v_sweeps,
const unsigned  K_size 
)

Copy the provided parameters to the members

Parameters:
tria: A reference to the mesh where the problem will live.
params: The collection of material parameters.
deg: The polynomial degree for the finite element spaces where the charge is defined.
u_prec: How often we update the preconditioner.
ee: The accuracy to which we iteratively solve the problem.
thres: The aggregation threshold for the AMG preconditioner.
v_sweeps: The number of smoothing steps the AMG preconditioner must perform.
K_size: The dimension of the Krylov subspace (used in GMRES).
template<int dim>
virtual Velocity< dim >::~Velocity ( ) [virtual]

Clear the data


Member Function Documentation

template<int dim>
void Velocity< dim >::advance ( )

Advance in time. CHNSE calls this function once the inner iteration loop has converged.

template<int dim>
double Velocity< dim >::it_error ( )

Compute the inner iteration error. CahnHilliard calls this function to verify wether the inner iteration has converged.

template<int dim>
virtual void Velocity< dim >::InitLADataSuffix ( ) [protected, virtual]

What to do after the Linear Algebra data has been initialized.

We set oldsol and olditsol to the correct size.

Reimplemented from Problem< dim >.

template<int dim>
virtual void Velocity< dim >::SetupDoFsSuffix ( ) [protected, virtual]

What to do after the DoFs are set.

We add the impermeability conditions as constraints to the ConstraintMatrix object. Notice that since this is applied not only at physical boundaries, but also at the interface between fluid and plates, special care must be taken. In general this seems like a difficult process, since the normals might not be well defined. What we do is assume that the interface is horizontal, so that by using the FESystem::system_to_component_index we know which local DOF corresponds to the vertical velocity.

Reimplemented from Problem< dim >.

template<int dim>
virtual void Velocity< dim >::AssembleSystem ( std::vector< AsFunction< dim > * > &  data,
const double  time,
const bool  m_threaded = true 
) [protected, virtual]

Assemble the system matrix and right hand side for the current status of the problem

Parameters:
data: The views of the other problems this one needs information from to be able to assemble.
  1. Phase field from CahnHilliard
  2. Chemical potential from CahnHilliard
  3. Charge
  4. Voltage
  5. Extrapolated pressure from Pressure
  6. Phase field at previous time step from CahnHilliard
time: The current time.
m_threaded: flag that indicates if the assembly should be done using multithreading.

Implements Problem< dim >.

template<int dim>
void Velocity< dim >::AssembleCell ( const Iterator  Its,
ScratchData scratch,
typename Problem< dim >::PerTaskData data 
) const [protected]

Assemble the local problems.

Parameters:
Its: The cell iterators for the current problem and the data.
scratch: The scratch data.
data: The per task data.
template<int dim>
virtual void Velocity< dim >::SolveSystemPreffix ( ) [protected, virtual]

What to do before solving the system. Since the solving is always part of an inner iteration loop, we save the previous iterate so we can check convergence.

Reimplemented from Problem< dim >.

template<int dim>
virtual void Velocity< dim >::ReinitPrec ( ) [protected, virtual]

Reinit the preconditioner

Implements Problem< dim >.

template<int dim>
virtual void Velocity< dim >::DoSolve ( SolverControl control) [protected, virtual]

Solve the system

Parameters:
control: Provides the maximal number of iterations and exit criterion.

Implements Problem< dim >.

template<int dim>
virtual void Velocity< dim >::SetInitialData ( ) [protected, virtual]

Set initial data. In this case $\mathbf{u}_0 = 0$.

Implements Problem< dim >.

template<int dim>
virtual void Velocity< dim >::PreRefinementPreffix ( ) [protected, virtual]

Used to handle changes in the triangulation: Copy to the solution transfer

Implements Problem< dim >.

template<int dim>
virtual void Velocity< dim >::PostRefinementSuffix ( const std::vector< TrilinosWrappers::Vector > &  sol_tmp) [protected, virtual]

Used to handle changes in the triangulation: Copy from the solution transfer

Implements Problem< dim >.


Member Data Documentation

template<int dim>
const material_function& Velocity< dim >::density [protected]

density

Definition at line 193 of file NSE.h.

template<int dim>
const material_function & Velocity< dim >::viscosity [protected]

viscosity

Definition at line 193 of file NSE.h.

template<int dim>
const material_function & Velocity< dim >::slip_coeff [protected]

slip coefficient

Definition at line 193 of file NSE.h.

template<int dim>
InterfaceEnergy& Velocity< dim >::gamma_fs [protected]

interface free energy

Definition at line 196 of file NSE.h.

template<int dim>
const double Velocity< dim >::alpha [protected]

surface transport parameter

Definition at line 197 of file NSE.h.

template<int dim>
const double Velocity< dim >::lambda [protected]

charge coupling parameter

Definition at line 197 of file NSE.h.

template<int dim>
const double Velocity< dim >::gamma [protected]

surface tension

Definition at line 197 of file NSE.h.

template<int dim>
const double Velocity< dim >::delta [protected]

interface thickness

Definition at line 197 of file NSE.h.

template<int dim>
const unsigned Velocity< dim >::Krylov_size [protected]

dimension of the Krylov subspace

Definition at line 201 of file NSE.h.

template<int dim>
TrilinosWrappers::PreconditionAMG Velocity< dim >::prec [protected]

Definition at line 202 of file NSE.h.

Preconditioner data.

Definition at line 203 of file NSE.h.

template<int dim>
TrilinosWrappers::Vector Velocity< dim >::oldsol [protected]

the old solution

Definition at line 204 of file NSE.h.

template<int dim>
TrilinosWrappers::Vector Velocity< dim >::olditsol [protected]

the old iterative solution

Definition at line 204 of file NSE.h.

template<int dim>
AsFunction<dim> Velocity< dim >::get_velocity

velocity

Definition at line 340 of file NSE.h.


The documentation for this class was generated from the following file: