[petsc-users] Norm_2 calculation
Fazlul Huq
huq2090 at gmail.com
Fri Jun 5 13:17:51 CDT 2020
The same is happening with "-pc_type hypre -pc_hypre_type boomerang" and
with "-pc_type ilu".
I tried with,
ierr =
KSPSetTolerances(ksp,1.e-10,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
this also.
Getting higher values of error even then.
Any suggestion to go over this issue?
Shall I send the code? It's not a large script (single routine, less than
100 lines except comment).
Thank you.
Sincerely,
Huq
On Fri, Jun 5, 2020 at 12:53 PM Matthew Knepley <knepley at gmail.com> wrote:
> On Fri, Jun 5, 2020 at 12:32 PM Fazlul Huq <huq2090 at gmail.com> wrote:
>
>> Input string:
>> mpiexec -np 1 ./poisson_m -n 10000000 -pc_type cholesky -ksp_view
>> -ksp_converged_reason -ksp_monitor_true_residual
>>
>> Output is attached herewith.
>>
>
> Solving the problem...
>
> 0 KSP preconditioned resid norm 9.741453255800e+07 true resid norm
> 1.802775638200e+01 ||r(i)||/||b|| 1.000000000000e+00
>
> 1 KSP preconditioned resid norm 5.796040438920e+00 true resid norm
> 3.146946078273e-05 ||r(i)||/||b|| 1.745611606675e-06
>
> Your Cholesky preconditioner is crap, probably because the system is
> either not symmetric or nearly singular. You can see this
> because the preconditioned residual is 6 orders of magnitude greater than
> the true residual. Also, if you want to evaluate convergence, you
> should probably use a lower tolerance like
>
> -ksp_rtol 1e-10
>
> Thanks,
>
> Matt
>
> Thank you.
>> Sincerely,
>> Huq
>>
>> On Fri, Jun 5, 2020 at 10:16 AM Matthew Knepley <knepley at gmail.com>
>> wrote:
>>
>>> On Fri, Jun 5, 2020 at 11:04 AM Fazlul Huq <huq2090 at gmail.com> wrote:
>>>
>>>> Hello All,
>>>> I'm trying to calculate the norm_2 error of a solution.
>>>> Here, vector s: Analytical solution
>>>> vector x: Numerical solution
>>>>
>>>> ierr = VecView(s,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
>>>> ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
>>>>
>>>> ierr = VecAXPY(x,-1.0,s);CHKERRQ(ierr);
>>>> ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
>>>> ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
>>>> ierr = KSPGetIterationNumber(ksp,&its);CHKERRQ(ierr);
>>>> if (norm > tol) {
>>>> ierr = PetscPrintf(PETSC_COMM_WORLD,"Second Norm of error %g\n",
>>>> double)norm);CHKERRQ(ierr);
>>>> ierr = PetscPrintf(PETSC_COMM_WORLD,"Iterations
>>>> %D\n",its);CHKERRQ(ierr);
>>>> }
>>>>
>>>> Am I calculating the "Norm_2" error correctly or making any mistake?
>>>> Unfortunately, for large sized matrix, say 10^6, I am getting very high
>>>> value of "Norm_2" error.
>>>>
>>>
>>> 1) I am guessing x comes from a KSPSolve(). It is only as accurate as
>>> your tolerance
>>>
>>> 2) This is the l_2 norm, not the L_2 norm, so if you are using a
>>> continuum method like FEM, this is likely wrong.
>>>
>>> Thanks,
>>>
>>> Matt
>>>
>>>
>>>> Thank you.
>>>> Sincerely,
>>>> Huq
>>>>
>>>> --
>>>>
>>>> Fazlul Huq
>>>> Graduate Research Assistant
>>>> Department of Nuclear, Plasma & Radiological Engineering (NPRE)
>>>> University of Illinois at Urbana-Champaign (UIUC)
>>>> E-mail: huq2090 at gmail.com
>>>>
>>>
>>>
>>> --
>>> What most experimenters take for granted before they begin their
>>> experiments is infinitely more interesting than any results to which their
>>> experiments lead.
>>> -- Norbert Wiener
>>>
>>> https://www.cse.buffalo.edu/~knepley/
>>> <http://www.cse.buffalo.edu/~knepley/>
>>>
>>
>>
>> --
>>
>> Fazlul Huq
>> Graduate Research Assistant
>> Department of Nuclear, Plasma & Radiological Engineering (NPRE)
>> University of Illinois at Urbana-Champaign (UIUC)
>> E-mail: huq2090 at gmail.com
>>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
> https://www.cse.buffalo.edu/~knepley/
> <http://www.cse.buffalo.edu/~knepley/>
>
--
Fazlul Huq
Graduate Research Assistant
Department of Nuclear, Plasma & Radiological Engineering (NPRE)
University of Illinois at Urbana-Champaign (UIUC)
E-mail: huq2090 at gmail.com
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200605/008157b0/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: output_boomeramg.log
Type: text/x-log
Size: 3983 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200605/008157b0/attachment-0002.bin>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: output_ilu.log
Type: text/x-log
Size: 3088 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200605/008157b0/attachment-0003.bin>
More information about the petsc-users
mailing list