Source

PetIGA / src / petigaaxis.c

Lisandro Dalcin f3d7eb8 











Lisandro Dalcin 8e13c37 

Lisandro Dalcin 1dae083 

Lisandro Dalcin d0cfa31 
Lisandro Dalcin f3d7eb8 



Lisandro Dalcin 8e13c37 
Lisandro Dalcin d0cfa31 
Lisandro Dalcin 8e13c37 



Lisandro Dalcin f3d7eb8 













Lisandro Dalcin 8e13c37 
Lisandro Dalcin af24fc2 
Lisandro Dalcin f3d7eb8 




Lisandro Dalcin 67ed8c1 






Lisandro Dalcin 8e13c37 
Lisandro Dalcin 1dae083 

Lisandro Dalcin 8e13c37 








Lisandro Dalcin af24fc2 



Lisandro Dalcin 8e13c37 

Lisandro Dalcin af24fc2 
Lisandro Dalcin 8e13c37 
Lisandro Dalcin 67ed8c1 



Lisandro Dalcin f3d7eb8 










Nathan Collier 9d5677a 

Lisandro Dalcin 4248d1e 
Nathan Collier 9d5677a 











Lisandro Dalcin f3d7eb8 





Lisandro Dalcin 3cfad3a 

Lisandro Dalcin f3d7eb8 
Lisandro Dalcin d0cfa31 
Lisandro Dalcin f3d7eb8 


Lisandro Dalcin d0cfa31 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin 3cfad3a 
Lisandro Dalcin d0cfa31 

Lisandro Dalcin 3cfad3a 
Lisandro Dalcin d0cfa31 

Lisandro Dalcin 3cfad3a 
Lisandro Dalcin f3d7eb8 

















Nathan Collier 9d5677a 

Lisandro Dalcin 4248d1e 
Nathan Collier 9d5677a 









Lisandro Dalcin f3d7eb8 



















Lisandro Dalcin 073d86a 
Nathan Collier 9d5677a 

Lisandro Dalcin 4248d1e 
Nathan Collier 9d5677a 









Lisandro Dalcin 073d86a 
Lisandro Dalcin f3d7eb8 




Lisandro Dalcin d0cfa31 



Lisandro Dalcin f3d7eb8 




Lisandro Dalcin 073d86a 

Lisandro Dalcin f3d7eb8 









Lisandro Dalcin 7b14e15 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin d0cfa31 
Lisandro Dalcin f3d7eb8 


Lisandro Dalcin 7b14e15 
Lisandro Dalcin cf4cf7d 
Lisandro Dalcin 5b18640 

Lisandro Dalcin 073d86a 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin 5b18640 
Lisandro Dalcin 7b14e15 

Lisandro Dalcin d0cfa31 






Lisandro Dalcin 7b14e15 



Lisandro Dalcin d0cfa31 






Lisandro Dalcin 5b18640 
Lisandro Dalcin 7b14e15 

Lisandro Dalcin d0cfa31 
Lisandro Dalcin f3d7eb8 



Lisandro Dalcin 7b14e15 
Lisandro Dalcin cf4cf7d 
Lisandro Dalcin d0cfa31 
Lisandro Dalcin 4248d1e 
Lisandro Dalcin d0cfa31 








Lisandro Dalcin 4248d1e 
Lisandro Dalcin f3d7eb8 
















Lisandro Dalcin 3a5d1ae 
























Lisandro Dalcin d0cfa31 







































Lisandro Dalcin 3a5d1ae 
Lisandro Dalcin 5b18640 
Lisandro Dalcin 5cb4f23 
Lisandro Dalcin 5b18640 

Lisandro Dalcin 3f4015e 
Lisandro Dalcin 5b18640 





Lisandro Dalcin 3f4015e 
Lisandro Dalcin 5b18640 


Lisandro Dalcin 073d86a 
Lisandro Dalcin 5b18640 








Lisandro Dalcin 3f4015e 


Lisandro Dalcin 5b18640 
Lisandro Dalcin 073d86a 
Lisandro Dalcin 3f4015e 
Lisandro Dalcin 5b18640 


Lisandro Dalcin cf4cf7d 





Lisandro Dalcin d0cfa31 

Lisandro Dalcin 5b18640 










