enabled alternative (and soon standard) keywords grid (resolution) and size (dimension).

removed errors for odd resolution
This commit is contained in:
Martin Diehl 2013-04-08 14:22:32 +00:00
parent ad879ebcf9
commit 7b5a202e8c
2 changed files with 53 additions and 60 deletions

View File

@ -1547,9 +1547,9 @@ subroutine IO_error(error_ID,e,i,g,ext_msg)
case (842_pInt) case (842_pInt)
msg = 'homogenization in spectral mesh' msg = 'homogenization in spectral mesh'
case (843_pInt) case (843_pInt)
msg = 'resolution in spectral mesh' msg = 'grid in spectral mesh'
case (844_pInt) case (844_pInt)
msg = 'dimension in spectral mesh' msg = 'size in spectral mesh'
case (845_pInt) case (845_pInt)
msg = 'incomplete information in spectral mesh header' msg = 'incomplete information in spectral mesh header'
case (846_pInt) case (846_pInt)

View File

@ -346,8 +346,8 @@ module mesh
private :: & private :: &
#ifdef Spectral #ifdef Spectral
mesh_spectral_getResolution, & mesh_spectral_getGrid, &
mesh_spectral_getDimension, & mesh_spectral_getSize, &
mesh_spectral_getHomogenization, & mesh_spectral_getHomogenization, &
mesh_spectral_count_nodesAndElements, & mesh_spectral_count_nodesAndElements, &
mesh_spectral_count_cpElements, & mesh_spectral_count_cpElements, &
@ -453,10 +453,10 @@ subroutine mesh_init(ip,el)
call mesh_build_FEdata ! get properties of the different types of elements call mesh_build_FEdata ! get properties of the different types of elements
#ifdef Spectral #ifdef Spectral
call IO_open_file(fileUnit,geometryFile) ! parse info from geometry file... call IO_open_file(fileUnit,geometryFile) ! parse info from geometry file...
res = mesh_spectral_getResolution(fileUnit) res = mesh_spectral_getGrid(fileUnit)
res1_red = res(1)/2_pInt + 1_pInt res1_red = res(1)/2_pInt + 1_pInt
wgt = 1.0/real(product(res),pReal) wgt = 1.0/real(product(res),pReal)
geomdim = mesh_spectral_getDimension(fileUnit) geomdim = mesh_spectral_getSize(fileUnit)
homog = mesh_spectral_getHomogenization(fileUnit) homog = mesh_spectral_getHomogenization(fileUnit)
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
@ -776,10 +776,10 @@ endfunction mesh_cellCenterCoordinates
#ifdef Spectral #ifdef Spectral
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Reads resolution information from geometry file. If fileUnit is given, !> @brief Reads grid information from geometry file. If fileUnit is given,
!! assumes an opened file, otherwise tries to open the one specified in geometryFile !! assumes an opened file, otherwise tries to open the one specified in geometryFile
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function mesh_spectral_getResolution(fileUnit) function mesh_spectral_getGrid(fileUnit)
use IO, only: & use IO, only: &
IO_checkAndRewind, & IO_checkAndRewind, &
IO_open_file, & IO_open_file, &
@ -793,17 +793,17 @@ function mesh_spectral_getResolution(fileUnit)
geometryFile geometryFile
implicit none implicit none
integer(pInt), dimension(1_pInt + 7_pInt*2_pInt) :: positions ! for a,b c + 3 values + keyword integer(pInt), dimension(3) :: mesh_spectral_getGrid
integer(pInt), intent(in), optional :: fileUnit integer(pInt), intent(in), optional :: fileUnit
integer(pInt) :: headerLength = 0_pInt integer(pInt), dimension(1_pInt + 7_pInt*2_pInt) :: positions ! for a,b,c + 3 values + keyword
integer(pInt), dimension(3) :: mesh_spectral_getResolution
integer(pInt) :: headerLength = 0_pInt
character(len=1024) :: line, & character(len=1024) :: line, &
keyword keyword
integer(pInt) :: i, j integer(pInt) :: i, j, myUnit
logical :: gotResolution = .false. logical :: gotGrid = .false.
integer(pInt) :: myUnit
mesh_spectral_getResolution = -1_pInt mesh_spectral_getGrid = -1_pInt
if(.not. present(fileUnit)) then if(.not. present(fileUnit)) then
myUnit = 289_pInt myUnit = 289_pInt
call IO_open_file(myUnit,trim(geometryFile)) call IO_open_file(myUnit,trim(geometryFile))
@ -819,23 +819,23 @@ function mesh_spectral_getResolution(fileUnit)
if (keyword(1:4) == 'head') then if (keyword(1:4) == 'head') then
headerLength = IO_intValue(line,positions,1_pInt) + 1_pInt headerLength = IO_intValue(line,positions,1_pInt) + 1_pInt
else else
call IO_error(error_ID=841_pInt, ext_msg='mesh_spectral_getResolution') call IO_error(error_ID=841_pInt, ext_msg='mesh_spectral_getGrid')
endif endif
rewind(myUnit) rewind(myUnit)
do i = 1_pInt, headerLength do i = 1_pInt, headerLength
read(myUnit,'(a1024)') line read(myUnit,'(a1024)') line
positions = IO_stringPos(line,7_pInt) positions = IO_stringPos(line,7_pInt)
select case ( IO_lc(IO_StringValue(line,positions,1_pInt,.true.)) ) select case ( IO_lc(IO_StringValue(line,positions,1_pInt,.true.)) )
case ('resolution') case ('resolution','grid')
gotResolution = .true. gotGrid = .true.
do j = 2_pInt,6_pInt,2_pInt do j = 2_pInt,6_pInt,2_pInt
select case (IO_lc(IO_stringValue(line,positions,j))) select case (IO_lc(IO_stringValue(line,positions,j)))
case('a') case('a')
mesh_spectral_getResolution(1) = IO_intValue(line,positions,j+1_pInt) mesh_spectral_getGrid(1) = IO_intValue(line,positions,j+1_pInt)
case('b') case('b')
mesh_spectral_getResolution(2) = IO_intValue(line,positions,j+1_pInt) mesh_spectral_getGrid(2) = IO_intValue(line,positions,j+1_pInt)
case('c') case('c')
mesh_spectral_getResolution(3) = IO_intValue(line,positions,j+1_pInt) mesh_spectral_getGrid(3) = IO_intValue(line,positions,j+1_pInt)
end select end select
enddo enddo
end select end select
@ -843,25 +843,21 @@ function mesh_spectral_getResolution(fileUnit)
if(.not. present(fileUnit)) close(myUnit) if(.not. present(fileUnit)) close(myUnit)
if (.not. gotResolution) & if (.not. gotGrid) &
call IO_error(error_ID = 845_pInt, ext_msg='resolution') call IO_error(error_ID = 845_pInt, ext_msg='grid')
if((mod(mesh_spectral_getResolution(1),2_pInt)/=0_pInt .or. & ! must be a even number if(mesh_spectral_getGrid(1) < 2_pInt .or. & ! must be at least 2
mesh_spectral_getResolution(1) < 2_pInt .or. & ! and larger than 1 mesh_spectral_getGrid(2) < 2_pInt .or. & ! -"-
mod(mesh_spectral_getResolution(2),2_pInt)/=0_pInt .or. & ! -"- mesh_spectral_getGrid(3) < 1_pInt) & ! must be at least 1
mesh_spectral_getResolution(2) < 2_pInt .or. & ! -"- call IO_error(error_ID = 843_pInt, ext_msg='mesh_spectral_getGrid')
(mod(mesh_spectral_getResolution(3),2_pInt)/=0_pInt .and. &
mesh_spectral_getResolution(3)/= 1_pInt)) .or. & ! third res might be 1
mesh_spectral_getResolution(3) < 1_pInt) &
call IO_error(error_ID = 843_pInt, ext_msg='mesh_spectral_getResolution')
end function mesh_spectral_getResolution end function mesh_spectral_getGrid
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Reads dimension information from geometry file. If fileUnit is given, !> @brief Reads size information from geometry file. If fileUnit is given,
!! assumes an opened file, otherwise tries to open the one specified in geometryFile !! assumes an opened file, otherwise tries to open the one specified in geometryFile
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function mesh_spectral_getDimension(fileUnit) function mesh_spectral_getSize(fileUnit)
use IO, only: & use IO, only: &
IO_checkAndRewind, & IO_checkAndRewind, &
IO_open_file, & IO_open_file, &
@ -875,17 +871,16 @@ function mesh_spectral_getDimension(fileUnit)
geometryFile geometryFile
implicit none implicit none
integer(pInt), dimension(1_pInt + 7_pInt*2_pInt) :: positions ! for a,b c + 3 values + keyword real(pReal), dimension(3) :: mesh_spectral_getSize
integer(pInt), intent(in), optional :: fileUnit integer(pInt), intent(in), optional :: fileUnit
integer(pInt) :: headerLength = 0_pInt integer(pInt), dimension(1_pInt + 7_pInt*2_pInt) :: positions ! for x,y,z + 3 values + keyword
real(pReal), dimension(3) :: mesh_spectral_getDimension integer(pInt) :: headerLength = 0_pInt
character(len=1024) :: line, & character(len=1024) :: line, &
keyword keyword
integer(pInt) :: i, j integer(pInt) :: i, j, myUnit
logical :: gotDimension = .false. logical :: gotSize = .false.
integer(pInt) :: myUnit
mesh_spectral_getDimension = -1.0_pReal mesh_spectral_getSize = -1.0_pReal
if(.not. present(fileUnit)) then if(.not. present(fileUnit)) then
myUnit = 289_pInt myUnit = 289_pInt
call IO_open_file(myUnit,trim(geometryFile)) call IO_open_file(myUnit,trim(geometryFile))
@ -901,23 +896,23 @@ function mesh_spectral_getDimension(fileUnit)
if (keyword(1:4) == 'head') then if (keyword(1:4) == 'head') then
headerLength = IO_intValue(line,positions,1_pInt) + 1_pInt headerLength = IO_intValue(line,positions,1_pInt) + 1_pInt
else else
call IO_error(error_ID=841_pInt, ext_msg='mesh_spectral_getDimension') call IO_error(error_ID=841_pInt, ext_msg='mesh_spectral_getSize')
endif endif
rewind(myUnit) rewind(myUnit)
do i = 1_pInt, headerLength do i = 1_pInt, headerLength
read(myUnit,'(a1024)') line read(myUnit,'(a1024)') line
positions = IO_stringPos(line,7_pInt) positions = IO_stringPos(line,7_pInt)
select case ( IO_lc(IO_StringValue(line,positions,1,.true.)) ) select case ( IO_lc(IO_StringValue(line,positions,1,.true.)) )
case ('dimension') case ('dimension', 'size')
gotDimension = .true. gotSize = .true.
do j = 2_pInt,6_pInt,2_pInt do j = 2_pInt,6_pInt,2_pInt
select case (IO_lc(IO_stringValue(line,positions,j))) select case (IO_lc(IO_stringValue(line,positions,j)))
case('x') case('x')
mesh_spectral_getDimension(1) = IO_floatValue(line,positions,j+1_pInt) mesh_spectral_getSize(1) = IO_floatValue(line,positions,j+1_pInt)
case('y') case('y')
mesh_spectral_getDimension(2) = IO_floatValue(line,positions,j+1_pInt) mesh_spectral_getSize(2) = IO_floatValue(line,positions,j+1_pInt)
case('z') case('z')
mesh_spectral_getDimension(3) = IO_floatValue(line,positions,j+1_pInt) mesh_spectral_getSize(3) = IO_floatValue(line,positions,j+1_pInt)
end select end select
enddo enddo
end select end select
@ -925,19 +920,19 @@ function mesh_spectral_getDimension(fileUnit)
if(.not. present(fileUnit)) close(myUnit) if(.not. present(fileUnit)) close(myUnit)
if (.not. gotDimension) & if (.not. gotSize) &
call IO_error(error_ID = 845_pInt, ext_msg='dimension') call IO_error(error_ID = 845_pInt, ext_msg='size')
if (any(mesh_spectral_getDimension<=0.0_pReal)) & if (any(mesh_spectral_getSize<=0.0_pReal)) &
call IO_error(error_ID = 844_pInt, ext_msg='mesh_spectral_getDimension') call IO_error(error_ID = 844_pInt, ext_msg='mesh_spectral_getSize')
end function mesh_spectral_getDimension end function mesh_spectral_getSize
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
!> @brief Reads homogenization information from geometry file. If fileUnit is given, !> @brief Reads homogenization information from geometry file. If fileUnit is given,
!! assumes an opened file, otherwise tries to open the one specified in geometryFile !! assumes an opened file, otherwise tries to open the one specified in geometryFile
!-------------------------------------------------------------------------------------------------- !--------------------------------------------------------------------------------------------------
function mesh_spectral_getHomogenization(fileUnit) integer(pInt) function mesh_spectral_getHomogenization(fileUnit)
use IO, only: & use IO, only: &
IO_checkAndRewind, & IO_checkAndRewind, &
IO_open_file, & IO_open_file, &
@ -950,15 +945,13 @@ function mesh_spectral_getHomogenization(fileUnit)
geometryFile geometryFile
implicit none implicit none
integer(pInt), intent(in), optional :: fileUnit
integer(pInt), dimension(1_pInt + 7_pInt*2_pInt) :: positions ! for a, b, c + 3 values + keyword integer(pInt), dimension(1_pInt + 7_pInt*2_pInt) :: positions ! for a, b, c + 3 values + keyword
integer(pInt), intent(in), optional :: fileUnit integer(pInt) :: headerLength = 0_pInt
integer(pInt) :: headerLength = 0_pInt
integer(pInt) :: mesh_spectral_getHomogenization
character(len=1024) :: line, & character(len=1024) :: line, &
keyword keyword
integer(pInt) :: i integer(pInt) :: i, myUnit
logical :: gotHomogenization = .false. logical :: gotHomogenization = .false.
integer(pInt) :: myUnit
mesh_spectral_getHomogenization = -1_pInt mesh_spectral_getHomogenization = -1_pInt
if(.not. present(fileUnit)) then if(.not. present(fileUnit)) then