PLUTO
AMRLevelPluto.cpp
Go to the documentation of this file.
1 #ifdef CH_LANG_CC
2 /*
3  * _______ __
4  * / ___/ / ___ __ _ / / ___
5  * / /__/ _ \/ _ \/ V \/ _ \/ _ \
6  * \___/_//_/\___/_/_/_/_.__/\___/
7  * Please refer to Copyright.txt, in Chombo's root directory.
8  */
9 #endif
10 
11 #include <iomanip>
12 
13 #include "parstream.H"
14 #include "ParmParse.H"
15 #include "LayoutIterator.H"
16 #include "CH_HDF5.H"
17 #include "SPMD.H"
18 #include "LoadBalance.H"
19 #include "BoxIterator.H"
20 #include "computeSum.H"
21 #include "AMRIO.H"
22 #include "AMRLevel.H"
23 
24 #include "AMRLevelPluto.H"
25 
26 #include "NamespaceHeader.H"
27 
28 // Constructor
29 AMRLevelPluto::AMRLevelPluto()
30 {
31  if (s_verbosity >= 3)
32  {
33  pout() << "AMRLevelPluto default constructor" << endl;
34  }
35 
36  m_patchPlutoFactory = NULL;
37  m_patchPluto = NULL;
38  m_paramsDefined = false;
39 }
40 
41 // Destructor
42 AMRLevelPluto::~AMRLevelPluto()
43 {
44  if (s_verbosity >= 3)
45  {
46  pout() << "AMRLevelPluto destructor" << endl;
47  }
48 
49  // Get rid of the patch integrator and its factory
50  if (m_patchPluto != NULL)
51  {
52  delete m_patchPluto;
53  }
54 
55  if (m_patchPlutoFactory != NULL)
56  {
57  delete m_patchPlutoFactory;
58  }
59 
60  m_paramsDefined = false;
61 }
62 
63 void AMRLevelPluto::defineParams(const Real& a_cfl,
64  const Real& a_domainLength,
65  const int& a_verbosity,
66  const Real& a_refineThresh,
67  const int& a_tagBufferSize,
68  const Real& a_initialDtMultiplier,
69  const PatchPluto* const a_patchPluto)
70 {
71  // Store the CFL number
72  m_cfl = a_cfl;
73 
74  // Store the physical dimension of the longest side of the domain
75  m_domainLength = a_domainLength;
76 
77  // Store the verbosity of the object
78  verbosity(a_verbosity);
79 
80  // Store the refinement threshold for gradient
81  m_refineThresh = a_refineThresh;
82 
83  // Store the tag buffer size
84  m_tagBufferSize = a_tagBufferSize;
85 
86  // Store the initial dt multiplier
87  initialDtMultiplier(a_initialDtMultiplier);
88 
89  // Delete any existing patch factory object
90  // (Factory has to be kept for the Patch integrator in LevelPluto ...
91  // ... go figure ...)
92  if (m_patchPlutoFactory != NULL)
93  {
94  delete m_patchPlutoFactory;
95  }
96 
97  m_patchPlutoFactory = a_patchPluto->new_patchPluto();
98 
99  // Delete any existing patch object
100  if (m_patchPlutoFactory != NULL)
101  {
102  delete m_patchPluto;
103  }
104 
105  m_patchPluto = a_patchPluto->new_patchPluto();
106 
107  m_paramsDefined = true;
108 }
109 
110 // This instance should never get called - historical
111 void AMRLevelPluto::define(AMRLevel* a_coarserLevelPtr,
112  const Box& a_problemDomain,
113  int a_level,
114  int a_refRatio)
115 {
116  ProblemDomain physdomain(a_problemDomain);
117 
118  MayDay::Error("AMRLevelPluto::define -\n\tShould never be called with a Box for a problem domain");
119 }
120 
121 // Define new AMR level
122 void AMRLevelPluto::define(AMRLevel* a_coarserLevelPtr,
123  const ProblemDomain& a_problemDomain,
124  int a_level,
125  int a_refRatio)
126 {
127  if (s_verbosity >= 3)
128  {
129  pout() << "AMRLevelPluto::define " << a_level << endl;
130  }
131 
132  // Call inherited define
133  AMRLevel::define(a_coarserLevelPtr,
134  a_problemDomain,
135  a_level,
136  a_refRatio);
137 
138  // Get setup information from the next coarser level
139  if (a_coarserLevelPtr != NULL)
140  {
141  AMRLevelPluto* amrGodPtr = dynamic_cast<AMRLevelPluto*>(a_coarserLevelPtr);
142 
143  if (amrGodPtr != NULL)
144  {
145  m_cfl = amrGodPtr->m_cfl;
146  m_domainLength = amrGodPtr->m_domainLength;
147  m_refineThresh = amrGodPtr->m_refineThresh;
148  m_tagBufferSize = amrGodPtr->m_tagBufferSize;
149  }
150  else
151  {
152  MayDay::Error("AMRLevelPluto::define: a_coarserLevelPtr is not castable to AMRLevelPluto*");
153  }
154  }
155 
156  // Compute the grid spacing
157  m_dx = m_domainLength / (a_problemDomain.domainBox().bigEnd(0)-a_problemDomain.domainBox().smallEnd(0)+1.);
158 
159  // Nominally, one layer of ghost cells is maintained permanently and
160  // individual computations may create local data with more
161  m_numGhost = 1;
162 
163  CH_assert(m_patchPluto != NULL);
164  CH_assert(isDefined());
165  m_patchPluto->define(m_problem_domain,m_dx,m_level,m_numGhost);
166 
167  // Get additional information from the patch integrator
168  m_numStates = m_patchPluto->numConserved();
169  m_ConsStateNames = m_patchPluto->ConsStateNames();
170  m_PrimStateNames = m_patchPluto->PrimStateNames();
171 }
172 
173 // Advance by one timestep
174 Real AMRLevelPluto::advance()
175 {
176  CH_assert(allDefined());
177 
178  if (s_verbosity >= 3)
179  {
180  pout() << "AMRLevelPluto::advance level " << m_level << " to time " << m_time + m_dt << endl;
181  }
182 
183  // Copy the new to the old
184 
185  DataIterator dit = m_UNew.dataIterator();
186  for(;dit.ok(); ++dit){
187  m_UOld[dit()].copy(m_UNew[dit()]);
188  }
189 
190  Real newDt = 0.0;
191 
192  // Set up arguments to LevelPluto::step based on whether there are
193  // coarser and finer levels
194 
195  // Undefined flux register in case we need it
196  LevelFluxRegister dummyFR;
197 
198  // Undefined leveldata in case we need it
199  const LevelData<FArrayBox> dummyData;
200 
201  // Set arguments to dummy values and then fix if real values are available
202  LevelFluxRegister* coarserFR = &dummyFR;
203  LevelFluxRegister* finerFR = &dummyFR;
204 
205  const LevelData<FArrayBox>* coarserDataOld = &dummyData;
206  const LevelData<FArrayBox>* coarserDataNew = &dummyData;
207 
208  Real tCoarserOld = 0.0;
209  Real tCoarserNew = 0.0;
210 
211  // A coarser level exists
212  if (m_hasCoarser)
213  {
214  AMRLevelPluto* coarserPtr = getCoarserLevel();
215 
216  // Recall that my flux register goes between my level and the next
217  // finer level
218  coarserFR = &coarserPtr->m_fluxRegister;
219 
220  coarserDataOld = &coarserPtr->m_UOld;
221  coarserDataNew = &coarserPtr->m_UNew;
222 
223  tCoarserNew = coarserPtr->m_time;
224  tCoarserOld = tCoarserNew - coarserPtr->m_dt;
225  }
226 
227  // A finer level exists
228  if (m_hasFiner)
229  {
230  // Recall that my flux register goes between my level and the next
231  // finer level
232  finerFR = &m_fluxRegister;
233  }
234 
235  // we don't need the flux in the simple hyperbolic case...
236  LevelData<FArrayBox> flux[SpaceDim];
237 
238 // Make sure m_time >= tCoarserOld to avoid negative alpha during levelPluto.step
239 // Thus if tm = tc*(1-epsilon) we set tm = tc
240 
241  if ( m_time < tCoarserOld &&
242  fabs(m_time - tCoarserOld) < 1.e-12*tCoarserOld){
243  m_time = tCoarserOld;
244  }
245 
246  g_intStage = 1;
247 
248  // Advance the solve one timestep
249  newDt = m_levelPluto.step(m_UNew,
250  flux,
251  *finerFR,
252  *coarserFR,
253  m_split_tags,
254  *coarserDataOld,
255  tCoarserOld,
256  *coarserDataNew,
257  tCoarserNew,
258  m_time,
259  m_dt,
260  m_cfl);
261 
262  #if (TIME_STEPPING == RK2)
263  g_intStage = 2;
264  Real DtCool; // The predictor returns the advective/diffusive timestep
265  // The corrector returns the cooling timestep
266  DtCool = m_levelPluto.step(m_UNew,
267  flux,
268  *finerFR,
269  *coarserFR,
270  m_split_tags,
271  *coarserDataOld,
272  tCoarserOld,
273  *coarserDataNew,
274  tCoarserNew,
275  m_time+m_dt,
276  m_dt,
277  m_cfl);
278 
279  newDt = Min(newDt,DtCool);
280  #endif
281 
282  // Update the time and store the new timestep
283  m_time += m_dt;
284  Real returnDt = m_cfl * newDt;
285 
286  m_dtNew = returnDt;
287 
288  return returnDt;
289 }
290 
291 Real AMRLevelPluto::getDlMin()
292 {
293 
294  Real dlMin = m_levelPluto.getDlMin();
295 
296  return dlMin;
297 
298 }
299 
300 // Things to do after a timestep
301 void AMRLevelPluto::postTimeStep()
302 {
303  CH_assert(allDefined());
304 
305  // Used for conservation tests
306  static Real orig_integral = 0.0;
307  static Real last_integral = 0.0;
308  static bool first = true;
309 
310  if (s_verbosity >= 3)
311  {
312  pout() << "AMRLevelPluto::postTimeStep " << m_level << endl;
313  }
314 
315  if (m_hasFiner)
316  {
317  // Reflux
318  Real scale = -1.0/m_dx;
319  m_fluxRegister.reflux(m_UNew,scale);
320 
321  // Average from finer level data
322  AMRLevelPluto* amrGodFinerPtr = getFinerLevel();
323 
324  amrGodFinerPtr->m_coarseAverage.averageToCoarse(m_UNew,
325  amrGodFinerPtr->m_UNew);
326  }
327 
328  if (s_verbosity >= 2 && m_level == 0)
329  {
330  int nRefFine = 1;
331 
332  pout() << "AMRLevelPluto::postTimeStep:" << endl;
333  pout() << " Sums:" << endl;
334  for (int comp = 0; comp < m_numStates; comp++)
335  {
336  Interval curComp(comp,comp);
337  Real integral = computeSum(m_UNew,NULL,nRefFine,m_dx,curComp);
338 
339  pout() << " " << setw(23)
340  << setprecision(16)
341  << setiosflags(ios::showpoint)
342  << setiosflags(ios::scientific)
343  << integral
344  << " --- " << m_ConsStateNames[comp];
345 
346  if (comp == 0 && !first) {
347  pout() << " (" << setw(23)
348  << setprecision(16)
349  << setiosflags(ios::showpoint)
350  << setiosflags(ios::scientific)
351  << (integral-last_integral)/last_integral
352  << " " << setw(23)
353  << setprecision(16)
354  << setiosflags(ios::showpoint)
355  << setiosflags(ios::scientific)
356  << (integral-orig_integral)/orig_integral
357  << ")";
358  }
359 
360  pout() << endl;
361 
362  if (comp == 0)
363  {
364  if (first)
365  {
366  orig_integral = integral;
367  first = false;
368  }
369 
370  last_integral = integral;
371  }
372  }
373  }
374 
375  if (s_verbosity >= 3)
376  {
377  pout() << "AMRLevelPluto::postTimeStep " << m_level << " finished" << endl;
378  }
379 }
380 
381 // Create tags for regridding
382 void AMRLevelPluto::tagCells(IntVectSet& a_tags)
383 {
384  CH_assert(allDefined());
385 
386  if (s_verbosity >= 3)
387  {
388  pout() << "AMRLevelPluto::tagCells " << m_level << endl;
389  }
390 
391  // Create tags based on undivided gradient of density
392  const DisjointBoxLayout& levelDomain = m_UNew.disjointBoxLayout();
393  IntVectSet localTags;
394  // If there is a coarser level interpolate undefined ghost cells
395  if (m_hasCoarser)
396  {
397  const AMRLevelPluto* amrGodCoarserPtr = getCoarserLevel();
398 
399  PiecewiseLinearFillPluto pwl;
400 
401  pwl.define(levelDomain,
402  amrGodCoarserPtr->m_UNew.disjointBoxLayout(),
403  m_numStates,
404  amrGodCoarserPtr->m_problem_domain,
405  amrGodCoarserPtr->m_ref_ratio,
406  m_dx,
407  1);
408 
409  pwl.fillInterp(m_UNew,
410  amrGodCoarserPtr->m_UNew,
411  amrGodCoarserPtr->m_UNew,
412  1.0,
413  0,
414  0,
415  m_numStates);
416  }
417  m_UNew.exchange(Interval(0,m_numStates-1));
418 
419  #if GEOMETRY != CARTESIAN
420  const LevelData<FArrayBox>& dV = m_levelPluto.getdV();
421  #endif
422 
423  // Compute relative gradient
424  DataIterator dit = levelDomain.dataIterator();
425  for (dit.begin(); dit.ok(); ++dit){
426  const Box& b = levelDomain[dit()];
427  FArrayBox gradFab(b,1);
428  FArrayBox& UFab = m_UNew[dit()];
429 
430  #if GEOMETRY != CARTESIAN
431  const FArrayBox& curdV = dV[dit()];
432  #else
433  const FArrayBox curdV;
434  #endif
435 
436  m_patchPluto->computeRefGradient(gradFab, UFab, curdV, b);
437 
438  // Tag where gradient exceeds threshold
439  BoxIterator bit(b);
440  for (bit.begin(); bit.ok(); ++bit){
441  const IntVect& iv = bit();
442  if (gradFab(iv) >= m_refineThresh) {
443  localTags |= iv;
444  }
445  }
446  }
447 
448  localTags.grow(m_tagBufferSize);
449 
450  // Need to do this in two steps unless a IntVectSet::operator &=
451  // (ProblemDomain) operator is defined
452  Box localTagsBox = localTags.minBox();
453  localTagsBox &= m_problem_domain;
454  localTags &= localTagsBox;
455 
456  a_tags = localTags;
457 }
458 
459 // Create tags at initialization
460 void AMRLevelPluto::tagCellsInit(IntVectSet& a_tags)
461 {
462  CH_assert(allDefined());
463 
464  if (s_verbosity >= 3)
465  {
466  pout() << "AMRLevelPolytropicGas::tagCellsInit " << m_level << endl;
467  }
468 
469  tagCells(a_tags);
470 }
471 
472 // Set up data on this level after regridding
473 void AMRLevelPluto::regrid(const Vector<Box>& a_newGrids)
474 {
475  CH_assert(allDefined());
476 
477  if (s_verbosity >= 3)
478  {
479  pout() << "AMRLevelPluto::regrid " << m_level << endl;
480  }
481 
482  // Save original grids and load balance
483  m_level_grids = a_newGrids;
484  m_grids = loadBalance(a_newGrids);
485 
486  if (s_verbosity >= 4)
487  {
488  // Indicate/guarantee that the indexing below is only for reading
489  // otherwise an error/assertion failure occurs
490  const DisjointBoxLayout& constGrids = m_grids;
491 
492  pout() << "new grids: " << endl;
493 
494  for (LayoutIterator lit = constGrids.layoutIterator(); lit.ok(); ++lit)
495  {
496  pout() << constGrids[lit()] << endl;
497  }
498  }
499 
500  // Save data for later
501  DataIterator dit = m_UNew.dataIterator();
502  for(;dit.ok(); ++dit){
503  m_UOld[dit()].copy(m_UNew[dit()]);
504  }
505 
506  // Reshape state with new grids
507  IntVect ivGhost = m_numGhost*IntVect::Unit;
508  m_UNew.define(m_grids,m_numStates,ivGhost);
509 
510 
511  // Set up data structures
512  levelSetup();
513 
514  // Interpolate from coarser level
515  if (m_hasCoarser)
516  {
517  AMRLevelPluto* amrGodCoarserPtr = getCoarserLevel();
518  m_fineInterp.interpToFine(m_UNew,amrGodCoarserPtr->m_UNew);
519  }
520 
521  // Copy from old state
522  m_UOld.copyTo(m_UOld.interval(),
523  m_UNew,
524  m_UNew.interval());
525 
526  m_UOld.define(m_grids,m_numStates,ivGhost);
527 }
528 
529 // Initialize grids
530 void AMRLevelPluto::initialGrid(const Vector<Box>& a_newGrids)
531 {
532  CH_assert(allDefined());
533 
534  if (s_verbosity >= 3)
535  {
536  pout() << "AMRLevelPluto::initialGrid " << m_level << endl;
537  }
538 
539  // Save original grids and load balance
540  m_level_grids = a_newGrids;
541  m_grids = loadBalance(a_newGrids);
542 
543  if (s_verbosity >= 4)
544  {
545  // Indicate/guarantee that the indexing below is only for reading
546  // otherwise an error/assertion failure occurs
547  const DisjointBoxLayout& constGrids = m_grids;
548 
549  pout() << "new grids: " << endl;
550  for (LayoutIterator lit = constGrids.layoutIterator(); lit.ok(); ++lit)
551  {
552  pout() << constGrids[lit()] << endl;
553  }
554  }
555 
556  // Define old and new state data structures
557  IntVect ivGhost = m_numGhost*IntVect::Unit;
558  m_UNew.define(m_grids,m_numStates,ivGhost);
559  m_UOld.define(m_grids,m_numStates,ivGhost);
560 
561  // Set up data structures
562  levelSetup();
563 }
564 
565 // Initialize data
566 void AMRLevelPluto::initialData()
567 {
568  if (s_verbosity >= 3)
569  {
570  pout() << "AMRLevelPluto::initialData " << m_level << endl;
571  }
572 
573  m_patchPluto->initialize(m_UNew);
574 
575 }
576 
577 // Things to do after initialization
578 void AMRLevelPluto::postInitialize()
579 {
580  CH_assert(allDefined());
581 
582  if (s_verbosity >= 3)
583  {
584  pout() << "AMRLevelPluto::postInitialize " << m_level << endl;
585  }
586 
587  if (m_hasFiner)
588  {
589  // Volume weighted average from finer level data
590  AMRLevelPluto* amrGodFinerPtr = getFinerLevel();
591 
592  amrGodFinerPtr->m_coarseAverage.averageToCoarse(m_UNew,
593  amrGodFinerPtr->m_UNew);
594  }
595 
596 }
597 
598 #ifdef CH_USE_HDF5
599 
600 // Write checkpoint header
601 void AMRLevelPluto::writeCheckpointHeader(HDF5Handle& a_handle) const
602 {
603  CH_assert(allDefined());
604 
605  if (s_verbosity >= 3)
606  {
607  pout() << "AMRLevelPluto::writeCheckpointHeader" << endl;
608  }
609 
610  // Setup the number of components
611  HDF5HeaderData header;
612  header.m_int["num_components"] = m_numStates;
613 
614  // Setup the component names
615  char compStr[30];
616  for (int comp = 0; comp < m_numStates; ++comp)
617  {
618  sprintf(compStr,"component_%d",comp);
619  header.m_string[compStr] = m_ConsStateNames[comp];
620  }
621 
622  // Write the header
623  header.writeToFile(a_handle);
624 
625  if (s_verbosity >= 3)
626  {
627  pout() << header << endl;
628  }
629 }
630 
631 // Write checkpoint data for this level
632 void AMRLevelPluto::writeCheckpointLevel(HDF5Handle& a_handle) const
633 {
634  CH_assert(allDefined());
635 
636  if (s_verbosity >= 3)
637  {
638  pout() << "AMRLevelPluto::writeCheckpointLevel" << endl;
639  }
640 
641  // Setup the level string
642  char levelStr[20];
643  sprintf(levelStr,"%d",m_level);
644  const std::string label = std::string("level_") + levelStr;
645 
646  a_handle.setGroup(label);
647 
648  // Setup the level header information
649  HDF5HeaderData header;
650 
651  header.m_int ["ref_ratio"] = m_ref_ratio;
652  header.m_int ["tag_buffer_size"] = m_tagBufferSize;
653  header.m_real["dx"] = m_dx;
654  header.m_int ["geometry"] = GEOMETRY;
655  header.m_int ["logr"] = CHOMBO_LOGR;
656  #if DIMENSIONS >= 2
657  header.m_real["g_x2stretch"] = g_x2stretch;
658  #endif
659  #if DIMENSIONS == 3
660  header.m_real["g_x3stretch"] = g_x3stretch;
661  #endif
662  header.m_real["domBeg1"] = g_domBeg[IDIR];
663  #if DIMENSIONS >= 2
664  header.m_real["domBeg2"] = g_domBeg[JDIR];
665  #endif
666  #if DIMENSIONS == 3
667  header.m_real["domBeg3"] = g_domBeg[KDIR];
668  #endif
669  header.m_real["dt"] = m_dt;
670  #ifdef GLM_MHD
671  header.m_real["ch"] = glm_ch_max;
672  #endif
673  header.m_real["time"] = m_time;
674  header.m_box ["prob_domain"] = m_problem_domain.domainBox();
675 
676  // Setup the periodicity info
677  D_TERM(if (m_problem_domain.isPeriodic(0))
678  {
679  header.m_int ["is_periodic_0"] = 1;
680  }
681  else
682  {
683  header.m_int ["is_periodic_0"] = 0;
684  } ,
685 
686  if (m_problem_domain.isPeriodic(1))
687  {
688  header.m_int ["is_periodic_1"] = 1;
689  }
690  else
691  {
692  header.m_int ["is_periodic_1"] = 0;
693  } ,
694 
695  if (m_problem_domain.isPeriodic(2))
696  {
697  header.m_int ["is_periodic_2"] = 1;
698  }
699  else
700  {
701  header.m_int ["is_periodic_2"] = 0;
702  } );
703 
704  // Write the header for this level
705  header.writeToFile(a_handle);
706 
707  if (s_verbosity >= 3)
708  {
709  pout() << header << endl;
710  }
711 
712  // Write the data for this level
713  #if GEOMETRY != CARTESIAN
714  LevelData<FArrayBox> tmp_U;
715  IntVect ivGhost = m_numGhost*IntVect::Unit;
716  tmp_U.define(m_grids,m_numStates,ivGhost);
717 
718  const LevelData<FArrayBox>& dV = m_levelPluto.getdV();
719 
720  DataIterator dit = m_UNew.dataIterator();
721  for(;dit.ok(); ++dit){
722  tmp_U[dit()].copy(m_UNew[dit()]);
723  FArrayBox& curU = tmp_U[dit()];
724  const FArrayBox& curdV = dV[dit()];
725  for (int ivar = 0; ivar < curU.nComp(); ivar++) curU.divide(curdV,0,ivar);
726  #if CHOMBO_CONS_AM == YES
727  #if ROTATING_FRAME == YES
728  Box curBox = curU.box();
729  for(BoxIterator bit(curBox); bit.ok(); ++bit) {
730  const IntVect& iv = bit();
731  curU(iv,iMPHI) /= curdV(iv,1);
732  curU(iv,iMPHI) -= curU(iv,RHO)*curdV(iv,1)*g_OmegaZ;
733  }
734  #else
735  curU.divide(curdV,1,iMPHI);
736  #endif
737  #endif
738  }
739  write(a_handle,tmp_U.boxLayout());
740  write(a_handle,tmp_U,"data");
741  #else
742  if (g_stretch_fact != 1.) {
743  LevelData<FArrayBox> tmp_U;
744  IntVect ivGhost = m_numGhost*IntVect::Unit;
745  tmp_U.define(m_grids,m_numStates,ivGhost);
746 
747  DataIterator dit = m_UNew.dataIterator();
748  for(;dit.ok(); ++dit){
749  tmp_U[dit()].copy(m_UNew[dit()]);
750  tmp_U[dit()] /= g_stretch_fact;
751  }
752  write(a_handle,tmp_U.boxLayout());
753  write(a_handle,tmp_U,"data");
754  } else {
755  write(a_handle,m_UNew.boxLayout());
756  write(a_handle,m_UNew,"data");
757  }
758  #endif
759 
760 }
761 
762 // Read checkpoint header
763 void AMRLevelPluto::readCheckpointHeader(HDF5Handle& a_handle)
764 {
765  if (s_verbosity >= 3)
766  {
767  pout() << "AMRLevelPluto::readCheckpointHeader" << endl;
768  }
769 
770  // Reader the header
771  HDF5HeaderData header;
772  header.readFromFile(a_handle);
773 
774  if (s_verbosity >= 3)
775  {
776  pout() << "hdf5 header data:" << endl;
777  pout() << header << endl;
778  }
779 
780  // Get the number of components
781  if (header.m_int.find("num_components") == header.m_int.end())
782  {
783  MayDay::Error("AMRLevelPluto::readCheckpointHeader: checkpoint file does not have num_components");
784  }
785 
786  int numStates = header.m_int["num_components"];
787  if (numStates != m_numStates)
788  {
789  MayDay::Error("AMRLevelPluto::readCheckpointHeader: num_components in checkpoint file does not match solver");
790  }
791 
792  // Get the component names
793  std::string stateName;
794  char compStr[60];
795  for (int comp = 0; comp < m_numStates; ++comp)
796  {
797  sprintf(compStr,"component_%d",comp);
798  if (header.m_string.find(compStr) == header.m_string.end())
799  {
800  MayDay::Error("AMRLevelPluto::readCheckpointHeader: checkpoint file does not have enough component names");
801  }
802 
803  stateName = header.m_string[compStr];
804  if (stateName != m_ConsStateNames[comp])
805  {
806  MayDay::Error("AMRLevelPluto::readCheckpointHeader: state_name in checkpoint does not match solver");
807  }
808  }
809 }
810 
811 // Read checkpoint data for this level
812 void AMRLevelPluto::readCheckpointLevel(HDF5Handle& a_handle)
813 {
814  if (s_verbosity >= 3)
815  {
816  pout() << "AMRLevelPluto::readCheckpointLevel" << endl;
817  }
818 
819  // Setup the level string
820  char levelStr[20];
821  sprintf(levelStr,"%d",m_level);
822  const std::string label = std::string("level_") + levelStr;
823 
824  // Read the header for this level
825  a_handle.setGroup(label);
826 
827  HDF5HeaderData header;
828  header.readFromFile(a_handle);
829 
830  if (s_verbosity >= 3)
831  {
832  pout() << "hdf5 header data:" << endl;
833  pout() << header << endl;
834  }
835 
836  // Get the refinement ratio
837  if (header.m_int.find("ref_ratio") == header.m_int.end())
838  {
839  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain ref_ratio");
840  }
841  m_ref_ratio = header.m_int["ref_ratio"];
842 
843  if (s_verbosity >= 2)
844  {
845  pout() << "read ref_ratio = " << m_ref_ratio << endl;
846  }
847 
848  // Get the tag buffer size
849  if (header.m_int.find("tag_buffer_size") == header.m_int.end())
850  {
851  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain tag_buffer_size");
852  }
853  m_tagBufferSize= header.m_int["tag_buffer_size"];
854 
855  if (s_verbosity >= 2)
856  {
857  pout() << "read tag_buffer_size = " << m_tagBufferSize << endl;
858  }
859 
860  // Get dx
861  if (header.m_real.find("dx") == header.m_real.end())
862  {
863  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain dx");
864  }
865  m_dx = header.m_real["dx"];
866 
867  if (s_verbosity >= 2)
868  {
869  pout() << "read dx = " << m_dx << endl;
870  }
871 
872  // Get dt
873  if (header.m_real.find("dt") == header.m_real.end())
874  {
875  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain dt");
876  }
877  m_dt = header.m_real["dt"];
878 
879  if (s_verbosity >= 2)
880  {
881  pout() << "read dt = " << m_dt << endl;
882  }
883 
884  // Get g_x2stretch and g_x3stretch factors
885 #if CH_SPACEDIM >= 2
886  if (header.m_real.find("g_x2stretch") == header.m_real.end())
887  {
888  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain g_x2stretch");
889  }
890  g_x2stretch = header.m_real["g_x2stretch"];
891 
892  if (s_verbosity >= 2)
893  {
894  pout() << "read g_x2stretch = " << g_x2stretch << endl;
895  }
896 #endif
897 #if CH_SPACEDIM == 3
898  if (header.m_real.find("g_x3stretch") == header.m_real.end())
899  {
900  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain g_x3stretch");
901  }
902  g_x3stretch = header.m_real["g_x3stretch"];
903 
904  if (s_verbosity >= 2)
905  {
906  pout() << "read g_x3stretch = " << g_x3stretch << endl;
907  }
908 #endif
909 
910 #if GEOMETRY == CARTESIAN
911  g_stretch_fact = g_x2stretch*g_x3stretch;
912 #endif
913 
914  #ifdef GLM_MHD
915  // Get glm_ch
916  if (header.m_real.find("ch") == header.m_real.end())
917  {
918  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain glm_ch");
919  }
920  glm_ch_max = header.m_real["ch"];
921 
922  if (s_verbosity >= 2)
923  {
924  pout() << "read ch max = " << glm_ch_max << endl;
925  }
926  #endif
927 
928  // Get the problem domain
929  if (header.m_box.find("prob_domain") == header.m_box.end())
930  {
931  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain prob_domain");
932  }
933 
934  Box domainBox = header.m_box["prob_domain"];
935 
936  // Get the periodicity info -- this is more complicated than it really
937  // needs to be in order to preserve backward compatibility
938  bool isPeriodic[SpaceDim];
939  D_TERM(if (!(header.m_int.find("is_periodic_0") == header.m_int.end()))
940  isPeriodic[0] = (header.m_int["is_periodic_0"] == 1);
941  else
942  isPeriodic[0] = false; ,
943 
944  if (!(header.m_int.find("is_periodic_1") == header.m_int.end()))
945  isPeriodic[1] = (header.m_int["is_periodic_1"] == 1);
946  else
947  isPeriodic[1] = false; ,
948 
949  if (!(header.m_int.find("is_periodic_2") == header.m_int.end()))
950  isPeriodic[2] = (header.m_int["is_periodic_2"] == 1);
951  else
952  isPeriodic[2] = false;);
953 
954  m_problem_domain = ProblemDomain(domainBox,isPeriodic);
955 
956  // Get the grids
957  Vector<Box> grids;
958  const int gridStatus = read(a_handle,grids);
959 
960  if (gridStatus != 0)
961  {
962  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain a Vector<Box>");
963  }
964 
965  // Create level domain
966  m_grids = loadBalance(grids);
967 
968  // Indicate/guarantee that the indexing below is only for reading
969  // otherwise an error/assertion failure occurs
970  const DisjointBoxLayout& constGrids = m_grids;
971 
972  LayoutIterator lit = constGrids.layoutIterator();
973  for (lit.begin(); lit.ok(); ++lit)
974  {
975  const Box& b = constGrids[lit()];
976  m_level_grids.push_back(b);
977  }
978 
979  if (s_verbosity >= 4)
980  {
981  pout() << "read level domain: " << endl;
982  LayoutIterator lit = m_grids.layoutIterator();
983  for (lit.begin(); lit.ok(); ++lit)
984  {
985  const Box& b = m_grids[lit()];
986  pout() << lit().intCode() << ": " << b << endl;
987  }
988  pout() << endl;
989  }
990 
991  // Reshape state with new grids
992  m_UNew.define(m_grids,m_numStates);
993  const int dataStatus = read<FArrayBox>(a_handle,
994  m_UNew,
995  "data",
996  m_grids);
997 
998  if (dataStatus != 0)
999  {
1000  MayDay::Error("AMRLevelPluto::readCheckpointLevel: file does not contain state data");
1001  }
1002 
1003  m_UOld.define(m_grids,m_numStates);
1004 
1005  // Set up data structures
1006  levelSetup();
1007 
1008  #if GEOMETRY != CARTESIAN
1009  const LevelData<FArrayBox>& dV = m_levelPluto.getdV();
1010  DataIterator dit = m_UNew.dataIterator();
1011  for(;dit.ok(); ++dit){
1012  FArrayBox& curU = m_UNew[dit()];
1013  const FArrayBox& curdV = dV[dit()];
1014  #if CHOMBO_CONS_AM == YES
1015  #if ROTATING_FRAME == YES
1016  Box curBox = curU.box();
1017  for(BoxIterator bit(curBox); bit.ok(); ++bit) {
1018  const IntVect& iv = bit();
1019  curU(iv,iMPHI) += curU(iv,RHO)*curdV(iv,1)*g_OmegaZ;
1020  curU(iv,iMPHI) *= curdV(iv,1);
1021  }
1022  #else
1023  curU.mult(curdV,1,iMPHI);
1024  #endif
1025  #endif
1026  for (int ivar = 0; ivar < curU.nComp(); ivar++) curU.mult(curdV,0,ivar);
1027  }
1028  #else
1029  if (g_stretch_fact != 1.) {
1030  DataIterator dit = m_UNew.dataIterator();
1031  for(;dit.ok(); ++dit) {
1032  FArrayBox& curU = m_UNew[dit()];
1033  curU *= g_stretch_fact;
1034  }
1035  }
1036  #endif
1037 }
1038 
1039 // Write plotfile header
1040 void AMRLevelPluto::writePlotHeader(HDF5Handle& a_handle) const
1041 {
1042  CH_assert(allDefined());
1043 
1044  if (s_verbosity >= 3)
1045  {
1046  pout() << "AMRLevelPluto::writePlotHeader" << endl;
1047  }
1048 
1049  // Setup the number of components
1050  HDF5HeaderData header;
1051  header.m_int["num_components"] = m_numStates;
1052 
1053  // Setup the component names
1054  char compStr[30];
1055  for (int comp = 0; comp < m_numStates; ++comp)
1056  {
1057  sprintf(compStr,"component_%d",comp);
1058  header.m_string[compStr] = m_PrimStateNames[comp];
1059  }
1060 
1061  // Write the header
1062  header.writeToFile(a_handle);
1063  a_handle.setGroup("/Expressions");
1064  HDF5HeaderData expressions;
1065  DefineExpressions(expressions);
1066  expressions.writeToFile(a_handle);
1067 
1068  if (s_verbosity >= 3)
1069  {
1070  pout() << header << endl;
1071  }
1072 }
1073 
1074 // Write plotfile data for this level
1075 void AMRLevelPluto::writePlotLevel(HDF5Handle& a_handle) const
1076 {
1077  CH_assert(allDefined());
1078 
1079  // Setup the level string
1080  char levelStr[20];
1081  sprintf(levelStr,"%d",m_level);
1082  const std::string label = std::string("level_") + levelStr;
1083 
1084  a_handle.setGroup(label);
1085 
1086  // Setup the level header information
1087  HDF5HeaderData header;
1088 
1089  header.m_int ["ref_ratio"] = m_ref_ratio;
1090  header.m_real["dx"] = m_dx;
1091  header.m_int ["geometry"] = GEOMETRY;
1092  header.m_int ["logr"] = CHOMBO_LOGR;
1093  #if DIMENSIONS >= 2
1094  header.m_real["g_x2stretch"] = g_x2stretch;
1095  #endif
1096  #if DIMENSIONS == 3
1097  header.m_real["g_x3stretch"] = g_x3stretch;
1098  #endif
1099  header.m_real["domBeg1"] = g_domBeg[IDIR];
1100  #if DIMENSIONS >= 2
1101  header.m_real["domBeg2"] = g_domBeg[JDIR];
1102  #endif
1103  #if DIMENSIONS == 3
1104  header.m_real["domBeg3"] = g_domBeg[KDIR];
1105  #endif
1106  header.m_real["dt"] = m_dt;
1107  header.m_real["time"] = m_time;
1108  header.m_box ["prob_domain"] = m_problem_domain.domainBox();
1109 
1110  // Write the header for this level
1111  header.writeToFile(a_handle);
1112 
1113  if (s_verbosity >= 3)
1114  {
1115  pout() << header << endl;
1116  }
1117 
1118  // Write the data for this level
1119  LevelData<FArrayBox> tmp_U;
1120  IntVect ivGhost = m_numGhost*IntVect::Unit;
1121  tmp_U.define(m_grids,m_numStates,ivGhost);
1122 
1123  #if GEOMETRY != CARTESIAN
1124  const LevelData<FArrayBox>& dV = m_levelPluto.getdV();
1125  #endif
1126 
1127  DataIterator dit = m_UNew.dataIterator();
1128  for(;dit.ok(); ++dit){
1129  tmp_U[dit()].copy(m_UNew[dit()]);
1130  FArrayBox& curU = tmp_U[dit()];
1131  #if GEOMETRY != CARTESIAN
1132  const FArrayBox& curdV = dV[dit()];
1133  for (int ivar = 0; ivar < curU.nComp(); ivar++) curU.divide(curdV,0,ivar);
1134  #if CHOMBO_CONS_AM == YES
1135  #if ROTATING_FRAME == YES
1136  Box curBox = curU.box();
1137  for(BoxIterator bit(curBox); bit.ok(); ++bit) {
1138  const IntVect& iv = bit();
1139  curU(iv,iMPHI) /= curdV(iv,1);
1140  curU(iv,iMPHI) -= curU(iv,RHO)*curdV(iv,1)*g_OmegaZ;
1141  }
1142  #else
1143  curU.divide(curdV,1,iMPHI);
1144  #endif
1145  #endif
1146  #else
1147  if (g_stretch_fact != 1.) curU /= g_stretch_fact;
1148  #endif
1149 
1150  m_patchPluto->convertFArrayBox(curU); /* -- convert data to primitive -- */
1151  }
1152  write(a_handle,tmp_U.boxLayout());
1153  write(a_handle,tmp_U,"data");
1154 
1155 }
1156 
1157 void AMRLevelPluto::DefineExpressions(HDF5HeaderData& a_expressions) const
1158 {
1159 
1160  char str_expr[128];
1161 
1162  // Grid/geometry expressions
1163 
1164  #if GEOMETRY == CARTESIAN
1165  sprintf (str_expr,"coords(Mesh)[0]+nodal_constant(Mesh,%12.6e)",g_domBeg[IDIR]);
1166  a_expressions.m_string["scalar X"] = str_expr;
1167  #if DIMENSIONS == 2
1168  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[1]+nodal_constant(Mesh,%12.6e)",g_x2stretch,g_domBeg[JDIR]);
1169  a_expressions.m_string["scalar Y"] = str_expr;
1170  a_expressions.m_string["vector Displacement"] = "{X,Y}-coords(Mesh)";
1171  #endif
1172  #if DIMENSIONS == 3
1173  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[1]+nodal_constant(Mesh,%12.6e)",g_x2stretch,g_domBeg[JDIR]);
1174  a_expressions.m_string["scalar Y"] = str_expr;
1175  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[2]+nodal_constant(Mesh,%12.6e)",g_x3stretch,g_domBeg[KDIR]);
1176  a_expressions.m_string["scalar Z"] = str_expr;
1177  a_expressions.m_string["vector Displacement"] = "{X,Y,Z}-coords(Mesh)";
1178  #endif
1179  #endif
1180 
1181  #if GEOMETRY == CYLINDRICAL
1182  sprintf (str_expr,"coords(Mesh)[0]+nodal_constant(Mesh,%12.6e)",g_domBeg[IDIR]);
1183  a_expressions.m_string["scalar R"] = str_expr;
1184  #if DIMENSIONS == 2
1185  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[1]+nodal_constant(Mesh,%12.6e)",g_x2stretch,g_domBeg[JDIR]);
1186  a_expressions.m_string["scalar Z"] = str_expr;
1187  a_expressions.m_string["vector Displacement"] = "{R,Z}-coords(Mesh)";
1188  #endif
1189  #endif
1190 
1191  #if GEOMETRY == SPHERICAL
1192  #if CHOMBO_LOGR == NO
1193  sprintf (str_expr,"coords(Mesh)[0]+nodal_constant(Mesh,%12.6e)",g_domBeg[IDIR]);
1194  a_expressions.m_string["scalar R"] = str_expr;
1195  #else
1196  sprintf (str_expr,"exp(coords(Mesh)[0])*nodal_constant(Mesh,%12.6e)",g_domBeg[IDIR]);
1197  a_expressions.m_string["scalar R"] = str_expr;
1198  #endif
1199  #if DIMENSIONS == 2
1200  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[1]+nodal_constant(Mesh,%12.6e)",g_x2stretch,g_domBeg[JDIR]);
1201  a_expressions.m_string["scalar Theta"] = str_expr;
1202  a_expressions.m_string["scalar X"] = "R*sin(Theta)";
1203  a_expressions.m_string["scalar Z"] = "R*cos(Theta)";
1204  a_expressions.m_string["vector Displacement"] = "{X,Z}-coords(Mesh)";
1205  #endif
1206  #if DIMENSIONS == 3
1207  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[1]+nodal_constant(Mesh,%12.6e)",g_x2stretch,g_domBeg[JDIR]);
1208  a_expressions.m_string["scalar Theta"] = str_expr;
1209  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[2]+nodal_constant(Mesh,%12.6e)",g_x3stretch,g_domBeg[KDIR]);
1210  a_expressions.m_string["scalar Phi"] = str_expr;
1211  a_expressions.m_string["scalar X"] = "R*sin(Theta)*cos(Phi)";
1212  a_expressions.m_string["scalar Y"] = "R*sin(Theta)*sin(Phi)";
1213  a_expressions.m_string["scalar Z"] = "R*cos(Theta)";
1214  a_expressions.m_string["vector Displacement"] = "{X,Y,Z}-coords(Mesh)";
1215  #endif
1216  #endif
1217 
1218  #if GEOMETRY == POLAR
1219  #if CHOMBO_LOGR == NO
1220  sprintf (str_expr,"coords(Mesh)[0]+nodal_constant(Mesh,%12.6e)",g_domBeg[IDIR]);
1221  a_expressions.m_string["scalar R"] = str_expr;
1222  #else
1223  sprintf (str_expr,"exp(coords(Mesh)[0])*nodal_constant(Mesh,%12.6e)",g_domBeg[IDIR]);
1224  a_expressions.m_string["scalar R"] = str_expr;
1225  #endif
1226  #if DIMENSIONS == 2
1227  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[1]+nodal_constant(Mesh,%12.6e)",g_x2stretch,g_domBeg[JDIR]);
1228  a_expressions.m_string["scalar Phi"] = str_expr;
1229  a_expressions.m_string["scalar X"] = "R*cos(Phi)";
1230  a_expressions.m_string["scalar Y"] = "R*sin(Phi)";
1231  a_expressions.m_string["vector Displacement"] = "{X,Y}-coords(Mesh)";
1232  #endif
1233  #if DIMENSIONS == 3
1234  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[1]+nodal_constant(Mesh,%12.6e)",g_x2stretch,g_domBeg[JDIR]);
1235  a_expressions.m_string["scalar Phi"] = str_expr;
1236  a_expressions.m_string["scalar X"] = "R*cos(Phi)";
1237  a_expressions.m_string["scalar Y"] = "R*sin(Phi)";
1238  sprintf (str_expr,"nodal_constant(Mesh,%12.6e)*coords(Mesh)[2]+nodal_constant(Mesh,%12.6e)",g_x3stretch,g_domBeg[KDIR]);
1239  a_expressions.m_string["scalar Z"] = str_expr;
1240  a_expressions.m_string["vector Displacement"] = "{X,Y,Z}-coords(Mesh)";
1241  #endif
1242  #endif
1243 
1244 }
1245 
1246 #endif
1247 
1248 // Returns the dt computed earlier for this level
1249 Real AMRLevelPluto::computeDt()
1250 {
1251  CH_assert(allDefined());
1252 
1253  if (s_verbosity >= 3)
1254  {
1255  pout() << "AMRLevelPluto::computeDt " << m_level << endl;
1256  }
1257 
1258  Real newDt;
1259  newDt = m_dtNew;
1260 
1261  return newDt;
1262 }
1263 
1264 // Compute dt using initial data
1265 Real AMRLevelPluto::computeInitialDt()
1266 {
1267  CH_assert(allDefined());
1268 
1269  if (s_verbosity >= 3)
1270  {
1271  pout() << "AMRLevelPluto::computeInitialDt " << m_level << endl;
1272  }
1273 
1274  Real newDT = m_initial_dt_multiplier*m_dx/m_domainLength;
1275 
1276  return newDT;
1277 }
1278 
1279 const LevelData<FArrayBox>& AMRLevelPluto::getStateNew() const
1280 {
1281  CH_assert(allDefined());
1282 
1283  return m_UNew;
1284 }
1285 
1286 const LevelData<FArrayBox>& AMRLevelPluto::getStateOld() const
1287 {
1288  CH_assert(allDefined());
1289 
1290  return m_UOld;
1291 }
1292 
1293 bool AMRLevelPluto::allDefined() const
1294 {
1295  return isDefined() &&
1296  m_paramsDefined ;
1297 }
1298 
1299 // Create a load-balanced DisjointBoxLayout from a collection of Boxes
1300 DisjointBoxLayout AMRLevelPluto::loadBalance(const Vector<Box>& a_grids)
1301 {
1302  CH_assert(allDefined());
1303 
1304  // Load balance and create boxlayout
1305  Vector<int> procMap;
1306 
1307 // use the four-argument call to LoadBalance to have more control over it.
1308 // this may be useful for dynamic load balancing (e.g. when cooling is
1309 // employed.
1310 
1311  Vector<long long> computeLoads(a_grids.size());
1312 
1313  for (int i = 0; i < a_grids.size(); ++i) {
1314 // if (i < 7) computeLoads[i] = a_grids[i].numPts()*(long long)100;
1315  computeLoads[i] = a_grids[i].numPts();
1316  }
1317 
1318  LoadBalance(procMap, computeLoads, a_grids, numProc());
1319 
1320  // appears to be faster for all procs to do the loadbalance (ndk)
1321 // LoadBalance(procMap,a_grids);
1322 
1323  if (s_verbosity >= 4)
1324  {
1325  pout() << "AMRLevelPluto::loadBalance (level "<<m_level<<"): procesor map: " << endl;
1326 
1327  for (int igrid = 0; igrid < a_grids.size(); ++igrid)
1328  {
1329  pout() << igrid << "(load = " << computeLoads[igrid] <<
1330  "): " << procMap[igrid] << " " << endl;
1331  }
1332  pout() << endl;
1333  }
1334 
1335  DisjointBoxLayout dbl(a_grids,procMap,m_problem_domain);
1336  dbl.close();
1337 
1338  return dbl;
1339 }
1340 
1341 // Setup menagerie of data structures
1342 void AMRLevelPluto::levelSetup()
1343 {
1344  CH_assert(allDefined());
1345 
1346  if (s_verbosity >= 3)
1347  {
1348  pout() << "AMRLevelPluto::levelSetup " << m_level << endl;
1349  }
1350 
1351  AMRLevelPluto* amrGodCoarserPtr = getCoarserLevel();
1352  AMRLevelPluto* amrGodFinerPtr = getFinerLevel();
1353 
1354  m_hasCoarser = (amrGodCoarserPtr != NULL);
1355  m_hasFiner = (amrGodFinerPtr != NULL);
1356 
1357 
1358 #ifdef SKIP_SPLIT_CELLS
1359  // Mark split/unsplit cells
1360  m_split_tags.define(m_grids,1);
1361  for (DataIterator dit = m_grids.dataIterator(); dit.ok(); ++dit) {
1362  m_split_tags[dit()].setVal(1.);
1363  }
1364 #endif
1365 
1366  if (m_hasCoarser)
1367  {
1368  int nRefCrse = m_coarser_level_ptr->refRatio();
1369 
1370  m_coarseAverage.define(m_grids,
1371  m_numStates,
1372  nRefCrse);
1373 
1374  m_fineInterp.define(m_grids,
1375  m_numStates,
1376  nRefCrse,
1377  m_dx,
1378  m_problem_domain);
1379 
1380  const DisjointBoxLayout& coarserLevelDomain = amrGodCoarserPtr->m_grids;
1381 
1382 #ifdef SKIP_SPLIT_CELLS
1383  // Mark split/unsplit cells of the coarser level
1384  amrGodCoarserPtr->mark_split(m_grids);
1385 #endif
1386 
1387  // Maintain levelPluto
1388  m_levelPluto.define(m_grids,
1389  coarserLevelDomain,
1390  m_problem_domain,
1391  nRefCrse,
1392  m_level,
1393  m_dx,
1394  m_patchPlutoFactory,
1395  m_hasCoarser,
1396  m_hasFiner);
1397 
1398  // This may look twisted but you have to do this this way because the
1399  // coarser levels get setup before the finer levels so, since a flux
1400  // register lives between this level and the next FINER level, the finer
1401  // level has to do the setup because it is the only one with the
1402  // information at the time of construction.
1403 
1404  // Maintain flux registers
1405  amrGodCoarserPtr->m_fluxRegister.define(m_grids,
1406  amrGodCoarserPtr->m_grids,
1407  m_problem_domain,
1408  amrGodCoarserPtr->m_ref_ratio,
1409  m_numStates);
1410  amrGodCoarserPtr->m_fluxRegister.setToZero();
1411  }
1412  else
1413  {
1414  m_levelPluto.define(m_grids,
1415  DisjointBoxLayout(),
1416  m_problem_domain,
1417  m_ref_ratio,
1418  m_level,
1419  m_dx,
1420  m_patchPlutoFactory,
1421  m_hasCoarser,
1422  m_hasFiner);
1423  }
1424 }
1425 
1426 // Get the next coarser level
1427 AMRLevelPluto* AMRLevelPluto::getCoarserLevel() const
1428 {
1429  CH_assert(allDefined());
1430 
1431  AMRLevelPluto* amrGodCoarserPtr = NULL;
1432 
1433  if (m_coarser_level_ptr != NULL)
1434  {
1435  amrGodCoarserPtr = dynamic_cast<AMRLevelPluto*>(m_coarser_level_ptr);
1436 
1437  if (amrGodCoarserPtr == NULL)
1438  {
1439  MayDay::Error("AMRLevelPluto::getCoarserLevel: dynamic cast failed");
1440  }
1441  }
1442 
1443  return amrGodCoarserPtr;
1444 }
1445 
1446 // Get the next finer level
1447 AMRLevelPluto* AMRLevelPluto::getFinerLevel() const
1448 {
1449  CH_assert(allDefined());
1450 
1451  AMRLevelPluto* amrGodFinerPtr = NULL;
1452 
1453  if (m_finer_level_ptr != NULL)
1454  {
1455  amrGodFinerPtr = dynamic_cast<AMRLevelPluto*>(m_finer_level_ptr);
1456 
1457  if (amrGodFinerPtr == NULL)
1458  {
1459  MayDay::Error("AMRLevelPluto::getFinerLevel: dynamic cast failed");
1460  }
1461  }
1462 
1463  return amrGodFinerPtr;
1464 }
1465 
1466 // Mark Split/unsplit cells of this level
1467 void AMRLevelPluto::mark_split(const DisjointBoxLayout& finerLevelDomain)
1468 {
1469  CH_assert(allDefined());
1470 
1471  DisjointBoxLayout finerDomain;
1472 
1473  coarsen(finerDomain,finerLevelDomain,m_ref_ratio);
1474 
1475  DataIterator dit = m_grids.dataIterator();
1476  LayoutIterator lit = finerDomain.layoutIterator();
1477  for (dit.begin(); dit.ok(); ++dit)
1478  {
1479  FArrayBox& splitU = m_split_tags[dit()];
1480  splitU.setVal(1.);
1481  Box splitBox = splitU.box();
1482 
1483  for (lit.begin(); lit.ok(); ++lit)
1484  {
1485  Box finBox = finerDomain.get(lit());
1486  const Box intBox = splitBox & finBox;
1487 
1488  BoxIterator bit(intBox);
1489  for (bit.begin(); bit.ok(); ++bit)
1490  {
1491  const IntVect& iv = bit();
1492  splitU(iv) = 0.0;
1493  }
1494  }
1495  }
1496 
1497 }
1498 
1499 #include "NamespaceFooter.H"
#define iMPHI
Definition: mod_defs.h:71
#define RHO
Definition: mod_defs.h:19
int g_intStage
Gives the current integration stage of the time stepping method (predictor = 0, 1st corrector = 1...
Definition: globals.h:98
#define g_OmegaZ
Definition: init.c:64
def read(fname)
Definition: file.py:29
#define KDIR
Definition: pluto.h:195
#define IDIR
Definition: pluto.h:193
def find(fname, str, action, want)
Definition: file.py:150
if(divB==NULL)
Definition: analysis.c:10
double g_domBeg[3]
Lower limits of the computational domain.
Definition: globals.h:125
#define GEOMETRY
Definition: definitions_01.h:4
int i
Definition: analysis.c:2
#define JDIR
Definition: pluto.h:194