Lisandro Dalcin e16ea00 
Lisandro Dalcin 5b18640 



Lisandro Dalcin d0cfa31 
Lisandro Dalcin 4248d1e 
Lisandro Dalcin d0cfa31 



Lisandro Dalcin 4248d1e 
Lisandro Dalcin 5b18640 



Lisandro Dalcin f3d7eb8 
Nathan Collier 9d5677a 

Lisandro Dalcin 4248d1e 
Nathan Collier 9d5677a 

















Lisandro Dalcin e1d8929 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin 5b18640 
Lisandro Dalcin 3f4015e 
Lisandro Dalcin f3d7eb8 




Lisandro Dalcin 3f4015e 
Lisandro Dalcin 67ed8c1 
Lisandro Dalcin e1d8929 

Lisandro Dalcin 073d86a 
Lisandro Dalcin e1d8929 
Lisandro Dalcin f3d7eb8 
Nathan Collier 9d5677a 
Lisandro Dalcin 5b18640 


Lisandro Dalcin e1d8929 


Lisandro Dalcin f3d7eb8 
Lisandro Dalcin cf4cf7d 
Lisandro Dalcin 5b18640 
Lisandro Dalcin 3f4015e 
Lisandro Dalcin e1d8929 
Lisandro Dalcin 5b18640 
Lisandro Dalcin cf4cf7d 





Lisandro Dalcin d0cfa31 

Lisandro Dalcin f3d7eb8 




Lisandro Dalcin 3f4015e 
Lisandro Dalcin 5b18640 
Lisandro Dalcin e1d8929 
Lisandro Dalcin f3d7eb8 


Lisandro Dalcin e16ea00 
Lisandro Dalcin 5b18640 
Lisandro Dalcin f3d7eb8 


Lisandro Dalcin d0cfa31 
Lisandro Dalcin 4248d1e 
Lisandro Dalcin d0cfa31 



Lisandro Dalcin 4248d1e 
Lisandro Dalcin f3d7eb8 



Lisandro Dalcin 392bc0f 

Lisandro Dalcin f3d7eb8 
Lisandro Dalcin cf4cf7d 

Lisandro Dalcin af24fc2 
Lisandro Dalcin f3d7eb8 


Lisandro Dalcin 5b18640 
Lisandro Dalcin 073d86a 
Lisandro Dalcin d0cfa31 
Lisandro Dalcin 5b18640 
Lisandro Dalcin f3d7eb8 
Lisandro Dalcin cf4cf7d 





Lisandro Dalcin d0cfa31 
Lisandro Dalcin 8e13c37 
Lisandro Dalcin cf4cf7d 
Lisandro Dalcin 8e13c37 

Lisandro Dalcin cf4cf7d 




Lisandro Dalcin f3d7eb8 

Lisandro Dalcin cf4cf7d 
Lisandro Dalcin 1dae083 
Lisandro Dalcin cf4cf7d 







Lisandro Dalcin 1dae083 
Lisandro Dalcin cf4cf7d 






Lisandro Dalcin 1dae083 


























  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
