xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 84df9cb40eca90ea9b18a456fab7a4ecc7f6c1a4)
1 
2 #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESReset_NRichardson"
6 PetscErrorCode SNESReset_NRichardson(SNES snes)
7 {
8   PetscErrorCode ierr;
9 
10   PetscFunctionBegin;
11   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
12   PetscFunctionReturn(0);
13 }
14 
15 /*
16   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
17 
18   Input Parameter:
19 . snes - the SNES context
20 
21   Application Interface Routine: SNESDestroy()
22 */
23 #undef __FUNCT__
24 #define __FUNCT__ "SNESDestroy_NRichardson"
25 PetscErrorCode SNESDestroy_NRichardson(SNES snes)
26 {
27   PetscErrorCode   ierr;
28   SNES_NRichardson *neP = (SNES_NRichardson*)snes->data;
29 
30   PetscFunctionBegin;
31   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
32   ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr);
33   ierr = PetscFree(snes->data);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 /*
38    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
39    of the SNESNRICHARDSON nonlinear solver.
40 
41    Input Parameters:
42 +  snes - the SNES context
43 -  x - the solution vector
44 
45    Application Interface Routine: SNESSetUp()
46  */
47 #undef __FUNCT__
48 #define __FUNCT__ "SNESSetUp_NRichardson"
49 PetscErrorCode SNESSetUp_NRichardson(SNES snes)
50 {
51   PetscErrorCode ierr;
52 
53   PetscFunctionBegin;
54   ierr = SNESDefaultGetWork(snes,1);CHKERRQ(ierr);
55   PetscFunctionReturn(0);
56 }
57 
58 EXTERN_C_BEGIN
59 #undef __FUNCT__
60 #define __FUNCT__ "SNESLineSearchSetMonitor_NRichardson"
61 PetscErrorCode  SNESLineSearchSetMonitor_NRichardson(SNES snes,PetscBool  flg)
62 {
63   SNES_NRichardson *neP = (SNES_NRichardson*)snes->data;
64   PetscErrorCode   ierr;
65 
66   PetscFunctionBegin;
67   if (flg && !neP->monitor) {
68     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&neP->monitor);CHKERRQ(ierr);
69   } else if (!flg && neP->monitor) {
70     ierr = PetscViewerDestroy(&neP->monitor);CHKERRQ(ierr);
71   }
72   PetscFunctionReturn(0);
73 }
74 EXTERN_C_END
75 
76 PetscErrorCode SNESLineSearchNoNorms_NRichardson(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
77 PetscErrorCode SNESLineSearchNo_NRichardson(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
78 PetscErrorCode SNESLineSearchQuadratic_NRichardson(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool*);
79 /*
80   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS method.
81 
82   Input Parameter:
83 . snes - the SNES context
84 
85   Application Interface Routine: SNESSetFromOptions()
86 */
87 #undef __FUNCT__
88 #define __FUNCT__ "SNESSetFromOptions_NRichardson"
89 static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes)
90 {
91   SNES_NRichardson   *ls = (SNES_NRichardson *)snes->data;
92   SNESLineSearchType indx;
93   PetscBool          flg,set;
94   PetscErrorCode     ierr;
95 
96   PetscFunctionBegin;
97   ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr);
98     ierr = PetscOptionsEnum("-snes_ls","Richardson line search type","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)ls->type,(PetscEnum*)&indx,&flg);CHKERRQ(ierr);
99     if (flg) {
100       ls->type = indx;
101       switch (indx) {
102       case SNES_LS_BASIC:
103         ierr = SNESLineSearchSet(snes,SNESLineSearchNo_NRichardson,PETSC_NULL);CHKERRQ(ierr);
104         break;
105       case SNES_LS_BASIC_NONORMS:
106         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms_NRichardson,PETSC_NULL);CHKERRQ(ierr);
107         break;
108       case SNES_LS_QUADRATIC:
109         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NRichardson,PETSC_NULL);CHKERRQ(ierr);
110         break;
111       case SNES_LS_CUBIC:
112         SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_SUP,"No support for cubic search");
113         break;
114       }
115     }
116     ierr = PetscOptionsReal("-snes_ls_damping","Damping parameter","SNES",ls->damping,&ls->damping,&flg);CHKERRQ(ierr);
117     ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
118     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
119   ierr = PetscOptionsTail();CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 /*
124   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
125 
126   Input Parameters:
127 + SNES - the SNES context
128 - viewer - visualization context
129 
130   Application Interface Routine: SNESView()
131 */
132 #undef __FUNCT__
133 #define __FUNCT__ "SNESView_NRichardson"
134 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
135 {
136   SNES_NRichardson *ls = (SNES_NRichardson *)snes->data;
137   PetscBool        iascii;
138   PetscErrorCode   ierr;
139 
140   PetscFunctionBegin;
141   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
142   if (iascii) {
143     ierr = PetscViewerASCIIPrintf(viewer,"  richardson variant: %s\n", SNESLineSearchTypeName(ls->type));CHKERRQ(ierr);
144   }
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "SNESLineSearchNo_NRichardson"
150 PetscErrorCode SNESLineSearchNo_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
151 {
152   PetscErrorCode   ierr;
153   SNES_NRichardson *neP = (SNES_NRichardson *) snes->data;
154 
155   PetscFunctionBegin;
156   ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr);
157   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
158   ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
159   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm");
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "SNESLineSearchNoNorms_NRichardson"
165 PetscErrorCode SNESLineSearchNoNorms_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
166 {
167   PetscErrorCode   ierr;
168   SNES_NRichardson *neP = (SNES_NRichardson *) snes->data;
169 
170   PetscFunctionBegin;
171   ierr = VecAXPY(X, neP->damping, Y);CHKERRQ(ierr);
172   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "SNESLineSearchQuadratic_NRichardson"
178 PetscErrorCode SNESLineSearchQuadratic_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
179 {
180   PetscInt         i;
181   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
182   PetscReal        norms[3];
183   PetscReal        alpha,a,b;
184   PetscErrorCode   ierr;
185   SNES_NRichardson *neP = (SNES_NRichardson *) snes->data;
186 
187   PetscFunctionBegin;
188   norms[0]  = fnorm;
189   for(i=1; i < 3; ++i) {
190     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
191     ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr);
192     ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr);
193   }
194   for(i = 0; i < 3; ++i) {
195     norms[i] = PetscSqr(norms[i]);
196   }
197   /* Fit a quadratic:
198        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
199        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
200        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
201        c = y_0
202        x_min = -b/2a
203 
204        If we let x_{0,1,2} = 0, 0.5, 1.0
205        a = 2 y_2 - 4 y_1 + 2 y_0
206        b =  -y_2 + 4 y_1 - 3 y_0
207        c =   y_0
208   */
209   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
210   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
211   /* Check for positive a (concave up) */
212   if (a >= 0.0) {
213     alpha = -b/(2.0*a);
214     alpha = PetscMin(alpha, alphas[2]);
215     alpha = PetscMax(alpha, alphas[0]);
216   } else {
217     alpha = 1.0;
218   }
219   if (neP->monitor) {
220     ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
221     ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr);
222     ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
223   }
224   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
225   if (alpha != 1.0) {
226     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
227     ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
228   } else {
229     *gnorm = PetscSqrtReal(norms[2]);
230   }
231   if (alpha == 0.0) *flag = PETSC_FALSE;
232   else              *flag = PETSC_TRUE;
233   PetscFunctionReturn(0);
234 }
235 
236 /*
237   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
238 
239   Input Parameters:
240 . snes - the SNES context
241 
242   Output Parameter:
243 . outits - number of iterations until termination
244 
245   Application Interface Routine: SNESSolve()
246 */
247 #undef __FUNCT__
248 #define __FUNCT__ "SNESSolve_NRichardson"
249 PetscErrorCode SNESSolve_NRichardson(SNES snes)
250 {
251   SNES_NRichardson    *neP = (SNES_NRichardson *) snes->data;
252   Vec                 X, Y, F, W;
253   PetscReal           fnorm;
254   PetscInt            maxits, i;
255   PetscErrorCode      ierr;
256   SNESConvergedReason reason;
257 
258   PetscFunctionBegin;
259   snes->reason = SNES_CONVERGED_ITERATING;
260 
261   maxits = snes->max_its;	     /* maximum number of iterations */
262   X      = snes->vec_sol;	     /* X^n */
263   Y      = snes->vec_sol_update; /* \tilde X */
264   F      = snes->vec_func;       /* residual vector */
265   W      = snes->work[0];        /* work vector */
266 
267   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
268   snes->iter = 0;
269   snes->norm = 0.;
270   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
271   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
272   if (snes->domainerror) {
273     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
274     PetscFunctionReturn(0);
275   }
276   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
277   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
278   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
279   snes->norm = fnorm;
280   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
281   SNESLogConvHistory(snes,fnorm,0);
282   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
283 
284   /* set parameter for default relative tolerance convergence test */
285   snes->ttol = fnorm*snes->rtol;
286   /* test convergence */
287   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
288   if (snes->reason) PetscFunctionReturn(0);
289 
290   for(i = 0; i < maxits; i++) {
291     PetscBool  lsSuccess = PETSC_TRUE;
292     PetscReal  dummyNorm;
293 
294     /* Call general purpose update function */
295     if (snes->ops->update) {
296       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
297     }
298     if (!snes->pc) {
299       ierr = VecCopy(F,Y);CHKERRQ(ierr);
300       ierr = VecScale(Y,-1.0);CHKERRQ(ierr);
301     } else {
302       ierr = VecCopy(X,Y);CHKERRQ(ierr);
303       ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr);
304       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
305       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
306         snes->reason = SNES_DIVERGED_INNER;
307         PetscFunctionReturn(0);
308       }
309       ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
310     }
311 
312       ierr = (*neP->LineSearch)(snes, neP->lsP, X, F, Y, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
313     if (!lsSuccess) {
314       if (++snes->numFailures >= snes->maxFailures) {
315         snes->reason = SNES_DIVERGED_LINE_SEARCH;
316         break;
317       }
318     }
319     if (snes->nfuncs >= snes->max_funcs) {
320       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
321       break;
322     }
323     if (snes->domainerror) {
324       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
325       PetscFunctionReturn(0);
326     }
327     /* Monitor convergence */
328     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
329     snes->iter = i+1;
330     snes->norm = fnorm;
331     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
332     SNESLogConvHistory(snes,snes->norm,0);
333     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
334     /* Test for convergence */
335     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
336     if (snes->reason) break;
337   }
338   if (i == maxits) {
339     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
340     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
341   }
342   PetscFunctionReturn(0);
343 }
344 
345 typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,void*,PetscBool *);                 /* force argument to next function to not be extern C*/
346 EXTERN_C_BEGIN
347 #undef __FUNCT__
348 #define __FUNCT__ "SNESLineSearchSetPreCheck_NRichardson"
349 PetscErrorCode  SNESLineSearchSetPreCheck_NRichardson(SNES snes, FCN1 func, void *checkctx)
350 {
351   PetscFunctionBegin;
352   ((SNES_NRichardson *)(snes->data))->precheckstep = func;
353   ((SNES_NRichardson *)(snes->data))->precheck     = checkctx;
354   PetscFunctionReturn(0);
355 }
356 EXTERN_C_END
357 
358 typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,PetscReal,PetscReal,Vec,Vec,PetscReal*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
359 EXTERN_C_BEGIN
360 #undef __FUNCT__
361 #define __FUNCT__ "SNESLineSearchSet_NRichardson"
362 PetscErrorCode  SNESLineSearchSet_NRichardson(SNES snes, FCN2 func, void *lsctx)
363 {
364   PetscFunctionBegin;
365   ((SNES_NRichardson *)(snes->data))->LineSearch = func;
366   ((SNES_NRichardson *)(snes->data))->lsP        = lsctx;
367   PetscFunctionReturn(0);
368 }
369 EXTERN_C_END
370 
371 typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/
372 EXTERN_C_BEGIN
373 #undef __FUNCT__
374 #define __FUNCT__ "SNESLineSearchSetPostCheck_NRichardson"
375 PetscErrorCode  SNESLineSearchSetPostCheck_NRichardson(SNES snes, FCN3 func, void *checkctx)
376 {
377   PetscFunctionBegin;
378   ((SNES_NRichardson *)(snes->data))->postcheckstep = func;
379   ((SNES_NRichardson *)(snes->data))->postcheck     = checkctx;
380   PetscFunctionReturn(0);
381 }
382 EXTERN_C_END
383 
384 /*MC
385   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
386 
387   Level: beginner
388 
389   Options Database:
390 +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
391 -   -snes_ls <basic,basicnormnorms,quadratic>
392 
393   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda (F(x^n) - b)
394             where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.
395          If an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner solver is called an initial solution x^n and the
396             nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
397 
398      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
399 
400 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
401 M*/
402 EXTERN_C_BEGIN
403 #undef __FUNCT__
404 #define __FUNCT__ "SNESCreate_NRichardson"
405 PetscErrorCode  SNESCreate_NRichardson(SNES snes)
406 {
407   SNES_NRichardson *neP;
408   PetscErrorCode   ierr;
409 
410   PetscFunctionBegin;
411   snes->ops->destroy	     = SNESDestroy_NRichardson;
412   snes->ops->setup	     = SNESSetUp_NRichardson;
413   snes->ops->setfromoptions  = SNESSetFromOptions_NRichardson;
414   snes->ops->view            = SNESView_NRichardson;
415   snes->ops->solve	     = SNESSolve_NRichardson;
416   snes->ops->reset           = SNESReset_NRichardson;
417 
418   snes->usesksp              = PETSC_FALSE;
419   snes->usespc               = PETSC_TRUE;
420 
421   ierr = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr);
422   snes->data = (void*) neP;
423   neP->maxstep	     = 1.e8;
424   neP->steptol       = 1.e-12;
425   neP->type          = SNES_LS_QUADRATIC;
426   neP->damping	     = 1.0;
427   neP->LineSearch    = SNESLineSearchQuadratic_NRichardson;
428   neP->lsP           = PETSC_NULL;
429   neP->postcheckstep = PETSC_NULL;
430   neP->postcheck     = PETSC_NULL;
431   neP->precheckstep  = PETSC_NULL;
432   neP->precheck      = PETSC_NULL;
433 
434   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_NRichardson",SNESLineSearchSetMonitor_NRichardson);CHKERRQ(ierr);
435   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_NRichardson",SNESLineSearchSet_NRichardson);CHKERRQ(ierr);
436   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_NRichardson",SNESLineSearchSetPostCheck_NRichardson);CHKERRQ(ierr);
437   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_NRichardson",SNESLineSearchSetPreCheck_NRichardson);CHKERRQ(ierr);
438   PetscFunctionReturn(0);
439 }
440 EXTERN_C_END
441