43 #ifndef MFEM_THREAD_SAFE
50 Diffusion2Integrator() { Q = NULL; }
65 #ifdef MFEM_THREAD_SAFE
78 if (el.
Space() == FunctionSpace::Pk)
87 if (el.
Space() == FunctionSpace::rQk)
109 w *= Q->Eval(Trans, ip);
112 for (
int j = 0; j < nd; j++)
114 for (
int i = 0; i < nd; i++)
116 elmat(i, j) += w*shape(i)*laplace(j);
124 int main(
int argc,
char *argv[])
127 const char *mesh_file =
"../../data/star.mesh";
128 const char *per_file =
"none";
132 bool static_cond =
false;
133 bool visualization = 1;
141 args.
AddOption(&mesh_file,
"-m",
"--mesh",
142 "Mesh file to use.");
143 args.
AddOption(&ref_levels,
"-r",
"--refine",
144 "Number of times to refine the mesh uniformly, -1 for auto.");
146 "Periodic BCS file.");
147 args.
AddOption(&master,
"-pm",
"--master",
148 "Master boundaries for periodic BCs");
150 "Slave boundaries for periodic BCs");
152 "Finite element order (polynomial degree) or -1 for"
153 " isoparametric space.");
154 args.
AddOption(&ibp,
"-ibp",
"--ibp",
"-no-ibp",
156 "Selects the standard weak form (IBP) or the nonstandard (NO-IBP).");
157 args.
AddOption(&strongBC,
"-sbc",
"--strong-bc",
"-wbc",
159 "Selects strong or weak enforcement of Dirichlet BCs.");
161 "Sets the SIPG penalty parameters, should be positive."
162 " Negative values are replaced with (order+1)^2.");
163 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
164 "--no-static-condensation",
"Enable static condensation.");
165 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
166 "--no-visualization",
167 "Enable or disable GLVis visualization.");
174 if (!strongBC & (kappa < 0))
176 kappa = 4*(order.
Max()+1)*(order.
Max()+1);
183 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
194 (int)floor(log(5000./mesh->
GetNE())/log(2.)/
dim);
197 for (
int l = 0; l < ref_levels; l++)
217 if (order.
Size() == 1)
224 if (order.
Size() != nkv ) {
mfem_error(
"Wrong number of orders set."); }
229 in.open(per_file, std::ifstream::in);
236 master.
Load(in, psize);
237 slave.
Load(in, psize);
244 else if (order[0] == -1)
250 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
254 cout <<
"Mesh does not have FEs --> Assume order 1.\n";
261 if (order.
Size() > 1) { cout <<
"Wrong number of orders set, needs one.\n"; }
267 cout <<
"Number of finite element unknowns: "
274 cout <<
"No integration by parts requires a NURBS mesh."<< endl;
279 cout <<
"No integration by parts requires a NURBS mesh, with only 1 patch."<<
281 cout <<
"A C_1 discretisation is required."<< endl;
282 cout <<
"Currently only C_0 multipatch coupling implemented."<< endl;
287 cout <<
"No integration by parts requires at least quadratic NURBS."<< endl;
288 cout <<
"A C_1 discretisation is required."<< endl;
311 for (
int i = 0; i < master.
Size(); i++)
313 ess_bdr[master[i]-1] = 0;
314 ess_bdr[slave[i]-1] = 0;
367 cout <<
"Size of linear system: " << A.
Height() << endl;
369 #ifndef MFEM_USE_SUITESPARSE
373 PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
377 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
379 umf_solver.
Mult(B, X);
387 ofstream mesh_ofs(
"refined.mesh");
388 mesh_ofs.precision(8);
389 mesh->
Print(mesh_ofs);
390 ofstream sol_ofs(
"sol.gf");
391 sol_ofs.precision(8);
400 sol_sock.precision(8);
401 sol_sock <<
"solution\n" << *mesh << x << flush;
413 if (own_fec) {
delete fec; }
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
int Size() const
Return the logical size of the array.
Class for domain integration L(v) := (f, v)
virtual void Print(std::ostream &out=mfem::out) const
int GetDim() const
Returns the reference space dimension for the finite element.
void Load(std::istream &in, int fmt=0)
Read an Array from the stream in using format fmt. The format fmt can be:
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
A coefficient that is constant across space and time.
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...
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Data type dense matrix using column-major storage.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
int GetNE() const
Returns number of elements.
virtual void Save()
Save the collection and a VisIt root file.
int Space() const
Returns the type of FunctionSpace on the element.
int main(int argc, char *argv[])
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
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...
virtual void PrintInfo(std::ostream &out=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
void PrintUsage(std::ostream &out) const
Print the usage message.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
double Control[UMFPACK_CONTROL]
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
virtual const char * Name() const
Class for integration point with weight.
void PrintOptions(std::ostream &out) const
Print the options.
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
bool Good() const
Return true if the command line options were parsed successfully.