aboutsummaryrefslogtreecommitdiff
path: root/src/setup_epsdis.c
diff options
context:
space:
mode:
authorknarf <knarf@850bcc8b-0e4f-0410-8c26-8d28fbf1eda9>2005-06-15 12:05:57 +0000
committerknarf <knarf@850bcc8b-0e4f-0410-8c26-8d28fbf1eda9>2005-06-15 12:05:57 +0000
commit63962de4f7a67068229d280924d74a2c8e58a1a4 (patch)
tree35712726bc19e42e63610bfc0a7965ab52cb73cf /src/setup_epsdis.c
parent006083638b7fcad16e5e2c6eed71f72f2c1a2c9c (diff)
remote the factor epsdisp from the slopes. This means that you have to change
the slope values in your old parameter files by dividing by the slope value used to get the old behaviour. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/Dissipation/trunk@11 850bcc8b-0e4f-0410-8c26-8d28fbf1eda9
Diffstat (limited to 'src/setup_epsdis.c')
-rw-r--r--src/setup_epsdis.c14
1 files changed, 7 insertions, 7 deletions
diff --git a/src/setup_epsdis.c b/src/setup_epsdis.c
index 9a3b2cf..5acb649 100644
--- a/src/setup_epsdis.c
+++ b/src/setup_epsdis.c
@@ -66,7 +66,7 @@ setup_epsdis (CCTK_ARGUMENTS)
for (i=0;i<outer_bound_npoints;i++) {
index = CCTK_GFINDEX3D(cctkGH,i,j,k);
indexP= CCTK_GFINDEX3D(cctkGH,outer_bound_npoints,j,k);
- epsdisA[index]=epsdis+epsdis*ob_slope*abs(x[index]-x[indexP]);
+ epsdisA[index]=epsdis+ob_slope*abs(x[index]-x[indexP]);
}
}
}
@@ -77,7 +77,7 @@ setup_epsdis (CCTK_ARGUMENTS)
for (i=ni;i>ni-outer_bound_npoints;i--) {
index = CCTK_GFINDEX3D(cctkGH,i,j,k);
indexP= CCTK_GFINDEX3D(cctkGH,ni-outer_bound_npoints,j,k);
- epsdisA[index]=epsdis+epsdis*ob_slope*abs(x[index]-x[indexP]);
+ epsdisA[index]=epsdis+ob_slope*abs(x[index]-x[indexP]);
}
}
}
@@ -89,7 +89,7 @@ setup_epsdis (CCTK_ARGUMENTS)
for (i=0;i<ni;i++) {
index = CCTK_GFINDEX3D(cctkGH,i,j,k);
indexP= CCTK_GFINDEX3D(cctkGH,i,outer_bound_npoints,k);
- epsdisA[index]=epsdis+epsdis*ob_slope*abs(y[index]-y[indexP]);
+ epsdisA[index]=epsdis+ob_slope*abs(y[index]-y[indexP]);
}
}
}
@@ -100,7 +100,7 @@ setup_epsdis (CCTK_ARGUMENTS)
for (i=0;i<ni;i++) {
index = CCTK_GFINDEX3D(cctkGH,i,j,k);
indexP= CCTK_GFINDEX3D(cctkGH,i,nj-outer_bound_npoints,k);
- epsdisA[index]=epsdis+epsdis*ob_slope*abs(y[index]-y[indexP]);
+ epsdisA[index]=epsdis+ob_slope*abs(y[index]-y[indexP]);
}
}
}
@@ -112,7 +112,7 @@ setup_epsdis (CCTK_ARGUMENTS)
for (i=0;i<ni;i++) {
index = CCTK_GFINDEX3D(cctkGH,i,j,k);
indexP= CCTK_GFINDEX3D(cctkGH,i,j,outer_bound_npoints);
- epsdisA[index]=epsdis+epsdis*ob_slope*abs(z[index]-z[indexP]);
+ epsdisA[index]=epsdis+ob_slope*abs(z[index]-z[indexP]);
}
}
}
@@ -123,7 +123,7 @@ setup_epsdis (CCTK_ARGUMENTS)
for (i=0;i<ni;i++) {
index = CCTK_GFINDEX3D(cctkGH,i,j,k);
indexP= CCTK_GFINDEX3D(cctkGH,i,j,nk-outer_bound_npoints);
- epsdisA[index]=epsdis+epsdis*ob_slope*abs(z[index]-z[indexP]);
+ epsdisA[index]=epsdis+ob_slope*abs(z[index]-z[indexP]);
}
}
}
@@ -219,7 +219,7 @@ setup_epsdis (CCTK_ARGUMENTS)
for (i=0;i<npts;i++) {
radp=sqrt(xa[i]*xa[i]+ya[i]*ya[i]+za[i]*za[i]);
if (radp<=rads[i]) {
- epsdisA[inds[i]]=epsdis+ ah_slope*(rads[i]-radp)*epsdis;
+ epsdisA[inds[i]]=epsdis+ ah_slope*(rads[i]-radp);
}
}