Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
519 dsic.upv.es!antodo 1
/*
2
     The basic EPS routines, Create, View, etc. are here.
3
*/
4
#include "src/eps/epsimpl.h"      /*I "slepceps.h" I*/
5
 
6
PetscTruth EPSRegisterAllCalled = PETSC_FALSE;
7
PetscFList EPSList = 0;
8
 
9
#undef __FUNCT__  
10
#define __FUNCT__ "EPSView"
11
/*@C
12
   EPSView - Prints the EPS data structure.
13
 
14
   Collective on EPS
15
 
16
   Input Parameters:
17
+  eps - the eigenproblem solver context
18
-  viewer - optional visualization context
19
 
20
   Options Database Key:
21
.  -eps_view -  Calls EPSView() at end of EPSSolve()
22
 
23
   Note:
24
   The available visualization contexts include
25
+     PETSC_VIEWER_STDOUT_SELF - standard output (default)
26
-     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
27
         output where only the first processor opens
28
         the file.  All other processors send their
29
         data to the first processor to print.
30
 
31
   The user can open an alternative visualization context with
32
   PetscViewerASCIIOpen() - output to a specified file.
33
 
34
   Level: beginner
35
 
36
.seealso: STView(), PetscViewerASCIIOpen()
37
@*/
38
PetscErrorCode EPSView(EPS eps,PetscViewer viewer)
39
{
40
  PetscErrorCode ierr;
41
  char        *type, *which;
42
  PetscTruth  isascii;
43
 
44
  PetscFunctionBegin;
45
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
46
  if (!viewer) viewer = PETSC_VIEWER_STDOUT_(eps->comm);
47
  PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
48
  PetscCheckSameComm(eps,1,viewer,2);
49
 
50
#if defined(PETSC_USE_COMPLEX)
51
#define HERM "hermitian"
52
#else
53
#define HERM "symmetric"
54
#endif
55
  ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
56
  if (isascii) {
57
    ierr = PetscViewerASCIIPrintf(viewer,"EPS Object:\n");CHKERRQ(ierr);
58
    switch (eps->problem_type) {
59
      case EPS_HEP:   type = HERM " eigenvalue problem"; break;
60
      case EPS_GHEP:  type = "generalized " HERM " eigenvalue problem"; break;
61
      case EPS_NHEP:  type = "non-" HERM " eigenvalue problem"; break;
62
      case EPS_GNHEP: type = "generalized non-" HERM " eigenvalue problem"; break;
63
      case 0:         type = "not yet set"; break;
64
      default: SETERRQ(1,"Wrong value of eps->problem_type");
65
    }
66
    ierr = PetscViewerASCIIPrintf(viewer,"  problem type: %s\n",type);CHKERRQ(ierr);
67
    ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
68
    if (type) {
780 dsic.upv.es!jroman 69
      ierr = PetscViewerASCIIPrintf(viewer,"  method: %s",type);CHKERRQ(ierr);
70
      switch (eps->solverclass) {
71
        case EPS_ONE_SIDE:
72
          ierr = PetscViewerASCIIPrintf(viewer,"\n",type);CHKERRQ(ierr); break;
73
        case EPS_TWO_SIDE:
74
          ierr = PetscViewerASCIIPrintf(viewer," (two-sided)\n",type);CHKERRQ(ierr); break;
75
        default: SETERRQ(1,"Wrong value of eps->solverclass");
76
      }
519 dsic.upv.es!antodo 77
    } else {
78
      ierr = PetscViewerASCIIPrintf(viewer,"  method: not yet set\n");CHKERRQ(ierr);
79
    }
80
    if (eps->ops->view) {
81
      ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
82
      ierr = (*eps->ops->view)(eps,viewer);CHKERRQ(ierr);
83
      ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
84
    }
85
    switch (eps->which) {
86
      case EPS_LARGEST_MAGNITUDE:  which = "largest eigenvalues in magnitude"; break;
87
      case EPS_SMALLEST_MAGNITUDE: which = "smallest eigenvalues in magnitude"; break;
88
      case EPS_LARGEST_REAL:       which = "largest real parts"; break;
89
      case EPS_SMALLEST_REAL:      which = "smallest real parts"; break;
90
      case EPS_LARGEST_IMAGINARY:  which = "largest imaginary parts"; break;
91
      case EPS_SMALLEST_IMAGINARY: which = "smallest imaginary parts"; break;
92
      default: SETERRQ(1,"Wrong value of eps->which");
93
    }
94
    ierr = PetscViewerASCIIPrintf(viewer,"  selected portion of the spectrum: %s\n",which);CHKERRQ(ierr);
95
    ierr = PetscViewerASCIIPrintf(viewer,"  number of eigenvalues (nev): %d\n",eps->nev);CHKERRQ(ierr);
600 dsic.upv.es!jroman 96
    ierr = PetscViewerASCIIPrintf(viewer,"  number of column vectors (ncv): %d\n",eps->ncv);CHKERRQ(ierr);
519 dsic.upv.es!antodo 97
    ierr = PetscViewerASCIIPrintf(viewer,"  maximum number of iterations: %d\n", eps->max_it);
98
    ierr = PetscViewerASCIIPrintf(viewer,"  tolerance: %g\n",eps->tol);CHKERRQ(ierr);
99
    ierr = PetscViewerASCIIPrintf(viewer,"  orthogonalization method: ");CHKERRQ(ierr);
100
    switch (eps->orthog_type) {
101
      case EPS_MGS_ORTH:
102
        ierr = PetscViewerASCIIPrintf(viewer,"modified Gram-Schmidt\n");CHKERRQ(ierr);
103
        break;
104
      case EPS_CGS_ORTH:
105
        ierr = PetscViewerASCIIPrintf(viewer,"classical Gram-Schmidt\n");CHKERRQ(ierr);
106
        break;
107
      default: SETERRQ(1,"Wrong value of eps->orth_type");
108
    }
109
    ierr = PetscViewerASCIIPrintf(viewer,"  orthogonalization refinement: ");CHKERRQ(ierr);
110
    switch (eps->orthog_ref) {
111
      case EPS_ORTH_REFINE_NEVER:
112
        ierr = PetscViewerASCIIPrintf(viewer,"never\n");CHKERRQ(ierr);
113
        break;
114
      case EPS_ORTH_REFINE_IFNEEDED:
115
        ierr = PetscViewerASCIIPrintf(viewer,"if needed (eta: %f)\n",eps->orthog_eta);CHKERRQ(ierr);
116
        break;
117
      case EPS_ORTH_REFINE_ALWAYS:
118
        ierr = PetscViewerASCIIPrintf(viewer,"always\n");CHKERRQ(ierr);
119
        break;
120
      default: SETERRQ(1,"Wrong value of eps->orth_ref");
121
    }
122
    ierr = PetscViewerASCIIPrintf(viewer,"  dimension of user-provided deflation space: %d\n",eps->nds);CHKERRQ(ierr);
123
    ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
124
    ierr = STView(eps->OP,viewer); CHKERRQ(ierr);
125
    ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
126
  } else {
127
    if (eps->ops->view) {
128
      ierr = (*eps->ops->view)(eps,viewer);CHKERRQ(ierr);
129
    }
130
    ierr = STView(eps->OP,viewer); CHKERRQ(ierr);
131
  }
132
  PetscFunctionReturn(0);
133
}
134
 
