Discontinuous Galerkin Library
#include "dg/algorithm.h"
Loading...
Searching...
No Matches
refined_gridX.h
Go to the documentation of this file.
1#pragma once
2
3#include "interpolation.h"
4#include "evaluationX.h"
5#include "weightsX.h"
6#include "gridX.h"
7#include "base_geometryX.h"
8#include "refined_grid.h"
9
10
11namespace dg
12{
15
18template<class real_type>
20{
29 void generate( const RealGridX2d<real_type>& g_old, thrust::host_vector<real_type>& weightsX, thrust::host_vector<real_type>& weightsY, thrust::host_vector<real_type>& abscissasX, thrust::host_vector<real_type>& abscissasY) const
30 {
31 thrust::host_vector<real_type> wx, ax, wy, ay;
32 RealGrid1d<real_type> gx( g_old.x0(), g_old.x1(), g_old.n(), g_old.Nx(), g_old.bcx());
33 do_generateX(gx,g_old.inner_Nx(), wx,ax);
34 RealGridX1d<real_type> gy( g_old.y0(), g_old.y1(), g_old.fy(), g_old.n(), g_old.Ny(), g_old.bcy());
35 do_generateY(gy,wy,ay);
36 unsigned size=wx.size()*wy.size();
37 weightsX.resize(size), weightsY.resize(size);
38 abscissasX.resize(size), abscissasY.resize(size);
39 //now make product space
40 for( unsigned i=0; i<wy.size(); i++)
41 for( unsigned j=0; j<wx.size(); j++)
42 {
43 weightsX[i*wx.size()+j] = wx[j];
44 weightsY[i*wx.size()+j] = wy[i];
45 abscissasX[i*wx.size()+j] = ax[j];
46 abscissasY[i*wx.size()+j] = ay[i];
47 }
48 }
49
54 unsigned nx_new( unsigned Nx_old, real_type fx_old) const
55 {
56 return do_Nx_new(Nx_old, fx_old);
57 }
58 unsigned ny_new( unsigned Ny_old, real_type fy_old) const
59 {
60 return do_Ny_new(Ny_old, fy_old);
61 }
62 real_type fx_new( unsigned Nx_old, real_type fx_old) const
63 {
64 return do_fx_new(Nx_old, fx_old);
65 }
66 real_type fy_new( unsigned Ny_old, real_type fy_old) const
67 {
68 return do_fy_new(Ny_old, fy_old);
69 }
70 virtual aRealRefinementX2d* clone()const=0;
71 virtual ~aRealRefinementX2d() = default;
72 protected:
73 aRealRefinementX2d() = default;
76 private:
77 virtual void do_generateX( const RealGrid1d<real_type>& gx, unsigned nodeXX, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const =0;
78 virtual void do_generateY( const RealGridX1d<real_type>& gy, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const =0;
79 virtual unsigned do_Nx_new( unsigned Nx_old, real_type fx) const =0;
80 virtual unsigned do_Ny_new( unsigned Ny_old, real_type fy) const =0;
81 virtual real_type do_fx_new( unsigned Nx_old, real_type fx) const =0;
82 virtual real_type do_fy_new( unsigned Ny_old, real_type fy) const =0;
83};
84
88template<class real_type>
90{
91 virtual RealIdentityXRefinement* clone()const override final{
92 return new RealIdentityXRefinement(*this);
93 }
94 private:
95 virtual void do_generateX( const RealGrid1d<real_type>& g, unsigned nodeXX, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const override final{
96 weights=dg::create::weights(g);
97 abscissas = g.abscissas()[0];
98 }
99 virtual void do_generateY( const RealGridX1d<real_type>& g, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const override final{
100 weights=dg::create::weights(g);
101 abscissas = g.abscissas()[0];
102 }
103 virtual unsigned do_Nx_new( unsigned Nx_old, real_type fx) const override final{ return Nx_old; }
104 virtual unsigned do_Ny_new( unsigned Ny_old, real_type fy) const override final{ return Ny_old; }
105 virtual real_type do_fx_new( unsigned Nx_old, real_type fx) const override final{ return fx; }
106 virtual real_type do_fy_new( unsigned Ny_old, real_type fy) const override final{ return fy; }
107};
108
112template<class real_type>
114{
115 RealEquidistXRefinement( unsigned add_x, unsigned add_y, unsigned howmanyX=1, unsigned howmanyY=1): add_x_(add_x), howm_x_(howmanyX), add_y_(add_y), howm_y_(howmanyY){ }
116 virtual RealEquidistXRefinement* clone()const override final{return new RealEquidistXRefinement(*this);}
117 private:
118 unsigned add_x_, howm_x_, add_y_, howm_y_;
119 virtual void do_generateX( const RealGrid1d<real_type>& gx, unsigned nodeXX, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const override final
120 {
121 RealEquidistRefinement<real_type> equi(add_x_, nodeXX, howm_x_);
122 equi.generate(gx,weights,abscissas);
123 }
124 virtual void do_generateY( const RealGridX1d<real_type>& g, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const override final
125 {
126 RealEquidistRefinement<real_type> equi0(add_y_,0,howm_y_);
127 if( add_y_ == 0 || howm_y_ == 0) {
128 equi0.generate(g.grid(), weights, abscissas);
129 return;
130 }
131 if( g.f() == 0) {
132 equi0.generate( RealGrid1d<real_type>( g.x0(), g.x1(), g.n(), g.N(), dg::PER), weights, abscissas);
133 return;
134 }
135 thrust::host_vector<real_type> w1, w2, w3;
136 thrust::host_vector<real_type> a1, a2, a3;
137 unsigned node1 = g.outer_N();
138 RealEquidistRefinement<real_type> equi1(add_y_, node1, howm_y_);
139 equi1.generate( RealGrid1d<real_type>(g.x0() ,g.x0()+g.outer_N()*g.h(), g.n(), g.outer_N(), dg::DIR), w1,a1); //left side
140 equi0.generate( RealGrid1d<real_type>(g.x0()+g.outer_N()*g.h(),g.x1()-g.outer_N()*g.h(), g.n(), g.inner_N(), dg::PER), w2,a2);//inner side
141 equi0.generate( RealGrid1d<real_type>(g.x1()-g.outer_N()*g.h(),g.x1() , g.n(), g.outer_N(), dg::DIR), w3,a3);//right side
142 //now combine and unnormalize weights
143 thrust::host_vector<real_type> wtot( w1.size() + w2.size() + w3.size()), atot( wtot);
144 for( unsigned i=0; i<w1.size() ; i++)
145 wtot[i] = w1[i]*equi1.N_new( g.outer_N(), dg::DIR)/g.outer_N();
146 for( unsigned i=0; i<w2.size(); i++)
147 wtot[w1.size()+i] = w2[i]*equi0.N_new( g.inner_N(), dg::PER)/g.inner_N();
148 for( unsigned i=0; i<w3.size(); i++)
149 wtot[w1.size()+w2.size()+i] = w3[i]*equi0.N_new( g.outer_N(), dg::DIR)/g.outer_N();
150 weights = wtot;
151 abscissas = detail::normalize_weights_and_compute_abscissas( g.grid(), weights);
152 }
153 virtual unsigned do_Ny_new( unsigned Ny, real_type fy) const override final{
154 if( fy==0 ) return Ny + 2*add_y_;
155 return Ny + 4*add_y_;
156 }
157 virtual unsigned do_Nx_new( unsigned Nx, real_type fx) const override final{
158 if( fx==0 ) return Nx + add_x_;
159 return Nx + 2*add_x_;
160 }
161 virtual real_type do_fx_new( unsigned Nx, real_type fx) const override final{
162 if( fx==0 ) return 0;
163 return (fx*(real_type)Nx + (real_type)add_x_)/(real_type)(Nx+2.*add_x_);
164 }
165 virtual real_type do_fy_new( unsigned Ny, real_type fy) const override final{
166 if( fy==0 ) return 0;
167 return (fy*(real_type)Ny + (real_type)add_y_)/(real_type)(Ny+4.*add_y_);
168 }
169};
170
174template<class real_type>
176{
177 RealExponentialXRefinement( unsigned add_x, unsigned add_y, unsigned howmanyX=1, unsigned howmanyY=1): add_x_(add_x), howm_x_(howmanyX), add_y_(add_y), howm_y_(howmanyY){ }
179 private:
180 unsigned add_x_, howm_x_, add_y_, howm_y_;
181 virtual void do_generateX( const RealGrid1d<real_type>& gx, unsigned nodeXX, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const override final
182 {
183 RealEquidistRefinement<real_type> equi(add_x_, nodeXX, howm_x_);
184 equi.generate(gx,weights,abscissas);
185 }
186 virtual void do_generateY( const RealGridX1d<real_type>& g, thrust::host_vector<real_type>& weights, thrust::host_vector<real_type>& abscissas) const override final
187 {
188 RealExponentialRefinement<real_type> expo0( add_x_, 0);
189 if( add_y_ == 0) { return expo0.generate( g.grid(), weights, abscissas); }
190 if( g.f() == 0) { return expo0.generate( RealGrid1d<real_type>( g.x0(), g.x1(), g.n(), g.N(), dg::PER), weights, abscissas); }
191 thrust::host_vector<real_type> w1, w2, w3;
192 thrust::host_vector<real_type> a1, a2, a3;
193 unsigned node1 = g.outer_N();
194 RealExponentialRefinement<real_type> expo1( add_y_, node1);
195 expo1.generate( RealGrid1d<real_type>(g.x0(),g.x1(),g.n(), g.outer_N(), dg::DIR), w1,a1); //left side
196 expo0.generate( RealGrid1d<real_type>(g.x0(),g.x1(),g.n(), g.inner_N(), dg::PER), w2,a2); //inner side
197 expo0.generate( RealGrid1d<real_type>(g.x0(),g.x1(),g.n(), g.outer_N(), dg::DIR), w3,a3); //right side
198 //now combine unnormalized weights
199 thrust::host_vector<real_type> wtot( g.size() + 4*g.n()*add_x_);
200 for( unsigned i=0; i<w1.size() ; i++)
201 wtot[i] = w1[i];
202 for( unsigned i=0; i<w2.size(); i++)
203 wtot[w1.size()+i] = w2[i];
204 for( unsigned i=0; i<w3.size(); i++)
205 wtot[w1.size()+w2.size()+i] = w3[i];
206 weights = wtot;
207
208 abscissas = detail::normalize_weights_and_compute_abscissas( g.grid(), weights);
209 }
210 virtual unsigned do_Ny_new( unsigned Ny, real_type fy) const override final{
211 if( fy==0 ) return Ny + 2*add_y_;
212 return Ny + 4*add_y_;
213 }
214 virtual unsigned do_Nx_new( unsigned Nx, real_type fx) const override final{
215 if( fx==0 ) return Nx + add_x_;
216 return Nx + 2*add_x_;
217 }
218 virtual real_type do_fx_new( unsigned Nx, real_type fx) const override final{
219 if( fx==0 ) return 0;
220 return (fx*(real_type)Nx + (real_type)add_x_)/(real_type)(Nx+2.*add_x_);
221 }
222 virtual real_type do_fy_new( unsigned Ny, real_type fy) const override final{
223 if( fy==0 ) return 0;
224 return (fy*(real_type)Ny + (real_type)add_y_)/(real_type)(Ny+4.*add_y_);
225 }
226};
227
233
234
239template<class real_type>
241{
243 real_type x0, real_type x1, real_type y0, real_type y1,
244 real_type fx, real_type fy,
245 unsigned n, unsigned Nx, unsigned Ny,
247 ref.fx_new(Nx, fx), ref.fy_new(Ny, fy), n, ref.nx_new(Nx, fx), ref.ny_new(Ny, fy), bcx, bcy), w_(2), abs_(2)
248 {
250 ref.generate(g,w_[0],w_[1],abs_[0],abs_[1]);
251 }
252
253 virtual RealCartesianRefinedGridX2d* clone()const override final{return new RealCartesianRefinedGridX2d(*this);}
254 private:
255 std::vector<thrust::host_vector<real_type> > w_,abs_;
256 virtual SparseTensor<thrust::host_vector<real_type> > do_compute_metric()const override final{
258 t.values().push_back( w_[0]);
259 t.values().push_back( w_[1]);
260 dg::blas1::pointwiseDot( w_[0], w_[0], t.values()[2]);
261 dg::blas1::pointwiseDot( w_[1], w_[1], t.values()[3]);
262 t.idx(0,0)=2, t.idx(1,1)=3;
263 return t;
264 }
265 virtual SparseTensor<thrust::host_vector<real_type> > do_compute_jacobian()const override final{
266 SparseTensor<thrust::host_vector<real_type> > t(*this);
267 t.values().push_back( w_[0]);
268 t.values().push_back( w_[1]);
269 t.idx(0,0)=2, t.idx(1,1)=3;
270 return t;
271 }
272 virtual std::vector<thrust::host_vector<real_type> > do_compute_map()const override final{
273 return abs_;
274 }
275};
276
281template<class real_type>
283{
285 real_type x0, real_type x1, real_type y0, real_type y1, real_type z0, real_type z1,
286 real_type fx, real_type fy,
287 unsigned n, unsigned Nx, unsigned Ny, unsigned Nz,
289 x0, x1, y0, y1,z0,z1,
290 ref.fx_new(Nx, fx), ref.fy_new(Ny, fy),
291 n, ref.nx_new(Nx, fx), ref.ny_new(Ny, fy), Nz,
292 bcx, bcy, bcz), w_(2), abs_(2)
293 {
295 ref.generate(g,w_[0],w_[1],abs_[0],abs_[1]);
296 //lift to 3d
297 w_[0].resize(this->size()), w_[1].resize(this->size()), abs_[0].resize(this->size()), abs_[1].resize(this->size());
298 unsigned size2d=this->n()*this->n()*this->Nx()*this->Ny();
299 for( unsigned i=1; i<Nz; i++)
300 for(unsigned k=0; k<size2d; k++)
301 {
302 w_[0][i*size2d+k]=w_[0][k];
303 w_[1][i*size2d+k]=w_[1][k];
304 abs_[0][i*size2d+k]=abs_[0][k];
305 abs_[1][i*size2d+k]=abs_[1][k];
306 }
307 }
308
310 private:
311 std::vector<thrust::host_vector<real_type> > w_,abs_;
312 virtual SparseTensor<thrust::host_vector<real_type> > do_compute_metric()const override final {
314 t.values().push_back( w_[0]);
315 t.values().push_back( w_[1]);
316 dg::blas1::pointwiseDot( w_[0], w_[0], t.values()[2]);
317 dg::blas1::pointwiseDot( w_[1], w_[1], t.values()[3]);
318 t.idx(0,0)=2, t.idx(1,1)=3;
319 return t;
320 }
321 virtual SparseTensor<thrust::host_vector<real_type> > do_compute_jacobian()const override final {
322 SparseTensor<thrust::host_vector<real_type> > t(*this);
323 t.values().push_back( w_[0]);
324 t.values().push_back( w_[1]);
325 t.idx(0,0)=2, t.idx(1,1)=3;
326 return t;
327 }
328 virtual std::vector<thrust::host_vector<real_type> > do_compute_map()const override final {
329 return abs_;
330 }
331};
332
338
339
340}//namespace dg
base X-point topology classes
void pointwiseDot(value_type alpha, const ContainerType1 &x1, const ContainerType2 &x2, value_type1 beta, ContainerType &y)
Definition blas1.h:406
bc
Switch between boundary conditions.
Definition enums.h:15
@ PER
periodic boundaries
Definition enums.h:16
@ DIR
homogeneous dirichlet boundaries
Definition enums.h:17
auto weights(const Topology &g)
Nodal weight coefficients.
Definition weights.h:62
dg::RealGrid< T, 1 > RealGrid1d
Definition grid.h:798
Interpolation matrix creation functions.
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin li...
Refined X-point grid.
Definition refined_gridX.h:241
RealCartesianRefinedGridX2d(const aRealRefinementX2d< real_type > &ref, real_type x0, real_type x1, real_type y0, real_type y1, real_type fx, real_type fy, unsigned n, unsigned Nx, unsigned Ny, bc bcx=dg::PER, bc bcy=dg::PER)
Definition refined_gridX.h:242
virtual RealCartesianRefinedGridX2d * clone() const override final
Geometries are cloneable.
Definition refined_gridX.h:253
Refined X-point grid.
Definition refined_gridX.h:283
RealCartesianRefinedGridX3d(const aRealRefinementX2d< real_type > &ref, real_type x0, real_type x1, real_type y0, real_type y1, real_type z0, real_type z1, real_type fx, real_type fy, unsigned n, unsigned Nx, unsigned Ny, unsigned Nz, bc bcx=dg::PER, bc bcy=dg::PER, bc bcz=dg::PER)
Definition refined_gridX.h:284
virtual RealCartesianRefinedGridX3d * clone() const
Geometries are cloneable.
Definition refined_gridX.h:309
Cell refinement around a given node.
Definition refined_grid.h:260
RealEquidistant cell refinement around the X-point.
Definition refined_gridX.h:114
virtual RealEquidistXRefinement * clone() const override final
Definition refined_gridX.h:116
RealEquidistXRefinement(unsigned add_x, unsigned add_y, unsigned howmanyX=1, unsigned howmanyY=1)
Definition refined_gridX.h:115
The exponential refinement around the X-point.
Definition refined_gridX.h:176
RealExponentialXRefinement(unsigned add_x, unsigned add_y, unsigned howmanyX=1, unsigned howmanyY=1)
Definition refined_gridX.h:177
virtual RealExponentialXRefinement * clone() const
Definition refined_gridX.h:178
The simplest implementation of aRealTopology.
Definition grid.h:710
1D grid for X-point topology
Definition gridX.h:68
The simplest implementation of aRealTopologyX2d.
Definition gridX.h:514
No refinement.
Definition refined_gridX.h:90
virtual RealIdentityXRefinement * clone() const override final
Definition refined_gridX.h:91
Class for 2x2 and 3x3 matrices sharing elements.
Definition tensor.h:51
This is the abstract interface class for a two-dimensional RealGeometryX.
Definition base_geometryX.h:18
This is the abstract interface class for a three-dimensional RealGeometryX.
Definition base_geometryX.h:63
Abstract base class for 2d grid refinement that increases the number of grid cells of a fixed basis g...
Definition refined_gridX.h:20
aRealRefinementX2d()=default
unsigned nx_new(unsigned Nx_old, real_type fx_old) const
the new number of cells
Definition refined_gridX.h:54
aRealRefinementX2d(const aRealRefinementX2d &src)=default
aRealRefinementX2d & operator=(const aRealRefinementX2d &src)=default
unsigned ny_new(unsigned Ny_old, real_type fy_old) const
Definition refined_gridX.h:58
real_type fy_new(unsigned Ny_old, real_type fy_old) const
Definition refined_gridX.h:66
void generate(const RealGridX2d< real_type > &g_old, thrust::host_vector< real_type > &weightsX, thrust::host_vector< real_type > &weightsY, thrust::host_vector< real_type > &abscissasX, thrust::host_vector< real_type > &abscissasY) const
Generate the grid transformation.
Definition refined_gridX.h:29
virtual ~aRealRefinementX2d()=default
real_type fx_new(unsigned Nx_old, real_type fx_old) const
Definition refined_gridX.h:62
virtual aRealRefinementX2d * clone() const =0
real_type y0() const
left boundary in y
Definition gridX.h:286
bc bcy() const
boundary conditions in y
Definition gridX.h:382
real_type y1() const
Right boundary in y.
Definition gridX.h:292
unsigned Nx() const
number of cells in x
Definition gridX.h:340
real_type x1() const
Right boundary in x.
Definition gridX.h:280
real_type x0() const
Left boundary in x.
Definition gridX.h:274
real_type fx() const
partition factor in x
Definition gridX.h:322
bc bcx() const
boundary conditions in x
Definition gridX.h:376
unsigned inner_Nx() const
number of topological cells in x
Definition gridX.h:346
unsigned Ny() const
number of cells in y
Definition gridX.h:358
real_type fy() const
partition factor in y
Definition gridX.h:328
unsigned n() const
number of polynomial coefficients in x and y
Definition gridX.h:334
bc bcx() const
boundary conditions in x
Definition gridX.h:693
unsigned n() const
number of polynomial coefficients in x and y
Definition gridX.h:645
real_type y1() const
right boundary in y
Definition gridX.h:576
real_type z1() const
right boundary in z
Definition gridX.h:589
real_type fx() const
partition factor in x
Definition gridX.h:633
unsigned Nz() const
number of points in z
Definition gridX.h:687
unsigned Ny() const
number of cells in y
Definition gridX.h:669
unsigned Nx() const
number of points in x
Definition gridX.h:651
bc bcy() const
boundary conditions in y
Definition gridX.h:699
unsigned size() const
real_typehe total number of points
Definition gridX.h:719
real_type z0() const
left boundary in z
Definition gridX.h:583
real_type y0() const
left boundary in y
Definition gridX.h:570
real_type x0() const
left boundary in x
Definition gridX.h:557
real_type x1() const
right boundary in x
Definition gridX.h:563
real_type fy() const
partition factor in y
Definition gridX.h:639
bc bcz() const
boundary conditions in z
Definition gridX.h:705