annotate src/assemble.cc @ 247:8ca45824938e

Add factories hierarchy for matrices' and vectors' assembly
author Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
date Sat, 12 Jul 2014 12:30:13 +0200
parents 66071811eef8
children 5e9b5bbdc56b
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 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"
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
20 #include "femfenics_factory.h"
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
21
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
22 DEFUN_DLD (assemble, args, nargout,
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
23 "-*- texinfo -*-\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
24 @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
25 assemble (@var{form_a}, @var{DirichletBC}) \n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
26 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
27 The input arguments are\n\
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
28 @itemize @bullet\n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
29 @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
30 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
31 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
32 of rank 0 (Functional).\n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
33 @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
34 the system. \n\
125
a80ac536c78a Fix bug in the texinfo code
gedeone-octave <marcovass89@hotmail.it>
parents: 116
diff changeset
35 @end itemize \n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
36 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
37 @itemize @bullet\n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
38 @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
39 @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
40 @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
41 @end itemize \n\
184
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
42 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
43 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
44 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
45 to the HyperElasticity example. \n\
66071811eef8 Improve the documentation for the pkg releae.
gedeone-octave <marcovass89@hotmail.it>
parents: 172
diff changeset
46 @seealso{BilinearForm, LinearForm, ResidualForm, JacobianForm, Functional}\n\
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
47 @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
48 {
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
49 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
50 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
51
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
52 if (! boundarycondition_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
53 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
54 boundarycondition::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
55 boundarycondition_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
56 mlock ();
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
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
59 if (! form_type_loaded)
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
60 {
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
61 form::register_type ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
62 form_type_loaded = true;
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
63 mlock ();
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
64 }
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
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 (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
67 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
68 if (nargin < 1)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
69 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
70 else
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 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
73 {
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
74 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
75
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
76 if (! error_state)
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
77 {
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
78 const dolfin::Form & a = frm.get_form ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
79 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
80
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
81 if (a.rank () == 2)
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
82 {
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
83 femfenics_factory factory;
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
84
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
85 dolfin::Matrix A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
86 dolfin::assemble (A, a);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
87
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
88 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
89 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
90 if (args(i).type_id () ==
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
91 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
92 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
93 const boundarycondition & bc
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
94 = static_cast<const boundarycondition&>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
95 (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
96
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
97 const
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
98 std::vector<boost::shared_ptr
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
99 <const dolfin::DirichletBC> >
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
100 & pbc = bc.get_bc ();
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
101
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
102 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
103 pbc[j]->apply(A);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
104 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
105 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
106 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
107 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
108
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
109 retval(0) = factory.matrix (A);
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
110 }
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
111
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
112 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
113 {
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
114 femfenics_factory factory;
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
115
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
116 dolfin::Vector A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
117 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
118
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
119 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
120 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
121 if (args(i).type_id () ==
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
122 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
123 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
124 const boundarycondition & bc
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
125 = static_cast<const boundarycondition&>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
126 (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
127
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
128 const std::vector<boost::shared_ptr
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
129 <const dolfin::DirichletBC> >
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
130 & 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
131
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
132 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
133 pbc[j]->apply(A);
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
134 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
135 else
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
136 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
137 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
138
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
139 retval(0) = factory.vector (A);
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
140 }
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
141
108
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
142 else if (a.rank () == 0)
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
143 {
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
144 double b = dolfin::assemble (a);
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
145 retval(0) = octave_value (b);
5cbc7341ded5 Added the possibility to assemble form of rank 0.
gedeone-octave <marcovass89@hotmail.it>
parents: 97
diff changeset
146 }
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
147
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
148 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
149 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
150 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
151 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
152 }
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
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 (nargout == 2)
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 if (nargin < 2)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
157 print_usage ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
158 else
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 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
161 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
162 const form & frm =
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
163 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
164 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
165
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
166 if (! error_state)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
167 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
168 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
169 a.check ();
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
170
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
171 if (a.rank () == 1)
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
172 {
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
173 femfenics_factory factory;
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
174
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
175 dolfin::Vector A;
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
176 dolfin::assemble (A, a);
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 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
179 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
180 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
181
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
182 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
183 {
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
184 if (args(i).type_id () ==
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
185 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
186 {
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
187 const boundarycondition & bc
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
188 = static_cast<const boundarycondition&>
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
189 (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
190
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
191 const std::vector<boost::shared_ptr
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
192 <const dolfin::DirichletBC> >
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
193 & 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
194
172
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
195 for (std::size_t j = 0;
8ea37cfc7a14 some formatting improvements
Carlo de Falco <cdf@users.sourceforge.net>
parents: 168
diff changeset
196 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
197 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
198 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
199 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
200 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
201 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
202
247
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
203 retval(0) = factory.vector (A);
8ca45824938e Add factories hierarchy for matrices' and vectors' assembly
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents: 184
diff changeset
204 retval(1) = factory.vector (x);
97
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
205 }
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
206
dc9987325fea Added the possibility to apply BC to a NL problem.
gedeone-octave <marcovass89@hotmail.it>
parents: 91
diff changeset
207 else
116
77eefe47f7ab Maint: improve the txt message
gedeone-octave <marcovass89@hotmail.it>
parents: 108
diff changeset
208 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
209 }
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 }
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
91
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
214 return retval;
51bfdf30b822 New DLD functions for assembling a matrix/vector from a bilinear/linear form
gedeone-octave <marcovass89@hotmail.it>
parents:
diff changeset
215 }