aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetRegrid/src/moving.cc
blob: 58a010e8ae8fae23281e15fcfc97275161fe2958 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
#include <cassert>
#include <cmath>

#include "cctk.h"
#include "cctk_Parameters.h"

#include "gh.hh"

#include "carpet.hh"
#include "regrid.hh"



namespace CarpetRegrid {
  
  using namespace std;
  using namespace Carpet;
  
  
  
  int Moving (cGH const * const cctkGH,
              gh<dim> const & hh,
              gh<dim>::rexts  & bbsss,
              gh<dim>::rbnds  & obss,
              gh<dim>::rprocs & pss)
  {
    DECLARE_CCTK_ARGUMENTS;
    DECLARE_CCTK_PARAMETERS;
    
    assert (refinement_levels >= 1);
    
    assert (bbsss.size() >= 1);
    
    bbsss.resize (refinement_levels);
    obss.resize (refinement_levels);
    pss.resize (refinement_levels);
    
    bvect const symmetric (symmetry_x, symmetry_y, symmetry_z);
    bbvect const ob (false);
    
    assert (! smart_outer_boundaries);
    
    for (size_t rl=1; rl<bbsss.size(); ++rl) {
      
      // calculate new extent
      CCTK_REAL const argument = 2*M_PI * moving_circle_frequency * cctk_time;
      rvect const pos
        (moving_centre_x + moving_circle_radius * cos(argument),
         moving_centre_y + moving_circle_radius * sin(argument),
         moving_centre_z);
      CCTK_REAL const radius = moving_region_radius / ipow(reffact, rl-1);
      
      rvect const rlb (symmetric.ifthen (rvect(0),      pos - rvect(radius)));
      rvect const rub (symmetric.ifthen (rvect(radius), pos + rvect(radius)));
      
      vector<ibbox> bbs;
      gh<dim>::cbnds obs;
      
      ManualCoordinates_OneLevel
        (cctkGH, hh, rl, refinement_levels, rlb, rub, ob, bbs, obs);
      
      // make multiprocessor aware
      gh<dim>::cprocs ps;
      SplitRegions (cctkGH, bbs, obs, ps);
      
      // make multigrid aware
      vector<vector<ibbox> > bbss;
      MakeMultigridBoxes (cctkGH, bbs, obs, bbss);
      
      bbsss.at(rl) = bbss;
      obss.at(rl) = obs;
      pss.at(rl) = ps;
      
    } // for rl
    
    return 1;
  }
  
} // namespace CarpetRegrid