Skip to content

Commit e2cccac

Browse files
committed
Add option to disable caching in FE::reinit()
1 parent 358b251 commit e2cccac

File tree

2 files changed

+15
-10
lines changed

2 files changed

+15
-10
lines changed

include/fe/fe.h

+7
Original file line numberDiff line numberDiff line change
@@ -793,6 +793,13 @@ class FE : public FEGenericBase<typename FEOutputType<T>::type>
793793
std::vector<Number> & nodal_soln_on_side,
794794
bool add_p_level = true);
795795

796+
/**
797+
* Whether we cache the node locations on the last
798+
* element we computed on to try to avoid calling
799+
* init_shape_functions and compute_shape_functions
800+
*/
801+
bool caching;
802+
796803
/**
797804
* An array of the node locations on the last
798805
* element we computed on

src/fe/fe.C

+8-10
Original file line numberDiff line numberDiff line change
@@ -47,6 +47,7 @@ namespace libMesh
4747
template <unsigned int Dim, FEFamily T>
4848
FE<Dim,T>::FE (const FEType & fet) :
4949
FEGenericBase<typename FEOutputType<T>::type> (Dim,fet),
50+
caching(!libMesh::on_command_line("--disable-caching")),
5051
last_side(INVALID_ELEM),
5152
last_edge(INVALID_ELEM)
5253
{
@@ -207,13 +208,6 @@ void FE<Dim,T>::reinit(const Elem * elem,
207208
this->_fe_map->template init_reference_to_physical_map<Dim>
208209
(this->qrule->get_points(), elem);
209210
this->init_shape_functions (this->qrule->get_points(), elem);
210-
211-
if (this->shapes_need_reinit())
212-
{
213-
cached_nodes.resize(elem->n_nodes());
214-
for (auto n : elem->node_index_range())
215-
cached_nodes[n] = elem->point(n);
216-
}
217211
}
218212
else
219213
{
@@ -238,12 +232,16 @@ void FE<Dim,T>::reinit(const Elem * elem,
238232
this->_fe_map->template init_reference_to_physical_map<Dim>
239233
(this->qrule->get_points(), elem);
240234
this->init_shape_functions (this->qrule->get_points(), elem);
241-
cached_nodes.resize(elem->n_nodes());
242-
for (auto n : elem->node_index_range())
243-
cached_nodes[n] = elem->point(n);
244235
}
245236
}
246237

238+
if (this->shapes_need_reinit() && caching && !cached_nodes_still_fit)
239+
{
240+
cached_nodes.resize(elem->n_nodes());
241+
for (auto n : elem->node_index_range())
242+
cached_nodes[n] = elem->point(n);
243+
}
244+
247245
// The shape functions correspond to the qrule
248246
this->shapes_on_quadrature = true;
249247
}

0 commit comments

Comments
 (0)