[Trilinos-Users] [EXTERNAL] question about MueLu (AMG) setting

Ouchi Hisanao ouchi at joe.co.jp
Sun Jul 12 00:52:53 EDT 2015


Hi Andrey,

Thank you very much for your really useful advise.
According to your advise, I provided physical coordinates information of
"every primary unknown" (in other words, "degrees of freedom of the
matrix") to MueLu  and add the following sentence in my input deck.

<Parameter name="repartition: enable" type="bool" value="true"/>

Then, parallel performance was drastically improved in many problems in my
PC.
I will check MueLu performance more on TACC.

Thank you very much again.

Hisanao


2015-07-10 17:26 GMT-05:00 Andrey Prokopenko <aprokop at sandia.gov>:

>  Hi Hisanao,
>
> >Do you have coordinates available?
>
>  Yes, I can use my physical coordinates data for repartitioning in my
> simulation.
>
>  In our 2D hydraulic fracturing simulation, each elements has four
> unknowns (displacement of x and y, pore pressure, fracture fluid pressure).
> At that case, if the number of elements is N, do I have to provide the
> element positions of N elements?  In addition, how do I give the
> information that each element has four unknowns.
>
> What are the degrees of freedom of the matrix that you pass to MueLu? The
> coordinates should correspond to those. So, for instance, if the dofs are
> pressure, you should pass pressure coordinates. For instance, if your
> pressure is associated with a center of element, you should pass
> coordinates of the center.
>
>  > Also, how do you construct MueLu preconditioner?
>
>  I don't use Tpetra because I am using Aztecoo.
> I construct MueLu preconditioner by MueLu::EpetraOperator.
> In this case, how do I provide the coordinate information?
> Could you show me an example code to provide the physical coordinates?
>
>
> If you are using Epetra, there is a similar function called
> CreateEpetraPreconditioner, which also takes coordinates as 3rd argument.
> This is the preferable way to construct the preconditioner. Otherwise, if I
> understand correctly what you mean, you first need to construct a MueLu
> hierarchy, and then construct MueLu::EpetraOperator.
> CreateEpetraPreconditioner, on the other hand, does that work for you.
>
> -Andrey
>
>
> On 07/10/2015 04:20 PM, Ouchi Hisanao wrote:
>
> Hi Andrey,
>
>  >Do you have coordinates available?
>
>  Yes, I can use my physical coordinates data for repartitioning in my
> simulation.
>
>  In our 2D hydraulic fracturing simulation, each elements has four
> unknowns (displacement of x and y, pore pressure, fracture fluid pressure).
> At that case, if the number of elements is N, do I have to provide the
> element positions of N elements?  In addition, how do I give the
> information that each element has four unknowns.
>
>  > Also, how do you construct MueLu preconditioner?
>
>  I don't use Tpetra because I am using Aztecoo.
> I construct MueLu preconditioner by MueLu::EpetraOperator.
> In this case, how do I provide the coordinate information?
> Could you show me an example code to provide the physical coordinates?
>
>  Hisanao
>
>
>
>
>
>
> 2015-07-10 16:52 GMT-05:00 Andrey Prokopenko <aprokop at sandia.gov>:
>
>>  Oh, I see. In order to use repartitioning at the moment you have to
>> provide the physical coordinates of nodes corresponding to degrees of
>> freedom of the matrix you pass to MueLu. We require that information as we
>> use geometric partitioning to do rebalancing.
>>
>> Do you have coordinates available? Also, how do you construct MueLu
>> preconditioner? If it is through CreateTpetraPreconditioner, you can simply
>> pass these coordinates as a 3rd argument.
>>
>> -Andrey
>>
>>
>> On 07/10/2015 03:46 PM, Ouchi Hisanao wrote:
>>
>>   Hi Andrey,
>>
>>  Thank you very much for your reply.
>>
>>  According to your advice, I modified my input deck based on a new style
>> (please see the bottom of my email) and simply add the following sentence
>> as you suggested.
>>
>>  <Parameter name="repartition: enable" type="bool" value="true"/>
>>
>>  However, I could not execute my code successfully since the following
>> error message appeared.
>>
>>  -------------- ---------------------------- error messeage
>> -------------------------------------------
>> terminate called after throwing an instance of
>> 'MueLu::Exceptions::RuntimeError'
>>   what():  CoordinatesTransferFactory::DeclareInput: (LevelID = 0) unable
>> to find or generate requested data "Coordinates" with generating factory
>> "null"
>>     during request for data "    Coordinates" on level 0 by factory
>> CoordinatesTransferFactory
>>     during request for data "    Coordinates" on level 1 by factory
>> RebalanceTransferFactory
>>     during request for data "              P" on level 1 by factory
>> NoFactory
>>
>> -----------------------------------------------------------------------------------------------------------
>>
>>  This error message appears not only in my program in my Ubuntu but also
>> in the sample code in the MueLu virtual machine image distributed by Sandia
>> as a tutorial. In the case of the MueLu virtual machine image, I just added
>> the following sentence in the simplest example input deck  (n1_easy.xml)
>> and executed the simplest example problem in the MueLu tutorial (Laplace 2D
>> 50x50).
>>
>>  <Parameter name="repartition: enable" type="bool" value="true"/>
>>
>>  Do you have any idea to resolve this error?
>>
>>  Best regards,
>>
>>  Hisanao
>>
>>  ----------------------------- modified input deck
>> -------------------------------------
>> <ParameterList name="MueLu">
>>
>>    <Parameter name="verbosity" type="string" value="high"/>
>>
>>    <Parameter name="max levels" type="int" value="10"/>
>>   <Parameter name="coarse: max size" type="int" value="10"/>
>>
>>    <Parameter name="multigrid algorithm" type="string" value="sa"/>
>>   <Parameter name="sa: damping factor"  type="double" value="1.0"/>
>>
>>    <!-- Smoothing -->
>>   <!-- Comment/uncomment different sections to try different smoothers -->
>>   <Parameter name="smoother: type" type="string" value="ILUT"/>
>>
>>  <!-- Rebalancing -->
>>   <Parameter name="repartition: enable" type="bool" value="true"/>
>>
>>    <!-- Aggregation -->
>>   <Parameter name="aggregation: type" type="string" value="uncoupled"/>
>>   <Parameter name="aggregation: min agg size" type="int" value="4"/>
>>
>>  </ParameterList>
>>
>>
>>
>>
>>
>>
>>
>>
>>
>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://trilinos.org/pipermail/trilinos-users/attachments/20150711/60f822b9/attachment.html>


More information about the Trilinos-Users mailing list