| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| integer, | public | :: | deg |
input parameter: degree of Spline/polynomial |
|||
| integer, | public | :: | degGP |
number of Gauss-points (degGP+1) per element >= deg |
|||
| integer, | public | :: | continuity |
input parameter: full spline (=deg-1) or discontinuous (=-1) |
|||
| integer, | public | :: | nGP |
global number of gausspoints = (degGP+1)*nElems |
|||
| integer, | public | :: | nGP_str |
local number of gausspoints = (degGP+1)*nElems per MPI subdomain |
|||
| integer, | public | :: | nGP_end |
local number of gausspoints = (degGP+1)*nElems per MPI subdomain |
|||
| integer, | public | :: | nbase |
total number of degree of freedom / global basis functions |
|||
| class(sll_c_spline_matrix), | public, | ALLOCATABLE | :: | mass | |||
| logical, | public | :: | initialized | = | .FALSE. |
set to true in init, set to false in free |
|
| class(t_sGrid), | public, | POINTER | :: | grid |
pointer to grid |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | xi_GP(:) |
element local gauss point positions for interval [-1,1], size(0:degGP) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | w_GPloc(:) |
element local gauss weights for interval [-1,1], size(0:degGP) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | w_GP(:) |
global radial integration weight size((degGP+1)*nElems) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | s_GP(:) |
global position of gauss points in s [0,1] , size((degGP+1)*nElems) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | s_IP(:) |
position of interpolation points for initialization, size(nBase) |
||
| integer, | public, | ALLOCATABLE | :: | base_offset(:) |
offset of 0:deg element local basis functions to global index of degree of freedom, allocated (1:nElems). iBase = offset(iElem)+j, j=0...deg |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | base_GP(:,:,:) |
basis functions, (0:degGP,0:deg,1:nElems), |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | base_ds_GP(:,:,:) |
s derivative of basis functions, (0:degGP,0:deg,1:nElems) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | base_dsAxis(:,:) |
all derivatives 1..deg of all basis functions at axis size(1:deg+1,0:deg) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | base_dsEdge(:,:) |
all derivatives 1..deg of all basis functions at edge size(nBase-deg:nBase,0:deg) |
||
| integer, | public, | ALLOCATABLE | :: | nDOF_BC(:) |
number of boundary dofs involved in bc of BC_TYPE, size(-(deg+1):6) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | A_Axis(:,:,:) |
matrix to apply boundary conditions after interpolation (direct) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | invA_Axis(:,:,:) |
inverse of A_Axis |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | R_Axis(:,:,:) |
matrix to apply boundary conditions for RHS (testfunction) size(1:deg+1,1:deg+1,-(deg+1):6) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | AR_Axis(:,:,:) |
matrix to apply boundary conditions via LGM with identity size(1:deg+1,1:deg+1,-(deg+1):6) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | A_Edge(:,:,:) |
matrix to apply boundary conditions after interpolation (direct) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | invA_Edge(:,:,:) |
inverse of A_Edge |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | R_Edge(:,:,:) |
matrix to apply boundary conditions for RHS size(nBase-deg:nBase,nBase-deg:nBase,6) |
||
| real(kind=wp), | public, | ALLOCATABLE | :: | AR_Edge(:,:,:) |
matrix to apply boundary conditions via LGM with identity size(nBase-deg:nBase,nBase-deg:nBase,6) possible Boundary conditions 1=1 : boundary left open 2=2 : 1st deriv fixed 3=3 : sol fixed (typically used at edge) 4=4 : all odd derivs = 0 5=5 : all even derivs & sol = 0 6=6 : all odd derivs & sol = 0 <=0: m-dependent BC: 0, m=0: 4 -1, m=1: 6 -2, m=2: 5 ... odd m<=deg: ANTISYMM + all derivatives 1,...,m-1 =0 ...even m<=deg: SYMMZERO + all derivatives 1,...,m-1 =0 |
||
| class(sll_c_bsplines), | public, | ALLOCATABLE | :: | bspl |
contains bspline functions |
||
| type(sll_t_spline_interpolator_1d), | public | :: | Interpol |
spline interpolator |
initialize the type sbase with polynomial degree, continuity ( -1: disc, 1: full) and number of gauss points per element
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(inout) | :: | sf |
self |
||
| integer, | intent(in) | :: | deg_in |
polynomial degree |
||
| integer, | intent(in) | :: | continuity_in | |||
| class(t_sGrid), | intent(in), | TARGET | :: | grid_in |
grid information |
|
| integer, | intent(in) | :: | degGP_in |
gauss quadrature points: nGP=degGP+1 per elements |
finalize the type sBase
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(inout) | :: | sf |
self |
copy onto sf <-- tocopy
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(inout) | :: | sf |
self |
||
| class(c_sbase), | intent(in) | :: | tocopy |
compare sf with input sbase
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| class(c_sbase), | intent(in) | :: | tocompare | |||
| logical, | intent(out), | optional | :: | is_same | ||
| logical, | intent(out), | optional | :: | cond_out(:) |
change data from old_sBase to self. using interpolations of the old data at the new interpolation points
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| class(c_sbase), | intent(in) | :: | old_sBase |
base of old_data |
||
| integer, | intent(in) | :: | iterDim |
iterate on first or second dimension or old_data/sf_data |
||
| real(kind=wp), | intent(in) | :: | old_data(:,:) | |||
| real(kind=wp), | intent(out) | :: | sf_data(:,:) |
evaluate sbase at position x [0,1], NOT EFFICIENT!!
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| real(kind=wp), | intent(in) | :: | x |
position [0,1] where to evaluate |
||
| integer, | intent(in) | :: | deriv | |||
| integer, | intent(out) | :: | iElem | |||
| real(kind=wp), | intent(out) | :: | base_x(:) |
all basis functions (0:deg) evaluated |
simply evaluate function or derivative at point x
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| real(kind=wp), | intent(in) | :: | x |
point positions in [0,1] |
||
| integer, | intent(in) | :: | deriv |
derivative (=0: solution) |
||
| real(kind=wp), | intent(in) | :: | DOFs(:) |
array of all degrees of freedom |
simply evaluate function or derivative at point x, for multiple DOF vectors
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| real(kind=wp), | intent(in) | :: | x |
point positions in [0,1] |
||
| integer, | intent(in) | :: | nd |
number of DOF vectors |
||
| integer, | intent(in) | :: | deriv |
derivative (=0: solution) |
||
| real(kind=wp), | intent(in) | :: | DOFs(1:sf%nBase,1:nd) |
simply evaluate function with a base or base derivative evaluated at a point and its corresponding iElem use together with sBase_eval(x) => iElem,base
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| integer, | intent(in) | :: | iElem |
element where evaluation point |
||
| real(kind=wp), | intent(in) | :: | base_x(0:sf%deg) |
evaluation of base or its derivative in element iElem at a point position |
||
| real(kind=wp), | intent(in) | :: | DOFs(:) |
degrees of freedom |
evaluate all degrees of freedom at all Gauss Points (deriv=0 solution, deriv=1 first derivative d/ds)
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| integer, | intent(in) | :: | deriv |
only 0 or 1 |
||
| real(kind=wp), | intent(in) | :: | DOFs(:) |
array of all degrees of freedom |
take values interpolated at sf%s_IP positions and give back the degrees of freedom
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| real(kind=wp), | intent(in) | :: | g_IP(:) |
interpolation values at s_IP positions [0,1] |
result of interpolation
apply boundary conditions at axis and edge, via solving the Lagrange multiplier problem: x_new=x_old & A*x_new = d
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| real(kind=wp), | intent(inout) | :: | DOFs(1:sf%nBase) |
DOFs with boundary conditions applied |
||
| integer, | intent(in) | :: | BC_Type(2) |
bc type on axis (1) and edge (2) |
||
| real(kind=wp), | intent(in) | :: | BC_Val(2) |
for dirichlet BC : value |
apply strong boundary conditions at axis and edge for solution update
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(t_sBase), | intent(in) | :: | sf |
self |
||
| real(kind=wp), | intent(inout) | :: | RHS(sf%nBase) |
DOFs with boundary conditions applied |
||
| integer, | intent(in) | :: | BC_Type(2) |
bc type on axis (1) and edge (2) |
TYPE,EXTENDS(t_sbase) :: t_sBase_spl CLASS(sll_c_bsplines),ALLOCATABLE :: bspl !! contains bspline functions TYPE(sll_t_spline_interpolator_1d):: Interpol !! spline interpolator END TYPE t_sBase_spl