#include "petiga.h"

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisCreate"
PetscErrorCode IGAAxisCreate(IGAAxis *_axis)
{
  IGAAxis        axis;
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(_axis,1);
  ierr = PetscNew(struct _n_IGAAxis,_axis);CHKERRQ(ierr);
  (*_axis)->refct = 1; axis = *_axis;

  /* */
  axis->periodic = PETSC_FALSE;
  /* */
  ierr = PetscMalloc1(2,PetscReal,&axis->U);CHKERRQ(ierr);
  axis->p = 0;
  axis->m = 1;
  axis->U[0] = -0.5;
  axis->U[1] = +0.5;
  /* */
  ierr = PetscMalloc1(1,PetscInt,&axis->span);CHKERRQ(ierr);
  axis->nnp = 1;
  axis->nel = 1;
  axis->span[0] = 0;

  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisDestroy"
PetscErrorCode IGAAxisDestroy(IGAAxis *_axis)
{
  IGAAxis        axis;
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(_axis,1);
  axis = *_axis; *_axis = 0;
  if (!axis) PetscFunctionReturn(0);
  if (--axis->refct > 0) PetscFunctionReturn(0);
  ierr = PetscFree(axis->U);CHKERRQ(ierr);
  ierr = PetscFree(axis->span);CHKERRQ(ierr);
  ierr = PetscFree(axis);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisReset"
PetscErrorCode IGAAxisReset(IGAAxis axis)
{
  PetscErrorCode ierr;
  PetscFunctionBegin;
  if (!axis) PetscFunctionReturn(0);
  PetscValidPointer(axis,1);

  axis->periodic = PETSC_FALSE;

  if (axis->m != 1) {
    ierr = PetscFree(axis->U);CHKERRQ(ierr);
    ierr = PetscMalloc1(2,PetscReal,&axis->U);CHKERRQ(ierr);
  }
  axis->p = 0;
  axis->m = 1;
  axis->U[0] = -0.5;
  axis->U[1] = +0.5;

  if (axis->nel != 1) {
    ierr = PetscFree(axis->span);CHKERRQ(ierr);
    ierr = PetscMalloc1(1,PetscInt,&axis->span);CHKERRQ(ierr);
  }
  axis->nnp = 1;
  axis->nel = 1;
  axis->span[0] = 0;

  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisReference"
PetscErrorCode IGAAxisReference(IGAAxis axis)
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  axis->refct++;
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisCopy"
/*@
   IGAAxisCopy - Copies an axis. axis <-- base

   Logically Collective on IGAAxis

   Input Parameters:
.  base - the axis

   Input Parameters:
.  axis - the copy

   Level: normal

.keywords: IGA, axis, copy
@*/
PetscErrorCode IGAAxisCopy(IGAAxis base,IGAAxis axis)
{
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(base,1);
  PetscValidPointer(axis,2);
  if(base == axis) PetscFunctionReturn(0);

  axis->periodic = base->periodic;

  axis->p = base->p;
  axis->m = base->m;
  ierr = PetscFree(axis->U);CHKERRQ(ierr);
  ierr = PetscMalloc1(axis->m+1,PetscReal,&axis->U);CHKERRQ(ierr);
  ierr = PetscMemcpy(axis->U,base->U,(axis->m+1)*sizeof(PetscReal));CHKERRQ(ierr);

  axis->nnp = base->nnp;
  axis->nel = base->nel;
  ierr = PetscFree(axis->span);CHKERRQ(ierr);
  ierr = PetscMalloc1(axis->nel,PetscInt,&axis->span);CHKERRQ(ierr);
  ierr = PetscMemcpy(axis->span,base->span,axis->nel*sizeof(PetscInt));CHKERRQ(ierr);

  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisDuplicate"
PetscErrorCode IGAAxisDuplicate(IGAAxis base,IGAAxis *axis)
{
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(base,1);
  PetscValidPointer(axis,2);
  ierr = PetscNew(struct _n_IGAAxis,axis);CHKERRQ(ierr);
  ierr = IGAAxisCopy(base,*axis);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisSetPeriodic"
/*@
   IGAAxisSetPeriodic - Sets the axis periodicity

   Logically Collective on IGAAxis

   Input Parameters:
+  axis - the IGAAxis context
-  periodic - TRUE or FALSE

   Level: normal

.keywords: IGA, axis, periodic
@*/
PetscErrorCode IGAAxisSetPeriodic(IGAAxis axis,PetscBool periodic)
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  axis->periodic = periodic ? PETSC_TRUE : PETSC_FALSE;
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisGetPeriodic"
PetscErrorCode IGAAxisGetPeriodic(IGAAxis axis,PetscBool *periodic)
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  PetscValidPointer(periodic,2);
  *periodic = axis->periodic;
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisSetDegree"
/*@
   IGAAxisSetDegree - Sets the axis degree

   Logically Collective on IGAAxis

   Input Parameters:
+  axis - the IGAAxis context
-  p - the polynomial degree

   Level: normal

.keywords: IGA, axis, degree
@*/
PetscErrorCode IGAAxisSetDegree(IGAAxis axis,PetscInt p)
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  if (p < 1)
    SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,
             "Polynomial degree must be greater than one, got %D",p);
  if (axis->p > 0 && axis->m > 1 && axis->p != p)
    SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,
             "Cannot change degree to %D after it was set to %D",p,axis->p);
  axis->p = p;
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisGetDegree"
PetscErrorCode IGAAxisGetDegree(IGAAxis axis,PetscInt *p)
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  PetscValidPointer(p,2);
  *p = axis->p;
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisSetKnots"
PetscErrorCode IGAAxisSetKnots(IGAAxis axis,PetscInt m,const PetscReal U[])
{
  PetscInt       p,n,k;
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  PetscValidPointer(U,3);

  if (axis->p < 1)
    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,
            "Must call IGAAxisSetDegree() first");
  if (m < 2*axis->p+1)
    SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,
             "Number of knots must be at least %D, got %D",
             2*(axis->p+1),m+1);

  p = axis->p;
  n = m -p - 1;

  for (k=1; k<=m;) {
    PetscInt i = k, s = 1;
    /* check increasing sequence */
    if (U[k-1] > U[k])
      SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,
               "Knot sequence must be increasing, got U[%D]=%G > U[%D]=%G",
               k-1,U[k-1],k,U[k]);
    /* check multiplicity */
    while (++k < m && U[i] == U[k]) s++;
    if (s > p)
      SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,
               "Knot U[%D]=%G has multiplicity %D greater than polynomial degree %D",
               i,U[i],s,p);
  }

  if (m != axis->m) {
    PetscReal *V;
    ierr = PetscMalloc1(m+1,PetscReal,&V);CHKERRQ(ierr);
    ierr = PetscFree(axis->U);CHKERRQ(ierr);
    axis->m = m;
    axis->U = V;
  }
  ierr = PetscMemcpy(axis->U,U,(m+1)*sizeof(PetscReal));CHKERRQ(ierr);

  axis->nel = 0;
  ierr = PetscFree(axis->span);CHKERRQ(ierr);
  ierr = IGAAxisGetSpans(axis,&axis->nel,&axis->span);CHKERRQ(ierr);

  if (axis->periodic) {
    PetscInt s = 1;
    while(s < p && U[m-p] == U[m-p+s]) s++;
    axis->nnp = n-p+s;
  } else {
    axis->nnp = n+1;
  }

  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisGetKnots"
