|
The Open FUSION Toolkit 1.0.0-beta5
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... | |
Functions/Subroutines | |
| subroutine | lag_base_pop (acors, afine) |
| Transfer a base level Lagrange scalar field to the next MPI level. | |
| subroutine | lag_base_push (afine, acors) |
| Transfer a MPI level Lagrange scalar field to the base level. | |
| subroutine | lag_div (a, reg) |
| Compute the divergence of a Lagrange vector field. | |
| subroutine | lag_getlop_pre (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_ginterpmatrix (mat) |
| Construct interpolation matrix for polynomial levels. | |
| subroutine | lag_inject (afine, acors) |
| Inject a fine level Lagrange scalar field to the next coarsest level. | |
| subroutine | lag_interp (acors, afine) |
| Interpolate a coarse level Lagrange scalar field to the next finest level. | |
| subroutine | lag_lop_eigs (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_pinterpmatrix (mat) |
| Construct interpolation matrix for polynomial levels. | |
| subroutine | lag_rinterp (self, cell, f, gop, val) |
| Reconstruct a Lagrange scalar field. | |
| subroutine | lag_rinterp_delete (self) |
| Destroy temporary internal storage. | |
| subroutine | lag_rinterp_setup (self) |
| Setup interpolator for Lagrange scalar fields. | |
| subroutine | lag_setup_interp (create_vec) |
| Construct interpolation matrices on each MG level. | |
| subroutine | lag_vbc_diag (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 (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_vinject (afine, acors) |
| Inject a fine level Lagrange scalar field to the next coarsest level. | |
| subroutine | lag_vinterp (acors, afine) |
| Interpolate a coarse level Lagrange scalar field to the next finest level. | |
| subroutine | lag_vrinterp (self, cell, f, gop, val) |
| Reconstruct a Lagrange vector field. | |
| subroutine | lag_vrinterp_delete (self) |
| Setup interpolator for Lagrange vector fields. | |
| subroutine | lag_vrinterp_setup (self) |
| Setup interpolator for Lagrange vector fields. | |
| subroutine | lag_vzerob (a) |
| Zero a surface Lagrange vector field at all edge nodes. | |
| subroutine | lag_vzeron (a) |
| Zero normal component of a Lagrange vector field at boundary nodes. | |
| subroutine | lag_vzerot (a) |
| Zero tangential component of a Lagrange vector field at boundary nodes. | |
| subroutine | lag_zerob (a) |
| Zero a Lagrange scalar field at all boundary nodes. | |
| subroutine | lag_zerogrnd (a) |
| Zero a Lagrange scalar field at the global grounding node. | |
| subroutine | oft_lag_getlop (mat, bc) |
| Construct laplacian matrix for Lagrange scalar representation. | |
| subroutine | oft_lag_getmop (mat, bc) |
| Construct mass matrix for Lagrange scalar representation. | |
| subroutine | oft_lag_getpdop (mat, field, bc, perp, be_flag) |
| Construct parallel diffusion matrix for Lagrange scalar representation. | |
| subroutine | oft_lag_project (field, x) |
| Project a scalar field onto a lagrange basis. | |
| subroutine | oft_lag_project_div (field, x) |
| Project the divergence of a scalar field onto a lagrange basis. | |
| subroutine | oft_lag_vgetmop (mat, bc) |
| Construct mass matrix for Lagrange vector representation. | |
| subroutine | oft_lag_vproject (field, x) |
| Project a vector field onto a lagrange basis. | |
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_vector), intent(inout) | acors, |
| class(oft_vector), intent(inout) | afine | ||
| ) |
Transfer a base level Lagrange scalar field to the next MPI level.
| [in] | acors | Vector to transfer |
| [in,out] | afine | Fine vector from transfer |
| subroutine lag_base_push | ( | class(oft_vector), intent(inout) | afine, |
| class(oft_vector), intent(inout) | acors | ||
| ) |
Transfer a MPI level Lagrange scalar field to the base level.
| [in] | afine | Vector to transfer |
| [in,out] | acors | Fine vector from transfer |
| subroutine lag_div | ( | 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 | ( | 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.
| 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 | Possition in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed gradient at f [3] |
| subroutine lag_ginterpmatrix | ( | class(oft_matrix), intent(inout), pointer | mat | ) |
Construct interpolation matrix for polynomial levels.
| subroutine lag_inject | ( | class(oft_vector), intent(inout) | afine, |
| class(oft_vector), intent(inout) | acors | ||
| ) |
Inject a fine level Lagrange scalar field to the next coarsest level.
| [in] | afine | Vector to inject |
| [in,out] | acors | Coarse vector from injection |
| subroutine lag_interp | ( | class(oft_vector), intent(inout) | acors, |
| class(oft_vector), intent(inout) | afine | ||
| ) |
Interpolate a coarse level Lagrange scalar field to the next finest level.
| [in] | acors | Vector to interpolate |
| [in,out] | afine | Fine vector from interpolation |
| subroutine lag_lop_eigs | ( | integer(i4), intent(in) | 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_pinterpmatrix | ( | class(oft_matrix), intent(inout), pointer | mat | ) |
Construct interpolation matrix for polynomial levels.
| 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 | Possition 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.
| subroutine lag_rinterp_setup | ( | class(oft_lag_rinterp), intent(inout) | self | ) |
Setup interpolator for Lagrange scalar fields.
Fetches local representation used for interpolation from vector object
| subroutine lag_setup_interp | ( | logical, intent(in), optional | create_vec | ) |
Construct interpolation matrices on each MG level.
| subroutine lag_vbc_diag | ( | 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 | ( | 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 | Possition in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed curl 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 | Possition in cell in logical coord [4] |
| [in] | gop | Logical gradient vectors at f [3,4] |
| [out] | val | Reconstructed tensor at f [9] |
| subroutine lag_vinject | ( | class(oft_vector), intent(inout) | afine, |
| class(oft_vector), intent(inout) | acors | ||
| ) |
Inject a fine level Lagrange scalar field to the next coarsest level.
| [in] | afine | Vector to inject |
| [in,out] | acors | Coarse vector from injection |
| subroutine lag_vinterp | ( | class(oft_vector), intent(inout) | acors, |
| class(oft_vector), intent(inout) | afine | ||
| ) |
Interpolate a coarse level Lagrange scalar field to the next finest level.
| [in] | acors | Vector to interpolate |
| [in,out] | afine | Fine vector from interpolation |
| 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 | Possition 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 | ) |
Setup interpolator for Lagrange vector fields.
Fetches local representation used for interpolation from vector object
| subroutine lag_vrinterp_setup | ( | class(oft_lag_vrinterp), intent(inout) | self | ) |
Setup interpolator for Lagrange vector fields.
Fetches local representation used for interpolation from vector object
| subroutine lag_vzerob | ( | class(oft_vector), intent(inout) | a | ) |
Zero a surface Lagrange vector field at all edge nodes.
| [in,out] | a | Field to be zeroed |
| subroutine lag_vzeron | ( | 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 lag_vzerot | ( | 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 lag_zerob | ( | class(oft_vector), intent(inout) | a | ) |
Zero a Lagrange scalar field at all boundary nodes.
| [in,out] | a | Field to be zeroed |
| subroutine lag_zerogrnd | ( | class(oft_vector), intent(inout) | a | ) |
Zero a Lagrange scalar field at the global grounding node.
| [in,out] | a | Field to be zeroed |
| subroutine oft_lag_getlop | ( | class(oft_matrix), intent(inout), pointer | mat, |
| character(len=*), intent(in) | bc | ||
| ) |
Construct laplacian matrix for Lagrange scalar representation.
Supported boundary conditions
'none' 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_matrix), intent(inout), pointer | mat, |
| character(len=*), intent(in) | bc | ||
| ) |
Construct mass matrix for Lagrange scalar representation.
Supported boundary conditions
'none' Full matrix'zerob' Dirichlet for all boundary DOF| [in,out] | mat | Matrix object |
| [in] | bc | Boundary condition |
| subroutine oft_lag_getpdop | ( | 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
'none' 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(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(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_matrix), intent(inout), pointer | mat, |
| character(len=*), intent(in) | bc | ||
| ) |
Construct mass matrix for Lagrange vector representation.
Supported boundary conditions
'none' 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(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 |
| 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() |