2016-10-24 14:03:01 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-30 19:36:22 +05:30
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Utilities used by the different spectral solver variants
!--------------------------------------------------------------------------------------------------
2015-12-16 02:15:54 +05:30
module spectral_utilities
2019-03-25 23:47:10 +05:30
use , intrinsic :: iso_c_binding
2018-05-17 15:34:21 +05:30
#include <petsc/finclude/petscsys.h>
use PETScSys
2019-06-10 00:50:38 +05:30
2019-06-06 21:58:10 +05:30
use prec
2020-03-09 18:35:49 +05:30
use DAMASK_interface
2020-09-13 13:49:38 +05:30
use parallelization
2019-06-06 21:58:10 +05:30
use math
2019-12-02 21:07:22 +05:30
use rotations
2019-06-06 21:58:10 +05:30
use IO
2020-03-20 11:28:55 +05:30
use discretization_grid
2019-06-06 21:58:10 +05:30
use discretization
use homogenization
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
implicit none
private
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
include 'fftw3-mpi.f03'
2020-02-22 20:18:40 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
! field labels information
2020-03-17 12:47:14 +05:30
enum , bind ( c ) ; enumerator :: &
FIELD_UNDEFINED_ID , &
FIELD_MECH_ID , &
FIELD_THERMAL_ID , &
FIELD_DAMAGE_ID
2019-03-25 23:47:10 +05:30
end enum
2015-10-15 00:06:19 +05:30
2013-05-08 21:22:29 +05:30
!--------------------------------------------------------------------------------------------------
! grid related information information
2019-10-24 02:08:17 +05:30
real ( pReal ) , protected , public :: wgt !< weighting factor 1/Nelems
integer , protected , public :: grid1Red !< grid(1)/2
real ( pReal ) , protected , public , dimension ( 3 ) :: scaledGeomSize !< scaled geometry size for calculation of divergence
2015-10-15 00:06:19 +05:30
2012-07-23 15:42:31 +05:30
!--------------------------------------------------------------------------------------------------
! variables storing information for spectral method and FFTW
2019-10-24 02:08:17 +05:30
2019-03-25 23:47:10 +05:30
real ( C_DOUBLE ) , public , dimension ( : , : , : , : , : ) , pointer :: tensorField_real !< real representation (some stress or deformation) of field_fourier
complex ( C_DOUBLE_COMPLEX ) , public , dimension ( : , : , : , : , : ) , pointer :: tensorField_fourier !< field on which the Fourier transform operates
real ( C_DOUBLE ) , public , dimension ( : , : , : , : ) , pointer :: vectorField_real !< vector field real representation for fftw
complex ( C_DOUBLE_COMPLEX ) , public , dimension ( : , : , : , : ) , pointer :: vectorField_fourier !< vector field fourier representation for fftw
real ( C_DOUBLE ) , public , dimension ( : , : , : ) , pointer :: scalarField_real !< scalar field real representation for fftw
complex ( C_DOUBLE_COMPLEX ) , public , dimension ( : , : , : ) , pointer :: scalarField_fourier !< scalar field fourier representation for fftw
complex ( pReal ) , private , dimension ( : , : , : , : , : , : , : ) , allocatable :: gamma_hat !< gamma operator (field) for spectral method
complex ( pReal ) , private , dimension ( : , : , : , : ) , allocatable :: xi1st !< wave vector field for first derivatives
complex ( pReal ) , private , dimension ( : , : , : , : ) , allocatable :: xi2nd !< wave vector field for second derivatives
real ( pReal ) , private , dimension ( 3 , 3 , 3 , 3 ) :: C_ref !< mechanic reference stiffness
2019-10-24 02:08:17 +05:30
2015-10-15 00:06:19 +05:30
2012-11-12 19:44:39 +05:30
!--------------------------------------------------------------------------------------------------
! plans for FFTW
2019-03-25 23:47:10 +05:30
type ( C_PTR ) , private :: &
planTensorForth , & !< FFTW MPI plan P(x) to P(k)
planTensorBack , & !< FFTW MPI plan F(k) to F(x)
planVectorForth , & !< FFTW MPI plan v(x) to v(k)
planVectorBack , & !< FFTW MPI plan v(k) to v(x)
planScalarForth , & !< FFTW MPI plan s(x) to s(k)
planScalarBack !< FFTW MPI plan s(k) to s(x)
2012-11-12 19:44:39 +05:30
2012-07-23 15:42:31 +05:30
!--------------------------------------------------------------------------------------------------
2012-10-24 17:01:40 +05:30
! variables controlling debugging
2019-03-25 23:47:10 +05:30
logical , private :: &
debugGeneral , & !< general debugging of spectral solver
debugRotation , & !< also printing out results in lab frame
debugPETSc !< use some in debug defined options for more verbose PETSc solution
2012-07-23 15:42:31 +05:30
2012-07-25 19:31:39 +05:30
!--------------------------------------------------------------------------------------------------
2012-08-03 14:55:48 +05:30
! derived types
2019-03-25 23:47:10 +05:30
type , public :: tSolutionState !< return type of solution from spectral solver variants
integer :: &
iterationsNeeded = 0
logical :: &
converged = . true . , &
stagConverged = . true . , &
termIll = . false .
end type tSolutionState
type , public :: tBoundaryCondition !< set of parameters defining a boundary condition
real ( pReal ) , dimension ( 3 , 3 ) :: values = 0.0_pReal , &
maskFloat = 0.0_pReal
logical , dimension ( 3 , 3 ) :: maskLogical = . false .
2020-01-26 16:54:35 +05:30
character ( len = pStringLen ) :: myType = 'None'
2019-03-25 23:47:10 +05:30
end type tBoundaryCondition
type , public :: tLoadCase
2019-12-02 21:07:22 +05:30
type ( rotation ) :: rot !< rotation of BC
2019-03-25 23:47:10 +05:30
type ( tBoundaryCondition ) :: stress , & !< stress BC
deformation !< deformation BC (Fdot or L)
real ( pReal ) :: time = 0.0_pReal !< length of increment
integer :: incs = 0 , & !< number of increments
outputfrequency = 1 , & !< frequency of result writes
2019-06-29 05:39:27 +05:30
restartfrequency = huge ( 0 ) , & !< frequency of restart writes
2019-03-25 23:47:10 +05:30
logscale = 0 !< linear/logarithmic time inc flag
logical :: followFormerTrajectory = . true . !< follow trajectory of former loadcase
integer ( kind ( FIELD_UNDEFINED_ID ) ) , allocatable :: ID ( : )
end type tLoadCase
type , public :: tSolutionParams !< @todo use here the type definition for a full loadcase
2019-12-03 00:36:58 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: stress_mask , stress_BC
type ( rotation ) :: rotation_BC
2019-03-25 23:47:10 +05:30
real ( pReal ) :: timeinc
real ( pReal ) :: timeincOld
end type tSolutionParams
2020-02-22 20:18:40 +05:30
type , private :: tNumerics
2019-03-25 23:47:10 +05:30
real ( pReal ) :: &
FFTW_timelimit !< timelimit for FFTW plan creation, see www.fftw.org
integer :: &
2019-03-26 12:06:55 +05:30
divergence_correction !< scale divergence/curl calculation: [0: no correction, 1: size scaled to 1, 2: size scaled to Npoints]
2019-03-25 23:47:10 +05:30
logical :: &
2019-03-26 12:06:55 +05:30
memory_efficient !< calculate gamma operator on the fly
2020-06-16 22:17:19 +05:30
character ( len = : ) , allocatable :: &
2019-03-26 12:06:55 +05:30
spectral_derivative , & !< approximation used for derivatives in Fourier space
2020-06-29 20:35:11 +05:30
FFTW_plan_mode !< FFTW plan mode, see www.fftw.org
2019-03-25 23:47:10 +05:30
end type tNumerics
2020-02-22 20:18:40 +05:30
2020-06-24 20:15:13 +05:30
type ( tNumerics ) , private :: num ! numerics parameters. Better name?
2019-03-25 23:47:10 +05:30
2020-03-17 12:47:14 +05:30
enum , bind ( c ) ; enumerator :: &
DERIVATIVE_CONTINUOUS_ID , &
DERIVATIVE_CENTRAL_DIFF_ID , &
DERIVATIVE_FWBW_DIFF_ID
2019-03-25 23:47:10 +05:30
end enum
integer ( kind ( DERIVATIVE_CONTINUOUS_ID ) ) :: &
spectral_derivative_ID
public :: &
2020-06-17 21:32:22 +05:30
spectral_utilities_init , &
2019-03-25 23:47:10 +05:30
utilities_updateGamma , &
utilities_FFTtensorForward , &
utilities_FFTtensorBackward , &
utilities_FFTvectorForward , &
utilities_FFTvectorBackward , &
utilities_FFTscalarForward , &
utilities_FFTscalarBackward , &
utilities_fourierGammaConvolution , &
utilities_fourierGreenConvolution , &
utilities_divergenceRMS , &
utilities_curlRMS , &
utilities_fourierScalarGradient , &
utilities_fourierVectorDivergence , &
utilities_fourierVectorGradient , &
utilities_fourierTensorDivergence , &
utilities_maskedCompliance , &
utilities_constitutiveResponse , &
utilities_calculateRate , &
utilities_forwardField , &
2019-09-28 03:32:36 +05:30
utilities_updateCoords , &
2019-10-24 02:20:01 +05:30
utilities_saveReferenceStiffness , &
2019-03-25 23:47:10 +05:30
FIELD_UNDEFINED_ID , &
FIELD_MECH_ID , &
FIELD_THERMAL_ID , &
FIELD_DAMAGE_ID
2012-08-03 14:55:48 +05:30
2015-10-15 00:06:19 +05:30
contains
2012-07-19 22:54:56 +05:30
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
2012-11-12 19:44:39 +05:30
!> @brief allocates all neccessary fields, sets debug flags, create plans for FFTW
2018-02-17 03:11:07 +05:30
!> @details Sets the debug levels for general, divergence, restart, and FFTW from the bitwise coding
2012-08-09 18:34:56 +05:30
!> provided by the debug module to logicals.
2018-02-17 03:11:07 +05:30
!> Allocate all fields used by FFTW and create the corresponding plans depending on the debug
2012-08-09 18:34:56 +05:30
!> level chosen.
!> Initializes FFTW.
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
2020-06-17 21:32:22 +05:30
subroutine spectral_utilities_init
2019-06-06 21:58:10 +05:30
2019-03-25 23:47:10 +05:30
PetscErrorCode :: ierr
integer :: i , j , k , &
FFTW_planner_flag
integer , dimension ( 3 ) :: k_s
type ( C_PTR ) :: &
tensorField , & !< field containing data for FFTW in real and fourier space (in place)
vectorField , & !< field containing data for FFTW in real space when debugging FFTW (no in place)
scalarField !< field containing data for FFTW in real space when debugging FFTW (no in place)
integer ( C_INTPTR_T ) , dimension ( 3 ) :: gridFFTW
integer ( C_INTPTR_T ) :: alloc_local , local_K , local_K_offset
integer ( C_INTPTR_T ) , parameter :: &
scalarSize = 1_C_INTPTR_T , &
vecSize = 3_C_INTPTR_T , &
tensorSize = 9_C_INTPTR_T
2020-07-03 20:15:11 +05:30
character ( len = * ) , parameter :: &
2020-06-18 21:44:53 +05:30
PETSCDEBUG = ' -snes_view -snes_monitor '
2020-07-03 20:15:11 +05:30
class ( tNode ) , pointer :: &
2020-06-18 02:30:03 +05:30
num_grid , &
2020-06-26 23:42:05 +05:30
debug_grid ! pointer to grid debug options
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
write ( 6 , '(/,a)' ) ' <<<+- spectral_utilities init -+>>>'
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
write ( 6 , '(/,a)' ) ' Diehl, Diploma Thesis TU München, 2010'
write ( 6 , '(a)' ) ' https://doi.org/10.13140/2.1.3234.3840'
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
write ( 6 , '(/,a)' ) ' Eisenlohr et al., International Journal of Plasticity 46:37– 53, 2013'
write ( 6 , '(a)' ) ' https://doi.org/10.1016/j.ijplas.2012.09.012'
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
write ( 6 , '(/,a)' ) ' Shanthraj et al., International Journal of Plasticity 66:31– 45, 2015'
write ( 6 , '(a)' ) ' https://doi.org/10.1016/j.ijplas.2014.02.006'
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
write ( 6 , '(/,a)' ) ' Shanthraj et al., Handbook of Mechanics of Materials, 2019'
write ( 6 , '(a)' ) ' https://doi.org/10.1007/978-981-10-6855-3_80'
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-23 15:42:31 +05:30
! set debugging parameters
2020-09-13 14:09:17 +05:30
debug_grid = > config_debug % get ( 'grid' , defaultVal = emptyList )
2020-06-18 21:13:25 +05:30
debugGeneral = debug_grid % contains ( 'basic' )
debugRotation = debug_grid % contains ( 'rotation' )
debugPETSc = debug_grid % contains ( 'petsc' )
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
if ( debugPETSc ) write ( 6 , '(3(/,a),/)' ) &
' Initializing PETSc with debug options: ' , &
trim ( PETScDebug ) , &
2020-06-29 20:35:11 +05:30
' add more using the PETSc_Options keyword in numerics.yaml ' ; flush ( 6 )
2020-02-22 20:18:40 +05:30
2020-09-13 14:09:17 +05:30
num_grid = > config_numerics % get ( 'grid' , defaultVal = emptyDict )
2020-06-18 02:30:03 +05:30
2019-03-25 23:47:10 +05:30
call PETScOptionsClear ( PETSC_NULL_OPTIONS , ierr )
CHKERRQ ( ierr )
if ( debugPETSc ) call PETScOptionsInsertString ( PETSC_NULL_OPTIONS , trim ( PETSCDEBUG ) , ierr )
CHKERRQ ( ierr )
2020-06-29 18:39:13 +05:30
call PETScOptionsInsertString ( PETSC_NULL_OPTIONS , &
num_grid % get_asString ( 'petsc_options' , defaultVal = '' ) , ierr )
2019-03-25 23:47:10 +05:30
CHKERRQ ( ierr )
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
grid1Red = grid ( 1 ) / 2 + 1
wgt = 1.0 / real ( product ( grid ) , pReal )
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
write ( 6 , '(/,a,3(i12 ))' ) ' grid a b c: ' , grid
write ( 6 , '(a,3(es12.5))' ) ' size x y z: ' , geomSize
2020-06-29 18:39:13 +05:30
2020-06-16 22:17:19 +05:30
num % memory_efficient = num_grid % get_asInt ( 'memory_efficient' , defaultVal = 1 ) > 0
num % FFTW_timelimit = num_grid % get_asFloat ( 'fftw_timelimit' , defaultVal = - 1.0_pReal )
num % divergence_correction = num_grid % get_asInt ( 'divergence_correction' , defaultVal = 2 )
num % spectral_derivative = num_grid % get_asString ( 'derivative' , defaultVal = 'continuous' )
num % FFTW_plan_mode = num_grid % get_asString ( 'fftw_plan_mode' , defaultVal = 'FFTW_MEASURE' )
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
if ( num % divergence_correction < 0 . or . num % divergence_correction > 2 ) &
call IO_error ( 301 , ext_msg = 'divergence_correction' )
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
select case ( num % spectral_derivative )
case ( 'continuous' )
spectral_derivative_ID = DERIVATIVE_CONTINUOUS_ID
case ( 'central_difference' )
spectral_derivative_ID = DERIVATIVE_CENTRAL_DIFF_ID
2020-06-16 22:17:19 +05:30
case ( 'FWBW_difference' )
2019-03-25 23:47:10 +05:30
spectral_derivative_ID = DERIVATIVE_FWBW_DIFF_ID
case default
call IO_error ( 892 , ext_msg = trim ( num % spectral_derivative ) )
end select
2015-12-14 23:42:09 +05:30
2013-05-08 21:22:29 +05:30
!--------------------------------------------------------------------------------------------------
2015-09-11 14:22:03 +05:30
! scale dimension to calculate either uncorrected, dimension-independent, or dimension- and
! resolution-independent divergence
2019-03-25 23:47:10 +05:30
if ( num % divergence_correction == 1 ) then
do j = 1 , 3
if ( j / = minloc ( geomSize , 1 ) . and . j / = maxloc ( geomSize , 1 ) ) &
scaledGeomSize = geomSize / geomSize ( j )
enddo
elseif ( num % divergence_correction == 2 ) then
do j = 1 , 3
if ( j / = int ( minloc ( geomSize / real ( grid , pReal ) , 1 ) ) &
. and . j / = int ( maxloc ( geomSize / real ( grid , pReal ) , 1 ) ) ) &
scaledGeomSize = geomSize / geomSize ( j ) * real ( grid ( j ) , pReal )
enddo
else
scaledGeomSize = geomSize
endif
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
select case ( IO_lc ( num % FFTW_plan_mode ) ) ! setting parameters for the plan creation of FFTW. Basically a translation from fftw3.f
2019-03-29 13:18:32 +05:30
case ( 'fftw_estimate' ) ! ordered from slow execution (but fast plan creation) to fast execution
FFTW_planner_flag = FFTW_ESTIMATE
case ( 'fftw_measure' )
FFTW_planner_flag = FFTW_MEASURE
case ( 'fftw_patient' )
FFTW_planner_flag = FFTW_PATIENT
case ( 'fftw_exhaustive' )
FFTW_planner_flag = FFTW_EXHAUSTIVE
2019-03-25 23:47:10 +05:30
case default
call IO_warning ( warning_ID = 47 , ext_msg = trim ( IO_lc ( num % FFTW_plan_mode ) ) )
2019-03-29 13:18:32 +05:30
FFTW_planner_flag = FFTW_MEASURE
2019-03-25 23:47:10 +05:30
end select
2019-03-25 20:24:51 +05:30
2019-03-25 22:58:07 +05:30
!--------------------------------------------------------------------------------------------------
! general initialization of FFTW (see manual on fftw.org for more details)
2020-09-13 15:50:44 +05:30
if ( pReal / = C_DOUBLE . or . kind ( 1 ) / = C_INT ) error stop 'C and Fortran datatypes do not match'
call fftw_set_timelimit ( num % FFTW_timelimit )
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
if ( debugGeneral ) write ( 6 , '(/,a)' ) ' FFTW initialized' ; flush ( 6 )
2019-03-25 22:58:07 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2015-03-25 21:36:19 +05:30
! MPI allocation
2019-03-25 23:47:10 +05:30
gridFFTW = int ( grid , C_INTPTR_T )
alloc_local = fftw_mpi_local_size_3d ( gridFFTW ( 3 ) , gridFFTW ( 2 ) , gridFFTW ( 1 ) / 2 + 1 , &
PETSC_COMM_WORLD , local_K , local_K_offset )
allocate ( xi1st ( 3 , grid1Red , grid ( 2 ) , grid3 ) , source = cmplx ( 0.0_pReal , 0.0_pReal , pReal ) ) ! frequencies for first derivatives, only half the size for first dimension
allocate ( xi2nd ( 3 , grid1Red , grid ( 2 ) , grid3 ) , source = cmplx ( 0.0_pReal , 0.0_pReal , pReal ) ) ! frequencies for second derivatives, only half the size for first dimension
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
tensorField = fftw_alloc_complex ( tensorSize * alloc_local )
call c_f_pointer ( tensorField , tensorField_real , [ 3_C_INTPTR_T , 3_C_INTPTR_T , &
2_C_INTPTR_T * ( gridFFTW ( 1 ) / 2_C_INTPTR_T + 1_C_INTPTR_T ) , gridFFTW ( 2 ) , local_K ] ) ! place a pointer for a real tensor representation
call c_f_pointer ( tensorField , tensorField_fourier , [ 3_C_INTPTR_T , 3_C_INTPTR_T , &
gridFFTW ( 1 ) / 2_C_INTPTR_T + 1_C_INTPTR_T , gridFFTW ( 2 ) , local_K ] ) ! place a pointer for a fourier tensor representation
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
vectorField = fftw_alloc_complex ( vecSize * alloc_local )
call c_f_pointer ( vectorField , vectorField_real , [ 3_C_INTPTR_T , &
2_C_INTPTR_T * ( gridFFTW ( 1 ) / 2_C_INTPTR_T + 1_C_INTPTR_T ) , gridFFTW ( 2 ) , local_K ] ) ! place a pointer for a real vector representation
call c_f_pointer ( vectorField , vectorField_fourier , [ 3_C_INTPTR_T , &
gridFFTW ( 1 ) / 2_C_INTPTR_T + 1_C_INTPTR_T , gridFFTW ( 2 ) , local_K ] ) ! place a pointer for a fourier vector representation
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
scalarField = fftw_alloc_complex ( scalarSize * alloc_local ) ! allocate data for real representation (no in place transform)
call c_f_pointer ( scalarField , scalarField_real , &
[ 2_C_INTPTR_T * ( gridFFTW ( 1 ) / 2_C_INTPTR_T + 1 ) , gridFFTW ( 2 ) , local_K ] ) ! place a pointer for a real scalar representation
call c_f_pointer ( scalarField , scalarField_fourier , &
[ gridFFTW ( 1 ) / 2_C_INTPTR_T + 1 , gridFFTW ( 2 ) , local_K ] ) ! place a pointer for a fourier scarlar representation
2015-10-15 00:06:19 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
! tensor MPI fftw plans
2019-03-25 23:47:10 +05:30
planTensorForth = fftw_mpi_plan_many_dft_r2c ( 3 , [ gridFFTW ( 3 ) , gridFFTW ( 2 ) , gridFFTW ( 1 ) ] , & ! dimension, logical length in each dimension in reversed order
tensorSize , FFTW_MPI_DEFAULT_BLOCK , FFTW_MPI_DEFAULT_BLOCK , & ! no. of transforms, default iblock and oblock
tensorField_real , tensorField_fourier , & ! input data, output data
PETSC_COMM_WORLD , FFTW_planner_flag ) ! use all processors, planer precision
if ( . not . C_ASSOCIATED ( planTensorForth ) ) call IO_error ( 810 , ext_msg = 'planTensorForth' )
planTensorBack = fftw_mpi_plan_many_dft_c2r ( 3 , [ gridFFTW ( 3 ) , gridFFTW ( 2 ) , gridFFTW ( 1 ) ] , & ! dimension, logical length in each dimension in reversed order
tensorSize , FFTW_MPI_DEFAULT_BLOCK , FFTW_MPI_DEFAULT_BLOCK , & ! no. of transforms, default iblock and oblock
tensorField_fourier , tensorField_real , & ! input data, output data
PETSC_COMM_WORLD , FFTW_planner_flag ) ! all processors, planer precision
if ( . not . C_ASSOCIATED ( planTensorBack ) ) call IO_error ( 810 , ext_msg = 'planTensorBack' )
2015-10-15 00:06:19 +05:30
2015-03-13 03:58:33 +05:30
!--------------------------------------------------------------------------------------------------
2015-06-03 23:00:31 +05:30
! vector MPI fftw plans
2019-03-25 23:47:10 +05:30
planVectorForth = fftw_mpi_plan_many_dft_r2c ( 3 , [ gridFFTW ( 3 ) , gridFFTW ( 2 ) , gridFFTW ( 1 ) ] , & ! dimension, logical length in each dimension in reversed order
vecSize , FFTW_MPI_DEFAULT_BLOCK , FFTW_MPI_DEFAULT_BLOCK , & ! no. of transforms, default iblock and oblock
vectorField_real , vectorField_fourier , & ! input data, output data
PETSC_COMM_WORLD , FFTW_planner_flag ) ! use all processors, planer precision
if ( . not . C_ASSOCIATED ( planVectorForth ) ) call IO_error ( 810 , ext_msg = 'planVectorForth' )
planVectorBack = fftw_mpi_plan_many_dft_c2r ( 3 , [ gridFFTW ( 3 ) , gridFFTW ( 2 ) , gridFFTW ( 1 ) ] , & ! dimension, logical length in each dimension in reversed order
vecSize , FFTW_MPI_DEFAULT_BLOCK , FFTW_MPI_DEFAULT_BLOCK , & ! no. of transforms, default iblock and oblock
vectorField_fourier , vectorField_real , & ! input data, output data
PETSC_COMM_WORLD , FFTW_planner_flag ) ! all processors, planer precision
if ( . not . C_ASSOCIATED ( planVectorBack ) ) call IO_error ( 810 , ext_msg = 'planVectorBack' )
2015-10-15 00:06:19 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-15 00:06:19 +05:30
! scalar MPI fftw plans
2019-03-25 23:47:10 +05:30
planScalarForth = fftw_mpi_plan_many_dft_r2c ( 3 , [ gridFFTW ( 3 ) , gridFFTW ( 2 ) , gridFFTW ( 1 ) ] , & ! dimension, logical length in each dimension in reversed order
scalarSize , FFTW_MPI_DEFAULT_BLOCK , FFTW_MPI_DEFAULT_BLOCK , & ! no. of transforms, default iblock and oblock
scalarField_real , scalarField_fourier , & ! input data, output data
PETSC_COMM_WORLD , FFTW_planner_flag ) ! use all processors, planer precision
if ( . not . C_ASSOCIATED ( planScalarForth ) ) call IO_error ( 810 , ext_msg = 'planScalarForth' )
planScalarBack = fftw_mpi_plan_many_dft_c2r ( 3 , [ gridFFTW ( 3 ) , gridFFTW ( 2 ) , gridFFTW ( 1 ) ] , & ! dimension, logical length in each dimension in reversed order, no. of transforms
scalarSize , FFTW_MPI_DEFAULT_BLOCK , FFTW_MPI_DEFAULT_BLOCK , & ! no. of transforms, default iblock and oblock
scalarField_fourier , scalarField_real , & ! input data, output data
PETSC_COMM_WORLD , FFTW_planner_flag ) ! use all processors, planer precision
if ( . not . C_ASSOCIATED ( planScalarBack ) ) call IO_error ( 810 , ext_msg = 'planScalarBack' )
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-23 15:42:31 +05:30
! calculation of discrete angular frequencies, ordered as in FFTW (wrap around)
2019-03-25 23:47:10 +05:30
do k = grid3Offset + 1 , grid3Offset + grid3
k_s ( 3 ) = k - 1
if ( k > grid ( 3 ) / 2 + 1 ) k_s ( 3 ) = k_s ( 3 ) - grid ( 3 ) ! running from 0,1,...,N/2,N/2+1,-N/2,-N/2+1,...,-1
do j = 1 , grid ( 2 )
k_s ( 2 ) = j - 1
if ( j > grid ( 2 ) / 2 + 1 ) k_s ( 2 ) = k_s ( 2 ) - grid ( 2 ) ! running from 0,1,...,N/2,N/2+1,-N/2,-N/2+1,...,-1
do i = 1 , grid1Red
k_s ( 1 ) = i - 1 ! symmetry, junst running from 0,1,...,N/2,N/2+1
xi2nd ( 1 : 3 , i , j , k - grid3Offset ) = utilities_getFreqDerivative ( k_s )
where ( mod ( grid , 2 ) == 0 . and . [ i , j , k ] == grid / 2 + 1 . and . &
spectral_derivative_ID == DERIVATIVE_CONTINUOUS_ID ) ! for even grids, set the Nyquist Freq component to 0.0
xi1st ( 1 : 3 , i , j , k - grid3Offset ) = cmplx ( 0.0_pReal , 0.0_pReal , pReal )
elsewhere
xi1st ( 1 : 3 , i , j , k - grid3Offset ) = xi2nd ( 1 : 3 , i , j , k - grid3Offset )
endwhere
enddo ; enddo ; enddo
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
if ( num % memory_efficient ) then ! allocate just single fourth order tensor
allocate ( gamma_hat ( 3 , 3 , 3 , 3 , 1 , 1 , 1 ) , source = cmplx ( 0.0_pReal , 0.0_pReal , pReal ) )
else ! precalculation of gamma_hat field
allocate ( gamma_hat ( 3 , 3 , 3 , 3 , grid1Red , grid ( 2 ) , grid3 ) , source = cmplx ( 0.0_pReal , 0.0_pReal , pReal ) )
endif
2012-08-03 14:55:48 +05:30
2020-06-17 21:32:22 +05:30
end subroutine spectral_utilities_init
2012-10-24 17:01:40 +05:30
2012-07-26 19:28:47 +05:30
2019-03-26 12:06:55 +05:30
!---------------------------------------------------------------------------------------------------
2018-02-17 03:11:07 +05:30
!> @brief updates reference stiffness and potentially precalculated gamma operator
2012-08-09 18:34:56 +05:30
!> @details Sets the current reference stiffness to the stiffness given as an argument.
!> If the gamma operator is precalculated, it is calculated with this stiffness.
2018-02-17 03:11:07 +05:30
!> In case of an on-the-fly calculation, only the reference stiffness is updated.
2019-03-26 12:06:55 +05:30
!---------------------------------------------------------------------------------------------------
2019-10-24 02:20:01 +05:30
subroutine utilities_updateGamma ( C )
2020-02-22 20:18:40 +05:30
2019-03-09 03:46:56 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 , 3 , 3 ) :: C !< input stiffness to store as reference stiffness
complex ( pReal ) , dimension ( 3 , 3 ) :: temp33_complex , xiDyad_cmplx
real ( pReal ) , dimension ( 6 , 6 ) :: A , A_inv
integer :: &
i , j , k , &
2019-10-24 02:20:01 +05:30
l , m , n , o
2019-03-09 03:46:56 +05:30
logical :: err
2020-02-22 20:18:40 +05:30
2019-03-09 03:46:56 +05:30
C_ref = C
2020-02-22 20:18:40 +05:30
2019-03-25 14:06:59 +05:30
if ( . not . num % memory_efficient ) then
2019-03-25 20:24:51 +05:30
gamma_hat = cmplx ( 0.0_pReal , 0.0_pReal , pReal ) ! for the singular point and any non invertible A
2019-03-09 03:46:56 +05:30
do k = grid3Offset + 1 , grid3Offset + grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
2019-03-25 20:24:51 +05:30
if ( any ( [ i , j , k ] / = 1 ) ) then ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1
2019-03-09 03:46:56 +05:30
forall ( l = 1 : 3 , m = 1 : 3 ) &
xiDyad_cmplx ( l , m ) = conjg ( - xi1st ( l , i , j , k - grid3Offset ) ) * xi1st ( m , i , j , k - grid3Offset )
forall ( l = 1 : 3 , m = 1 : 3 ) &
temp33_complex ( l , m ) = sum ( cmplx ( C_ref ( l , 1 : 3 , m , 1 : 3 ) , 0.0_pReal ) * xiDyad_cmplx )
A ( 1 : 3 , 1 : 3 ) = real ( temp33_complex ) ; A ( 4 : 6 , 4 : 6 ) = real ( temp33_complex )
A ( 1 : 3 , 4 : 6 ) = aimag ( temp33_complex ) ; A ( 4 : 6 , 1 : 3 ) = - aimag ( temp33_complex )
if ( abs ( math_det33 ( A ( 1 : 3 , 1 : 3 ) ) ) > 1e-16 ) then
2019-09-21 06:50:33 +05:30
call math_invert ( A_inv , err , A )
2019-03-09 03:46:56 +05:30
temp33_complex = cmplx ( A_inv ( 1 : 3 , 1 : 3 ) , A_inv ( 1 : 3 , 4 : 6 ) , pReal )
forall ( l = 1 : 3 , m = 1 : 3 , n = 1 : 3 , o = 1 : 3 ) &
gamma_hat ( l , m , n , o , i , j , k - grid3Offset ) = temp33_complex ( l , n ) * &
2015-12-14 23:42:09 +05:30
conjg ( - xi1st ( o , i , j , k - grid3Offset ) ) * xi1st ( m , i , j , k - grid3Offset )
2019-03-09 03:46:56 +05:30
endif
endif
enddo ; enddo ; enddo
endif
2020-02-22 20:18:40 +05:30
2012-10-24 17:01:40 +05:30
end subroutine utilities_updateGamma
2019-03-09 03:46:56 +05:30
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
2015-09-11 18:35:46 +05:30
!> @brief forward FFT of data in field_real to field_fourier
2019-10-29 20:48:58 +05:30
!> @details Does an unweighted FFT transform from real to complex. Extra padding entries are set
! to 0.0
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
2019-03-25 14:06:59 +05:30
subroutine utilities_FFTtensorForward
2013-03-07 01:04:30 +05:30
2020-03-09 18:37:31 +05:30
tensorField_real ( 1 : 3 , 1 : 3 , grid ( 1 ) + 1 : grid1Red * 2 , : , : ) = 0.0_pReal
2019-03-25 23:47:10 +05:30
call fftw_mpi_execute_dft_r2c ( planTensorForth , tensorField_real , tensorField_fourier )
2015-10-15 00:06:19 +05:30
2015-06-03 23:00:31 +05:30
end subroutine utilities_FFTtensorForward
2012-07-26 19:28:47 +05:30
2012-08-09 18:34:56 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief backward FFT of data in field_fourier to field_real
2015-09-11 14:22:03 +05:30
!> @details Does an weighted inverse FFT transform from complex to real
2012-08-09 18:34:56 +05:30
!--------------------------------------------------------------------------------------------------
2019-03-25 14:06:59 +05:30
subroutine utilities_FFTtensorBackward
2013-04-10 15:49:16 +05:30
2019-03-25 23:47:10 +05:30
call fftw_mpi_execute_dft_c2r ( planTensorBack , tensorField_fourier , tensorField_real )
tensorField_real = tensorField_real * wgt ! normalize the result by number of elements
2012-08-03 14:55:48 +05:30
2015-06-03 23:00:31 +05:30
end subroutine utilities_FFTtensorBackward
2012-07-26 19:28:47 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
2015-09-11 14:22:03 +05:30
!> @brief forward FFT of data in scalarField_real to scalarField_fourier
2019-10-29 20:48:58 +05:30
!> @details Does an unweighted FFT transform from real to complex. Extra padding entries are set
! to 0.0
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
2019-03-25 14:06:59 +05:30
subroutine utilities_FFTscalarForward
2015-10-15 00:06:19 +05:30
2020-03-09 18:37:31 +05:30
scalarField_real ( grid ( 1 ) + 1 : grid1Red * 2 , : , : ) = 0.0_pReal
2019-03-25 23:47:10 +05:30
call fftw_mpi_execute_dft_r2c ( planScalarForth , scalarField_real , scalarField_fourier )
2015-06-03 23:00:31 +05:30
end subroutine utilities_FFTscalarForward
2015-09-24 23:08:49 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
2015-09-11 14:22:03 +05:30
!> @brief backward FFT of data in scalarField_fourier to scalarField_real
!> @details Does an weighted inverse FFT transform from complex to real
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
2019-03-25 14:06:59 +05:30
subroutine utilities_FFTscalarBackward
2015-09-11 14:22:03 +05:30
2019-03-25 23:47:10 +05:30
call fftw_mpi_execute_dft_c2r ( planScalarBack , scalarField_fourier , scalarField_real )
scalarField_real = scalarField_real * wgt ! normalize the result by number of elements
2015-06-03 23:00:31 +05:30
end subroutine utilities_FFTscalarBackward
2015-09-24 23:08:49 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief forward FFT of data in field_real to field_fourier with highest freqs. removed
2019-10-29 20:48:58 +05:30
!> @details Does an unweighted FFT transform from real to complex. Extra padding entries are set
! to 0.0
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
2019-03-25 14:06:59 +05:30
subroutine utilities_FFTvectorForward
2015-09-11 18:35:46 +05:30
2020-03-09 18:37:31 +05:30
vectorField_real ( 1 : 3 , grid ( 1 ) + 1 : grid1Red * 2 , : , : ) = 0.0_pReal
2019-03-25 23:47:10 +05:30
call fftw_mpi_execute_dft_r2c ( planVectorForth , vectorField_real , vectorField_fourier )
2015-06-03 23:00:31 +05:30
end subroutine utilities_FFTvectorForward
2015-09-24 23:08:49 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief backward FFT of data in field_fourier to field_real
2015-09-11 14:22:03 +05:30
!> @details Does an weighted inverse FFT transform from complex to real
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
2019-03-25 14:06:59 +05:30
subroutine utilities_FFTvectorBackward
2015-09-11 14:22:03 +05:30
2019-03-25 23:47:10 +05:30
call fftw_mpi_execute_dft_c2r ( planVectorBack , vectorField_fourier , vectorField_real )
vectorField_real = vectorField_real * wgt ! normalize the result by number of elements
2015-06-03 23:00:31 +05:30
end subroutine utilities_FFTvectorBackward
2012-07-26 19:28:47 +05:30
2013-07-26 21:55:37 +05:30
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
2012-10-24 17:01:40 +05:30
!> @brief doing convolution gamma_hat * field_real, ensuring that average value = fieldAim
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
2015-06-03 23:00:31 +05:30
subroutine utilities_fourierGammaConvolution ( fieldAim )
2020-02-22 20:18:40 +05:30
2019-03-09 12:17:01 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: fieldAim !< desired average value of the field after convolution
complex ( pReal ) , dimension ( 3 , 3 ) :: temp33_complex , xiDyad_cmplx
real ( pReal ) , dimension ( 6 , 6 ) :: A , A_inv
2020-02-22 20:18:40 +05:30
2019-03-09 12:17:01 +05:30
integer :: &
i , j , k , &
l , m , n , o
logical :: err
2020-02-22 20:18:40 +05:30
2019-03-09 12:17:01 +05:30
write ( 6 , '(/,a)' ) ' ... doing gamma convolution ...............................................'
flush ( 6 )
2015-10-15 00:06:19 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2013-04-18 22:10:49 +05:30
! do the actual spectral method calculation (mechanical equilibrium)
2019-03-25 14:06:59 +05:30
memoryEfficient : if ( num % memory_efficient ) then
2019-03-09 12:17:01 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
2019-03-25 23:47:10 +05:30
if ( any ( [ i , j , k + grid3Offset ] / = 1 ) ) then ! singular point at xi=(0.0,0.0,0.0) i.e. i=j=k=1
2019-03-09 12:17:01 +05:30
forall ( l = 1 : 3 , m = 1 : 3 ) &
xiDyad_cmplx ( l , m ) = conjg ( - xi1st ( l , i , j , k ) ) * xi1st ( m , i , j , k )
forall ( l = 1 : 3 , m = 1 : 3 ) &
temp33_complex ( l , m ) = sum ( cmplx ( C_ref ( l , 1 : 3 , m , 1 : 3 ) , 0.0_pReal ) * xiDyad_cmplx )
A ( 1 : 3 , 1 : 3 ) = real ( temp33_complex ) ; A ( 4 : 6 , 4 : 6 ) = real ( temp33_complex )
A ( 1 : 3 , 4 : 6 ) = aimag ( temp33_complex ) ; A ( 4 : 6 , 1 : 3 ) = - aimag ( temp33_complex )
if ( abs ( math_det33 ( A ( 1 : 3 , 1 : 3 ) ) ) > 1e-16 ) then
2019-09-21 06:50:33 +05:30
call math_invert ( A_inv , err , A )
2019-03-09 12:17:01 +05:30
temp33_complex = cmplx ( A_inv ( 1 : 3 , 1 : 3 ) , A_inv ( 1 : 3 , 4 : 6 ) , pReal )
forall ( l = 1 : 3 , m = 1 : 3 , n = 1 : 3 , o = 1 : 3 ) &
gamma_hat ( l , m , n , o , 1 , 1 , 1 ) = temp33_complex ( l , n ) * conjg ( - xi1st ( o , i , j , k ) ) * xi1st ( m , i , j , k )
2020-02-22 20:18:40 +05:30
else
2019-03-09 12:17:01 +05:30
gamma_hat ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , 1 , 1 ) = cmplx ( 0.0_pReal , 0.0_pReal , pReal )
endif
forall ( l = 1 : 3 , m = 1 : 3 ) &
temp33_Complex ( l , m ) = sum ( gamma_hat ( l , m , 1 : 3 , 1 : 3 , 1 , 1 , 1 ) * tensorField_fourier ( 1 : 3 , 1 : 3 , i , j , k ) )
tensorField_fourier ( 1 : 3 , 1 : 3 , i , j , k ) = temp33_Complex
endif
enddo ; enddo ; enddo
else memoryEfficient
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
forall ( l = 1 : 3 , m = 1 : 3 ) &
temp33_Complex ( l , m ) = sum ( gamma_hat ( l , m , 1 : 3 , 1 : 3 , i , j , k ) * tensorField_fourier ( 1 : 3 , 1 : 3 , i , j , k ) )
tensorField_fourier ( 1 : 3 , 1 : 3 , i , j , k ) = temp33_Complex
enddo ; enddo ; enddo
endif memoryEfficient
2020-02-22 20:18:40 +05:30
2019-03-09 12:17:01 +05:30
if ( grid3Offset == 0 ) tensorField_fourier ( 1 : 3 , 1 : 3 , 1 , 1 , 1 ) = cmplx ( fieldAim / wgt , 0.0_pReal , pReal )
2012-07-30 19:36:22 +05:30
2015-06-03 23:00:31 +05:30
end subroutine utilities_fourierGammaConvolution
2015-10-15 00:06:19 +05:30
2015-09-24 23:08:49 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief doing convolution DamageGreenOp_hat * field_real
!--------------------------------------------------------------------------------------------------
subroutine utilities_fourierGreenConvolution ( D_ref , mobility_ref , deltaT )
2019-06-06 21:58:10 +05:30
2019-03-07 11:31:09 +05:30
real ( pReal ) , dimension ( 3 , 3 ) , intent ( in ) :: D_ref
real ( pReal ) , intent ( in ) :: mobility_ref , deltaT
complex ( pReal ) :: GreenOp_hat
2019-03-25 22:58:07 +05:30
integer :: i , j , k
2020-02-22 20:18:40 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
! do the actual spectral method calculation
2019-03-25 22:58:07 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
2019-03-07 11:31:09 +05:30
GreenOp_hat = cmplx ( 1.0_pReal , 0.0_pReal , pReal ) / &
( cmplx ( mobility_ref , 0.0_pReal , pReal ) + cmplx ( deltaT , 0.0_pReal ) * &
2019-03-25 23:47:10 +05:30
sum ( conjg ( xi1st ( 1 : 3 , i , j , k ) ) * matmul ( cmplx ( D_ref , 0.0_pReal ) , xi1st ( 1 : 3 , i , j , k ) ) ) )
2019-03-07 11:31:09 +05:30
scalarField_fourier ( i , j , k ) = scalarField_fourier ( i , j , k ) * GreenOp_hat
enddo ; enddo ; enddo
2015-06-03 23:00:31 +05:30
end subroutine utilities_fourierGreenConvolution
2012-07-30 19:36:22 +05:30
2015-09-24 23:08:49 +05:30
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculate root mean square of divergence of field_fourier
!--------------------------------------------------------------------------------------------------
2015-10-15 00:06:19 +05:30
real ( pReal ) function utilities_divergenceRMS ( )
2012-10-24 17:01:40 +05:30
2019-03-25 22:58:07 +05:30
integer :: i , j , k , ierr
2019-03-07 11:31:09 +05:30
complex ( pReal ) , dimension ( 3 ) :: rescaledGeom
2012-07-30 19:36:22 +05:30
2019-03-07 11:31:09 +05:30
write ( 6 , '(/,a)' ) ' ... calculating divergence ................................................'
flush ( 6 )
2016-09-03 17:57:56 +05:30
2019-03-07 11:31:09 +05:30
rescaledGeom = cmplx ( geomSize / scaledGeomSize , 0.0_pReal )
2013-04-10 15:49:16 +05:30
2012-07-19 22:54:56 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-26 19:28:47 +05:30
! calculating RMS divergence criterion in Fourier space
2019-03-07 11:31:09 +05:30
utilities_divergenceRMS = 0.0_pReal
2019-03-25 22:58:07 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 )
do i = 2 , grid1Red - 1 ! Has somewhere a conj. complex counterpart. Therefore count it twice.
2019-03-07 11:31:09 +05:30
utilities_divergenceRMS = utilities_divergenceRMS &
+ 2.0_pReal * ( sum ( real ( matmul ( tensorField_fourier ( 1 : 3 , 1 : 3 , i , j , k ) , & ! (sqrt(real(a)**2 + aimag(a)**2))**2 = real(a)**2 + aimag(a)**2. do not take square root and square again
conjg ( - xi1st ( 1 : 3 , i , j , k ) ) * rescaledGeom ) ) ** 2.0_pReal ) & ! --> sum squared L_2 norm of vector
+ sum ( aimag ( matmul ( tensorField_fourier ( 1 : 3 , 1 : 3 , i , j , k ) , &
conjg ( - xi1st ( 1 : 3 , i , j , k ) ) * rescaledGeom ) ) ** 2.0_pReal ) )
enddo
utilities_divergenceRMS = utilities_divergenceRMS & ! these two layers (DC and Nyquist) do not have a conjugate complex counterpart (if grid(1) /= 1)
+ sum ( real ( matmul ( tensorField_fourier ( 1 : 3 , 1 : 3 , 1 , j , k ) , &
conjg ( - xi1st ( 1 : 3 , 1 , j , k ) ) * rescaledGeom ) ) ** 2.0_pReal ) &
+ sum ( aimag ( matmul ( tensorField_fourier ( 1 : 3 , 1 : 3 , 1 , j , k ) , &
conjg ( - xi1st ( 1 : 3 , 1 , j , k ) ) * rescaledGeom ) ) ** 2.0_pReal ) &
+ sum ( real ( matmul ( tensorField_fourier ( 1 : 3 , 1 : 3 , grid1Red , j , k ) , &
conjg ( - xi1st ( 1 : 3 , grid1Red , j , k ) ) * rescaledGeom ) ) ** 2.0_pReal ) &
+ sum ( aimag ( matmul ( tensorField_fourier ( 1 : 3 , 1 : 3 , grid1Red , j , k ) , &
conjg ( - xi1st ( 1 : 3 , grid1Red , j , k ) ) * rescaledGeom ) ) ** 2.0_pReal )
enddo ; enddo
2019-03-25 22:58:07 +05:30
if ( grid ( 1 ) == 1 ) utilities_divergenceRMS = utilities_divergenceRMS * 0.5_pReal ! counted twice in case of grid(1) == 1
2019-03-07 11:31:09 +05:30
call MPI_Allreduce ( MPI_IN_PLACE , utilities_divergenceRMS , 1 , MPI_DOUBLE , MPI_SUM , PETSC_COMM_WORLD , ierr )
2019-03-25 22:58:07 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'utilities_divergenceRMS' )
2019-03-07 11:31:09 +05:30
utilities_divergenceRMS = sqrt ( utilities_divergenceRMS ) * wgt ! RMS in real space calculated with Parsevals theorem from Fourier space
2012-10-24 17:01:40 +05:30
end function utilities_divergenceRMS
2015-10-15 00:06:19 +05:30
2013-03-06 20:01:13 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculate max of curl of field_fourier
!--------------------------------------------------------------------------------------------------
real ( pReal ) function utilities_curlRMS ( )
2019-06-06 21:58:10 +05:30
2019-03-25 22:58:07 +05:30
integer :: i , j , k , l , ierr
2019-03-07 11:31:09 +05:30
complex ( pReal ) , dimension ( 3 , 3 ) :: curl_fourier
complex ( pReal ) , dimension ( 3 ) :: rescaledGeom
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
write ( 6 , '(/,a)' ) ' ... calculating curl ......................................................'
flush ( 6 )
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
rescaledGeom = cmplx ( geomSize / scaledGeomSize , 0.0_pReal )
2020-02-22 20:18:40 +05:30
2016-05-27 15:16:34 +05:30
!--------------------------------------------------------------------------------------------------
2013-03-06 20:01:13 +05:30
! calculating max curl criterion in Fourier space
2019-03-07 11:31:09 +05:30
utilities_curlRMS = 0.0_pReal
2020-02-22 20:18:40 +05:30
2019-03-25 22:58:07 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ;
do i = 2 , grid1Red - 1
do l = 1 , 3
2019-03-07 11:31:09 +05:30
curl_fourier ( l , 1 ) = ( + tensorField_fourier ( l , 3 , i , j , k ) * xi1st ( 2 , i , j , k ) * rescaledGeom ( 2 ) &
- tensorField_fourier ( l , 2 , i , j , k ) * xi1st ( 3 , i , j , k ) * rescaledGeom ( 3 ) )
curl_fourier ( l , 2 ) = ( + tensorField_fourier ( l , 1 , i , j , k ) * xi1st ( 3 , i , j , k ) * rescaledGeom ( 3 ) &
- tensorField_fourier ( l , 3 , i , j , k ) * xi1st ( 1 , i , j , k ) * rescaledGeom ( 1 ) )
curl_fourier ( l , 3 ) = ( + tensorField_fourier ( l , 2 , i , j , k ) * xi1st ( 1 , i , j , k ) * rescaledGeom ( 1 ) &
- tensorField_fourier ( l , 1 , i , j , k ) * xi1st ( 2 , i , j , k ) * rescaledGeom ( 2 ) )
enddo
utilities_curlRMS = utilities_curlRMS &
+ 2.0_pReal * sum ( real ( curl_fourier ) ** 2.0_pReal + aimag ( curl_fourier ) ** 2.0_pReal ) ! Has somewhere a conj. complex counterpart. Therefore count it twice.
enddo
2019-03-25 22:58:07 +05:30
do l = 1 , 3
2019-03-07 11:31:09 +05:30
curl_fourier = ( + tensorField_fourier ( l , 3 , 1 , j , k ) * xi1st ( 2 , 1 , j , k ) * rescaledGeom ( 2 ) &
- tensorField_fourier ( l , 2 , 1 , j , k ) * xi1st ( 3 , 1 , j , k ) * rescaledGeom ( 3 ) )
curl_fourier = ( + tensorField_fourier ( l , 1 , 1 , j , k ) * xi1st ( 3 , 1 , j , k ) * rescaledGeom ( 3 ) &
- tensorField_fourier ( l , 3 , 1 , j , k ) * xi1st ( 1 , 1 , j , k ) * rescaledGeom ( 1 ) )
curl_fourier = ( + tensorField_fourier ( l , 2 , 1 , j , k ) * xi1st ( 1 , 1 , j , k ) * rescaledGeom ( 1 ) &
- tensorField_fourier ( l , 1 , 1 , j , k ) * xi1st ( 2 , 1 , j , k ) * rescaledGeom ( 2 ) )
enddo
utilities_curlRMS = utilities_curlRMS &
+ sum ( real ( curl_fourier ) ** 2.0_pReal + aimag ( curl_fourier ) ** 2.0_pReal ) ! this layer (DC) does not have a conjugate complex counterpart (if grid(1) /= 1)
2019-03-25 22:58:07 +05:30
do l = 1 , 3
2019-03-07 11:31:09 +05:30
curl_fourier = ( + tensorField_fourier ( l , 3 , grid1Red , j , k ) * xi1st ( 2 , grid1Red , j , k ) * rescaledGeom ( 2 ) &
- tensorField_fourier ( l , 2 , grid1Red , j , k ) * xi1st ( 3 , grid1Red , j , k ) * rescaledGeom ( 3 ) )
curl_fourier = ( + tensorField_fourier ( l , 1 , grid1Red , j , k ) * xi1st ( 3 , grid1Red , j , k ) * rescaledGeom ( 3 ) &
- tensorField_fourier ( l , 3 , grid1Red , j , k ) * xi1st ( 1 , grid1Red , j , k ) * rescaledGeom ( 1 ) )
curl_fourier = ( + tensorField_fourier ( l , 2 , grid1Red , j , k ) * xi1st ( 1 , grid1Red , j , k ) * rescaledGeom ( 1 ) &
- tensorField_fourier ( l , 1 , grid1Red , j , k ) * xi1st ( 2 , grid1Red , j , k ) * rescaledGeom ( 2 ) )
enddo
utilities_curlRMS = utilities_curlRMS &
+ sum ( real ( curl_fourier ) ** 2.0_pReal + aimag ( curl_fourier ) ** 2.0_pReal ) ! this layer (Nyquist) does not have a conjugate complex counterpart (if grid(1) /= 1)
enddo ; enddo
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
call MPI_Allreduce ( MPI_IN_PLACE , utilities_curlRMS , 1 , MPI_DOUBLE , MPI_SUM , PETSC_COMM_WORLD , ierr )
2019-03-25 22:58:07 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'utilities_curlRMS' )
2019-03-07 11:31:09 +05:30
utilities_curlRMS = sqrt ( utilities_curlRMS ) * wgt
2019-03-25 22:58:07 +05:30
if ( grid ( 1 ) == 1 ) utilities_curlRMS = utilities_curlRMS * 0.5_pReal ! counted twice in case of grid(1) == 1
2013-03-06 20:01:13 +05:30
end function utilities_curlRMS
2012-07-24 22:37:10 +05:30
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
2013-04-18 22:10:49 +05:30
!> @brief calculates mask compliance tensor used to adjust F to fullfill stress BC
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
2012-10-24 17:01:40 +05:30
function utilities_maskedCompliance ( rot_BC , mask_stress , C )
2019-06-06 21:58:10 +05:30
2019-12-03 00:53:50 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: utilities_maskedCompliance !< masked compliance
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 , 3 , 3 ) :: C !< current average stiffness
type ( rotation ) , intent ( in ) :: rot_BC !< rotation of load frame
logical , intent ( in ) , dimension ( 3 , 3 ) :: mask_stress !< mask of stress BC
2020-02-22 20:18:40 +05:30
integer :: i , j
logical , dimension ( 9 ) :: mask_stressVector
logical , dimension ( 9 , 9 ) :: mask
real ( pReal ) , dimension ( 9 , 9 ) :: temp99_real
2019-03-25 23:47:10 +05:30
integer :: size_reduced = 0
real ( pReal ) , dimension ( : , : ) , allocatable :: &
s_reduced , & !< reduced compliance matrix (depending on number of stress BC)
c_reduced , & !< reduced stiffness (depending on number of stress BC)
sTimesC !< temp variable to check inversion
logical :: errmatinv
2020-01-04 23:31:36 +05:30
character ( len = pStringLen ) :: formatString
2020-02-22 20:18:40 +05:30
2019-03-25 23:47:10 +05:30
mask_stressVector = reshape ( transpose ( mask_stress ) , [ 9 ] )
size_reduced = count ( mask_stressVector )
2020-02-22 20:18:40 +05:30
if ( size_reduced > 0 ) then
temp99_real = math_3333to99 ( rot_BC % rotate ( C ) )
2019-03-25 23:47:10 +05:30
if ( debugGeneral ) then
write ( 6 , '(/,a)' ) ' ... updating masked compliance ............................................'
write ( 6 , '(/,a,/,9(9(2x,f12.7,1x)/))' , advance = 'no' ) ' Stiffness C (load) / GPa =' , &
transpose ( temp99_Real ) * 1.0e-9_pReal
flush ( 6 )
endif
2020-02-22 20:18:40 +05:30
do i = 1 , 9 ; do j = 1 , 9
mask ( i , j ) = mask_stressVector ( i ) . and . mask_stressVector ( j )
enddo ; enddo
c_reduced = reshape ( pack ( temp99_Real , mask ) , [ size_reduced , size_reduced ] )
allocate ( s_reduced , mold = c_reduced )
2019-09-21 06:50:33 +05:30
call math_invert ( s_reduced , errmatinv , c_reduced ) ! invert reduced stiffness
2019-03-25 23:47:10 +05:30
if ( any ( IEEE_is_NaN ( s_reduced ) ) ) errmatinv = . true .
if ( errmatinv ) call IO_error ( error_ID = 400 , ext_msg = 'utilities_maskedCompliance' )
2015-10-15 00:06:19 +05:30
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
2013-03-07 01:04:30 +05:30
! check if inversion was successful
2019-03-25 23:47:10 +05:30
sTimesC = matmul ( c_reduced , s_reduced )
2020-09-13 14:48:09 +05:30
errmatinv = errmatinv . or . any ( dNeq ( sTimesC , math_eye ( size_reduced ) , 1.0e-12_pReal ) )
2019-03-25 23:47:10 +05:30
if ( debugGeneral . or . errmatinv ) then
write ( formatString , '(i2)' ) size_reduced
formatString = '(/,a,/,' / / trim ( formatString ) / / '(' / / trim ( formatString ) / / '(2x,es9.2,1x)/))'
write ( 6 , trim ( formatString ) , advance = 'no' ) ' C * S (load) ' , &
transpose ( matmul ( c_reduced , s_reduced ) )
write ( 6 , trim ( formatString ) , advance = 'no' ) ' S (load) ' , transpose ( s_reduced )
if ( errmatinv ) call IO_error ( error_ID = 400 , ext_msg = 'utilities_maskedCompliance' )
endif
2020-02-22 20:18:40 +05:30
temp99_real = reshape ( unpack ( reshape ( s_reduced , [ size_reduced ** 2 ] ) , reshape ( mask , [ 81 ] ) , 0.0_pReal ) , [ 9 , 9 ] )
2019-03-25 23:47:10 +05:30
else
temp99_real = 0.0_pReal
endif
2020-02-22 20:18:40 +05:30
utilities_maskedCompliance = math_99to3333 ( temp99_Real )
2019-03-25 23:47:10 +05:30
if ( debugGeneral ) then
write ( 6 , '(/,a,/,9(9(2x,f10.5,1x)/),/)' , advance = 'no' ) &
' Masked Compliance (load) * GPa =' , transpose ( temp99_Real ) * 1.0e9_pReal
flush ( 6 )
endif
2012-10-24 17:01:40 +05:30
2015-10-15 00:06:19 +05:30
end function utilities_maskedCompliance
2012-10-24 17:01:40 +05:30
2015-06-03 23:00:31 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculate scalar gradient in fourier field
!--------------------------------------------------------------------------------------------------
subroutine utilities_fourierScalarGradient ( )
2020-02-22 20:18:40 +05:30
2019-03-25 22:58:07 +05:30
integer :: i , j , k
2020-02-22 20:18:40 +05:30
2019-10-29 01:57:57 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
vectorField_fourier ( 1 : 3 , i , j , k ) = scalarField_fourier ( i , j , k ) * xi1st ( 1 : 3 , i , j , k ) ! ToDo: no -conjg?
enddo ; enddo ; enddo
2015-06-03 23:00:31 +05:30
end subroutine utilities_fourierScalarGradient
!--------------------------------------------------------------------------------------------------
!> @brief calculate vector divergence in fourier field
!--------------------------------------------------------------------------------------------------
subroutine utilities_fourierVectorDivergence ( )
2020-02-22 20:18:40 +05:30
2019-10-29 01:57:57 +05:30
integer :: i , j , k
2020-02-22 20:18:40 +05:30
2019-10-29 01:57:57 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
scalarField_fourier ( i , j , k ) = sum ( vectorField_fourier ( 1 : 3 , i , j , k ) * conjg ( - xi1st ( 1 : 3 , i , j , k ) ) )
enddo ; enddo ; enddo
2015-06-03 23:00:31 +05:30
end subroutine utilities_fourierVectorDivergence
2015-09-11 18:35:46 +05:30
2015-12-14 23:42:09 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculate vector gradient in fourier field
!--------------------------------------------------------------------------------------------------
subroutine utilities_fourierVectorGradient ( )
2020-02-22 20:18:40 +05:30
2019-03-25 22:58:07 +05:30
integer :: i , j , k , m , n
2020-02-22 20:18:40 +05:30
2019-03-25 22:58:07 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
do m = 1 , 3 ; do n = 1 , 3
2019-03-07 11:31:09 +05:30
tensorField_fourier ( m , n , i , j , k ) = vectorField_fourier ( m , i , j , k ) * xi1st ( n , i , j , k )
enddo ; enddo
enddo ; enddo ; enddo
2015-12-14 23:42:09 +05:30
end subroutine utilities_fourierVectorGradient
!--------------------------------------------------------------------------------------------------
!> @brief calculate tensor divergence in fourier field
!--------------------------------------------------------------------------------------------------
subroutine utilities_fourierTensorDivergence ( )
2019-06-06 21:58:10 +05:30
2019-10-29 01:57:57 +05:30
integer :: i , j , k
2020-02-22 20:18:40 +05:30
2019-03-25 22:58:07 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
2019-10-29 01:57:57 +05:30
vectorField_fourier ( : , i , j , k ) = matmul ( tensorField_fourier ( : , : , i , j , k ) , conjg ( - xi1st ( : , i , j , k ) ) )
2019-03-07 11:31:09 +05:30
enddo ; enddo ; enddo
2015-12-14 23:42:09 +05:30
end subroutine utilities_fourierTensorDivergence
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
2018-02-16 20:06:18 +05:30
!> @brief calculate constitutive response from materialpoint_F0 to F during timeinc
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
2018-02-16 20:06:18 +05:30
subroutine utilities_constitutiveResponse ( P , P_av , C_volAvg , C_minmaxAvg , &
F , timeinc , rotation_BC )
2020-02-22 20:18:40 +05:30
2019-12-03 00:53:50 +05:30
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 , 3 , 3 ) :: C_volAvg , C_minmaxAvg !< average stiffness
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 ) :: P_av !< average PK stress
real ( pReal ) , intent ( out ) , dimension ( 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) :: P !< PK stress
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) :: F !< deformation gradient target
real ( pReal ) , intent ( in ) :: timeinc !< loading time
type ( rotation ) , intent ( in ) , optional :: rotation_BC !< rotation of load frame
2020-02-22 20:18:40 +05:30
2019-03-25 22:58:07 +05:30
integer :: &
2019-03-07 11:31:09 +05:30
i , ierr
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: dPdF_max , dPdF_min
real ( pReal ) :: dPdF_norm_max , dPdF_norm_min
real ( pReal ) , dimension ( 2 ) :: valueAndRank !< pair of min/max norm of dPdF to synchronize min/max of dPdF
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
write ( 6 , '(/,a)' ) ' ... evaluating constitutive response ......................................'
flush ( 6 )
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
materialpoint_F = reshape ( F , [ 3 , 3 , 1 , product ( grid ( 1 : 2 ) ) * grid3 ] ) ! set materialpoint target F to estimated field
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
call materialpoint_stressAndItsTangent ( . true . , timeinc ) ! calculate P field
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
P = reshape ( materialpoint_P , [ 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ] )
P_av = sum ( sum ( sum ( P , dim = 5 ) , dim = 4 ) , dim = 3 ) * wgt ! average of P
call MPI_Allreduce ( MPI_IN_PLACE , P_av , 9 , MPI_DOUBLE , MPI_SUM , PETSC_COMM_WORLD , ierr )
if ( debugRotation ) &
write ( 6 , '(/,a,/,3(3(2x,f12.4,1x)/))' , advance = 'no' ) ' Piola--Kirchhoff stress (lab) / MPa =' , &
transpose ( P_av ) * 1.e-6_pReal
2019-12-03 00:53:50 +05:30
if ( present ( rotation_BC ) ) &
2020-02-22 20:18:40 +05:30
P_av = rotation_BC % rotate ( P_av )
2019-03-07 11:31:09 +05:30
write ( 6 , '(/,a,/,3(3(2x,f12.4,1x)/))' , advance = 'no' ) ' Piola--Kirchhoff stress / MPa =' , &
transpose ( P_av ) * 1.e-6_pReal
flush ( 6 )
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
dPdF_max = 0.0_pReal
dPdF_norm_max = 0.0_pReal
dPdF_min = huge ( 1.0_pReal )
dPdF_norm_min = huge ( 1.0_pReal )
2019-03-25 22:58:07 +05:30
do i = 1 , product ( grid ( 1 : 2 ) ) * grid3
2019-03-07 11:31:09 +05:30
if ( dPdF_norm_max < sum ( materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , i ) ** 2.0_pReal ) ) then
dPdF_max = materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , i )
dPdF_norm_max = sum ( materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , i ) ** 2.0_pReal )
endif
if ( dPdF_norm_min > sum ( materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , i ) ** 2.0_pReal ) ) then
dPdF_min = materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , i )
dPdF_norm_min = sum ( materialpoint_dPdF ( 1 : 3 , 1 : 3 , 1 : 3 , 1 : 3 , 1 , i ) ** 2.0_pReal )
endif
end do
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
valueAndRank = [ dPdF_norm_max , real ( worldrank , pReal ) ]
call MPI_Allreduce ( MPI_IN_PLACE , valueAndRank , 1 , MPI_2DOUBLE_PRECISION , MPI_MAXLOC , PETSC_COMM_WORLD , ierr )
2019-03-25 22:58:07 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'MPI_Allreduce max' )
2019-03-07 11:31:09 +05:30
call MPI_Bcast ( dPdF_max , 81 , MPI_DOUBLE , int ( valueAndRank ( 2 ) ) , PETSC_COMM_WORLD , ierr )
2019-03-25 22:58:07 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'MPI_Bcast max' )
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
valueAndRank = [ dPdF_norm_min , real ( worldrank , pReal ) ]
call MPI_Allreduce ( MPI_IN_PLACE , valueAndRank , 1 , MPI_2DOUBLE_PRECISION , MPI_MINLOC , PETSC_COMM_WORLD , ierr )
2019-03-25 22:58:07 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'MPI_Allreduce min' )
2019-03-07 11:31:09 +05:30
call MPI_Bcast ( dPdF_min , 81 , MPI_DOUBLE , int ( valueAndRank ( 2 ) ) , PETSC_COMM_WORLD , ierr )
2019-03-25 22:58:07 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'MPI_Bcast min' )
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
C_minmaxAvg = 0.5_pReal * ( dPdF_max + dPdF_min )
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
C_volAvg = sum ( sum ( materialpoint_dPdF , dim = 6 ) , dim = 5 )
call MPI_Allreduce ( MPI_IN_PLACE , C_volAvg , 81 , MPI_DOUBLE , MPI_SUM , PETSC_COMM_WORLD , ierr )
C_volAvg = C_volAvg * wgt
2019-05-12 16:41:30 +05:30
2015-10-15 00:06:19 +05:30
2012-10-24 17:01:40 +05:30
end subroutine utilities_constitutiveResponse
2012-07-30 19:36:22 +05:30
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates forward rate, either guessing or just add delta/timeinc
!--------------------------------------------------------------------------------------------------
2018-02-16 20:06:18 +05:30
pure function utilities_calculateRate ( heterogeneous , field0 , field , dt , avRate )
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 ) :: &
avRate !< homogeneous addon
real ( pReal ) , intent ( in ) :: &
dt !< timeinc between field0 and field
logical , intent ( in ) :: &
heterogeneous !< calculate field of rates
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) :: &
field0 , & !< data of previous step
field !< data of current step
real ( pReal ) , dimension ( 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) :: &
utilities_calculateRate
2020-02-22 20:18:40 +05:30
2019-03-07 11:31:09 +05:30
if ( heterogeneous ) then
utilities_calculateRate = ( field - field0 ) / dt
else
utilities_calculateRate = spread ( spread ( spread ( avRate , 3 , grid ( 1 ) ) , 4 , grid ( 2 ) ) , 5 , grid3 )
endif
2012-08-09 18:34:56 +05:30
2012-10-24 17:01:40 +05:30
end function utilities_calculateRate
2012-10-02 20:56:56 +05:30
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
2015-10-15 00:06:19 +05:30
!> @brief forwards a field with a pointwise given rate, if aim is given,
2013-01-08 15:42:03 +05:30
!> ensures that the average matches the aim
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
2015-03-25 21:36:19 +05:30
function utilities_forwardField ( timeinc , field_lastInc , rate , aim )
2015-10-15 00:06:19 +05:30
2019-03-07 11:42:35 +05:30
real ( pReal ) , intent ( in ) :: &
timeinc !< timeinc of current step
real ( pReal ) , intent ( in ) , dimension ( 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) :: &
field_lastInc , & !< initial field
rate !< rate by which to forward
real ( pReal ) , intent ( in ) , optional , dimension ( 3 , 3 ) :: &
aim !< average field value aim
real ( pReal ) , dimension ( 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) :: &
utilities_forwardField
real ( pReal ) , dimension ( 3 , 3 ) :: fieldDiff !< <a + adot*t> - aim
PetscErrorCode :: ierr
utilities_forwardField = field_lastInc + rate * timeinc
if ( present ( aim ) ) then !< correct to match average
fieldDiff = sum ( sum ( sum ( utilities_forwardField , dim = 5 ) , dim = 4 ) , dim = 3 ) * wgt
call MPI_Allreduce ( MPI_IN_PLACE , fieldDiff , 9 , MPI_DOUBLE , MPI_SUM , PETSC_COMM_WORLD , ierr )
fieldDiff = fieldDiff - aim
utilities_forwardField = utilities_forwardField - &
spread ( spread ( spread ( fieldDiff , 3 , grid ( 1 ) ) , 4 , grid ( 2 ) ) , 5 , grid3 )
endif
2012-10-02 20:56:56 +05:30
2012-10-24 17:01:40 +05:30
end function utilities_forwardField
2012-10-02 20:56:56 +05:30
2012-08-07 22:53:13 +05:30
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculates filter for fourier convolution depending on type given in numerics.config
2016-02-08 22:03:17 +05:30
!> @details this is the full operator to calculate derivatives, i.e. 2 \pi i k for the
! standard approach
2012-10-24 17:01:40 +05:30
!--------------------------------------------------------------------------------------------------
2015-12-14 23:42:09 +05:30
pure function utilities_getFreqDerivative ( k_s )
2019-03-25 22:58:07 +05:30
integer , intent ( in ) , dimension ( 3 ) :: k_s !< indices of frequency
complex ( pReal ) , dimension ( 3 ) :: utilities_getFreqDerivative
2019-03-07 11:42:35 +05:30
select case ( spectral_derivative_ID )
2020-02-22 20:18:40 +05:30
case ( DERIVATIVE_CONTINUOUS_ID )
2019-03-07 11:42:35 +05:30
utilities_getFreqDerivative = cmplx ( 0.0_pReal , 2.0_pReal * PI * real ( k_s , pReal ) / geomSize , pReal )
2020-02-22 20:18:40 +05:30
case ( DERIVATIVE_CENTRAL_DIFF_ID )
2019-03-07 11:42:35 +05:30
utilities_getFreqDerivative = cmplx ( 0.0_pReal , sin ( 2.0_pReal * PI * real ( k_s , pReal ) / real ( grid , pReal ) ) , pReal ) / &
2020-02-22 20:18:40 +05:30
cmplx ( 2.0_pReal * geomSize / real ( grid , pReal ) , 0.0_pReal , pReal )
case ( DERIVATIVE_FWBW_DIFF_ID )
2019-03-07 11:42:35 +05:30
utilities_getFreqDerivative ( 1 ) = &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 1 ) , pReal ) / real ( grid ( 1 ) , pReal ) ) - 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 1 ) , pReal ) / real ( grid ( 1 ) , pReal ) ) , pReal ) * &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 2 ) , pReal ) / real ( grid ( 2 ) , pReal ) ) + 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 2 ) , pReal ) / real ( grid ( 2 ) , pReal ) ) , pReal ) * &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 3 ) , pReal ) / real ( grid ( 3 ) , pReal ) ) + 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 3 ) , pReal ) / real ( grid ( 3 ) , pReal ) ) , pReal ) / &
2020-02-22 20:18:40 +05:30
cmplx ( 4.0_pReal * geomSize ( 1 ) / real ( grid ( 1 ) , pReal ) , 0.0_pReal , pReal )
2019-03-07 11:42:35 +05:30
utilities_getFreqDerivative ( 2 ) = &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 1 ) , pReal ) / real ( grid ( 1 ) , pReal ) ) + 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 1 ) , pReal ) / real ( grid ( 1 ) , pReal ) ) , pReal ) * &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 2 ) , pReal ) / real ( grid ( 2 ) , pReal ) ) - 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 2 ) , pReal ) / real ( grid ( 2 ) , pReal ) ) , pReal ) * &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 3 ) , pReal ) / real ( grid ( 3 ) , pReal ) ) + 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 3 ) , pReal ) / real ( grid ( 3 ) , pReal ) ) , pReal ) / &
2020-02-22 20:18:40 +05:30
cmplx ( 4.0_pReal * geomSize ( 2 ) / real ( grid ( 2 ) , pReal ) , 0.0_pReal , pReal )
2019-03-07 11:42:35 +05:30
utilities_getFreqDerivative ( 3 ) = &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 1 ) , pReal ) / real ( grid ( 1 ) , pReal ) ) + 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 1 ) , pReal ) / real ( grid ( 1 ) , pReal ) ) , pReal ) * &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 2 ) , pReal ) / real ( grid ( 2 ) , pReal ) ) + 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 2 ) , pReal ) / real ( grid ( 2 ) , pReal ) ) , pReal ) * &
cmplx ( cos ( 2.0_pReal * PI * real ( k_s ( 3 ) , pReal ) / real ( grid ( 3 ) , pReal ) ) - 1.0_pReal , &
sin ( 2.0_pReal * PI * real ( k_s ( 3 ) , pReal ) / real ( grid ( 3 ) , pReal ) ) , pReal ) / &
2020-02-22 20:18:40 +05:30
cmplx ( 4.0_pReal * geomSize ( 3 ) / real ( grid ( 3 ) , pReal ) , 0.0_pReal , pReal )
2019-03-07 11:42:35 +05:30
end select
2012-07-25 19:31:39 +05:30
2015-12-14 23:42:09 +05:30
end function utilities_getFreqDerivative
2012-10-02 20:56:56 +05:30
2015-03-13 03:58:33 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief calculate coordinates in current configuration for given defgrad field
! using integration in Fourier space. Similar as in mesh.f90, but using data already defined for
! convolution
!--------------------------------------------------------------------------------------------------
2019-09-28 03:32:36 +05:30
subroutine utilities_updateCoords ( F )
2020-02-22 20:18:40 +05:30
2019-03-07 11:42:35 +05:30
real ( pReal ) , dimension ( 3 , 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) , intent ( in ) :: F
2019-09-28 03:14:28 +05:30
real ( pReal ) , dimension ( 3 , grid ( 1 ) , grid ( 2 ) , grid3 ) :: IPcoords
2019-09-29 22:24:20 +05:30
real ( pReal ) , dimension ( 3 , grid ( 1 ) , grid ( 2 ) , grid3 + 2 ) :: IPfluct_padded ! Fluctuations of cell center displacement (padded along z for MPI)
real ( pReal ) , dimension ( 3 , grid ( 1 ) + 1 , grid ( 2 ) + 1 , grid3 + 1 ) :: nodeCoords
integer :: &
i , j , k , n , &
rank_t , &
rank_b , &
c , r , &
ierr
integer , dimension ( MPI_STATUS_SIZE ) :: &
s
2019-10-29 01:46:43 +05:30
real ( pReal ) , dimension ( 3 ) :: step
2019-03-07 11:42:35 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: Favg
2019-09-29 22:24:20 +05:30
integer , dimension ( 3 ) :: me
integer , dimension ( 3 , 8 ) :: &
neighbor = reshape ( [ &
0 , 0 , 0 , &
1 , 0 , 0 , &
1 , 1 , 0 , &
0 , 1 , 0 , &
0 , 0 , 1 , &
1 , 0 , 1 , &
1 , 1 , 1 , &
0 , 1 , 1 ] , [ 3 , 8 ] )
2020-02-22 20:18:40 +05:30
2019-09-29 23:56:57 +05:30
step = geomSize / real ( grid , pReal )
2019-03-07 11:42:35 +05:30
!--------------------------------------------------------------------------------------------------
2019-09-29 23:56:57 +05:30
! integration in Fourier space to get fluctuations of cell center discplacements
2019-10-29 20:48:58 +05:30
tensorField_real ( 1 : 3 , 1 : 3 , 1 : grid ( 1 ) , 1 : grid ( 2 ) , 1 : grid3 ) = F
2019-03-07 11:42:35 +05:30
call utilities_FFTtensorForward ( )
2019-10-29 13:45:35 +05:30
2019-03-25 22:58:07 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid1Red
2019-10-29 20:48:58 +05:30
if ( any ( [ i , j , k + grid3Offset ] / = 1 ) ) then
vectorField_fourier ( 1 : 3 , i , j , k ) = matmul ( tensorField_fourier ( 1 : 3 , 1 : 3 , i , j , k ) , xi2nd ( 1 : 3 , i , j , k ) ) &
2020-03-09 18:37:31 +05:30
/ sum ( conjg ( - xi2nd ( 1 : 3 , i , j , k ) ) * xi2nd ( 1 : 3 , i , j , k ) ) * cmplx ( wgt , 0.0 , pReal )
2019-10-29 20:48:58 +05:30
else
vectorField_fourier ( 1 : 3 , i , j , k ) = cmplx ( 0.0 , 0.0 , pReal )
endif
2019-03-07 11:42:35 +05:30
enddo ; enddo ; enddo
2019-10-29 20:48:58 +05:30
2019-03-07 11:42:35 +05:30
call fftw_mpi_execute_dft_c2r ( planVectorBack , vectorField_fourier , vectorField_real )
2020-02-22 20:18:40 +05:30
2019-03-07 11:42:35 +05:30
!--------------------------------------------------------------------------------------------------
! average F
2019-03-25 22:58:07 +05:30
if ( grid3Offset == 0 ) Favg = real ( tensorField_fourier ( 1 : 3 , 1 : 3 , 1 , 1 , 1 ) , pReal ) * wgt
2019-03-07 11:42:35 +05:30
call MPI_Bcast ( Favg , 9 , MPI_DOUBLE , 0 , PETSC_COMM_WORLD , ierr )
2019-09-29 23:56:57 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'update_IPcoords/MPI_Bcast' )
2020-02-22 20:18:40 +05:30
2019-09-29 22:24:20 +05:30
!--------------------------------------------------------------------------------------------------
2019-09-29 23:56:57 +05:30
! pad cell center fluctuations along z-direction (needed when running MPI simulation)
2019-09-29 22:24:20 +05:30
IPfluct_padded ( 1 : 3 , 1 : grid ( 1 ) , 1 : grid ( 2 ) , 2 : grid3 + 1 ) = vectorField_real ( 1 : 3 , 1 : grid ( 1 ) , 1 : grid ( 2 ) , 1 : grid3 )
2019-09-29 23:56:57 +05:30
c = product ( shape ( IPfluct_padded ( : , : , : , 1 ) ) ) !< amount of data to transfer
2019-09-29 22:24:20 +05:30
rank_t = modulo ( worldrank + 1 , worldsize )
rank_b = modulo ( worldrank - 1 , worldsize )
2020-02-22 20:18:40 +05:30
2019-09-29 23:56:57 +05:30
! send bottom layer to process below
2019-09-29 22:24:20 +05:30
call MPI_Isend ( IPfluct_padded ( : , : , : , 2 ) , c , MPI_DOUBLE , rank_b , 0 , PETSC_COMM_WORLD , r , ierr )
2019-09-29 23:56:57 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'update_IPcoords/MPI_Isend' )
2019-09-29 22:24:20 +05:30
call MPI_Irecv ( IPfluct_padded ( : , : , : , grid3 + 2 ) , c , MPI_DOUBLE , rank_t , 0 , PETSC_COMM_WORLD , r , ierr )
2019-09-29 23:56:57 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'update_IPcoords/MPI_Irecv' )
2019-09-29 22:24:20 +05:30
call MPI_Wait ( r , s , ierr )
2020-03-29 22:48:06 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'update_IPcoords/MPI_Wait' )
2020-02-22 20:18:40 +05:30
2019-09-29 23:56:57 +05:30
! send top layer to process above
2019-09-29 22:24:20 +05:30
call MPI_Isend ( IPfluct_padded ( : , : , : , grid3 + 1 ) , c , MPI_DOUBLE , rank_t , 0 , PETSC_COMM_WORLD , r , ierr )
2019-09-29 23:56:57 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'update_IPcoords/MPI_Isend' )
2019-09-29 22:24:20 +05:30
call MPI_Irecv ( IPfluct_padded ( : , : , : , 1 ) , c , MPI_DOUBLE , rank_b , 0 , PETSC_COMM_WORLD , r , ierr )
2019-09-29 23:56:57 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'update_IPcoords/MPI_Irecv' )
2019-09-29 22:24:20 +05:30
call MPI_Wait ( r , s , ierr )
2020-03-29 22:48:06 +05:30
if ( ierr / = 0 ) call IO_error ( 894 , ext_msg = 'update_IPcoords/MPI_Wait' )
2020-02-22 20:18:40 +05:30
2019-09-29 23:56:57 +05:30
!--------------------------------------------------------------------------------------------------
2020-02-22 20:18:40 +05:30
! calculate nodal displacements
2019-09-29 22:24:20 +05:30
nodeCoords = 0.0_pReal
do k = 0 , grid3 ; do j = 0 , grid ( 2 ) ; do i = 0 , grid ( 1 )
2019-10-30 19:19:08 +05:30
nodeCoords ( 1 : 3 , i + 1 , j + 1 , k + 1 ) = matmul ( Favg , step * ( real ( [ i , j , k + grid3Offset ] , pReal ) ) )
averageFluct : do n = 1 , 8
2019-09-29 22:24:20 +05:30
me = [ i + neighbor ( 1 , n ) , j + neighbor ( 2 , n ) , k + neighbor ( 3 , n ) ]
2019-09-29 23:56:57 +05:30
nodeCoords ( 1 : 3 , i + 1 , j + 1 , k + 1 ) = nodeCoords ( 1 : 3 , i + 1 , j + 1 , k + 1 ) &
2019-10-30 19:19:08 +05:30
+ IPfluct_padded ( 1 : 3 , modulo ( me ( 1 ) - 1 , grid ( 1 ) ) + 1 , modulo ( me ( 2 ) - 1 , grid ( 2 ) ) + 1 , me ( 3 ) + 1 ) * 0.125_pReal
enddo averageFluct
2019-09-29 22:24:20 +05:30
enddo ; enddo ; enddo
2020-02-22 20:18:40 +05:30
2019-09-29 22:24:20 +05:30
!--------------------------------------------------------------------------------------------------
! calculate cell center displacements
2019-03-25 22:58:07 +05:30
do k = 1 , grid3 ; do j = 1 , grid ( 2 ) ; do i = 1 , grid ( 1 )
2019-09-28 03:14:28 +05:30
IPcoords ( 1 : 3 , i , j , k ) = vectorField_real ( 1 : 3 , i , j , k ) &
2020-03-09 18:37:31 +05:30
+ matmul ( Favg , step * ( real ( [ i , j , k + grid3Offset ] , pReal ) - 0.5_pReal ) )
2019-03-07 11:42:35 +05:30
enddo ; enddo ; enddo
2020-02-22 20:18:40 +05:30
2019-09-29 23:56:57 +05:30
call discretization_setNodeCoords ( reshape ( NodeCoords , [ 3 , ( grid ( 1 ) + 1 ) * ( grid ( 2 ) + 1 ) * ( grid3 + 1 ) ] ) )
call discretization_setIPcoords ( reshape ( IPcoords , [ 3 , grid ( 1 ) * grid ( 2 ) * grid3 ] ) )
2020-02-22 20:18:40 +05:30
2019-09-28 03:32:36 +05:30
end subroutine utilities_updateCoords
2015-03-13 03:58:33 +05:30
2019-10-24 02:20:01 +05:30
!---------------------------------------------------------------------------------------------------
!> @brief Write out the current reference stiffness for restart.
!---------------------------------------------------------------------------------------------------
subroutine utilities_saveReferenceStiffness
2020-02-22 20:18:40 +05:30
2019-10-24 02:20:01 +05:30
integer :: &
2020-09-13 14:58:48 +05:30
fileUnit , ierr
2019-10-24 02:20:01 +05:30
if ( worldrank == 0 ) then
2019-12-07 19:54:45 +05:30
write ( 6 , '(a)' ) ' writing reference stiffness data required for restart to file' ; flush ( 6 )
2020-09-13 14:58:48 +05:30
open ( newunit = fileUnit , file = getSolverJobName ( ) / / '.C_ref' , &
status = 'replace' , access = 'stream' , action = 'write' , iostat = ierr )
if ( ierr / = 0 ) call IO_error ( 100 , ext_msg = 'could not open file ' / / getSolverJobName ( ) / / '.C_ref' )
2019-10-24 16:36:42 +05:30
write ( fileUnit ) C_ref
close ( fileUnit )
2019-10-24 02:20:01 +05:30
endif
2020-02-22 20:18:40 +05:30
2019-10-24 02:20:01 +05:30
end subroutine utilities_saveReferenceStiffness
2015-12-16 02:15:54 +05:30
end module spectral_utilities