41 #ifndef MFEM_THREAD_SAFE
48 Diffusion2Integrator() { Q = NULL; }
63 #ifdef MFEM_THREAD_SAFE
76 if (el.
Space() == FunctionSpace::Pk)
85 if (el.
Space() == FunctionSpace::rQk)
107 w *= Q->Eval(Trans, ip);
110 for (
int j = 0; j < nd; j++)
112 for (
int i = 0; i < nd; i++)
114 elmat(i, j) += w*shape(i)*laplace(j);
122 int main(
int argc,
char *argv[])
125 const char *mesh_file =
"../../data/star.mesh";
126 const char *per_file =
"none";
129 bool static_cond =
false;
130 bool visualization = 1;
136 args.
AddOption(&mesh_file,
"-m",
"--mesh",
137 "Mesh file to use.");
139 "Periodic BCS file.");
140 args.
AddOption(&master,
"-pm",
"--master",
141 "Master boundaries for periodic BCs");
143 "Slave boundaries for periodic BCs");
145 "Finite element order (polynomial degree) or -1 for"
146 " isoparametric space.");
147 args.
AddOption(&ibp,
"-ibp",
"--ibp",
"-no-ibp",
149 "Selects the standard weak form (IBP) or the nonstandard (NO-IBP).");
150 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
151 "--no-static-condensation",
"Enable static condensation.");
152 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
153 "--no-visualization",
154 "Enable or disable GLVis visualization.");
166 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
175 (int)floor(log(5000./mesh->
GetNE())/log(2.)/
dim);
176 for (
int l = 0; l < ref_levels; l++)
195 if (order.
Size() == 1)
202 if (order.
Size() != nkv ) {
mfem_error(
"Wrong number of orders set."); }
207 in.open(per_file, std::ifstream::in);
214 master.
Load(in, psize);
215 slave.
Load(in, psize);
222 else if (order[0] == -1)
228 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
232 cout <<
"Mesh does not have FEs --> Assume order 1.\n";
239 if (order.
Size() > 1) { cout <<
"Wrong number of orders set, needs one.\n"; }
245 cout <<
"Number of finite element unknowns: "
252 cout <<
"No integration by parts requires a NURBS mesh."<< endl;
257 cout <<
"No integration by parts requires a NURBS mesh, with only 1 patch."<<
259 cout <<
"A C_1 discretisation is required."<< endl;
260 cout <<
"Currently only C_0 multipatch coupling implemented."<< endl;
265 cout <<
"No integration by parts requires at least quadratic NURBS."<< endl;
266 cout <<
"A C_1 discretisation is required."<< endl;
283 for (
int i = 0; i < master.
Size(); i++)
285 ess_bdr[master[i]-1] = 0;
286 ess_bdr[slave[i]-1] = 0;
329 cout <<
"Size of linear system: " << A.
Height() << endl;
331 #ifndef MFEM_USE_SUITESPARSE
335 PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
339 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
341 umf_solver.
Mult(B, X);
349 ofstream mesh_ofs(
"refined.mesh");
350 mesh_ofs.precision(8);
351 mesh->
Print(mesh_ofs);
352 ofstream sol_ofs(
"sol.gf");
353 sol_ofs.precision(8);
359 char vishost[] =
"localhost";
362 sol_sock.precision(8);
363 sol_sock <<
"solution\n" << *mesh << x << flush;
375 if (own_fec) {
delete fec; }
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for Finite Elements.
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
int Size() const
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.
Subclass constant coefficient.
void SetSize(int s)
Resize the vector to size s.
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)
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 space on each 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.
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...
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
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 Coefficient that may optionally depend on time.
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)
void SetSize(int nsize)
Change 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
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.