PetscErrorCode IGAAxisGetKnots(IGAAxis axis,PetscInt *m,PetscReal *U[])
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  if (m) PetscValidPointer(m,2);
  if (U) PetscValidPointer(U,3);
  if (m) *m = axis->m;
  if (U) *U = axis->U;
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisGetLimits"
PetscErrorCode IGAAxisGetLimits(IGAAxis axis,PetscReal *Ui,PetscReal *Uf)
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  if (Ui) PetscValidRealPointer(Ui,2);
  if (Uf) PetscValidRealPointer(Uf,3);
  if (Ui) *Ui = axis->U[axis->p];
  if (Uf) *Uf = axis->U[axis->m-axis->p];
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisGetSizes"
PetscErrorCode IGAAxisGetSizes(IGAAxis axis,PetscInt *nel,PetscInt *nnp)
{
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  if (nel) PetscValidIntPointer(nel,2);
  if (nnp) PetscValidIntPointer(nnp,3);
  if (nel) *nel = axis->nel;
  if (nnp) *nnp = axis->nnp;
  PetscFunctionReturn(0);
}

EXTERN_C_BEGIN
extern PetscInt IGA_SpanCount(PetscInt n,PetscInt p,const PetscReal U[]);
extern PetscInt IGA_SpanIndex(PetscInt n,PetscInt p,const PetscReal U[],PetscInt index[]);
EXTERN_C_END

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisGetSpans"
PetscErrorCode IGAAxisGetSpans(IGAAxis axis,PetscInt *nel,PetscInt *span[])
{
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  if (nel)  PetscValidIntPointer(nel,2);
  if (span) PetscValidPointer(span,3);
  if (!axis->span) {
    PetscInt p = axis->p;
    PetscInt m = axis->m;
    PetscInt n = m - p - 1;
    axis->nel = IGA_SpanCount(n,p,axis->U);
    ierr = PetscMalloc1(axis->nel,PetscInt,&axis->span);CHKERRQ(ierr);
    (void)IGA_SpanIndex(n,p,axis->U,axis->span);
  }
  if (nel)  *nel  = axis->nel;
  if (span) *span = axis->span;
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisInit"
PetscErrorCode IGAAxisInit(IGAAxis axis,PetscInt p,PetscInt m,const PetscReal U[])
{
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  axis->p = 0;
  ierr = IGAAxisSetDegree(axis,p);CHKERRQ(ierr);
  ierr = IGAAxisSetKnots(axis,m,U);CHKERRQ(ierr);
  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisInitBreaks"
PetscErrorCode IGAAxisInitBreaks(IGAAxis axis,PetscInt nu,const PetscReal u[],PetscInt C)
{
  PetscInt       i,j,k;
  PetscInt       p,s,n,m,r;
  PetscReal      *U;
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  PetscValidPointer(u,3);

  if (C == PETSC_DECIDE) C = axis->p-1;

  if (axis->p < 1)
    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,
            "Must call IGAAxisSetDegree() first");
  if (nu < 2)
    SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,
             "Number of breaks must be at least two, got %D",nu);
  for (i=1; i<nu; i++)
    if (u[i-1] >= u[i])
      SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,
               "Break sequence must be strictly increasing, "
               "got u[%D]=%G %s u[%D]=%G",
               i-1,u[i-1],i,u[i],(u[i]==u[i-1])?"==":">");
  if (C < 0 || C >= axis->p)
    SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,
             "Continuity must be in range [0,%D], got %D",axis->p-1,C);

  p = axis->p; /* polynomial degree */
  s = p - C; /* multiplicity */
  r = nu - 1; /* last break index */
  m = 2*(p+1) + (r-1)*s - 1; /* last knot index */
  n = m - p - 1; /* last basis function index */

  if (m != axis->m) {
    ierr = PetscMalloc1(m+1,PetscReal,&U);CHKERRQ(ierr);
    ierr = PetscFree(axis->U);CHKERRQ(ierr);
    axis->m = m;
    axis->U = U;
  } else
    U = axis->U;

  for(k=0; k<=p; k++) { /* open part */
    U[k]   = u[0];
    U[m-k] = u[r];
  }
  for(i=1; i<=r-1; i++) { /* r-1 breaks */
    for(j=0; j<s; j++)    /* s times */
      U[k++] = u[i];
  }
  if (axis->periodic) {
    for(k=0; k<=C; k++) { /* periodic part */
      U[C-k]   = U[p] - U[m-p] + U[n-k];
      U[m-C+k] = U[m-p] - U[p] + U[p+1+k];
    }
  }

  axis->nel = r;
  ierr = PetscFree(axis->span);CHKERRQ(ierr);
  ierr = PetscMalloc1(axis->nel,PetscInt,&axis->span);CHKERRQ(ierr);
  for(i=0; i<axis->nel; i++) axis->span[i] = p + i*s;

  axis->nnp = axis->periodic ? n-C : n+1;

  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisInitUniform"
