First step of performance improvement

Finalysing first step of performance improvement focusing on reducing
iteration CPU time by improving calculation of basic element functions,
which took a lot of the CPU time
This commit is contained in:
Jorge Gonzalez 2023-01-06 21:02:54 +01:00
commit 746c5bea09
13 changed files with 260 additions and 252 deletions

View file

@ -47,7 +47,6 @@ MODULE moduleMesh2DCyl
CLASS(meshNode), POINTER:: n1 => NULL(), n2 => NULL(), n3 => NULL(), n4 => NULL()
!Connectivity to adjacent elements
CLASS(meshElement), POINTER:: e1 => NULL(), e2 => NULL(), e3 => NULL(), e4 => NULL()
REAL(8):: arNodes(1:4) = 0.D0
CONTAINS
!meshCell DEFERRED PROCEDURES
@ -67,7 +66,7 @@ MODULE moduleMesh2DCyl
PROCEDURE, PASS:: phy2log => phy2logQuad
PROCEDURE, PASS:: neighbourElement => neighbourElementQuad
!PARTICLUAR PROCEDURES
PROCEDURE, PASS, PRIVATE:: area => areaQuad
PROCEDURE, PASS, PRIVATE:: vol => volumeQuad
END TYPE meshCell2DCylQuad
@ -79,7 +78,6 @@ MODULE moduleMesh2DCyl
CLASS(meshNode), POINTER:: n1 => NULL(), n2 => NULL(), n3 => NULL()
!Connectivity to adjacent elements
CLASS(meshElement), POINTER:: e1 => NULL(), e2 => NULL(), e3 => NULL()
REAL(8):: arNodes(1:3) = 0.D0
CONTAINS
!meshCell DEFERRED PROCEDURES
@ -99,13 +97,13 @@ MODULE moduleMesh2DCyl
PROCEDURE, PASS:: phy2log => phy2logTria
PROCEDURE, PASS:: neighbourElement => neighbourElementTria
!PARTICULAR PROCEDURES
PROCEDURE, PASS, PRIVATE:: area => areaTria
PROCEDURE, PASS, PRIVATE:: vol => volumeTria
END TYPE meshCell2DCylTria
CONTAINS
!NODE FUNCTIONS
!Inits node element
!Init node element
SUBROUTINE initNode2DCyl(self, n, r)
USE moduleSpecies
USE moduleRefParam
@ -141,7 +139,7 @@ MODULE moduleMesh2DCyl
END FUNCTION getCoord2DCyl
!EDGE FUNCTIONS
!Inits edge element
!Init edge element
SUBROUTINE initEdge2DCyl(self, n, p, bt, physicalSurface)
USE moduleSpecies
USE moduleBoundary
@ -198,6 +196,7 @@ MODULE moduleMesh2DCyl
END FUNCTION getNodes2DCyl
!Calculate intersection between position and edge
PURE FUNCTION intersection2DCylEdge(self, r0) RESULT(r)
IMPLICIT NONE
@ -216,15 +215,15 @@ MODULE moduleMesh2DCyl
END FUNCTION intersection2DCylEdge
!Calculates a random position in edge
!Calculate a random position in edge
FUNCTION randPosEdge(self) RESULT(r)
USE moduleRandom
IMPLICIT NONE
CLASS(meshEdge2DCyl), INTENT(in):: self
REAL(8):: rnd
REAL(8):: dr, dz
REAL(8):: r(1:3)
REAL(8):: dr, dz
rnd = random()
dr = self%r(2) - self%r(1)
@ -245,7 +244,7 @@ MODULE moduleMesh2DCyl
!VOLUME FUNCTIONS
!QUAD FUNCTIONS
!Inits quadrilateral element
!Init element
SUBROUTINE initCellQuad2DCyl(self, n, p, nodes)
USE moduleRefParam
IMPLICIT NONE
@ -276,11 +275,7 @@ MODULE moduleMesh2DCyl
self%r = (/r1(2), r2(2), r3(2), r4(2)/)
!Assign node volume
CALL self%area()
self%n1%v = self%n1%v + self%arNodes(1)
self%n2%v = self%n2%v + self%arNodes(2)
self%n3%v = self%n3%v + self%arNodes(3)
self%n4%v = self%n4%v + self%arNodes(4)
CALL self%vol()
CALL OMP_INIT_LOCK(self%lock)
@ -289,7 +284,7 @@ MODULE moduleMesh2DCyl
END SUBROUTINE initCellQuad2DCyl
!Gets nodes from quadrilateral element
!Get nodes from quadrilateral element
PURE FUNCTION getNodesQuad(self, nNodes) RESULT(n)
IMPLICIT NONE
@ -297,7 +292,7 @@ MODULE moduleMesh2DCyl
INTEGER, INTENT(in):: nNodes
INTEGER:: n(1:nNodes)
n = (/self%n1%n, self%n2%n, self%n3%n, self%n4%n /)
n = (/ self%n1%n, self%n2%n, self%n3%n, self%n4%n /)
END FUNCTION getNodesQuad
@ -331,12 +326,12 @@ MODULE moduleMesh2DCyl
INTEGER, INTENT(in):: nNodes
REAL(8):: fPsi(1:nNodes)
fPsi = (/ (1.D0-Xi(1)) * (1.D0-Xi(2)), &
(1.D0+Xi(1)) * (1.D0-Xi(2)), &
(1.D0+Xi(1)) * (1.D0+Xi(2)), &
(1.D0-Xi(1)) * (1.D0+Xi(2)) /)
fPsi = (/ (1.D0 - Xi(1)) * (1.D0 - Xi(2)), &
(1.D0 + Xi(1)) * (1.D0 - Xi(2)), &
(1.D0 + Xi(1)) * (1.D0 + Xi(2)), &
(1.D0 - Xi(1)) * (1.D0 + Xi(2)) /)
fPsi = fPsi*0.25D0
fPsi = fPsi * 0.25D0
END FUNCTION fPsiQuad
@ -350,15 +345,15 @@ MODULE moduleMesh2DCyl
dPsi = 0.D0
dPsi(1,:) = (/ -(1.D0 - Xi(2)), &
(1.D0 - Xi(2)), &
(1.D0 + Xi(2)), &
-(1.D0 + Xi(2)) /)
dPsi(1, 1:4) = (/ -(1.D0 - Xi(2)), &
(1.D0 - Xi(2)), &
(1.D0 + Xi(2)), &
-(1.D0 + Xi(2)) /)
dPsi(2,:) = (/ -(1.D0 - Xi(1)), &
-(1.D0 + Xi(1)), &
(1.D0 + Xi(1)), &
(1.D0 - Xi(1)) /)
dPsi(2, 1:4) = (/ -(1.D0 - Xi(1)), &
-(1.D0 + Xi(1)), &
(1.D0 + Xi(1)), &
(1.D0 - Xi(1)) /)
dPsi = dPsi * 0.25D0
@ -379,6 +374,7 @@ MODULE moduleMesh2DCyl
DOT_PRODUCT(dPsi(2,1:4),self%z(1:4)) /)
pDer(2, 1:2) = (/ DOT_PRODUCT(dPsi(1,1:4),self%r(1:4)), &
DOT_PRODUCT(dPsi(2,1:4),self%r(1:4)) /)
pDer(3,3) = 1.D0
END FUNCTION partialDerQuad
@ -439,10 +435,11 @@ MODULE moduleMesh2DCyl
REAL(8):: invJ(1:3,1:3), detJ
INTEGER:: l, m
localK=0.D0
Xi=0.D0
localK = 0.D0
Xi = 0.D0
!Start 2D Gauss Quad Integral
DO l=1, 3
DO l = 1, 3
Xi(2) = corQuad(l)
DO m = 1, 3
Xi(1) = corQuad(m)
@ -479,8 +476,9 @@ MODULE moduleMesh2DCyl
INTEGER:: l, m
localF = 0.D0
Xi = 0.D0
DO l=1, 3
Xi = 0.D0
DO l = 1, 3
Xi(1) = corQuad(l)
DO m = 1, 3
Xi(2) = corQuad(m)
@ -489,7 +487,7 @@ MODULE moduleMesh2DCyl
detJ = self%detJac(pDer)
fPsi = self%fPsi(Xi, 4)
r = DOT_PRODUCT(fPsi, self%r)
f = DOT_PRODUCT(fPsi,source)
f = DOT_PRODUCT(fPsi, source)
localF = localF + r*f*fPsi*wQuad(l)*wQuad(m)*detJ
END DO
@ -498,7 +496,7 @@ MODULE moduleMesh2DCyl
END FUNCTION elemFQuad
!Checks if a particle is inside a quad element
!Checks if Xi is inside the element
PURE FUNCTION insideQuad(Xi) RESULT(ins)
IMPLICIT NONE
@ -510,7 +508,7 @@ MODULE moduleMesh2DCyl
END FUNCTION insideQuad
!Transforms physical coordinates to element coordinates
!Transform physical coordinates to element coordinates
PURE FUNCTION phy2logQuad(self,r) RESULT(Xi)
IMPLICIT NONE
@ -544,7 +542,7 @@ MODULE moduleMesh2DCyl
END FUNCTION phy2logQuad
!Get the next element for a logical position Xi
!Get the neighbour element for a logical position Xi
SUBROUTINE neighbourElementQuad(self, Xi, neighbourElement)
IMPLICIT NONE
@ -571,8 +569,8 @@ MODULE moduleMesh2DCyl
END SUBROUTINE neighbourElementQuad
!Computes element area
PURE SUBROUTINE areaQuad(self)
!Compute element volume
PURE SUBROUTINE volumeQuad(self)
USE moduleConstParam, ONLY: PI8
IMPLICIT NONE
@ -584,34 +582,33 @@ MODULE moduleMesh2DCyl
REAL(8):: dPsi(1:3, 1:4), pDer(1:3, 1:3)
self%volume = 0.D0
self%arNodes = 0.D0
!2D 1 point Gauss Quad Integral
Xi = 0.D0
dPsi = self%dPsi(Xi, 4)
pDer = self%partialDer(4, dPsi)
detJ = self%detJac(pDer)
fPsi = self%fPsi(Xi, 4)
!Computes total volume of the cell
r = DOT_PRODUCT(fPsi,self%r)
!Computes total volume of the cell
self%volume = r*detJ*PI8 !4*2*pi
!Computes volume per node
Xi = (/-5.D-1, -5.D-1, 0.D0/)
r = self%gatherF(Xi, 4, self%r)
self%arNodes(1) = fPsi(1)*self%volume
self%n1%v = self%n1%v + fPsi(1)*r*detJ*PI8
Xi = (/ 5.D-1, -5.D-1, 0.D0/)
r = self%gatherF(Xi, 4, self%r)
self%arNodes(2) = fPsi(2)*self%volume
self%n2%v = self%n2%v + fPsi(2)*r*detJ*PI8
Xi = (/ 5.D-1, 5.D-1, 0.D0/)
r = self%gatherF(Xi, 4, self%r)
self%arNodes(3) = fPsi(3)*self%volume
self%n3%v = self%n3%v + fPsi(3)*r*detJ*PI8
Xi = (/-5.D-1, 5.D-1, 0.D0/)
r = self%gatherF(Xi, 4, self%r)
self%arNodes(4) = fPsi(4)*self%volume
self%n4%v = self%n4%v + fPsi(4)*r*detJ*PI8
END SUBROUTINE areaQuad
END SUBROUTINE volumeQuad
!TRIA ELEMENT
!Init tria element
!TRIA FUNCTIONS
!Init element
SUBROUTINE initCellTria2DCyl(self, n, p, nodes)
USE moduleRefParam
IMPLICIT NONE
@ -639,10 +636,7 @@ MODULE moduleMesh2DCyl
self%z = (/r1(1), r2(1), r3(1)/)
self%r = (/r1(2), r2(2), r3(2)/)
!Assign node volume
CALL self%area()
self%n1%v = self%n1%v + self%arNodes(1)
self%n2%v = self%n2%v + self%arNodes(2)
self%n3%v = self%n3%v + self%arNodes(3)
CALL self%vol()
CALL OMP_INIT_LOCK(self%lock)
@ -651,7 +645,7 @@ MODULE moduleMesh2DCyl
END SUBROUTINE initCellTria2DCyl
!Gets node indexes from triangular element
!Random position in cell
PURE FUNCTION getNodesTria(self, nNodes) RESULT(n)
IMPLICIT NONE
@ -663,7 +657,7 @@ MODULE moduleMesh2DCyl
END FUNCTION getNodesTria
!Random position in quadrilateral volume
!Random position in cell
FUNCTION randPosCellTria(self) RESULT(r)
USE moduleRandom
IMPLICIT NONE
@ -685,7 +679,7 @@ MODULE moduleMesh2DCyl
END FUNCTION randPosCellTria
!Shape functions for triangular element
!Compute element functions in point Xi
PURE FUNCTION fPsiTria(Xi, nNodes) RESULT(fPsi)
IMPLICIT NONE
@ -699,7 +693,7 @@ MODULE moduleMesh2DCyl
END FUNCTION fPsiTria
!Derivative element function at coordinates Xi
!Compute element derivative functions in point Xi
PURE FUNCTION dPsiTria(Xi, nNodes) RESULT(dPsi)
IMPLICIT NONE
@ -714,6 +708,7 @@ MODULE moduleMesh2DCyl
END FUNCTION dPsiTria
!Compute the derivatives in global coordinates
PURE FUNCTION partialDerTria(self, nNodes, dPsi) RESULT(pDer)
IMPLICIT NONE
@ -731,6 +726,7 @@ MODULE moduleMesh2DCyl
END FUNCTION partialDerTria
!Gather electric field at position Xi
PURE FUNCTION gatherEFTria(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
@ -746,6 +742,7 @@ MODULE moduleMesh2DCyl
END FUNCTION gatherEFTria
!Gather magnetic field at position Xi
PURE FUNCTION gatherMFTria(self, Xi) RESULT(array)
IMPLICIT NONE
CLASS(meshCell2DCylTria), INTENT(in):: self
@ -769,7 +766,7 @@ MODULE moduleMesh2DCyl
END FUNCTION gatherMFTria
!Computes element local stiffness matrix
!Compute cell local stiffness matrix
PURE FUNCTION elemKTria(self, nNodes) RESULT(localK)
USE moduleConstParam, ONLY: PI2
IMPLICIT NONE
@ -784,7 +781,8 @@ MODULE moduleMesh2DCyl
REAL(8):: invJ(1:3,1:3), detJ
INTEGER:: l
localK=0.D0
localK = 0.D0
Xi=0.D0
!Start 2D Gauss Quad Integral
DO l=1, 4
@ -802,7 +800,7 @@ MODULE moduleMesh2DCyl
END FUNCTION elemKTria
!Computes element local source vector
!Compute element local source vector
PURE FUNCTION elemFTria(self, nNodes, source) RESULT(localF)
USE moduleConstParam, ONLY: PI2
IMPLICIT NONE
@ -838,6 +836,7 @@ MODULE moduleMesh2DCyl
END FUNCTION elemFTria
!Check if Xi is inside the element
PURE FUNCTION insideTria(Xi) RESULT(ins)
IMPLICIT NONE
@ -850,7 +849,7 @@ MODULE moduleMesh2DCyl
END FUNCTION insideTria
!Transforms physical coordinates to element coordinates
!Transform physical coordinates to element coordinates
PURE FUNCTION phy2logTria(self,r) RESULT(Xi)
IMPLICIT NONE
@ -873,6 +872,7 @@ MODULE moduleMesh2DCyl
END FUNCTION phy2logTria
!Get the neighbour cell for a logical position Xi
SUBROUTINE neighbourElementTria(self, Xi, neighbourElement)
IMPLICIT NONE
@ -896,8 +896,8 @@ MODULE moduleMesh2DCyl
END SUBROUTINE neighbourElementTria
!Calculates area for triangular element
PURE SUBROUTINE areaTria(self)
!Calculate volume for triangular element
PURE SUBROUTINE volumeTria(self)
USE moduleConstParam, ONLY: PI
IMPLICIT NONE
@ -909,23 +909,24 @@ MODULE moduleMesh2DCyl
REAL(8):: r
self%volume = 0.D0
self%arNodes = 0.D0
!2D 1 point Gauss Quad Integral
Xi = (/ 1.D0/3.D0, 1.D0/3.D0, 0.D0 /)
dPsi = self%dPsi(Xi, 3)
pDer = self%partialDer(3, dPsi)
detJ = self%detJac(pDer)
fPsi = self%fPsi(Xi, 3)
!Computes total volume of the cell
r = DOT_PRODUCT(fPsi, self%r)
!Computes total volume of the cell
self%volume = r*detJ*PI !2PI*1/2
!Computes volume per node
self%arNodes = fPsi*self%volume
self%n1%v = self%n1%v + fPsi(1)*self%volume
self%n2%v = self%n2%v + fPsi(2)*self%volume
self%n3%v = self%n3%v + fPsi(3)*self%volume
END SUBROUTINE areaTria
END SUBROUTINE volumeTria
!COMMON FUNCTIONS FOR CYLINDRICAL VOLUME ELEMENTS
!Computes element Jacobian determinant
!Compute element Jacobian determinant
PURE FUNCTION detJ2DCyl(pDer) RESULT(dJ)
IMPLICIT NONE
@ -936,7 +937,7 @@ MODULE moduleMesh2DCyl
END FUNCTION detJ2DCyl
!Computes element Jacobian inverse matrix (without determinant)
!Compute element Jacobian inverse matrix (without determinant)
PURE FUNCTION invJ2DCyl(pDer) RESULT(invJ)
IMPLICIT NONE