Electrowetting on dielectric
/export/home/gkanscha/ANS/submissions/Salgado/include/NSE.h
Go to the documentation of this file.
00001 #ifndef _NSE_H_
00002 #define _NSE_H_
00003 
00004 #include "Problem.h"
00005 
00024 template<int dim> class Pressure: public Problem<dim>{
00025   public:
00038     Pressure( const Triangulation<dim> &tria, Material_Parameters &params, const unsigned deg,
00039               const double ee, const double thres, const unsigned sweeps );
00043     virtual ~Pressure();
00044   protected:
00045     const double rho_min; 
00046     TrilinosWrappers::PreconditionAMG prec; 
00047     TrilinosWrappers::PreconditionAMG::AdditionalData prec_data; 
00048     TrilinosWrappers::Vector p, 
00049                              pres_extr; 
00050 
00055     virtual void SetupDoFsSuffix();
00059     virtual void InitLADataSuffix();
00068     typedef IteratorGroup<2, typename hp::DoFHandler<dim>::active_cell_iterator> Iterator;
00077     struct ScratchData{
00078       hp::QCollection<dim> quad; 
00079       UpdateFlags flags; 
00080       hp::FEValues<dim> fe_val; 
00081       AsFunction<dim> velocity; 
00082       std::vector<double> loc_div_vel;
00083       ScratchData( const hp::FECollection<dim> &fe, std::vector<AsFunction<dim> *> &data,
00084                    const hp::QCollection<dim> &q, const UpdateFlags u_flags ) :
00085           quad(q), flags( u_flags ), fe_val( fe, quad, flags ), velocity( *( data[0] ) ),
00086           loc_div_vel( quad.max_n_quadrature_points() ) {}
00087       ScratchData( const ScratchData &scratch ) :
00088           quad( scratch.quad ), flags( scratch.flags ),
00089           fe_val( scratch.fe_val.get_fe_collection(), quad, flags ), velocity( scratch.velocity ),
00090           loc_div_vel( quad.max_n_quadrature_points() ) {}
00091     };
00102     virtual void AssembleSystem( std::vector< AsFunction<dim> *> &data, const double time,
00103                                  const bool m_threaded = true );
00110     void AssembleCell( const Iterator Its, ScratchData &scratch,
00111                        typename Problem<dim>::PerTaskData &data );
00115     virtual void ReinitPrec();
00119     virtual void DoSolve( SolverControl &control );
00125     virtual void SolveSystemSuffix();
00129     virtual void SetInitialData();
00133     virtual void PreRefinementPreffix();
00137     virtual void PostRefinementSuffix( const std::vector<TrilinosWrappers::Vector> &sol_tmp );
00138   public:
00139     AsFunction<dim> get_pressure, 
00140                     get_extrapolated_pressure; 
00141 };
00142 
00164 template<int dim> class Velocity: public Problem<dim>{
00165   public:
00177     Velocity( const Triangulation<dim> &tria, Material_Parameters &params, const unsigned deg, const unsigned u_prec,
00178               const double ee, const double thres, const unsigned v_sweeps, const unsigned K_size );
00182     virtual ~Velocity();
00186     void advance();
00191     double it_error();
00192   protected:
00193     const material_function &density, 
00194                             &viscosity, 
00195                             &slip_coeff; 
00196     InterfaceEnergy &gamma_fs; 
00197     const double alpha, 
00198                  lambda, 
00199                  gamma, 
00200                  delta; 
00201     const unsigned Krylov_size; 
00202     TrilinosWrappers::PreconditionAMG prec; // Preconditioner (AMG)
00203     TrilinosWrappers::PreconditionAMG::AdditionalData prec_data; 
00204     TrilinosWrappers::Vector oldsol, 
00205                              olditsol; 
00206 
00211     virtual void InitLADataSuffix();
00222     virtual void SetupDoFsSuffix();
00234     typedef IteratorGroup<7,typename hp::DoFHandler<dim>::active_cell_iterator> Iterator;
00243     struct ScratchData{
00244       hp::QCollection<dim> quad; 
00245       UpdateFlags flags; 
00246       hp::FEValues<dim> fe_val; 
00247       hp::QCollection<dim-1> fquad; 
00248       UpdateFlags f_flags; 
00249       hp::FEFaceValues<dim> fe_face_val; 
00250       AsFunction<dim> phase, 
00251                       mu, 
00252                       q, 
00253                       V, 
00254                       psharp, 
00255                       oldphase; 
00256       std::vector<double> loc_phase, loc_old_phase, face_phase, face_old_phase, loc_q,
00257                           loc_psharp, loc_divu, loc_mu;
00258       std::vector< Tensor<1,dim> > loc_grad_q, loc_grad_V, loc_grad_phase, face_grad_phase, loc_u, face_base,
00259                                    loc_face_u;
00260       Point<dim> normal;
00261       ScratchData( const hp::FECollection<dim> &fe, std::vector<AsFunction<dim> *> &data,
00262                    const hp::QCollection<dim> &q, const UpdateFlags u_flags,
00263                    const hp::QCollection<dim-1> &fq, const UpdateFlags u_f_flags ) :
00264           quad(q), flags( u_flags ), fe_val( fe, quad, flags ),
00265           fquad( fq ), f_flags( u_f_flags ), fe_face_val( fe, fquad, f_flags ),
00266           phase( *( data[0] ) ), mu( *( data[1] ) ), q( *( data[2] ) ), V( *( data[3] ) ),
00267           psharp( *( data[4] ) ), oldphase( *( data[5] ) ),
00268           loc_phase( quad.max_n_quadrature_points() ), loc_old_phase( quad.max_n_quadrature_points() ),
00269           face_phase( fquad.max_n_quadrature_points() ), face_old_phase( fquad.max_n_quadrature_points() ),
00270           loc_q( quad.max_n_quadrature_points() ), loc_psharp( quad.max_n_quadrature_points() ),
00271           loc_divu( quad.max_n_quadrature_points() ), loc_mu( quad.max_n_quadrature_points() ),
00272           loc_grad_q( quad.max_n_quadrature_points() ), loc_grad_V( quad.max_n_quadrature_points() ),
00273           loc_grad_phase( quad.max_n_quadrature_points() ), face_grad_phase( fquad.max_n_quadrature_points() ),
00274           loc_u( quad.max_n_quadrature_points() ), face_base( fe.max_dofs_per_cell() ),
00275           loc_face_u( fquad.max_n_quadrature_points() ) {}
00276       ScratchData( const ScratchData &scratch ) :
00277           quad( scratch.quad ), flags( scratch.flags ), fe_val( scratch.fe_val.get_fe_collection(), quad, flags ),
00278           fquad( scratch.fquad ), f_flags( scratch.f_flags ),
00279           fe_face_val( scratch.fe_face_val.get_fe_collection(), fquad, f_flags ),
00280           phase( scratch.phase ), mu( scratch.mu ), q( scratch.q ), V( scratch.V ),
00281           psharp( scratch.psharp ), oldphase( scratch.oldphase ),
00282           loc_phase( quad.max_n_quadrature_points() ), loc_old_phase( quad.max_n_quadrature_points() ),
00283           face_phase( fquad.max_n_quadrature_points() ), face_old_phase( fquad.max_n_quadrature_points() ),
00284           loc_q( quad.max_n_quadrature_points() ), loc_psharp( quad.max_n_quadrature_points() ),
00285           loc_divu( quad.max_n_quadrature_points() ), loc_mu( quad.max_n_quadrature_points() ),
00286           loc_grad_q( quad.max_n_quadrature_points() ), loc_grad_V( quad.max_n_quadrature_points() ),
00287           loc_grad_phase( quad.max_n_quadrature_points() ), face_grad_phase( fquad.max_n_quadrature_points() ),
00288           loc_u( quad.max_n_quadrature_points() ), face_base( fe_val.get_fe_collection().max_dofs_per_cell() ),
00289           loc_face_u( fquad.max_n_quadrature_points() ) {}
00290     };
00305     virtual void AssembleSystem( std::vector< AsFunction<dim> *> &data, const double time, const bool m_threaded = true );
00312     void AssembleCell( const Iterator Its, ScratchData &scratch, typename Problem<dim>::PerTaskData &data ) const;
00317     virtual void SolveSystemPreffix();
00321     virtual void ReinitPrec();
00326     virtual void DoSolve( SolverControl &control );
00330     virtual void SetInitialData();
00334     virtual void PreRefinementPreffix();
00338     virtual void PostRefinementSuffix( const std::vector<TrilinosWrappers::Vector> &sol_tmp );
00339   public:
00340     AsFunction<dim> get_velocity; 
00341 };
00342 
00343 #endif // _NSE_H_
00344