annotate src/assemble_system.cc @ 247:8ca45824938e

Add factories hierarchy for matrices' and vectors' assembly
author Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
date Sat, 12 Jul 2014 12:30:13 +0200
parents 66071811eef8
children 5e9b5bbdc56b
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
1 /*
151
5fe2a157f4eb Update GPL to v3
gedeone-octave <marcovass89@hotmail.it>
parents: 143
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
5fe2a157f4eb Update GPL to v3
gedeone-octave <marcovass89@hotmail.it>
parents: 143
diff changeset
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"
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
21
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
22 DEFUN_DLD (assemble_system, args, nargout,
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
23 "-*- texinfo -*-\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
24 @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
25 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
26 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
27 The input arguments are\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
28 @itemize @bullet\n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 173
diff changeset
29 @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
30 @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
31 @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
32 the system. \n\
125
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
33 @end itemize \n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 173
diff changeset
34 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
35 @var{b} represents @var{form_L}. \n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 173
diff changeset
36 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
37 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
38 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
39 to the HyperElasticity example. \n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 173
diff changeset
40 @seealso{BilinearForm, LinearForm, ResidualForm, JacobianForm, Functional}\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
41 @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
42 {
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
43 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
44 octave_value_list retval;
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
45 if (! form_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
46 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
47 form::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
48 form_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
49 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
50 }
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 if (! boundarycondition_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
53 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
54 boundarycondition::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
55 boundarycondition_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
56 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
57 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
58
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
59 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
60 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
61 if (nargin < 2)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
62 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
63 else
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
64 {
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
65
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
66 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
67 && 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
68 {
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
69 const form & frm1 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
70 static_cast<const form&> (args(0).get_rep ());
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
71 const form & frm2 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
72 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
73
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
74 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
75 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
76 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
77 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
78 a.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
79 b.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
80
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
81 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
82 {
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
83 femfenics_factory factory;
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
84
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
85 dolfin::Matrix A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
86 dolfin::assemble (A, a);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
87 dolfin::Vector B;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
88 dolfin::assemble (B, b);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
89
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
90 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
91 {
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
92 if (args(i).type_id () ==
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
93 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
94 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
95 const boundarycondition & bc
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
96 = static_cast<const boundarycondition&>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
97 (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
98
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
99 const std::vector<boost::shared_ptr
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
100 <const dolfin::DirichletBC> >
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
101 & 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
102
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
103 for (std::size_t j = 0; j < pbc.size (); ++j)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
104 pbc[j]->apply(A, B);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
105 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
106 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
107 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
108 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
109
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
110 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
111 retval(1) = factory.vector (B);
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
112 }
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 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
115 error ("assemble_system: unknown size");
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
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 else if (nargout == 3)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
120 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
121 if (nargin < 3)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
122 print_usage ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
123 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
124 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
125 if (args(0).type_id () == form::static_type_id ()
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
126 && args(1).type_id () == form::static_type_id ())
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
127 {
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
128 const form & frm1 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
129 static_cast<const form&> (args(0).get_rep ());
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
130 const form & frm2 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
131 static_cast<const form&> (args(1).get_rep ());
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
132 const Array<double> myx = args(2).array_value ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
133
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
134 if (! error_state)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
135 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
136 const dolfin::Form & a = frm1.get_form ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
137 const dolfin::Form & b = frm2.get_form ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
138 a.check ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
139 b.check ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
140
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
141 if (a.rank () == 2 && b.rank () == 1)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
142 {
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
143 femfenics_factory factory;
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
144
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
145 dolfin::Matrix A;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
146 dolfin::assemble (A, a);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
147 dolfin::Vector B;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
148 dolfin::assemble (B, b);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
149 dolfin::Vector x (myx.length ());
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
150
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
151 for (std::size_t i = 0; i < myx.length (); ++i)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
152 x.setitem (i, myx.xelem (i));
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
153
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
154 for (std::size_t i = 3; i < nargin; ++i)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
155 {
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
156 if (args(i).type_id () ==
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
157 boundarycondition::static_type_id ())
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 const boundarycondition & bc
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
160 = static_cast<const boundarycondition&>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
161 (args(i).get_rep ());
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
162
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
163 const std::vector<boost::shared_ptr
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
164 <const dolfin::DirichletBC> >
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
165 & pbc = bc.get_bc ();
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
166
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
167 for (std::size_t j = 0;
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
168 j < pbc.size (); ++j)
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
169 pbc[j]->apply(A, B, x);
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
170
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
171 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
172 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
173 error ("assemble_system: unknown argument type");
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
174 }
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
175
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
176 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
177 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
178 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
179 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
180 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
181 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
182 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
183 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
184 }
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 return retval;
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
187 }