135
#undef __FUNCT__  
136
#define __FUNCT__ "EPSPublish_Petsc"
137
static PetscErrorCode EPSPublish_Petsc(PetscObject object)
138
{
139
  PetscFunctionBegin;
140
  PetscFunctionReturn(0);
141
}
142
 
143
#undef __FUNCT__  
144
#define __FUNCT__ "EPSCreate"
145
/*@C
146
   EPSCreate - Creates the default EPS context.
147
 
148
   Collective on MPI_Comm
149
 
150
   Input Parameter:
151
.  comm - MPI communicator
152
 
153
   Output Parameter:
154
.  eps - location to put the EPS context
155
 
156
   Note:
157
   The default EPS type is EPSARNOLDI
158
 
159
   Level: beginner
160
 
161
.seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
162
@*/
163
PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
164
{
165
  PetscErrorCode ierr;
166
  EPS            eps;
167
 
168
  PetscFunctionBegin;
169
  PetscValidPointer(outeps,2);
170
  *outeps = 0;
171
#ifndef PETSC_USE_DYNAMIC_LIBRARIES
172
  ierr = EPSRegisterAll(PETSC_NULL);CHKERRQ(ierr);
173
#endif
174
 
175
  PetscHeaderCreate(eps,_p_EPS,struct _EPSOps,EPS_COOKIE,-1,"EPS",comm,EPSDestroy,EPSView);
176
  PetscLogObjectCreate(eps);
177
  *outeps = eps;
178
 
179
  eps->bops->publish   = EPSPublish_Petsc;
180
  ierr = PetscMemzero(eps->ops,sizeof(struct _EPSOps));CHKERRQ(ierr);
181
 
182
  eps->type            = -1;
183
  eps->max_it          = 0;
184
  eps->nev             = 1;
185
  eps->ncv             = 0;
186
  eps->allocated_ncv   = 0;
187
  eps->nds             = 0;
188
  eps->tol             = 0.0;
189
  eps->which           = EPS_LARGEST_MAGNITUDE;
190
  eps->evecsavailable  = PETSC_FALSE;
191
  eps->problem_type    = (EPSProblemType)0;
780 dsic.upv.es!jroman 192
  eps->solverclass     = (EPSClass)0;
519 dsic.upv.es!antodo 193
 
194
  eps->vec_initial     = 0;
195
  eps->vec_initial_set = PETSC_FALSE;
780 dsic.upv.es!jroman 196
  eps->vec_initial_left= 0;
519 dsic.upv.es!antodo 197
  eps->V               = 0;
198
  eps->AV              = 0;
780 dsic.upv.es!jroman 199
  eps->W               = 0;
200
  eps->AW              = 0;
519 dsic.upv.es!antodo 201
  eps->T               = 0;
202
  eps->DS              = 0;
203
  eps->ds_ortho        = PETSC_TRUE;
204
  eps->eigr            = 0;
205
  eps->eigi            = 0;
206
  eps->errest          = 0;
780 dsic.upv.es!jroman 207
  eps->errest_left     = 0;
519 dsic.upv.es!antodo 208
  eps->OP              = 0;
209
  eps->data            = 0;
210
  eps->nconv           = 0;
211
  eps->its             = 0;
212
  eps->perm            = PETSC_NULL;
213
 
214
  eps->nwork           = 0;
215
  eps->work            = 0;
216
  eps->isgeneralized   = PETSC_FALSE;
217
  eps->ishermitian     = PETSC_FALSE;
218
  eps->setupcalled     = 0;
219
  eps->reason          = EPS_CONVERGED_ITERATING;
220
 
221
  eps->numbermonitors  = 0;
222
 
223
  eps->orthog_type     = EPS_CGS_ORTH;
224
  eps->orthog_ref      = EPS_ORTH_REFINE_IFNEEDED;
691 dsic.upv.es!antodo 225
  eps->orthog_eta      = PETSC_DEFAULT;
519 dsic.upv.es!antodo 226
 
227
  ierr = STCreate(comm,&eps->OP); CHKERRQ(ierr);
228
  PetscLogObjectParent(eps,eps->OP);
229
  ierr = PetscPublishAll(eps);CHKERRQ(ierr);
230
  PetscFunctionReturn(0);
231
}
232
 