/*@
   IGAAxisInitUniform - Initializes an axis with uniformly spaces knots.

   Logically Collective on IGAAxis

   Input Parameters:
+  axis - the IGAAxis context
.  N - the number of equally-spaced, nonzero spans (elements)
.  Ui - the initial knot value
.  Uf - the final knot value
-  C - the global continuity order

   Notes: You must have called IGAAxisSetDegree() prior to this
   command. Creates a function space which consists of N spans of
   piecewise polynomials of the degree set with continuity order C at
   element interfaces.

   Level: normal

.keywords: IGA, axis, initialize, uniform
@*/
PetscErrorCode IGAAxisInitUniform(IGAAxis axis,PetscInt N,PetscReal Ui,PetscReal Uf,PetscInt C)
{
  PetscInt       i,j,k;
  PetscInt       p,s,n,m,r;
  PetscReal      *U;
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(axis,1);

  if (C == PETSC_DECIDE) C = axis->p-1;

  if (axis->p < 1)
    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,
            "Must call IGAAxisSetDegree() first");
  if (N < 1)
    SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,
             "Number of elements must be greater than zero, got %D",N);
  if (Ui >= Uf)
    SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,
             "Initial value %G must be less than final value %G",Ui,Uf);
  if (C < 0 || C >= axis->p)
    SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,
             "Continuity must be in range [0,%D], got %D",axis->p-1,C);

  p = axis->p; /* polynomial degree */
  s = p - C; /* multiplicity */
  r = N ; /* last break index */
  m = 2*(p+1) + (N-1)*s - 1; /* last knot index */
  n = m - p - 1; /* last basis function index */

  if (m != axis->m) {
    ierr = PetscMalloc1(m+1,PetscReal,&U);CHKERRQ(ierr);
    ierr = PetscFree(axis->U);CHKERRQ(ierr);
    axis->m = m;
    axis->U = U;
  } else
    U = axis->U;

  for(k=0; k<=p; k++) { /* open part */
    U[k]   = Ui;
    U[m-k] = Uf;
  }
  for(i=1; i<=r-1; i++) { /* (N-1) breaks */
    for(j=1; j<=s; j++)     /* s times */
      U[k++] = Ui + i * ((Uf-Ui)/N);
  }
  if (axis->periodic) {
    for(k=0; k<=C; k++) { /* periodic part */
      U[C-k]   = U[p] - U[m-p] + U[n-k];
      U[m-C+k] = U[m-p] - U[p] + U[p+1+k];
    }
  }

  axis->nel = r;
  ierr = PetscFree(axis->span);CHKERRQ(ierr);
  ierr = PetscMalloc1(axis->nel,PetscInt,&axis->span);CHKERRQ(ierr);
  for(i=0; i<axis->nel; i++) axis->span[i] = p + i*s;

  axis->nnp = axis->periodic ? n-C : n+1;

  PetscFunctionReturn(0);
}

