diff --git a/src/modules/BaseContinuity/src/BaseContinuity_Method.F90 b/src/modules/BaseContinuity/src/BaseContinuity_Method.F90 index 6a1b4b190..ce28c82a9 100644 --- a/src/modules/BaseContinuity/src/BaseContinuity_Method.F90 +++ b/src/modules/BaseContinuity/src/BaseContinuity_Method.F90 @@ -129,20 +129,40 @@ END SUBROUTINE BaseContinuity_Copy ! date: 2023-08-09 ! summary: Returns a string name of base interpolation type -FUNCTION BaseContinuity_ToString(obj) RESULT(ans) +FUNCTION BaseContinuity_ToString(obj, isUpper) RESULT(ans) CLASS(BaseContinuity_), INTENT(IN) :: obj + LOGICAL(LGT), OPTIONAL, INTENT(IN) :: isUpper TYPE(String) :: ans + ! internal variables + LOGICAL(LGT) :: isUpper0 + + isUpper0 = .FALSE. + IF (PRESENT(isUpper)) isUpper0 = isUpper + SELECT TYPE (obj) CLASS IS (H1_) ans = "H1" + CLASS IS (HCurl_) - ans = "HCurl" + IF (isUpper0) THEN + ans = "HCURL" + ELSE + ans = "HCurl" + END IF + CLASS IS (HDiv_) - ans = "HDiv" + IF (isUpper0) THEN + ans = "HDIV" + ELSE + ans = "HDiv" + END IF + CLASS IS (DG_) ans = "DG" + CLASS DEFAULT + CALL ErrorMsg(msg="NO CASE FOUND for type of obj", & routine="BaseContinuity_toString()", & line=__LINE__, unitno=stderr, file=__FILE__) diff --git a/src/modules/BaseInterpolation/src/BaseInterpolation_Method.F90 b/src/modules/BaseInterpolation/src/BaseInterpolation_Method.F90 index 0b77c81b1..7e424ed84 100644 --- a/src/modules/BaseInterpolation/src/BaseInterpolation_Method.F90 +++ b/src/modules/BaseInterpolation/src/BaseInterpolation_Method.F90 @@ -401,25 +401,52 @@ END SUBROUTINE BaseInterpolation_FromString ! date: 2023-08-09 ! summary: Returns a string name of base interpolation type -FUNCTION BaseInterpolation_ToString1(obj) RESULT(ans) +FUNCTION BaseInterpolation_ToString1(obj, isUpper) RESULT(ans) CLASS(BaseInterpolation_), INTENT(IN) :: obj + LOGICAL(LGT), OPTIONAL, INTENT(IN) :: isUpper TYPE(String) :: ans + ! internal variables + LOGICAL(LGT) :: isUpper0 + + isUpper0 = .FALSE. + IF (PRESENT(isUpper)) isUpper0 = isUpper + SELECT TYPE (obj) CLASS IS (LagrangeInterpolation_) - ans = "LagrangeInterpolation" + IF (isUpper0) THEN + ans = "LAGRANGEINTERPOLATION" + ELSE + ans = "LagrangeInterpolation" + END IF CLASS IS (SerendipityInterpolation_) - ans = "SerendipityInterpolation" + IF (isUpper0) THEN + ans = "SERENDIPITYINTERPOLATION" + ELSE + ans = "SerendipityInterpolation" + END IF CLASS IS (HermitInterpolation_) - ans = "HermitInterpolation" + IF (isUpper0) THEN + ans = "HERMITINTERPOLATION" + ELSE + ans = "HermitInterpolation" + END IF CLASS IS (HierarchyInterpolation_) - ans = "HierarchyInterpolation" + IF (isUpper0) THEN + ans = "HIERARCHYINTERPOLATION" + ELSE + ans = "HierarchyInterpolation" + END IF CLASS IS (OrthogonalInterpolation_) - ans = "OrthogonalInterpolation" + IF (isUpper0) THEN + ans = "ORTHOGONALINTERPOLATION" + ELSE + ans = "OrthogonalInterpolation" + END IF CLASS DEFAULT CALL ErrorMsg(msg="No Case Found For Type of obj2", & @@ -434,40 +461,87 @@ END FUNCTION BaseInterpolation_ToString1 ! BaseType_ToChar !---------------------------------------------------------------------------- -FUNCTION BaseType_ToChar(name) RESULT(ans) +FUNCTION BaseType_ToChar(name, isUpper) RESULT(ans) INTEGER(I4B), INTENT(IN) :: name + LOGICAL(LGT), OPTIONAL, INTENT(IN) :: isUpper CHARACTER(:), ALLOCATABLE :: ans + ! internal variable + LOGICAL(LGT) :: isUpper0 + + isUpper0 = .FALSE. + IF (PRESENT(isUpper)) isUpper0 = isUpper + SELECT CASE (name) CASE (poly%monomial) - ans = "Monomial" + IF (isUpper0) THEN + ans = "MONOMIAL" + ELSE + ans = "Monomial" + END IF CASE (poly%lagrange) - ans = "LagrangeInterpolation" + IF (isUpper0) THEN + ans = "LAGRANGEINTERPOLATION" + ELSE + ans = "LagrangeInterpolation" + END IF CASE (poly%serendipity) - ans = "SerendipityInterpolation" + IF (isUpper0) THEN + ans = "SERENDIPITYINTERPOLATION" + ELSE + ans = "SerendipityInterpolation" + END IF CASE (poly%hermit) - ans = "HermitInterpolation" + IF (isUpper0) THEN + ans = "HERMITINTERPOLATION" + ELSE + ans = "HermitInterpolation" + END IF CASE (poly%hierarchical) - ans = "HierarchyInterpolation" + IF (isUpper0) THEN + ans = "HIERARCHYINTERPOLATION" + ELSE + ans = "HierarchyInterpolation" + END IF CASE (poly%orthogonal) - ans = "OrthogonalInterpolation" + IF (isUpper0) THEN + ans = "ORTHOGONALINTERPOLATION" + ELSE + ans = "OrthogonalInterpolation" + END IF CASE (poly%legendre) - ans = "LegendreInterpolation" + IF (isUpper0) THEN + ans = "LEGENDREINTERPOLATION" + ELSE + ans = "LegendreInterpolation" + END IF CASE (poly%jacobi) - ans = "JacobiInterpolation" + IF (isUpper0) THEN + ans = "JACOBIINTERPOLATION" + ELSE + ans = "JacobiInterpolation" + END IF CASE (poly%ultraspherical) - ans = "UltrasphericalInterpolation" + IF (isUpper0) THEN + ans = "ULTRASPHERICALINTERPOLATION" + ELSE + ans = "UltrasphericalInterpolation" + END IF CASE (poly%chebyshev) - ans = "ChebyshevInterpolation" + IF (isUpper0) THEN + ans = "CHEBYSHEVINTERPOLATION" + ELSE + ans = "ChebyshevInterpolation" + END IF CASE DEFAULT CALL ErrorMsg(msg="No Case Found For name "//tostring(name), & @@ -482,83 +556,175 @@ END FUNCTION BaseType_ToChar ! QuadraturePointIDToName !---------------------------------------------------------------------------- -FUNCTION BaseInterpolation_ToString2(name) RESULT(ans) +FUNCTION BaseInterpolation_ToString2(name, isUpper) RESULT(ans) INTEGER(I4B), INTENT(IN) :: name + LOGICAL(LGT), OPTIONAL, INTENT(IN) :: isUpper TYPE(String) :: ans - ans = BaseInterpolation_ToChar(name) + ans = BaseInterpolation_ToChar(name=name, isUpper=isUpper) END FUNCTION BaseInterpolation_ToString2 !---------------------------------------------------------------------------- ! BaseInterpolation_ToChar !---------------------------------------------------------------------------- -FUNCTION BaseInterpolation_ToChar(name) RESULT(ans) +FUNCTION BaseInterpolation_ToChar(name, isUpper) RESULT(ans) INTEGER(I4B), INTENT(IN) :: name + LOGICAL(LGT), OPTIONAL, INTENT(IN) :: isUpper CHARACTER(:), ALLOCATABLE :: ans + ! internal varibles + LOGICAL(LGT) :: isUpper0 + + isUpper0 = .FALSE. + IF (PRESENT(isUpper)) isUpper0 = isUpper + SELECT CASE (name) CASE (ip%equidistance) - ans = "Equidistance" + IF (isUpper0) THEN + ans = "EQUIDISTANCE" + ELSE + ans = "Equidistance" + END IF CASE (ip%GaussLegendre) - ans = "GaussLegendre" + IF (isUpper0) THEN + ans = "GAUSSLEGENDRE" + ELSE + ans = "GaussLegendre" + END IF CASE (ip%GaussLegendreLobatto) - ans = "GaussLegendreLobatto" + IF (isUpper0) THEN + ans = "GAUSSLEGENDRELOBATTO" + ELSE + ans = "GaussLegendreLobatto" + END IF CASE (ip%GaussLegendreRadau) - ans = "GaussLegendreRadau" + IF (isUpper0) THEN + ans = "GAUSSLEGENDRERADAU" + ELSE + ans = "GaussLegendreRadau" + END IF CASE (ip%GaussLegendreRadauLeft) - ans = "GaussLegendreRadauLeft" + IF (isUpper0) THEN + ans = "GAUSSLEGENDRERADAULEFT" + ELSE + ans = "GaussLegendreRadauLeft" + END IF CASE (ip%GaussLegendreRadauRight) - ans = "GaussLegendreRadauRight" + IF (isUpper0) THEN + ans = "GAUSSLEGENDRERADAURIGHT" + ELSE + ans = "GaussLegendreRadauRight" + END IF CASE (ip%GaussChebyshev) - ans = "GaussChebyshev" + IF (isUpper0) THEN + ans = "GAUSSCHEBYSHEV" + ELSE + ans = "GaussChebyshev" + END IF CASE (ip%GaussChebyshevLobatto) - ans = "GaussChebyshevLobatto" + IF (isUpper0) THEN + ans = "GAUSSCHEBYSHEVLOBATTO" + ELSE + ans = "GaussChebyshevLobatto" + END IF CASE (ip%GaussChebyshevRadau) - ans = "GaussChebyshevRadau" + IF (isUpper0) THEN + ans = "GAUSSCHEBYSHEVRADAU" + ELSE + ans = "GaussChebyshevRadau" + END IF CASE (ip%GaussChebyshevRadauLeft) - ans = "GaussChebyshevRadauLeft" + IF (isUpper0) THEN + ans = "GAUSSCHEBYSHEVRADAULEFT" + ELSE + ans = "GaussChebyshevRadauLeft" + END IF CASE (ip%GaussChebyshevRadauRight) - ans = "GaussChebyshevRadauRight" + IF (isUpper0) THEN + ans = "GAUSSCHEBYSHEVRADAURIGHT" + ELSE + ans = "GaussChebyshevRadauRight" + END IF CASE (ip%GaussJacobi) - ans = "GaussJacobi" + IF (isUpper0) THEN + ans = "GAUSSJACOBI" + ELSE + ans = "GaussJacobi" + END IF CASE (ip%GaussJacobiLobatto) - ans = "GaussJacobiLobatto" + IF (isUpper0) THEN + ans = "GAUSSJACOBILOBATTO" + ELSE + ans = "GaussJacobiLobatto" + END IF CASE (ip%GaussJacobiRadau) - ans = "GaussJacobiRadau" + IF (isUpper0) THEN + ans = "GAUSSJACOBIRADAU" + ELSE + ans = "GaussJacobiRadau" + END IF CASE (ip%GaussJacobiRadauLeft) - ans = "GaussJacobiRadauLeft" + IF (isUpper0) THEN + ans = "GAUSSJACOBIRADAULEFT" + ELSE + ans = "GaussJacobiRadauLeft" + END IF CASE (ip%GaussJacobiRadauRight) - ans = "GaussJacobiRadauRight" + IF (isUpper0) THEN + ans = "GAUSSJACOBIRADAURIGHT" + ELSE + ans = "GaussJacobiRadauRight" + END IF CASE (ip%GaussUltraspherical) - ans = "GaussUltraspherical" + IF (isUpper0) THEN + ans = "GAUSSULTRASPHERICAL" + ELSE + ans = "GaussUltraspherical" + END IF CASE (ip%GaussUltrasphericalLobatto) - ans = "GaussUltrasphericalLobatto" + IF (isUpper0) THEN + ans = "GAUSSULTRASPHERICALLOBATTO" + ELSE + ans = "GaussUltrasphericalLobatto" + END IF CASE (ip%GaussUltrasphericalRadau) - ans = "GaussUltrasphericalRadau" + IF (isUpper0) THEN + ans = "GAUSSULTRASPHERICALRADAU" + ELSE + ans = "GaussUltrasphericalRadau" + END IF CASE (ip%GaussUltrasphericalRadauLeft) - ans = "GaussUltrasphericalRadauLeft" + IF (isUpper0) THEN + ans = "GAUSSULTRASPHERICALRADAULEFT" + ELSE + ans = "GaussUltrasphericalRadauLeft" + END IF CASE (ip%GaussUltrasphericalRadauRight) - ans = "GaussUltrasphericalRadauRight" + IF (isUpper0) THEN + ans = "GAUSSULTRASPHERICALRADAURIGHT" + ELSE + ans = "GaussUltrasphericalRadauRight" + END IF CASE DEFAULT CALL Errormsg(msg="No case found for given quadratureType name", & @@ -570,4 +736,8 @@ FUNCTION BaseInterpolation_ToChar(name) RESULT(ans) END FUNCTION BaseInterpolation_ToChar +!---------------------------------------------------------------------------- +! +!---------------------------------------------------------------------------- + END MODULE BaseInterpolation_Method diff --git a/src/modules/ElemshapeData/src/ElemshapeData_Method.F90 b/src/modules/ElemshapeData/src/ElemshapeData_Method.F90 index faff724e8..841d55eda 100644 --- a/src/modules/ElemshapeData/src/ElemshapeData_Method.F90 +++ b/src/modules/ElemshapeData/src/ElemshapeData_Method.F90 @@ -20,6 +20,15 @@ MODULE ElemshapeData_Method USE ElemshapeData_GetMethods USE ElemshapeData_GradientMethods +! USE ElemshapeData_H1Methods +! USE ElemshapeData_HCurlMethods +! USE ElemshapeData_HDivMethods +! USE ElemshapeData_DGMethods + +USE ElemshapeData_Lagrange +USE ElemshapeData_Hierarchical +USE ElemshapeData_Orthogonal + USE ElemshapeData_HRGNParamMethods USE ElemshapeData_HRQIParamMethods USE ElemshapeData_HminHmaxMethods @@ -32,14 +41,4 @@ MODULE ElemshapeData_Method USE ElemshapeData_StabilizationParamMethods USE ElemshapeData_UnitNormalMethods -! USE ElemshapeData_H1Methods -! USE ElemshapeData_HCurlMethods -! USE ElemshapeData_HDivMethods -! USE ElemshapeData_DGMethods - -USE ElemshapeData_Lagrange -USE ElemshapeData_Hierarchical -USE ElemshapeData_Orthogonal - - END MODULE ElemshapeData_Method diff --git a/src/modules/IntVector/src/IntVector_ConstructorMethod.F90 b/src/modules/IntVector/src/IntVector_ConstructorMethod.F90 index cd3af48cd..d96efafe3 100644 --- a/src/modules/IntVector/src/IntVector_ConstructorMethod.F90 +++ b/src/modules/IntVector/src/IntVector_ConstructorMethod.F90 @@ -31,6 +31,8 @@ MODULE IntVector_ConstructorMethod PUBLIC :: IntVector PUBLIC :: IntVector_Pointer PUBLIC :: Convert +PUBLIC :: Copy +PUBLIC :: Copy_ !---------------------------------------------------------------------------- ! Shape@Constructor @@ -76,10 +78,10 @@ END FUNCTION intVec_Size ! This function returns the total dimension (or rank) of an array, INTERFACE GetTotalDimension - MODULE PURE FUNCTION IntVec_getTotalDimension(obj) RESULT(Ans) + MODULE PURE FUNCTION intVec_getTotalDimension(obj) RESULT(Ans) TYPE(IntVector_), INTENT(IN) :: obj INTEGER(I4B) :: ans - END FUNCTION IntVec_getTotalDimension + END FUNCTION intVec_getTotalDimension END INTERFACE GetTotalDimension !---------------------------------------------------------------------------- @@ -371,4 +373,121 @@ MODULE PURE SUBROUTINE obj_convert_int(From, To) END SUBROUTINE obj_convert_int END INTERFACE Convert +!---------------------------------------------------------------------------- +! Copy@Constructor +!---------------------------------------------------------------------------- + +!> author: Vikas Sharma, Ph. D. +! date: 2025-06-22 +! summary: Copy y into x, x will be reallocated +! +!# Introduction +! +! Get the size of y and reallocate x to the same size. +! If x is already allocated, it will be reallocated to the size of y. + +INTERFACE Copy + MODULE PURE SUBROUTINE obj_Copy_Int8(x, y) + INTEGER(INT8), INTENT(INOUT), ALLOCATABLE :: x(:) + INTEGER(INT8), INTENT(IN) :: y(:) + END SUBROUTINE obj_Copy_Int8 +END INTERFACE Copy + +!---------------------------------------------------------------------------- +! Copy@Constructor +!---------------------------------------------------------------------------- + +!> author: Vikas Sharma, Ph. D. +! date: 2025-06-22 +! summary: Copy y into x, x will be reallocated +! +!# Introduction +! +! Get the size of y and reallocate x to the same size. +! If x is already allocated, it will be reallocated to the size of y. + +INTERFACE Copy + MODULE PURE SUBROUTINE obj_Copy_Int16(x, y) + INTEGER(INT16), INTENT(INOUT), ALLOCATABLE :: x(:) + INTEGER(INT16), INTENT(IN) :: y(:) + END SUBROUTINE obj_Copy_Int16 +END INTERFACE Copy + +!---------------------------------------------------------------------------- +! Copy@Constructor +!---------------------------------------------------------------------------- + +!> author: Vikas Sharma, Ph. D. +! date: 2025-06-22 +! summary: Copy y into x, x will be reallocated +! +! Introduction +! +! Get the size of y and reallocate x to the same size. +! If x is already allocated, it will be reallocated to the size of y. + +INTERFACE Copy + MODULE PURE SUBROUTINE obj_Copy_Int32(x, y) + INTEGER(INT32), INTENT(INOUT), ALLOCATABLE :: x(:) + INTEGER(INT32), INTENT(IN) :: y(:) + END SUBROUTINE obj_Copy_Int32 +END INTERFACE Copy + +!---------------------------------------------------------------------------- +! Copy@Constructor +!---------------------------------------------------------------------------- + +!> author: Vikas Sharma, Ph. D. +! date: 2025-06-22 +! summary: Copy y into x, x will be reallocated +! +!# Introduction +! +! Get the size of y and reallocate x to the same size. +! If x is already allocated, it will be reallocated to the size of y. + +INTERFACE Copy + MODULE PURE SUBROUTINE obj_Copy_Int64(x, y) + INTEGER(INT64), INTENT(INOUT), ALLOCATABLE :: x(:) + INTEGER(INT64), INTENT(IN) :: y(:) + END SUBROUTINE obj_Copy_Int64 +END INTERFACE Copy + +!---------------------------------------------------------------------------- +! Copy@Constructor +!---------------------------------------------------------------------------- + +!> author: Vikas Sharma, Ph. D. +! date: 2025-06-22 +! summary: Copy portion of y into x + +INTERFACE Copy_ + MODULE PURE SUBROUTINE obj_Copy1_(x, x_start, y, y_start, y_end) + INTEGER(I4B), INTENT(INOUT) :: x(:) + !! x vector should be allocated + INTEGER(I4B), INTENT(IN) :: y(:) + INTEGER(I4B), INTENT(IN) :: x_start + INTEGER(I4B), INTENT(IN) :: y_start, y_end + END SUBROUTINE obj_Copy1_ +END INTERFACE Copy_ + +!---------------------------------------------------------------------------- +! Copy@Constructor +!---------------------------------------------------------------------------- + +!> author: Vikas Sharma, Ph. D. +! date: 2025-06-22 +! summary: Copy y into x + +INTERFACE Copy_ + MODULE PURE SUBROUTINE obj_Copy2_(x, y) + INTEGER(I4B), INTENT(INOUT) :: x(:) + INTEGER(I4B), INTENT(IN) :: y(:) + END SUBROUTINE obj_Copy2_ +END INTERFACE Copy_ + +!---------------------------------------------------------------------------- +! +!---------------------------------------------------------------------------- + END MODULE IntVector_ConstructorMethod diff --git a/src/modules/QuadraturePoint/src/QuadraturePoint_Method.F90 b/src/modules/QuadraturePoint/src/QuadraturePoint_Method.F90 index f3d5492ff..2d288dd6d 100755 --- a/src/modules/QuadraturePoint/src/QuadraturePoint_Method.F90 +++ b/src/modules/QuadraturePoint/src/QuadraturePoint_Method.F90 @@ -28,6 +28,8 @@ MODULE QuadraturePoint_Method PRIVATE PUBLIC :: Initiate +PUBLIC :: Copy +PUBLIC :: ASSIGNMENT(=) PUBLIC :: QuadraturePoint PUBLIC :: QuadraturePoint_Pointer PUBLIC :: DEALLOCATE @@ -67,12 +69,13 @@ END FUNCTION QuadraturePointNameToId !> author: Vikas Sharma, Ph. D. ! date: 2023-08-06 -! summary: Quadrature point name to quadrature point id +! summary: Convert Quadrature point from int id to string name INTERFACE - MODULE FUNCTION QuadraturePointIdToName(name) RESULT(ans) + MODULE FUNCTION QuadraturePointIdToName(name, isUpper) RESULT(ans) INTEGER(I4B), INTENT(IN) :: name TYPE(String) :: ans + LOGICAL, INTENT(IN), OPTIONAL :: isUpper END FUNCTION QuadraturePointIdToName END INTERFACE @@ -80,10 +83,15 @@ END FUNCTION QuadraturePointIdToName ! QuadraturePoint_ToChar@ConstructorMethods !---------------------------------------------------------------------------- +!> author: Vikas Sharma, Ph. D. +! date: 2025-06-18 +! summary: Convert Quadrature poitn from int id to char name + INTERFACE - MODULE FUNCTION QuadraturePoint_ToChar(name) RESULT(ans) + MODULE FUNCTION QuadraturePoint_ToChar(name, isUpper) RESULT(ans) INTEGER(I4B), INTENT(IN) :: name - TYPE(String) :: ans + LOGICAL(LGT), OPTIONAL, INTENT(IN) :: isUpper + CHARACTER(:), ALLOCATABLE :: ans END FUNCTION QuadraturePoint_ToChar END INTERFACE @@ -110,6 +118,29 @@ MODULE FUNCTION obj_QuadratureNumber1(topo, order, quadratureType) & END FUNCTION obj_QuadratureNumber1 END INTERFACE QuadratureNumber +!---------------------------------------------------------------------------- +! Copy@ConstructorMethods +!---------------------------------------------------------------------------- + +!> author: Vikas Sharma, Ph. D. +! date: 23 July 2021 +! summary: This routine Initiates the quadrature points + +INTERFACE Initiate + MODULE PURE SUBROUTINE obj_Copy(obj, obj2) + TYPE(QuadraturePoint_), INTENT(INOUT) :: obj + TYPE(QuadraturePoint_), INTENT(IN) :: obj2 + END SUBROUTINE obj_Copy +END INTERFACE Initiate + +INTERFACE Copy + MODULE PROCEDURE obj_Copy +END INTERFACE Copy + +INTERFACE ASSIGNMENT(=) + MODULE PROCEDURE obj_Copy +END INTERFACE ASSIGNMENT(=) + !---------------------------------------------------------------------------- ! Initiate@ConstructorMethods !---------------------------------------------------------------------------- diff --git a/src/submodules/ElemshapeData/src/ElemshapeData_SetMethods@Methods.F90 b/src/submodules/ElemshapeData/src/ElemshapeData_SetMethods@Methods.F90 index c81920b61..085d4e2ca 100644 --- a/src/submodules/ElemshapeData/src/ElemshapeData_SetMethods@Methods.F90 +++ b/src/submodules/ElemshapeData/src/ElemshapeData_SetMethods@Methods.F90 @@ -50,9 +50,8 @@ MODULE PROCEDURE elemsd_SetBarycentricCoord INTEGER(I4B) :: nns -nns = SIZE(N, 1) -obj%coord(1:obj%nsd, 1:obj%nips) = MATMUL(val(1:obj%nsd, 1:nns), & - N(1:nns, 1:obj%nips)) +obj%coord(1:obj%nsd, 1:obj%nips) = MATMUL(val(1:obj%nsd, :), & + N(:, 1:obj%nips)) END PROCEDURE elemsd_SetBarycentricCoord !---------------------------------------------------------------------------- @@ -142,7 +141,7 @@ MODULE PROCEDURE elemsd_SetJacobian obj%jacobian(1:obj%nsd, 1:obj%xidim, 1:obj%nips) = & - MATMUL(val(1:obj%nsd, 1:obj%nns), dNdXi(1:obj%nns, 1:obj%xidim, 1:obj%nips)) + MATMUL(val(1:obj%nsd, :), dNdXi(:, 1:obj%xidim, 1:obj%nips)) END PROCEDURE elemsd_SetJacobian !---------------------------------------------------------------------------- @@ -151,8 +150,8 @@ MODULE PROCEDURE stsd_SetJacobian obj%jacobian(1:obj%nsd, 1:obj%xidim, 1:obj%nips) = & - MATMUL(MATMUL(val(1:obj%nsd, 1:obj%nns, :), T), & - dNdXi(1:obj%nns, 1:obj%xidim, 1:obj%nips)) + MATMUL(MATMUL(val(1:obj%nsd, :, :), T), & + dNdXi(:, 1:obj%xidim, 1:obj%nips)) END PROCEDURE stsd_SetJacobian !---------------------------------------------------------------------------- diff --git a/src/submodules/IntVector/src/IntVector_ConstructorMethod@Methods.F90 b/src/submodules/IntVector/src/IntVector_ConstructorMethod@Methods.F90 index bcbeb6ae0..cefd52e88 100644 --- a/src/submodules/IntVector/src/IntVector_ConstructorMethod@Methods.F90 +++ b/src/submodules/IntVector/src/IntVector_ConstructorMethod@Methods.F90 @@ -241,4 +241,80 @@ END IF END PROCEDURE obj_convert_int +!---------------------------------------------------------------------------- +! Copy +!---------------------------------------------------------------------------- + +MODULE PROCEDURE obj_Copy_Int8 +INTEGER(I4B) :: tsize, ii +tsize = SIZE(y) +CALL Reallocate(x, tsize) +DO ii = 1, tsize + x(ii) = y(ii) +END DO + +END PROCEDURE obj_Copy_Int8 + +!---------------------------------------------------------------------------- +! Copy +!---------------------------------------------------------------------------- + +MODULE PROCEDURE obj_Copy_Int16 +INTEGER(I4B) :: tsize, ii +tsize = SIZE(y) +CALL Reallocate(x, tsize) +DO ii = 1, tsize + x(ii) = y(ii) +END DO +END PROCEDURE obj_Copy_Int16 + +!---------------------------------------------------------------------------- +! Copy +!---------------------------------------------------------------------------- + +MODULE PROCEDURE obj_Copy_Int32 +INTEGER(I4B) :: tsize, ii +tsize = SIZE(y) +CALL Reallocate(x, tsize) +DO ii = 1, tsize + x(ii) = y(ii) +END DO +END PROCEDURE obj_Copy_Int32 + +!---------------------------------------------------------------------------- +! Copy +!---------------------------------------------------------------------------- + +MODULE PROCEDURE obj_Copy_Int64 +INTEGER(I4B) :: tsize, ii +tsize = SIZE(y) +CALL Reallocate(x, tsize) +DO ii = 1, tsize + x(ii) = y(ii) +END DO +END PROCEDURE obj_Copy_Int64 + +!---------------------------------------------------------------------------- +! Copy +!---------------------------------------------------------------------------- + +MODULE PROCEDURE obj_Copy1_ +INTEGER(I4B) :: xx, yy + +DO yy = y_start, y_end + xx = x_start + yy - y_start + x(xx) = y(yy) +END DO +END PROCEDURE obj_Copy1_ + +!---------------------------------------------------------------------------- +! Copy +!---------------------------------------------------------------------------- + +MODULE PROCEDURE obj_Copy2_ +INTEGER(I4B) :: tsize +tsize = SIZE(y) +CALL obj_Copy1_(x=x, y=y, x_start=1, y_start=1, y_end=tsize) +END PROCEDURE obj_Copy2_ + END SUBMODULE Methods diff --git a/src/submodules/QuadraturePoint/src/QuadraturePoint_Method@ConstructorMethods.F90 b/src/submodules/QuadraturePoint/src/QuadraturePoint_Method@ConstructorMethods.F90 index 38be36089..cf90c2a59 100755 --- a/src/submodules/QuadraturePoint/src/QuadraturePoint_Method@ConstructorMethods.F90 +++ b/src/submodules/QuadraturePoint/src/QuadraturePoint_Method@ConstructorMethods.F90 @@ -145,6 +145,20 @@ END PROCEDURE obj_QuadratureNumber1 +!---------------------------------------------------------------------------- +! Copy +!---------------------------------------------------------------------------- + +MODULE PROCEDURE obj_Copy +INTEGER(I4B) :: s(2) +obj%txi = obj2%txi +IF (ALLOCATED(obj2%points)) THEN + s = SHAPE(obj2%points) + CALL Reallocate(obj%points, s(1), s(2)) + obj%points(1:s(1), 1:s(2)) = obj2%points(1:s(1), 1:s(2)) +END IF +END PROCEDURE obj_Copy + !---------------------------------------------------------------------------- ! Initiate !----------------------------------------------------------------------------