New physical surfaces implemented

This commit is contained in:
Jorge Gonzalez 2026-02-20 09:50:42 +01:00
commit ce6b6a41a6
11 changed files with 122 additions and 99 deletions

View file

@ -98,7 +98,7 @@ MODULE moduleMesh1DCart
!EDGE FUNCTIONS
!Init edge element
SUBROUTINE initEdge1DCart(self, n, p, btPart, btEM)
SUBROUTINE initEdge1DCart(self, n, p, ps)
USE moduleSpecies, only: nSpecies
USE moduleErrors
IMPLICIT NONE
@ -106,8 +106,8 @@ MODULE moduleMesh1DCart
CLASS(meshEdge1DCart), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
INTEGER, INTENT(in):: btPart(1:nSpecies)
integer, intent(in):: btEM
INTEGER, INTENT(in):: ps
integer:: ps_index
REAL(8), DIMENSION(1:3):: r1
INTEGER:: s
@ -123,16 +123,12 @@ MODULE moduleMesh1DCart
self%normal = (/ 1.D0, 0.D0, 0.D0 /)
!Boundary particle linking
allocate(self%boundariesParticle(1:nSpecies))
!Assign functions to boundary
do s = 1, nSpecies
self%boundariesParticle(s)%obj => boundariesParticle(btPart(s))%obj
! Get Physical Surface index based on input numering
ps_index = physicalSurface_to_index(ps)
end do
! Add nodes to EM boundary
call boundariesEMLinking(btEM)%model%addNodes(self%getNodes(self%nNodes))
! Add elements to physical surface
call physicalSurfaces(ps_index)%addEdge(self%n)
call physicalSurfaces(ps_index)%addNode(self%n1%n)
end subroutine initEdge1DCart

View file

@ -98,7 +98,7 @@ MODULE moduleMesh1DRad
!EDGE FUNCTIONS
!Init edge element
SUBROUTINE initEdge1DRad(self, n, p, btPart, btEM)
SUBROUTINE initEdge1DRad(self, n, p, ps)
USE moduleSpecies, only:nSpecies
USE moduleErrors
IMPLICIT NONE
@ -106,8 +106,8 @@ MODULE moduleMesh1DRad
CLASS(meshEdge1DRad), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
INTEGER, INTENT(in):: btPart(1:nSpecies)
integer, intent(in):: btEM
INTEGER, INTENT(in):: ps
integer:: ps_index
REAL(8), DIMENSION(1:3):: r1
INTEGER:: s
@ -123,16 +123,12 @@ MODULE moduleMesh1DRad
self%normal = (/ 1.D0, 0.D0, 0.D0 /)
!Boundary particle linking
allocate(self%boundariesParticle(1:nSpecies))
!Assign functions to boundary
do s = 1, nSpecies
self%boundariesParticle(s)%obj => boundariesParticle(btPart(s))%obj
! Get Physical Surface index based on input numering
ps_index = physicalSurface_to_index(ps)
end do
! Add nodes to EM boundary
call boundariesEMLinking(btEM)%model%addNodes(self%getNodes(self%nNodes))
! Add elements to physical surface
call physicalSurfaces(ps_index)%addEdge(self%n)
call physicalSurfaces(ps_index)%addNode(self%n1%n)
end subroutine initEdge1DRad

View file

@ -133,7 +133,7 @@ MODULE moduleMesh2DCart
!EDGE FUNCTIONS
!Init edge element
SUBROUTINE initEdge2DCart(self, n, p, btPart, btEm)
SUBROUTINE initEdge2DCart(self, n, p, ps)
USE moduleSpecies, only:nSpecies
USE moduleErrors
IMPLICIT NONE
@ -141,8 +141,8 @@ MODULE moduleMesh2DCart
CLASS(meshEdge2DCart), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
INTEGER, INTENT(in):: btPart(1:nSpecies)
integer, intent(in):: btEM
INTEGER, INTENT(in):: ps
integer:: ps_index
REAL(8), DIMENSION(1:3):: r1, r2
INTEGER:: s
@ -162,16 +162,13 @@ MODULE moduleMesh2DCart
0.D0 /)
self%normal = self%normal/NORM2(self%normal)
!Boundary particle linking
allocate(self%boundariesParticle(1:nSpecies))
!Assign functions to boundary
do s = 1, nSpecies
self%boundariesParticle(s)%obj => boundariesParticle(btPart(s))%obj
! Get Physical Surface index based on input numering
ps_index = physicalSurface_to_index(ps)
end do
! Add nodes to EM boundary
call boundariesEMLinking(btEM)%model%addNodes(self%getNodes(self%nNodes))
! Add elements to physical surface
call physicalSurfaces(ps_index)%addEdge(self%n)
call physicalSurfaces(ps_index)%addNode(self%n1%n)
call physicalSurfaces(ps_index)%addNode(self%n2%n)
END SUBROUTINE initEdge2DCart

