aboutsummaryrefslogtreecommitdiff
path: root/Carpet/LoopControl/src/loopcontrol.c
blob: dd7d3f57fe2afbc55299fc8f8f800b9a3dd90964 (plain)
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
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
#include <assert.h>
#include <errno.h>
#include <float.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>

#include <sys/time.h>

#ifdef _OPENMP
#  include <omp.h>
#endif

#include <cctk_Parameters.h>

#include "loopcontrol.h"



#ifndef _OPENMP
/* Replacements for some OpenMP routines if OpenMP is not available */

static inline
int
omp_get_thread_num (void)
{
  return 0;
}

static inline
int
omp_get_num_threads (void)
{
  return 1;
}

static inline
double
omp_get_wtime (void)
{
  struct timeval tv;
  gettimeofday (& tv, NULL);
  return tv.tv_sec + 1.0e-6 * tv.tv_usec;
}

#endif



/* Linked list of all loop statistics structures */
struct lc_statmap_t * lc_statmap_list = NULL;



/* Find all possible thread topologies */
static
void
find_thread_topologies (struct lc_topology_t * restrict const topologies,
                        const int maxntopologies,
                        int * restrict const ntopologies,
                        int const nthreads)
{
  * ntopologies = 0;
  for (int nk=1; nk<=nthreads; ++nk) {
    if (nthreads % nk == 0) {
      for (int nj=1; nj<=nthreads/nk; ++nj) {
        if (nthreads % (nj*nk) == 0) {
          int const ni = nthreads/(nj*nk);
          if (nthreads == ni*nj*nk) {
            assert (* ntopologies < maxntopologies);
            topologies[* ntopologies].ni = ni;
            topologies[* ntopologies].nj = nj;
            topologies[* ntopologies].nk = nk;
            ++ * ntopologies;
          }
        }
      }
    }
  }
}



/* Find "good" tiling specifications */
static
int tiling_compare (const void * const a, const void * const b)
{
  struct lc_tiling_t const * const aa = a;
  struct lc_tiling_t const * const bb = b;
  return aa->npoints - bb->npoints;
}

static
void
find_tiling_specifications (struct lc_tiling_t * restrict const tilings,
                            const int maxntilings,
                            int * restrict const ntilings,
                            int const npoints)
{
  /* In order to reduce the number of possible tilings, require that
     the step sizes differ by more than 10% */
  double const distance_factor = 1.1;
  /* Determine the "good" step sizes in two passes: first small step
     sizes from 1 up to snpoints, then large step sizes from npoints
     down to snpoints+1 */
  int const snpoints = floor (sqrt (npoints));
  /* For N grid points and a minimum spacing factor F, there are at
     most log(N) / log(F) possible tilings.  There will be fewer,
     since the actual spacings will be rounded up to integers.  */
  
  * ntilings = 0;
  
  /* Small step sizes */
  int minnpoints = 0;
  for (int n=1; n<=snpoints; ++n) {
    if ((double) n > minnpoints * distance_factor) {
      assert (* ntilings < maxntilings);
      tilings[* ntilings].npoints = n;
      minnpoints = n;
      ++ * ntilings;
    }
  }
  
  /* Large step sizes */
  int maxnpoints = 1000000000;
  for (int n=npoints; n>snpoints; --n) {
    if (n * distance_factor < (double) maxnpoints) {
      assert (* ntilings < maxntilings);
      tilings[* ntilings].npoints = n;
      maxnpoints = n;
      ++ * ntilings;
    }
  }
  
  /* Sort */
  qsort (tilings, * ntilings, sizeof * tilings, tiling_compare);
}



