Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion .fortls
Original file line number Diff line number Diff line change
Expand Up @@ -11,4 +11,4 @@
"lowercase_intrinsics": false,
"debug_log": false,
"disable_diagnostics": false
}
}
4 changes: 2 additions & 2 deletions src/modules/BaseType/src/BaseType.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1161,9 +1161,9 @@ END SUBROUTINE lag_elem_refelem

!> authors: Vikas Sharma, Ph. D.
! date: 4 March 2021
! summary: This data type contains shapefunction related data defined at all gauss points of an elements
! summary: Datatype for data defined at all gauss points of an elements
!
!{!pages/ElemshapeData.md}
!{!pages/ElemshapeData_.md}
TYPE :: ElemShapeData_
REAL( DFP ), ALLOCATABLE :: N( :, : )
!! Shape function value `N(I,ips)`
Expand Down
144 changes: 105 additions & 39 deletions src/modules/ElemshapeData/src/ElemshapeData_Method.F90
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,7 @@

!> authors: Vikas Sharma, Ph. D.
! date: 1 March 2021
! summary: This module implements methods related to [[elemShapeData_]] datatype
! summary: Methods related to [[elemShapeData_]] datatype

MODULE ElemshapeData_Method
USE BaseType
Expand All @@ -30,18 +30,39 @@ MODULE ElemshapeData_Method

!> authors: Vikas Sharma, Ph. D.
! date: 4 March 2021
! summary: This subroutine allocate the memory for various matrices in the object
! summary: Allocate the memory for various matrices in the object
!
!# Introduction
! This subroutine allocates the memory for various matrices in the object. This subroutine belongs to the generic interface called `AllocateData()`.
!
! - This routine also belongs to generic routien called `initiate`
! This subroutine allocates the memory for various matrices in the object.
! This subroutine belongs to the generic interface called `AllocateData()`.
!
!@note
! This routine also belongs to generic routien called `initiate`
!@endnote
!
!### Usage
!
!```fortran
!
! PROGRAM main
! USE easifemBase
! IMPLICIT NONE
! TYPE( ElemshapeData_ ) :: obj
! TYPE( QuadraturePoint_ ) :: quad
! TYPE( ReferenceLine_ ) :: refelem
! INTEGER( I4B ), PARAMETER :: nsd=1, order=1
! !> main
! ! #ReferenceLine/ReferenceLine
! ! #ReferenceLine
! refelem = ReferenceLine( nsd = nsd )
! ! #GaussLegendreQuadrature
! quad = GaussLegendreQuadrature( refelem = refelem, order = order )
! CALL AllocateData( obj = obj, nsd = refelem%nsd, &
! & xidim = refelem%xidimension, nns = 2, nips = 1 )
! CALL Initiate( obj = obj, quad = quad, refelem = refelem, &
! & ContinuityType= typeH1, InterpolType = TypeLagrangeInterpolation )
! CALL Display( obj, "obj" )
! END PROGRAM main
!```

INTERFACE
Expand All @@ -65,6 +86,12 @@ END SUBROUTINE elemsd_AllocateData

PUBLIC :: AllocateData

INTERFACE Allocate
MODULE PROCEDURE elemsd_AllocateData
END INTERFACE Allocate

PUBLIC :: Allocate

!----------------------------------------------------------------------------
! Initiate@Constructor
!----------------------------------------------------------------------------
Expand All @@ -77,6 +104,29 @@ END SUBROUTINE elemsd_AllocateData
!
! This routine initiates the shape function related data inside the element.
!
!## Usage
!
!```fortran
! PROGRAM main
! USE easifemBase
! IMPLICIT NONE
! TYPE( ElemshapeData_ ) :: obj
! TYPE( QuadraturePoint_ ) :: quad
! TYPE( ReferenceLine_ ) :: refelem
! INTEGER( I4B ), PARAMETER :: nsd=1, order=1
! !> main
! ! #ReferenceLine/ReferenceLine
! ! #ReferenceLine
! refelem = ReferenceLine( nsd = nsd )
! ! #GaussLegendreQuadrature
! quad = GaussLegendreQuadrature( refelem = refelem, order = order )
! CALL AllocateData( obj = obj, nsd = refelem%nsd, &
! & xidim = refelem%xidimension, nns = 2, nips = 1 )
! CALL Initiate( obj = obj, quad = quad, refelem = refelem, &
! & ContinuityType= typeH1, InterpolType = TypeLagrangeInterpolation )
! CALL Display( obj, "obj" )
! END PROGRAM main
!```