#undef  __FUNCT__
#define __FUNCT__ "IGAAxisSetUp"
PetscErrorCode IGAAxisSetUp(IGAAxis axis)
{
  PetscInt p,m,n;
  const PetscReal *U;
  PetscErrorCode ierr;
  PetscFunctionBegin;
  PetscValidPointer(axis,1);
  if (axis->p < 1)
    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,
            "Must call IGAAxisSetDegree() first");
  if (axis->m < 2*axis->p+1)
    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,
            "Must call IGAAxisSetKnots() first");

  p = axis->p;
  m = axis->m;
  n = m - p - 1;
  U = axis->U;

  ierr = IGAAxisGetSpans(axis,&axis->nel,&axis->span);CHKERRQ(ierr);

  if (axis->periodic) {
    PetscInt s = 1;
    while(s < p && U[m-p] == U[m-p+s]) s++;
    axis->nnp = n-p+s;
  } else {
    axis->nnp = n+1;
  }

  PetscFunctionReturn(0);
}

PetscInt IGA_SpanCount(PetscInt n,PetscInt p,const PetscReal U[])
{
  PetscInt i, span = 0;
  for (i=p; i<=n; i++)
    if (U[i] != U[i+1])
      span++;
  return span;
}

PetscInt IGA_SpanIndex(PetscInt n,PetscInt p,const PetscReal U[],PetscInt index[])
{
  PetscInt i, span = 0;
  for (i=p; i<=n; i++)
    if (U[i] != U[i+1])
      index[span++] = i;
  return span;
}

/*
void IGA_Stencil(PetscInt n,PetscInt p,const PetscReal U[],
                 PetscBool periodic,PetscInt first[],PetscInt last[])
{
  PetscInt i,k;
  for (i=0; i<=n; i++) {
    k = i;
    while (U[k]==U[k+1]) k++;
    first[i] = k - p;
    k = i + p + 1;
    while (U[k]==U[k-1]) k--;
    last[i] = k-1;
  }
  if (!periodic)
    for (i=0; i<=p; i++) {
      first[i]  = 0;
      last[n-i] = n;
    }
  else {
    PetscInt s = 1;
    while (s < p && U[m-p]==U[m-p+s]) s++;
    k = n - p + s;
    while (U[k]==U[k+1]) k++;
    first[0] = k - s - n;
  }
}
*/
Tip: Filter by directory path e.g. /media app.js to search for public/media/app.js.
Tip: Use camelCasing e.g. ProjME to search for ProjectModifiedEvent.java.
Tip: Filter by extension type e.g. /repo .js to search for all .js files in the /repo directory.
Tip: Separate your search with spaces e.g. /ssh pom.xml to search for src/ssh/pom.xml.
Tip: Use ↑ and ↓ arrow keys to navigate and return to view the file.
Tip: You can also navigate files with Ctrl+j (next) and Ctrl+k (previous) and view the file with Ctrl+o.
Tip: You can also navigate files with Alt+j (next) and Alt+k (previous) and view the file with Alt+o.