annotate src/assemble_system.cc @ 143:9486cbdc0a2e

Maint: update the documentation
author gedeone-octave <marcovass89@hotmail.it>
date Mon, 09 Sep 2013 22:40:46 +0200
parents a80ac536c78a
children 5fe2a157f4eb
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 /*
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
2 Copyright (C) 2013 Marco Vassallo
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
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
6 Foundation; either version 2 of the License, or (at your option) any later
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
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
21 DEFUN_DLD (assemble_system, args, nargout, "-*- texinfo -*-\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
22 @deftypefn {Function File} {[@var{A}], [@var{b}], [@var{x}(Optional)]} = \
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
23 assemble_system (@var{form a}, @var{form L}, @var{DirichletBC}(Optional), @var{...}) \n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
24 The input arguments are\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
25 @itemize @bullet\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
26 @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
27 @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
28 @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
29 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
30 @end itemize \n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
31 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
32 @itemize @bullet\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
33 @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
34 @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
35 @end itemize \n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
36 If you need to apply boundary condition to a system for a nonlinear problem \n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
37 then you should provide as 3rd argument the vector and you will receive it back\n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
38 as the third output argument. For an example of this situation, you can look\n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 117
diff changeset
39 the example HyperElasticity.m\n\
143
9486cbdc0a2e Maint: update the documentation
gedeone-octave <marcovass89@hotmail.it>
parents: 125
diff changeset
40 @seealso{BilinearForm, LinearForm, ResidualForm, JacobianForm}\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 {
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
69 const form & frm1 = static_cast<const form&> (args(0).get_rep ());
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
70 const form & frm2 = 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
71
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
72 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
73 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
74 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
75 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
76 a.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
77 b.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
78
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
79 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
80 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
81 dolfin::Matrix A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
82 dolfin::assemble (A, a);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
83 dolfin::Vector B;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
84 dolfin::assemble (B, b);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
85
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
86 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
87 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
88 if (args(i).type_id () == boundarycondition::static_type_id ())
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 const boundarycondition & bc
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
91 = static_cast<const boundarycondition&> (args(i).get_rep ());
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
92
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
93 const std::vector<boost::shared_ptr <const dolfin::DirichletBC> > & pbc
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
94 = 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
95
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
96 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
97 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
98 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
99 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
100 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
101 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
102
117
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
103 std::cout << "Assembling Matrix from the bilinear form..."
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
104 << std::endl;
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
105 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
106 std::vector<double> data_tmp;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
107 std::vector<std::size_t> cidx_tmp;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
108
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
109 octave_idx_type dims = A.size (0), nz = 0, ii = 0;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
110 ColumnVector ridx (dims), cidx (dims), data (dims);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
111
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
112 for (std::size_t i = 0; i < nr; ++i)
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 A.getrow (i, cidx_tmp, data_tmp);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
115 nz += cidx_tmp.size ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
116 if (dims < nz)
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 dims = 1.2 * ((nr * nz) / (i + 1));;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
119 ridx.resize (dims);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
120 cidx.resize (dims);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
121 data.resize (dims);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
122 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
123
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
124 for (octave_idx_type j = 0; j < cidx_tmp.size (); ++j)
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 ridx.xelem (ii + j) = i + 1;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
127 cidx.xelem (ii + j) = cidx_tmp [j] + 1;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
128 data.xelem (ii + j) = data_tmp [j];
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
129 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
130 ii = nz;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
131 }
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 ridx.resize (ii);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
134 cidx.resize (ii);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
135 data.resize (ii);
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 SparseMatrix sm (data, ridx, cidx, nr, nc);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
138 retval(0) = sm;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
139
117
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
140 std::cout << "Assembling Vector from the linear form..."
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
141 << std::endl;
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
142 dim_vector dim;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
143 dim.resize (2);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
144 dim(0) = B.size ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
145 dim(1) = 1;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
146 Array<double> myb (dim);
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 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
149 myb.xelem (i) = B[i];
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 retval(1) = myb;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
152 }
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 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
155 error ("assemble_system: unknown size");
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 }
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 else if (nargout == 3)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
160 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
161 std::cout << "Assemble_system: apply boundary condition to a vector for a nonlinear problem..."
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
162 << std::endl;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
163 if (nargin < 3)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
164 print_usage ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
165 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
166 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
167 if (args(0).type_id () == form::static_type_id ()
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
168 && args(1).type_id () == form::static_type_id ())
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
169 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
170 const form & frm1 = static_cast<const form&> (args(0).get_rep ());
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
171 const form & frm2 = static_cast<const form&> (args(1).get_rep ());
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
172 const Array<double> myx = args(2).array_value ();
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 if (! error_state)
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 const dolfin::Form & a = frm1.get_form ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
177 const dolfin::Form & b = frm2.get_form ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
178 a.check ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
179 b.check ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
180
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
181 if (a.rank () == 2 && b.rank () == 1)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
182 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
183 dolfin::Matrix A;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
184 dolfin::assemble (A, a);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
185 dolfin::Vector B;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
186 dolfin::assemble (B, b);
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
187 dolfin::Vector x (myx.length ());
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 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
190 x.setitem (i, myx.xelem (i));
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
191
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
192 for (std::size_t i = 3; i < nargin; ++i)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
193 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
194 if (args(i).type_id () == boundarycondition::static_type_id ())
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 const boundarycondition & bc
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
197 = static_cast<const boundarycondition&> (args(i).get_rep ());
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
198
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
199 const std::vector<boost::shared_ptr <const dolfin::DirichletBC> > & pbc
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
200 = bc.get_bc ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
201
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
202 for (std::size_t j = 0; j < pbc.size (); ++j)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
203 pbc[j]->apply(A, B, x);
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 else
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
206 error ("assemble_system: unknown argument type");
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
207 }
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
208
117
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
209 std::cout << "Assembling Matrix from the bilinear form..."
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
210 << std::endl;
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
211 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
212 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
213 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
214
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
215 octave_idx_type dims = A.size (0), nz = 0, ii = 0;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
216 ColumnVector ridx (dims), cidx (dims), data (dims);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
217
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
218 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
219 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
220 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
221 nz += cidx_tmp.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
222 if (dims < nz)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
223 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
224 dims = 1.2 * ((nr * nz) / (i + 1));;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
225 ridx.resize (dims);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
226 cidx.resize (dims);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
227 data.resize (dims);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
228 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
229
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
230 for (octave_idx_type j = 0; j < cidx_tmp.size (); ++j)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
231 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
232 ridx.xelem (ii + j) = i + 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
233 cidx.xelem (ii + j) = cidx_tmp [j] + 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
234 data.xelem (ii + j) = data_tmp [j];
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
235 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
236
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
237 ii = nz;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
238 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
239
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
240 ridx.resize (ii);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
241 cidx.resize (ii);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
242 data.resize (ii);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
243
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
244 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
245 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
246
117
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
247 std::cout << "Assembling Vector from the linear form..."
70322135b25f Maint: Add copyright notice
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
248 << std::endl;
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
249 dim_vector dim;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
250 dim.resize (2);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
251 dim(0) = B.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
252 dim(1) = 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
253 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
254
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
255 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
256 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
257 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
258 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
259 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
260
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
261 retval(1) = myb;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
262 retval(2) = myc;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
263 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
264 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
265 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
266 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
267 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
268 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
269 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
270 return retval;
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
271 }