Mercurial > fem-fenics-eugenio
annotate src/assemble_system.cc @ 268:61830a4f9ab9
Improve formatting
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Thu, 14 Aug 2014 12:26:55 +0200 |
parents | 5e9b5bbdc56b |
children | f4d6ae912a08 |
rev | line source |
---|---|
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
1 /* |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
2 Copyright (C) 2013 Marco Vassallo <gedeone-octave@users.sourceforge.net> |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
3 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
4 This program is free software; you can redistribute it and/or modify it under |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
5 the terms of the GNU General Public License as published by the Free Software |
151 | 6 Foundation; either version 3 of the License, or (at your option) any later |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
7 version. |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
8 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
9 This program is distributed in the hope that it will be useful, but WITHOUT |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
10 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
11 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
12 details. |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
13 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
14 You should have received a copy of the GNU General Public License along with |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
15 this program; if not, see <http://www.gnu.org/licenses/>. |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
16 */ |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
17 |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
18 #include "form.h" |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
19 #include "boundarycondition.h" |
247
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
20 #include "femfenics_factory.h" |
253
5e9b5bbdc56b
Support both DOLFIN 1.3.0 and 1.4.0
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
247
diff
changeset
|
21 #include "dolfin_compat.h" |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
22 |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
23 DEFUN_DLD (assemble_system, args, nargout, |
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
24 "-*- texinfo -*-\n\ |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
25 @deftypefn {Function File} {[@var{A}], [@var{b}], [@var{x}(Optional)]} = \ |
184
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
26 assemble_system (@var{form_a}, @var{form_L}, @var{DirichletBC})\n\ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
27 Construct the discretization of a system and apply essential BC.\n\ |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
28 The input arguments are\n\ |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
29 @itemize @bullet\n\ |
184
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
30 @item @var{form_a} which is the BilinearForm to assemble.\n\ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
31 @item @var{form_L} which is the LinearForm to assemble.\n\ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
32 @item @var{DirichletBC} represents the optional BC applied to \ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
33 the system. \n\ |
125
a80ac536c78a
Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents:
117
diff
changeset
|
34 @end itemize \n\ |
184
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
35 The output @var{A} is a matrix representing the @var{form_a} while \ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
36 @var{b} represents @var{form_L}. \n\ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
37 If boundary conditions have to be applied to a vector for a nonlinear problem \ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
38 then it should be provide as 3rd argument and it will be given back \ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
39 as the 3rd output argument. For an example of this situation, please refer \ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
40 to the HyperElasticity example. \n\ |
66071811eef8
Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents:
173
diff
changeset
|
41 @seealso{BilinearForm, LinearForm, ResidualForm, JacobianForm, Functional}\n\ |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
42 @end deftypefn") |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
43 { |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
44 int nargin = args.length (); |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
45 octave_value_list retval; |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
46 if (! form_type_loaded) |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
47 { |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
48 form::register_type (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
49 form_type_loaded = true; |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
50 mlock (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
51 } |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
52 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
53 if (! boundarycondition_type_loaded) |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
54 { |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
55 boundarycondition::register_type (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
56 boundarycondition_type_loaded = true; |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
57 mlock (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
58 } |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
59 |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
60 if (nargout == 2) |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
61 { |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
62 if (nargin < 2) |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
63 { print_usage (); } |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
64 else |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
65 { |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
66 |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
67 if (args(0).type_id () == form::static_type_id () |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
68 && args(1).type_id () == form::static_type_id ()) |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
69 { |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
70 const form & frm1 = |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
71 static_cast<const form &> (args(0).get_rep ()); |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
72 const form & frm2 = |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
73 static_cast<const form &> (args(1).get_rep ()); |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
74 |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
75 if (! error_state) |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
76 { |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
77 const dolfin::Form & a = frm1.get_form (); |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
78 const dolfin::Form & b = frm2.get_form (); |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
79 a.check (); |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
80 b.check (); |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
81 |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
82 if (a.rank () == 2 && b.rank () == 1) |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
83 { |
247
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
84 femfenics_factory factory; |
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
85 |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
86 dolfin::Matrix A; |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
87 dolfin::assemble (A, a); |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
88 dolfin::Vector B; |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
89 dolfin::assemble (B, b); |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
90 |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
91 for (std::size_t i = 2; i < nargin; ++i) |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
92 { |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
93 if (args(i).type_id () == |
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
94 boundarycondition::static_type_id ()) |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
95 { |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
96 const boundarycondition & bc |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
97 = static_cast<const boundarycondition &> |
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
98 (args(i).get_rep ()); |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
99 |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
100 const std::vector<SHARED_PTR |
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
101 <const dolfin::DirichletBC> > |
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
102 & pbc = bc.get_bc (); |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
103 |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
104 for (std::size_t j = 0; j < pbc.size (); ++j) |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
105 { pbc[j]->apply (A, B); } |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
106 } |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
107 else |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
108 { error ("assemble_system: unknown argument type"); } |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
109 } |
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
110 |
247
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
111 retval(0) = factory.matrix (A); |
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
112 retval(1) = factory.vector (B); |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
113 } |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
114 } |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
115 else |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
116 { error ("assemble_system: unknown size"); } |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
117 } |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
118 } |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
119 } |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
120 else if (nargout == 3) |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
121 { |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
122 if (nargin < 3) |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
123 { print_usage (); } |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
124 else |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
125 { |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
126 if (args(0).type_id () == form::static_type_id () |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
127 && args(1).type_id () == form::static_type_id ()) |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
128 { |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
129 const form & frm1 = |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
130 static_cast<const form &> (args(0).get_rep ()); |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
131 const form & frm2 = |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
132 static_cast<const form &> (args(1).get_rep ()); |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
133 const Array<double> myx = args(2).array_value (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
134 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
135 if (! error_state) |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
136 { |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
137 const dolfin::Form & a = frm1.get_form (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
138 const dolfin::Form & b = frm2.get_form (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
139 a.check (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
140 b.check (); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
141 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
142 if (a.rank () == 2 && b.rank () == 1) |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
143 { |
247
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
144 femfenics_factory factory; |
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
145 |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
146 dolfin::Matrix A; |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
147 dolfin::assemble (A, a); |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
148 dolfin::Vector B; |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
149 dolfin::assemble (B, b); |
253
5e9b5bbdc56b
Support both DOLFIN 1.3.0 and 1.4.0
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
247
diff
changeset
|
150 #ifdef LATEST_DOLFIN |
5e9b5bbdc56b
Support both DOLFIN 1.3.0 and 1.4.0
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
247
diff
changeset
|
151 dolfin::Vector x (MPI_COMM_WORLD, myx.length ()); |
5e9b5bbdc56b
Support both DOLFIN 1.3.0 and 1.4.0
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
247
diff
changeset
|
152 #else |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
153 dolfin::Vector x (myx.length ()); |
253
5e9b5bbdc56b
Support both DOLFIN 1.3.0 and 1.4.0
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
247
diff
changeset
|
154 #endif |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
155 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
156 for (std::size_t i = 0; i < myx.length (); ++i) |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
157 { x.setitem (i, myx.xelem (i)); } |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
158 |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
159 for (std::size_t i = 3; i < nargin; ++i) |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
160 { |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
161 if (args(i).type_id () == |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
162 boundarycondition::static_type_id ()) |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
163 { |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
164 const boundarycondition & bc |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
165 = static_cast<const boundarycondition &> |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
166 (args(i).get_rep ()); |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
167 |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
168 const std::vector<SHARED_PTR |
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
169 <const dolfin::DirichletBC> > |
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
170 & pbc = bc.get_bc (); |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
171 |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
172 for (std::size_t j = 0; |
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
173 j < pbc.size (); ++j) |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
174 { pbc[j]->apply (A, B, x); } |
173
9e944b0d0fc8
Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents:
169
diff
changeset
|
175 |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
176 } |
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
177 else |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
178 { error ("assemble_system: unknown argument type"); } |
116
77eefe47f7ab
Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents:
97
diff
changeset
|
179 } |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
180 |
247
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
181 retval(0) = factory.matrix (A); |
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
182 retval(1) = factory.vector (B); |
8ca45824938e
Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
184
diff
changeset
|
183 retval(2) = factory.vector (x); |
97
dc9987325fea
Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents:
91
diff
changeset
|
184 } |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
185 } |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
186 else |
268
61830a4f9ab9
Improve formatting
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
253
diff
changeset
|
187 { error ("assemble_system: unknown size"); } |
91
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
188 } |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
189 } |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
190 } |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
191 return retval; |
51bfdf30b822
New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
192 } |