233
#undef __FUNCT__  
234
#define __FUNCT__ "EPSSetType"
235
/*@C
236
   EPSSetType - Selects the particular solver to be used in the EPS object.
237
 
238
   Collective on EPS
239
 
240
   Input Parameters:
241
+  eps      - the eigensolver context
242
-  type     - a known method
243
 
244
   Options Database Key:
245
.  -eps_type <method> - Sets the method; use -help for a list
246
    of available methods
247
 
248
   Notes:  
545 dsic.upv.es!jroman 249
   See "slepc/include/slepceps.h" for available methods. The default
250
   is EPSARNOLDI.
519 dsic.upv.es!antodo 251
 
252
   Normally, it is best to use the EPSSetFromOptions() command and
253
   then set the EPS type from the options database rather than by using
254
   this routine.  Using the options database provides the user with
255
   maximum flexibility in evaluating the different available methods.
256
   The EPSSetType() routine is provided for those situations where it
257
   is necessary to set the iterative solver independently of the command
258
   line or options database.
259
 
260
   Level: intermediate
261
 
262
.seealso: STSetType(), EPSType
263
@*/
264
PetscErrorCode EPSSetType(EPS eps,EPSType type)
265
{
266
  PetscErrorCode ierr,(*r)(EPS);
267
  PetscTruth match;
268
 
269
  PetscFunctionBegin;
270
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
271
  PetscValidCharPointer(type,2);
272
 
273
  ierr = PetscTypeCompare((PetscObject)eps,type,&match);CHKERRQ(ierr);
274
  if (match) PetscFunctionReturn(0);
275
 
276
  if (eps->data) {
277
    /* destroy the old private EPS context */
278
    ierr = (*eps->ops->destroy)(eps); CHKERRQ(ierr);
279
    eps->data = 0;
280
  }
281
  /* Get the function pointers for the iterative method requested */
282
  if (!EPSRegisterAllCalled) {ierr = EPSRegisterAll(PETSC_NULL); CHKERRQ(ierr);}
283
 
284
  ierr = PetscFListFind(eps->comm,EPSList,type,(void (**)(void)) &r);CHKERRQ(ierr);
285
 
286
  if (!r) SETERRQ1(1,"Unknown EPS type given: %s",type);
287
 
288
  eps->setupcalled = 0;
289
  ierr = PetscMemzero(eps->ops,sizeof(struct _EPSOps));CHKERRQ(ierr);
290
  ierr = (*r)(eps); CHKERRQ(ierr);
291
 
292
  ierr = PetscObjectChangeTypeName((PetscObject)eps,type);CHKERRQ(ierr);
293
  PetscFunctionReturn(0);
294
}
295
 
