57#ifndef MFEM_THREAD_SAFE
64 Diffusion2Integrator() { Q = NULL; }
79#ifdef MFEM_THREAD_SAFE
123 w *= Q->
Eval(Trans, ip);
126 for (
int jj = 0; jj < nd; jj++)
128 for (
int ii = 0; ii < nd; ii++)
130 elmat(ii, jj) += w*shape(ii)*laplace(jj);
138int main(
int argc,
char *argv[])
147 const char *mesh_file =
"../../data/star.mesh";
151 bool static_cond =
false;
152 bool visualization = 1;
160 args.
AddOption(&mesh_file,
"-m",
"--mesh",
161 "Mesh file to use.");
162 args.
AddOption(&ref_levels,
"-r",
"--refine",
163 "Number of times to refine the mesh uniformly, -1 for auto.");
164 args.
AddOption(&master,
"-pm",
"--master",
165 "Master boundaries for periodic BCs");
167 "Slave boundaries for periodic BCs");
169 "Finite element order (polynomial degree) or -1 for"
170 " isoparametric space.");
171 args.
AddOption(&ibp,
"-ibp",
"--ibp",
"-no-ibp",
173 "Selects the standard weak form (IBP) or the nonstandard (NO-IBP).");
174 args.
AddOption(&strongBC,
"-sbc",
"--strong-bc",
"-wbc",
176 "Selects strong or weak enforcement of Dirichlet BCs.");
178 "Sets the SIPG penalty parameters, should be positive."
179 " Negative values are replaced with (order+1)^2.");
180 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
181 "--no-static-condensation",
"Enable static condensation.");
182 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
183 "--no-visualization",
184 "Enable or disable GLVis visualization.");
194 if (!strongBC & (
kappa < 0))
206 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
217 (int)floor(log(5000./mesh->
GetNE())/log(2.)/
dim);
220 for (
int l = 0; l < ref_levels; l++)
237 int par_ref_levels = 2;
238 for (
int l = 0; l < par_ref_levels; l++)
257 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
261 cout <<
"Mesh does not have FEs --> Assume order 1.\n";
266 else if (pmesh->
NURBSext && (order[0] > 0) )
272 if (order.
Size() == 1)
278 if (order.
Size() != nkv ) {
mfem_error(
"Wrong number of orders set."); }
282 if (master.
Size() > 0)
286 cout<<
"Connecting boundaries"<<endl;
287 cout<<
" - master : "; master.
Print();
288 cout<<
" - slave : "; slave.
Print();
296 if (order.
Size() > 1) { cout <<
"Wrong number of orders set, needs one.\n"; }
304 cout <<
"Number of finite element unknowns: " << size << endl;
311 cout <<
"No integration by parts requires a NURBS mesh."<< endl;
316 cout <<
"No integration by parts requires a NURBS mesh, with only 1 patch."<<
318 cout <<
"A C_1 discretisation is required."<< endl;
319 cout <<
"Currently only C_0 multipatch coupling implemented."<< endl;
324 cout <<
"No integration by parts requires at least quadratic NURBS."<< endl;
325 cout <<
"A C_1 discretisation is required."<< endl;
348 for (
int i = 0; i < master.
Size(); i++)
350 ess_bdr[master[i]-1] = 0;
351 ess_bdr[slave[i]-1] = 0;
367 b->AddBdrFaceIntegrator(
387 a->AddDomainIntegrator(
new Diffusion2Integrator(one));
398 if (static_cond) {
a->EnableStaticCondensation(); }
403 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
422 a->RecoverFEMSolution(X, *
b, x);
427 ostringstream mesh_name, sol_name;
428 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
429 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
431 ofstream mesh_ofs(mesh_name.str().c_str());
432 mesh_ofs.precision(8);
433 pmesh->
Print(mesh_ofs);
435 ofstream sol_ofs(sol_name.str().c_str());
436 sol_ofs.precision(8);
446 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
447 sol_sock.precision(8);
448 sol_sock <<
"solution\n" << *pmesh << x << flush;
462 if (own_fec) {
delete fec; }
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip)=0
Evaluate the coefficient in the element described by T at the point ip.
A coefficient that is constant across space and time.
Data type dense matrix using column-major storage.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Class for domain integration .
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
virtual const char * Name() const
Abstract class for all finite elements.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
int GetDim() const
Returns the reference space dimension for the finite element.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
int Space() const
Returns the type of FunctionSpace on the element.
virtual void CalcPhysLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Evaluate the Laplacian of all shape functions of a scalar finite element in reference space at the gi...
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
@ Pk
Polynomials of order k.
@ rQk
Refined tensor products of polynomials of order k.
Arbitrary order H1-conforming (continuous) finite elements.
The BoomerAMG solver in hypre.
void SetPrintLevel(int print_lvl)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
void SetMaxIter(int max_iter)
Wrapper for hypre's ParCSR matrix class.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Abstract class for hypre's solvers and preconditioners.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetNodes(Vector &node_coord) const
virtual void PrintInfo(std::ostream &os=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
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).
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
Class for parallel grid function.
void Save(std::ostream &out) const override
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
void SetSize(int s)
Resize the vector to size s.
Data collection with VisIt I/O routines.
virtual void Save() override
Save the collection and a VisIt root file.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
void mfem_error(const char *msg)
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)