/* Initialise control parameter set statistics */
static
void
lc_stattime_init (struct lc_stattime_t * restrict const lt,
                  struct lc_statset_t * restrict const ls,
                  int const topology,
                  int const tiling[3])
{
  DECLARE_CCTK_PARAMETERS;
  
  /* Check arguments */
  assert (lt);
  assert (ls);
  
  /*** Topology ***************************************************************/
  
  lt->topology = topology;
  
  if (topology == -1) {
    
    /* User-specified topology */
    lt->inthreads = -1;
    lt->jnthreads = -1;
    lt->knthreads = -1;
    
  } else {
    
    assert (topology >= 0 && topology < ls->ntopologies);
    
    lt->inthreads = ls->topologies[lt->topology].ni;
    lt->jnthreads = ls->topologies[lt->topology].nj;
    lt->knthreads = ls->topologies[lt->topology].nk;
    
  }
  
  if (debug) {
    printf ("Thread topology #%d [%d,%d,%d]\n",
            lt->topology, lt->inthreads, lt->jnthreads, lt->knthreads);
  }
  
  /* Assert thread topology consistency */
  if (lt->topology != -1) {
    assert (lt->inthreads >= 1);
    assert (lt->jnthreads >= 1);
    assert (lt->knthreads >= 1);
    assert (lt->inthreads * lt->jnthreads * lt->knthreads == ls->num_threads);
  }
  
  /*** Tilings ****************************************************************/
  
  for (int d=0; d<3; ++d) {
    lt->tiling[d] = tiling[d];
    if (tiling[d] != -1) {
      assert (tiling[d] >= 0 && tiling[d] < ls->ntilings[d]);
    }
  }
  
  if (tiling[0] != -1) {
    lt->inpoints = ls->tilings[0][lt->tiling[0]].npoints;
  }
  if (tiling[1] != -1) {
    lt->jnpoints = ls->tilings[1][lt->tiling[1]].npoints;
  }
  if (tiling[2] != -1) {
    lt->knpoints = ls->tilings[2][lt->tiling[2]].npoints;
  }
  
  if (debug) {
    printf ("Tiling stride [%d,%d,%d]\n",
            lt->inpoints, lt->jnpoints, lt->knpoints);
  }
  
  /* Assert tiling specification consistency */
  if (tiling[0] != -1) {
    assert (lt->inpoints > 0);
  }
  if (tiling[1] != -1) {
    assert (lt->jnpoints > 0);
  }
  if (tiling[2] != -1) {
    assert (lt->knpoints > 0);
  }
  
  
  
  /* Initialise statistics */
  lt->time_count      = 0.0;
  lt->time_setup_sum  = 0.0;
  lt->time_setup_sum2 = 0.0;
  lt->time_calc_sum   = 0.0;
  lt->time_calc_sum2  = 0.0;
  
  /* Append to loop statistics list */
/*   _Pragma ("omp critical") { */
    lt->next = ls->stattime_list;
    ls->stattime_list = lt;
/*   } */
}

static
struct lc_stattime_t *
lc_stattime_find (struct lc_statset_t * restrict const ls,
                  int const topology,
                  int const tiling[3])
{
  assert (ls);
  
  struct lc_stattime_t * lt;
  
  for (lt = ls->stattime_list; lt; lt = lt->next) {
    if (lt->topology == topology &&
        lt->tiling[0] == tiling[0] &&
        lt->tiling[1] == tiling[1] &&
        lt->tiling[2] == tiling[2])
    {
      break;
    }
  }
  
  if (! lt) {
    lt = malloc (sizeof * lt);
    lc_stattime_init (lt, ls, topology, tiling);
  }
  
  return lt;
}



/* Initialise user parameter set statistics */
static
void
lc_statset_init (struct lc_statset_t * restrict const ls,
                 struct lc_statmap_t * restrict const lm,
                 int const num_threads,
                 int const npoints[3])
{
  DECLARE_CCTK_PARAMETERS;
  
  /* Check arguments */
  assert (ls);
  assert (lm);
  assert (num_threads >= 1);
  for (int d=0; d<3; ++d) {
    assert (npoints[d] >= 0);
  }
  
  /*** Threads ****************************************************************/
  
  ls->num_threads = num_threads;
  
  /* For up to 1024 threads, there are at most 270 possible
     topologies */
  int const maxntopologies = 1000;
  if (debug) {
    printf ("Running on %d threads\n", ls->num_threads);
  }
  ls->topologies = malloc (maxntopologies * sizeof * ls->topologies);
  find_thread_topologies
    (ls->topologies, maxntopologies, & ls->ntopologies, ls->num_threads);
#if 0
  ls->topologies =
    realloc (ls->topologies, ls->ntopologies * sizeof * ls->topologies);
#endif
  if (debug) {
    printf ("Found %d possible thread topologies\n", ls->ntopologies);
    for (int n = 0; n < ls->ntopologies; ++n) {
      printf ("   %2d: %2d %2d %2d\n",
              n,
              ls->topologies[n].ni, ls->topologies[n].nj, ls->topologies[n].nk);
    }
  }
  
  /*** Tilings ****************************************************************/
  
  for (int d=0; d<3; ++d) {
    ls->npoints[d] = npoints[d];
  }
  
  /* For up to 1000000 grid points, there are at most 126 possible
     tilings (assuming a minimum distance of 10%) */
  int const maxntilings = 1000;
  for (int d=0; d<3; ++d) {
    if (debug) {
      printf ("Dimension %d: %d points\n", d, ls->npoints[d]);
    }
    ls->tilings[d] = malloc (maxntilings * sizeof * ls->tilings[d]);
    find_tiling_specifications
      (ls->tilings[d], maxntilings, & ls->ntilings[d], ls->npoints[d]);
#if 0
    ls->tilings[d] =
      realloc (ls->tilings[d], ls->ntilings[d] * sizeof * ls->tilings[d]);
#endif
    if (debug) {
      printf ("   Found %d possible tilings\n", ls->ntilings[d]);
      printf ("     ");
      for (int n = 0; n < ls->ntilings[d]; ++n) {
        printf (" %d", ls->tilings[d][n].npoints);
      }
      printf ("\n");
    }
  }
  
  
  
  /* Initialise list */
  ls->stattime_list = NULL;
  
  /* Append to loop statistics list */
/*   _Pragma ("omp critical") { */
    ls->next = lm->statset_list;
    lm->statset_list = ls;
/*   } */
}