296
#undef __FUNCT__  
297
#define __FUNCT__ "EPSRegisterDestroy"
298
/*@C
299
   EPSRegisterDestroy - Frees the list of EPS methods that were
300
   registered by EPSRegisterDynamic().
301
 
302
   Not Collective
303
 
304
   Level: advanced
305
 
306
.seealso: EPSRegisterDynamic(), EPSRegisterAll()
307
@*/
308
PetscErrorCode EPSRegisterDestroy(void)
309
{
310
  PetscErrorCode ierr;
311
 
312
  PetscFunctionBegin;
313
  if (EPSList) {
314
    ierr = PetscFListDestroy(&EPSList);CHKERRQ(ierr);
315
    EPSList = 0;
316
  }
317
  EPSRegisterAllCalled = PETSC_FALSE;
318
  PetscFunctionReturn(0);
319
}
320
 
321
#undef __FUNCT__  
322
#define __FUNCT__ "EPSGetType"
323
/*@C
324
   EPSGetType - Gets the EPS type as a string from the EPS object.
325
 
326
   Not Collective
327
 
328
   Input Parameter:
329
.  eps - the eigensolver context
330
 
331
   Output Parameter:
332
.  name - name of EPS method
333
 
334
   Level: intermediate
335
 
336
.seealso: EPSSetType()
337
@*/
338
PetscErrorCode EPSGetType(EPS eps,EPSType *type)
339
{
340
  PetscFunctionBegin;
341
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
342
  *type = eps->type_name;
343
  PetscFunctionReturn(0);
344
}
345
 
