161 lines
4.5 KiB
Fortran
161 lines
4.5 KiB
Fortran
submodule(moduleMesh) boundaryEM
|
|
CONTAINS
|
|
module function boundaryEMName_to_Index(boundaryName) result(bp)
|
|
use moduleErrors
|
|
implicit none
|
|
|
|
character(:), allocatable:: boundaryName
|
|
integer:: bp
|
|
integer:: b
|
|
|
|
bp = 0
|
|
do b = 1, nBoundariesEM
|
|
if (boundaryName == boundariesEM(b)%obj%name) then
|
|
bp = boundariesEM(b)%obj%n
|
|
|
|
end if
|
|
|
|
end do
|
|
|
|
if (bp == 0) then
|
|
call criticalError('Boundary ' // boundaryName // ' not found', 'boundaryEMName_to_Index')
|
|
|
|
end if
|
|
|
|
end function boundaryEMName_to_Index
|
|
|
|
! Initialize Dirichlet boundary condition
|
|
module SUBROUTINE initDirichlet(self, config, object)
|
|
use json_module
|
|
USE moduleRefParam, ONLY: Volt_ref
|
|
use moduleErrors
|
|
IMPLICIT NONE
|
|
|
|
CLASS(boundaryEMGeneric), ALLOCATABLE, INTENT(inout):: self
|
|
type(json_file), intent(inout):: config
|
|
character(:), allocatable, intent(in):: object
|
|
REAL(8):: potential
|
|
logical:: found
|
|
|
|
select type(self)
|
|
type is(boundaryEMDirichlet)
|
|
call config%get(object // '.potential', potential, found)
|
|
if (.not. found) then
|
|
call criticalError('Required parameter "potential" for Dirichlet boundary condition not found', 'readBoundaryEM')
|
|
|
|
end if
|
|
|
|
self%potential = potential / Volt_ref
|
|
|
|
end select
|
|
|
|
end subroutine initDirichlet
|
|
|
|
! Initialize Dirichlet boundary condition
|
|
module subroutine initDirichletTime(self, config, object)
|
|
use json_module
|
|
use moduleRefParam, ONLY: Volt_ref, ti_ref
|
|
use moduleErrors
|
|
implicit none
|
|
|
|
class(boundaryEMGeneric), allocatable, intent(inout):: self
|
|
type(json_file), intent(inout):: config
|
|
character(:), allocatable, intent(in):: object
|
|
real(8):: potential
|
|
character(:), allocatable:: temporalProfile
|
|
logical:: found
|
|
character(:), allocatable:: temporalProfilePath
|
|
|
|
select type(self)
|
|
type is(boundaryEMDirichletTime)
|
|
call config%get(object // '.potential', potential, found)
|
|
if (.not. found) then
|
|
call criticalError('Required parameter "potential" for Dirichlet Time boundary condition not found', &
|
|
'readBoundaryEM')
|
|
|
|
end if
|
|
|
|
call config%get(object // '.temporalProfile', temporalProfile, found)
|
|
if (.not. found) then
|
|
call criticalError('Required parameter "temporalProfile" for Dirichlet Time boundary condition not found', &
|
|
'readBoundaryEM')
|
|
|
|
end if
|
|
temporalProfilePath = path // temporalProfile
|
|
|
|
self%potential = potential / Volt_ref
|
|
|
|
call self%temporalProfile%init(temporalProfile)
|
|
|
|
call self%temporalProfile%convert(1.D0/ti_ref, 1.D0)
|
|
|
|
self%update => updateDirichletTime
|
|
|
|
end select
|
|
|
|
END SUBROUTINE initDirichletTime
|
|
|
|
!Apply Dirichlet boundary condition to the poisson equation
|
|
module SUBROUTINE applyDirichlet(self, vectorF)
|
|
USE moduleMesh
|
|
IMPLICIT NONE
|
|
|
|
CLASS(boundaryEMDirichlet), INTENT(in):: self
|
|
REAL(8), INTENT(inout):: vectorF(:)
|
|
integer:: n
|
|
|
|
DO n = 1, self%nNodes
|
|
self%nodes(n)%obj%emData%phi = self%potential
|
|
vectorF(self%nodes(n)%obj%n) = self%nodes(n)%obj%emData%phi
|
|
|
|
END DO
|
|
|
|
END SUBROUTINE applyDirichlet
|
|
|
|
!Apply Dirichlet boundary condition with time temporal profile
|
|
module SUBROUTINE applyDirichletTime(self, vectorF)
|
|
USE moduleMesh
|
|
USE moduleCaseParam, ONLY: timeStep, tauMin
|
|
IMPLICIT NONE
|
|
|
|
CLASS(boundaryEMDirichletTime), INTENT(in):: self
|
|
REAL(8), INTENT(inout):: vectorF(:)
|
|
integer:: n
|
|
|
|
DO n = 1, self%nNodes
|
|
self%nodes(n)%obj%emData%phi = self%potential * self%timeFactor
|
|
vectorF(self%nodes(n)%obj%n) = self%nodes(n)%obj%emData%phi
|
|
|
|
END DO
|
|
|
|
END SUBROUTINE applyDirichletTime
|
|
|
|
module subroutine updateDirichletTime(self)
|
|
implicit none
|
|
|
|
class(boundaryEMGeneric), intent(inout):: self
|
|
|
|
select type(self)
|
|
type is(boundaryEMDirichletTime)
|
|
self%timeFactor = self%temporalProfile%get(DBLE(timeStep)*tauMin)
|
|
|
|
end select
|
|
|
|
end subroutine updateDirichletTime
|
|
|
|
module subroutine boundariesEM_update()
|
|
implicit none
|
|
|
|
integer:: b
|
|
|
|
do b = 1, nBoundariesEM
|
|
if (associated(boundariesEM(b)%obj%update)) then
|
|
call boundariesEM(b)%obj%update
|
|
|
|
end if
|
|
|
|
end do
|
|
|
|
end subroutine boundariesEM_update
|
|
|
|
end submodule boundaryEM
|