static
struct lc_statset_t *
lc_statset_find (struct lc_statmap_t * restrict const lm,
                 int const num_threads,
                 int const npoints[3])
{
  assert (lm);
  
  struct lc_statset_t * ls;
  
  for (ls = lm->statset_list; ls; ls = ls->next) {
    if (ls->num_threads == num_threads &&
        ls->npoints[0] == npoints[0] &&
        ls->npoints[1] == npoints[1] &&
        ls->npoints[2] == npoints[2])
    {
      break;
    }
  }
  
  if (! ls) {
    ls = malloc (sizeof * ls);
    lc_statset_init (ls, lm, num_threads, npoints);
  }
  
  return ls;
}
                         


/* Initialise loop statistics */
void
lc_statmap_init (struct lc_statmap_t * restrict const lm,
                 char const * restrict const name)
{
  /* Check arguments */
  assert (lm);
  
  /* Set name */
  lm->name = strdup (name);
  
  /* Initialise list */
  lm->statset_list = NULL;
  
  /* Append to loop statistics list */
  lm->next = lc_statmap_list;
  lc_statmap_list = lm;
}



void
lc_control_init (struct lc_control_t * restrict const lc,
                 struct lc_statmap_t * restrict const lm,
                 int const imin, int const jmin, int const kmin,
                 int const imax, int const jmax, int const kmax,
                 int const ilsh, int const jlsh, int const klsh)
{
  DECLARE_CCTK_PARAMETERS;
  
  /* Check arguments */
  assert (lc);
  
  /* Timer */
  lc->time_setup_begin = omp_get_wtime();
  
  /* Check arguments */
  assert (imin >= 0 && imax <= ilsh && ilsh >= 0);
  assert (jmin >= 0 && jmax <= jlsh && jlsh >= 0);
  assert (kmin >= 0 && kmax <= klsh && klsh >= 0);
  
  /* Copy arguments */
  lc->imin = imin;
  lc->jmin = jmin;
  lc->kmin = kmin;
  lc->imax = imax;
  lc->jmax = jmax;
  lc->kmax = kmax;
  lc->ilsh = ilsh;
  lc->jlsh = jlsh;
  lc->klsh = klsh;
  
  
  
  struct lc_statset_t * restrict ls;
  _Pragma ("omp single copyprivate (ls)") {
    /* Get number of threads */
    int const num_threads = omp_get_num_threads();
    
    /* Calculate number of points */
    int npoints[3];
    npoints[0] = lc_max (imax - imin, 0);
    npoints[1] = lc_max (jmax - jmin, 0);
    npoints[2] = lc_max (kmax - kmin, 0);
    
    ls = lc_statset_find (lm, num_threads, npoints);
  }
  
  
  
  struct lc_stattime_t * restrict lt;
  _Pragma ("omp single copyprivate (lt)") {
    
    /* Select topology */
    
    int topology;
    
    if (lc_inthreads != -1 || lc_jnthreads != -1 || lc_knthreads != -1)
    {
      /* User-specified thread topology */
      
      if (lc_inthreads == -1 || lc_jnthreads == -1 || lc_knthreads == -1) {
        CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "Illegal thread topology [%d,%d,%d] specified\n",
                    (int)lc_inthreads, (int)lc_jnthreads, (int)lc_knthreads);
      }
      if (lc_inthreads * lc_jnthreads * lc_knthreads != ls->num_threads) {
        CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "Specified thread topology [%d,%d,%d] is not compatible with the number of threads %d\n",
                    (int)lc_inthreads, (int)lc_jnthreads, (int)lc_knthreads,
                    ls->num_threads);
      }
      
      topology = -1;
      
    } else {
      
      /* Split in the k direction */
      topology = ls->ntopologies - 1;
      
    }
    
    /* Select tiling */
    
    int tiling[3];
    
    if (lc_inpoints != -1) {
      /* User-specified tiling */
      tiling[0] = -1;
    } else {
      tiling[0] = ls->ntilings[0] - 1; /* as many points as possible */
    }
    
    if (lc_jnpoints != -1) {
      /* User-specified tiling */
      tiling[1] = -1;
    } else {
      if (cycle_j_tilings) {
        /* cycle through all tilings */
        static int count = 0;
        tiling[1] = (count ++) % ls->ntilings[1];
      } else {
        /* as few points as possible */
        tiling[1] = 0;
      }
    }
    
    if (lc_knpoints != -1) {
      /* User-specified tiling */
      tiling[2] = -1;
    } else {
      tiling[2] = ls->ntilings[2] - 1; /* as many points as possible */
    }
    
    /* Find or create database entry */
    
    lt = lc_stattime_find (ls, topology, tiling);
    
    /* Topology */
    
    if (topology == -1) {
      /* User-specified topology */
      lt->inthreads = lc_inthreads;
      lt->jnthreads = lc_jnthreads;
      lt->knthreads = lc_knthreads;
    }
    
    /* Assert thread topology consistency */
    assert (lt->inthreads >= 1);
    assert (lt->jnthreads >= 1);
    assert (lt->knthreads >= 1);
    assert (lt->inthreads * lt->jnthreads * lt->knthreads == ls->num_threads);
    
    /* Tilings */
    
    if (tiling[0] == -1) {
      /* User-specified tiling */
      lt->inpoints = lc_inpoints;
    }
    if (tiling[1] == -1) {
      /* User-specified tiling */
      lt->jnpoints = lc_jnpoints;
    }
    if (tiling[2] == -1) {
      /* User-specified tiling */
      lt->knpoints = lc_knpoints;
    }
    
    /* Assert tiling specification consistency */
    assert (lt->inpoints > 0);
    assert (lt->jnpoints > 0);
    assert (lt->knpoints > 0);
    
  } /* omp single */
  
  
  
  lc->statmap  = lm;
  lc->statset  = ls;
  lc->stattime = lt;
  
  

  /*** Threads ****************************************************************/
  
  
  /* Thread loop settings */
  lc->iiimin = imin;
  lc->jjjmin = jmin;
  lc->kkkmin = kmin;
  lc->iiimax = imax;
  lc->jjjmax = jmax;
  lc->kkkmax = kmax;
  lc->iiistep = (lc->iiimax - lc->iiimin + lt->inthreads-1) / lt->inthreads;
  lc->jjjstep = (lc->jjjmax - lc->jjjmin + lt->jnthreads-1) / lt->jnthreads;
  lc->kkkstep = (lc->kkkmax - lc->kkkmin + lt->knthreads-1) / lt->knthreads;
  
  /* Find location of current thread */
  lc->thread_num  = omp_get_thread_num();
  int c = lc->thread_num;
  int const ci = c % lt->inthreads; c /= lt->inthreads;
  int const cj = c % lt->jnthreads; c /= lt->jnthreads;
  int const ck = c % lt->knthreads; c /= lt->knthreads;
  assert (c == 0);
  lc->iii = lc->iiimin + ci * lc->iiistep;
  lc->jjj = lc->jjjmin + cj * lc->jjjstep;
  lc->kkk = lc->kkkmin + ck * lc->kkkstep;
  
  
  
  /*** Tilings ****************************************************************/
  
  /* Tiling loop settings */
  lc->iimin = lc->iii;
  lc->jjmin = lc->jjj;
  lc->kkmin = lc->kkk;
  lc->iimax = lc_min (lc->iii + lc->iiistep, lc->iiimax);
  lc->jjmax = lc_min (lc->jjj + lc->jjjstep, lc->jjjmax);
  lc->kkmax = lc_min (lc->kkk + lc->kkkstep, lc->kkkmax);
  lc->iistep = lt->inpoints;
  lc->jjstep = lt->jnpoints;
  lc->kkstep = lt->knpoints;
  
  
  
  /****************************************************************************/
  
  /* Timer */
  lc->time_calc_begin = omp_get_wtime();
}



