diff options
author | hawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b> | 2003-05-21 09:12:14 +0000 |
---|---|---|
committer | hawke <hawke@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b> | 2003-05-21 09:12:14 +0000 |
commit | f3a2829b5b32b6fad8383ac16c16ac6c2f58b696 (patch) | |
tree | 19a58318d2bcb7f95ecf044e5713453d9bff0bea /src/RK2.c | |
parent | 7e9dff2f70bd9950b27092e816f9a3cc14ff3582 (diff) |
Add support for evolving complex GFs and (real and complex) GAs.
Only works with ICN or RK2 for now - in fact this commit may break the generic RK methods temporarily.
Note the documentation isn't quite right - there's no longer a seperate function for each different type...
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@12 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
Diffstat (limited to 'src/RK2.c')
-rw-r--r-- | src/RK2.c | 100 |
1 files changed, 100 insertions, 0 deletions
@@ -67,12 +67,22 @@ void MoL_RK2Add(CCTK_ARGUMENTS) DECLARE_CCTK_ARGUMENTS; DECLARE_CCTK_PARAMETERS; + + cGroupDynamicData arraydata; + CCTK_INT groupindex, ierr; + CCTK_INT arraytotalsize, arraydim; CCTK_INT index, var; CCTK_INT totalsize; CCTK_REAL *OldVar; CCTK_REAL *UpdateVar; CCTK_REAL *RHSVar; + + CCTK_COMPLEX *OldComplexVar; + CCTK_COMPLEX *UpdateComplexVar; + CCTK_COMPLEX *RHSComplexVar; + CCTK_COMPLEX Complex_Delta_Time = CCTK_Cmplx(CCTK_DELTA_TIME, 0); + CCTK_COMPLEX Complex_Half = CCTK_Cmplx(0.5, 0); #ifdef MOLDEBUG printf("Inside RK2.\nStep %d.\nRefinement %d.\nTimestep %g.\nSpacestep %g.\nTime %g\n", @@ -102,6 +112,48 @@ void MoL_RK2Add(CCTK_ARGUMENTS) UpdateVar[index] += CCTK_DELTA_TIME * RHSVar[index]; } } + + for (var = 0; var < MoLNumEvolvedArrayVariables; var++) + { + UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0, + EvolvedArrayVariableIndex[var]); + RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0, + RHSArrayVariableIndex[var]); + groupindex = CCTK_GroupIndexFromVarI(EvolvedArrayVariableIndex[var]); + ierr = CCTK_GroupDynamicData(cctkGH, groupindex, + &arraydata); + if (ierr) + { + CCTK_VWarn(0, __LINE__, __FILE__, "MoL", + "The driver does not return group information for group '%s'.", + CCTK_GroupName(groupindex)); + } + arraytotalsize = 1; + for (arraydim = 0; arraydim < arraydata.dim; arraydim++) + { + arraytotalsize *= arraydata.lsh[arraydim]; + } + + for (index = 0; index < arraytotalsize; index++) + { + UpdateVar[index] += CCTK_DELTA_TIME * RHSVar[index]; + } + } + + for (var = 0; var < MoLNumEvolvedComplexVariables; var++) + { + UpdateComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0, + EvolvedComplexVariableIndex[var]); + RHSComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0, + RHSComplexVariableIndex[var]); + + for (index = 0; index < totalsize; index++) + { + UpdateComplexVar[index] = CCTK_CmplxAdd(UpdateComplexVar[index], + CCTK_CmplxMul(Complex_Delta_Time, + RHSComplexVar[index])); + } + } break; } case 1: @@ -121,6 +173,54 @@ void MoL_RK2Add(CCTK_ARGUMENTS) CCTK_DELTA_TIME * RHSVar[index]; } } + + for (var = 0; var < MoLNumEvolvedArrayVariables; var++) + { + OldVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 1, + EvolvedVariableIndex[var]); + UpdateVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0, + EvolvedVariableIndex[var]); + RHSVar = (CCTK_REAL*)CCTK_VarDataPtrI(cctkGH, 0, + RHSVariableIndex[var]); + + groupindex = CCTK_GroupIndexFromVarI(EvolvedArrayVariableIndex[var]); + ierr = CCTK_GroupDynamicData(cctkGH, groupindex, + &arraydata); + if (ierr) + { + CCTK_VWarn(0, __LINE__, __FILE__, "MoL", + "The driver does not return group information for group '%s'.", + CCTK_GroupName(groupindex)); + } + arraytotalsize = 1; + for (arraydim = 0; arraydim < arraydata.dim; arraydim++) + { + arraytotalsize *= arraydata.lsh[arraydim]; + } + + for (index = 0; index < arraytotalsize; index++) + { + UpdateVar[index] = 0.5 * (OldVar[index] + UpdateVar[index]) + + CCTK_DELTA_TIME * RHSVar[index]; + } + } + + for (var = 0; var < MoLNumEvolvedComplexVariables; var++) + { + OldComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 1, + EvolvedComplexVariableIndex[var]); + UpdateComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0, + EvolvedComplexVariableIndex[var]); + RHSComplexVar = (CCTK_COMPLEX*)CCTK_VarDataPtrI(cctkGH, 0, + RHSComplexVariableIndex[var]); + + for (index = 0; index < totalsize; index++) + { + UpdateComplexVar[index] = + CCTK_CmplxAdd(CCTK_CmplxMul(Complex_Half, + (CCTK_CmplxAdd(OldComplexVar[index], UpdateComplexVar[index]))), CCTK_CmplxMul(Complex_Delta_Time, RHSComplexVar[index])); + } + } break; } default: |