2011-04-07 12:50:28 +05:30
! Copyright 2011 Max-Planck-Institut für Eisenforschung GmbH
2011-04-04 19:39:54 +05:30
!
! This file is part of DAMASK,
2011-04-07 12:50:28 +05:30
! the Düsseldorf Advanced MAterial Simulation Kit.
2011-04-04 19:39:54 +05:30
!
! DAMASK is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! DAMASK is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with DAMASK. If not, see <http://www.gnu.org/licenses/>.
!
!##############################################################
2009-08-31 20:39:15 +05:30
!* $Id$
2009-06-15 18:41:21 +05:30
!##############################################################
MODULE numerics
!##############################################################
use prec , only : pInt , pReal
2011-12-01 17:31:13 +05:30
use IO , only : IO_warning
2009-06-15 18:41:21 +05:30
implicit none
2012-01-31 01:46:19 +05:30
character ( len = 64 ) , parameter :: numerics_configFile = 'numerics.config' ! name of configuration file
integer ( pInt ) :: iJacoStiffness = 1_pInt , & ! frequency of stiffness update
iJacoLpresiduum = 1_pInt , & ! frequency of Jacobian update of residuum in Lp
nHomog = 20_pInt , & ! homogenization loop limit (only for debugging info, loop limit is determined by "subStepMinHomog")
nMPstate = 10_pInt , & ! materialpoint state loop limit
nCryst = 20_pInt , & ! crystallite loop limit (only for debugging info, loop limit is determined by "subStepMinCryst")
nState = 10_pInt , & ! state loop limit
nStress = 40_pInt , & ! stress loop limit
pert_method = 1_pInt , & ! method used in perturbation technique for tangent
numerics_integrationMode = 0_pInt ! integrationMode 1 = central solution ; integrationMode 2 = perturbation, Default 0: undefined, is not read from file
integer ( pInt ) , dimension ( 2 ) :: numerics_integrator = 1_pInt ! method used for state integration (central & perturbed state), Default 1: fix-point iteration for both states
real ( pReal ) :: relevantStrain = 1.0e-7_pReal , & ! strain increment considered significant (used by crystallite to determine whether strain inc is considered significant)
defgradTolerance = 1.0e-7_pReal , & ! deviation of deformation gradient that is still allowed (used by CPFEM to determine outdated ffn1)
pert_Fg = 1.0e-7_pReal , & ! strain perturbation for FEM Jacobi
subStepMinCryst = 1.0e-3_pReal , & ! minimum (relative) size of sub-step allowed during cutback in crystallite
subStepMinHomog = 1.0e-3_pReal , & ! minimum (relative) size of sub-step allowed during cutback in homogenization
subStepSizeCryst = 0.25_pReal , & ! size of first substep when cutback in crystallite
subStepSizeHomog = 0.25_pReal , & ! size of first substep when cutback in homogenization
stepIncreaseCryst = 1.5_pReal , & ! increase of next substep size when previous substep converged in crystallite
stepIncreaseHomog = 1.5_pReal , & ! increase of next substep size when previous substep converged in homogenization
rTol_crystalliteState = 1.0e-6_pReal , & ! relative tolerance in crystallite state loop
rTol_crystalliteTemperature = 1.0e-6_pReal , & ! relative tolerance in crystallite temperature loop
rTol_crystalliteStress = 1.0e-6_pReal , & ! relative tolerance in crystallite stress loop
aTol_crystalliteStress = 1.0e-8_pReal , & ! absolute tolerance in crystallite stress loop, Default 1.0e-8: residuum is in Lp and hence strain is on this order
2009-07-31 17:32:20 +05:30
2012-01-31 01:46:19 +05:30
absTol_RGC = 1.0e+4_pReal , & ! absolute tolerance of RGC residuum
relTol_RGC = 1.0e-3_pReal , & ! relative tolerance of RGC residuum
absMax_RGC = 1.0e+10_pReal , & ! absolute maximum of RGC residuum
relMax_RGC = 1.0e+2_pReal , & ! relative maximum of RGC residuum
pPert_RGC = 1.0e-7_pReal , & ! perturbation for computing RGC penalty tangent
xSmoo_RGC = 1.0e-5_pReal , & ! RGC penalty smoothing parameter (hyperbolic tangent)
viscPower_RGC = 1.0e+0_pReal , & ! power (sensitivity rate) of numerical viscosity in RGC scheme, Default 1.0e0: Newton viscosity (linear model)
viscModus_RGC = 0.0e+0_pReal , & ! stress modulus of RGC numerical viscosity, Default 0.0e0: No viscosity is applied
refRelaxRate_RGC = 1.0e-3_pReal , & ! reference relaxation rate in RGC viscosity
maxdRelax_RGC = 1.0e+0_pReal , & ! threshold of maximum relaxation vector increment (if exceed this then cutback)
maxVolDiscr_RGC = 1.0e-5_pReal , & ! threshold of maximum volume discrepancy allowed
volDiscrMod_RGC = 1.0e+12_pReal , & ! stiffness of RGC volume discrepancy (zero = without volume discrepancy constraint)
volDiscrPow_RGC = 5.0_pReal , & ! powerlaw penalty for volume discrepancy
2010-10-13 21:34:44 +05:30
!* spectral parameters:
2012-01-31 01:46:19 +05:30
err_div_tol = 1.0e-4_pReal , & ! error of divergence in fourier space, Default 1.0e-4: proposed by Suquet
err_stress_tolrel = 0.01_pReal , & ! relative tolerance for fullfillment of stress BC, Default: 0.01 allowing deviation of 1% of maximum stress
fftw_timelimit = - 1.0_pReal , & ! sets the timelimit of plan creation for FFTW, see manual on www.fftw.org, Default -1.0: disable timelimit
rotation_tol = 1.0e-12_pReal ! tolerance of rotation specified in loadcase, Default 1.0e-12: first guess
character ( len = 64 ) :: fftw_planner_string = 'FFTW_PATIENT' ! reads the planing-rigor flag, see manual on www.fftw.org, Default FFTW_PATIENT: use patiant planner flag
2012-02-02 01:50:05 +05:30
integer ( pInt ) :: fftw_planner_flag = - 1_pInt ! conversion of fftw_planner_string to integer, basically what is usually done in the include file of fftw
2012-01-31 01:46:19 +05:30
logical :: memory_efficient = . true . , & ! for fast execution (pre calculation of gamma_hat), Default .true.: do not precalculate
divergence_correction = . false . , & ! correct divergence calculation in fourier space, Default .false.: no correction
update_gamma = . false . , & ! update gamma operator with current stiffness, Default .false.: use initial stiffness
simplified_algorithm = . true . ! use short algorithm without fluctuation field, Default .true.: use simplified algorithm
real ( pReal ) :: cut_off_value = 0.0_pReal ! percentage of frequencies to cut away, Default 0.0: use all frequencies
integer ( pInt ) :: itmax = 20_pInt , & ! maximum number of iterations
2009-06-15 18:41:21 +05:30
2011-02-07 20:05:42 +05:30
2011-01-31 22:37:42 +05:30
!* Random seeding parameters
2012-01-31 01:46:19 +05:30
fixedSeed = 0_pInt ! fixed seeding for pseudo-random number generator, Default 0: use random seed
2011-01-31 22:37:42 +05:30
!* OpenMP variable
2012-01-31 01:46:19 +05:30
integer ( pInt ) :: DAMASK_NumThreadsInt = 0_pInt ! value stored in environment variable DAMASK_NUM_THREADS, set to zero if no OpenMP directive
2010-12-02 16:34:29 +05:30
2009-11-10 19:06:27 +05:30
2009-06-15 18:41:21 +05:30
CONTAINS
!*******************************************
! initialization subroutine
!*******************************************
subroutine numerics_init ( )
2012-02-13 19:38:07 +05:30
use , intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran 4.6 at the moment)
2009-06-15 18:41:21 +05:30
!*** variables and functions from other modules ***!
use prec , only : pInt , &
pReal
use IO , only : IO_error , &
2012-02-13 23:11:27 +05:30
IO_open_file_stat , &
2009-06-15 18:41:21 +05:30
IO_isBlank , &
IO_stringPos , &
IO_stringValue , &
IO_lc , &
IO_floatValue , &
IO_intValue
2010-12-02 16:34:29 +05:30
!$ use OMP_LIB ! the openMP function library
2009-06-15 18:41:21 +05:30
implicit none
!*** local variables ***!
2012-01-30 19:22:41 +05:30
integer ( pInt ) , parameter :: fileunit = 300_pInt
integer ( pInt ) , parameter :: maxNchunks = 2_pInt
integer ( pInt ) :: gotDAMASK_NUM_THREADS = 1_pInt
2009-06-15 18:41:21 +05:30
integer ( pInt ) , dimension ( 1 + 2 * maxNchunks ) :: positions
2012-01-31 01:46:19 +05:30
character ( len = 64 ) :: tag
character ( len = 1024 ) :: line
2009-06-15 18:41:21 +05:30
2010-12-02 16:34:29 +05:30
! OpenMP variable
2012-01-31 01:46:19 +05:30
!$ character(len=6) DAMASK_NumThreadsString !environment variable DAMASK_NUM_THREADS
2009-08-27 21:00:40 +05:30
2012-01-31 01:46:19 +05:30
!$OMP CRITICAL (write2out)
write ( 6 , * )
write ( 6 , * ) '<<<+- numerics init -+>>>'
write ( 6 , * ) '$Id$'
2012-02-01 00:48:55 +05:30
#include "compilation_info.f90"
2012-01-31 01:46:19 +05:30
!$OMP END CRITICAL (write2out)
2010-12-02 16:34:29 +05:30
2012-01-31 01:46:19 +05:30
!$ call GET_ENVIRONMENT_VARIABLE(NAME='DAMASK_NUM_THREADS',VALUE=DAMASK_NumThreadsString,STATUS=gotDAMASK_NUM_THREADS) ! get environment variable DAMASK_NUM_THREADS...
2012-02-10 16:54:53 +05:30
!$ if(gotDAMASK_NUM_THREADS /= 0_pInt) call IO_warning(47_pInt,ext_msg=DAMASK_NumThreadsString)
2012-01-31 01:46:19 +05:30
!$ read(DAMASK_NumThreadsString,'(i6)') DAMASK_NumThreadsInt ! ...convert it to integer...
2012-01-30 19:22:41 +05:30
!$ if (DAMASK_NumThreadsInt < 1) DAMASK_NumThreadsInt = 1 ! ...ensure that its at least one...
!$ call omp_set_num_threads(DAMASK_NumThreadsInt) ! ...and use it as number of threads for parallel execution
2010-12-02 16:34:29 +05:30
2009-06-18 19:58:02 +05:30
! try to open the config file
2012-02-13 23:11:27 +05:30
if ( IO_open_file_stat ( fileunit , numerics_configFile ) ) then
2009-06-18 19:58:02 +05:30
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
!$OMP CRITICAL (write2out)
write ( 6 , * ) ' ... using values from config file'
write ( 6 , * )
!$OMP END CRITICAL (write2out)
2009-06-18 19:58:02 +05:30
2012-02-13 23:11:27 +05:30
!* read variables from config file and overwrite parameters
2009-06-18 19:58:02 +05:30
line = ''
do
read ( fileunit , '(a1024)' , END = 100 ) line
if ( IO_isBlank ( line ) ) cycle ! skip empty lines
positions = IO_stringPos ( line , maxNchunks )
2012-02-09 18:05:55 +05:30
tag = IO_lc ( IO_stringValue ( line , positions , 1_pInt ) ) ! extract key
2009-06-18 19:58:02 +05:30
select case ( tag )
case ( 'relevantstrain' )
2012-02-09 18:05:55 +05:30
relevantStrain = IO_floatValue ( line , positions , 2_pInt )
2010-05-20 20:25:11 +05:30
case ( 'defgradtolerance' )
2012-02-09 18:05:55 +05:30
defgradTolerance = IO_floatValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'ijacostiffness' )
2012-02-09 18:05:55 +05:30
iJacoStiffness = IO_intValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'ijacolpresiduum' )
2012-02-09 18:05:55 +05:30
iJacoLpresiduum = IO_intValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'pert_fg' )
2012-02-09 18:05:55 +05:30
pert_Fg = IO_floatValue ( line , positions , 2_pInt )
2009-11-10 19:06:27 +05:30
case ( 'pert_method' )
2012-02-09 18:05:55 +05:30
pert_method = IO_intValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'nhomog' )
2012-02-09 18:05:55 +05:30
nHomog = IO_intValue ( line , positions , 2_pInt )
2009-08-11 22:01:57 +05:30
case ( 'nmpstate' )
2012-02-09 18:05:55 +05:30
nMPstate = IO_intValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'ncryst' )
2012-02-09 18:05:55 +05:30
nCryst = IO_intValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'nstate' )
2012-02-09 18:05:55 +05:30
nState = IO_intValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'nstress' )
2012-02-09 18:05:55 +05:30
nStress = IO_intValue ( line , positions , 2_pInt )
2009-10-26 22:13:43 +05:30
case ( 'substepmincryst' )
2012-02-09 18:05:55 +05:30
subStepMinCryst = IO_floatValue ( line , positions , 2_pInt )
2009-11-10 19:06:27 +05:30
case ( 'substepsizecryst' )
2012-02-09 18:05:55 +05:30
subStepSizeCryst = IO_floatValue ( line , positions , 2_pInt )
2009-11-10 19:06:27 +05:30
case ( 'stepincreasecryst' )
2012-02-09 18:05:55 +05:30
stepIncreaseCryst = IO_floatValue ( line , positions , 2_pInt )
2009-10-26 22:13:43 +05:30
case ( 'substepminhomog' )
2012-02-09 18:05:55 +05:30
subStepMinHomog = IO_floatValue ( line , positions , 2_pInt )
2009-11-10 19:06:27 +05:30
case ( 'substepsizehomog' )
2012-02-09 18:05:55 +05:30
subStepSizeHomog = IO_floatValue ( line , positions , 2_pInt )
2009-11-10 19:06:27 +05:30
case ( 'stepincreasehomog' )
2012-02-09 18:05:55 +05:30
stepIncreaseHomog = IO_floatValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'rtol_crystallitestate' )
2012-02-09 18:05:55 +05:30
rTol_crystalliteState = IO_floatValue ( line , positions , 2_pInt )
2009-07-22 21:37:19 +05:30
case ( 'rtol_crystallitetemperature' )
2012-02-09 18:05:55 +05:30
rTol_crystalliteTemperature = IO_floatValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'rtol_crystallitestress' )
2012-02-09 18:05:55 +05:30
rTol_crystalliteStress = IO_floatValue ( line , positions , 2_pInt )
2009-06-18 19:58:02 +05:30
case ( 'atol_crystallitestress' )
2012-02-09 18:05:55 +05:30
aTol_crystalliteStress = IO_floatValue ( line , positions , 2_pInt )
2010-10-01 17:48:49 +05:30
case ( 'integrator' )
2012-02-09 18:05:55 +05:30
numerics_integrator ( 1 ) = IO_intValue ( line , positions , 2_pInt )
2010-10-01 17:48:49 +05:30
case ( 'integratorstiffness' )
2012-02-09 18:05:55 +05:30
numerics_integrator ( 2 ) = IO_intValue ( line , positions , 2_pInt )
2009-07-31 17:32:20 +05:30
2012-02-13 23:11:27 +05:30
!* RGC parameters:
2009-07-31 17:32:20 +05:30
case ( 'atol_rgc' )
2012-02-09 18:05:55 +05:30
absTol_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-07-31 17:32:20 +05:30
case ( 'rtol_rgc' )
2012-02-09 18:05:55 +05:30
relTol_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-07-31 17:32:20 +05:30
case ( 'amax_rgc' )
2012-02-09 18:05:55 +05:30
absMax_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-07-31 17:32:20 +05:30
case ( 'rmax_rgc' )
2012-02-09 18:05:55 +05:30
relMax_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-07-31 17:32:20 +05:30
case ( 'perturbpenalty_rgc' )
2012-02-09 18:05:55 +05:30
pPert_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-07-31 17:32:20 +05:30
case ( 'relevantmismatch_rgc' )
2012-02-09 18:05:55 +05:30
xSmoo_RGC = IO_floatValue ( line , positions , 2_pInt )
2010-03-24 18:50:12 +05:30
case ( 'viscositypower_rgc' )
2012-02-09 18:05:55 +05:30
viscPower_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-11-17 19:12:38 +05:30
case ( 'viscositymodulus_rgc' )
2012-02-09 18:05:55 +05:30
viscModus_RGC = IO_floatValue ( line , positions , 2_pInt )
2010-03-24 18:50:12 +05:30
case ( 'refrelaxationrate_rgc' )
2012-02-09 18:05:55 +05:30
refRelaxRate_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-11-17 19:12:38 +05:30
case ( 'maxrelaxation_rgc' )
2012-02-09 18:05:55 +05:30
maxdRelax_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-12-16 21:50:53 +05:30
case ( 'maxvoldiscrepancy_rgc' )
2012-02-09 18:05:55 +05:30
maxVolDiscr_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-12-16 21:50:53 +05:30
case ( 'voldiscrepancymod_rgc' )
2012-02-09 18:05:55 +05:30
volDiscrMod_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-12-16 21:50:53 +05:30
case ( 'discrepancypower_rgc' )
2012-02-09 18:05:55 +05:30
volDiscrPow_RGC = IO_floatValue ( line , positions , 2_pInt )
2009-07-31 17:32:20 +05:30
2012-02-13 23:11:27 +05:30
!* spectral parameters
2011-01-31 22:37:42 +05:30
case ( 'err_div_tol' )
2012-02-09 18:05:55 +05:30
err_div_tol = IO_floatValue ( line , positions , 2_pInt )
2011-01-31 22:37:42 +05:30
case ( 'err_stress_tolrel' )
2012-02-09 18:05:55 +05:30
err_stress_tolrel = IO_floatValue ( line , positions , 2_pInt )
2011-01-31 22:37:42 +05:30
case ( 'itmax' )
2012-02-09 18:05:55 +05:30
itmax = IO_intValue ( line , positions , 2_pInt )
2011-02-21 20:07:38 +05:30
case ( 'memory_efficient' )
2012-02-09 18:05:55 +05:30
memory_efficient = IO_intValue ( line , positions , 2_pInt ) > 0_pInt
2011-10-18 14:46:18 +05:30
case ( 'fftw_timelimit' )
2012-02-09 18:05:55 +05:30
fftw_timelimit = IO_floatValue ( line , positions , 2_pInt )
2011-12-01 17:31:13 +05:30
case ( 'fftw_planner_string' )
2012-02-09 18:05:55 +05:30
fftw_planner_string = IO_stringValue ( line , positions , 2_pInt )
2011-10-25 19:08:24 +05:30
case ( 'rotation_tol' )
2012-02-09 18:05:55 +05:30
rotation_tol = IO_floatValue ( line , positions , 2_pInt )
2011-11-15 23:24:18 +05:30
case ( 'divergence_correction' )
2012-02-09 18:05:55 +05:30
divergence_correction = IO_intValue ( line , positions , 2_pInt ) > 0_pInt
2012-01-25 14:26:46 +05:30
case ( 'update_gamma' )
2012-02-09 18:05:55 +05:30
update_gamma = IO_intValue ( line , positions , 2_pInt ) > 0_pInt
2012-01-25 14:26:46 +05:30
case ( 'simplified_algorithm' )
2012-02-09 18:05:55 +05:30
simplified_algorithm = IO_intValue ( line , positions , 2_pInt ) > 0_pInt
2012-01-25 14:26:46 +05:30
case ( 'cut_off_value' )
2012-02-09 18:05:55 +05:30
cut_off_value = IO_floatValue ( line , positions , 2_pInt )
2010-10-13 21:34:44 +05:30
2012-02-13 23:11:27 +05:30
!* Random seeding parameters
2009-08-27 21:00:40 +05:30
case ( 'fixed_seed' )
2012-02-09 18:05:55 +05:30
fixedSeed = IO_intValue ( line , positions , 2_pInt )
2012-02-13 23:11:27 +05:30
case default
call IO_error ( 300_pInt , ext_msg = tag )
2009-06-18 19:58:02 +05:30
endselect
enddo
100 close ( fileunit )
! no config file, so we use standard values
else
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
!$OMP CRITICAL (write2out)
write ( 6 , * ) ' ... using standard values'
write ( 6 , * )
!$OMP END CRITICAL (write2out)
2009-06-18 19:58:02 +05:30
endif
2012-02-13 23:11:27 +05:30
2011-12-01 17:31:13 +05:30
select case ( IO_lc ( fftw_planner_string ) ) ! setting parameters for the plan creation of FFTW. Basically a translation from fftw3.f
case ( 'estimate' , 'fftw_estimate' ) ! ordered from slow execution (but fast plan creation) to fast execution
2012-02-02 01:50:05 +05:30
fftw_planner_flag = 64_pInt
2011-12-01 17:31:13 +05:30
case ( 'measure' , 'fftw_measure' )
2012-02-02 01:50:05 +05:30
fftw_planner_flag = 0_pInt
2011-12-01 17:31:13 +05:30
case ( 'patient' , 'fftw_patient' )
2012-02-02 01:50:05 +05:30
fftw_planner_flag = 32_pInt
2011-12-01 17:31:13 +05:30
case ( 'exhaustive' , 'fftw_exhaustive' )
2012-02-02 01:50:05 +05:30
fftw_planner_flag = 8_pInt
2011-12-01 17:31:13 +05:30
case default
call IO_warning ( warning_ID = 47_pInt , ext_msg = trim ( IO_lc ( fftw_planner_string ) ) )
2012-02-02 01:50:05 +05:30
fftw_planner_flag = 32_pInt
2011-12-01 17:31:13 +05:30
end select
2009-06-15 18:41:21 +05:30
2012-02-13 23:11:27 +05:30
!* writing parameters to output file
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
!$OMP CRITICAL (write2out)
2012-02-13 23:11:27 +05:30
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,e8.1)' ) ' relevantStrain: ' , relevantStrain
write ( 6 , '(a24,1x,e8.1)' ) ' defgradTolerance: ' , defgradTolerance
write ( 6 , '(a24,1x,i8)' ) ' iJacoStiffness: ' , iJacoStiffness
write ( 6 , '(a24,1x,i8)' ) ' iJacoLpresiduum: ' , iJacoLpresiduum
write ( 6 , '(a24,1x,e8.1)' ) ' pert_Fg: ' , pert_Fg
write ( 6 , '(a24,1x,i8)' ) ' pert_method: ' , pert_method
write ( 6 , '(a24,1x,i8)' ) ' nCryst: ' , nCryst
write ( 6 , '(a24,1x,e8.1)' ) ' subStepMinCryst: ' , subStepMinCryst
write ( 6 , '(a24,1x,e8.1)' ) ' subStepSizeCryst: ' , subStepSizeCryst
write ( 6 , '(a24,1x,e8.1)' ) ' stepIncreaseCryst: ' , stepIncreaseCryst
write ( 6 , '(a24,1x,i8)' ) ' nState: ' , nState
write ( 6 , '(a24,1x,i8)' ) ' nStress: ' , nStress
write ( 6 , '(a24,1x,e8.1)' ) ' rTol_crystalliteState: ' , rTol_crystalliteState
write ( 6 , '(a24,1x,e8.1)' ) ' rTol_crystalliteTemp: ' , rTol_crystalliteTemperature
write ( 6 , '(a24,1x,e8.1)' ) ' rTol_crystalliteStress: ' , rTol_crystalliteStress
write ( 6 , '(a24,1x,e8.1)' ) ' aTol_crystalliteStress: ' , aTol_crystalliteStress
write ( 6 , '(a24,2(1x,i8),/)' ) ' integrator: ' , numerics_integrator
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,i8)' ) ' nHomog: ' , nHomog
write ( 6 , '(a24,1x,e8.1)' ) ' subStepMinHomog: ' , subStepMinHomog
write ( 6 , '(a24,1x,e8.1)' ) ' subStepSizeHomog: ' , subStepSizeHomog
write ( 6 , '(a24,1x,e8.1)' ) ' stepIncreaseHomog: ' , stepIncreaseHomog
write ( 6 , '(a24,1x,i8,/)' ) ' nMPstate: ' , nMPstate
2009-07-31 17:32:20 +05:30
2012-02-13 23:11:27 +05:30
!* RGC parameters
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,e8.1)' ) ' aTol_RGC: ' , absTol_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' rTol_RGC: ' , relTol_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' aMax_RGC: ' , absMax_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' rMax_RGC: ' , relMax_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' perturbPenalty_RGC: ' , pPert_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' relevantMismatch_RGC: ' , xSmoo_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' viscosityrate_RGC: ' , viscPower_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' viscositymodulus_RGC: ' , viscModus_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' maxrelaxation_RGC: ' , maxdRelax_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' maxVolDiscrepancy_RGC: ' , maxVolDiscr_RGC
write ( 6 , '(a24,1x,e8.1)' ) ' volDiscrepancyMod_RGC: ' , volDiscrMod_RGC
write ( 6 , '(a24,1x,e8.1,/)' ) ' discrepancyPower_RGC: ' , volDiscrPow_RGC
2010-10-13 21:34:44 +05:30
2012-02-13 23:11:27 +05:30
!* spectral parameters
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,e8.1)' ) ' err_div_tol: ' , err_div_tol
write ( 6 , '(a24,1x,e8.1)' ) ' err_stress_tolrel: ' , err_stress_tolrel
write ( 6 , '(a24,1x,i8)' ) ' itmax: ' , itmax
write ( 6 , '(a24,1x,L8)' ) ' memory_efficient: ' , memory_efficient
2012-01-25 14:26:46 +05:30
if ( fftw_timelimit < 0.0_pReal ) then
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,L8)' ) ' fftw_timelimit: ' , . false .
2011-10-18 14:46:18 +05:30
else
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,e8.1)' ) ' fftw_timelimit: ' , fftw_timelimit
2011-10-18 14:46:18 +05:30
endif
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,a)' ) ' fftw_planner_string: ' , trim ( fftw_planner_string )
write ( 6 , '(a24,1x,i8)' ) ' fftw_planner_flag: ' , fftw_planner_flag
write ( 6 , '(a24,1x,e8.1)' ) ' rotation_tol: ' , rotation_tol
write ( 6 , '(a24,1x,L8,/)' ) ' divergence_correction: ' , divergence_correction
write ( 6 , '(a24,1x,L8,/)' ) ' update_gamma: ' , update_gamma
write ( 6 , '(a24,1x,L8,/)' ) ' simplified_algorithm: ' , simplified_algorithm
write ( 6 , '(a24,1x,e8.1)' ) ' cut_off_value: ' , cut_off_value
2012-02-13 23:11:27 +05:30
!* Random seeding parameters
2012-02-01 00:48:55 +05:30
write ( 6 , '(a24,1x,i16,/)' ) ' fixed_seed: ' , fixedSeed
2012-02-13 23:11:27 +05:30
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
!$OMP END CRITICAL (write2out)
2010-12-02 16:34:29 +05:30
2012-02-13 23:11:27 +05:30
!* openMP parameter
2012-02-01 00:48:55 +05:30
!$ write(6,'(a24,1x,i8,/)') ' number of threads: ',DAMASK_NumThreadsInt
2009-06-18 19:58:02 +05:30
2012-02-13 23:11:27 +05:30
!* sanity check
if ( relevantStrain < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'relevantStrain' )
if ( defgradTolerance < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'defgradTolerance' )
if ( iJacoStiffness < 1_pInt ) call IO_error ( 301_pInt , ext_msg = 'iJacoStiffness' )
if ( iJacoLpresiduum < 1_pInt ) call IO_error ( 301_pInt , ext_msg = 'iJacoLpresiduum' )
if ( pert_Fg < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'pert_Fg' )
2009-11-10 19:06:27 +05:30
if ( pert_method < = 0_pInt . or . pert_method > = 4_pInt ) &
2012-02-13 23:11:27 +05:30
call IO_error ( 301_pInt , ext_msg = 'pert_method' )
if ( nHomog < 1_pInt ) call IO_error ( 301_pInt , ext_msg = 'nHomog' )
if ( nMPstate < 1_pInt ) call IO_error ( 301_pInt , ext_msg = 'nMPstate' )
if ( nCryst < 1_pInt ) call IO_error ( 301_pInt , ext_msg = 'nCryst' )
if ( nState < 1_pInt ) call IO_error ( 301_pInt , ext_msg = 'nState' )
if ( nStress < 1_pInt ) call IO_error ( 301_pInt , ext_msg = 'nStress' )
if ( subStepMinCryst < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'subStepMinCryst' )
if ( subStepSizeCryst < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'subStepSizeCryst' )
if ( stepIncreaseCryst < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'stepIncreaseCryst' )
if ( subStepMinHomog < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'subStepMinHomog' )
if ( subStepSizeHomog < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'subStepSizeHomog' )
if ( stepIncreaseHomog < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'stepIncreaseHomog' )
if ( rTol_crystalliteState < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'rTol_crystalliteState' )
if ( rTol_crystalliteTemperature < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'rTol_crystalliteTemperature' )
if ( rTol_crystalliteStress < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'rTol_crystalliteStress' )
if ( aTol_crystalliteStress < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'aTol_crystalliteStress' )
2011-02-23 18:00:52 +05:30
if ( any ( numerics_integrator < = 0_pInt ) . or . any ( numerics_integrator > = 6_pInt ) ) &
2012-02-13 23:11:27 +05:30
call IO_error ( 301_pInt , ext_msg = 'integrator' )
2009-07-31 17:32:20 +05:30
2012-02-13 23:11:27 +05:30
!* RGC parameters
2012-02-09 18:05:55 +05:30
2012-02-13 23:11:27 +05:30
if ( absTol_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'absTol_RGC' )
if ( relTol_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'relTol_RGC' )
if ( absMax_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'absMax_RGC' )
if ( relMax_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'relMax_RGC' )
if ( pPert_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'pPert_RGC' )
if ( xSmoo_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'xSmoo_RGC' )
if ( viscPower_RGC < 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'viscPower_RGC' )
if ( viscModus_RGC < 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'viscModus_RGC' )
if ( refRelaxRate_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'refRelaxRate_RGC' )
if ( maxdRelax_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'maxdRelax_RGC' )
if ( maxVolDiscr_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'maxVolDiscr_RGC' )
if ( volDiscrMod_RGC < 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'volDiscrMod_RGC' )
if ( volDiscrPow_RGC < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'volDiscrPw_RGC' )
2010-10-13 21:34:44 +05:30
2012-02-13 23:11:27 +05:30
!* spectral parameters
if ( err_div_tol < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'err_div_tol' )
if ( err_stress_tolrel < = 0.0_pReal ) call IO_error ( 301_pInt , ext_msg = 'err_stress_tolrel' )
if ( itmax < = 1.0_pInt ) call IO_error ( 301_pInt , ext_msg = 'itmax' )
2010-10-13 21:34:44 +05:30
openmp parallelization working again (at least for j2 and nonlocal constitutive model).
In order to keep it like that, please follow these simple rules:
DON'T use implicit array subscripts:
example: real, dimension(3,3) :: A,B
A(:,2) = B(:,1) <--- DON'T USE
A(1:3,2) = B(1:3,1) <--- BETTER USE
In many cases the use of explicit array subscripts is inevitable for parallelization. Additionally, it is an easy means to prevent memory leaks.
Enclose all write statements with the following:
!$OMP CRITICAL (write2out)
<your write statement>
!$OMP END CRITICAL (write2out)
Whenever you change something in the code and are not sure if it affects parallelization and leads to nonconforming behavior, please ask me and/or Franz to check this.
2011-03-17 16:16:17 +05:30
if ( fixedSeed < = 0_pInt ) then
!$OMP CRITICAL (write2out)
write ( 6 , '(a)' ) 'Random is random!'
!$OMP END CRITICAL (write2out)
endif
2009-06-15 18:41:21 +05:30
endsubroutine
2009-07-31 17:32:20 +05:30
END MODULE numerics