This file is indexed.

/usr/include/dune/pdelab/multistep/method.hh is in libdune-pdelab-dev 2.0.0-1.

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
// -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=8 sw=2 sts=2:

#ifndef DUNE_PDELAB_MULTISTEP_METHOD_HH
#define DUNE_PDELAB_MULTISTEP_METHOD_HH

#include <iomanip>
#include <iostream>

#include <dune/common/ios_state.hh>
#include <dune/common/shared_ptr.hh>
#include <dune/common/timer.hh>

#include <dune/pdelab/multistep/parameter.hh>

namespace Dune {
  namespace PDELab {

    //! \addtogroup MultiStepMethods
    //! \{

    //! Do one step of a multi-step time-stepping scheme
    /**
     * \tparam T          Type to represent time values
     * \tparam MGOS       Assembler for multi-step instationary problems
     * \tparam PDESOLVER  Solver problem in each step (typically Newton)
     * \tparam TrialV     Vector type to represent coefficients of solutions
     * \tparam TestV      Vector type to represent residuals
     */
    template<class T, class MGOS, class PDESolver, class TrialV,
             class TestV = TrialV>
    class MultiStepMethod {
      typedef MultiStepParameterInterface<T, MGOS::order> Parameters;

      const Parameters* parameters;
      MGOS& mgos;
      PDESolver& pdeSolver;
      unsigned verbosity;
      unsigned step;

    public:
      //! construct a new multi-step scheme
      /**
       * \param parameters_ Parameter object.
       * \param mgos_       Assembler object (MultiStepGridOperatorSpace).
       * \param pdesolver_  Solver object (typically Newton).
       *
       * The contructed method object stores references to the object it is
       * constructed with, so these objects should be valid for as long as the
       * constructed object is used.
       */
      MultiStepMethod(const Parameters& parameters_,
                      MGOS& mgos_, PDESolver& pdeSolver_)
        : parameters(&parameters_), mgos(mgos_), pdeSolver(pdeSolver_),
          verbosity(1), step(1)
      {
        if(mgos.trialGridFunctionSpace().gridView().comm().rank() > 0)
          verbosity = 0;
      }

      //! change verbosity level; 0 means completely quiet
      void setVerbosityLevel (int level)
      {
        if(mgos.trialGridFunctionSpace().gridView().comm().rank() == 0)
          verbosity = level;
      }

      //! redefine the method to be used; can be done before every step
      void setMethod(const Parameters &parameters_) {
        parameters = &parameters_;
      }

      //! do one step;
      /*
       * \param time      Start of time step
       * \param dt        Suggested time step size
       * \param oldValues Vector of pointers to the old values.  Must support
       *                  the expression *oldvalues[i], which should yield a
       *                  reference to the old value at time time-i*dt.
       * \param xnew      Where to store the new value.
       *
       * \return Actual time step size
       */
      template<class OldValues>
      T apply(T time, T dt, const OldValues& oldValues, TrialV& xnew)
      {
        Timer allTimer;
        Timer subTimer;

        // save formatting attributes
        ios_base_all_saver format_attribute_saver(std::cout);

        if(verbosity >= 1)
          std::cout << "TIME STEP [" << parameters->name() << "] "
                    << std::setw(6) << step
                    << " time (from): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time
                    << " dt: "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << dt
                    << " time (to): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time+dt
                    << std::endl;

        // prepare assembler
        if(verbosity >= 2) {
          std::cout << "== prepare assembler" << std::endl;
          subTimer.reset();
        }
        mgos.preStep(time,dt, oldValues);
        if(verbosity >= 2)
          std::cout << "== prepare assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // solve
        if(verbosity >= 2) {
          std::cout << "== apply solver" << std::endl;
          subTimer.reset();
        }
        pdeSolver.apply(xnew);
        if(verbosity >= 2)
          std::cout << "== apply solver (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // postprocessing in the assembler
        if(verbosity >= 2) {
          std::cout << "== cleanup assembler" << std::endl;
          subTimer.reset();
        }
        mgos.postStep();
        if(verbosity >= 2)
          std::cout << "== cleanup assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        ++step;

        if(verbosity >= 2)
          std::cout << "== time step done (" << allTimer.elapsed() << "s)"
                    << std::endl;

        return dt;
      }

      //! do one step;
      /* This is a version which interpolates constraints at the start of each
       * stage
       *
       * \param time      Start of time step
       * \param dt        Suggested time step size
       * \param oldValues Vector of pointers to the old values.  Must support
       *                  the expression *oldValues[i], which should yield a
       *                  reference to the old value at time time-i*dt.
       * \param f         Function to interpolate boundary conditions from.
       *                  Should support the method setTime().
       * \param xnew      Where to store the new value.
       *
       * \return Actual time step size
       */
      template<typename OldValues, typename F>
      T apply(T time, T dt, const OldValues& oldValues, F& f, TrialV& xnew)
      {
        Timer allTimer;
        Timer subTimer;

        // save formatting attributes
        ios_base_all_saver format_attribute_saver(std::cout);

        if(verbosity >= 1)
          std::cout << "TIME STEP [" << parameters->name() << "] "
                    << std::setw(6) << step
                    << " time (from): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time
                    << " dt: "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << dt
                    << " time (to): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time+dt
                    << std::endl;

        // prepare assembler
        if(verbosity >= 2) {
          std::cout << "== prepare assembler" << std::endl;
          subTimer.reset();
        }
        mgos.preStep(time, dt, oldValues);
        if(verbosity >= 2)
          std::cout << "== prepare assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // set boundary conditions and initial value
        if(verbosity >= 2) {
          std::cout << "== setup result vector" << std::endl;
          subTimer.reset();
        }
        f.setTime(time+dt);
        mgos.interpolate(*oldValues[0],f,xnew);
        if(verbosity >= 2)
          std::cout << "== setup result vector (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // solve stage
        if(verbosity >= 2) {
          std::cout << "== apply solver" << std::endl;
          subTimer.reset();
        }
        pdeSolver.apply(xnew);
        if(verbosity >= 2)
          std::cout << "== apply solver (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // postprocessing in the assembler
        if(verbosity >= 2) {
          std::cout << "== cleanup assembler" << std::endl;
          subTimer.reset();
        }
        mgos.postStep();
        if(verbosity >= 2)
          std::cout << "== cleanup assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        step++;

        if(verbosity >= 2)
          std::cout << "== time step done (" << allTimer.elapsed() << "s)"
                    << std::endl;

        return dt;
      }

