[petsc-users] Preconditioning systems of equations with complex numbers
Mark Adams
mfadams at lbl.gov
Fri Feb 1 08:12:52 CST 2019
>
> Both GAMG and ILU are nice and dandy for this,
>
I would test Richardson/SOR and Chebyshev/Jacobi on the tiny system and
converge it way down, say rtol = 1.e-12. See which one is better in the
early iteration and pick it. It would be nice to check that it solves the
problem ...
The residual drops about 5 orders in the first iteration and then
flatlines. That is very bad. Check that the smoother can actually solve the
problem.
> but as soon as I look at a bigger system, like a network with 8500 buses,
> the out-of-box gamg craps out. I am not sure where to start when it comes
> to tuning GAMG.
>
First try -pc_gamg_nsmooths 0
You can run with -info and grep on GAMG to see some diagnostic output.
Eigen estimates are fragile in practice and with this parameter and SOR
there are no eigen estimates needed. The max eigen values for all levels
should be between say 2 (or less) and say 3-4. Much higher is a sign of a
problem.
> Attached is the ksp monitor/view output for gamg on the unsuccessful solve
>
> I'm also attaching a zip file which contains the simple PETSc script that
> loads the binary matrix/vector as well as two test cases, if you guys are
> interested in trying it out. It only works if you have PETSc configured
> with complex numbers.
>
> Thanks
>
> Justin
>
> PS - A couple years ago I had asked if there was a paper/tutorial on
> using/tuning GAMG. Does such a thing exist today?
>
There is a write up in the manual that is tutorial like.
>
> On Thu, Jan 31, 2019 at 5:00 PM Matthew Knepley <knepley at gmail.com> wrote:
>
>> On Thu, Jan 31, 2019 at 6:22 PM Justin Chang <jychang48 at gmail.com> wrote:
>>
>>> Here's IMHO the simplest explanation of the equations I'm trying to
>>> solve:
>>>
>>> http://home.eng.iastate.edu/~jdm/ee458_2011/PowerFlowEquations.pdf
>>>
>>> Right now we're just trying to solve eq(5) (in section 1), inverting the
>>> linear Y-bus matrix. Eventually we have to be able to solve equations like
>>> those in the next section.
>>>
>>
>> Maybe I am reading this wrong, but the Y-bus matrix looks like an
>> M-matrix to me (if all the y's are positive). This means
>> that it should be really easy to solve, and I think GAMG should do it.
>> You can start out just doing relaxation, like SOR, on
>> small examples.
>>
>> Thanks,
>>
>> Matt
>>
>>
>>> On Thu, Jan 31, 2019 at 1:47 PM Matthew Knepley <knepley at gmail.com>
>>> wrote:
>>>
>>>> On Thu, Jan 31, 2019 at 3:20 PM Justin Chang via petsc-users <
>>>> petsc-users at mcs.anl.gov> wrote:
>>>>
>>>>> Hi all,
>>>>>
>>>>> I'm working with some folks to extract a linear system of equations
>>>>> from an external software package that solves power flow equations in
>>>>> complex form. Since that external package uses serial direct solvers like
>>>>> KLU from suitesparse, I want a proof-of-concept where the same matrix can
>>>>> be solved in PETSc using its parallel solvers.
>>>>>
>>>>> I got mumps to achieve a very minor speedup across two MPI processes
>>>>> on a single node (went from solving a 300k dog system in 1.8 seconds to 1.5
>>>>> seconds). However I want to use iterative solvers and preconditioners but I
>>>>> have never worked with complex numbers so I am not sure what the "best"
>>>>> options are given PETSc's capabilities.
>>>>>
>>>>> So far I tried GMRES/BJACOBI and it craps out (unsurprisingly). I
>>>>> believe I also tried BICG with BJACOBI and while it did converge it
>>>>> converged slowly. Does anyone have recommendations on how one would go
>>>>> about preconditioning PETSc matrices with complex numbers? I was originally
>>>>> thinking about converting it to cartesian form: Declaring all voltages =
>>>>> sqrt(real^2+imaginary^2) and all angles to be something like a conditional
>>>>> arctan(imaginary/real) because all the papers I've seen in literature that
>>>>> claim to successfully precondition power flow equations operate in this
>>>>> form.
>>>>>
>>>>
>>>> 1) We really need to see the (simplified) equations
>>>>
>>>> 2) All complex equations can be converted to a system of real equations
>>>> twice as large, but this is not necessarily the best way to go
>>>>
>>>> Thanks,
>>>>
>>>> Matt
>>>>
>>>>
>>>>> Justin
>>>>>
>>>>
>>>>
>>>> --
>>>> 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/>
>>>>
>>>
>>
>> --
>> 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/>
>>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20190201/0febddff/attachment.html>
More information about the petsc-users
mailing list