annotate src/assemble.cc @ 184:66071811eef8

Improve the documentation for the pkg releae.
author gedeone-octave <marcovass89@hotmail.it>
date Sat, 09 Nov 2013 10:44:51 +0000
parents 8ea37cfc7a14
children 8ca45824938e
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: 108
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: 108
diff changeset
3
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
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: 108
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: 108
diff changeset
7 version.
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
8
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
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: 108
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: 108
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: 108
diff changeset
12 details.
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
13
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
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: 108
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: 108
diff changeset
16 */
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
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
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
21 DEFUN_DLD (assemble, args, nargout,
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
22 "-*- texinfo -*-\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
23 @deftypefn {Function File} {[@var{A}], [@var{x}(Optional)]} = \
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
24 assemble (@var{form_a}, @var{DirichletBC}) \n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
25 Construct the discretization of a Form and apply essential BC.\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
26 The input arguments are\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
27 @itemize @bullet\n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
28 @item @var{form_a} which is the form to assemble.\n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
29 It can be a form of rank 2 (BilinearForm or JacobianForm), \
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
30 a form of rank 1 (LinearForm or ResidualForm) or a form \
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
31 of rank 0 (Functional).\n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
32 @item @var{DirichletBC} represents the optional BC applied to \
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
33 the system. \n\
125
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
34 @end itemize \n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
35 The output @var{A} is a discretized representation of the @var{form_a}:\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
36 @itemize @bullet\n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
37 @item @var{A} is a sparse Matrix if @var{form_a} is a bilinear form\n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
38 @item @var{A} is a Vector if @var{form_a} is a linear form\n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
39 @item @var{A} is a Double if @var{form_a} is a functional\n\
125
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
40 @end itemize \n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
41 If boundary condition has to be applied to a vector for a nonlinear problem \
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
42 then it should be provide as 2nd argument and it will be given back \
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
43 as the second output argument. For an example of this situation, please refer \
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
44 to the HyperElasticity example. \n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
45 @seealso{BilinearForm, LinearForm, ResidualForm, JacobianForm, Functional}\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
46 @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
47 {
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
48 int nargin = args.length ();
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
49 octave_value_list retval;
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
50
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
51 if (! boundarycondition_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
52 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
53 boundarycondition::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
54 boundarycondition_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
55 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
56 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
57
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
58 if (! form_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
59 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
60 form::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
61 form_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
62 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
63 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
64
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
65 if (nargout == 1)
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
66 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
67 if (nargin < 1)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
68 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
69 else
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
70 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
71 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
72 {
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
73 const form & frm = static_cast<const form&> (args(0).get_rep ());
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
74
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
75 if (! error_state)
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
76 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
77 const dolfin::Form & a = frm.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 ();
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
79
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
80 if (a.rank () == 2)
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
81 {
166
40ae9e5dfb93 Use the uBLAS for an estimation of the nnz elements for excess.
gedeone-octave <marcovass89@hotmail.it>
parents: 151
diff changeset
82 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
83 dolfin::Matrix A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
84 dolfin::assemble (A, a);
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 = 1; i < nargin; ++i)
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
87 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
88 if (args(i).type_id () ==
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
89 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
90 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
91 const boundarycondition & bc
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
92 = static_cast<const boundarycondition&>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
93 (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
94
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
95 const
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
96 std::vector<boost::shared_ptr
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
97 <const dolfin::DirichletBC> >
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
98 & pbc = bc.get_bc ();
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
99
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
100 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
101 pbc[j]->apply(A);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
102 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
103 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
104 error ("assemble: unknown argument type");
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
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
107 // Get capacity of the dolfin sparse matrix
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
108 boost::tuples::tuple<const std::size_t*,
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
109 const std::size_t*,
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
110 const double*, int>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
111 aa = A.data ();
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
112
166
40ae9e5dfb93 Use the uBLAS for an estimation of the nnz elements for excess.
gedeone-octave <marcovass89@hotmail.it>
parents: 151
diff changeset
113 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
114 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
115 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
116 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
117
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
118 dim_vector dims (nnz, 1);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
119 octave_idx_type nz = 0, ii = 0;
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
120 Array<octave_idx_type>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
121 ridx (dims, 0),
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
122 cidx (dims, 0);
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
123 Array<double> data (dims, 0);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
124
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
125 octave_idx_type* orow = ridx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
126 octave_idx_type* oc = cidx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
127 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
128
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
129 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
130 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
131 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
132 nz += cidx_tmp.size ();
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
133
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
134 for (octave_idx_type j = 0;
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
135 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
136 {
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
137 orow [ii + j] = i;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
138 oc [ii + j] = cidx_tmp [j];
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
139 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
140 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
141
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
142 ii = nz;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
143 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
144
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
145 dims(0) = ii;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
146 ridx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
147 cidx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
148 data.resize (dims);
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
149
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
150 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
151 retval(0) = sm;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
152 }
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
153
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
154 else if (a.rank () == 1)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
155 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
156 dolfin::Vector A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
157 dolfin::assemble (A, a);
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
158
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
159 for (std::size_t i = 1; i < nargin; ++i)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
160 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
161 if (args(i).type_id () ==
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
162 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
163 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
164 const boundarycondition & bc
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
165 = static_cast<const boundarycondition&>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
166 (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
167
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
168 const std::vector<boost::shared_ptr
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
169 <const dolfin::DirichletBC> >
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
170 & pbc = bc.get_bc ();
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
171
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
172 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
173 pbc[j]->apply(A);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
174 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
175 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
176 error ("assemble: unknown argument type");
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
177 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
178
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
179 dim_vector dims;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
180 dims.resize (2);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
181 dims(0) = A.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
182 dims(1) = 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
183 Array<double> myb (dims);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
184
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
185 for (std::size_t i = 0; i < A.size (); ++i)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
186 myb.xelem (i) = A[i];
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
187
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
188 retval(0) = myb;
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
189 }
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
190
108
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
191 else if (a.rank () == 0)
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
192 {
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
193 double b = dolfin::assemble (a);
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
194 retval(0) = octave_value (b);
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
195 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
196
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
197 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
198 error ("assemble: unknown form size");
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
199 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
200 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
201 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
202 }
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
203 else if (nargout == 2)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
204 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
205 if (nargin < 2)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
206 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
207 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
208 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
209 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
210 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
211 const form & frm =
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
212 static_cast<const form &> (args(0).get_rep ());
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
213 const Array<double> myx = args(1).array_value ();
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 if (! error_state)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
216 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
217 const dolfin::Form & a = frm.get_form ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
218 a.check ();
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 if (a.rank () == 1)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
221 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
222 dolfin::Vector A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
223 dolfin::assemble (A, a);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
224
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
225 dolfin::Vector x (myx.length ());
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
226 for (std::size_t i = 0; i < myx.length (); ++i)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
227 x.setitem (i, myx.xelem (i));
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 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
230 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
231 if (args(i).type_id () ==
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
232 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
233 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
234 const boundarycondition & bc
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
235 = static_cast<const boundarycondition&>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
236 (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
237
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
238 const std::vector<boost::shared_ptr
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
239 <const dolfin::DirichletBC> >
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
240 & pbc = bc.get_bc ();
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
241
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
242 for (std::size_t j = 0;
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
243 j < pbc.size (); ++j)
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
244 pbc[j]->apply(A, x);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
245 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
246 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
247 error ("assemble: unknown argument type");
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
248 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
249
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
250 dim_vector dims;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
251 dims.resize (2);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
252 dims(0) = A.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
253 dims(1) = 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
254 Array<double> myb (dims), myc (dims);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
255
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
256 for (std::size_t i = 0; i < A.size (); ++i)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
257 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
258 myb.xelem (i) = A[i];
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
259 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
260 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
261
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
262 retval(0) = myb;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
263 retval(1) = myc;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
264 }
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 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
267 error ("assemble: unknown size");
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
268 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
269 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
270 }
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
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
273 return retval;
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
274 }