mus_derQuanNernstPlanck_module Module

This module provides the MUSUBI specific functions for calculating macroscopic quantities from the state variables. The depending common interface between MUSUBI and ATELES is defined in the tem_derived_module. The functionality for accessing a variable from the state and evaluating a lua function are also provided in the tem_derived module. A Novel lattice boltzmann model for nernstPlanck equation author> Zhenhua Chai , Baochang Shi A Coupled Lattice Boltzmann method to solve Nernst -Planck Model for simulating Electro-Osmotic Flows public :: deriveSrc_electricField


Uses

Used by

  • module~~mus_derquannernstplanck_module~~UsedByGraph module~mus_derquannernstplanck_module mus_derQuanNernstPlanck_module module~mus_variable_module mus_variable_module module~mus_variable_module->module~mus_derquannernstplanck_module module~mus_scheme_module mus_scheme_module module~mus_scheme_module->module~mus_variable_module module~mus_hvs_config_module mus_hvs_config_module module~mus_hvs_config_module->module~mus_scheme_module module~mus_config_module mus_config_module module~mus_hvs_config_module->module~mus_config_module module~mus_config_module->module~mus_scheme_module module~mus_tools_module mus_tools_module module~mus_config_module->module~mus_tools_module program~mus_harvesting mus_harvesting program~mus_harvesting->module~mus_scheme_module program~mus_harvesting->module~mus_hvs_config_module module~mus_dynloadbal_module mus_dynLoadBal_module module~mus_dynloadbal_module->module~mus_scheme_module module~mus_dynloadbal_module->module~mus_tools_module module~mus_program_module mus_program_module module~mus_program_module->module~mus_scheme_module module~mus_program_module->module~mus_dynloadbal_module module~mus_program_module->module~mus_tools_module module~mus_tools_module->module~mus_scheme_module program~musubi musubi program~musubi->module~mus_config_module program~musubi->module~mus_program_module module~mus_hvs_aux_module mus_hvs_aux_module module~mus_hvs_aux_module->module~mus_tools_module module~mus_tracking_module mus_tracking_module module~mus_tracking_module->module~mus_tools_module module~mus_aux_module mus_aux_module module~mus_aux_module->module~mus_tools_module module~mus_interpolate_verify_module mus_interpolate_verify_module module~mus_interpolate_verify_module->module~mus_config_module

Contents


Subroutines

public subroutine deriveAuxNP_fromState(derVarPos, state, neigh, iField, nElems, nSize, iLevel, stencil, varSys, auxField)

subroutine to add derive variables for weakly compressible PB (schemekind = 'nernstPlanck') to the varsys. A Coupled Lattice Boltzmann Method to Solve Nernst-Planck Model for Simulating Electro-Osmotic flows author> Xuguang yang This routine computes auxField 'mole_density' from state array

Read more…

Arguments

TypeIntentOptionalAttributesName
class(mus_derVarPos_type), intent(in) :: derVarPos

Position of derive variable in variable system

real(kind=rk), intent(in) :: state(:)

Array of state n * layout%stencil(1)%QQ * nFields

integer, intent(in) :: neigh(:)

connectivity vector

integer, intent(in) :: iField

Current field

integer, intent(in) :: nElems

number of elements

integer, intent(in) :: nSize

number of elements in state array

integer, intent(in) :: iLevel

current level

type(tem_stencilHeader_type), intent(in) :: stencil

stencil header contains discrete velocity vectors

type(tem_varSys_type), intent(in) :: varSys

variable system which is required to access fieldProp information via variable method data c_ptr

real(kind=rk), intent(inout) :: auxField(:)

Output of this routine Size: nElems*nAuxScalars

public subroutine deriveEquilNP_fromAux(derVarPos, auxField, iField, nElems, varSys, layout, fEq)

This routine computes equilbrium from auxField

Read more…

Arguments

TypeIntentOptionalAttributesName
class(mus_derVarPos_type), intent(in) :: derVarPos

Position of derive variable in variable system

real(kind=rk), intent(in) :: auxField(:)

Array of auxField. Single species: dens_1, vel_1, dens_2, vel_2, .. dens_n, vel_n multi-species: dens_1_sp1, vel_1_spc1, dens_1_sp2, vel_1_spc2, dens_2_sp1, vel_2_spc2, dens_2_sp2, vel_2_spc2 ... dens_n_sp1, vel_n_sp1, dens_n_sp2, vel_n_spc2 Access: (iElem-1)*nAuxScalars + auxField_varPos

integer, intent(in) :: iField

Current field

integer, intent(in) :: nElems

number of elements

type(tem_varSys_type), intent(in) :: varSys

