Alamo
Stencil.H
Go to the documentation of this file.
1#ifndef TEST_NUMERIC_STENCIL
2#define TEST_NUMERIC_STENCIL
3
4#include <AMReX.H>
5
6#include "Set/Set.H"
7#include "IC/Trig2.H"
8#include "Numeric/Stencil.H"
9
10namespace Test
11{
12/// Tests for the Numeric namespace classes
13namespace Numeric
14{
16{
17
18public:
19 Stencil() {};
21
22 void Define(amrex::IntVect _ncells)
23 {
24 ncells = _ncells;
25 amrex::RealBox rb({AMREX_D_DECL(0.,0.,0.)}, {AMREX_D_DECL(L,L,L)});
26 amrex::Geometry::Setup(&rb, 0);
27
28 amrex::Box domain(amrex::IntVect{AMREX_D_DECL(0,0,0)}, ncells,
29 amrex::IntVect::TheCellVector());
30 geom.resize(1);
31 grids.resize(1);
32 dmap.resize(1);
33 phi.resize(1);
34 DphiNumeric.resize(1);
35 DphiExact.resize(1);
36
37 geom[0].define(domain);
38 grids[0].define(domain);
39 grids[0].maxSize(ncells[0]/2);
40 dmap[0].define(grids[0]);
41
42 phi .Define(0,grids[0],dmap[0],1,nghost);
44 DphiExact .Define(0,grids[0],dmap[0],1,nghost);
45
46 phi[0]->setVal(0.0);
47 // Initial function is
48 // cos(2*pi*x/L) * cos(2*pi*y/L) * cos(2*pi*z / L)
49 IC::Trig2 ic(geom,1.0,AMREX_D_DECL(0.0,0.0,0.0),AMREX_D_DECL(1,1,1));
50 ic.Add(0,phi,0.0);
51 }
52 void Define(int _ncells)
53 {
54 Define(amrex::IntVect{AMREX_D_DECL(_ncells,_ncells,_ncells)});
55 }
56
57 template<int i,int j,int k>
58 bool Derivative(int verbose,std::string plotfile = "")
59 {
60 const Set::Scalar tolerance = 1E-2;
61
62 DphiExact[0]->setVal(0.0);
63 DphiNumeric[0]->setVal(0.0);
64
65 //
66 // Compute exact derivative
67 //
68 Set::Scalar AMREX_D_DECL(theta1 = 0.0, theta2 = 0.0, theta3 = 0.0);
69 Set::Scalar fac = 1.0;
70 AMREX_D_TERM( for (int p = 0; p < i; p++) theta1 += 0.5*Set::Constant::Pi;,
71 for (int q = 0; q < j; q++) theta2 += 0.5*Set::Constant::Pi;,
72 for (int r = 0; r < k; r++) theta3 += 0.5*Set::Constant::Pi; );
73
74 AMREX_D_TERM( for (int p = 0; p < i; p++) fac *= 1.0*Set::Constant::Pi / L;,
75 for (int q = 0; q < j; q++) fac *= 1.0*Set::Constant::Pi / L;,
76 for (int r = 0; r < k; r++) fac *= 1.0*Set::Constant::Pi / L;);
77 IC::Trig2 ic(geom,fac,AMREX_D_DECL(theta1,theta2,theta3),AMREX_D_DECL(1,1,1));
78 ic.Add(0,DphiExact,0.0);
79
80
81 //
82 // Compute numeric derivative
83 //
84 const amrex::Real* DX = geom[0].CellSize();
85 for ( amrex::MFIter mfi(*phi[0],amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi )
86 {
87 const amrex::Box& bx = mfi.tilebox();
88 amrex::Array4<const amrex::Real> const& Phi = phi[0]->array(mfi);
89 amrex::Array4<amrex::Real> const& DPhi = DphiNumeric[0]->array(mfi);
90 amrex::ParallelFor (bx,[=] AMREX_GPU_DEVICE(int p, int q, int r){ DPhi(p,q,r) = ::Numeric::Stencil<Set::Scalar,i,j,k>::D(Phi,p,q,r,0,DX); });
91 }
92
93 amrex::MultiFab diff(grids[0],dmap[0],ncomp,nghost);
94 amrex::MultiFab::Copy (diff,*DphiExact[0] ,0,0,ncomp,nghost);
95 amrex::MultiFab::Subtract(diff,*DphiNumeric[0],0,0,ncomp,nghost);
96
97 Set::Scalar error = diff.norm0(0,0,false) / DphiExact[0]->norm0(0,0,false);
98
99 if (verbose) Util::Message(INFO,"Infinity norm of error = ", 100*error,"%");
100
101 if (plotfile != "")
102 WritePlotFile(plotfile);
103
104 if (error < tolerance) return 0;
105 else return 1;
106 }
107
108
109 void WritePlotFile(std::string plotfile)
110 {
111 amrex::Vector<amrex::MultiFab> plotmf(1);
112 plotmf[0].define(grids[0], dmap[0], 3, 0);
113 for (amrex::MFIter mfi(plotmf[0],false); mfi.isValid(); ++mfi)
114 {
115 amrex::Box bx = mfi.validbox();
116 plotmf[0][mfi].copy((*phi[0])[mfi] ,bx, 0, bx, 0, 1);
117 plotmf[0][mfi].copy((*DphiExact[0])[mfi] ,bx, 0, bx, 1, 1);
118 plotmf[0][mfi].copy((*DphiNumeric[0])[mfi],bx, 0, bx, 2, 1);
119 }
120 amrex::WriteSingleLevelPlotfile(plotfile,plotmf[0],varnames,geom[0],0.0,1);
121 }
122
123
124private:
125 amrex::IntVect ncells;
126 const Set::Scalar L = 1.0;
127 const int ncomp = 1, nghost = 2;
128 amrex::Vector<amrex::Geometry> geom;
129 amrex::Vector<amrex::BoxArray> grids;
130 amrex::Vector<amrex::DistributionMapping> dmap;
131
135 amrex::Vector<std::string> varnames = {"phi","dphi_exact","dphi_numeric"};
136};
137}
138}
139
140
141
142#endif
#define INFO
Definition Util.H:20
Initialize using a trigonometric series.
Definition Trig2.H:17
virtual void Add(const int &lev, Set::Field< Set::Scalar > &field, Set::Scalar)
Definition Trig2.H:53
void Define(int a_levs, const amrex::Vector< amrex::BoxArray > &a_grids, const amrex::Vector< amrex::DistributionMapping > &a_dmap, int a_ncomp, int a_nghost)
Definition Set.H:241
const Set::Scalar L
Definition Stencil.H:126
amrex::Vector< amrex::DistributionMapping > dmap
Definition Stencil.H:130
Set::Field< Set::Scalar > DphiNumeric
Definition Stencil.H:134
Set::Field< Set::Scalar > phi
Definition Stencil.H:132
void WritePlotFile(std::string plotfile)
Definition Stencil.H:109
void Define(int _ncells)
Definition Stencil.H:52
Set::Field< Set::Scalar > DphiExact
Definition Stencil.H:133
amrex::Vector< amrex::BoxArray > grids
Definition Stencil.H:129
bool Derivative(int verbose, std::string plotfile="")
Definition Stencil.H:58
amrex::IntVect ncells
Definition Stencil.H:125
void Define(amrex::IntVect _ncells)
Definition Stencil.H:22
amrex::Vector< std::string > varnames
Definition Stencil.H:135
amrex::Vector< amrex::Geometry > geom
Definition Stencil.H:128
This namespace contains some numerical tools.
Definition Function.H:6
static const Set::Scalar Pi
Definition Set.H:286
amrex::Real Scalar
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:20
Definition GB.H:8
void Message(std::string file, std::string func, int line, Args const &... args)
Definition Util.H:141