44 int main(
int argc,
char *argv[])
48 MPI_Init(&argc, &argv);
49 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
50 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
53 const char *mesh_file =
"../data/beam-tri.mesh";
55 bool visualization = 1;
58 args.
AddOption(&mesh_file,
"-m",
"--mesh",
61 "Finite element order (polynomial degree).");
62 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
64 "Enable or disable GLVis visualization.");
80 ifstream imesh(mesh_file);
84 cerr <<
"\nCan not open mesh file: " << mesh_file <<
'\n' << endl;
88 mesh =
new Mesh(imesh, 1, 1);
95 cerr <<
"\nInput mesh should have at least two materials and "
96 <<
"two boundary attributes! (See schematic in ex2.cpp)\n"
113 (int)floor(log(1000./mesh->
GetNE())/log(2.)/dim);
114 for (
int l = 0; l < ref_levels; l++)
124 int par_ref_levels = 1;
125 for (
int l = 0; l < par_ref_levels; l++)
149 cout <<
"Number of unknowns: " << size << endl
150 <<
"Assembling: " << flush;
161 for (
int i = 0; i < dim-1; i++)
166 pull_force(1) = -1.0e-2;
173 cout <<
"r.h.s. ... " << flush;
189 lambda(0) = lambda(1)*50;
199 cout <<
"matrix ... " << flush;
207 cout <<
"done." << endl;
251 ostringstream mesh_name, sol_name;
252 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
253 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
255 ofstream mesh_ofs(mesh_name.str().c_str());
256 mesh_ofs.precision(8);
257 pmesh->
Print(mesh_ofs);
259 ofstream sol_ofs(sol_name.str().c_str());
260 sol_ofs.precision(8);
268 char vishost[] =
"localhost";
271 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
272 sol_sock.precision(8);
273 sol_sock <<
"solution\n" << *pmesh << x << flush;
Vector coefficient defined by an array of scalar coefficients.
Class for grid function - Vector with associated FE space.
Subclass constant coefficient.
void DegreeElevate(int t)
int GetNE() const
Returns number of elements.
virtual void Save(std::ostream &out) const
Abstract parallel finite element space.
void SetPrintLevel(int print_lvl)
void SetSystemsOptions(int dim)
The BoomerAMG solver in hypre.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetNodalFESpace(FiniteElementSpace *nfes)
void PrintUsage(std::ostream &out) const
void SetMaxIter(int max_iter)
Wrapper for hypre's parallel vector class.
Array< int > bdr_attributes
int main(int argc, char *argv[])
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)
NURBSExtension * NURBSext
class for piecewise constant coefficient
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void PrintOptions(std::ostream &out) const
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
Class for parallel grid function.
Wrapper for hypre's ParCSR matrix class.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Class for parallel meshes.
void ParallelAverage(Vector &tv) const
Returns the vector averaged on the true dofs.
void Set(int i, Coefficient *c)
Sets coefficient in the vector.
virtual void Print(std::ostream &out=std::cout) const