variable system which is required to access fieldProp information via variable method data c_ptr

type(mus_scheme_layout_type), intent(in) :: layout

scheme layout contains stencil definition and lattice weights

real(kind=rk), intent(out) :: fEq(:)

Output of this routine Dimension: n*QQ of res

public recursive subroutine mus_deriveMoleDensity(fun, varsys, stencil, iLevel, posInState, pdf, res, nVals)

Calculate the potential of a given set of pdfs of elements

Read more…

Arguments

TypeIntentOptionalAttributesName
class(tem_varSys_op_type), intent(in) :: fun

description of the method to obtain the variables, here some preset values might be stored, like the space time function to use or the required variables.

type(tem_varSys_type), intent(in) :: varsys

the variable system to obtain the variable from.

type(tem_stencilHeader_type), intent(in) :: stencil

fluid stencil defintion

integer, intent(in) :: iLevel

current Level

integer, intent(in) :: posInState(:)

Position of element in levelwise state array

real(kind=rk), intent(in) :: pdf(:)

pdf array

real(kind=rk), intent(out) :: res(:)

results

integer, intent(in) :: nVals

nVals to get

public recursive subroutine deriveMoleDensity_forElement(fun, varSys, elempos, time, tree, nElems, nDofs, res)

Calculate the potential of a given set of elements (sum up all links). This routine is used to compute potential for all scheme kinds

Read more…

Arguments

TypeIntentOptionalAttributesName
class(tem_varSys_op_type), intent(in) :: fun

Description of the method to obtain the variables, here some preset values might be stored, like the space time function to use or the required variables.

type(tem_varSys_type), intent(in) :: varSys

The variable system to obtain the variable from.

integer, intent(in) :: elempos(:)

Position of the TreeID of the element to get the variable for in the global treeID list.

type(tem_time_type), intent(in) :: time

Point in time at which to evaluate the variable.

type(treelmesh_type), intent(in) :: tree

global treelm mesh info

integer, intent(in) :: nElems

Number of values to obtain for this variable (vectorized access).

integer, intent(in) :: nDofs

Number of degrees of freedom within an element.

real(kind=rk), intent(out) :: res(:)

Resulting values for the requested variable.

Linearized array dimension: (n requested entries) x (nComponents of this variable) x (nDegrees of freedom) Access: (iElem-1)fun%nComponentsnDofs + (iDof-1)*fun%nComponents + iComp

public recursive subroutine deriveMoleDensity_fromIndex(fun, varSys, time, iLevel, idx, idxLen, nVals, res)

Calculate the potential of a given set of elements (sum up all links). This routine is used to compute potential for all scheme kinds

Read more…

Arguments

TypeIntentOptionalAttributesName
class(tem_varSys_op_type), intent(in) :: fun

Description of the method to obtain the variables, here some preset values might be stored, like the space time function to use or the required variables.

type(tem_varSys_type), intent(in) :: varSys

The variable system to obtain the variable from.

type(tem_time_type), intent(in) :: time

Point in time at which to evaluate the variable.

integer, intent(in) :: iLevel

Level on which values are requested

integer, intent(in) :: idx(:)

Index of points in the growing array and variable val array to return. Size: nVals

integer, intent(in), optional :: idxLen(:)

With idx as start index in contiguous memory, idxLength defines length of each contiguous memory Size: nVals

integer, intent(in) :: nVals

Number of values to obtain for this variable (vectorized access).

real(kind=rk), intent(out) :: res(:)

Resulting values for the requested variable.

Dimension: n requested entries x nComponents of this variable Access: (iElem-1)*fun%nComponents + iComp

public subroutine applySrc_electricFieldNP(fun, inState, outState, neigh, auxField, nPdfSize, iLevel, varSys, time, phyConvFac, derVarPos)

Update state with source variable "electric field"

Read more…

Arguments

TypeIntentOptionalAttributesName
class(mus_source_op_type), intent(in) :: fun

Description of method to apply source terms

real(kind=rk), intent(in) :: inState(:)

input pdf vector

real(kind=rk), intent(inout) :: outState(:)

output pdf vector

integer, intent(in) :: neigh(:)

connectivity Array corresponding to state vector

real(kind=rk), intent(in) :: auxField(:)

auxField array

integer, intent(in) :: nPdfSize

number of elements in state Array

integer, intent(in) :: iLevel

current level

type(tem_varSys_type), intent(in) :: varSys

variable system

type(tem_time_type), intent(in) :: time

Point in time at which to evaluate the variable.

type(mus_convertFac_type), intent(in) :: phyConvFac

Physics conversion factor for current level

type(mus_derVarPos_type), intent(in) :: derVarPos(:)

position of derived quantities in varsys