View file

@ -133,7 +133,7 @@ MODULE moduleMesh2DCyl
!EDGE FUNCTIONS
!Init edge element
SUBROUTINE initEdge2DCyl(self, n, p, btPart, btEM)
SUBROUTINE initEdge2DCyl(self, n, p, ps)
USE moduleSpecies
USE moduleErrors
USE moduleConstParam, ONLY: PI
@ -142,8 +142,8 @@ MODULE moduleMesh2DCyl
CLASS(meshEdge2DCyl), INTENT(out):: self
INTEGER, INTENT(in):: n
INTEGER, INTENT(in):: p(:)
INTEGER, INTENT(in):: btPart(1:nsPecies)
integer, intent(in):: btEM
INTEGER, INTENT(in):: ps
integer:: ps_index
REAL(8), DIMENSION(1:3):: r1, r2
INTEGER:: s
@ -171,16 +171,13 @@ MODULE moduleMesh2DCyl
0.D0 /)
self%normal = self%normal/NORM2(self%normal)
!Boundary particle linking
allocate(self%boundariesParticle(1:nSpecies))
!Assign functions to boundary
do s = 1, nSpecies
self%boundariesParticle(s)%obj => boundariesParticle(btPart(s))%obj
! Get Physical Surface index based on input numering
ps_index = physicalSurface_to_index(ps)
end do
! Add nodes to EM boundary
call boundariesEMLinking(btEM)%model%addNodes(self%getNodes(self%nNodes))
! Add elements to physical surface
call physicalSurfaces(ps_index)%addEdge(self%n)
call physicalSurfaces(ps_index)%addNode(self%n1%n)
call physicalSurfaces(ps_index)%addNode(self%n2%n)
END SUBROUTINE initEdge2DCyl

View file

@ -144,8 +144,14 @@ MODULE moduleMesh3DCart
self%surface = 1.D0/L_ref**2 !TODO: FIX THIS WHEN MOVING TO 3D
! Get Physical Surface index based on input numering
ps_index = physicalSurface_to_index(ps)
! Add elements to physical surface
call physicalSurfaces(ps_index)%addEdge(self%n)
call physicalSurfaces(ps_index)%addNode(self%n1%n)
call physicalSurfaces(ps_index)%addNode(self%n2%n)
call physicalSurfaces(ps_index)%addNode(self%n3%n)
END SUBROUTINE initEdge3DCartTria

View file

@ -37,7 +37,6 @@ MODULE moduleMeshInputGmsh2
REAL(8):: r(1:3) !3 generic coordinates
INTEGER, ALLOCATABLE:: p(:) !Array for nodes
INTEGER:: e = 0, n = 0, eTemp = 0, elemType = 0
integer:: btPart = 0, btEM = 0
INTEGER:: totalNumElem
INTEGER:: numEdges
INTEGER:: boundaryType
@ -199,11 +198,8 @@ MODULE moduleMeshInputGmsh2
END SELECT
!Associate boundary condition procedure.
btPart = physicalSurface_to_id(boundariesParticleLinking,boundaryType)
btEM = physicalSurface_to_id(boundariesEMLinking, boundaryType)
CALL self%edges(e)%obj%init(n, p, boundaryType)
CALL self%edges(e)%obj%init(n, p, boundariesParticleLinking(btPart)%speciesIndex, btEM)
DEALLOCATE(p)
END DO

View file

@ -47,7 +47,6 @@ module moduleMeshInputText
integer:: physicalID
integer:: n, c
integer, allocatable:: p(:)
integer:: btPart, btEM
fileID = 10
@ -141,10 +140,8 @@ module moduleMeshInputText
allocate(p(1))
p(1) = n
!Associate boundary condition procedure.
btPart = physicalSurface_to_id(boundariesParticleLinking, physicalID)
btEM = physicalSurface_to_id(boundariesEMLinking, physicalID)
call self%edges(physicalID)%obj%init(physicalID, p, boundariesParticleLinking(btPart)%speciesIndex, btEM)
call self%edges(physicalID)%obj%init(physicalID, p, physicalID)
deallocate(p)