void
lc_control_finish (struct lc_control_t * restrict const lc)
{
  /* Timer */
  double const time_calc_end   = omp_get_wtime();
  double const time_calc_begin = lc->time_calc_begin;
  
  double const time_setup_end   = time_calc_begin;
  double const time_setup_begin = lc->time_setup_begin;
  
  double const time_setup_sum  = time_setup_end - time_setup_begin;
  double const time_setup_sum2 = pow (time_setup_sum, 2);
  
  double const time_calc_sum  = time_calc_end - time_calc_begin;
  double const time_calc_sum2 = pow (time_calc_sum, 2);
  
  /* Update statistics */
  struct lc_stattime_t * restrict const lt = lc->stattime;
  _Pragma ("omp critical") {
    lt->time_count += 1.0;
    
    lt->time_setup_sum  += time_setup_sum;
    lt->time_setup_sum2 += time_setup_sum2;
    
    lt->time_calc_sum  += time_calc_sum;
    lt->time_calc_sum2 += time_calc_sum2;
  }
}



static
double
avg (double const c, double const s)
{
  if (c == 0.0) return 0.0;
  return s / c;
}

static
double
stddev (double const c, double const s, double const s2)
{
  if (c == 0.0) return 0.0;
  return sqrt (s2 / c - pow (s / c, 2));
}



