aboutsummaryrefslogtreecommitdiff
path: root/Carpet
diff options
context:
space:
mode:
authorschnetter <>2004-08-07 18:07:00 +0000
committerschnetter <>2004-08-07 18:07:00 +0000
commit69b18dea57e95ec720b259d1cc78db3cf0d8d2d4 (patch)
tree8cd237c355055054bae8e7cb2560a957298a1dc2 /Carpet
parentf1c0efba773d01fa53f2f48a9d16b0f5c3a62d37 (diff)
Code cleanup:
Code cleanup: Rename variables that shadow other variables. Properly cast all size_t to int when comparing. darcs-hash:20040807180727-07bb3-283dc6358053e4cdbd88dff0760b7f0b772701a5.gz
Diffstat (limited to 'Carpet')
-rw-r--r--Carpet/Carpet/src/Evolve.cc12
-rw-r--r--Carpet/Carpet/src/Recompose.cc42
-rw-r--r--Carpet/Carpet/src/SetupGH.cc17
-rw-r--r--Carpet/Carpet/src/modes.cc11
4 files changed, 41 insertions, 41 deletions
diff --git a/Carpet/Carpet/src/Evolve.cc b/Carpet/Carpet/src/Evolve.cc
index 4bdcffdc1..2121f08b7 100644
--- a/Carpet/Carpet/src/Evolve.cc
+++ b/Carpet/Carpet/src/Evolve.cc
@@ -31,7 +31,7 @@
#include "carpet.hh"
extern "C" {
- static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Evolve.cc,v 1.52 2004/08/07 19:47:12 schnetter Exp $";
+ static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Evolve.cc,v 1.53 2004/08/07 20:07:27 schnetter Exp $";
CCTK_FILEVERSION(Carpet_Carpet_Evolve_cc);
}
@@ -284,11 +284,11 @@ namespace Carpet {
int finest_active_reflevel = -1;
{
- for (int rl=0; rl<reflevels; ++rl) {
- const int do_every
- = ipow(mgfact, ml) * (maxreflevelfact / ipow(reffact, rl));
- if (cgh->cctk_iteration % do_every == 0) {
- finest_active_reflevel = rl;
+ for (int rl_=0; rl_<reflevels; ++rl_) {
+ const int do_every_
+ = ipow(mgfact, ml) * (maxreflevelfact / ipow(reffact, rl_));
+ if (cgh->cctk_iteration % do_every_ == 0) {
+ finest_active_reflevel = rl_;
}
}
assert (finest_active_reflevel >= 0);
diff --git a/Carpet/Carpet/src/Recompose.cc b/Carpet/Carpet/src/Recompose.cc
index a6e9c301a..8a649fd2c 100644
--- a/Carpet/Carpet/src/Recompose.cc
+++ b/Carpet/Carpet/src/Recompose.cc
@@ -27,7 +27,7 @@
#include "modes.hh"
extern "C" {
- static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Recompose.cc,v 1.70 2004/08/07 19:47:11 schnetter Exp $";
+ static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/Recompose.cc,v 1.71 2004/08/07 20:07:27 schnetter Exp $";
CCTK_FILEVERSION(Carpet_Carpet_Recompose_cc);
}
@@ -289,11 +289,11 @@ namespace Carpet {
file << "iteration " << cgh->cctk_iteration << endl;
file << "maps " << maps << endl;
file << m << " reflevels " << bbsss.size() << endl;
- for (size_t rl=0; rl<bbsss.size(); ++rl) {
+ for (int rl=0; rl<(int)bbsss.size(); ++rl) {
file << m << " " << rl << " components " << bbsss.at(rl).size() << endl;
- for (size_t c=0; c<bbsss.at(rl).size(); ++c) {
+ for (int c=0; c<(int)bbsss.at(rl).size(); ++c) {
file << m << " " << rl << " " << c << " mglevels " << bbsss.at(rl).at(c).size() << endl;
- for (size_t ml=0; ml<bbsss.at(rl).at(c).size(); ++ml) {
+ for (int ml=0; ml<(int)bbsss.at(rl).at(c).size(); ++ml) {
file << m << " " << rl << " " << c << " " << ml << " " << pss.at(rl).at(c) << " " << bbsss.at(rl).at(c).at(ml) << obss.at(rl).at(c) << endl;
}
}
@@ -386,7 +386,7 @@ namespace Carpet {
if (c<nprocs-1) obs.at(c)[dir][1] = false;
}
- for (size_t n=0; n<ps.size(); ++n) {
+ for (int n=0; n<(int)ps.size(); ++n) {
assert (ps.at(n) == n);
}
}
@@ -468,12 +468,12 @@ namespace Carpet {
// store
bbs.insert (bbs.end(), nprocs, newbb);
obs.insert (obs.end(), nprocs, newob);
- for (int p=0; p<nprocs; ++p) ps.insert (ps.end(), 1, newp+p);
+ for (int pp=0; pp<nprocs; ++pp) ps.insert (ps.end(), 1, newp+pp);
// check postconditions
- assert (bbs.size() == nprocs);
- assert (obs.size() == nprocs);
- assert (ps.size() == nprocs);
+ assert ((int)bbs.size() == nprocs);
+ assert ((int)obs.size() == nprocs);
+ assert ((int)ps.size() == nprocs);
if (DEBUG) cout << "SRAR exit" << endl;
return;
}
@@ -566,20 +566,20 @@ namespace Carpet {
if (DEBUG) cout << "SRAR " << mydim << " newps " << newps << endl;
// store
- assert (newbbs.size() == mynprocs.at(n));
- assert (newobs.size() == mynprocs.at(n));
- assert (newps.size() == mynprocs.at(n));
+ assert ((int)newbbs.size() == mynprocs.at(n));
+ assert ((int)newobs.size() == mynprocs.at(n));
+ assert ((int)newps.size() == mynprocs.at(n));
bbs.insert (bbs.end(), newbbs.begin(), newbbs.end());
obs.insert (obs.end(), newobs.begin(), newobs.end());
ps.insert (ps.end(), newps.begin(), newps.end());
}
// check postconditions
- assert (bbs.size() == nprocs);
- assert (obs.size() == nprocs);
- assert (ps.size() == nprocs);
- for (size_t n=0; n<ps.size(); ++n) {
- assert (ps.at(n) == p+n);
+ assert ((int)bbs.size() == nprocs);
+ assert ((int)obs.size() == nprocs);
+ assert ((int)ps.size() == nprocs);
+ for (int n=0; n<(int)ps.size(); ++n) {
+ assert ((int)ps.at(n) == p+n);
}
if (DEBUG) cout << "SRAR exit" << endl;
}
@@ -695,7 +695,7 @@ namespace Carpet {
bbs = allbbs;
obs = allobs;
ps = allps;
- for (size_t n=0; n<ps.size(); ++n) {
+ for (int n=0; n<(int)ps.size(); ++n) {
ps.at(n) /= ncomps;
assert (ps.at(n) >= 0 && ps.at(n) < nprocs);
}
@@ -783,7 +783,7 @@ namespace Carpet {
}
}
- for (size_t n=0; n<ps.size(); ++n) {
+ for (int n=0; n<(int)ps.size(); ++n) {
assert (ps.at(n) == n);
}
}
@@ -848,11 +848,11 @@ namespace Carpet {
{
assert (bbs.size() == obs.size());
ibbox base;
- for (size_t c=0; c<bbs.size(); ++c) {
+ for (int c=0; c<(int)bbs.size(); ++c) {
base = base.expanded_containing(bbs.at(c));
}
bbss.resize(bbs.size());
- for (size_t c=0; c<bbs.size(); ++c) {
+ for (int c=0; c<(int)bbs.size(); ++c) {
MakeMultigridBoxes (cgh, base, bbs.at(c), obs.at(c), bbss.at(c));
}
}
diff --git a/Carpet/Carpet/src/SetupGH.cc b/Carpet/Carpet/src/SetupGH.cc
index 8bccd5015..d19ab3dc3 100644
--- a/Carpet/Carpet/src/SetupGH.cc
+++ b/Carpet/Carpet/src/SetupGH.cc
@@ -24,7 +24,7 @@
#include "carpet.hh"
extern "C" {
- static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/SetupGH.cc,v 1.84 2004/08/07 19:47:11 schnetter Exp $";
+ static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/SetupGH.cc,v 1.85 2004/08/07 20:07:27 schnetter Exp $";
CCTK_FILEVERSION(Carpet_Carpet_SetupGH_cc);
}
@@ -144,6 +144,8 @@ namespace Carpet {
{
assert (group>=0 && group<CCTK_NumGroups());
+ int ierr;
+
if (CCTK_GroupTypeI(group) != CCTK_GF) {
// Ignore everything but true grid functions
return op_error;
@@ -152,7 +154,7 @@ namespace Carpet {
const bool can_transfer = CanTransferVariableType (cgh, group);
cGroup gp;
- const int ierr = CCTK_GroupData (group, &gp);
+ ierr = CCTK_GroupData (group, &gp);
assert (!ierr);
// Get prolongation method
@@ -199,8 +201,7 @@ namespace Carpet {
if (have_prolong_param_string) {
char * thorn;
char * name;
- int const ierr
- = CCTK_DecomposeName (prolong_param_string, &thorn, &name);
+ ierr = CCTK_DecomposeName (prolong_param_string, &thorn, &name);
if (ierr < 0) {
char * const groupname = CCTK_GroupName (group);
CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
@@ -714,7 +715,7 @@ namespace Carpet {
groupname, basemglevel, convergence_factor);
free (groupname);
}
- if (any(abs(sizes - real_sizes) > 0.001)) {
+ if (any(abs(rvect(sizes) - real_sizes) > 0.001)) {
char * const groupname = CCTK_GroupName(group);
CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
"The shape of group \"%s\" scaled for convergence level %d with convergence factor %d is not integer",
@@ -902,7 +903,7 @@ namespace Carpet {
for (int group=0; group<CCTK_NumGroups(); ++group) {
cGroup data;
- const int ierr = CCTK_GroupData (group, &data);
+ ierr = CCTK_GroupData (group, &data);
assert (!ierr);
switch (data.grouptype) {
@@ -927,10 +928,10 @@ namespace Carpet {
CCTK_VInfo (CCTK_THORNSTRING,
" There are %d grid scalars in %d groups",
num_array_vars.at(0), num_array_groups.at(0));
- for (int dim=1; dim<=3; ++dim) {
+ for (int d=1; d<=3; ++d) {
CCTK_VInfo (CCTK_THORNSTRING,
" There are %d %d-dimensional grid arrays in %d groups",
- num_array_vars.at(dim), dim, num_array_groups.at(dim));
+ num_array_vars.at(d), d, num_array_groups.at(d));
}
CCTK_VInfo (CCTK_THORNSTRING,
" (The number of variables counts all time levels)");
diff --git a/Carpet/Carpet/src/modes.cc b/Carpet/Carpet/src/modes.cc
index c9f6b9bbf..e831b2dd8 100644
--- a/Carpet/Carpet/src/modes.cc
+++ b/Carpet/Carpet/src/modes.cc
@@ -12,7 +12,7 @@
#include "carpet.hh"
extern "C" {
- static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/modes.cc,v 1.8 2004/08/03 19:39:07 schnetter Exp $";
+ static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/Carpet/src/modes.cc,v 1.9 2004/08/07 20:07:27 schnetter Exp $";
CCTK_FILEVERSION(Carpet_Carpet_modes_cc);
}
@@ -90,7 +90,6 @@ namespace Carpet {
for (int group=0; group<CCTK_NumGroups(); ++group) {
if (CCTK_GroupTypeI(group) != CCTK_GF) {
- const int ml = mglevel;
const int rl = 0;
const int m = 0;
const int c = CCTK_MyProc(cgh);
@@ -238,8 +237,8 @@ namespace Carpet {
cgh->cctk_timefac = reflevelfact;
// Set current time
- assert (mglevel>=0 && mglevel<leveltimes.size());
- assert (reflevel>=0 && reflevel<leveltimes.at(mglevel).size());
+ assert (mglevel>=0 && mglevel<(int)leveltimes.size());
+ assert (reflevel>=0 && reflevel<(int)leveltimes.at(mglevel).size());
cgh->cctk_time = leveltimes.at(mglevel).at(reflevel);
assert (is_level_mode());
@@ -253,8 +252,8 @@ namespace Carpet {
if (reflevel == -1) return; // early return
// Save and unset current time
- assert (mglevel>=0 && mglevel<leveltimes.size());
- assert (reflevel>=0 && reflevel<leveltimes.at(mglevel).size());
+ assert (mglevel>=0 && mglevel<(int)leveltimes.size());
+ assert (reflevel>=0 && reflevel<(int)leveltimes.at(mglevel).size());
leveltimes.at(mglevel).at(reflevel) = cgh->cctk_time;
cgh->cctk_time = global_time;