Skip to content

Commit d25f5ad

Browse files
committed
Combination of MoveInsertable fixes (#2170)
The master hashes of the commits cherry-picked and squashed to create this commit: PR #2720 783091b PR #2731 b67b3e4 0a1dbe4 73d58f2 3a61395 a6f14fe 00cd55d 1310eac A similar "combined" commit should also be cherry-picked onto the 1.4.x and 1.5.x series and new patch-level releases should be made for them, since this is an issue which prevents certain compilers from being able to compile the library. Refs #2710
1 parent 595b980 commit d25f5ad

File tree

4 files changed

+164
-104
lines changed

4 files changed

+164
-104
lines changed

include/numerics/parsed_fem_function.h

+70-33
Original file line numberDiff line numberDiff line change
@@ -40,6 +40,7 @@
4040
#include <sstream>
4141
#include <string>
4242
#include <vector>
43+
#include <memory>
4344

4445
namespace libMesh
4546
{
@@ -68,21 +69,16 @@ class ParsedFEMFunction : public FEMFunctionBase<Output>
6869
const std::vector<Output> * initial_vals=nullptr);
6970

7071
/**
71-
* This class contains a const reference so it can't be copy or move-assigned.
72+
* Constructors
73+
* - This class contains a const reference so it can't be default
74+
* copy or move-assigned. We manually implement the former
75+
* - This class contains unique_ptrs so it can't be default copy
76+
* constructed or assigned.
77+
* - It can still be default moved and deleted.
7278
*/
73-
ParsedFEMFunction & operator= (const ParsedFEMFunction &) = delete;
79+
ParsedFEMFunction & operator= (const ParsedFEMFunction &);
7480
ParsedFEMFunction & operator= (ParsedFEMFunction &&) = delete;
75-
76-
/**
77-
* The remaining 5 special functions can be safely defaulted.
78-
*
79-
* \note The underlying FunctionParserBase class has a copy
80-
* constructor, so this class should be default-constructible. And,
81-
* although FunctionParserBase's move constructor is deleted, _this_
82-
* class should still be move-constructible because
83-
* FunctionParserBase only appears in a vector.
84-
*/
85-
ParsedFEMFunction (const ParsedFEMFunction &) = default;
81+
ParsedFEMFunction (const ParsedFEMFunction &);
8682
ParsedFEMFunction (ParsedFEMFunction &&) = default;
8783
virtual ~ParsedFEMFunction () = default;
8884

@@ -168,9 +164,9 @@ class ParsedFEMFunction : public FEMFunctionBase<Output>
168164
_n_requested_hess_components;
169165
bool _requested_normals;
170166
#ifdef LIBMESH_HAVE_FPARSER
171-
std::vector<FunctionParserBase<Output>> parsers;
167+
std::vector<std::unique_ptr<FunctionParserBase<Output>>> parsers;
172168
#else
173-
std::vector<char> parsers;
169+
std::vector<char*> parsers;
174170
#endif
175171
std::vector<Output> _spacetime;
176172

@@ -222,6 +218,47 @@ ParsedFEMFunction<Output>::ParsedFEMFunction (const System & sys,
222218
}
223219

224220

221+
template <typename Output>
222+
inline
223+
ParsedFEMFunction<Output>::ParsedFEMFunction (const ParsedFEMFunction<Output> & other) :
224+
FEMFunctionBase<Output>(other),
225+
_sys(other._sys),
226+
_expression(other._expression),
227+
_subexpressions(other._subexpressions),
228+
_n_vars(other._n_vars),
229+
_n_requested_vars(other._n_requested_vars),
230+
_n_requested_grad_components(other._n_requested_grad_components),
231+
_n_requested_hess_components(other._n_requested_hess_components),
232+
_requested_normals(other._requested_normals),
233+
_spacetime(other._spacetime),
234+
_need_var(other._need_var),
235+
_need_var_grad(other._need_var_grad),
236+
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
237+
_need_var_hess(other._need_var_hess),
238+
#endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
239+
variables(other.variables),
240+
_additional_vars(other._additional_vars),
241+
_initial_vals(other._initial_vals)
242+
{
243+
this->reparse(_expression);
244+
}
245+
246+
247+
template <typename Output>
248+
inline
249+
ParsedFEMFunction<Output> &
250+
ParsedFEMFunction<Output>::operator= (const ParsedFEMFunction<Output> & other)
251+
{
252+
// We can only be assigned another ParsedFEMFunction defined on the same System
253+
libmesh_assert(&_sys == &other._sys);
254+
255+
// Use copy-and-swap idiom
256+
ParsedFEMFunction<Output> tmp(other);
257+
std::swap(tmp, *this);
258+
return *this;
259+
}
260+
261+
225262
template <typename Output>
226263
inline
227264
void
@@ -401,7 +438,7 @@ ParsedFEMFunction<Output>::operator() (const FEMContext & c,
401438
{
402439
eval_args(c, p, time);
403440

404-
return eval(parsers[0], "f", 0);
441+
return eval(*parsers[0], "f", 0);
405442
}
406443

407444

@@ -421,7 +458,7 @@ ParsedFEMFunction<Output>::operator() (const FEMContext & c,
421458
libmesh_assert_equal_to (size, parsers.size());
422459

423460
for (unsigned int i=0; i != size; ++i)
424-
output(i) = eval(parsers[i], "f", i);
461+
output(i) = eval(*parsers[i], "f", i);
425462
}
426463

427464

@@ -436,7 +473,7 @@ ParsedFEMFunction<Output>::component (const FEMContext & c,
436473
eval_args(c, p, time);
437474

438475
libmesh_assert_less (i, parsers.size());
439-
return eval(parsers[i], "f", i);
476+
return eval(*parsers[i], "f", i);
440477
}
441478

442479
template <typename Output>
@@ -479,16 +516,16 @@ ParsedFEMFunction<Output>::get_inline_value(const std::string & inline_var_name)
479516
#ifdef LIBMESH_HAVE_FPARSER
480517
// Parse and evaluate the new subexpression.
481518
// Add the same constants as we used originally.
482-
FunctionParserBase<Output> fp;
483-
fp.AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
484-
fp.AddConstant("pi", std::acos(Real(-1)));
485-
fp.AddConstant("e", std::exp(Real(1)));
519+
auto fp = libmesh_make_unique<FunctionParserBase<Output>>();
520+
fp->AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
521+
fp->AddConstant("pi", std::acos(Real(-1)));
522+
fp->AddConstant("e", std::exp(Real(1)));
486523
libmesh_error_msg_if
487-
(fp.Parse(new_subexpression, variables) != -1, // -1 for success
524+
(fp->Parse(new_subexpression, variables) != -1, // -1 for success
488525
"ERROR: FunctionParser is unable to parse modified expression: "
489-
<< new_subexpression << '\n' << fp.ErrorMsg());
526+
<< new_subexpression << '\n' << fp->ErrorMsg());
490527

491-
Output new_var_value = this->eval(fp, new_subexpression, 0);
528+
Output new_var_value = this->eval(*fp, new_subexpression, 0);
492529
#ifdef NDEBUG
493530
return new_var_value;
494531
#else
@@ -621,16 +658,16 @@ ParsedFEMFunction<Output>::partial_reparse (const std::string & expression)
621658
#ifdef LIBMESH_HAVE_FPARSER
622659
// Parse (and optimize if possible) the subexpression.
623660
// Add some basic constants, to Real precision.
624-
FunctionParserBase<Output> fp;
625-
fp.AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
626-
fp.AddConstant("pi", std::acos(Real(-1)));
627-
fp.AddConstant("e", std::exp(Real(1)));
661+
auto fp = libmesh_make_unique<FunctionParserBase<Output>>();
662+
fp->AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
663+
fp->AddConstant("pi", std::acos(Real(-1)));
664+
fp->AddConstant("e", std::exp(Real(1)));
628665
libmesh_error_msg_if
629-
(fp.Parse(_subexpressions.back(), variables) != -1, // -1 for success
666+
(fp->Parse(_subexpressions.back(), variables) != -1, // -1 for success
630667
"ERROR: FunctionParser is unable to parse expression: "
631-
<< _subexpressions.back() << '\n' << fp.ErrorMsg());
632-
fp.Optimize();
633-
parsers.push_back(fp);
668+
<< _subexpressions.back() << '\n' << fp->ErrorMsg());
669+
fp->Optimize();
670+
parsers.push_back(std::move(fp));
634671
#else
635672
libmesh_error_msg("ERROR: This functionality requires fparser!");
636673
#endif

include/numerics/parsed_function.h

+71-52
Original file line numberDiff line numberDiff line change
@@ -66,27 +66,13 @@ class ParsedFunction : public FunctionBase<Output>
6666
const std::vector<Output> * initial_vals=nullptr);
6767

6868
/**
69-
* This class cannot be (default) copy assigned because the
70-
* underlying FunctionParserADBase class does not define a custom
71-
* copy assignment operator, and manually manages memory.
69+
* Constructors
70+
* - This class contains unique_ptrs so it can't be default copy
71+
* constructed or assigned, only default moved and deleted.
7272
*/
73-
ParsedFunction & operator= (const ParsedFunction &) = delete;
73+
ParsedFunction (const ParsedFunction &);
74+
ParsedFunction & operator= (const ParsedFunction &);
7475

75-
/**
76-
* The remaining special functions can be defaulted for this class.
77-
*
78-
* \note Despite the fact that the underlying FunctionParserADBase
79-
* class is not move-assignable or move-constructible, it is still
80-
* possible for _this_ class to be move-assigned and
81-
* move-constructed, because FunctionParserADBase objects only
82-
* appear within std::vectors in this class, and std::vectors can
83-
* generally still be move-assigned and move-constructed even when
84-
* their contents cannot. There are some allocator-specific
85-
* exceptions to this, but it should be guaranteed to work for
86-
* std::allocator in C++14 and beyond. See also:
87-
* https://stackoverflow.com/q/42051917/659433
88-
*/
89-
ParsedFunction (const ParsedFunction &) = default;
9076
ParsedFunction (ParsedFunction &&) = default;
9177
ParsedFunction & operator= (ParsedFunction &&) = default;
9278
virtual ~ParsedFunction () = default;
@@ -183,18 +169,18 @@ class ParsedFunction : public FunctionBase<Output>
183169

184170
std::string _expression;
185171
std::vector<std::string> _subexpressions;
186-
std::vector<FunctionParserADBase<Output>> parsers;
172+
std::vector<std::unique_ptr<FunctionParserADBase<Output>>> parsers;
187173
std::vector<Output> _spacetime;
188174

189175
// derivative functions
190-
std::vector<FunctionParserADBase<Output>> dx_parsers;
176+
std::vector<std::unique_ptr<FunctionParserADBase<Output>>> dx_parsers;
191177
#if LIBMESH_DIM > 1
192-
std::vector<FunctionParserADBase<Output>> dy_parsers;
178+
std::vector<std::unique_ptr<FunctionParserADBase<Output>>> dy_parsers;
193179
#endif
194180
#if LIBMESH_DIM > 2
195-
std::vector<FunctionParserADBase<Output>> dz_parsers;
181+
std::vector<std::unique_ptr<FunctionParserADBase<Output>>> dz_parsers;
196182
#endif
197-
std::vector<FunctionParserADBase<Output>> dt_parsers;
183+
std::vector<std::unique_ptr<FunctionParserADBase<Output>>> dt_parsers;
198184
bool _valid_derivatives;
199185

200186
// Variables/values that can be parsed and handled by the function parser
@@ -225,6 +211,37 @@ ParsedFunction<Output,OutputGradient>::ParsedFunction (const std::string & expre
225211
}
226212

227213

214+
template <typename Output, typename OutputGradient>
215+
inline
216+
ParsedFunction<Output,OutputGradient>::ParsedFunction (const ParsedFunction<Output,OutputGradient> & other) :
217+
FunctionBase<Output>(other),
218+
_expression(other._expression),
219+
_subexpressions(other._subexpressions),
220+
_spacetime(other._spacetime),
221+
_valid_derivatives(other._valid_derivatives),
222+
variables(other.variables),
223+
_additional_vars(other._additional_vars),
224+
_initial_vals(other._initial_vals)
225+
{
226+
// parsers can be generated from scratch by reparsing expression
227+
this->reparse(this->_expression);
228+
this->_initialized = true;
229+
}
230+
231+
232+
233+
template <typename Output, typename OutputGradient>
234+
inline
235+
ParsedFunction<Output,OutputGradient> &
236+
ParsedFunction<Output,OutputGradient>::operator= (const ParsedFunction<Output,OutputGradient> & other)
237+
{
238+
// Use copy-and-swap idiom
239+
ParsedFunction<Output,OutputGradient> tmp(other);
240+
std::swap(tmp, *this);
241+
return *this;
242+
}
243+
244+
228245
template <typename Output, typename OutputGradient>
229246
inline
230247
void
@@ -262,7 +279,7 @@ Output
262279
ParsedFunction<Output,OutputGradient>::operator() (const Point & p, const Real time)
263280
{
264281
set_spacetime(p, time);
265-
return eval(parsers[0], "f", 0);
282+
return eval(*parsers[0], "f", 0);
266283
}
267284

268285
template <typename Output, typename OutputGradient>
@@ -271,7 +288,7 @@ Output
271288
ParsedFunction<Output,OutputGradient>::dot (const Point & p, const Real time)
272289
{
273290
set_spacetime(p, time);
274-
return eval(dt_parsers[0], "df/dt", 0);
291+
return eval(*dt_parsers[0], "df/dt", 0);
275292
}
276293

277294
template <typename Output, typename OutputGradient>
@@ -282,12 +299,12 @@ ParsedFunction<Output,OutputGradient>::gradient (const Point & p, const Real tim
282299
OutputGradient grad;
283300
set_spacetime(p, time);
284301

285-
grad(0) = eval(dx_parsers[0], "df/dx", 0);
302+
grad(0) = eval(*dx_parsers[0], "df/dx", 0);
286303
#if LIBMESH_DIM > 1
287-
grad(1) = eval(dy_parsers[0], "df/dy", 0);
304+
grad(1) = eval(*dy_parsers[0], "df/dy", 0);
288305
#endif
289306
#if LIBMESH_DIM > 2
290-
grad(2) = eval(dz_parsers[0], "df/dz", 0);
307+
grad(2) = eval(*dz_parsers[0], "df/dz", 0);
291308
#endif
292309

293310
return grad;
@@ -310,7 +327,7 @@ ParsedFunction<Output,OutputGradient>::operator()
310327
// The remaining locations in _spacetime are currently fixed at construction
311328
// but could potentially be made dynamic
312329
for (unsigned int i=0; i != size; ++i)
313-
output(i) = eval(parsers[i], "f", i);
330+
output(i) = eval(*parsers[i], "f", i);
314331
}
315332

316333
/**
@@ -330,7 +347,7 @@ ParsedFunction<Output,OutputGradient>::component (unsigned int i,
330347
// The remaining locations in _spacetime are currently fixed at construction
331348
// but could potentially be made dynamic
332349
libmesh_assert_less(i, parsers.size());
333-
return eval(parsers[i], "f", i);
350+
return eval(*parsers[i], "f", i);
334351
}
335352

336353
/**
@@ -541,51 +558,53 @@ ParsedFunction<Output,OutputGradient>::partial_reparse (const std::string & expr
541558

542559
// Parse (and optimize if possible) the subexpression.
543560
// Add some basic constants, to Real precision.
544-
FunctionParserADBase<Output> fp;
545-
fp.AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
546-
fp.AddConstant("pi", std::acos(Real(-1)));
547-
fp.AddConstant("e", std::exp(Real(1)));
561+
auto fp = libmesh_make_unique<FunctionParserADBase<Output>>();
562+
fp->AddConstant("NaN", std::numeric_limits<Real>::quiet_NaN());
563+
fp->AddConstant("pi", std::acos(Real(-1)));
564+
fp->AddConstant("e", std::exp(Real(1)));
548565
libmesh_error_msg_if
549-
(fp.Parse(_subexpressions.back(), variables) != -1, // -1 for success
566+
(fp->Parse(_subexpressions.back(), variables) != -1, // -1 for success
550567
"ERROR: FunctionParser is unable to parse expression: "
551-
<< _subexpressions.back() << '\n' << fp.ErrorMsg());
568+
<< _subexpressions.back() << '\n' << fp->ErrorMsg());
552569

553570
// use of derivatives is optional. suppress error output on the console
554571
// use the has_derivatives() method to check if AutoDiff was successful.
555572
// also enable immediate optimization
556-
fp.SetADFlags(FunctionParserADBase<Output>::ADSilenceErrors |
573+
fp->SetADFlags(FunctionParserADBase<Output>::ADSilenceErrors |
557574
FunctionParserADBase<Output>::ADAutoOptimize);
558575

559576
// optimize original function
560-
fp.Optimize();
561-
parsers.push_back(fp);
577+
fp->Optimize();
562578

563579
// generate derivatives through automatic differentiation
564-
FunctionParserADBase<Output> dx_fp(fp);
565-
if (dx_fp.AutoDiff("x") != -1) // -1 for success
580+
auto dx_fp = libmesh_make_unique<FunctionParserADBase<Output>>(*fp);
581+
if (dx_fp->AutoDiff("x") != -1) // -1 for success
566582
_valid_derivatives = false;
567-
dx_parsers.push_back(dx_fp);
583+
dx_parsers.push_back(std::move(dx_fp));
568584
#if LIBMESH_DIM > 1
569-
FunctionParserADBase<Output> dy_fp(fp);
570-
if (dy_fp.AutoDiff("y") != -1) // -1 for success
585+
auto dy_fp = libmesh_make_unique<FunctionParserADBase<Output>>(*fp);
586+
if (dy_fp->AutoDiff("y") != -1) // -1 for success
571587
_valid_derivatives = false;
572-
dy_parsers.push_back(dy_fp);
588+
dy_parsers.push_back(std::move(dy_fp));
573589
#endif
574590
#if LIBMESH_DIM > 2
575-
FunctionParserADBase<Output> dz_fp(fp);
576-
if (dz_fp.AutoDiff("z") != -1) // -1 for success
591+
auto dz_fp = libmesh_make_unique<FunctionParserADBase<Output>>(*fp);
592+
if (dz_fp->AutoDiff("z") != -1) // -1 for success
577593
_valid_derivatives = false;
578-
dz_parsers.push_back(dz_fp);
594+
dz_parsers.push_back(std::move(dz_fp));
579595
#endif
580-
FunctionParserADBase<Output> dt_fp(fp);
581-
if (dt_fp.AutoDiff("t") != -1) // -1 for success
596+
auto dt_fp = libmesh_make_unique<FunctionParserADBase<Output>>(*fp);
597+
if (dt_fp->AutoDiff("t") != -1) // -1 for success
582598
_valid_derivatives = false;
583-
dt_parsers.push_back(dt_fp);
599+
dt_parsers.push_back(std::move(dt_fp));
584600

585601
// If at end, use nextstart=maxSize. Else start at next
586602
// character.
587603
nextstart = (end == std::string::npos) ?
588604
std::string::npos : end + 1;
605+
606+
// Store fp for later use
607+
parsers.push_back(std::move(fp));
589608
}
590609
}
591610

0 commit comments

Comments
 (0)