aboutsummaryrefslogtreecommitdiff
path: root/src/Registration.c
blob: 715230443bffc7f53d20a1eae67525de9d97829d (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
 /*@@
   @file      Registration.c
   @date      Thu May 30 11:21:44 2002
   @author    Ian Hawke
   @desc 
   The external functions called (via function aliasing) by physics
   thorns to tell MoL that they want these GFs to be treated as a
   given type.
   @enddesc 
   @version   $Header$
 @@*/

#include "cctk.h"
#include "cctk_Parameters.h"

#include "ExternalVariables.h"

static const char *rcsid = "$Header$";

CCTK_FILEVERSION(CactusMoL2_MoL2_Registration_c);

/********************************************************************
 *********************     Local Data Types   ***********************
 ********************************************************************/

/********************************************************************
 ********************* Local Routine Prototypes *********************
 ********************************************************************/

/********************************************************************
 ***************** Scheduled Routine Prototypes *********************
 ********************************************************************/

/********************************************************************
 ********************* Other Routine Prototypes *********************
 ********************************************************************/

CCTK_INT MoL_RegisterEvolved(CCTK_INT EvolvedIndex, CCTK_INT RHSIndex);

CCTK_INT MoL_RegisterConstrained(CCTK_INT ConstrainedIndex);

CCTK_INT MoL_RegisterSaveAndRestore(CCTK_INT SandRIndex);

CCTK_INT MoL_RegisterEvolvedGroup(CCTK_INT EvolvedGroupIndex, 
                                 CCTK_INT RHSGroupIndex);

CCTK_INT MoL_RegisterConstrainedGroup(CCTK_INT ConstrainedGroupIndex);

CCTK_INT MoL_RegisterSaveAndRestoreGroup(CCTK_INT SandRGroupIndex);

/********************************************************************
 *********************     Local Data   *****************************
 ********************************************************************/

/********************************************************************
 *********************     External Routines   **********************
 ********************************************************************/

 /*@@
   @routine    MoL_RegisterEvolved
   @date       Thu May 30 11:36:59 2002
   @author     Ian Hawke
   @desc 
   Given the index of the GF to be evolved and the RHS GF, it stores 
   the indexes for later use together with various error checking.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

CCTK_INT MoL_RegisterEvolved(CCTK_INT EvolvedIndex, CCTK_INT RHSIndex)
{

  DECLARE_CCTK_PARAMETERS;

  CCTK_INT ierr, index, varused, numtimelevs1, numtimelevs2;
  
  /*
  printf("Arrived in MoLRegisterEvolved \n");
  printf("The indexes are %d and %d.\n",EvolvedIndex, RHSIndex);
  printf("These correspond to variables %s and %s.\n",
         CCTK_VarName(EvolvedIndex),CCTK_VarName(RHSIndex));
  printf("The pointer to EvolvedVariableIndex: %p\n",
         EvolvedVariableIndex);
  */

  numtimelevs1 = CCTK_NumTimeLevelsFromVarI(EvolvedIndex);
  numtimelevs2 = CCTK_NumTimeLevelsFromVarI(RHSIndex);
  
  if ( (numtimelevs1 < 0) || (numtimelevs2 < 0) ) 
  {
    CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for variable index %i", EvolvedIndex); 
    CCTK_WARN(0, "The index passed does not correspond to a GF.");
  }

  if (numtimelevs1 < 2)
  {
    CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for variable index %i name %s", EvolvedIndex, CCTK_VarName(EvolvedIndex));
    CCTK_WARN(0, "The GF passed only has one timelevel. It must have at least two.");
  }

  varused = 0;

  for (index = 0; (index < MoLNumEvolvedVariables)&&(!varused); index++)
  {
    varused = (EvolvedIndex == EvolvedVariableIndex[index]);
    /*
    printf("Registering %d. Checking index %d which is %d\n",EvolvedIndex,
           index,EvolvedVariableIndex[index]);
    */
  }

  if (varused)
  {

    CCTK_VWarn(2,__LINE__,__FILE__,"MoL",
               "The GF %s has already been registered as an evolved variable with RHS GF %s. The attempt to register with RHS GF %s will be ignored",
               CCTK_VarName(EvolvedIndex),
               CCTK_VarName(RHSVariableIndex[index-1]),
               CCTK_VarName(RHSIndex));

  }
  else
  {
  
    EvolvedVariableIndex[MoLNumEvolvedVariables] = EvolvedIndex;
    RHSVariableIndex[MoLNumEvolvedVariables] = RHSIndex;
  
    MoLNumEvolvedVariables++;
    
    if (MoLNumEvolvedVariables > MoL_Num_Evolved_Vars)
    {
      CCTK_WARN(0,"You have tried to register more evolved variables than the accumulator parameter MoL_Num_Evolved_Variables allows. Check that you are accumulating onto this parameter correctly");
    }

  }

  varused = 0;

  for (index = 0; (index < MoLNumConstrainedVariables)&&(!varused); index++)
  {
    varused = (EvolvedIndex == ConstrainedVariableIndex[index]);
  }
  
  if (varused)
  {
    for (index = varused; index < MoLNumConstrainedVariables-1; index++)
    {
      ConstrainedVariableIndex[index] = ConstrainedVariableIndex[index+1];
    }
    MoLNumConstrainedVariables--;
  }
    
  varused = 0;

  for (index = 0; (index < MoLNumSandRVariables)&&(!varused); index++)
  {
    varused = (EvolvedIndex == SandRVariableIndex[index]);
  }
  
  if (varused)
  {
    for (index = varused; index < MoLNumSandRVariables-1; index++)
    {
      SandRVariableIndex[index] = SandRVariableIndex[index+1];
    }
    MoLNumSandRVariables--;
  }
  
  return 0;

}

 /*@@
   @routine    MoL_RegisterConstrained
   @date       Thu May 30 12:35:58 2002
   @author     Ian Hawke
   @desc 
   Given the index of the GF, register it as a constrained variable.
   If there's only one timelevel then ignore it as there will be no
   rotation and so MoL doesn't have to do anything.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

CCTK_INT MoL_RegisterConstrained(CCTK_INT ConstrainedIndex)
{

  DECLARE_CCTK_PARAMETERS;
  
  CCTK_INT numtimelevs, varused, evolved, index;
  
  numtimelevs = CCTK_NumTimeLevelsFromVarI(ConstrainedIndex);

  if (numtimelevs < 1) {

    CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for constrained variable index %i", ConstrainedIndex);
    CCTK_WARN(0, "The index passed does not correspond to a GF.");

  }
  else if (numtimelevs > 1) {
    
    varused = 0;
    
    for (evolved = 0; (evolved < MoLNumEvolvedVariables)&&(!varused); evolved++)
    {
      varused = (EvolvedVariableIndex[evolved] == ConstrainedIndex);
    }
    
    for (evolved = 0; (evolved < MoLNumConstrainedVariables)&&(!varused); evolved++)
    {
      varused = (ConstrainedVariableIndex[evolved] == ConstrainedIndex);
    }
    
    if (!varused)
    {
      ConstrainedVariableIndex[MoLNumConstrainedVariables] = ConstrainedIndex;
      MoLNumConstrainedVariables++;

      if (MoLNumConstrainedVariables > MoL_Num_Constrained_Vars)
      {
        CCTK_WARN(0,"You have tried to register more evolved variables than the accumulator parameter MoL_Num_Evolved_Variables allows. Check that you are accumulating onto this parameter correctly");
      }

    }
    
    varused = 0;
    
    for (evolved = 0; (evolved < MoLNumSandRVariables)&&(!varused); evolved++)
    {
      varused = (SandRVariableIndex[evolved] == ConstrainedIndex);
    }
    
    if (varused)
    {
      for (index = evolved; index < MoLNumSandRVariables-1; index++)
      {
        SandRVariableIndex[index] = SandRVariableIndex[index+1];
      }
      MoLNumSandRVariables--;
    }
    
  }
  else
  {
    
    CCTK_VInfo(CCTK_THORNSTRING,
               "MoL will not treat variable %s as a constrained variable at it has only one timelevel. This should not cause problems with the evolution.\n", 
               CCTK_VarName(ConstrainedIndex));

  }
  
  return 0;
  
}

 /*@@
   @routine    MoL_RegisterSaveAndRestore
   @date       Thu May 30 12:37:40 2002
   @author     Ian Hawke
   @desc 
   Given a GF index store it for later use as a save and restore type.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

CCTK_INT MoL_RegisterSaveAndRestore(CCTK_INT SandRIndex)
{

  DECLARE_CCTK_PARAMETERS;

  CCTK_INT numtimelevs, varused, evolved;
  
  numtimelevs = CCTK_NumTimeLevelsFromVarI(SandRIndex);

  if (numtimelevs < 1) {

    CCTK_VWarn(1,__LINE__,__FILE__,"MoL","Warning for save and restore variable index %i", SandRIndex);
    CCTK_WARN(0, "The index passed does not correspond to a GF.");

  }
  else if (numtimelevs > 1) {
    
    varused = 0;
    
    for (evolved = 0; (evolved < MoLNumEvolvedVariables)&&(!varused); evolved++)
    {
      varused = (EvolvedVariableIndex[evolved] == SandRIndex);
    }
    
    for (evolved = 0; (evolved < MoLNumConstrainedVariables)&&(!varused); evolved++)
    {
      varused = (ConstrainedVariableIndex[evolved] == SandRIndex);
    }
    
    for (evolved = 0; (evolved < MoLNumSandRVariables)&&(!varused); evolved++)
    {
      varused = (SandRVariableIndex[evolved] == SandRIndex);
    }

    if (!varused)
    {
      SandRVariableIndex[MoLNumSandRVariables] = SandRIndex;
      MoLNumSandRVariables++;

      if (MoLNumSandRVariables > MoL_Num_SaveAndRestore_Vars)
      {
        CCTK_WARN(0,"You have tried to register more evolved variables than the accumulator parameter MoL_Num_Evolved_Variables allows. Check that you are accumulating onto this parameter correctly");
      }

    }
        
  }
  else
  {
    
    CCTK_VInfo(CCTK_THORNSTRING,
               "MoL will not treat variable %s as a save and restore variable at it has only one timelevel. This should not cause problems with the evolution.\n", 
               CCTK_VarName(SandRIndex));

  }
  
  return 0;
  
}

CCTK_INT MoL_RegisterEvolvedGroup(CCTK_INT EvolvedGroupIndex, 
                                 CCTK_INT RHSGroupIndex)
{

  CCTK_INT EvolvedGroupFirstVar, RHSGroupFirstVar, GroupNumVars, retval;
  CCTK_INT EvolvedVar, RHSVar;
  
  EvolvedGroupFirstVar = CCTK_FirstVarIndexI(EvolvedGroupIndex);
  if (EvolvedGroupFirstVar < 0)
  {
    CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
               "Evolved group index %i is not a real group index.",
              EvolvedGroupIndex);
  }
  
  RHSGroupFirstVar = CCTK_FirstVarIndexI(RHSGroupIndex);
  if (RHSGroupFirstVar < 0)
  {
    CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
               "RHS group index %d is not a real group index.",
              RHSGroupIndex);
  }

  GroupNumVars = CCTK_NumVarsInGroupI(EvolvedGroupIndex);
  if (CCTK_NumVarsInGroupI(RHSGroupIndex) != GroupNumVars)
  {
    CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
               "There are a different number of variables in evolved group %d and RHS group %d.", EvolvedGroupIndex, RHSGroupIndex);
  }
  
  retval = 0;
  
  for (EvolvedVar = EvolvedGroupFirstVar, RHSVar = RHSGroupFirstVar;
       EvolvedVar < EvolvedGroupFirstVar + GroupNumVars;
       EvolvedVar++, RHSVar++)
  {
    retval += MoLRegisterEvolved(EvolvedVar, RHSVar);
  }
  
  return retval;
}


CCTK_INT MoL_RegisterConstrainedGroup(CCTK_INT ConstrainedGroupIndex)
{
  
  CCTK_INT ConstrainedGroupFirstVar, GroupNumVars, retval;
  CCTK_INT ConstrainedVar;
  
  ConstrainedGroupFirstVar = CCTK_FirstVarIndexI(ConstrainedGroupIndex);
  if (ConstrainedGroupFirstVar < 0)
  {
    CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
               "Constrained group index %i is not a real group index.",
              ConstrainedGroupIndex);
  }
  
  GroupNumVars = CCTK_NumVarsInGroupI(ConstrainedGroupIndex);
  
  retval = 0;
  
  for (ConstrainedVar = ConstrainedGroupFirstVar;
       ConstrainedVar < ConstrainedGroupFirstVar + GroupNumVars;
       ConstrainedVar++)
  {
    retval += MoLRegisterConstrained(ConstrainedVar);
  }
  
  return retval;
}

CCTK_INT MoL_RegisterSaveAndRestoreGroup(CCTK_INT SandRGroupIndex)
{
  
  CCTK_INT SandRGroupFirstVar, GroupNumVars, retval;
  CCTK_INT SandRVar;
  
  SandRGroupFirstVar = CCTK_FirstVarIndexI(SandRGroupIndex);
  if (SandRGroupFirstVar < 0)
  {
    CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
               "Save and Restore group index %i is not a real group index.",
              SandRGroupIndex);
  }
  
  GroupNumVars = CCTK_NumVarsInGroupI(SandRGroupIndex);
  
  retval = 0;
  
  for (SandRVar = SandRGroupFirstVar;
       SandRVar < SandRGroupFirstVar + GroupNumVars;
       SandRVar++)
  {
    retval += MoLRegisterSaveAndRestore(SandRVar);
  }
  
  return retval;
}

/*
  Old function names. Just calls the new version. 
  Included for compatibility
*/

