00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052 #ifndef _AMR_H_
00053 #define _AMR_H_
00054
00055 #include <iostream>
00056 #include <string>
00057 #include <ctime>
00058
00059 #include "REAL.H"
00060 #include "Vector.H"
00061 #include "AMRLevel.H"
00062 #include "AMRLevelFactory.H"
00063 #include "BRMeshRefine.H"
00064 #include "ProblemDomain.H"
00065 #include "Box.H"
00066 #include "CH_HDF5.H"
00067 #ifdef CH_USE_TIMER
00068 #include "CH_Timer.H"
00069 #endif
00070
00072
00093 class AMR
00094 {
00095 public:
00097
00101 AMR();
00102
00104
00107 ~AMR();
00108
00110
00130 void define(int a_max_level,
00131 const Vector<int>& a_ref_ratios,
00132 const Box& a_prob_domain,
00133 const AMRLevelFactory* const a_amrLevelFact);
00134
00136
00156 void define(int a_max_level,
00157 const Vector<int>& a_ref_ratios,
00158 const ProblemDomain& a_prob_domain,
00159 const AMRLevelFactory* const a_amrLevelFact);
00160
00162
00166 bool isDefined() const;
00167
00169
00173 bool isSetUp() const;
00174
00175 #ifdef CH_USE_HDF5
00176
00177
00183 void setupForRestart(HDF5Handle& a_handle);
00184 #endif
00185
00187
00193 void setupForNewAMRRun();
00194
00200 void setupForFixedHierarchyRun(const Vector<Vector<Box> >& a_amr_grids,
00201 int a_proper_nest = 1);
00202
00204
00209 void conclude() const;
00210
00212
00216 void run(Real a_max_time, int a_max_step);
00217
00219
00224 void plotPrefix(const std::string& a_plotfile_prefix);
00225
00227
00232 void checkpointPrefix(const std::string& a_checkpointfile_prefix);
00233
00235
00239 void plotInterval(int a_plot_interval);
00240
00242
00246 void checkpointInterval(int a_checkpoint_interval);
00247
00249
00253 void maxGridSize(int a_max_grid_size);
00254
00256
00261 void maxBaseGridSize(int a_max_base_grid_size);
00262
00264
00269 void dtToleranceFactor(Real a_dt_tolerance_factor);
00270
00272
00276 void fillRatio(Real a_fillRat);
00277
00279
00283 void blockFactor(int a_blockFactor);
00284
00286
00290 void gridBufferSize(int a_grid_buffer_size);
00291
00293
00296 int maxGridSize() const;
00297
00299
00302 int maxBaseGridSize() const;
00303
00305
00316 void verbosity (int a_verbosity);
00317
00319
00323 void regridIntervals(const Vector<int>& a_regridIntervals);
00324
00326
00335 int verbosity () const;
00336
00338
00339 void maxDtGrow(Real a_dtGrowFactor);
00340
00342
00347 void timeEps(Real a_timeEps);
00348
00350
00351 Real maxDtGrow() const;
00352
00354
00357 Real timeEps() const;
00358
00360
00361 void fixedDt(Real a_dt);
00362
00364
00365 Real fixedDt() const;
00366
00368
00371 Vector<AMRLevel*> getAMRLevels();
00372
00374
00380 void initialTime(Real a_initialTime);
00381
00383
00386 Real getCurrentTime() const;
00387
00388 #ifdef CH_USE_TIMER
00389
00390
00393 Chombo::Timer * timer(Chombo::Timer *a_timer = NULL );
00394 #endif
00395
00396 protected:
00397
00398
00399 int timeStep(int a_level, int a_stepsLeft, bool a_coarseTimeBoundary);
00400
00401
00402 void clearMemory();
00403
00404
00405 void regrid(int a_base_level);
00406
00407
00408
00409 void initialGrid();
00410
00411 void writePlotFile() const;
00412
00413 void writeCheckpointFile() const;
00414
00415
00416
00417 void assignDt();
00418
00419
00420 bool needToRegrid(int a_level, int a_numStepsLeft) const;
00421
00422 void makeBaseLevelMesh (Vector<Box>& a_grids) const;
00423
00424 void setDefaultValues();
00425
00426 int m_blockFactor;
00427 Real m_fillRatio;
00428 int m_max_level;
00429 int m_finest_level_old;
00430 int m_finest_level;
00431 int m_checkpoint_interval;
00432 int m_plot_interval;
00433 int m_max_grid_size;
00434 int m_max_base_grid_size;
00435 Real m_dt_tolerance_factor;
00436 Real m_fixedDt;
00437
00438 bool m_isDefined;
00439 bool m_isSetUp;
00440 Vector<AMRLevel*> m_amrlevels;
00441 Vector<int> m_ref_ratios;
00442 Vector<int> m_reduction_factor;
00443 Vector<int> m_regrid_intervals;
00444
00445 Real m_dt_base;
00446
00447 Vector<Real> m_dt_new;
00448
00449 Vector<Real> m_dt_cur;
00450 Real m_maxDtGrow;
00451 Real m_time_eps;
00452
00453 BRMeshRefine m_mesh_refine;
00454 bool m_use_meshrefine;
00455
00456 Vector<Vector<Box> > m_amr_grids;
00457
00458 int m_cur_step;
00459 int m_restart_step;
00460 int m_lastcheck_step;
00461
00462 Real m_cur_time;
00463 Vector<int> m_steps_since_regrid;
00464 Vector<long> m_cell_updates;
00465
00466 std::string m_plotfile_prefix;
00467 std::string m_checkpointfile_prefix;
00468
00469 int m_verbosity;
00470 #ifdef CH_USE_TIMER
00471 Chombo::Timer *m_timer;
00472 #endif
00473 };
00474
00475 #endif