Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
implicit.h File Reference
#include "poisson.h"
Include dependency graph for implicit.h:

Go to the source code of this file.

Functions

void event_defaults (void)
 Event: defaults (i = 0)
 
void event_tracer_diffusion (void)
 Event: tracer_diffusion (i++)
 

Variables

double dt
 
scalar phi []
 The electric potential and the volume charge density are scalars while the permittivity and conductivity are face vectors.
 
scalar rhoe []
 
vector epsilon []
 
vector K []
 
mgstats mgphi
 

Function Documentation

◆ event_defaults()

void event_defaults ( void  )

Event: defaults (i = 0)

Boundary conditions for VOF-advected tracers usually depend on boundary conditions for the VOF field.

Event: defaults (i = 0)

Event: defaults (i = 0)

Initialisation

We set the default values.

The restriction/refine attributes of the charge density are those of a tracer otherwise the conservation is not guaranteed.

By default the permittivity is unity and other quantities are zero.

Definition at line 39 of file implicit.h.

References _i, epsilon, fm, K, refine_linear(), restriction_volume_average(), rhoe, vector::x, and x.

Here is the call graph for this function:

◆ event_tracer_diffusion()

void event_tracer_diffusion ( void  )

Event: tracer_diffusion (i++)

The r.h.s of the poisson equation is set to \(-\rho_e^n\),

We compute the coefficients of the Laplacian: $(K dt +\epsilon)$.

The poisson equation is solved to obtain \(\phi^{n+1}\).

Finally \(\rho_e^{n+1}\) is computed as \(\rho_e^{n+1} = -\nabla \cdot (\epsilon \nabla \phi^{n+1})\).

In the absence of conductivity, the electric potential only varies if the electrical permittivity varies,

Definition at line 59 of file implicit.h.

References _i, cm, dimension, dt, epsilon, f, scalar::i, K, mgphi, phi, poisson, rhoe, sq(), vector::x, and x.

Here is the call graph for this function:

Variable Documentation

◆ dt

double dt
extern

Ohmic conduction

This function approximates implicitly the EHD equation set given by the electric Poisson equation and the ohmic conduction term of the charge density conservation (the advection term is computed elsewhere using a tracer advection scheme),

\[ \nabla \cdot( \epsilon \nabla \phi^{n+1}) = -\rho_e^{n+1} \quad \mathrm{and} \quad (\rho_e^{n+1}-\rho_e^n) = \Delta t\nabla \cdot (K\nabla \phi^{n+1}) \]

where \(\rho_e\) is the charge density, and \(\phi\) is the electric potential. \(K\) and \(\epsilon\) are the conductivity and permittivity respectively.

Substituting the Poisson equation into the conservation equation gives the following time-implicit scheme,

\[ \nabla \cdot [(K \, \Delta t + \epsilon) \nabla \phi^{n+1}] = -\rho_e^{n} \]

We need the Poisson solver and the timestep dt.

Generic time loop

The run() function below implements a generic time loop which executes events until termination.

The timestep dt can be accessed as a global variable.

Definition at line 18 of file predictor-corrector.h.

Referenced by advance_saint_venant(), advect(), advection(), airsea_integrated_fluxes(), correct_qz(), correction(), dphidt(), dtnext(), event_acceleration(), event_advance(), event_advection_term(), event_end_timestep(), event_face_fields(), event_half_advection(), event_perfs(), event_pressure(), event_projection(), event_stability(), event_tracer_advection(), event_tracer_diffusion(), event_velocity(), event_viscous_term(), gotm_deprecated_output(), gotm_step(), h_relax(), h_residual(), hllc(), horizontal_diffusion(), implicit_horizontal_diffusion(), kinetic(), kurganov(), meanflow_buoyancy(), meanflow_coriolis(), meanflow_salinity(), meanflow_stratification(), meanflow_temperature(), meanflow_uequation(), meanflow_updategrid(), meanflow_vequation(), meanflow_wequation(), ohmic_flux(), prediction(), project(), relax_diffusion(), relax_nh(), relax_viscosity(), residual_diffusion(), residual_nh(), residual_viscosity(), riemann(), run(), runge_kutta(), timestep(), timestep(), tracer_fluxes(), turbulence_dissipationeq(), turbulence_do_epsb(), turbulence_do_kb(), turbulence_do_lengthscale(), turbulence_do_tke(), turbulence_do_turbulence(), turbulence_genericeq(), turbulence_kbeq(), turbulence_lengthscaleeq(), turbulence_q2over2eq(), turbulence_tkeeq(), update(), update_green_naghdi(), update_tracer(), util_adv_center(), util_diff_center(), util_diff_face(), util_findp_bisection(), util_lagrange(), vertical_diffusion(), vertical_viscosity(), viscosity(), and viscosity_explicit().

◆ epsilon

vector epsilon[]

Definition at line 35 of file implicit.h.

Referenced by event_defaults(), and event_tracer_diffusion().

◆ K

vector K[]

Definition at line 35 of file implicit.h.

Referenced by event_defaults(), event_tracer_diffusion(), and ohmic_flux().

◆ mgphi

mgstats mgphi

Definition at line 36 of file implicit.h.

Referenced by event_tracer_diffusion().

◆ phi

scalar phi[]

The electric potential and the volume charge density are scalars while the permittivity and conductivity are face vectors.

In mgphi we will store the statistics for the multigrid resolution of the electric poisson equation.

Definition at line 34 of file implicit.h.

Referenced by correct_qz(), dphidt(), event_acceleration(), event_defaults(), event_init(), event_metric(), event_pressure(), event_tracer_diffusion(), fraction(), if(), if(), levelset_to_vof(), matrix(), matrix1(), refine_cm_lonlat(), refine_face_y_lonlat(), relax_nh(), relax_nh1(), residual_nh(), residual_nh2(), residual_nh3(), and solid().

◆ rhoe

scalar rhoe[]

Definition at line 34 of file implicit.h.

Referenced by event_defaults(), and event_tracer_diffusion().