INTERFACE
MODULE SUBROUTINE elemsd_initiate( obj, quad, refElem, continuityType, &
Expand All @@ -99,13 +149,17 @@ END SUBROUTINE elemsd_initiate

!> authors: Vikas Sharma, Ph. D.
! date: 4 March 2021
! summary: This subroutine initiate time shape function data in [[stelemshapedata_]]
! summary: Initiate time shape function data in [[stelemshapedata_]]
!
!# Introduction
!
! * This subroutine initiate the time shape function data in [[stelemshapedata_]].
! * For the effeciency purpose, user should supply an instance of [[Elemshapedata_]] on time element. This object will have information of location time shape function data such as `T, dTdtheta` etc.
! * This routine uses `elemsd` to set `obj%T`, `obj%dTdTheta`, `obj%Jt`, `obj%Wt`, `obj%Theta`.
! * This subroutine initiate the time shape function data in
! [[stelemshapedata_]].
! * For the effeciency purpose, user should supply an instance of
! [[Elemshapedata_]] on time element. This object will have information of
! location time shape function data such as `T, dTdtheta` etc.
! * This routine uses `elemsd` to set `obj%T`, `obj%dTdTheta`, `obj%Jt`,
! `obj%Wt`, `obj%Theta`.
! * The following examples shows how to use it.
!
!### Usage
Expand Down Expand Up @@ -134,10 +188,11 @@ END SUBROUTINE stsd_initiate

!> authors: Vikas Sharma, Ph. D.
! date: 4 March 2021
! summary: This subroutine deallocates the data stored inside [[elemshapedata_]]
! summary: Deallocates the data stored inside [[elemshapedata_]]
!
!# Introduction
! This routine deallocates the data stored inside [[elemshapedata_]]. This routine belongs to `AllocateData()`
! This routine deallocates the data stored inside [[elemshapedata_]]. This
! routine belongs to `AllocateData()`
!
!
!### Usage
Expand All @@ -153,7 +208,7 @@ END SUBROUTINE elemsd_DeallocateData
END INTERFACE

INTERFACE DeallocateData
MODULE PROCEDURE elemsd_DeallocateData
MODULE PROCEDURE elemsd_DeallocateData
END INTERFACE DeallocateData

PUBLIC :: DeallocateData
Expand Down Expand Up @@ -206,10 +261,11 @@ END FUNCTION elemsd_BaseContinuity

!> authors: Vikas Sharma, Ph. D.
! date: 4 March 2021
! summary: This subroutine display the content of [[elemshapedata_]] and [[stelemshapedata_]]
! summary: Display the content of [[elemshapedata_]] and [[stelemshapedata_]]
!
!# Introduction
! This subroutine displays the content of [[elemshapedata_]] and [[stelemshapedata_]] on screen. this routine belongs to `Display()`.
! This subroutine displays the content of [[elemshapedata_]] and
! [[stelemshapedata_]] on screen. this routine belongs to `Display()`.
!
!### Usage
!
Expand Down Expand Up @@ -345,7 +401,7 @@ END SUBROUTINE H1_Hierarchy


!----------------------------------------------------------------------------
! Initiate@H1DivLagrange
! Initiate@H1DivLagrange
!----------------------------------------------------------------------------

!> authors: Vikas Sharma, Ph. D.
Expand Down Expand Up @@ -373,7 +429,7 @@ END SUBROUTINE H1Div_Lagrange
END INTERFACE Initiate

!----------------------------------------------------------------------------
! Initiate@H1DivHermit
! Initiate@H1DivHermit
!----------------------------------------------------------------------------

!> authors: Vikas Sharma, Ph. D.
Expand Down Expand Up @@ -401,7 +457,7 @@ END SUBROUTINE H1Div_Hermit
END INTERFACE Initiate

!----------------------------------------------------------------------------
! Initiate@H1DivSerendipity
! Initiate@H1DivSerendipity
!----------------------------------------------------------------------------

!> authors: Vikas Sharma, Ph. D.
Expand Down Expand Up @@ -985,7 +1041,7 @@ END SUBROUTINE stsd_set_dNTdXt_internally
PUBLIC :: setdNTdXt

!----------------------------------------------------------------------------
! setValue@setMethod
! Set@setMethod
!----------------------------------------------------------------------------

!> authors: Vikas Sharma, Ph. D.
Expand All @@ -1011,19 +1067,19 @@ END SUBROUTINE stsd_set_dNTdXt_internally
!@endnote

INTERFACE
MODULE PURE SUBROUTINE set_value( obj, Val, N, dNdXi )
MODULE PURE SUBROUTINE elemsd_set1( obj, Val, N, dNdXi )
CLASS( ElemshapeData_ ), INTENT( INOUT ) :: obj
REAL( DFP ), INTENT( IN ) :: Val( :, : )
!! Spatial nodal coordinates
REAL( DFP ), INTENT( IN ) :: N(:, :)
!! Shape function for geometry
REAL( DFP ), INTENT( IN ) :: dNdXi( :, :, : )
!! Local derivative of shape functions for geometry
END SUBROUTINE set_value
END SUBROUTINE elemsd_set1
END INTERFACE

!----------------------------------------------------------------------------
! setValue@setMethod
! set@setMethod
!----------------------------------------------------------------------------

!> authors: Vikas Sharma, Ph. D.
Expand Down Expand Up @@ -1051,21 +1107,21 @@ END SUBROUTINE set_value
!@endnote

INTERFACE
MODULE PURE SUBROUTINE stsd_set_value( obj, Val, N, T, dNdXi )
MODULE PURE SUBROUTINE stelemsd_set1( obj, Val, N, T, dNdXi )
CLASS( STElemshapeData_ ), INTENT( INOUT ) :: obj
REAL( DFP ), INTENT( IN ) :: Val( :, :, : )
!! Spatial nodal coordinates
REAL( DFP ), INTENT( IN ) :: N( :, : )
REAL( DFP ), INTENT( IN ) :: T( : )
REAL( DFP ), INTENT( IN ) :: dNdXi( :, :, : )
END SUBROUTINE stsd_set_value
END SUBROUTINE stelemsd_set1
END INTERFACE

INTERFACE setValue
MODULE PROCEDURE set_value, stsd_set_value
END INTERFACE setValue
INTERFACE set
MODULE PROCEDURE elemsd_set1, stelemsd_set1
END INTERFACE set

PUBLIC :: setValue
PUBLIC :: set

!----------------------------------------------------------------------------
! setNormal@setMethod
Expand Down Expand Up @@ -1299,14 +1355,17 @@ END FUNCTION interpol_matrix
! getSTInterpolation@getMethod
!----------------------------------------------------------------------------

INTERFACE
!! This subroutine performs interpolations of scalar

!> authors: Dr. Vikas Sharma
!> authors: Vikas Sharma, Ph. D.
! date: 1 Nov 2021
! summary: This subroutine performs interpolations of scalar
!
!# Introduction
!
! This subroutine performs interpolation of a scalar from its space-time nodal
! values.
! $$u=u^{a}_{I}N^{I}T_{a}$$

INTERFACE
MODULE PURE SUBROUTINE stsd_get_interpol_scalar( obj, Interpol, Val )
CLASS( STElemshapeData_ ), INTENT( IN ) :: obj
REAL( DFP ), INTENT( INOUT ), ALLOCATABLE :: Interpol( : )
Expand All @@ -1332,15 +1391,18 @@ END SUBROUTINE stsd_get_interpol_fevar_scalar
! getSTInterpolation@getMethod
!----------------------------------------------------------------------------

INTERFACE
!! This subroutine performs interpolation of a vector

!> authors: Dr. Vikas Sharma
!> authors: Vikas Sharma, Ph. D.
! date: 1 Nov 2021
! summary: This subroutine performs interpolation of a vector
!
!# Introduction
!
! This subroutine performs interpolation of a vector from its space-time
! nodal values
!
! $$u_{i}=u^{a}_{iI}N^{I}T_{a}$$

INTERFACE
MODULE PURE SUBROUTINE stsd_get_interpol_vector( obj, Interpol, Val )
CLASS( STElemshapeData_ ), INTENT( IN ) :: obj
REAL( DFP ), INTENT( INOUT ), ALLOCATABLE :: Interpol( :, : )
Expand All @@ -1354,13 +1416,16 @@ END SUBROUTINE stsd_get_interpol_vector
! getSTInterpolation@getMethod
!----------------------------------------------------------------------------

INTERFACE
!! This subroutine performs interpolation of matrix

!> authors: Dr. Vikas Sharma
!> authors: Vikas Sharma, Ph. D.
! date: 1 Nov 2021
! summary: This subroutine performs interpolation of matrix
!
!# Introduction
!
! This subroutine performs interpolation of matrix from its space-time
! nodal values

INTERFACE
MODULE PURE SUBROUTINE stsd_get_interpol_matrix( obj, Interpol, Val )
CLASS( STElemshapeData_ ), INTENT( IN ) :: obj
REAL( DFP ), INTENT( INOUT ), ALLOCATABLE :: Interpol( :, :, : )
Expand Down Expand Up @@ -1517,6 +1582,7 @@ END SUBROUTINE getLocalGradient_vector
!
! $$\frac{\partial \phi }{\partial \xi_{i} } =\phi^{a}_{I} T_{a}\frac
!{\partial N^{I}}{\partial \xi_{i} }$$

MODULE PURE SUBROUTINE stsd_getLocalGradient_scalar( obj, dPhidXi, Val )
CLASS( STElemShapeData_ ), INTENT( IN ) :: obj
REAL( DFP ), ALLOCATABLE, INTENT( INOUT ) :: dPhidXi( :, : )
Expand Down
6 changes: 4 additions & 2 deletions src/modules/FEMatrix/src/FEMatrix_Method.F90
Original file line number Diff line number Diff line change
Expand Up @@ -70,11 +70,13 @@ END FUNCTION Space_MassMatrix
!
!# Introduction
!
! This subroutine makes space matrix in space domain, Here Rho $\rho$ is a finite element variable. Following expression can be evaluated
! This subroutine makes space matrix in space domain, Here Rho $\rho$ is a
! finite element variable. Following expression can be evaluated
!
! $$\int_{\Omega } N^{I}T_{a}\rho N^{J}T_{b}d\Omega$$
! $$\iint \frac{\partial N^{I}T_{a}}{\partial t} \rho N^{J}T_{b}d\Omega dt$$
! $$\iint \frac{\partial N^{I}T_{a}}{\partial t} \rho \frac{\partial N^{J}T_{b}}{\partial t} d\Omega dt$$
! $$\iint \frac{\partial N^{I}T_{a}}{\partial t} \rho
! \frac{\partial N^{J}T_{b}}{\partial t} d\Omega dt$$
! $$\iint N^{I}T_{a}\rho \frac{\partial N^{J}T_{b}}{\partial t} d\Omega dt$$
!
!### Usage
Expand Down
1 change: 1 addition & 0 deletions src/modules/GlobalData/src/GlobalData.F90
Original file line number Diff line number Diff line change
Expand Up @@ -479,6 +479,7 @@ MODULE GlobalData
INTEGER( I4B ), PARAMETER, PUBLIC :: Time= 3
INTEGER( I4B ), PARAMETER, PUBLIC :: SpaceTime = 4
INTEGER( I4B ), PARAMETER, PUBLIC :: SolutionDependent = 5
INTEGER( I4B ), PARAMETER, PUBLIC :: RandomSpace = 6
!>
INTEGER( I4B ), PARAMETER, PUBLIC :: Scalar = 1
INTEGER( I4B ), PARAMETER, PUBLIC :: Vector = 2
Expand Down
Loading