346
/*MC
347
   EPSRegisterDynamic - Adds a method to the eigenproblem solver package.
348
 
349
   Synopsis:
350
   EPSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(EPS))
351
 
352
   Not Collective
353
 
354
   Input Parameters:
355
+  name_solver - name of a new user-defined solver
356
.  path - path (either absolute or relative) the library containing this solver
357
.  name_create - name of routine to create the solver context
358
-  routine_create - routine to create the solver context
359
 
360
   Notes:
361
   EPSRegisterDynamic() may be called multiple times to add several user-defined solvers.
362
 
363
   If dynamic libraries are used, then the fourth input argument (routine_create)
364
   is ignored.
365
 
366
   Sample usage:
367
.vb
368
   EPSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
369
               "MySolverCreate",MySolverCreate);
370
.ve
371
 
372
   Then, your solver can be chosen with the procedural interface via
373
$     EPSSetType(eps,"my_solver")
374
   or at runtime via the option
375
$     -eps_type my_solver
376
 
377
   Level: advanced
378
 
379
   Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR}, ${BOPT},
380
   and others of the form ${any_environmental_variable} occuring in pathname will be
381
   replaced with appropriate values.
382
 
383
.seealso: EPSRegisterAll(), EPSRegisterDestroy()
384
 
385
M*/
386
 