View file

@ -173,7 +173,6 @@ MODULE moduleMeshInputVTU
REAL(8), ALLOCATABLE, DIMENSION(:):: coordinates
INTEGER:: n, e, c
INTEGER, ALLOCATABLE:: p(:)
INTEGER:: btPart, btEM
fileID = 10

View file

@ -71,7 +71,7 @@ MODULE moduleMesh
! Array of pointers to nodes.
TYPE:: meshNodePointer
CLASS(meshNode), POINTER:: obj
CONTAINS
contains
procedure, pass:: meshNodePointer_equal_type_type
procedure, pass:: meshNodePointer_equal_type_int
generic:: operator(==) => meshNodePointer_equal_type_type, meshNodePointer_equal_type_int
@ -136,14 +136,14 @@ MODULE moduleMesh
ABSTRACT INTERFACE
!Inits the edge parameters
subroutine initEdge_interface(self, n, p, physicalSurface)
subroutine initEdge_interface(self, n, p, ps)
use moduleSpecies, only:nSpecies
import:: meshEdge
class(meshEdge), intent(out):: self
integer, intent(in):: n
integer, intent(in):: p(:)
integer, intent(in):: physicalSurface
integer, intent(in):: ps
end subroutine initEdge_interface
@ -184,9 +184,29 @@ MODULE moduleMesh
! Array of pointers to edges.
type:: meshEdgePointer
class(meshEdge), pointer:: obj
contains
procedure, pass:: meshEdgePointer_equal_type_type
procedure, pass:: meshEdgePointer_equal_type_int
generic:: operator(==) => meshEdgePointer_equal_type_type, meshEdgePointer_equal_type_int
end type meshEdgePointer
interface
module function meshEdgePointer_equal_type_type(self, other) result(isEqual)
class(meshEdgePointer), intent(in):: self, other
logical:: isEqual
end function meshEdgePointer_equal_type_type
module function meshEdgePointer_equal_type_int(self, other) result(isEqual)
class(meshEdgePointer), intent(in):: self
integer, intent(in):: other
logical:: isEqual
end function meshEdgePointer_equal_type_int
end interface
!Parent of cell element
TYPE, PUBLIC, ABSTRACT, EXTENDS(meshElement):: meshCell
!Number of nodes in the cell

View file

@ -39,6 +39,27 @@ submodule(moduleMesh) elements
end function meshNodePointer_equal_type_int
function meshEdgePointer_equal_type_type(self, other) result(isEqual)
implicit none
class(meshEdgePointer), intent(in):: self, other
logical:: isEqual
isEqual = self%obj%n == other%obj%n
end function meshEdgePointer_equal_type_type
function meshEdgePointer_equal_type_int(self, other) result(isEqual)
implicit none
class(meshEdgePointer), intent(in):: self
integer, intent(in):: other
logical:: isEqual
isEqual = self%obj%n == other
end function meshEdgePointer_equal_type_int
!Constructs the global K matrix
SUBROUTINE constructGlobalK(self)
IMPLICIT NONE

View file

@ -24,15 +24,16 @@ submodule(moduleMesh) surfaces
class(physicalSurface), intent(inout):: self
integer, intent(in):: node
integer:: n, nn
integer:: nn
logical:: inArray
type(meshNodePointer):: nodeToAdd
associate(nodeToAdd => mesh%nodes(node)%obj)
nodeToAdd%obj => mesh%nodes(node)%obj
inArray = .false.
! I cannot use the procedure 'any' for this
do nn = 1 , size(self%nodes)
IF (self%nodes(nn)%obj == nodeToAdd) THEN
IF (self%nodes(nn) == nodeToAdd) THEN
! Node already in array
inArray = .true.
exit
@ -47,8 +48,6 @@ submodule(moduleMesh) surfaces
end if
end associate
end subroutine addNode
subroutine addEdge(self, edge)
@ -56,11 +55,12 @@ submodule(moduleMesh) surfaces
class(physicalSurface), intent(inout):: self
integer, intent(in):: edge
integer:: n, nn
integer:: nn
logical:: inArray
type(meshEdgePointer):: edgeToAdd
associate(edgeToAdd => mesh%edges(edge)%obj)
edgeToAdd %obj => mesh%edges(edge)%obj
inArray = .false.
! I cannot use the procedure 'any' for this
do nn = 1 , size(self%edges)
@ -79,8 +79,6 @@ submodule(moduleMesh) surfaces
end if
end associate
end subroutine addEdge
end submodule surfaces