CCTK_INT MoL_RegisterVar(CCTK_INT molvarindex,CCTK_INT molrhsvarindex)
{
  return MoL_RegisterEvolved(molvarindex, molrhsvarindex);
}

CCTK_INT MoL_RegisterPrimitive(CCTK_INT primitiveindex)
{
  return MoL_RegisterConstrained(primitiveindex);
}

CCTK_INT MoL_RegisterDepends(CCTK_INT dependsindex)
{
  return MoL_RegisterSaveAndRestore(dependsindex);
}

CCTK_INT MoL_RegisterVarGroup(CCTK_INT groupindex,CCTK_INT rhsgroupindex)
{
  return MoL_RegisterEvolvedGroup(groupindex, rhsgroupindex);
}

CCTK_INT MoL_RegisterPrimitiveGroup(CCTK_INT groupindex)
{
  return MoL_RegisterConstrainedGroup(groupindex);
}

CCTK_INT MoL_RegisterDependsGroup(CCTK_INT groupindex)
{
  return MoL_RegisterConstrainedGroup(groupindex);
}

CCTK_INT MoL_ChangeVarToEvolved(CCTK_INT varindex, CCTK_INT rhsindex)
{
  return MoL_ChangeToEvolved(varindex, rhsindex);
}

CCTK_INT MoL_ChangeVarToDependent(CCTK_INT dependsindex)
{
  return MoL_ChangeToSaveAndRestore(dependsindex);
}