void
lc_printstats (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_PARAMETERS;
  if (! verbose) return;
  
  int nmaps = 0;
  for (struct lc_statmap_t * lm = lc_statmap_list; lm; lm = lm->next) {
    printf ("statmap #%d \"%s\":\n",
            nmaps,
            lm->name);
    int nsets = 0;
    for (struct lc_statset_t * ls = lm->statset_list; ls; ls = ls->next) {
      printf ("   statset #%d nthreads=%d npoints=[%d,%d,%d]\n",
              nsets,
              ls->num_threads, ls->npoints[0], ls->npoints[1], ls->npoints[2]);
      double sum_setup = 0.0;
      double sum_calc  = 0.0;
      double min_calc = DBL_MAX;
      int imin_calc = -1;
      int ntimes = 0;
      for (struct lc_stattime_t * lt = ls->stattime_list; lt; lt = lt->next) {
        printf ("      stattime #%d topology=%d [%d,%d,%d] tiling=[%d,%d,%d]\n",
                ntimes,
                lt->topology, lt->inthreads, lt->jnthreads, lt->knthreads,
                lt->inpoints, lt->jnpoints, lt->knpoints);
        printf ("         count: %g   setup: %g   calc: %g\n",
                lt->time_count, lt->time_setup_sum, lt->time_calc_sum);
        sum_setup += lt->time_setup_sum;
        sum_calc  += lt->time_calc_sum;
        if (lt->time_calc_sum < min_calc) {
          min_calc = lt->time_calc_sum;
          imin_calc = ntimes;
        }
        ++ ntimes;
      }
      printf ("      total setup: %g   total calc: %g   min calc: %g (#%d)\n",
              sum_setup, sum_calc, min_calc, imin_calc);
      ++ nsets;
    }
    ++ nmaps;
  }
}



CCTK_FCALL
void
CCTK_FNAME (lc_statmap_init) (struct lc_statmap_t * restrict const lm,
                              ONE_FORTSTRING_ARG)
{
  ONE_FORTSTRING_CREATE (name);
  lc_statmap_init (lm, name);
  free (name);
}

CCTK_FCALL
void
CCTK_FNAME (lc_control_init) (struct lc_control_t * restrict const lc,
                              struct lc_statmap_t * restrict const lm,
                              int const * restrict const imin,
                              int const * restrict const jmin,
                              int const * restrict const kmin,
                              int const * restrict const imax,
                              int const * restrict const jmax,
                              int const * restrict const kmax,
                              int const * restrict const ilsh,
                              int const * restrict const jlsh,
                              int const * restrict const klsh)
{
  lc_control_init (lc, lm,
                   * imin - 1, * jmin - 1, * kmin - 1,
                   * imax, * jmax, * kmax,
                   * ilsh, * jlsh, * klsh);
}

CCTK_FCALL
void
CCTK_FNAME (lc_control_finish) (struct lc_control_t * restrict const lc)
{
  lc_control_finish (lc);
}