      //! do one step (with caching)
      /**
       * \param time Start of time step
       * \param dt   Time step size
       *
       * \return A shared_ptr to the new value
       *
       * The old values are expected in the cache of the GridOperatorSpace.
       * The computed value is store in the cache as well.
       */
      shared_ptr<const TrialV> apply(T time, T dt)
      {
        Timer allTimer;
        Timer subTimer;

        // save formatting attributes
        ios_base_all_saver format_attribute_saver(std::cout);

        if(verbosity >= 1)
          std::cout << "TIME STEP [" << parameters->name() << "] "
                    << std::setw(6) << step
                    << " time (from): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time
                    << " dt: "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << dt
                    << " time (to): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time+dt
                    << std::endl;

        // prepare assembler
        if(verbosity >= 2) {
          std::cout << "== prepare assembler" << std::endl;
          subTimer.reset();
        }
        mgos.preStep(step, time, dt);
        if(verbosity >= 2)
          std::cout << "== prepare assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // create vector, using last time step as start
        if(verbosity >= 2) {
          std::cout << "== setup result vector" << std::endl;
          subTimer.reset();
        }
        shared_ptr<TrialV>
          xnew(new TrialV(*mgos.getCache()->getUnknowns(step-1)));
        if(verbosity >= 2)
          std::cout << "== setup result vector (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // solve
        if(verbosity >= 2) {
          std::cout << "== apply solver" << std::endl;
          subTimer.reset();
        }
        pdeSolver.apply(*xnew);
        if(verbosity >= 2)
          std::cout << "== apply solver (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // postprocessing in the assembler
        if(verbosity >= 2) {
          std::cout << "== cleanup assembler" << std::endl;
          subTimer.reset();
        }
        mgos.postStep();
        if(verbosity >= 2)
          std::cout << "== cleanup assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // store result for next step
        mgos.getCache()->setUnknowns(step, xnew);

        ++step;

        if(verbosity >= 2)
          std::cout << "== time step done (" << allTimer.elapsed() << "s)"
                    << std::endl;

        return xnew;
      }

      //! do one step (with caching)
      /**
       * This is a version which interpolates constraints at the start of each
       * stage
       *
       * \param time Start of time step
       * \param dt   Time step size
       * \param f    Function to interpolate boundary conditions from.
       *             Should support the method setTime().
       *
       * \return A shared_ptr to the new value
       *
       * The old values are expected in the cache of the GridOperatorSpace.
       * The computed value is store in the cache as well.
       */
      template<typename F>
      shared_ptr<const TrialV> apply(T time, T dt, F& f)
      {
        Timer allTimer;
        Timer subTimer;

        // save formatting attributes
        ios_base_all_saver format_attribute_saver(std::cout);

        if(verbosity >= 1)
          std::cout << "TIME STEP [" << parameters->name() << "] "
                    << std::setw(6) << step
                    << " time (from): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time
                    << " dt: "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << dt
                    << " time (to): "
                    << std::setw(12) << std::setprecision(4) << std::scientific
                    << time+dt
                    << std::endl;

        // prepare assembler
        if(verbosity >= 2) {
          std::cout << "== prepare assembler" << std::endl;
          subTimer.reset();
        }
        mgos.preStep(step, time, dt);
        if(verbosity >= 2)
          std::cout << "== prepare assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // setup vector
        if(verbosity >= 2) {
          std::cout << "== setup result vector" << std::endl;
          subTimer.reset();
        }
        shared_ptr<TrialV> xnew(new TrialV(mgos.trialGridFunctionSpace()));
        // set boundary conditions and initial value
        f.setTime(time+dt);
        mgos.interpolate(*mgos.getCache()->getUnknowns(step-1),f,*xnew);
        if(verbosity >= 2)
          std::cout << "== setup result vector (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // solve stage
        if(verbosity >= 2) {
          std::cout << "== apply solver" << std::endl;
          subTimer.reset();
        }
        pdeSolver.apply(*xnew);
        if(verbosity >= 2)
          std::cout << "== apply solver (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // postprocessing in the assembler
        if(verbosity >= 2) {
          std::cout << "== cleanup assembler" << std::endl;
          subTimer.reset();
        }
        mgos.postStep();
        if(verbosity >= 2)
          std::cout << "== cleanup assembler (" << subTimer.elapsed() << "s)"
                    << std::endl;

        // store result for next step
        mgos.getCache()->setUnknowns(step, xnew);

        ++step;

        if(verbosity >= 2)
          std::cout << "== time step done (" << allTimer.elapsed() << "s)"
                    << std::endl;

        return xnew;
      }
    };

    //! \} group MultiStepMethods
  } // namespace PDELab
} // namespace Dune

#endif // DUNE_PDELAB_MULTISTEP_METHOD_HH