42 char *title,
int pos_x,
int pos_y);
44int main (
int argc,
char *argv[])
51 const char *mesh_file =
"../gslib/triple-pt-1.mesh";
52 const char *sltn_file =
"../gslib/triple-pt-1.gf";
54 bool visualization =
true;
57 bool continuous =
true;
62 args.
AddOption(&mesh_file,
"-m",
"--mesh",
64 args.
AddOption(&sltn_file,
"-s",
"--sltn",
65 "Solution file to use.");
66 args.
AddOption(&ref,
"-ref",
"--piecewise-linear-ref-factor",
67 "Scaling factor for resolution of piecewise linear bounds."
68 " If less than 2, the resolution is picked automatically");
69 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
71 "Enable or disable GLVis visualization.");
72 args.
AddOption(&visit,
"-visit",
"--visit",
"-no-visit",
74 "Enable or disable VisIt output.");
75 args.
AddOption(&b_type,
"-bt",
"--basis-type",
76 "Project input function to a different bases. "
77 "-1 = don't project (default)."
78 "0 = Gauss-Legendre nodes. "
79 "1 = Gauss-Lobatto nodes. "
80 "2 = uniformly spaced nodes. ");
81 args.
AddOption(&continuous,
"-h1",
"--h1",
"-l2",
"--l2",
82 "Use continuous or discontinuous space.");
83 args.
AddOption(&nbrute,
"-nb",
"--nbrute",
84 "Brute force search for minimum in an array of nxnxn points "
88 Mesh mesh(mesh_file, 1, 1,
false);
90 if (continuous && b_type != -1)
92 MFEM_VERIFY(b_type > 0,
"Continuous space do not support GL nodes. "
93 "Please use basis type: 1 for Lagrange interpolants on GLL "
94 " nodes 2 for positive bases on uniformly spaced nodes.");
97 std::unique_ptr<int[]> partition(
101 ifstream mat_stream_1(sltn_file);
102 std::unique_ptr<GridFunction> func(
new GridFunction(&mesh, mat_stream_1));
104 ParMesh pmesh(MPI_COMM_WORLD, mesh, partition.get());
106 int func_order = func->FESpace()->GetMaxElementOrder();
108 int nel = pmesh.
GetNE();
137 cout <<
"fec name: " << fec->
Name() << endl;
156 Vector bound_min(vdim), bound_max(vdim);
157 for (
int d = 0; d < vdim; d++)
159 Vector lowerT(lowerb.GetData() + d*nel, nel);
161 bound_min(d) = lowerT.
Min();
162 bound_max(d) = upperT.
Max();
165 MPI_Allreduce(MPI_IN_PLACE, bound_min.GetData(), vdim,
167 MPI_Allreduce(MPI_IN_PLACE, bound_max.
GetData(), vdim,
173 char title1[] =
"Input gridfunction";
177 char title1p[] =
"Projected gridfunction";
180 char title2[] =
"Element-wise lower bound";
182 char title3[] =
"Element-wise upper bound";
203 Vector global_min(vdim), global_max(vdim);
204 global_min = numeric_limits<real_t>::max();
205 global_max = numeric_limits<real_t>::min();
208 for (
int e = 0; e < pmesh.
GetNE(); e++)
211 for (
int k = 0; k < (
dim > 2 ? nbrute : 1); k++)
213 ip.
z = k/(nbrute-1.0);
214 for (
int j = 0; j < (
dim > 1 ? nbrute : 1); j++)
216 ip.
y = j/(nbrute-1.0);
217 for (
int i = 0; i < nbrute; i++)
219 ip.
x = i/(nbrute-1.0);
220 for (
int d = 0; d < vdim; d++)
223 global_min(d) = min(global_min(d), val);
224 global_max(d) = max(global_max(d), val);
231 MPI_Allreduce(MPI_IN_PLACE, global_min.GetData(), vdim,
233 MPI_Allreduce(MPI_IN_PLACE, global_max.
GetData(), vdim,
237 for (
int d = 0; d < vdim; d++)
239 cout <<
"Brute force and bounding comparison for component " <<
241 cout <<
"Brute force minimum and minimum bound: " << global_min(d)
242 <<
" " << bound_min(d) << endl;
244 cout <<
"Brute force maximum and maximum bound: " << global_max(d)
245 <<
" " << bound_max(d) << endl;
247 cout <<
"The difference in bounds is: " <<
248 global_min(d)-bound_min(d) <<
" " <<
249 bound_max(d)-global_max(d) << endl;
256 for (
int d = 0; d < vdim; d++)
258 cout <<
"Minimum bound for component " << d <<
" is " <<
259 bound_min(d) << endl;
260 cout <<
"Maximum bound for component " << d <<
" is " <<
261 bound_max(d) << endl;
273 char *title,
int pos_x,
int pos_y)
278 sock.
open(
"localhost", 19916);
279 sock <<
"solution\n";
285 sock <<
"window_title '"<< title <<
"'\n"
286 <<
"window_geometry "
287 << pos_x <<
" " << pos_y <<
" " << 400 <<
" " << 400 <<
"\n"
288 <<
"keys jRmclApppppppppppp//]]]]]]]]" << endl;
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
virtual const char * Name() const
Ordering::Type GetOrdering() const
Return the ordering method.
const FiniteElementCollection * FEColl() const
int GetVDim() const
Returns the vector dimension of the finite element space.
Class for grid function - Vector with associated FE space.
PLBound GetElementBounds(Vector &lower, Vector &upper, const int ref_factor=1, const int vdim=-1) const
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec_owned and fes.
FiniteElementSpace * FESpace()
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh.
Arbitrary order H1-conforming (continuous) finite elements.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Class for integration point with weight.
Arbitrary order "L2-conforming" discontinuous finite elements.
void Clear()
Clear the contents of the Mesh.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
int * GeneratePartitioning(int nparts, int part_method=1)
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
void ParseCheck(std::ostream &out=mfem::out)
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
Abstract parallel finite element space.
Class for parallel grid function.
real_t GetValue(int i, const IntegrationPoint &ip, int vdim=1) const override
void SaveAsOne(const char *fname, int precision=16) const
Class for parallel meshes.
void PrintAsOne(std::ostream &out=mfem::out, const std::string &comments="") const
real_t Max() const
Returns the maximal element of the vector.
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
real_t Min() const
Returns the minimal element of the vector.
Data collection with VisIt I/O routines.
void Save() override
Save the collection and a VisIt root file.
void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
Helper struct to convert a C++ type to an MPI type.