grid solver related changes
This commit is contained in:
parent
61f1cdf70d
commit
12cb043d6d
2
PRIVATE
2
PRIVATE
|
@ -1 +1 @@
|
|||
Subproject commit d61d62667fb683a61dcc41cd90194a2d9b279879
|
||||
Subproject commit 6c1a6993d31c62711553a5ba5d6b0cf5b6620634
|
34
src/IO.f90
34
src/IO.f90
|
@ -22,6 +22,8 @@ module IO
|
|||
character, parameter, public :: &
|
||||
IO_EOL = new_line('DAMASK'), & !< end of line character
|
||||
IO_COMMENT = '#'
|
||||
character(len=*), parameter :: LOWER = 'abcdefghijklmnopqrstuvwxyz'
|
||||
character(len=len(LOWER)), parameter :: UPPER = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
|
||||
character, parameter :: &
|
||||
CR = achar(13), &
|
||||
LF = IO_EOL
|
||||
|
@ -32,6 +34,7 @@ module IO
|
|||
IO_readlines, &
|
||||
IO_isBlank, &
|
||||
IO_wrapLines, &
|
||||
IO_postfix, &
|
||||
IO_strPos, &
|
||||
IO_strValue, &
|
||||
IO_intValue, &
|
||||
|
@ -238,6 +241,30 @@ pure function IO_strPos(str)
|
|||
end function IO_strPos
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief Append postfix to each indicator character that is followed by a lowercase letter.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
function IO_postfix(string,indicator,postfix)
|
||||
|
||||
character(len=*), intent(in) :: string
|
||||
character, intent(in) :: indicator
|
||||
character(len=*), intent(in) :: postfix
|
||||
character(len=:), allocatable :: IO_postfix
|
||||
|
||||
integer :: i,N
|
||||
|
||||
|
||||
IO_postfix = ''
|
||||
N = len(string)
|
||||
do i = 1, N
|
||||
IO_postfix = IO_postfix//string(i:i)
|
||||
if (string(i:i) == indicator .and. verify(IO_lc(string(min(i+1,N):min(i+1,N))),LOWER) == 0) &
|
||||
IO_postfix = IO_postfix//postfix
|
||||
end do
|
||||
|
||||
end function IO_postfix
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief Read string value at myChunk from string.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
@ -294,9 +321,6 @@ pure function IO_lc(str)
|
|||
character(len=*), intent(in) :: str !< string to convert
|
||||
character(len=len(str)) :: IO_lc
|
||||
|
||||
character(len=*), parameter :: LOWER = 'abcdefghijklmnopqrstuvwxyz'
|
||||
character(len=len(LOWER)), parameter :: UPPER = 'ABCDEFGHIJKLMNOPQRSTUVWXYZ'
|
||||
|
||||
integer :: i,n
|
||||
|
||||
|
||||
|
@ -832,6 +856,10 @@ subroutine selfTest()
|
|||
if ('abc,'//IO_EOL//'xxdefg,'//IO_EOL//'xxhij' /= IO_wrapLines('abc,defg, hij',filler='xx',length=4)) &
|
||||
error stop 'IO_wrapLines/7'
|
||||
|
||||
str='-a -1 -more 123 -flag -'
|
||||
out=IO_postfix(str,'-+','p_')
|
||||
if (out /= '-p_a -1 -p_more 123 -p_flag -') error stop 'IO_postfix'
|
||||
|
||||
end subroutine selfTest
|
||||
|
||||
end module IO
|
||||
|
|
|
@ -108,6 +108,7 @@ program DAMASK_grid
|
|||
quit
|
||||
type(tDict), pointer :: &
|
||||
config_load, &
|
||||
num_solver, &
|
||||
num_grid, &
|
||||
load_step, &
|
||||
solver, &
|
||||
|
@ -122,6 +123,7 @@ program DAMASK_grid
|
|||
character(len=:), allocatable :: &
|
||||
fileContent, fname
|
||||
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! init DAMASK (all modules)
|
||||
|
||||
|
@ -151,9 +153,12 @@ program DAMASK_grid
|
|||
end if
|
||||
|
||||
call parallelization_bcast_str(fileContent)
|
||||
config_load => YAML_parse_str_asDict(fileContent)
|
||||
config_load => YAML_parse_str_asDict(fileContent) !ToDo: misleading prefix (overlaps with entities from config module)
|
||||
solver => config_load%get_dict('solver')
|
||||
|
||||
num_solver => config_numerics%get_dict('solver',defaultVal=emptyDict)
|
||||
num_grid => num_solver%get_dict('grid',defaultVal=emptyDict)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! assign mechanics solver depending on selected type
|
||||
|
||||
|
@ -313,21 +318,21 @@ program DAMASK_grid
|
|||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! doing initialization depending on active solvers
|
||||
call spectral_Utilities_init()
|
||||
call spectral_utilities_init()
|
||||
|
||||
do field = 2, nActiveFields
|
||||
select case (ID(field))
|
||||
|
||||
case (FIELD_THERMAL_ID)
|
||||
call grid_thermal_spectral_init()
|
||||
call grid_thermal_spectral_init(num_grid)
|
||||
|
||||
case (FIELD_DAMAGE_ID)
|
||||
call grid_damage_spectral_init()
|
||||
call grid_damage_spectral_init(num_grid)
|
||||
|
||||
end select
|
||||
end do
|
||||
|
||||
call mechanical_init()
|
||||
call mechanical_init(num_grid)
|
||||
call config_numerics_deallocate()
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -70,7 +70,9 @@ contains
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief allocates all neccessary fields and fills them with data
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine grid_damage_spectral_init()
|
||||
subroutine grid_damage_spectral_init(num_grid)
|
||||
|
||||
type(tDict), pointer, intent(in) :: num_grid
|
||||
|
||||
PetscInt, dimension(0:worldsize-1) :: localK
|
||||
integer :: i, j, k, ce
|
||||
|
@ -82,10 +84,12 @@ subroutine grid_damage_spectral_init()
|
|||
integer(HID_T) :: fileHandle, groupHandle
|
||||
real(pREAL), dimension(1,product(cells(1:2))*cells3) :: tempN
|
||||
type(tDict), pointer :: &
|
||||
num_grid, &
|
||||
num_generic
|
||||
num_grid_damage
|
||||
character(len=pSTRLEN) :: &
|
||||
snes_type
|
||||
character(len=:), allocatable :: &
|
||||
petsc_options
|
||||
|
||||
|
||||
print'(/,1x,a)', '<<<+- grid_spectral_damage init -+>>>'
|
||||
|
||||
|
@ -96,25 +100,24 @@ subroutine grid_damage_spectral_init()
|
|||
|
||||
!-------------------------------------------------------------------------------------------------
|
||||
! read numerical parameters and do sanity checks
|
||||
num_grid => config_numerics%get_dict('grid',defaultVal=emptyDict)
|
||||
num%itmax = num_grid%get_asInt ('itmax',defaultVal=250)
|
||||
num%eps_damage_atol = num_grid%get_asReal ('eps_damage_atol',defaultVal=1.0e-2_pREAL)
|
||||
num%eps_damage_rtol = num_grid%get_asReal ('eps_damage_rtol',defaultVal=1.0e-6_pREAL)
|
||||
num_grid_damage => num_grid%get_dict('damage',defaultVal=emptyDict)
|
||||
|
||||
num_generic => config_numerics%get_dict('generic',defaultVal=emptyDict)
|
||||
num%phi_min = num_generic%get_asReal('phi_min', defaultVal=1.0e-6_pREAL)
|
||||
num%itmax = num_grid_damage%get_asInt('N_iter_max', defaultVal=250)
|
||||
|
||||
if (num%phi_min < 0.0_pREAL) call IO_error(301,ext_msg='phi_min')
|
||||
if (num%itmax <= 1) call IO_error(301,ext_msg='itmax')
|
||||
if (num%eps_damage_atol <= 0.0_pREAL) call IO_error(301,ext_msg='eps_damage_atol')
|
||||
if (num%eps_damage_rtol <= 0.0_pREAL) call IO_error(301,ext_msg='eps_damage_rtol')
|
||||
num%eps_damage_atol = num_grid_damage%get_asReal('eps_abs_phi',defaultVal=1.0e-2_pReal)
|
||||
num%eps_damage_rtol = num_grid_damage%get_asReal('eps_rel_phi',defaultVal=1.0e-6_pReal)
|
||||
num%phi_min = num_grid_damage%get_asReal('phi_min', defaultVal=1.0e-6_pReal)
|
||||
|
||||
if (num%phi_min < 0.0_pReal) call IO_error(301,ext_msg='phi_min')
|
||||
if (num%itmax <= 1) call IO_error(301,ext_msg='N_iter_max')
|
||||
if (num%eps_damage_atol <= 0.0_pReal) call IO_error(301,ext_msg='eps_abs_phi')
|
||||
if (num%eps_damage_rtol <= 0.0_pReal) call IO_error(301,ext_msg='eps_rel_phi')
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! set default and user defined options for PETSc
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,'-damage_snes_type newtonls -damage_snes_mf &
|
||||
&-damage_snes_ksp_ew -damage_ksp_type fgmres',err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,num_grid%get_asStr('petsc_options',defaultVal=''),err_PETSc)
|
||||
petsc_options = IO_postfix('-snes_type newtonls -snes_mf -snes_ksp_ew -ksp_type fgmres '// &
|
||||
num_grid_damage%get_asStr('PETSc_options',defaultVal=''), '-','damage_')
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,petsc_options,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -94,9 +94,11 @@ module grid_mechanical_FEM
|
|||
contains
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief allocates all necessary fields and fills them with data, potentially from restart info
|
||||
!> @brief Allocate all necessary fields and fills them with data, potentially from restart info.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine grid_mechanical_FEM_init
|
||||
subroutine grid_mechanical_FEM_init(num_grid)
|
||||
|
||||
type(tDict), pointer, intent(in) :: num_grid
|
||||
|
||||
real(pREAL), parameter :: HGCoeff = 0.0e-2_pREAL
|
||||
real(pREAL), parameter, dimension(4,8) :: &
|
||||
|
@ -118,41 +120,42 @@ subroutine grid_mechanical_FEM_init
|
|||
PetscInt, dimension(0:worldsize-1) :: localK
|
||||
integer(HID_T) :: fileHandle, groupHandle
|
||||
type(tDict), pointer :: &
|
||||
num_grid
|
||||
num_grid_mech
|
||||
character(len=pSTRLEN) :: &
|
||||
extmsg = ''
|
||||
character(len=:), allocatable :: &
|
||||
petsc_options
|
||||
|
||||
|
||||
print'(/,1x,a)', '<<<+- grid_mechanical_FEM init -+>>>'; flush(IO_STDOUT)
|
||||
|
||||
!-------------------------------------------------------------------------------------------------
|
||||
! read numerical parameters and do sanity checks
|
||||
num_grid => config_numerics%get_dict('grid',defaultVal=emptyDict)
|
||||
num_grid_mech => num_grid%get_dict('mechanical',defaultVal=emptyDict)
|
||||
|
||||
num%eps_div_atol = num_grid%get_asReal('eps_div_atol', defaultVal=1.0e-4_pREAL)
|
||||
num%eps_div_rtol = num_grid%get_asReal('eps_div_rtol', defaultVal=5.0e-4_pREAL)
|
||||
num%eps_stress_atol = num_grid%get_asReal('eps_stress_atol',defaultVal=1.0e3_pREAL)
|
||||
num%eps_stress_rtol = num_grid%get_asReal('eps_stress_rtol',defaultVal=1.0e-3_pREAL)
|
||||
num%itmin = num_grid%get_asInt ('itmin',defaultVal=1)
|
||||
num%itmax = num_grid%get_asInt ('itmax',defaultVal=250)
|
||||
num%eps_div_atol = num_grid_mech%get_asReal('eps_abs_div(P)',defaultVal=1.0e-4_pReal)
|
||||
num%eps_div_rtol = num_grid_mech%get_asReal('eps_rel_div(P)',defaultVal=5.0e-4_pReal)
|
||||
num%eps_stress_atol = num_grid_mech%get_asReal('eps_abs_P', defaultVal=1.0e3_pReal)
|
||||
num%eps_stress_rtol = num_grid_mech%get_asReal('eps_rel_P', defaultVal=1.0e-3_pReal)
|
||||
|
||||
if (num%eps_div_atol <= 0.0_pREAL) extmsg = trim(extmsg)//' eps_div_atol'
|
||||
if (num%eps_div_rtol < 0.0_pREAL) extmsg = trim(extmsg)//' eps_div_rtol'
|
||||
if (num%eps_stress_atol <= 0.0_pREAL) extmsg = trim(extmsg)//' eps_stress_atol'
|
||||
if (num%eps_stress_rtol < 0.0_pREAL) extmsg = trim(extmsg)//' eps_stress_rtol'
|
||||
if (num%itmax <= 1) extmsg = trim(extmsg)//' itmax'
|
||||
if (num%itmin > num%itmax .or. num%itmin < 1) extmsg = trim(extmsg)//' itmin'
|
||||
num%itmin = num_grid_mech%get_asInt('N_iter_min',defaultVal=1)
|
||||
num%itmax = num_grid_mech%get_asInt('N_iter_max',defaultVal=250)
|
||||
|
||||
if (num%eps_div_atol <= 0.0_pReal) extmsg = trim(extmsg)//' eps_abs_div(P)'
|
||||
if (num%eps_div_rtol < 0.0_pReal) extmsg = trim(extmsg)//' eps_rel_div(P)'
|
||||
if (num%eps_stress_atol <= 0.0_pReal) extmsg = trim(extmsg)//' eps_abs_P'
|
||||
if (num%eps_stress_rtol < 0.0_pReal) extmsg = trim(extmsg)//' eps_rel_P'
|
||||
if (num%itmax <= 1) extmsg = trim(extmsg)//' N_iter_max'
|
||||
if (num%itmin > num%itmax .or. num%itmin < 1) extmsg = trim(extmsg)//' N_iter_min'
|
||||
|
||||
if (extmsg /= '') call IO_error(301,ext_msg=trim(extmsg))
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! set default and user defined options for PETSc
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS, &
|
||||
'-mechanical_snes_type newtonls -mechanical_ksp_type fgmres &
|
||||
&-mechanical_ksp_max_it 25', &
|
||||
err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,num_grid%get_asStr('petsc_options',defaultVal=''),err_PETSc)
|
||||
|
||||
petsc_options = IO_postfix('-snes_type newtonls -ksp_type fgmres -ksp_max_it 25 '// &
|
||||
num_grid_mech%get_asStr('PETSc_options',defaultVal=''), '-','mechanical_')
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,petsc_options,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -22,7 +22,6 @@ module grid_mechanical_spectral_basic
|
|||
use math
|
||||
use rotations
|
||||
use spectral_utilities
|
||||
use config
|
||||
use homogenization
|
||||
use discretization_grid
|
||||
|
||||
|
@ -103,7 +102,9 @@ contains
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief allocates all necessary fields and fills them with data, potentially from restart info
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine grid_mechanical_spectral_basic_init()
|
||||
subroutine grid_mechanical_spectral_basic_init(num_grid)
|
||||
|
||||
type(tDict), pointer, intent(in) :: num_grid
|
||||
|
||||
real(pREAL), dimension(3,3,cells(1),cells(2),cells3) :: P
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
@ -114,9 +115,12 @@ subroutine grid_mechanical_spectral_basic_init()
|
|||
real(pREAL), dimension(3,3,product(cells(1:2))*cells3) :: temp33n
|
||||
integer(HID_T) :: fileHandle, groupHandle
|
||||
type(tDict), pointer :: &
|
||||
num_grid
|
||||
num_grid_fft, &
|
||||
num_grid_mech
|
||||
character(len=pSTRLEN) :: &
|
||||
extmsg = ''
|
||||
character(len=:), allocatable :: &
|
||||
petsc_options
|
||||
|
||||
|
||||
print'(/,1x,a)', '<<<+- grid_mechanical_spectral_basic init -+>>>'; flush(IO_STDOUT)
|
||||
|
@ -129,30 +133,32 @@ subroutine grid_mechanical_spectral_basic_init()
|
|||
|
||||
!-------------------------------------------------------------------------------------------------
|
||||
! read numerical parameters and do sanity checks
|
||||
num_grid => config_numerics%get_dict('grid',defaultVal=emptyDict)
|
||||
num_grid_fft => num_grid%get_dict('FFT',defaultVal=emptyDict)
|
||||
num_grid_mech => num_grid%get_dict('mechanical',defaultVal=emptyDict)
|
||||
|
||||
num%update_gamma = num_grid%get_asBool('update_gamma', defaultVal=.false.)
|
||||
num%eps_div_atol = num_grid%get_asReal('eps_div_atol', defaultVal=1.0e-4_pREAL)
|
||||
num%eps_div_rtol = num_grid%get_asReal('eps_div_rtol', defaultVal=5.0e-4_pREAL)
|
||||
num%eps_stress_atol = num_grid%get_asReal('eps_stress_atol',defaultVal=1.0e3_pREAL)
|
||||
num%eps_stress_rtol = num_grid%get_asReal('eps_stress_rtol',defaultVal=1.0e-3_pREAL)
|
||||
num%itmin = num_grid%get_asInt ('itmin',defaultVal=1)
|
||||
num%itmax = num_grid%get_asInt ('itmax',defaultVal=250)
|
||||
num%itmin = num_grid_mech%get_asInt ('N_iter_min',defaultVal=1)
|
||||
num%itmax = num_grid_mech%get_asInt ('N_iter_max',defaultVal=250)
|
||||
|
||||
if (num%eps_div_atol <= 0.0_pREAL) extmsg = trim(extmsg)//' eps_div_atol'
|
||||
if (num%eps_div_rtol < 0.0_pREAL) extmsg = trim(extmsg)//' eps_div_rtol'
|
||||
if (num%eps_stress_atol <= 0.0_pREAL) extmsg = trim(extmsg)//' eps_stress_atol'
|
||||
if (num%eps_stress_rtol < 0.0_pREAL) extmsg = trim(extmsg)//' eps_stress_rtol'
|
||||
if (num%itmax <= 1) extmsg = trim(extmsg)//' itmax'
|
||||
if (num%itmin > num%itmax .or. num%itmin < 1) extmsg = trim(extmsg)//' itmin'
|
||||
num%update_gamma = num_grid_fft%get_asBool ('update_gamma',defaultVal=.false.)
|
||||
|
||||
num%eps_div_atol = num_grid_mech%get_asReal('eps_abs_div(P)', defaultVal=1.0e-4_pReal)
|
||||
num%eps_div_rtol = num_grid_mech%get_asReal('eps_rel_div(P)', defaultVal=5.0e-4_pReal)
|
||||
num%eps_stress_atol = num_grid_mech%get_asReal('eps_abs_P', defaultVal=1.0e3_pReal)
|
||||
num%eps_stress_rtol = num_grid_mech%get_asReal('eps_rel_P', defaultVal=1.0e-3_pReal)
|
||||
|
||||
if (num%eps_div_atol <= 0.0_pReal) extmsg = trim(extmsg)//' eps_abs_div(P)'
|
||||
if (num%eps_div_rtol < 0.0_pReal) extmsg = trim(extmsg)//' eps_rel_div(P)'
|
||||
if (num%eps_stress_atol <= 0.0_pReal) extmsg = trim(extmsg)//' eps_abs_P'
|
||||
if (num%eps_stress_rtol < 0.0_pReal) extmsg = trim(extmsg)//' eps_rel_P'
|
||||
if (num%itmax <= 1) extmsg = trim(extmsg)//' N_iter_max'
|
||||
if (num%itmin > num%itmax .or. num%itmin < 1) extmsg = trim(extmsg)//' N_iter_min'
|
||||
|
||||
if (extmsg /= '') call IO_error(301,ext_msg=trim(extmsg))
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! set default and user defined options for PETSc
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,'-mechanical_snes_type ngmres',err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,num_grid%get_asStr('petsc_options',defaultVal=''),err_PETSc)
|
||||
petsc_options = IO_postfix('-snes_type ngmres '//num_grid_mech%get_asStr('PETSc_options',defaultVal=''), '-','mechanical_')
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,petsc_options,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -112,9 +112,11 @@ module grid_mechanical_spectral_polarisation
|
|||
contains
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief allocates all necessary fields and fills them with data, potentially from restart info
|
||||
!> @brief Allocate all necessary fields and fills them with data, potentially from restart info.
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine grid_mechanical_spectral_polarisation_init()
|
||||
subroutine grid_mechanical_spectral_polarisation_init(num_grid)
|
||||
|
||||
type(tDict), pointer, intent(in) :: num_grid
|
||||
|
||||
real(pREAL), dimension(3,3,cells(1),cells(2),cells3) :: P
|
||||
PetscErrorCode :: err_PETSc
|
||||
|
@ -127,9 +129,12 @@ subroutine grid_mechanical_spectral_polarisation_init()
|
|||
real(pREAL), dimension(3,3,product(cells(1:2))*cells3) :: temp33n
|
||||
integer(HID_T) :: fileHandle, groupHandle
|
||||
type(tDict), pointer :: &
|
||||
num_grid
|
||||
num_grid_fft,&
|
||||
num_grid_mech
|
||||
character(len=pSTRLEN) :: &
|
||||
extmsg = ''
|
||||
character(len=:), allocatable :: &
|
||||
petsc_options
|
||||
|
||||
|
||||
print '(/,1x,a)', '<<<+- grid_mechanical_spectral_polarization init -+>>>'; flush(IO_STDOUT)
|
||||
|
@ -137,41 +142,42 @@ subroutine grid_mechanical_spectral_polarisation_init()
|
|||
print '(/,1x,a)', 'P. Shanthraj et al., International Journal of Plasticity 66:31–45, 2015'
|
||||
print '( 1x,a)', 'https://doi.org/10.1016/j.ijplas.2014.02.006'
|
||||
|
||||
|
||||
!-------------------------------------------------------------------------------------------------
|
||||
! read numerical parameters and do sanity checks
|
||||
num_grid => config_numerics%get_dict('grid',defaultVal=emptyDict)
|
||||
num_grid_fft => num_grid%get_dict('FFT',defaultVal=emptyDict)
|
||||
num_grid_mech => num_grid%get_dict('mechanical',defaultVal=emptyDict)
|
||||
|
||||
num%update_gamma = num_grid%get_asBool('update_gamma', defaultVal=.false.)
|
||||
num%eps_div_atol = num_grid%get_asReal('eps_div_atol', defaultVal=1.0e-4_pREAL)
|
||||
num%eps_div_rtol = num_grid%get_asReal('eps_div_rtol', defaultVal=5.0e-4_pREAL)
|
||||
num%eps_curl_atol = num_grid%get_asReal('eps_curl_atol', defaultVal=1.0e-10_pREAL)
|
||||
num%eps_curl_rtol = num_grid%get_asReal('eps_curl_rtol', defaultVal=5.0e-4_pREAL)
|
||||
num%eps_stress_atol = num_grid%get_asReal('eps_stress_atol',defaultVal=1.0e3_pREAL)
|
||||
num%eps_stress_rtol = num_grid%get_asReal('eps_stress_rtol',defaultVal=1.0e-3_pREAL)
|
||||
num%itmin = num_grid%get_asInt ('itmin', defaultVal=1)
|
||||
num%itmax = num_grid%get_asInt ('itmax', defaultVal=250)
|
||||
num%alpha = num_grid%get_asReal('alpha', defaultVal=1.0_pREAL)
|
||||
num%beta = num_grid%get_asReal('beta', defaultVal=1.0_pREAL)
|
||||
num%update_gamma = num_grid_fft%get_asBool ('update_gamma',defaultVal=.false.)
|
||||
|
||||
if (num%eps_div_atol <= 0.0_pREAL) extmsg = trim(extmsg)//' eps_div_atol'
|
||||
if (num%eps_div_rtol < 0.0_pREAL) extmsg = trim(extmsg)//' eps_div_rtol'
|
||||
if (num%eps_curl_atol <= 0.0_pREAL) extmsg = trim(extmsg)//' eps_curl_atol'
|
||||
if (num%eps_curl_rtol < 0.0_pREAL) extmsg = trim(extmsg)//' eps_curl_rtol'
|
||||
if (num%eps_stress_atol <= 0.0_pREAL) extmsg = trim(extmsg)//' eps_stress_atol'
|
||||
if (num%eps_stress_rtol < 0.0_pREAL) extmsg = trim(extmsg)//' eps_stress_rtol'
|
||||
if (num%itmax <= 1) extmsg = trim(extmsg)//' itmax'
|
||||
if (num%itmin > num%itmax .or. num%itmin < 1) extmsg = trim(extmsg)//' itmin'
|
||||
if (num%alpha <= 0.0_pREAL .or. num%alpha > 2.0_pREAL) extmsg = trim(extmsg)//' alpha'
|
||||
if (num%beta < 0.0_pREAL .or. num%beta > 2.0_pREAL) extmsg = trim(extmsg)//' beta'
|
||||
num%eps_div_atol = num_grid_mech%get_asReal('eps_abs_div(P)', defaultVal=1.0e-4_pReal)
|
||||
num%eps_div_rtol = num_grid_mech%get_asReal('eps_rel_div(P)', defaultVal=5.0e-4_pReal)
|
||||
num%eps_curl_atol = num_grid_mech%get_asReal('eps_abs_curl(F)',defaultVal=1.0e-10_pReal)
|
||||
num%eps_curl_rtol = num_grid_mech%get_asReal('eps_rel_curl(F)',defaultVal=5.0e-4_pReal)
|
||||
num%eps_stress_atol = num_grid_mech%get_asReal('eps_abs_P', defaultVal=1.0e3_pReal)
|
||||
num%eps_stress_rtol = num_grid_mech%get_asReal('eps_rel_P', defaultVal=1.0e-3_pReal)
|
||||
num%alpha = num_grid_mech%get_asReal('alpha', defaultVal=1.0_pReal)
|
||||
num%beta = num_grid_mech%get_asReal('beta', defaultVal=1.0_pReal)
|
||||
|
||||
num%itmin = num_grid_mech%get_asInt ('N_iter_min',defaultVal=1)
|
||||
num%itmax = num_grid_mech%get_asInt ('N_iter_max',defaultVal=250)
|
||||
|
||||
if (num%eps_div_atol <= 0.0_pReal) extmsg = trim(extmsg)//' eps_abs_div(P)'
|
||||
if (num%eps_div_rtol < 0.0_pReal) extmsg = trim(extmsg)//' eps_rel_div(P)'
|
||||
if (num%eps_curl_atol <= 0.0_pReal) extmsg = trim(extmsg)//' eps_abs_curl(F)'
|
||||
if (num%eps_curl_rtol < 0.0_pReal) extmsg = trim(extmsg)//' eps_rel_curl(F)'
|
||||
if (num%eps_stress_atol <= 0.0_pReal) extmsg = trim(extmsg)//' eps_abs_P'
|
||||
if (num%eps_stress_rtol < 0.0_pReal) extmsg = trim(extmsg)//' eps_rel_P'
|
||||
if (num%itmax <= 1) extmsg = trim(extmsg)//' N_iter_max'
|
||||
if (num%itmin > num%itmax .or. num%itmin < 1) extmsg = trim(extmsg)//' N_iter_min'
|
||||
if (num%alpha <= 0.0_pReal .or. num%alpha > 2.0_pReal) extmsg = trim(extmsg)//' alpha'
|
||||
if (num%beta < 0.0_pReal .or. num%beta > 2.0_pReal) extmsg = trim(extmsg)//' beta'
|
||||
|
||||
if (extmsg /= '') call IO_error(301,ext_msg=trim(extmsg))
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! set default and user defined options for PETSc
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,'-mechanical_snes_type ngmres',err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,num_grid%get_asStr('petsc_options',defaultVal=''),err_PETSc)
|
||||
petsc_options = IO_postfix('-snes_type ngmres '//num_grid_mech%get_asStr('PETSc_options',defaultVal=''), '-','mechanical_')
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,petsc_options,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -69,7 +69,9 @@ contains
|
|||
!--------------------------------------------------------------------------------------------------
|
||||
!> @brief allocates all neccessary fields and fills them with data
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
subroutine grid_thermal_spectral_init()
|
||||
subroutine grid_thermal_spectral_init(num_grid)
|
||||
|
||||
type(tDict), pointer, intent(in) :: num_grid
|
||||
|
||||
PetscInt, dimension(0:worldsize-1) :: localK
|
||||
integer :: i, j, k, ce
|
||||
|
@ -80,7 +82,10 @@ subroutine grid_thermal_spectral_init()
|
|||
integer(HID_T) :: fileHandle, groupHandle
|
||||
real(pREAL), dimension(1,product(cells(1:2))*cells3) :: tempN
|
||||
type(tDict), pointer :: &
|
||||
num_grid
|
||||
num_grid_thermal
|
||||
character(len=:), allocatable :: &
|
||||
petsc_options
|
||||
|
||||
|
||||
print'(/,1x,a)', '<<<+- grid_thermal_spectral init -+>>>'
|
||||
|
||||
|
@ -91,21 +96,22 @@ subroutine grid_thermal_spectral_init()
|
|||
|
||||
!-------------------------------------------------------------------------------------------------
|
||||
! read numerical parameters and do sanity checks
|
||||
num_grid => config_numerics%get_dict('grid',defaultVal=emptyDict)
|
||||
num%itmax = num_grid%get_asInt ('itmax', defaultVal=250)
|
||||
num%eps_thermal_atol = num_grid%get_asReal('eps_thermal_atol',defaultVal=1.0e-2_pREAL)
|
||||
num%eps_thermal_rtol = num_grid%get_asReal('eps_thermal_rtol',defaultVal=1.0e-6_pREAL)
|
||||
num_grid_thermal => num_grid%get_dict('thermal',defaultVal=emptyDict)
|
||||
|
||||
if (num%itmax <= 1) call IO_error(301,ext_msg='itmax')
|
||||
if (num%eps_thermal_atol <= 0.0_pREAL) call IO_error(301,ext_msg='eps_thermal_atol')
|
||||
if (num%eps_thermal_rtol <= 0.0_pREAL) call IO_error(301,ext_msg='eps_thermal_rtol')
|
||||
num%itmax = num_grid_thermal%get_asInt('N_iter_max', defaultVal=250)
|
||||
|
||||
num%eps_thermal_atol = num_grid_thermal%get_asReal('eps_abs_T', defaultVal=1.0e-2_pReal)
|
||||
num%eps_thermal_rtol = num_grid_thermal%get_asReal('eps_rel_T', defaultVal=1.0e-6_pReal)
|
||||
|
||||
if (num%itmax <= 1) call IO_error(301,ext_msg='N_iter_max')
|
||||
if (num%eps_thermal_atol <= 0.0_pReal) call IO_error(301,ext_msg='eps_abs_T')
|
||||
if (num%eps_thermal_rtol <= 0.0_pReal) call IO_error(301,ext_msg='eps_rel_T')
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! set default and user defined options for PETSc
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,'-thermal_snes_type newtonls -thermal_snes_mf &
|
||||
&-thermal_snes_ksp_ew -thermal_ksp_type fgmres',err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,num_grid%get_asStr('petsc_options',defaultVal=''),err_PETSc)
|
||||
petsc_options = IO_postfix('-snes_type newtonls -snes_mf -snes_ksp_ew -ksp_type fgmres '// &
|
||||
num_grid_thermal%get_asStr('PETSc_options',defaultVal=''), '-','thermal_')
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,petsc_options,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
|
|
@ -146,8 +146,9 @@ subroutine spectral_utilities_init()
|
|||
vectorSize = 3_C_INTPTR_T, &
|
||||
tensorSize = 9_C_INTPTR_T
|
||||
type(tDict) , pointer :: &
|
||||
num_grid
|
||||
|
||||
num_solver, &
|
||||
num_grid, &
|
||||
num_grid_fft
|
||||
|
||||
print'(/,1x,a)', '<<<+- spectral_utilities init -+>>>'
|
||||
|
||||
|
@ -163,8 +164,10 @@ subroutine spectral_utilities_init()
|
|||
print'( 1x,a)', 'P. Shanthraj et al., Handbook of Mechanics of Materials, 2019'
|
||||
print'( 1x,a)', 'https://doi.org/10.1007/978-981-10-6855-3_80'
|
||||
|
||||
num_solver => config_numerics%get_dict('solver',defaultVal=emptyDict)
|
||||
num_grid => num_solver%get_dict('grid',defaultVal=emptyDict)
|
||||
num_grid_fft => num_grid%get_dict('FFT',defaultVal=emptyDict)
|
||||
|
||||
num_grid => config_numerics%get_dict('grid',defaultVal=emptyDict)
|
||||
call PetscOptionsClear(PETSC_NULL_OPTIONS,err_PETSc)
|
||||
CHKERRQ(err_PETSc)
|
||||
call PetscOptionsInsertString(PETSC_NULL_OPTIONS,&
|
||||
|
@ -174,13 +177,13 @@ subroutine spectral_utilities_init()
|
|||
cells1Red = cells(1)/2 + 1
|
||||
wgt = real(product(cells),pREAL)**(-1)
|
||||
|
||||
num%memory_efficient = num_grid%get_asInt('memory_efficient', defaultVal=1) > 0 ! ToDo: should be logical in YAML file
|
||||
num%divergence_correction = num_grid%get_asInt('divergence_correction', defaultVal=2)
|
||||
num%memory_efficient = num_grid_fft%get_asBool('memory_efficient', defaultVal=.true.)
|
||||
num%divergence_correction = num_grid_fft%get_asInt('divergence_correction', defaultVal=2)
|
||||
|
||||
if (num%divergence_correction < 0 .or. num%divergence_correction > 2) &
|
||||
call IO_error(301,ext_msg='divergence_correction')
|
||||
|
||||
select case (num_grid%get_asStr('derivative',defaultVal='continuous'))
|
||||
select case (num_grid_fft%get_asStr('derivative',defaultVal='continuous'))
|
||||
case ('continuous')
|
||||
spectral_derivative_ID = DERIVATIVE_CONTINUOUS_ID
|
||||
case ('central_difference')
|
||||
|
@ -188,7 +191,7 @@ subroutine spectral_utilities_init()
|
|||
case ('FWBW_difference')
|
||||
spectral_derivative_ID = DERIVATIVE_FWBW_DIFF_ID
|
||||
case default
|
||||
call IO_error(892,ext_msg=trim(num_grid%get_asStr('derivative')))
|
||||
call IO_error(892,ext_msg=trim(num_grid_fft%get_asStr('derivative')))
|
||||
end select
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
|
@ -209,7 +212,7 @@ subroutine spectral_utilities_init()
|
|||
scaledGeomSize = geomSize
|
||||
end if
|
||||
|
||||
select case(IO_lc(num_grid%get_asStr('fftw_plan_mode',defaultVal='FFTW_MEASURE')))
|
||||
select case(IO_lc(num_grid_fft%get_asStr('plan_mode',defaultVal='FFTW_MEASURE')))
|
||||
case('fftw_estimate') ! ordered from slow execution (but fast plan creation) to fast execution
|
||||
FFTW_planner_flag = FFTW_ESTIMATE
|
||||
case('fftw_measure')
|
||||
|
@ -219,14 +222,14 @@ subroutine spectral_utilities_init()
|
|||
case('fftw_exhaustive')
|
||||
FFTW_planner_flag = FFTW_EXHAUSTIVE
|
||||
case default
|
||||
call IO_warning(47,'using default FFTW_MEASURE instead of "'//trim(num_grid%get_asStr('fftw_plan_mode'))//'"')
|
||||
call IO_warning(47,'using default FFTW_MEASURE instead of "'//trim(num_grid_fft%get_asStr('plan_mode'))//'"')
|
||||
FFTW_planner_flag = FFTW_MEASURE
|
||||
end select
|
||||
|
||||
!--------------------------------------------------------------------------------------------------
|
||||
! general initialization of FFTW (see manual on fftw.org for more details)
|
||||
if (pREAL /= C_DOUBLE .or. kind(1) /= C_INT) error stop 'C and Fortran datatypes do not match'
|
||||
call fftw_set_timelimit(num_grid%get_asReal('fftw_timelimit',defaultVal=300.0_pREAL))
|
||||
call fftw_set_timelimit(num_grid_fft%get_asReal('fftw_timelimit',defaultVal=300.0_pREAL))
|
||||
|
||||
print'(/,1x,a)', 'FFTW initialized'; flush(IO_STDOUT)
|
||||
|
||||
|
|
Loading…
Reference in New Issue