annotate src/assemble_system.cc @ 173:9e944b0d0fc8

Some Formatting improvements (?)
author gedeone-octave <marcovass89@hotmail.it>
date Sat, 12 Oct 2013 16:06:00 +0100
parents f65252c56853
children 66071811eef8
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"
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
20
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
21 DEFUN_DLD (assemble_system, args, nargout,
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
22 "-*- texinfo -*-\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
23 @deftypefn {Function File} {[@var{A}], [@var{b}], [@var{x}(Optional)]} = \
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
24 assemble_system (@var{form a}, @var{form L}, @var{DirichletBC}(Optional), \
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
25 @var{...}) \n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
26 The input arguments are\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
27 @itemize @bullet\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
28 @item @var{form a} the bilinear form to assemble.\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
29 @item @var{form a} the linear form to assemble.\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
30 @item @var{DirichletBC} represents the optional BC that you wish to apply to\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
31 the system. If more than one BC has to be applied, just list them.\n\
125
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
32 @end itemize \n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
33 The output @var{A} is a discretized representation of the system:\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
34 @itemize @bullet\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
35 @item @var{A} is the sparse Matrix corresponding to the @var{form a}\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
36 @item @var{A} is the Vector corresponding to the @var{form L}\n\
125
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
37 @end itemize \n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
38 If you need to apply boundary condition to a system for a nonlinear problem \n\
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
39 then you should provide as 3rd argument the vector and you will receive it \n\
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
40 back as the third output argument.\n\
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
41 For an example of this situation, you can look the example HyperElasticity.m\n\
143
9486cbdc0a2e Maint: update the documentation
gedeone-octave <marcovass89@hotmail.it>
parents: 125
diff changeset
42 @seealso{BilinearForm, LinearForm, ResidualForm, JacobianForm}\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
43 @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
44 {
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
45 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
46 octave_value_list retval;
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
47 if (! form_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
48 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
49 form::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
50 form_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
51 mlock ();
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
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
54 if (! boundarycondition_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
55 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
56 boundarycondition::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
57 boundarycondition_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
58 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
59 }
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 (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
62 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
63 if (nargin < 2)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
64 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
65 else
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
66 {
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
67
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
68 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
69 && 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
70 {
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
71 const form & frm1 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
72 static_cast<const form&> (args(0).get_rep ());
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
73 const form & frm2 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
74 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
75
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
76 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
77 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
78 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
79 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
80 a.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
81 b.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
82
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
83 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
84 {
169
f65252c56853 Add the linear algebra backend also for assemble_system.cc
gedeone-octave <marcovass89@hotmail.it>
parents: 168
diff changeset
85 dolfin::parameters["linear_algebra_backend"] = "uBLAS";
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
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
97 = static_cast<const boundarycondition&>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
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
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
100 const std::vector<boost::shared_ptr
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
101 <const dolfin::DirichletBC> >
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
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)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
105 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
106 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
107 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
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
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
111 // Get capacity of the dolfin sparse matrix
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
112 boost::tuples::tuple<const std::size_t*,
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
113 const std::size_t*,
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
114 const double*, int>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
115 aa = A.data ();
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
116
166
40ae9e5dfb93 Use the uBLAS for an estimation of the nnz elements for excess.
gedeone-octave <marcovass89@hotmail.it>
parents: 151
diff changeset
117 int nnz = aa.get<3> ();
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
118 std::size_t nr = A.size (0), nc = A.size (1);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
119 std::vector<double> data_tmp;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
120 std::vector<std::size_t> cidx_tmp;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
121
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
122 dim_vector dims (nnz, 1);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
123 octave_idx_type nz = 0, ii = 0;
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
124 Array<octave_idx_type>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
125 ridx (dims, 0),
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
126 cidx (dims, 0);
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
127 Array<double> data (dims, 0);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
128
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
129 octave_idx_type* orow = ridx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
130 octave_idx_type* oc = cidx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
131 double* ov = data.fortran_vec ();
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
132
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
133 for (std::size_t i = 0; i < nr; ++i)
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 A.getrow (i, cidx_tmp, data_tmp);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
136 nz += cidx_tmp.size ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
137
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
138 for (octave_idx_type j = 0;
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
139 j < cidx_tmp.size (); ++j)
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
140 {
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
141 orow [ii + j] = i;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
142 oc [ii + j] = cidx_tmp [j];
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
143 ov [ii + j] = data_tmp [j];
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
144 }
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
145
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
146 ii = nz;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
147 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
148
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
149 dims(0) = ii;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
150 ridx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
151 cidx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
152 data.resize (dims);
116
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 SparseMatrix sm (data, ridx, cidx, nr, nc);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
155 retval(0) = sm;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
156
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
157 dim_vector dim;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
158 dim.resize (2);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
159 dim(0) = B.size ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
160 dim(1) = 1;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
161 Array<double> myb (dim);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
162
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
163 for (std::size_t i = 0; i < B.size (); ++i)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
164 myb.xelem (i) = B[i];
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
165
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
166 retval(1) = myb;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
167 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
168 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
169 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
170 error ("assemble_system: unknown size");
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 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
173 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
174 else if (nargout == 3)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
175 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
176 if (nargin < 3)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
177 print_usage ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
178 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
179 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
180 if (args(0).type_id () == form::static_type_id ()
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
181 && args(1).type_id () == form::static_type_id ())
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
182 {
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
183 const form & frm1 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
184 static_cast<const form&> (args(0).get_rep ());
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
185 const form & frm2 =
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
186 static_cast<const form&> (args(1).get_rep ());
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
187 const Array<double> myx = args(2).array_value ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
188
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
189 if (! error_state)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
190 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
191 const dolfin::Form & a = frm1.get_form ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
192 const dolfin::Form & b = frm2.get_form ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
193 a.check ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
194 b.check ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
195
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
196 if (a.rank () == 2 && b.rank () == 1)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
197 {
169
f65252c56853 Add the linear algebra backend also for assemble_system.cc
gedeone-octave <marcovass89@hotmail.it>
parents: 168
diff changeset
198 dolfin::parameters["linear_algebra_backend"] = "uBLAS";
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
199 dolfin::Matrix A;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
200 dolfin::assemble (A, a);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
201 dolfin::Vector B;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
202 dolfin::assemble (B, b);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
203 dolfin::Vector x (myx.length ());
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
204
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
205 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
206 x.setitem (i, myx.xelem (i));
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
207
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
208 for (std::size_t i = 3; i < nargin; ++i)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
209 {
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
210 if (args(i).type_id () ==
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
211 boundarycondition::static_type_id ())
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
212 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
213 const boundarycondition & bc
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
214 = static_cast<const boundarycondition&>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
215 (args(i).get_rep ());
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
216
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
217 const std::vector<boost::shared_ptr
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
218 <const dolfin::DirichletBC> >
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
219 & pbc = bc.get_bc ();
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
220
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
221 for (std::size_t j = 0;
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
222 j < pbc.size (); ++j)
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
223 pbc[j]->apply(A, B, x);
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
224
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
225 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
226 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
227 error ("assemble_system: unknown argument type");
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
228 }
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
229
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
230 // Get capacity of the dolfin sparse matrix
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
231 boost::tuples::tuple<const std::size_t*,
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
232 const std::size_t*,
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
233 const double*, int>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
234 aa = A.data ();
166
40ae9e5dfb93 Use the uBLAS for an estimation of the nnz elements for excess.
gedeone-octave <marcovass89@hotmail.it>
parents: 151
diff changeset
235 int nnz = aa.get<3> ();
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
236 std::size_t nr = A.size (0), nc = A.size (1);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
237 std::vector<double> data_tmp;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
238 std::vector<std::size_t> cidx_tmp;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
239
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
240 dim_vector dims (nnz, 1);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
241 octave_idx_type nz = 0, ii = 0;
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
242 Array<octave_idx_type>
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
243 ridx (dims, 0),
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
244 cidx (dims, 0);
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
245 Array<double> data (dims, 0);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
246
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
247 octave_idx_type* orow = ridx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
248 octave_idx_type* oc = cidx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
249 double* ov = data.fortran_vec ();
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
250
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
251 for (std::size_t i = 0; i < nr; ++i)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
252 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
253 A.getrow (i, cidx_tmp, data_tmp);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
254 nz += cidx_tmp.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
255
173
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
256 for (octave_idx_type j = 0;
9e944b0d0fc8 Some Formatting improvements (?)
gedeone-octave <marcovass89@hotmail.it>
parents: 169
diff changeset
257 j < cidx_tmp.size (); ++j)
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
258 {
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
259 orow [ii + j] = i;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
260 oc [ii + j] = cidx_tmp [j];
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
261 ov [ii + j] = data_tmp [j];
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
262 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
263
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
264 ii = nz;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
265 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
266
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
267 dims(0) = ii;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
268 ridx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
269 cidx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
270 data.resize (dims);
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
271
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
272 SparseMatrix sm (data, ridx, cidx, nr, nc);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
273 retval(0) = sm;
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
274
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
275 dim_vector dim;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
276 dim.resize (2);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
277 dim(0) = B.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
278 dim(1) = 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
279 Array<double> myb (dim), myc (dim);
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
280
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
281 for (std::size_t i = 0; i < B.size (); ++i)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
282 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
283 myb.xelem (i) = B[i];
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
284 myc.xelem (i) = x[i];
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
285 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
286
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
287 retval(1) = myb;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
288 retval(2) = myc;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
289 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
290 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
291 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
292 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
293 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
294 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
295 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
296 return retval;
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
297 }