/usr/include/rheolef/field.h is in librheolef-dev 6.6-1build2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 | # ifndef _RHEOLEF_FIELD_H
# define _RHEOLEF_FIELD_H
///
/// This file is part of Rheolef.
///
/// Copyright (C) 2000-2009 Pierre Saramito <Pierre.Saramito@imag.fr>
///
/// Rheolef is free software; you can redistribute it and/or modify
/// it under the terms of the GNU General Public License as published by
/// the Free Software Foundation; either version 2 of the License, or
/// (at your option) any later version.
///
/// Rheolef is distributed in the hope that it will be useful,
/// but WITHOUT ANY WARRANTY; without even the implied warranty of
/// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
/// GNU General Public License for more details.
///
/// You should have received a copy of the GNU General Public License
/// along with Rheolef; if not, write to the Free Software
/// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
///
/// =========================================================================
#include "rheolef/skit.h"
#include "rheolef/space.h"
#include "rheolef/undeterminated.h"
#include "rheolef/misc_algo.h" // copy_n
#include "rheolef/quadrature.h"
namespace rheolef {
// forward declaration:
template <class Expr> struct field_expr;
template <class T, class M> class field_indirect;
template <class T, class M> class field_indirect_const;
template <class T, class M> class field_component;
template <class T, class M> class field_component_const;
template <class T, class M> class field_concat_value;
template <class T, class M> class band_basic;
//<field:
/*Class:field
NAME: @code{field} - piecewise polynomial finite element field
DESCRIPTION:
@noindent
Store degrees of freedom associated to a mesh and
a piecewise polynomial approximation, with respect
to the numbering defined by the underlying @ref{space class}.
@noindent
This class contains two vectors, namely unknown and blocked
degrees of freedoms, and the associated finite element space.
Blocked and unknown degrees of freedom can be set by using
domain name indexation:
@example
geo omega ("circle");
space Xh (omega, "P1");
Xh.block ("boundary");
field uh (Xh);
uh ["boundary"] = 0;
@end example
@findex interpolate
INTERPOLATION:
Interpolation of a function @code{u} in a field @code{uh} with respect to
the interpolation writes:
@example
Float u (const point& x) @{ return x[0]*x[1]; @}
...
field uh = interpolate (Xh, u);
@end example
@findex dual
LINEAR ALGEBRA:
Linear algebra, such as @code{uh+vh}, @code{uh-vh} and @code{lambda*uh + mu*vh},
where @code{lambda} and @code{mu} are of type @code{Float}, are supported.
The duality product between two fields @code{lh} and @code{vh}
writes simply @code{dual(lh,vh)}: for discrete fields, it corresponds to a simple Euclidian dot
product in @code{IR^n}.
The application of a bilinear form (@pxref{form class})
writes @code{m(uh,vh)} and is equivalent to @code{dual(m*uh,vh)}.
NON-LINEAR ALGEBRA:
Non-linear operations, such as @code{sqrt(uh)} or @code{1/uh} are also available.
Notice that non-linear operations do not returns in general picewise polynomials:
the value returned by @code{sqrt(uh)} may be filtered by @code{interpolate},
@example
field vh = interpolate (Xh, sqrt(uh));
@end example
the Lagrange interpolant, to becomes a piecewise polynomial:
All standard unary and binary math functions @code{abs, cos, sin}... are available
on fields. Also @code{sqr(uh)}, the square of a field, and @code{min(uh,vh)}, @code{max(uh,vh)}
are provided. Binary functions can be used also with a scalar, as in
@example
field vh = interpolate (Xh, max (abs(uh), 0));
field wh = interpolate (Xh, pow (abs(uh), 1./3));
@end example
For applying a user-provided function to a field, use
the @code{compose} function:
@example
field vh = interpolate(Xh, compose(f, uh));
field wh = interpolate(Xh, compose(f, uh, vh));
@end example
The composition supports also general unary and binary class-functions.
Also, the multiplication @code{uh*vh} and the division @code{uh/vh}
returns a result that is not in the same discrete finite element space:
its result may be filtered by the @code{interpolate} operator:
@example
field wh = interpolate(Xh, uh*vh);
@end example
Any function or class function can be used in nonlinear expressions:
the function is interpolated in the specified finite element space.
@findex normal
There is a special predefined class-function named @code{normal} that represents the
outer unnit normal vector on a boundary domain or surfacic mesh:
@example
size_t k = omega.order();
string n_approx = "P" + itos(k-1) + "d";
space Nh (omega["boundary"], n_approx, "vector");
field nh = interpolate(Nh, normal());
@end example
The normal() function could appear in any nonlinear field expression:
it is evaluated on the fly, based on the current mesh.
Notice that when using isoparametric elements, the normal vector is no more
constant along any face of the mesh.
Also, on general curved domains, the unit normal vector is discontinuous
accross boundary element interfaces.
ACCESS BY DOMAIN:
The restriction of a field to a geometric domain, says @code{"boundary"}
writes @code{uh["boundary"]}: it represents the trace of the field on
the boundary:
@example
space Xh (omega, "P1");
uh["boundary"] = 0;
@end example
Extraction of the trace as a field is also possible:
@example
field wh = uh["boundary"];
@end example
The space associated to the trace writes @code{wh.get_space()} and
is equivalent to @code{space(omega["boundary"], "P1")}.
See @pxref{space class}.
VECTOR VALUED FIELD:
A vector-valued field contains several components, as:
@example
space Xh (omega, "P2", "vector");
field uh (Xh);
field vh = uh[0] - uh[1];
field nh = norm (uh);
@end example
The @code{norm} function returns the euclidian norm of the vector-valuated
field at each degree of freedom: its result is a scalar field.
TENSOR VALUED FIELD:
A tensor-valued field can be constructed and used as:
@example
space Th (omega, "P1d", "tensor");
field sigma_h (Xh);
field trace_h = sigma(0,0) + sigma_h(1,1);
field nh = norm (sigma_h);
@end example
The @code{norm} function returns the euclidian norm of the tensor-valuated
field at each degree of freedom: its result is a scalar field.
Notice that, as tensor-valued fields are symmetric, extra-diagonals are counted twice.
GENERAL MULTI-COMPONENT INTERFACE:
A general multi-component field writes:
@example
space Th (omega, "P1d", "tensor");
space Vh (omega, "P2", "vector");
space Qh (omega, "P1");
space Xh = Th*Vh*Qh;
field xh (Xh);
field tau_h = xh[0]; // tensor-valued
field uh = xh[1]; // vector-valued
field qh = xh[2]; // scalar
@end example
Remark the hierarchical multi-component field structure:
the first-component is tensor-valued and the second-one is vector-valued.
There is no limitation upon the hierarchical number of levels in use.
For any field @code{xh}, the string @code{xh.valued()} returns @code{"scalar"}
for a scalar field and @code{"vector"} for a vector-valued one. Other possible
valued are @code{"tensor"} and @code{"other"}.
The @code{xh.size()} returns the number of field components.
When the field is scalar, it returns zero by convention, and @code{xh[0]} is undefined.
A vector-valued field has @code{d} components, where @code{d=omega.dimension()}.
A tensor-valued field has @code{d*(d+1)/2} components, where @code{d=omega.dimension()}.
BLOCKED AND UNBLOCKED ARRAYS:
The field class contains two arrays of degrees-of-freedom (dof) associated
respectively to blocked and unknown dofs. Blocked dofs corresponds to
Dirichlet boundary conditions, as specified by space (See @pxref{space class}).
For simpliity, direct public access to these array is allowed, as @code{uh.b} and @code{uh.u}:
see @pxref{vec class}.
LOW-LEVEL DEGREE-OF-FREEDOM ACCESS:
The field class provides a STL-like container interface for accessing the degrees-of-freedom
(dof) of a finite element field @code{uh}. The number of dofs is @code{uh.ndof()}
and any dof can be accessed via @code{uh.dof(idof)}.
A non-local dof at the partition interface can be obtain via @code{uh.dis_dof(dis_idof)}
where @code{dis_idof} is the (global) distribued index assoiated to the
distribution @code{uh.ownership()}.
For performances, a STL-like iterator interface is available, with
@code{uh.begin_dof()} and @code{uh.end_dof()} returns iterators to the
arrays of dofs on the current processor.
See @pxref{array class} for more about distributed arrays.
For convenience, @code{uh.max()}, @code{uh.min()} and @code{uh.max_abs()}
retuns respectively the maximum, minimum and maximum of the absolute value
of the degrees of freedom.
FILE FORMAT:
TODO
IMPLEMENTATION NOTE:
The field expression use the expression template technics in
order to avoid temporaries when evaluating complex expressions.
AUTHOR: Pierre.Saramito@imag.fr
DATE: 23 february 2011
End:
*/
//<field:
template <class T, class M = rheo_default_memory_model>
class field_basic : public std::unary_function<point_basic<typename scalar_traits<T>::type>,T> {
public :
// typedefs:
typedef typename std::size_t size_type;
typedef M memory_type;
typedef T scalar_type;
typedef typename float_traits<T>::type float_type;
// typedef undeterminated_basic<T> value_type; // TODO
typedef T value_type; // TO_CLEAN
typedef space_constant::valued_type valued_type;
typedef geo_basic <float_type,M> geo_type;
typedef space_basic<float_type,M> space_type;
class iterator;
class const_iterator;
// allocator/deallocator:
field_basic();
explicit field_basic (
const space_type& V,
const T& init_value = std::numeric_limits<T>::max());
void resize (
const space_type& V,
const T& init_value = std::numeric_limits<T>::max());
field_basic (const field_indirect<T,M>&);
field_basic<T, M>& operator= (const field_indirect<T,M>&);
field_basic (const field_indirect_const<T,M>&);
field_basic<T, M>& operator= (const field_indirect_const<T,M>&);
field_basic (const field_component<T,M>&);
field_basic<T, M>& operator= (const field_component<T,M>&);
field_basic (const field_component_const<T,M>&);
field_basic<T, M>& operator= (const field_component_const<T,M>&);
template <class Expr> field_basic (const field_expr<Expr>&);
template <class Expr> field_basic<T, M>& operator= (const field_expr<Expr>&);
field_basic<T, M>& operator= (const T&);
// initializer list (c++ 2011):
#ifdef _RHEOLEF_HAVE_STD_INITIALIZER_LIST
field_basic (const std::initializer_list<field_concat_value<T,M> >& init_list);
field_basic<T,M>& operator= (const std::initializer_list<field_concat_value<T,M> >& init_list);
#endif // _RHEOLEF_HAVE_STD_INITIALIZER_LIST
// accessors:
const space_type& get_space() const { return _V; }
const geo_type& get_geo() const { return _V.get_geo(); }
std::string stamp() const { return _V.stamp(); }
std::string get_approx() const { return _V.get_approx(); }
valued_type valued_tag() const { return _V.valued_tag(); }
const std::string& valued() const { return _V.valued(); }
// accessors & modifiers to unknown & blocked parts:
const vec<T,M>& u() const { dis_dof_update_needed(); return _u; }
const vec<T,M>& b() const { dis_dof_update_needed(); return _b; }
vec<T,M>& set_u() { return _u; }
vec<T,M>& set_b() { return _b; }
// accessors to extremas:
T min() const;
T max() const;
T max_abs() const;
T min_abs() const;
// accessors by domains:
field_indirect<T,M> operator[] (const geo_basic<T,M>& dom);
field_indirect_const<T,M> operator[] (const geo_basic<T,M>& dom) const;
field_indirect<T,M> operator[] (std::string dom_name);
field_indirect_const<T,M> operator[] (std::string dom_name) const;
// accessors by components:
size_type size() const { return _V.size(); }
field_component<T,M> operator[] (size_type i_comp);
field_component_const<T,M> operator[] (size_type i_comp) const;
field_component<T,M> operator() (size_type i_comp, size_type j_comp);
field_component_const<T,M> operator() (size_type i_comp, size_type j_comp) const;
// accessors by degrees-of-freedom (dof):
const distributor& ownership() const { return get_space().ownership(); }
const communicator& comm() const { return ownership().comm(); }
size_type ndof() const { return ownership().size(); }
size_type dis_ndof() const { return ownership().dis_size(); }
T& dof (size_type idof);
const T& dof (size_type idof) const;
const T& dis_dof (size_type dis_idof) const;
iterator begin_dof();
iterator end_dof();
const_iterator begin_dof() const;
const_iterator end_dof() const;
// input/output:
idiststream& get (idiststream& ips);
odiststream& put (odiststream& ops) const;
odiststream& put_field (odiststream& ops) const;
// evaluate uh(x) where x is given locally as hat_x in K:
T dis_evaluate (const point_basic<T>& x, size_type i_comp = 0) const;
T operator() (const point_basic<T>& x) const { return dis_evaluate (x,0); }
point_basic<T> dis_vector_evaluate (const point_basic<T>& x) const;
// internals:
public:
// evaluate uh(x) where x is given locally as hat_x in K:
// requiers to call field::dis_dof_upgrade() before.
T evaluate (const geo_element& K, const point_basic<T>& hat_xq, size_type i_comp = 0) const;
// propagate changed values shared at partition boundaries to others procs
void dis_dof_update() const;
template <class Expr>
void assembly_internal (
const geo_basic<T,M>& dom,
const geo_basic<T,M>& band,
const band_basic<T,M>& gh,
const Expr& expr,
const quadrature_option_type& qopt,
bool is_on_band);
template <class Expr>
void assembly (
const geo_basic<T,M>& domain,
const Expr& expr,
const quadrature_option_type& qopt);
template <class Expr>
void assembly (
const band_basic<T,M>& gh,
const Expr& expr,
const quadrature_option_type& qopt);
protected:
void dis_dof_update_internal() const;
void dis_dof_update_needed() const;
// data:
space_type _V;
vec<T,M> _u;
vec<T,M> _b;
mutable bool _dis_dof_update_needed;
};
template <class T, class M>
idiststream& operator >> (odiststream& ips, field_basic<T,M>& u);
template <class T, class M>
odiststream& operator << (odiststream& ops, const field_basic<T,M>& uh);
typedef field_basic<Float> field;
typedef field_basic<Float,sequential> field_sequential;
//>field:
// =================================================================================
// field::iterator
// =================================================================================
template <class T, class M>
class field_basic<T,M>::iterator {
protected:
typedef typename vec<T,M>::iterator data_t; // random acess
typedef typename array<space_pair_type,M>::const_iterator iter_t; // forward access
public:
// typedefs:
typedef std::forward_iterator_tag iterator_category; // TODO: not fully random yet
typedef typename vec<T,M>::size_type size_type;
typedef T value_type;
typedef T& reference;
typedef T* pointer;
typedef std::ptrdiff_t difference_type;
// allocators:
iterator ()
: _blk_iub_iter(), _u(), _b()
{}
iterator (iter_t blk_iub_iter, data_t u, data_t b)
: _blk_iub_iter(blk_iub_iter), _u(u), _b(b)
{}
// accessors & modifiers:
reference operator* () const {
size_type iub = (*_blk_iub_iter).iub();
size_type blk = (*_blk_iub_iter).is_blocked();
if (!blk) return _u[iub]; else return _b[iub];
}
iterator& operator++ () { _blk_iub_iter++; return *this; }
iterator operator++ (int) { iterator tmp = *this; operator++(); return *this; }
iterator& operator+= (difference_type n) { _blk_iub_iter += n; return *this; }
iterator operator+ (difference_type n) const { iterator tmp = *this; return tmp += n; }
reference operator[] (size_type n) const { return *(*this + n); }
// comparators:
bool operator== (const iterator& j) const { return _blk_iub_iter == j._blk_iub_iter; }
bool operator!= (const iterator& j) const { return ! operator== (j); }
protected:
template <class T1, class M1> friend class field_basic<T1,M1>::const_iterator;
// data:
iter_t _blk_iub_iter;
data_t _u;
data_t _b;
};
template <class T, class M>
inline
typename field_basic<T,M>::iterator
field_basic<T,M>::begin_dof ()
{
dis_dof_update_needed();
return iterator (_V.data()._idof2blk_iub.begin(), _u.begin(), _b.begin());
}
template <class T, class M>
inline
typename field_basic<T,M>::iterator
field_basic<T,M>::end_dof ()
{
dis_dof_update_needed();
return iterator (_V.data()._idof2blk_iub.end(), _u.begin(), _b.begin());
}
// =================================================================================
// field::const_iterator
// =================================================================================
template <class T, class M>
class field_basic<T,M>::const_iterator {
protected:
typedef typename vec<T,M>::const_iterator data_t;
typedef typename array<space_pair_type,M>::const_iterator iter_t;
public:
// typedefs:
typedef std::forward_iterator_tag iterator_category;
typedef typename vec<T,M>::size_type size_type;
typedef T value_type;
typedef const T& reference;
typedef const T* pointer;
typedef std::ptrdiff_t difference_type;
// allocators:
const_iterator ()
: _blk_iub_iter(), _u(), _b()
{}
const_iterator (iter_t blk_iub_iter, data_t u, data_t b)
: _blk_iub_iter(blk_iub_iter), _u(u), _b(b)
{}
const_iterator (iterator i)
: _blk_iub_iter(i._blk_iub_iter), _u(i._u), _b(i._b)
{}
// accessors & modifiers:
reference operator* () const {
size_type iub = (*_blk_iub_iter).iub();
size_type blk = (*_blk_iub_iter).is_blocked();
if (!blk) return _u[iub]; else return _b[iub];
}
const_iterator& operator++ () { _blk_iub_iter++; return *this; }
const_iterator operator++ (int) { const_iterator tmp = *this; operator++(); return *this; }
const_iterator& operator+= (difference_type n) { _blk_iub_iter += n; return *this; }
const_iterator operator+ (difference_type n) const { const_iterator tmp = *this; return tmp += n; }
reference operator[] (size_type n) const { return *(*this + n); }
// comparators:
bool operator== (const const_iterator& j) const { return _blk_iub_iter == j._blk_iub_iter; }
bool operator!= (const const_iterator& j) const { return ! operator== (j); }
protected:
// data:
iter_t _blk_iub_iter;
data_t _u;
data_t _b;
};
template <class T, class M>
inline
typename field_basic<T,M>::const_iterator
field_basic<T,M>::begin_dof () const
{
return const_iterator (_V.data()._idof2blk_iub.begin(), _u.begin(), _b.begin());
}
template <class T, class M>
inline
typename field_basic<T,M>::const_iterator
field_basic<T,M>::end_dof () const
{
return const_iterator (_V.data()._idof2blk_iub.end(), _u.begin(), _b.begin());
}
// =================================================================================
// field: inlined
// =================================================================================
template <class T, class M>
inline
field_basic<T,M>::field_basic ()
: _V(),
_u(),
_b(),
_dis_dof_update_needed(true)
{
}
template <class T, class M>
void
field_basic<T,M>::dis_dof_update_needed() const
{
_dis_dof_update_needed = true;
}
template <class T, class M>
void
field_basic<T,M>::dis_dof_update() const
{
if (!_dis_dof_update_needed) return;
_dis_dof_update_needed = false;
if (is_distributed<M>::value) dis_dof_update_internal();
}
template<class T, class M>
inline
field_basic<T,M>&
field_basic<T,M>::operator= (const T& value)
{
check_macro (stamp() != "", "field=constant : uninitialized field in affectation");
std::fill (begin_dof(), end_dof(), value);
dis_dof_update_needed();
return *this;
}
template <class T, class M>
inline
T&
field_basic<T,M>::dof (size_type idof)
{
dis_dof_update_needed();
if (! _V.is_blocked(idof)) {
return _u [_V.iub(idof)];
} else {
return _b [_V.iub(idof)];
}
}
template <class T, class M>
inline
const T&
field_basic<T,M>::dof (size_type idof) const
{
if (! _V.is_blocked(idof)) {
return _u [_V.iub(idof)];
} else {
return _b [_V.iub(idof)];
}
}
template <class T, class M>
inline
T
field_basic<T,M>::min () const
{
T val = std::numeric_limits<T>::max();
for (const_iterator iter = begin_dof(), last = end_dof(); iter != last; iter++) {
val = std::min(val, *iter);
}
#ifdef _RHEOLEF_HAVE_MPI
if (is_distributed<M>::value) {
val = mpi::all_reduce (comm(), val, mpi::minimum<T>());
}
#endif // _RHEOLEF_HAVE_MPI
return val;
}
template <class T, class M>
inline
T
field_basic<T,M>::max () const
{
T val = std::numeric_limits<T>::min();
for (const_iterator iter = begin_dof(), last = end_dof(); iter != last; iter++) {
val = std::max(val, *iter);
}
#ifdef _RHEOLEF_HAVE_MPI
if (is_distributed<M>::value) {
val = mpi::all_reduce (comm(), val, mpi::maximum<T>());
}
#endif // _RHEOLEF_HAVE_MPI
return val;
}
template <class T, class M>
inline
T
field_basic<T,M>::min_abs () const
{
T val = std::numeric_limits<T>::max();
for (const_iterator iter = begin_dof(), last = end_dof(); iter != last; iter++) {
val = std::min(val, abs(*iter));
}
#ifdef _RHEOLEF_HAVE_MPI
if (is_distributed<M>::value) {
val = mpi::all_reduce (comm(), val, mpi::minimum<T>());
}
#endif // _RHEOLEF_HAVE_MPI
return val;
}
template <class T, class M>
inline
T
field_basic<T,M>::max_abs () const
{
T val = 0;
for (const_iterator iter = begin_dof(), last = end_dof(); iter != last; iter++) {
val = std::max(val, abs(*iter));
}
#ifdef _RHEOLEF_HAVE_MPI
if (is_distributed<M>::value) {
val = mpi::all_reduce (comm(), val, mpi::maximum<T>());
}
#endif // _RHEOLEF_HAVE_MPI
return val;
}
template <class T, class M>
inline
idiststream&
operator >> (idiststream& ips, field_basic<T,M>& uh)
{
return uh.get (ips);
}
template <class T, class M>
inline
odiststream&
operator << (odiststream& ops, const field_basic<T,M>& uh)
{
return uh.put (ops);
}
#ifdef TO_CLEAN
template <class T, class M>
inline
field_basic<T,M>
norm (const field_basic<T,M>& uh)
{
return sqrt(norm2(uh));
}
#endif // TO_CLEAN
}// namespace rheolef
# endif /* _RHEOLEF_FIELD_H */
|