Electrowetting on dielectric
|
Classes | |
struct | ScratchData |
Public Member Functions | |
Velocity (const Triangulation< dim > &tria, Material_Parameters ¶ms, 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_function & | density |
density | |
const material_function & | viscosity |
viscosity | |
const material_function & | slip_coeff |
slip coefficient | |
InterfaceEnergy & | gamma_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 |
The vector valued convection diffusion problem that defines the evolution of the velocity
where
This system is supplemented with the boundary conditions
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:
Iterator
[0]: Phase field from CahnHilliard Iterator
[1]: Chemical Potential from CahnHilliard Iterator
[2]: Charge Iterator
[3]: Voltage Iterator
[4]: Pressure extrapolation from Pressure Iterator
[5]: Value of the phase field variable at previous time step from CahnHilliard 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
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). |
Advance in time. CHNSE calls this function once the inner iteration loop has converged.
Compute the inner iteration error. CahnHilliard calls this function to verify wether the inner iteration has converged.
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 >.
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 >.
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
data | : The views of the other problems this one needs information from to be able to assemble.
|
time | : The current time. |
m_threaded | : flag that indicates if the assembly should be done using multithreading. |
Implements Problem< dim >.
void Velocity< dim >::AssembleCell | ( | const Iterator | Its, |
ScratchData & | scratch, | ||
typename Problem< dim >::PerTaskData & | data | ||
) | const [protected] |
Assemble the local problems.
Its | : The cell iterators for the current problem and the data. |
scratch | : The scratch data. |
data | : The per task data. |
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 >.
virtual void Velocity< dim >::ReinitPrec | ( | ) | [protected, virtual] |
Reinit the preconditioner
Implements Problem< dim >.
virtual void Velocity< dim >::DoSolve | ( | SolverControl & | control | ) | [protected, virtual] |
Solve the system
control | : Provides the maximal number of iterations and exit criterion. |
Implements Problem< dim >.
virtual void Velocity< dim >::SetInitialData | ( | ) | [protected, virtual] |
Set initial data. In this case .
Implements Problem< dim >.
virtual void Velocity< dim >::PreRefinementPreffix | ( | ) | [protected, virtual] |
Used to handle changes in the triangulation: Copy to the solution transfer
Implements Problem< 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 >.
const material_function& Velocity< dim >::density [protected] |
const material_function & Velocity< dim >::viscosity [protected] |
const material_function & Velocity< dim >::slip_coeff [protected] |
InterfaceEnergy& Velocity< dim >::gamma_fs [protected] |
const unsigned Velocity< dim >::Krylov_size [protected] |
TrilinosWrappers::PreconditionAMG Velocity< dim >::prec [protected] |
TrilinosWrappers::PreconditionAMG::AdditionalData Velocity< dim >::prec_data [protected] |
TrilinosWrappers::Vector Velocity< dim >::oldsol [protected] |
TrilinosWrappers::Vector Velocity< dim >::olditsol [protected] |
AsFunction<dim> Velocity< dim >::get_velocity |