387
#undef __FUNCT__  
388
#define __FUNCT__ "EPSRegister"
389
PetscErrorCode EPSRegister(char *sname,char *path,char *name,int (*function)(EPS))
390
{
391
  PetscErrorCode ierr;
392
  char           fullname[256];
393
 
394
  PetscFunctionBegin;
395
  ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
396
  ierr = PetscFListAdd(&EPSList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
397
  PetscFunctionReturn(0);
398
}
399
 
400
#undef __FUNCT__  
401
#define __FUNCT__ "EPSDestroy"
402
/*@
403
   EPSDestroy - Destroys the EPS context.
404
 
405
   Collective on EPS
406
 
407
   Input Parameter:
408
.  eps - eigensolver context obtained from EPSCreate()
409
 
410
   Level: beginner
411
 
412
.seealso: EPSCreate(), EPSSetUp(), EPSSolve()
413
@*/
414
PetscErrorCode EPSDestroy(EPS eps)
415
{
416
  PetscErrorCode ierr;
417
 
418
  PetscFunctionBegin;
419
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
420
  if (--eps->refct > 0) PetscFunctionReturn(0);
421
 
422
  /* if memory was published with AMS then destroy it */
423
  ierr = PetscObjectDepublish(eps);CHKERRQ(ierr);
424
 
425
  ierr = STDestroy(eps->OP);CHKERRQ(ierr);
426
 
427
  if (eps->ops->destroy) {
428
    ierr = (*eps->ops->destroy)(eps); CHKERRQ(ierr);
429
  }
430
 
431
  if (eps->T) {
432
    ierr = PetscFree(eps->T);CHKERRQ(ierr);
433
  }
434
 
435
  if (eps->perm) {
436
    ierr = PetscFree(eps->perm);CHKERRQ(ierr);
437
  }
438
 
439
  if (eps->vec_initial) {
440
    ierr = VecDestroy(eps->vec_initial);CHKERRQ(ierr);
441
  }
442
 
780 dsic.upv.es!jroman 443
  if (eps->vec_initial_left) {
444
    ierr = VecDestroy(eps->vec_initial_left);CHKERRQ(ierr);
445
  }
446
 
519 dsic.upv.es!antodo 447
  if (eps->nds > 0) {
448
    ierr = VecDestroyVecs(eps->DS, eps->nds);
449
  }
450
 
451
  if (eps->DSV) {
452
    ierr = PetscFree(eps->DSV);CHKERRQ(ierr);
453
  }
454
 
455
  PetscLogObjectDestroy(eps);
456
  PetscHeaderDestroy(eps);
457
  PetscFunctionReturn(0);
458
}
459
 
460
#undef __FUNCT__  
461
#define __FUNCT__ "EPSSetST"
462
/*@
463
   EPSSetST - Associates a spectral transformation object to the
464
   eigensolver.
465
 
466
   Collective on EPS
467
 
468
   Input Parameters:
469
+  eps - eigensolver context obtained from EPSCreate()
470
-  st   - the spectral transformation object
471
 
472
   Note:
473
   Use EPSGetST() to retrieve the spectral transformation context (for example,
474
   to free it at the end of the computations).
475
 
476
   Level: advanced
477
 
478
.seealso: EPSGetST()
479
@*/
480
PetscErrorCode EPSSetST(EPS eps,ST st)
481
{
482
  PetscErrorCode ierr;
483
 
484
  PetscFunctionBegin;
485
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
486
  PetscValidHeaderSpecific(st,ST_COOKIE,2);
487
  PetscCheckSameComm(eps,1,st,2);
488
  ierr = STDestroy(eps->OP); CHKERRQ(ierr);
489
  eps->OP = st;
490
  PetscObjectReference((PetscObject)eps->OP);
491
  PetscFunctionReturn(0);
492
}
493
 
494
#undef __FUNCT__  
495
#define __FUNCT__ "EPSGetST"
707 dsic.upv.es!antodo 496
/*@C
519 dsic.upv.es!antodo 497
   EPSGetST - Obtain the spectral transformation (ST) object associated
498
   to the eigensolver object.
499
 
500
   Not Collective
501
 
502
   Input Parameters:
503
.  eps - eigensolver context obtained from EPSCreate()
504
 
505
   Output Parameter:
506
.  st - spectral transformation context
507
 
508
   Level: beginner
509
 
510
.seealso: EPSSetST()
511
@*/
512
PetscErrorCode EPSGetST(EPS eps, ST *st)
513
{
514
  PetscFunctionBegin;
515
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
516
  *st = eps->OP;
517
  PetscFunctionReturn(0);
518
}
519
 
520
#undef __FUNCT__  
521
#define __FUNCT__ "EPSIsGeneralized"
522
/*@
523
   EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
524
   eigenvalue problem.
525
 
526
   Not collective
527
 
528
   Input Parameter:
529
.  eps - the eigenproblem solver context
530
 
531
   Output Parameter:
532
.  is - the answer
533
 
534
   Level: intermediate
535
 
536
@*/
537
PetscErrorCode EPSIsGeneralized(EPS eps,PetscTruth* is)
538
{
539
  PetscErrorCode ierr;
540
  Mat            B;
541
 
542
  PetscFunctionBegin;
543
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
544
  ierr = STGetOperators(eps->OP,PETSC_NULL,&B);CHKERRQ(ierr);
545
  if( B ) *is = PETSC_TRUE;
546
  else *is = PETSC_FALSE;
547
  if( eps->setupcalled ) {
548
    if( eps->isgeneralized != *is ) {
549
      SETERRQ(0,"Warning: Inconsistent EPS state");
550
    }
551
  }
552
  PetscFunctionReturn(0);
553
}
554
 
555
#undef __FUNCT__  
556
#define __FUNCT__ "EPSIsHermitian"
557
/*@
558
   EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
559
   eigenvalue problem.
560
 
561
   Not collective
562
 
563
   Input Parameter:
564
.  eps - the eigenproblem solver context
565
 
566
   Output Parameter:
567
.  is - the answer
568
 
569
   Level: intermediate
570
 
571
@*/
572
PetscErrorCode EPSIsHermitian(EPS eps,PetscTruth* is)
573
{
574
  PetscFunctionBegin;
575
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
576
  if( eps->ishermitian ) *is = PETSC_TRUE;
577
  else *is = PETSC_FALSE;
578
  PetscFunctionReturn(0);
579
}