annotate src/assemble.cc @ 168:67944f307560

Remove the verbose debug output
author gedeone-octave <marcovass89@hotmail.it>
date Sun, 06 Oct 2013 22:31:48 +0100
parents 3895c8d4c066
children 8ea37cfc7a14
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 */
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
17 #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
18 #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
19
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
20 DEFUN_DLD (assemble, args, nargout, "-*- texinfo -*-\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
21 @deftypefn {Function File} {[@var{A}], [@var{x}(Optional)]} = \
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
22 assemble (@var{form a}, @var{DirichletBC}(Optional), @var{...}) \n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
23 The input arguments are\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
24 @itemize @bullet\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
25 @item @var{form a} which is the form to assemble. It can be a form of rank 2\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
26 (bilinear), a form of rank 1 (linear) or a form of rank 0 (functional).\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
27 @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: 108
diff changeset
28 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: 116
diff changeset
29 @end itemize \n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
30 The output @var{A} is a discretized representation of the @var{form a}:\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
31 @itemize @bullet\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
32 @item @var{A} is a sparse Matrix if @var{form a} is a bilinear form\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
33 @item @var{A} is a Vector if @var{form a} is a linear form\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
34 @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
35 @end itemize \n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
36 If you need to apply boundary condition to a vector for a nonlinear problem \n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
37 then you should provide as 2nd argument the vector and you will receive it back\n\
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
38 as the second 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: 116
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: 108
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 ();
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
44 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
45
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
46 if (! boundarycondition_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
47 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
48 boundarycondition::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
49 boundarycondition_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
50 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
51 }
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 if (! form_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
54 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
55 form::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
56 form_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
57 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
58 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
59
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
60 if (nargout == 1)
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
61 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
62 if (nargin < 1)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
63 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
64 else
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
65 {
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 {
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
68 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
69
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
70 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
71 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
72 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
73 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
74
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
75 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
76 {
166
40ae9e5dfb93 Use the uBLAS for an estimation of the nnz elements for excess.
gedeone-octave <marcovass89@hotmail.it>
parents: 151
diff changeset
77 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
78 dolfin::Matrix A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
79 dolfin::assemble (A, a);
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 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
82 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
83 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
84 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
85 const boundarycondition & bc
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
86 = 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
87
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
88 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
89 = bc.get_bc ();
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 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
92 pbc[j]->apply(A);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
93 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
94 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
95 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
96 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
97
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
98 // It provides an upper boung for the nnz elements
166
40ae9e5dfb93 Use the uBLAS for an estimation of the nnz elements for excess.
gedeone-octave <marcovass89@hotmail.it>
parents: 151
diff changeset
99 boost::tuples::tuple<const std::size_t*, const std::size_t*, const double*, int> aa = A.data ();
40ae9e5dfb93 Use the uBLAS for an estimation of the nnz elements for excess.
gedeone-octave <marcovass89@hotmail.it>
parents: 151
diff changeset
100 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
101 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
102 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
103 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
104
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
105 dim_vector dims (nnz, 1);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
106 octave_idx_type nz = 0, ii = 0;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
107 Array<octave_idx_type> ridx (dims, 0), cidx (dims, 0);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
108 Array<double> data (dims, 0);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
109
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
110 octave_idx_type* orow = ridx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
111 octave_idx_type* oc = cidx.fortran_vec ();
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
112 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
113
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
114 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
115 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
116 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
117 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
118
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
119 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
120 {
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
121 orow [ii + j] = i;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
122 oc [ii + j] = cidx_tmp [j];
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
123 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
124 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
125
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
126 ii = nz;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
127 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
128
167
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
129 dims(0) = ii;
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
130 ridx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
131 cidx.resize (dims);
3895c8d4c066 Improved assembling of the matrix.
gedeone-octave <marcovass89@hotmail.it>
parents: 166
diff changeset
132 data.resize (dims);
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
133
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
134 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
135 retval(0) = sm;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
136 }
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
137
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
138 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
139 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
140 dolfin::Vector A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
141 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
142
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
143 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
144 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
145 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
146 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
147 const boundarycondition & bc
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
148 = 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
149
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
150 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
151 = bc.get_bc ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
152
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
153 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
154 pbc[j]->apply(A);
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 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
157 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
158 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
159
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
160 dim_vector dims;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
161 dims.resize (2);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
162 dims(0) = A.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
163 dims(1) = 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
164 Array<double> myb (dims);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
165
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
166 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
167 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
168
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
169 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
170 }
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
171
108
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
172 else if (a.rank () == 0)
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
173 {
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
174 double b = dolfin::assemble (a);
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
175 retval(0) = octave_value (b);
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
176 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
177
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
178 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
179 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
180 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
181 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
182 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
183 }
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
184 else if (nargout == 2)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
185 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
186 if (nargin < 2)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
187 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
188 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
189 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
190 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
191 {
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
192 const form & frm = 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
193 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
194
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
195 if (! error_state)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
196 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
197 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
198 a.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
199
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
200 if (a.rank () == 1)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
201 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
202 dolfin::Vector A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
203 dolfin::assemble (A, a);
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 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
206 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
207 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
208
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
209 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
210 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
211 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
212 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
213 const boundarycondition & bc
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
214 = 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
215
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
216 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
217 = bc.get_bc ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
218
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
219 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
220 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
221 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
222 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
223 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
224 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
225
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
226 dim_vector dims;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
227 dims.resize (2);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
228 dims(0) = A.size ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
229 dims(1) = 1;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
230 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
231
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
232 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
233 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
234 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
235 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
236 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
237
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
238 retval(0) = myb;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
239 retval(1) = myc;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
240 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
241
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
242 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
243 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
244 }
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 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
247 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
248
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
249 return retval;
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
250 }