CCTK_INT MoL_ChangeVarToPrimitive(CCTK_INT primitiveindex)
{
  return MoL_ChangeToConstrained(primitiveindex);
}

CCTK_INT MoL_ChangeVarToNone(CCTK_INT removeindex)
{
  return MoL_ChangeToNone(removeindex);
}

/* 
  Fortran wrappers for the above functions. 
  Should be replaced by using function aliasing eventually.
*/

void CCTK_FCALL CCTK_FNAME(MoL_RegisterEvolved)(int *ierr, 
                                            CCTK_INT *EvolvedIndex,
                                            CCTK_INT *RHSIndex)
{
  *ierr = MoL_RegisterEvolved(*EvolvedIndex, *RHSIndex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_RegisterConstrained)(int *ierr, 
                                            CCTK_INT *EvolvedIndex)
{
  *ierr = MoL_RegisterConstrained(*EvolvedIndex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_RegisterSaveAndRestore)(int *ierr, 
                                            CCTK_INT *EvolvedIndex)
{
  *ierr = MoL_RegisterSaveAndRestore(*EvolvedIndex);
  return;
}

/* And for MoL compatibility... */

void CCTK_FCALL CCTK_FNAME(MoL_RegisterVar)(int *ierr, 
                                            CCTK_INT *molvarindex,
                                            CCTK_INT *molrhsvarindex)
{
  *ierr = MoL_RegisterVar(*molvarindex, *molrhsvarindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_RegisterDepends)(int *ierr, 
                                                CCTK_INT *moldependsindex)
{
  *ierr = MoL_RegisterDepends(*moldependsindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_RegisterPrimitive)(int *ierr, 
                                                CCTK_INT *molprimitiveindex)
{
  *ierr = MoL_RegisterPrimitive(*molprimitiveindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_RegisterVarGroup)(int *ierr, 
                                                 CCTK_INT *groupindex,
                                                 CCTK_INT *rhsgroupindex)
{
  *ierr = MoL_RegisterVarGroup(*groupindex, *rhsgroupindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_RegisterPrimitiveGroup)(int *ierr, 
                                                       CCTK_INT *groupindex)
{
  *ierr = MoL_RegisterPrimitiveGroup(*groupindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_RegisterDependsGroup)(int *ierr, 
                                                       CCTK_INT *groupindex)
{
  *ierr = MoL_RegisterDependsGroup(*groupindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_ChangeVarToEvolved)(int *ierr, 
                                                   CCTK_INT *varindex,
                                                   CCTK_INT *rhsindex)
{
  *ierr = MoL_ChangeVarToEvolved(*varindex, *rhsindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_ChangeVarToDependent)(int *ierr, 
                                                     CCTK_INT *dependsindex)
{
  *ierr = MoL_ChangeVarToDependent(*dependsindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_ChangeVarToPrimitive)(int *ierr, 
                                                     CCTK_INT *primitiveindex)
{
  *ierr = MoL_ChangeVarToPrimitive(*primitiveindex);
  return;
}

void CCTK_FCALL CCTK_FNAME(MoL_ChangeVarToNone)(int *ierr, 
                                                CCTK_INT *removeindex)
{
  *ierr = MoL_ChangeVarToNone(*removeindex);
  return;
}

/********************************************************************
 *********************     Local Routines   *************************
 ********************************************************************/