38 using std::endl;
using std::cout;
using std::cerr;
39 using std::ends;
using std::cin;
44 using bgeot::scalar_type;
51 typedef gmm::row_matrix<sparse_vector_type> sparse_matrix_type;
52 typedef gmm::col_matrix<sparse_vector_type> col_sparse_matrix_type;
53 typedef std::vector<scalar_type> plain_vector;
59 base_small_vector sol_K;
61 scalar_type sol_u(
const base_node &x) {
return sin(gmm::vect_sp(sol_K, x)); }
63 scalar_type sol_f(
const base_node &x)
64 {
return gmm::vect_sp(sol_K, sol_K) * sin(gmm::vect_sp(sol_K, x)); }
66 base_small_vector sol_grad(
const base_node &x)
67 {
return sol_K * cos(gmm::vect_sp(sol_K, x)); }
72 struct laplacian_problem {
74 enum { DIRICHLET_BOUNDARY_NUM = 0, NEUMANN_BOUNDARY_NUM = 1};
85 sparse_matrix_type SM;
86 std::vector<scalar_type> U, B;
88 std::vector<scalar_type> Ud;
89 col_sparse_matrix_type NS;
92 std::string datafilename;
93 bgeot::md_param PARAM;
99 laplacian_problem(
void) : mim(mesh), mf_u(mesh), mf_rhs(mesh),
106 void laplacian_problem::init(
void) {
108 std::string MESH_TYPE = PARAM.string_value(
"MESH_TYPE",
"Mesh type ");
109 std::string FEM_TYPE = PARAM.string_value(
"FEM_TYPE",
"FEM name");
110 std::string INTEGRATION = PARAM.string_value(
"INTEGRATION",
111 "Name of integration method");
113 cout <<
"MESH_TYPE=" << MESH_TYPE <<
"\n";
114 cout <<
"FEM_TYPE=" << FEM_TYPE <<
"\n";
115 cout <<
"INTEGRATION=" << INTEGRATION <<
"\n";
121 std::vector<size_type> nsubdiv(N);
122 std::fill(nsubdiv.begin(),nsubdiv.end(),
123 PARAM.int_value(
"NX",
"Nomber of space steps "));
125 PARAM.int_value(
"MESH_NOISED") != 0);
127 bgeot::base_matrix M(N,N);
129 static const char *t[] = {
"LX",
"LY",
"LZ"};
130 M(i,i) = (i<3) ? PARAM.real_value(t[i],t[i]) : 1.0;
132 if (N>1) { M(0,1) = PARAM.real_value(
"INCLINE") * PARAM.real_value(
"LY"); }
135 mesh.transformation(M);
137 datafilename = PARAM.string_value(
"ROOTFILENAME",
"Base name of data files.");
138 scalar_type FT = PARAM.real_value(
"FT",
"parameter for exact solution");
139 residual = PARAM.real_value(
"RESIDUAL");
140 if (residual == 0.) residual = 1e-10;
143 sol_K[j] = ((j & 1) == 0) ? FT : -FT;
149 mim.set_integration_method(mesh.convex_index(), ppi);
150 mf_u.set_finite_element(mesh.convex_index(), pf_u);
154 std::string data_fem_name = PARAM.string_value(
"DATA_FEM_TYPE");
155 if (data_fem_name.size() == 0) {
156 GMM_ASSERT1(pf_u->is_lagrange(),
"You are using a non-lagrange FEM. "
157 <<
"In that case you need to set "
158 <<
"DATA_FEM_TYPE in the .param file");
159 mf_rhs.set_finite_element(mesh.convex_index(), pf_u);
161 mf_rhs.set_finite_element(mesh.convex_index(),
168 mf_coef.set_finite_element(mesh.convex_index(),
173 gen_dirichlet = PARAM.int_value(
"GENERIC_DIRICHLET");
175 if (!pf_u->is_lagrange() && !gen_dirichlet)
176 GMM_WARNING2(
"With non lagrange fem prefer the generic "
177 "Dirichlet condition option");
179 cout <<
"Selecting Neumann and Dirichlet boundaries\n";
184 base_node un = mesh.normal_of_face_of_convex(i.cv(), i.f());
186 if (gmm::abs(un[N-1] - 1.0) < 1.0E-7) {
187 mesh.region(NEUMANN_BOUNDARY_NUM).add(i.cv(), i.f());
189 mesh.region(DIRICHLET_BOUNDARY_NUM).add(i.cv(), i.f());
194 void laplacian_problem::assembly(
void) {
202 cout <<
"Number of dof : " << nb_dof << endl;
203 cout <<
"Assembling stiffness matrix" << endl;
205 std::vector<scalar_type>(mf_coef.nb_dof(), 1.0));
207 cout <<
"Assembling source term" << endl;
208 std::vector<scalar_type> F(nb_dof_rhs);
212 cout <<
"Assembling Neumann condition" << endl;
216 NEUMANN_BOUNDARY_NUM);
218 cout <<
"take Dirichlet condition into account" << endl;
219 if (!gen_dirichlet) {
220 std::vector<scalar_type> D(nb_dof);
222 getfem::assembling_Dirichlet_condition(SM, B, mf_u,
223 DIRICHLET_BOUNDARY_NUM, D);
230 col_sparse_matrix_type H(nb_dof_rhs, nb_dof);
231 std::vector<scalar_type> R(nb_dof_rhs);
232 std::vector<scalar_type> RHaux(nb_dof);
236 mf_rhs, F, DIRICHLET_BOUNDARY_NUM);
245 gmm::mult(SM, Ud, gmm::scaled(B, -1.0), RHaux);
248 gmm::mult(gmm::transposed(NS), gmm::scaled(RHaux, -1.0), B);
249 sparse_matrix_type SMaux(nbcols, nb_dof);
250 gmm::mult(gmm::transposed(NS), SM, SMaux);
253 sparse_matrix_type NSaux(nb_dof, nbcols);
gmm::copy(NS, NSaux);
259 bool laplacian_problem::solve(
void) {
263 cout <<
"Compute preconditionner\n";
265 double time = gmm::uclock_sec();
275 cout <<
"Time to compute preconditionner : "
276 << gmm::uclock_sec() - time <<
" seconds\n";
285 #if defined(GMM_USES_SUPERLU)
286 gmm::SuperLU_solve(SM, U, B, rcond);
288 cout <<
"cond = " << 1/rcond <<
"\n";
291 cout <<
"Total time to solve : "
292 << gmm::uclock_sec() - time <<
" seconds\n";
295 std::vector<scalar_type> Uaux(mf_u.nb_dof());
301 return (iter.converged());
305 void laplacian_problem::compute_error() {
306 std::vector<scalar_type> V(mf_rhs.nb_basic_dof());
308 for (
size_type i = 0; i < mf_rhs.nb_basic_dof(); ++i)
309 V[i] -= sol_u(mf_rhs.point_of_basic_dof(i));
320 int main(
int argc,
char *argv[]) {
322 GETFEM_MPI_INIT(argc, argv);
323 GMM_SET_EXCEPTION_DEBUG;
328 p.PARAM.read_command_line(argc, argv);
330 p.mesh.write_to_file(p.datafilename +
".mesh");
332 if (!p.solve()) GMM_ASSERT1(
false,
"Solve procedure has failed");
335 GMM_STANDARD_CATCH_ERROR;
Describe a finite element method linked to a mesh.
Describe an integration method linked to a mesh.
"iterator" class for regions.
structure used to hold a set of convexes and/or convex faces.
Describe a mesh (collection of convexes (elements) and points).
Incomplete LU with threshold and K fill-in Preconditioner.
The Iteration object calculates whether the solution has reached the desired accuracy,...
sparse vector built upon std::vector.
Miscelleanous assembly routines for common terms. Use the low-level generic assembly....
Compute the gradient of a field on a getfem::mesh_fem.
Export solutions to various formats.
Include common gmm files.
void copy(const L1 &l1, L2 &l2)
*/
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norminf(const V &v)
Infinity norm of a vector.
void clear(L &l)
clear (fill with zeros) a vector or matrix.
void resize(V &v, size_type n)
*/
void clean(L &l, double threshold)
Clean a vector or matrix (replace near-zero entries with zeroes).
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*/
strongest_value_type< V1, V2 >::value_type vect_sp(const V1 &v1, const V2 &v2)
*/
void gmres(const Mat &A, Vec &x, const VecB &b, const Precond &M, int restart, iteration &outer, Basis &KS)
Generalized Minimum Residual.
void asm_stiffness_matrix_for_laplacian(MAT &M, const mesh_im &mim, const mesh_fem &mf, const mesh_fem &mf_data, const VECT &A, const mesh_region &rg=mesh_region::all_convexes())
assembly of , where is scalar.
size_type Dirichlet_nullspace(const MAT1 &H, MAT2 &NS, const VECT1 &R, VECT2 &U0)
Build an orthogonal basis of the kernel of H in NS, gives the solution of minimal norm of H*U = R in ...
void asm_source_term(const VECT1 &B, const mesh_im &mim, const mesh_fem &mf, const mesh_fem &mf_data, const VECT2 &F, const mesh_region &rg=mesh_region::all_convexes())
source term (for both volumic sources and boundary (Neumann) sources).
scalar_type asm_L2_norm(const mesh_im &mim, const mesh_fem &mf, const VEC &U, const mesh_region &rg=mesh_region::all_convexes())
compute , U might be real or complex
scalar_type asm_H1_norm(const mesh_im &mim, const mesh_fem &mf, const VEC &U, const mesh_region &rg=mesh_region::all_convexes())
compute the H1 norm of U.
void asm_dirichlet_constraints(MAT &H, VECT1 &R, const mesh_im &mim, const mesh_fem &mf_u, const mesh_fem &mf_mult, const mesh_fem &mf_r, const VECT2 &r_data, const mesh_region ®ion, int version=ASMDIR_BUILDALL)
Assembly of Dirichlet constraints in a weak form.
void asm_normal_source_term(VECT1 &B, const mesh_im &mim, const mesh_fem &mf, const mesh_fem &mf_data, const VECT2 &F, const mesh_region &rg)
Normal source term (for boundary (Neumann) condition).
void APIDECL outer_faces_of_mesh(const mesh &m, const dal::bit_vector &cvlst, convex_face_ct &flist)
returns a list of "exterior" faces of a mesh (i.e.
std::shared_ptr< const getfem::virtual_fem > pfem
type of pointer on a fem description
pfem fem_descriptor(const std::string &name)
get a fem descriptor from its string name.
pfem classical_fem(bgeot::pgeometric_trans pgt, short_type k, bool complete=false)
Give a pointer on the structures describing the classical polynomial fem of degree k on a given conve...
pgeometric_trans geometric_trans_descriptor(std::string name)
Get the geometric transformation from its string name.
size_t size_type
used as the common size type in the library
std::shared_ptr< const bgeot::geometric_trans > pgeometric_trans
pointer type for a geometric transformation
void interpolation_function(mesh_fem &mf_target, const VECT &VV, F &f, mesh_region rg=mesh_region::all_convexes())
interpolation of a function f on mf_target.
void regular_unit_mesh(mesh &m, std::vector< size_type > nsubdiv, bgeot::pgeometric_trans pgt, bool noised=false)
Build a regular mesh of the unit square/cube/, etc.
void interpolation(const mesh_fem &mf_source, const mesh_fem &mf_target, const VECTU &U, VECTV &V, int extrapolation=0, double EPS=1E-10, mesh_region rg_source=mesh_region::all_convexes(), mesh_region rg_target=mesh_region::all_convexes())
interpolation/extrapolation of (mf_source, U) on mf_target.
pintegration_method int_method_descriptor(std::string name, bool throw_if_not_found=true)
Get an integration method from its name .