Mercurial > fem-fenics-eugenio
comparison src/assemble_system.cc @ 253:5e9b5bbdc56b
Support both DOLFIN 1.3.0 and 1.4.0
* src/dolfin_compat.h: use a macro to set the correct shared_ptr (std or boost)
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Tue, 29 Jul 2014 18:05:56 +0200 |
parents | 8ca45824938e |
children | 61830a4f9ab9 |
comparison
equal
deleted
inserted
replaced
252:7f33554e439a | 253:5e9b5bbdc56b |
---|---|
16 */ | 16 */ |
17 | 17 |
18 #include "form.h" | 18 #include "form.h" |
19 #include "boundarycondition.h" | 19 #include "boundarycondition.h" |
20 #include "femfenics_factory.h" | 20 #include "femfenics_factory.h" |
21 #include "dolfin_compat.h" | |
21 | 22 |
22 DEFUN_DLD (assemble_system, args, nargout, | 23 DEFUN_DLD (assemble_system, args, nargout, |
23 "-*- texinfo -*-\n\ | 24 "-*- texinfo -*-\n\ |
24 @deftypefn {Function File} {[@var{A}], [@var{b}], [@var{x}(Optional)]} = \ | 25 @deftypefn {Function File} {[@var{A}], [@var{b}], [@var{x}(Optional)]} = \ |
25 assemble_system (@var{form_a}, @var{form_L}, @var{DirichletBC})\n\ | 26 assemble_system (@var{form_a}, @var{form_L}, @var{DirichletBC})\n\ |
94 { | 95 { |
95 const boundarycondition & bc | 96 const boundarycondition & bc |
96 = static_cast<const boundarycondition&> | 97 = static_cast<const boundarycondition&> |
97 (args(i).get_rep ()); | 98 (args(i).get_rep ()); |
98 | 99 |
99 const std::vector<boost::shared_ptr | 100 const std::vector<SHARED_PTR |
100 <const dolfin::DirichletBC> > | 101 <const dolfin::DirichletBC> > |
101 & pbc = bc.get_bc (); | 102 & pbc = bc.get_bc (); |
102 | 103 |
103 for (std::size_t j = 0; j < pbc.size (); ++j) | 104 for (std::size_t j = 0; j < pbc.size (); ++j) |
104 pbc[j]->apply(A, B); | 105 pbc[j]->apply(A, B); |
144 | 145 |
145 dolfin::Matrix A; | 146 dolfin::Matrix A; |
146 dolfin::assemble (A, a); | 147 dolfin::assemble (A, a); |
147 dolfin::Vector B; | 148 dolfin::Vector B; |
148 dolfin::assemble (B, b); | 149 dolfin::assemble (B, b); |
150 #ifdef LATEST_DOLFIN | |
151 dolfin::Vector x (MPI_COMM_WORLD, myx.length ()); | |
152 #else | |
149 dolfin::Vector x (myx.length ()); | 153 dolfin::Vector x (myx.length ()); |
154 #endif | |
150 | 155 |
151 for (std::size_t i = 0; i < myx.length (); ++i) | 156 for (std::size_t i = 0; i < myx.length (); ++i) |
152 x.setitem (i, myx.xelem (i)); | 157 x.setitem (i, myx.xelem (i)); |
153 | 158 |
154 for (std::size_t i = 3; i < nargin; ++i) | 159 for (std::size_t i = 3; i < nargin; ++i) |
158 { | 163 { |
159 const boundarycondition & bc | 164 const boundarycondition & bc |
160 = static_cast<const boundarycondition&> | 165 = static_cast<const boundarycondition&> |
161 (args(i).get_rep ()); | 166 (args(i).get_rep ()); |
162 | 167 |
163 const std::vector<boost::shared_ptr | 168 const std::vector<SHARED_PTR |
164 <const dolfin::DirichletBC> > | 169 <const dolfin::DirichletBC> > |
165 & pbc = bc.get_bc (); | 170 & pbc = bc.get_bc (); |
166 | 171 |
167 for (std::size_t j = 0; | 172 for (std::size_t j = 0; |
168 j < pbc.size (); ++j) | 173 j < pbc.size (); ++j) |