diff options
author | diener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a> | 2007-10-16 14:45:23 +0000 |
---|---|---|
committer | diener <diener@f69c4107-0314-4c4f-9ad4-17e986b73f4a> | 2007-10-16 14:45:23 +0000 |
commit | f4cc81a52f102c61b07178b88d4a560b307d928f (patch) | |
tree | ae8d9b3dcbf84f6c2180622f0e7da3661d4fdcd4 /src/Coefficients_6_5.F90 | |
parent | bb568eb19331f243f79abc3a5b3e9331959c8773 (diff) |
Moved definitions of coefficients to a single module in order to only have
them in one place. Added routines to return coefficients in all the
missing cases. Also added routines to return 2 derivatives (only 2-1 and
4-2 so far). The routines calculating derivatives for whole grid functions
have not been tested yet, so if you are using them, you will probably want
to hold off on updating until I give the clear signal.
git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/SummationByParts/trunk@85 f69c4107-0314-4c4f-9ad4-17e986b73f4a
Diffstat (limited to 'src/Coefficients_6_5.F90')
-rw-r--r-- | src/Coefficients_6_5.F90 | 76 |
1 files changed, 3 insertions, 73 deletions
diff --git a/src/Coefficients_6_5.F90 b/src/Coefficients_6_5.F90 index 361ab5b..22677f7 100644 --- a/src/Coefficients_6_5.F90 +++ b/src/Coefficients_6_5.F90 @@ -4,13 +4,13 @@ subroutine set_coeff_6_5 ( nsize, loc_order, bb, gsize, imin, imax, dd ) + use All_Coeffs_mod + implicit none DECLARE_CCTK_FUNCTIONS DECLARE_CCTK_PARAMETERS - CCTK_REAL, parameter :: zero = 0.0 - integer, parameter :: wp = kind(zero) CCTK_INT, intent(IN) :: nsize, loc_order CCTK_INT, dimension(2), intent(IN) :: bb CCTK_INT, intent(IN) :: gsize @@ -25,77 +25,7 @@ subroutine set_coeff_6_5 ( nsize, loc_order, bb, gsize, imin, imax, dd ) logical, save :: first = .true. if ( first ) then - a(1) = 3.0_wp/4.0_wp; a(2) = -3.0_wp/20.0_wp; a(3) = 1.0_wp/60.0_wp - - q(1,1) = -137.0_wp/60.0_wp; q(2,1) = 5.0_wp; q(3,1) = -5.0_wp; - q(4,1) = 10.0_wp/3.0_wp; q(5,1) = -5.0_wp/4.0_wp; q(6,1) = 1.0_wp/5.0_wp; - q(7,1) = zero; q(8,1) = zero; q(9,1) = zero; q(10,1) = zero; - - q(1,2) = -0.3209098489859668118347_wp - q(2,2) = -0.359441670715467075894_wp - q(3,2) = 0.195756852998105503889_wp - q(4,2) = 1.394685510250317033169_wp - q(5,2) = -1.448965775497476572821_wp - q(6,2) = 0.6519476244467816674835_wp - q(7,2) = -0.1115052611983591304249_wp - q(8,2) = -0.00156743129793461356831_wp - q(9,2) = zero; q(10,2) = zero - - q(1,3) = 0.0980796259621024870281_wp - q(2,3) = -0.784752101846956387607_wp - q(3,3) = 0.364206221878666751300_wp - q(4,3) = 0.102097753636344358379_wp - q(5,3) = 0.377167650934576412896_wp - q(6,3) = -0.173241400242683302020_wp - q(7,3) = 0.0062120424243610783650_wp - q(8,3) = 0.01153111791917461507898_wp - q(9,3) = -0.00130091066558601341956_wp - q(10,3) = zero - - q(1,4) = 0.0771947815148234674634_wp - q(2,4) = -0.4124719432622107102831_wp - q(3,4) = 0.654798717867622763906_wp - q(4,4) = -1.873474657556531216053_wp - q(5,4) = 2.161961221935152129759_wp - q(6,4) = -0.729147814207596325277_wp - q(7,4) = 0.1292509053479867634674_wp - q(8,4) = -0.01092898521375837831427_wp - q(9,4) = 0.00316983426828354261009_wp - q(10,4) = -0.000352060693772037278516_wp - - q(1,5) = -0.01135836398682271847110_wp - q(2,5) = 0.0114974985688122284189_wp - q(3,5) = 0.228765565489867047381_wp - q(4,5) = -1.179122428698534075609_wp - q(5,5) = 0.774618306056840861215_wp - q(6,5) = 0.016612403904600041087_wp - q(7,5) = 0.2399304273448163041220_wp - q(8,5) = -0.0959180672407791626759_wp - q(9,5) = 0.01612460763754670343830_wp - q(10,5) = -0.001149949076347228906108_wp - - q(1,6) = -0.01912046962023373131259_wp - q(2,6) = 0.1569245474161812198493_wp - q(3,6) = -0.567250839387340928718_wp - q(4,6) = 1.230411310316952671353_wp - q(5,6) = -2.048795717924254691221_wp - q(6,6) = 0.982715833107347313532_wp - q(7,6) = 0.2876584784911885215053_wp - q(8,6) = -0.0207657038198929292487_wp - q(9,6) = -0.00335290957957120811713_wp - q(10,6) = 0.001575470999623762377395_wp - - q(1,7) = -0.0076072322555736921424_wp - q(2,7) = 0.034503414345314189270_wp - q(3,7) = -0.044377233426986476366_wp - q(4,7) = -0.049408595803628143797_wp - q(5,7) = 0.304693537653175477283_wp - q(6,7) = -0.935863321236636595622_wp - q(7,7) = 0.111375967229143365899_wp - q(8,7) = 0.7149322477536284813995_wp - q(9,7) = -0.1444768161656941519630_wp - q(10,7) = 0.01622803190725754603783_wp - + call coeffs_1_6_5 ( a, q ) first = .false. end if |