|
The Open FUSION Toolkit 1.0.0-beta6
Modeling tools for plasma and fusion research and engineering
|
Lagrange FE operator definitions.
Data Types | |
| type | oft_lag_ginterp |
| Interpolate \( \nabla \) of a Lagrange field. More... | |
| type | oft_lag_rinterp |
| Interpolate a Lagrange field. More... | |
| type | oft_lag_vcinterp |
| Interpolate \( \nabla \times \) of a Lagrange vector field. More... | |
| type | oft_lag_vdinterp |
| Interpolate \( \nabla \) of a Lagrange vector field. More... | |
| type | oft_lag_vrinterp |
| Interpolate a Lagrange vector field. More... | |
| type | oft_lag_zerob |
| Zero Lagrange FE vector at all global boundary nodes. More... | |
| type | oft_lag_zerogrnd |
| Zero Lagrange FE vector at all global grounding node(s) More... | |
| type | oft_vlag_zerob |
| Zero all components of vector Lagrange FE vector at all global boundary nodes. More... | |
| type | oft_vlag_zeron |
| Zero normal component of vector Lagrange FE vector at all global boundary nodes. More... | |
| type | oft_vlag_zerot |
| Zero tangential component of vector Lagrange FE vector at all global boundary nodes. More... | |
Functions/Subroutines | |
| subroutine | lag_base_pop (self, acors, afine) |
| Transfer a base level Lagrange scalar field to the next MPI level. | |
| subroutine | lag_base_push (self, afine, acors) |
| Transfer a MPI level Lagrange scalar field to the base level. | |
| subroutine | lag_div (fe_rep, a, reg) |
| Compute the divergence of a Lagrange vector field. | |
| subroutine | lag_getlop_pre (ml_lag_rep, pre, mats, level, nlevels) |
| Construct default MG preconditioner for LAG::LOP. | |
| subroutine | lag_ginterp_apply (self, cell, f, gop, val) |
| Reconstruct the gradient of a Lagrange scalar field. | |
| subroutine | lag_lop_eigs (ml_lag_rep, minlev) |
| Compute eigenvalues and smoothing coefficients for the operator LAG::LOP. | |
| subroutine | lag_mloptions () |
| Read-in options for the basic Lagrange ML preconditioners. | |
| subroutine | lag_rinterp (self, cell, f, gop, val) |
| Reconstruct a Lagrange scalar field. | |
| subroutine | lag_rinterp_delete (self) |
| Destroy temporary internal storage and nullify references. | |
| subroutine | lag_rinterp_setup (self, lag_rep) |
| Setup interpolator for Lagrange scalar fields. | |
| subroutine | lag_setup_interp (ml_lag_rep, ml_vlag_rep) |
| Construct interpolation matrices on each MG level. | |
| subroutine | lag_vbc_diag (lag_rep, j_lag, bc_type, nn) |
| Get diagonal entries for a given boundary condition and desired node in a Lagrange vector field. | |
| subroutine | lag_vbc_tensor (lag_rep, j_lag, bc_type, nn) |
| Get boundary condition tensor for desired node in a Lagrange vector field. | |
| subroutine | lag_vcinterp (self, cell, f, gop, val) |
| Reconstruct the curl of a Lagrange vector field. | |
| subroutine | lag_vdinterp (self, cell, f, gop, val) |
| Reconstruct \( \nabla v \) for a Lagrange vector field. | |
| subroutine | lag_vrinterp (self, cell, f, gop, val) |
| Reconstruct a Lagrange vector field. | |
| subroutine | lag_vrinterp_delete (self) |
| Destroy temporary internal storage and nullify references. | |
| subroutine | lag_vrinterp_setup (self, lag_rep) |
| Setup interpolator for Lagrange vector fields. | |
| subroutine | oft_lag_getlop (fe_rep, mat, bc) |
| Construct laplacian matrix for Lagrange scalar representation. | |
| subroutine | oft_lag_getmop (fe_rep, mat, bc) |
| Construct mass matrix for Lagrange scalar representation. | |
| subroutine | oft_lag_getpdop (fe_rep, mat, field, bc, perp, be_flag) |
| Construct parallel diffusion matrix for Lagrange scalar representation. | |
| subroutine | oft_lag_project (fe_rep, field, x) |
| Project a scalar field onto a lagrange basis. | |
| subroutine | oft_lag_project_div (fe_rep, field, x) |
| Project the divergence of a scalar field onto a lagrange basis. | |
| subroutine | oft_lag_vgetmop (vlag_rep, mat, bc) |
| Construct mass matrix for Lagrange vector representation. | |
| subroutine | oft_lag_vproject (fe_rep, field, x) |
| Project a vector field onto a lagrange basis. | |
| subroutine | vzerob_apply (self, a) |
| Zero a surface Lagrange vector field at all edge nodes. | |
| subroutine | vzerob_delete (self) |
| Destroy temporary internal storage and nullify references. | |
| subroutine | vzeron_apply (self, a) |
| Zero normal component of a Lagrange vector field at boundary nodes. | |
| subroutine | vzerot_apply (self, a) |
| Zero tangential component of a Lagrange vector field at boundary nodes. | |
| subroutine | zerob_apply (self, a) |
| Zero a Lagrange scalar field at all boundary nodes. | |
| subroutine | zerob_delete (self) |
| Destroy temporary internal storage and nullify references. | |
| subroutine | zerogrnd_apply (self, a) |
| Zero a Lagrange scalar field at the global grounding node. | |
Variables | |
| real(r8), dimension(fem_max_levels) | df_lop =-1.d99 |
| real(r8), dimension(fem_max_levels), private | df_pdop =-1.d99 |
| integer(i4), dimension(fem_max_levels) | nu_lop =0 |
| integer(i4), dimension(fem_max_levels), private | nu_pdop =0 |
| real(r8), dimension(:,:), pointer | oft_lag_gop => NULL() |
| real(r8), dimension(:), pointer | oft_lag_rop => NULL() |
| subroutine lag_base_pop | ( | class(oft_ml_fe_vecspace), intent(inout) | self, |
| class(oft_vector), intent(inout) | acors, | ||
| class(oft_vector), intent(inout) | afine | ||
| ) |
Transfer a base level Lagrange scalar field to the next MPI level.
| [in,out] | acors | Vector to transfer |
| [in,out] | afine | Fine vector from transfer |
| subroutine lag_base_push | ( | class(oft_ml_fe_vecspace), intent(inout) | self, |
| class(oft_vector), intent(inout) | afine, | ||
| class(oft_vector), intent(inout) | acors | ||
| ) |
Transfer a MPI level Lagrange scalar field to the base level.
| [in,out] | afine | Needs docs |
| [in,out] | acors | Needs docs |
| subroutine lag_div | ( | class(oft_afem_type), intent(inout), target | fe_rep, |
| class(oft_vector), intent(inout) | a, | ||
| real(r8), intent(out) | reg | ||
| ) |
Compute the divergence of a Lagrange vector field.
| [in,out] | a | Input field |
| [out] | reg | \( \int_v \nabla \cdot a \; dV \) |
| subroutine lag_getlop_pre | ( | type(oft_ml_fem_type), intent(inout), target | ml_lag_rep, |
| class(oft_solver), intent(out), pointer | pre, | ||
| type(oft_matrix_ptr), dimension(:), intent(inout), pointer | mats, | ||
| integer(i4), intent(in), optional | level, | ||
| integer(i4), intent(in), optional | nlevels | ||
| ) |
Construct default MG preconditioner for LAG::LOP.
| [out] | pre | Needs docs |
| [in,out] | mats | Needs docs |
| [in] | level | Needs docs |
| [in] | nlevels | Needs docs |
| subroutine lag_ginterp_apply | ( | class(oft_lag_ginterp), intent(inout) | self, |
| integer(i4), intent(in) | cell, | ||
| real(r8), dimension(:), intent(in) | f, | ||
| real(r8), dimension(3,4), intent(in) | gop, | ||
| real(r8), dimension(:), intent(out) | val | ||
| ) |
Reconstruct the gradient of a Lagrange scalar field.
| [in] | cell | Cell for interpolation |
| [in] | f | Position in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed gradient at f [3] |
| subroutine lag_lop_eigs | ( | type(oft_ml_fem_type), intent(inout), target | ml_lag_rep, |
| integer(i4), intent(in) | minlev | ||
| ) |
Compute eigenvalues and smoothing coefficients for the operator LAG::LOP.
| [in] | minlev | Needs docs |
| subroutine lag_mloptions |
Read-in options for the basic Lagrange ML preconditioners.
| subroutine lag_rinterp | ( | class(oft_lag_rinterp), intent(inout) | self, |
| integer(i4), intent(in) | cell, | ||
| real(r8), dimension(:), intent(in) | f, | ||
| real(r8), dimension(3,4), intent(in) | gop, | ||
| real(r8), dimension(:), intent(out) | val | ||
| ) |
Reconstruct a Lagrange scalar field.
| [in] | cell | Cell for interpolation |
| [in] | f | Position in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed field at f [1] |
| subroutine lag_rinterp_delete | ( | class(oft_lag_rinterp), intent(inout) | self | ) |
Destroy temporary internal storage and nullify references.
| subroutine lag_rinterp_setup | ( | class(oft_lag_rinterp), intent(inout) | self, |
| class(oft_afem_type), intent(inout), target | lag_rep | ||
| ) |
Setup interpolator for Lagrange scalar fields.
Fetches local representation used for interpolation from vector object
| subroutine lag_setup_interp | ( | class(oft_ml_fem_type), intent(inout) | ml_lag_rep, |
| class(oft_ml_fem_comp_type), intent(inout), optional | ml_vlag_rep | ||
| ) |
Construct interpolation matrices on each MG level.
| subroutine lag_vbc_diag | ( | class(oft_scalar_fem), intent(inout), target | lag_rep, |
| integer(i4), intent(in) | j_lag, | ||
| integer(i4), intent(in) | bc_type, | ||
| real(r8), dimension(3,3), intent(out) | nn | ||
| ) |
Get diagonal entries for a given boundary condition and desired node in a Lagrange vector field.
| [in] | j_lag | Local index of Lagrange node |
| [in] | bc_type | Desired BC (1 -> norm, 2-> tang) |
| [out] | nn | Diagonal entries (3,3) |
| subroutine lag_vbc_tensor | ( | class(oft_scalar_fem), intent(inout), target | lag_rep, |
| integer(i4), intent(in) | j_lag, | ||
| integer(i4), intent(in) | bc_type, | ||
| real(r8), dimension(3,3), intent(out) | nn | ||
| ) |
Get boundary condition tensor for desired node in a Lagrange vector field.
| [in] | j_lag | Local index of Lagrange node |
| [in] | bc_type | Desired BC (1 -> norm, 2-> tang) |
| [out] | nn | BC tensor (3,3) |
| subroutine lag_vcinterp | ( | class(oft_lag_vcinterp), intent(inout) | self, |
| integer(i4), intent(in) | cell, | ||
| real(r8), dimension(:), intent(in) | f, | ||
| real(r8), dimension(3,4), intent(in) | gop, | ||
| real(r8), dimension(:), intent(out) | val | ||
| ) |
Reconstruct the curl of a Lagrange vector field.
| [in] | cell | Cell for interpolation |
| [in] | f | Position in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed field at f [3] |
| subroutine lag_vdinterp | ( | class(oft_lag_vdinterp), intent(inout) | self, |
| integer(i4), intent(in) | cell, | ||
| real(r8), dimension(:), intent(in) | f, | ||
| real(r8), dimension(3,4), intent(in) | gop, | ||
| real(r8), dimension(:), intent(out) | val | ||
| ) |
Reconstruct \( \nabla v \) for a Lagrange vector field.
The tensor is packed using reshape, to retrieve use
| [in] | cell | Cell for interpolation |
| [in] | f | Position in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed field at f [3] |
| subroutine lag_vrinterp | ( | class(oft_lag_vrinterp), intent(inout) | self, |
| integer(i4), intent(in) | cell, | ||
| real(r8), dimension(:), intent(in) | f, | ||
| real(r8), dimension(3,4), intent(in) | gop, | ||
| real(r8), dimension(:), intent(out) | val | ||
| ) |
Reconstruct a Lagrange vector field.
| [in] | cell | Cell for interpolation |
| [in] | f | Position in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed field at f [3] |
| subroutine lag_vrinterp_delete | ( | class(oft_lag_vrinterp), intent(inout) | self | ) |
Destroy temporary internal storage and nullify references.
| subroutine lag_vrinterp_setup | ( | class(oft_lag_vrinterp), intent(inout) | self, |
| class(oft_afem_type), intent(inout), target | lag_rep | ||
| ) |
Setup interpolator for Lagrange vector fields.
Fetches local representation used for interpolation from vector object
| subroutine oft_lag_getlop | ( | class(oft_afem_type), intent(inout), target | fe_rep, |
| class(oft_matrix), intent(inout), pointer | mat, | ||
| character(len=*), intent(in) | bc | ||
| ) |
Construct laplacian matrix for Lagrange scalar representation.
Supported boundary conditions
Full matrix -'zerob'Dirichlet for all boundary DOF -'grnd'` Dirichlet for only groundin point | [in,out] | mat | Matrix object |
| [in] | bc | Boundary condition |
| subroutine oft_lag_getmop | ( | class(oft_afem_type), intent(inout), target | fe_rep, |
| class(oft_matrix), intent(inout), pointer | mat, | ||
| character(len=*), intent(in) | bc | ||
| ) |
Construct mass matrix for Lagrange scalar representation.
Supported boundary conditions
Full matrix -'zerob'` Dirichlet for all boundary DOF | [in,out] | mat | Matrix object |
| [in] | bc | Boundary condition |
| subroutine oft_lag_getpdop | ( | class(oft_afem_type), intent(inout), target | fe_rep, |
| class(oft_matrix), intent(inout), pointer | mat, | ||
| class(fem_interp), intent(inout) | field, | ||
| character(len=*), intent(in) | bc, | ||
| real(r8), intent(in), optional | perp, | ||
| logical, dimension(:), intent(in), optional | be_flag | ||
| ) |
Construct parallel diffusion matrix for Lagrange scalar representation.
Supported boundary conditions
Full matrix -'zerob'Dirichlet for all boundary DOF -'grnd'` Dirichlet for only groundin point | [in,out] | mat | Matrix object |
| [in,out] | field | Vector field defining \( \hat{b} \) |
| [in] | bc | Boundary condition |
| [in] | perp | Value of perpendicular conductivity (optional) |
| [in] | be_flag | Flag for dirichlet nodes if different from boundary [ne] (optional) |
| subroutine oft_lag_project | ( | class(oft_afem_type), intent(inout), target | fe_rep, |
| class(fem_interp), intent(inout) | field, | ||
| class(oft_vector), intent(inout) | x | ||
| ) |
Project a scalar field onto a lagrange basis.
| [in,out] | field | Scalar field for projection |
| [in,out] | x | Field projected onto Lagrange basis |
| subroutine oft_lag_project_div | ( | class(oft_afem_type), intent(inout), target | fe_rep, |
| class(fem_interp), intent(inout) | field, | ||
| class(oft_vector), intent(inout) | x | ||
| ) |
Project the divergence of a scalar field onto a lagrange basis.
| [in,out] | field | Scalar field for projection |
| [in,out] | x | Field projected onto Lagrange basis |
| subroutine oft_lag_vgetmop | ( | class(oft_fem_comp_type), intent(inout), target | vlag_rep, |
| class(oft_matrix), intent(inout), pointer | mat, | ||
| character(len=*), intent(in) | bc | ||
| ) |
Construct mass matrix for Lagrange vector representation.
Supported boundary conditions
Full matrix -'all'Dirichlet for all components at boundary -'norm'Dirichlet for normal component at boundary -'tang'` Dirichlet for tangential component at boundary | [in,out] | mat | Matrix object |
| [in] | bc | Boundary condition |
| subroutine oft_lag_vproject | ( | class(oft_afem_type), intent(inout), target | fe_rep, |
| class(fem_interp), intent(inout) | field, | ||
| class(oft_vector), intent(inout) | x | ||
| ) |
Project a vector field onto a lagrange basis.
| [in,out] | field | Vector field for projection |
| [in,out] | x | Field projected onto Lagrange basis |
| subroutine vzerob_apply | ( | class(oft_vlag_zerob), intent(inout) | self, |
| class(oft_vector), intent(inout) | a | ||
| ) |
Zero a surface Lagrange vector field at all edge nodes.
| [in,out] | a | Field to be zeroed |
| subroutine vzerob_delete | ( | class(oft_vlag_zerob), intent(inout) | self | ) |
Destroy temporary internal storage and nullify references.
| subroutine vzeron_apply | ( | class(oft_vlag_zeron), intent(inout) | self, |
| class(oft_vector), intent(inout) | a | ||
| ) |
Zero normal component of a Lagrange vector field at boundary nodes.
| [in,out] | a | Field to be zeroed |
| subroutine vzerot_apply | ( | class(oft_vlag_zerot), intent(inout) | self, |
| class(oft_vector), intent(inout) | a | ||
| ) |
Zero tangential component of a Lagrange vector field at boundary nodes.
| [in,out] | a | Field to be zeroed |
| subroutine zerob_apply | ( | class(oft_lag_zerob), intent(inout) | self, |
| class(oft_vector), intent(inout) | a | ||
| ) |
Zero a Lagrange scalar field at all boundary nodes.
| [in,out] | a | Field to be zeroed |
| subroutine zerob_delete | ( | class(oft_lag_zerob), intent(inout) | self | ) |
Destroy temporary internal storage and nullify references.
| subroutine zerogrnd_apply | ( | class(oft_lag_zerogrnd), intent(inout) | self, |
| class(oft_vector), intent(inout) | a | ||
| ) |
Zero a Lagrange scalar field at the global grounding node.
| [in,out] | a | Field to be zeroed |
| real(r8), dimension(fem_max_levels) df_lop =-1.d99 |
|
private |
| integer(i4), dimension(fem_max_levels) nu_lop =0 |
|
private |
| real(r8), dimension(:,:), pointer oft_lag_gop => NULL() |
| real(r8), dimension(:), pointer oft_lag_rop => NULL() |