Overture
Version 25
Main Page
Namespaces
Classes
Files
File List
File Members
Overture.v25.d
include
Regrid.h
Go to the documentation of this file.
1
#ifndef REGRID_H
2
#define REGRID_H
3
4
#include "
Overture.h
"
5
#include "BoxLib.H"
6
#include "BoxList.H"
7
#include "
LoadBalancer.h
"
8
9
class
RotatedBox
;
10
class
ListOfRotatedBox
;
11
class
GenericGraphicsInterface
;
12
13
class
Regrid
14
{
15
public
:
16
17
enum
GridAdditionOption
18
{
19
addGridsAsRefinementGrids
,
20
addGridsAsBaseGrids
21
};
22
23
enum
GridAlgorithmOption
24
{
25
aligned
,
26
rotated
27
};
28
29
Regrid
();
30
~Regrid
();
31
32
int
getDefaultNumberOfRefinementLevels
()
const
;
33
int
getRefinementRatio
()
const
;
34
35
int
displayParameters
(FILE *file = stdout )
const
;
36
37
bool
loadBalancingIsOn
()
const
;
// return true is load balancing is turned on
38
39
// regrid based on an error mask
40
int
regrid
(
GridCollection
& gc,
// grid to regrid.
41
GridCollection
& gcNew,
// put new grid here (must be different from gc)
42
intGridCollectionFunction
& errorMask,
// =1 at points to refine
43
int
refinementLevel = 1,
// highest level to refine
44
int
baseLevel = -1 );
// keep this level and below fixed, by default baseLevel=refinementLevel-1.
45
46
// regrid based on an error function and error tolerance
47
int
regrid
(
GridCollection
& gc,
// grid to regrid.
48
GridCollection
& gcNew,
// put new grid here (must be different from gc)
49
realGridCollectionFunction
& error,
50
real
errorThreshhold,
51
int
refinementLevel = 1,
// highest level to refine
52
int
baseLevel = -1 );
// keep this level and below fixed, by default baseLevel=refinementLevel-1.
53
54
static
int
outputRefinementInfo
(
GridCollection
& gc,
55
const
aString
& gridFileName,
56
const
aString
& fileName );
57
58
int
printStatistics
(
GridCollection
& gc, FILE *file =
NULL
,
59
int
*numberOfGridPoints=
NULL
);
60
61
void
setEfficiency
(
real
efficiency
);
// gridding efficiency, 0 < efficiency < 1
62
63
void
setIndexCoarseningFactor
(
int
factor);
64
65
void
setNumberOfBufferZones
(
int
numberOfBufferZones
);
// exapansion of tagged error points.
66
67
void
setMinimumBoxSize
(
int
numberOfGridPoints);
68
69
void
setMinimumBoxWidth
(
int
numberOfGridPoints);
70
71
void
setWidthOfProperNesting
(
int
widthOfProperNesting
);
// distance between levels
72
73
void
setRefinementRatio
(
int
refinementRatio
);
74
75
void
setUseSmartBisection
(
bool
trueOrFalse=
true
);
76
77
void
setGridAdditionOption
(
GridAdditionOption
gridAdditionOption
);
78
GridAdditionOption
getGridAdditionOption
()
const
;
79
80
void
setGridAlgorithmOption
(
GridAlgorithmOption
gridAlgorithmOption
);
81
82
void
setMaximumNumberOfSplits
(
int
num );
83
84
void
setMergeBoxes
(
bool
trueOrFalse=
true
);
// allow boxes to be merged?
85
86
void
turnOnLoadBalacing
(
bool
trueOrFalse=
true
);
87
88
LoadBalancer
&
getLoadBalancer
();
// here is the Loadbalancer used by Regrid
89
90
int
update
(
GenericGraphicsInterface
& gi );
91
92
int
get
(
const
GenericDataBase
& dir,
const
aString
& name);
93
94
int
put
(
GenericDataBase
& dir,
const
aString
& name)
const
;
95
96
int
debug
;
97
98
protected
:
99
100
enum
CutStatus
101
{
102
invalidCut
,
103
holeCut
,
104
steepCut
,
105
bisectCut
106
};
107
108
// int addRefinementsAsBaseGrids(GridCollection & gc, int level0, int numberOfRefinementLevels0,
109
// IntegerArray **gridInfo );
110
111
int
buildProperNestingDomains
(
GridCollection
& gc,
112
int
baseGrid,
113
int
refinementLevel,
114
int
baseLevel,
115
int
numberOfRefinementLevels );
116
117
int
buildTaggedCells
(
MappedGrid
& mg,
118
intMappedGridFunction
& tag,
119
const
realArray
& error,
120
real
errorThreshhold,
121
bool
useErrorMask,
122
bool
cellCentred =
true
);
123
124
Box
cellCenteredBox
(
MappedGrid
& mg,
int
ratio=1 );
125
Box
cellCenteredBaseBox
(
MappedGrid
& mg );
126
127
int
findCut
(
int
*hist,
int
lo,
int
hi,
CutStatus
&status);
128
int
findCutPoint
( BOX & box,
const
intSerialArray & ia,
int
& cutDirection,
int
& cutPoint );
129
// int fixPeriodicBox( MappedGrid & mg, BOX & mainBox, const intArray & ia, int level );
130
131
BOX
getBox
(
const
intArray & ia );
132
BOX
buildBox
(Index Iv[3] );
133
BOX
getBoundedBox
(
const
intSerialArray & ia,
const
Box & boundingBox );
134
135
#ifdef USE_PPP
136
BOX
getBox
(
const
intSerialArray & ia );
137
#endif
138
139
real
getEfficiency
(
const
intSerialArray & ia,
const
BOX & box );
140
141
int
buildGrids
(
GridCollection
& gc,
142
GridCollection
& gcNew,
143
int
baseGrid,
int
baseLevel,
int
refinementLevel, BoxList *refinementBoxList,
144
IntegerArray
**gridInfo);
145
146
int
regridAligned
(
GridCollection
& gc,
147
GridCollection
& gcNew,
148
bool
useErrorFunction,
149
realGridCollectionFunction
*pError,
150
real
errorThreshhold,
151
intGridCollectionFunction
& tagCollection,
152
int
refinementLevel = 1,
153
int
baseLevel = -1 );
154
155
int
regridRotated
(
GridCollection
& gc,
156
GridCollection
& gcNew,
157
bool
useErrorFunction,
158
realGridCollectionFunction
*pError,
159
real
errorThreshhold,
160
intGridCollectionFunction
& tagCollection,
161
int
refinementLevel = 1,
162
int
baseLevel = -1 );
163
164
int
splitBox
( BOX & box,
const
intSerialArray & ia, BoxList & boxList,
int
refinementLevel );
165
int
splitBoxRotated
(
RotatedBox
& box,
ListOfRotatedBox
& boxList,
166
realArray
& xa,
int
refinementLevel );
167
168
int
merge
(
ListOfRotatedBox
& boxList );
169
170
inline
int
coarsenIndexLower
(
int
i,
int
dir )
const
;
171
inline
int
coarsenIndexUpper
(
int
i,
int
dir )
const
;
172
inline
int
refineIndex
(
int
i,
int
dir )
const
;
173
174
void
setupCoarseIndexSpace
(
GridCollection
& gc,
int
baseGrid,
int
level );
175
176
int
defaultNumberOfRefinementLevels
;
177
real
efficiency
;
178
int
refinementRatio
;
179
180
int
numberOfBufferZones
;
// increase tagged cells by this many zones.
181
int
widthOfProperNesting
;
// number of cells between grids at level l with those at l-1
182
183
bool
useCoarsenedIndexSpace
;
// set to true if we are using a coarsened index-space
184
int
indexCoarseningFactor
;
// build amr grids on an index space that is coarsened by this amount.
185
int
piab
[6];
// holds index bounds for use with the indexCoarseningFactor
186
187
int
maximumNumberOfSplits
;
// limit max no of splits for testing
188
int
splitNumber
;
189
int
numberOfDimensions
;
190
191
int
minimumBoxSize
;
192
int
minimumBoxWidth
;
193
bool
useSmartBisection
;
194
bool
mergeBoxes
;
195
196
int
myid
;
197
198
bool
loadBalance
;
199
LoadBalancer
loadBalancer
;
200
201
GridAdditionOption
gridAdditionOption
;
202
GridAlgorithmOption
gridAlgorithmOption
;
203
204
BoxList *
properNestingDomain
;
205
BoxList *
complementOfProperNestingDomain
;
206
207
real
timeForRegrid
;
208
real
timeForBuildGrids
;
209
real
timeForBuildTaggedCells
;
210
real
timeForSomethingElse1
;
211
real
timeForSomethingElse2
;
212
real
timeForSomethingElse3
;
213
214
};
215
216
217
#endif
Generated on Fri Jan 4 2013 10:17:58 for Overture by
1.8.3