2013-03-22 23:05:05 +05:30
! Copyright 2011-13 Max-Planck-Institut für Eisenforschung GmbH
!
! This file is part of DAMASK,
! the Düsseldorf Advanced MAterial Simulation Kit.
!
! 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/>.
!
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
2012-08-03 14:55:48 +05:30
! $Id$
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
!> @author Pratheek Shanthraj, Max-Planck-Institut für Eisenforschung GmbH
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @author Philip Eisenlohr, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Basic scheme solver
2012-08-29 00:49:47 +05:30
!> @details this solver follows closely the original large strain formulation presented by
2013-01-29 15:58:01 +05:30
!> Suquet. The iterative procedure is solved using a fix-point iteration.
2012-07-30 19:36:22 +05:30
!--------------------------------------------------------------------------------------------------
2012-07-24 22:37:10 +05:30
module DAMASK_spectral_SolverBasic
2012-08-03 14:55:48 +05:30
use prec , only : &
pInt , &
pReal
2012-07-30 19:36:22 +05:30
use math , only : &
2012-08-03 14:55:48 +05:30
math_I3
use DAMASK_spectral_Utilities , only : &
2012-10-02 20:56:56 +05:30
tSolutionState
2012-07-25 19:31:39 +05:30
2012-08-03 14:55:48 +05:30
implicit none
2013-03-28 16:07:00 +05:30
private
2012-07-25 19:31:39 +05:30
character ( len = * ) , parameter , public :: &
DAMASK_spectral_SolverBasic_label = 'basic'
!--------------------------------------------------------------------------------------------------
2012-11-12 19:44:39 +05:30
! pointwise global data
real ( pReal ) , private , dimension ( : , : , : , : , : ) , allocatable :: &
F , & !< deformation gradient field
F_lastInc , & !< deformation gradient field last increment
Fdot !< assumed rate for F n to F n+1
2013-01-09 23:38:08 +05:30
2012-07-25 19:31:39 +05:30
!--------------------------------------------------------------------------------------------------
! stress, stiffness and compliance average etc.
2012-08-03 14:55:48 +05:30
real ( pReal ) , private , dimension ( 3 , 3 ) :: &
2012-11-12 19:44:39 +05:30
F_aim = math_I3 , & !< deformation gradient aim
F_aim_lastInc = math_I3 , & !< deformation gradient aim last increment
F_aimDot = 0.0_pReal !< assumed rate
2012-08-03 14:55:48 +05:30
real ( pReal ) , private , dimension ( 3 , 3 , 3 , 3 ) :: &
2013-03-06 20:01:13 +05:30
C = 0.0_pReal , C_minmaxAvg = 0.0_pReal , & !< average stiffness
2012-11-12 19:44:39 +05:30
C_lastInc = 0.0_pReal !< average stiffness last increment
2013-03-28 13:10:30 +05:30
public :: &
basic_init , &
basic_solution , &
basic_destroy
2012-08-29 00:49:47 +05:30
contains
2012-07-24 22:37:10 +05:30
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief allocates all neccessary fields and fills them with data, potentially from restart info
!--------------------------------------------------------------------------------------------------
2013-01-09 23:38:08 +05:30
subroutine basic_init ( temperature )
2012-08-27 13:34:47 +05:30
use , intrinsic :: iso_fortran_env ! to get compiler_version and compiler_options (at least for gfortran >4.6 at the moment)
2012-08-03 14:55:48 +05:30
use IO , only : &
IO_read_JobBinaryFile , &
2012-08-31 01:56:28 +05:30
IO_write_JobBinaryFile , &
2013-02-25 22:04:59 +05:30
IO_intOut , &
IO_timeStamp
2012-08-03 14:55:48 +05:30
use FEsolving , only : &
restartInc
use DAMASK_interface , only : &
getSolverJobName
use DAMASK_spectral_Utilities , only : &
Utilities_init , &
Utilities_constitutiveResponse , &
Utilities_updateGamma , &
2012-10-02 20:56:56 +05:30
debugRestart
2012-08-03 14:55:48 +05:30
use mesh , only : &
res , &
2012-12-14 01:50:04 +05:30
wgt , &
2012-11-14 20:08:10 +05:30
geomdim , &
2013-01-16 16:10:53 +05:30
scaledDim , &
2012-11-14 20:08:10 +05:30
mesh_ipCoordinates , &
mesh_NcpElems , &
mesh_deformedCoordsFFT
2012-08-03 14:55:48 +05:30
implicit none
2013-01-09 23:38:08 +05:30
real ( pReal ) , intent ( inout ) :: &
temperature
2012-11-12 19:44:39 +05:30
real ( pReal ) , dimension ( 3 , 3 , res ( 1 ) , res ( 2 ) , res ( 3 ) ) :: P
real ( pReal ) , dimension ( 3 , 3 ) :: &
2012-12-14 20:48:04 +05:30
temp33_Real = 0.0_pReal
2012-11-12 19:44:39 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
temp3333_Real
2012-07-30 19:36:22 +05:30
2012-08-03 14:55:48 +05:30
call Utilities_Init ( )
2013-03-27 17:58:55 +05:30
write ( 6 , '(/,a)' ) ' <<<+- DAMASK_spectral_solverBasic init -+>>>'
write ( 6 , '(a)' ) ' $Id$'
write ( 6 , '(a15,a)' ) ' Current time: ' , IO_timeStamp ( )
2012-07-30 19:36:22 +05:30
#include "compilation_info.f90"
2013-01-16 16:10:53 +05:30
write ( 6 , '(a,3(f12.5)/)' ) ' scaledDim x y z:' , scaledDim
2013-01-10 03:49:32 +05:30
2012-11-12 19:44:39 +05:30
!--------------------------------------------------------------------------------------------------
! allocate global fields
2012-12-14 20:48:04 +05:30
allocate ( F ( 3 , 3 , res ( 1 ) , res ( 2 ) , res ( 3 ) ) , source = 0.0_pReal )
allocate ( F_lastInc ( 3 , 3 , res ( 1 ) , res ( 2 ) , res ( 3 ) ) , source = 0.0_pReal )
allocate ( Fdot ( 3 , 3 , res ( 1 ) , res ( 2 ) , res ( 3 ) ) , source = 0.0_pReal )
2012-07-25 19:31:39 +05:30
!--------------------------------------------------------------------------------------------------
2012-11-12 19:44:39 +05:30
! init fields and average quantities
2012-08-03 14:55:48 +05:30
if ( restartInc == 1_pInt ) then ! no deformation (no restart)
2012-08-09 18:34:56 +05:30
F = spread ( spread ( spread ( math_I3 , 3 , res ( 1 ) ) , 4 , res ( 2 ) ) , 5 , res ( 3 ) ) ! initialize to identity
F_lastInc = F
2012-08-03 14:55:48 +05:30
elseif ( restartInc > 1_pInt ) then ! using old values from file
2013-01-10 03:49:32 +05:30
if ( debugRestart ) write ( 6 , '(/,a,' / / IO_intOut ( restartInc - 1_pInt ) / / ',a)' ) &
'reading values of increment' , restartInc - 1_pInt , 'from file'
flush ( 6 )
2012-12-14 23:00:22 +05:30
call IO_read_jobBinaryFile ( 777 , 'F' , &
2012-07-25 19:31:39 +05:30
trim ( getSolverJobName ( ) ) , size ( F ) )
2012-08-03 14:55:48 +05:30
read ( 777 , rec = 1 ) F
close ( 777 )
2012-12-14 23:00:22 +05:30
call IO_read_jobBinaryFile ( 777 , 'F_lastInc' , &
2012-07-25 19:31:39 +05:30
trim ( getSolverJobName ( ) ) , size ( F_lastInc ) )
2012-08-03 14:55:48 +05:30
read ( 777 , rec = 1 ) F_lastInc
close ( 777 )
2012-12-14 01:50:04 +05:30
F_aim = sum ( sum ( sum ( F , dim = 5 ) , dim = 4 ) , dim = 3 ) * wgt ! average of F
F_aim_lastInc = sum ( sum ( sum ( F_lastInc , dim = 5 ) , dim = 4 ) , dim = 3 ) * wgt ! average of F_lastInc
2013-01-24 00:03:46 +05:30
call IO_read_jobBinaryFile ( 777 , 'F_aimDot' , trim ( getSolverJobName ( ) ) , size ( f_aimDot ) )
read ( 777 , rec = 1 ) f_aimDot
2012-10-24 17:01:40 +05:30
close ( 777 )
call IO_read_jobBinaryFile ( 777 , 'C' , trim ( getSolverJobName ( ) ) , size ( C ) )
read ( 777 , rec = 1 ) C
close ( 777 )
2013-01-24 00:03:46 +05:30
call IO_read_jobBinaryFile ( 777 , 'C_lastInc' , trim ( getSolverJobName ( ) ) , size ( C_lastInc ) )
read ( 777 , rec = 1 ) C_lastInc
2012-10-24 17:01:40 +05:30
close ( 777 )
call IO_read_jobBinaryFile ( 777 , 'C_ref' , trim ( getSolverJobName ( ) ) , size ( temp3333_Real ) )
read ( 777 , rec = 1 ) temp3333_Real
close ( 777 )
2012-08-03 14:55:48 +05:30
endif
2013-02-05 18:01:44 +05:30
mesh_ipCoordinates = reshape ( mesh_deformedCoordsFFT ( geomdim , F ) , [ 3 , 1 , mesh_NcpElems ] )
2013-03-06 20:01:13 +05:30
call Utilities_constitutiveResponse ( F , F , temperature , 0.0_pReal , P , C , C_minmaxAvg , temp33_Real , . false . , math_I3 ) ! constitutive response with no deformation in no time to get reference stiffness
2012-10-24 17:01:40 +05:30
if ( restartInc == 1_pInt ) then ! use initial stiffness as reference stiffness
2013-03-06 20:01:13 +05:30
temp3333_Real = C_minmaxAvg
2012-10-24 17:01:40 +05:30
endif
2012-07-25 19:31:39 +05:30
2012-10-24 17:01:40 +05:30
call Utilities_updateGamma ( temp3333_Real , . True . )
2012-07-24 22:37:10 +05:30
2012-08-03 14:55:48 +05:30
end subroutine basic_init
!--------------------------------------------------------------------------------------------------
!> @brief solution for the basic scheme with internal iterations
!--------------------------------------------------------------------------------------------------
2012-10-02 20:56:56 +05:30
type ( tSolutionState ) function &
2012-11-09 01:02:00 +05:30
basic_solution ( incInfo , guess , timeinc , timeinc_old , P_BC , F_BC , temperature_bc , rotation_BC )
2012-07-24 22:37:10 +05:30
use numerics , only : &
2012-07-26 19:28:47 +05:30
itmax , &
itmin , &
update_gamma
2012-08-03 14:55:48 +05:30
use math , only : &
math_mul33x33 , &
math_rotate_backward33 , &
math_transpose33 , &
2012-08-27 13:34:47 +05:30
math_mul3333xx33
2012-08-03 14:55:48 +05:30
use mesh , only : &
res , &
2012-08-27 13:34:47 +05:30
geomdim , &
2012-11-14 20:08:10 +05:30
wgt , &
mesh_ipCoordinates , &
mesh_NcpElems , &
mesh_deformedCoordsFFT
2012-07-24 22:37:10 +05:30
use IO , only : &
2012-08-31 01:56:28 +05:30
IO_write_JobBinaryFile , &
IO_intOut
2012-08-03 14:55:48 +05:30
use DAMASK_spectral_Utilities , only : &
2012-10-02 20:56:56 +05:30
tBoundaryCondition , &
2012-08-03 14:55:48 +05:30
field_real , &
Utilities_forwardField , &
Utilities_maskedCompliance , &
2012-10-02 20:56:56 +05:30
Utilities_FFTforward , &
2012-08-03 14:55:48 +05:30
Utilities_divergenceRMS , &
Utilities_fourierConvolution , &
2012-10-02 20:56:56 +05:30
Utilities_FFTbackward , &
2012-08-03 14:55:48 +05:30
Utilities_updateGamma , &
2012-10-02 20:56:56 +05:30
Utilities_constitutiveResponse , &
2012-11-28 20:34:05 +05:30
Utilities_calculateRate , &
debugRotation
2012-07-24 22:37:10 +05:30
use FEsolving , only : &
2012-08-29 11:20:42 +05:30
restartWrite , &
2012-10-24 17:01:40 +05:30
restartRead , &
2012-08-29 11:20:42 +05:30
terminallyIll
2012-10-24 17:01:40 +05:30
use DAMASK_spectral_Utilities , only : &
cutBack
2012-07-24 22:37:10 +05:30
implicit none
!--------------------------------------------------------------------------------------------------
! input data for solution
2012-11-12 19:44:39 +05:30
real ( pReal ) , intent ( in ) :: &
timeinc , & !< increment in time for current solution
2013-01-09 23:38:08 +05:30
timeinc_old !< increment in time of last increment
2012-11-12 19:44:39 +05:30
logical , intent ( in ) :: &
guess !< if .false., assume homogeneous addon
type ( tBoundaryCondition ) , intent ( in ) :: &
P_BC , & !< stress boundary conditions
F_BC !< deformation boundary conditions
character ( len = * ) , intent ( in ) :: &
incInfo !< string with information of current increment for output to screen
real ( pReal ) , dimension ( 3 , 3 ) , intent ( in ) :: &
rotation_BC !< rotation load to lab
2013-01-09 23:38:08 +05:30
real ( pReal ) , intent ( inout ) :: &
temperature_bc
2012-11-12 19:44:39 +05:30
real ( pReal ) , dimension ( 3 , 3 , 3 , 3 ) :: &
S !< current average compliance
2012-11-28 20:34:05 +05:30
real ( pReal ) , dimension ( 3 , 3 ) :: &
F_aim_lastIter , & !< aim of last iteration
P_av
2012-11-12 19:44:39 +05:30
real ( pReal ) , dimension ( 3 , 3 , res ( 1 ) , res ( 2 ) , res ( 3 ) ) :: P
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
! loop variables, convergence etc.
2012-07-25 19:31:39 +05:30
real ( pReal ) :: err_div , err_stress
2013-03-13 00:18:28 +05:30
integer ( pInt ) :: iter , row , column
2012-07-25 19:31:39 +05:30
logical :: ForwardData
2012-07-24 22:37:10 +05:30
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
2012-11-12 19:44:39 +05:30
! write restart information for spectral solver
2012-07-24 22:37:10 +05:30
if ( restartWrite ) then
2013-01-10 03:49:32 +05:30
write ( 6 , '(/,a)' ) ' writing converged results for restart'
flush ( 6 )
2012-12-14 23:00:22 +05:30
call IO_write_jobBinaryFile ( 777 , 'F' , size ( F ) ) ! writing deformation gradient field to file
2012-10-24 17:01:40 +05:30
write ( 777 , rec = 1 ) F
close ( 777 )
2012-12-14 23:00:22 +05:30
call IO_write_jobBinaryFile ( 777 , 'F_lastInc' , size ( F_lastInc ) ) ! writing F_lastInc field to file
2012-10-24 17:01:40 +05:30
write ( 777 , rec = 1 ) F_lastInc
2012-07-24 22:37:10 +05:30
close ( 777 )
2013-01-18 17:00:52 +05:30
call IO_write_jobBinaryFile ( 777 , 'F_aimDot' , size ( f_aimDot ) )
write ( 777 , rec = 1 ) f_aimDot
close ( 777 )
2012-07-24 22:37:10 +05:30
call IO_write_jobBinaryFile ( 777 , 'C' , size ( C ) )
write ( 777 , rec = 1 ) C
close ( 777 )
2012-10-24 17:01:40 +05:30
call IO_write_jobBinaryFile ( 777 , 'C_lastInc' , size ( C_lastInc ) )
write ( 777 , rec = 1 ) C_lastInc
close ( 777 )
2012-07-24 22:37:10 +05:30
endif
2012-11-12 19:44:39 +05:30
!--------------------------------------------------------------------------------------------------
! forward data
if ( cutBack ) then
2012-10-24 17:01:40 +05:30
F_aim = F_aim_lastInc
F = F_lastInc
C = C_lastInc
else
2012-11-12 19:44:39 +05:30
C_lastInc = C
2013-02-05 18:01:44 +05:30
mesh_ipCoordinates = reshape ( mesh_deformedCoordsFFT ( geomdim , F ) , [ 3 , 1 , mesh_NcpElems ] )
2012-12-16 05:22:06 +05:30
!--------------------------------------------------------------------------------------------------
! calculate rate for aim
2012-10-02 20:56:56 +05:30
if ( F_BC % myType == 'l' ) then ! calculate f_aimDot from given L and current F
f_aimDot = F_BC % maskFloat * math_mul33x33 ( F_BC % values , F_aim )
elseif ( F_BC % myType == 'fdot' ) then ! f_aimDot is prescribed
f_aimDot = F_BC % maskFloat * F_BC % values
endif
2012-11-09 01:02:00 +05:30
if ( guess ) f_aimDot = f_aimDot + P_BC % maskFloat * ( F_aim - F_aim_lastInc ) / timeinc_old
2012-10-02 20:56:56 +05:30
F_aim_lastInc = F_aim
2012-12-16 05:22:06 +05:30
!--------------------------------------------------------------------------------------------------
! update rate and forward last inc
2012-11-28 20:34:05 +05:30
Fdot = utilities_calculateRate ( math_rotate_backward33 ( f_aimDot , rotation_BC ) , &
timeinc_old , guess , F_lastInc , F )
2012-10-02 20:56:56 +05:30
F_lastInc = F
endif
F_aim = F_aim + f_aimDot * timeinc
2013-01-08 15:42:03 +05:30
F = Utilities_forwardField ( timeinc , F_lastInc , Fdot , math_rotate_backward33 ( F_aim , rotation_BC ) )
2012-10-19 14:14:21 +05:30
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
! update stiffness (and gamma operator)
S = Utilities_maskedCompliance ( rotation_BC , P_BC % maskLogical , C )
2013-03-06 20:01:13 +05:30
if ( update_gamma ) call Utilities_updateGamma ( C_minmaxAvg , restartWrite )
2012-08-03 14:55:48 +05:30
2012-11-12 19:44:39 +05:30
!--------------------------------------------------------------------------------------------------
! iteration till converged
2012-10-24 17:01:40 +05:30
if ( . not . restartRead ) ForwardData = . True .
2012-08-03 14:55:48 +05:30
iter = 0_pInt
2013-01-11 16:10:16 +05:30
basic_solution % iterationsNeeded = itmax
2012-08-03 14:55:48 +05:30
convergenceLoop : do while ( iter < itmax )
2012-07-24 22:37:10 +05:30
iter = iter + 1_pInt
!--------------------------------------------------------------------------------------------------
! report begin of new iteration
2013-01-10 03:49:32 +05:30
write ( 6 , '(1x,a,3(a,' / / IO_intOut ( itmax ) / / '))' ) trim ( incInfo ) , &
' @ Iteration ' , itmin , '≤' , iter , '≤' , itmax
2012-11-28 20:34:05 +05:30
if ( debugRotation ) &
2013-01-10 03:49:32 +05:30
write ( 6 , '(/,a,/,3(3(f12.7,1x)/))' , advance = 'no' ) ' deformation gradient aim (lab)=' , &
2012-11-28 20:34:05 +05:30
math_transpose33 ( math_rotate_backward33 ( F_aim , rotation_BC ) )
2013-01-10 03:49:32 +05:30
write ( 6 , '(/,a,/,3(3(f12.7,1x)/))' , advance = 'no' ) ' deformation gradient aim =' , &
2012-11-28 20:34:05 +05:30
math_transpose33 ( F_aim )
2013-01-10 03:49:32 +05:30
flush ( 6 )
2012-07-24 22:37:10 +05:30
!--------------------------------------------------------------------------------------------------
! evaluate constitutive response
2012-11-28 20:34:05 +05:30
F_aim_lastIter = F_aim
2012-10-02 20:56:56 +05:30
basic_solution % termIll = . false .
2013-01-09 23:38:08 +05:30
call Utilities_constitutiveResponse ( F_lastInc , F , temperature_bc , timeinc , &
2013-03-06 20:01:13 +05:30
P , C , C_minmaxAvg , P_av , ForwardData , rotation_BC )
2012-08-29 11:20:42 +05:30
basic_solution % termIll = terminallyIll
2012-11-12 19:44:39 +05:30
terminallyIll = . false .
ForwardData = . false .
2012-07-25 19:31:39 +05:30
2012-07-24 22:37:10 +05:30
!--------------------------------------------------------------------------------------------------
! stress BC handling
2012-10-02 20:56:56 +05:30
F_aim = F_aim - math_mul3333xx33 ( S , ( ( P_av - P_BC % values ) ) ) ! S = 0.0 for no bc
err_stress = maxval ( abs ( P_BC % maskFloat * ( P_av - P_BC % values ) ) ) ! mask = 0.0 for no bc
2012-11-28 20:34:05 +05:30
2012-07-24 22:37:10 +05:30
!--------------------------------------------------------------------------------------------------
2012-08-03 14:55:48 +05:30
! updated deformation gradient using fix point algorithm of basic scheme
field_real = 0.0_pReal
2012-08-09 18:34:56 +05:30
field_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) = reshape ( P , [ res ( 1 ) , res ( 2 ) , res ( 3 ) , 3 , 3 ] , &
order = [ 4 , 5 , 1 , 2 , 3 ] ) ! field real has a different order
2012-10-02 20:56:56 +05:30
call Utilities_FFTforward ( )
2012-08-03 14:55:48 +05:30
err_div = Utilities_divergenceRMS ( )
2012-11-28 20:34:05 +05:30
call Utilities_fourierConvolution ( math_rotate_backward33 ( F_aim_lastIter - F_aim , rotation_BC ) )
2012-10-02 20:56:56 +05:30
call Utilities_FFTbackward ( )
2012-08-09 18:34:56 +05:30
F = F - reshape ( field_real ( 1 : res ( 1 ) , 1 : res ( 2 ) , 1 : res ( 3 ) , 1 : 3 , 1 : 3 ) , shape ( F ) , order = [ 3 , 4 , 5 , 1 , 2 ] ) ! F(x)^(n+1) = F(x)^(n) + correction; *wgt: correcting for missing normalization
2012-08-06 14:23:12 +05:30
basic_solution % converged = basic_Converged ( err_div , P_av , err_stress , P_av )
2013-01-10 03:49:32 +05:30
write ( 6 , '(/,a)' ) ' =========================================================================='
flush ( 6 )
2012-10-19 14:14:21 +05:30
if ( ( basic_solution % converged . and . iter > = itmin ) . or . basic_solution % termIll ) then
basic_solution % iterationsNeeded = iter
exit
endif
2012-07-24 22:37:10 +05:30
enddo convergenceLoop
end function basic_solution
2012-08-03 14:55:48 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief convergence check for basic scheme based on div of P and deviation from stress aim
!--------------------------------------------------------------------------------------------------
2012-08-06 14:23:12 +05:30
logical function basic_Converged ( err_div , pAvgDiv , err_stress , pAvgStress )
2012-07-30 19:36:22 +05:30
use numerics , only : &
2012-08-03 14:55:48 +05:30
itmin , &
err_div_tol , &
err_stress_tolrel , &
err_stress_tolabs
use math , only : &
math_mul33x33 , &
math_eigenvalues33 , &
math_transpose33
2012-07-26 19:28:47 +05:30
2012-07-30 19:36:22 +05:30
implicit none
2012-08-03 14:55:48 +05:30
real ( pReal ) , dimension ( 3 , 3 ) , intent ( in ) :: &
pAvgDiv , &
pAvgStress
real ( pReal ) , intent ( in ) :: &
err_div , &
err_stress
real ( pReal ) :: &
err_stress_tol , &
pAvgDivL2
2012-07-26 19:28:47 +05:30
2012-08-03 14:55:48 +05:30
pAvgDivL2 = sqrt ( maxval ( math_eigenvalues33 ( math_mul33x33 ( pAvgDiv , math_transpose33 ( pAvgDiv ) ) ) ) ) ! L_2 norm of average stress (http://mathworld.wolfram.com/SpectralNorm.html)
2013-03-22 20:16:55 +05:30
err_stress_tol = max ( maxval ( abs ( pAvgStress ) ) * err_stress_tolrel , err_stress_tolabs )
2012-08-03 14:55:48 +05:30
2012-08-06 14:23:12 +05:30
basic_Converged = all ( [ err_div / pAvgDivL2 / err_div_tol , &
2012-08-03 14:55:48 +05:30
err_stress / err_stress_tol ] < 1.0_pReal )
2012-07-26 19:28:47 +05:30
2013-03-13 00:18:28 +05:30
write ( 6 , '(/,a,f10.2,a,es11.5,a,es11.4,a)' ) ' error divergence = ' , &
2013-01-10 03:49:32 +05:30
err_div / pAvgDivL2 / err_div_tol , ' (' , err_div / pAvgDivL2 , ' / m, tol =' , err_div_tol , ')'
write ( 6 , '(a,f8.2,a,es11.5,a,es11.4,a)' ) ' error stress BC = ' , &
err_stress / err_stress_tol , ' (' , err_stress , ' Pa , tol =' , err_stress_tol , ')'
flush ( 6 )
2012-07-26 19:28:47 +05:30
2012-08-06 14:23:12 +05:30
end function basic_Converged
2012-07-25 19:31:39 +05:30
2012-10-02 20:56:56 +05:30
2012-12-14 20:48:04 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief does the cleaning up after the simulation has finished
!--------------------------------------------------------------------------------------------------
2012-08-03 14:55:48 +05:30
subroutine basic_destroy ( )
2012-11-12 19:44:39 +05:30
use DAMASK_spectral_Utilities , only : &
Utilities_destroy
2012-08-03 14:55:48 +05:30
implicit none
call Utilities_destroy ( )
2012-07-25 19:31:39 +05:30
end subroutine basic_destroy
2012-10-02 20:56:56 +05:30
2012-07-25 19:31:39 +05:30
end module DAMASK_spectral_SolverBasic