DMMG Time dependent matrix

jens.madsen at risoe.dk jens.madsen at risoe.dk
Mon Sep 18 15:40:24 CDT 2006


Hi Again 

Really hope that you have an answer :-) My problem is that the MG code gets slower and slower the longer it runs. For example lets say that that it took 10 sec to get from iteration 1000 to iteration 2000. If I use the fields from iteration 1000 as initial condition and start the code it takes only 7 sec to calculate 1000 iterations. Note that the solutions are identical. I do not think that this is caused by a memory leaks; I have logged the memory use and there is no increase in memory use when using "top" or -malloc_info/log. 

I think that the "slowing down" is happening because I am using the MG matrices on the coarser levels from the first iteration for all iterations. More specificly the amat and pmat belonging to mg[i]->smoothd(u) are not being updated ? When I am testing with a time independent matrix the code does not slow down at all.  

I am sure that the coarser level matrices are calculated/restricted for each timestep. I have tested this simply by subtracting the coarser matrices for different timesteps. 

Do you know whether this: ierr = KSPGetPC(dmmg[i]->ksp,&pc);
				    ierr = PCMGGetSmoother(pc,i,&lksp);CHKERRQ(ierr);
				    KSPSetOperators(lksp,dmmg[i]->B,dmmg[i]->B,DIFFERENT_NONZERO_PATTERN);	

will do the job ? I have tried a lot of different things.. I have noticed that there is a comment in mg.c in PCSetUp_MG. Do not see whether there is a connection ? 

Cheers Jens
 


-----Original Message-----
From: Barry Smith [mailto:bsmith at mcs.anl.gov]
Sent: Mon 9/11/2006 1:50 AM
To: petsc-users at mcs.anl.gov; jens.madsen at risoe.dk
Subject: Re: DMMG Time dependent matrix
 

   Jens,

    You are correct; functionality to do this does not exist in DMMG. So you
have to write it yourself as you have, let us know if you have any trouble
with it.

    Barry


On Sun, 10 Sep 2006, jens.madsen at risoe.dk wrote:

> Hi again
>
> My linear problem Ax = b. In my case the matrix A is timepdependent. In my KSP code I calculate A and call KSPSetOperators() in each timestep. When using DMMG I have not been able to find such a functionallity ? I am using galerkin matrices on all coarser MG levels. Here is the code that I am currently using
>
> /*Matrix is time dependent*/
>   			if(para.PolEq==GLOBAL)
> 			{
> 				ierr = MGComputePolaMatrix(*dmmg,dmmg[para.MG_levels-1]->J,dmmg[para.MG_levels-1]->B);
> 				/*ierr = MatView(dmmg[DMMGGetLevels(dmmg)-1]->B,PETSC_VIEWER_STDOUT_WORLD);*/
>    			for (i=DMMGGetLevels(dmmg)-2; i>-1; i--)
>    			{
> 	       			if (dmmg[i]->galerkin)
> 	       			{
>         				MatPtAP(dmmg[i+1]->B,dmmg[i+1]->R,MAT_REUSE_MATRIX,1.0,&dmmg[i]->B);
> 	         			if (!dmmg[i]->J)
>    	     			{
>  			         		dmmg[i]->J = dmmg[i]->B;
>         				}
>       				}
> 					/*ierr = KSPSetOperators(dmmg[i]->ksp,dmmg[i]->B,dmmg[i]->B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);*/
					dmmg[i]->matricesset = PETSC_TRUE;
> 				}
                            for(i=0;i<dmmg[0]->nlevels;i++)
			    {
			    	ierr = KSPGetPC(dmmg[i]->ksp,&pc);
				    ierr = PCMGGetSmoother(pc,i,&lksp);CHKERRQ(ierr);
				    KSPSetOperators(lksp,dmmg[i]->B,dmmg[i]->B,DIFFERENT_NONZERO_PATTERN);				    
			    }
> 		    }
> 	     	else
> 	     	{
> 				ierr = KSPSetOperators(dmmg[i]->ksp,dmmg[i]->B,dmmg[i]->B,SAME_PRECONDITIONER);CHKERRQ(ierr);
> 		    }

> 		 	ierr = DMMGSolve(dmmg);CHKERRQ(ierr);
>
> I this the way to do it or have missed some functionallity in DMMG  ?
>
> Cheers Jens
>
> -----Original Message-----
> From: Barry Smith [mailto:bsmith at mcs.anl.gov]
> Sent: Thu 8/31/2006 10:12 PM
> To: jens.madsen at risoe.dk
> Subject: Re: SV: SV: SV: [PETSC #14613] DMMG question; more stupid danish questions :-)
>
>
>   Yes, one of the arguments to DACreate is which directions you want
> periodicity.
>
>     Barry
>
>
> On Thu, 31 Aug 2006, jens.madsen at risoe.dk wrote:
>
>> Ok I had not figured that out. Think I misunderstood a response from you that you mailed in the winter :-)
>>
>> I am now trying to implement DMMG in my code. Just one question: is DMMG able to handle periodic boundary conditions ? I think I have made it work but I have not performed intensive testing ....
>>
>> thx Jens :-)
>>
>> ________________________________
>>
>> Fra: Barry Smith [mailto:bsmith at mcs.anl.gov]
>> Sendt: on 23-08-2006 00:06
>> Til: jens.madsen at risoe.dk
>> Cc: petsc-maint at mcs.anl.gov
>> Emne: Re: SV: SV: [PETSC #14613] DMMG question
>>
>>
>>
>>
>>  Actually DMMG does not always start at the coarsest level.
>> It will only do that if you pass in the -dmmg_grid_sequence option.
>> You can run with -pc_mg_type multiplicative -pc_mg_cycles 1 for
>> V cycle and -pc_mg_cycles 2 for W cycles.
>>
>>   There may be other reasons you might not be able to use DMMG like
>> it only works on logically rectangular grids, etc.
>>
>> The only example for the mg directly is src/ksp/ksp/examples/tests/ex19.c
>> it is terribly ugly.
>>
>>
>>  Barry
>>
>> It is a different ex19.c
>>
>>
>>
>> On Tue, 22 Aug 2006 jens.madsen at risoe.dk wrote:
>>
>>> Hi again
>>>
>>> As far as I can see petsc/src/snes/examples/tutorials/ex19.c uses DMMG ?
>>>
>>> I am trying to use the preconditioner itself. I do not want to use DMMG because DMMG always starts at the coarsest grid in order to improve the initial guess. I want to start on the finest grid in each timestep and do a V or W...
>>>
>>>
>>>
>>> Kind regards
>>>
>>>
>>>
>>> JEns
>>>
>>>
>>>
>>>
>>>
>>>
>>> ________________________________
>>>
>>> Fra: Hong Zhang [mailto:petsc-maint at mcs.anl.gov]
>>> Sendt: lø 19-08-2006 03:04
>>> Til: jens.madsen at risoe.dk
>>> Cc: petsc-maint at mcs.anl.gov
>>> Emne: RE: SV: [PETSC #14613] DMMG question
>>>
>>>
>>>
>>>
>>> See
>>> ~petsc/src/snes/examples/tutorials/ex19.c
>>>
>>> Hong
>>>
>>> On Fri, 18 Aug 2006 jens.madsen at risoe.dk wrote:
>>>
>>>> Hi again
>>>>
>>>> First I will thank you guys for developing such a great product. Since i wrote the mails below I have developed a code simulating 2D plasmaphysics. At the moment I am using the Krylov subspace methods, KSP. Now that this is runnning I would like to see which improvements a MG precoditioner would give. Before starting doing so I would ask you guys whether you have any examples ? The documentation on the PCMG is a bit sparse (:-))
>>>>
>>>> cheers Jens
>>>>
>>>> PS. I will not forget to site you @:-)
>>>>
>>>>
>>>> -----Original Message-----
>>>> From: Barry Smith [mailto:petsc-maint at mcs.anl.gov]
>>>> Sent: Thu 3/30/2006 2:35 AM
>>>> To: jens.madsen at risoe.dk
>>>> Cc: petsc-maint at mcs.anl.gov
>>>> Subject: Re: SV: [PETSC #14613] DMMG question
>>>>
>>>>
>>>>
>>>> On Wed, 29 Mar 2006, jens.madsen at risoe.dk wrote:
>>>>
>>>>> Hi Barry and Matt
>>>>>
>>>>> Thank you for you quick response :-)
>>>>>
>>>>> I have actually tried those command line options that you suggest. It is probably me not using the right terminology; I am not that experienced. I try to rephrase my problem.
>>>>>
>>>>> As far as I understand the full multigrid always starts on the coarsest grid G_0. On the coarsest grid you calculate an approximate solution v_0 (using jacobi, GaussSeiddel, Krylov etc). This approximate solution is interpolated to a fine grid and used as an initial guess on this finer grid G_1. Now you make a few iterative sweeps on G_1 smoothening out the high k modes and get v_1. Restrict this approximate solution to G_0. Relax on G_0. Correct the solution v_1 and use this as an initial guess on an even finer grid G_2 etc. etc. In other words it combines "nested iteration" and "coarse grid correction".
>>>>
>>>>    You are correct this is exactly traditional full multigrid.  The "PETSc full multigrid" is slightly
>>>> different. We start with a right hand side (and initial guess) on the finest grid, restrict the residual
>>>> to the coarsest grid and then start up the grids with nested iteration AND coarse grid correction.
>>>> Thus unlike "traditional" full multigrid you only need to define your problem on the finest grid.
>>>>
>>>> We've found that this usually works better than just using regular V or W cycles. (BTW: PETSc full multigrid
>>>> can, of course, use either V or W cycles.
>>>>
>>>>>
>>>>> My "algorithm" is as follows:
>>>>>
>>>>> 1) apply initial conditions to Density,n, and Temperature,T.
>>>>> 2) find \phi solving a Poisson like equation using a multigrid scheme. Use \phi from previous timestep as an initial guess. n and T are variables in this equation.
>>>>> 3) Step n and T forward in time using the "stiffly stable" time stepping scheme. This is to be done using a pre-LU-factorized matrix.
>>>>> 4) goto 2)
>>>>>
>>>>> I would like to make a V (or W) cycle starting on the finest grid instead. On the finest grid I would like to apply boundary conditions, provide the previous time step as an initial guess and use "coarse grid correction" only. My own code is (probably) full of errors :-) so I have tried running the multigrid examples under KSP and SNES with -pc_type richardson and/or, -pc_mg_type additive and/or -ksp_type preonly etc.What I can see using -ksp/snes_view is that the starting point is always the coarsest grid with dimensions given by DACreateNd ? The next MG-Grids are always finer. Can I make DMMG start on the finest grid with dimensions given in DACreateNd ?
>>>>
>>>>    You need to create a DA with a coarser size so that after it is refined the number of times it gives you
>>>> the grid you want. That is if you want 5 grid points with two levels you would create a DA with 3 grid
>>>> points and pass that into the DMMG. Sorry there is no way to start with a DA on the finest (but you get the
>>>> same effect by starting with a coarser DA.
>>>>
>>>>     Barry
>>>>>
>>>>> Cheers and thanks Jens
>>>>>
>>>>> ________________________________
>>>>>
>>>>> Fra: Barry Smith [mailto:petsc-maint at mcs.anl.gov]
>>>>> Sendt: on 29-03-2006 17:20
>>>>> Til: jens.madsen at risoe.dk
>>>>> Cc: petsc-maint at mcs.anl.gov
>>>>> Emne: Re: [PETSC #14613] DMMG question
>>>>>
>>>>>
>>>>>
>>>>>
>>>>>   Jens,
>>>>>
>>>>>    YOu can access any of the low level PCMG options from
>>>>> DMMG. For example, -pc_type richardson gives you "multigrid
>>>>> as a solve". -pc_mg_type additive gives you additive
>>>>> -pc_mg_type multiplicative gives you standard v or w cycel
>>>>> -pc_mg_type full gives "full" multigrid".  -pg_mg_cycles 2
>>>>> gives W cycle. -pc_mg_levels_pc_type sor gives SOR as the smoother
>>>>> etc etc etc. Run with -help to see all the choices for the parts
>>>>> of the multigrid process.
>>>>>
>>>>>    Barry
>>>>>
>>>>> As Matt noted, using the PCMG directly requires YOU provide
>>>>> a mesh infrastructure that manages the meshes, thus it is not
>>>>> realistic for us to provide this whole infrastructure in an
>>>>> example. We currently only provide the full infrastructure
>>>>> for structured grids (in DMMG).
>>>>>
>>>>> On Wed, 29 Mar 2006, jens.madsen at risoe.dk wrote:
>>>>>
>>>>>> Hi Petsc :-)
>>>>>>
>>>>>>
>>>>>>
>>>>>> I am trying to use Petsc for solving plasma fluid equations. Is it
>>>>>> possible to use the DMMG with multiplicative or additive multigrid
>>>>>> schemes ? Also can I use multigrid as a solver not only as a
>>>>>> preconditioner in DMMG ?
>>>>>>
>>>>>>
>>>>>>
>>>>>> Also I cannot find any examples using the "low-level"  PCMG multigrid
>>>>>> interface ? Only a testprogram in the PC/TEST directory ?
>>>>>>
>>>>>>
>>>>>>
>>>>>> Cheers Jens Madsen
>>>>>>
>>>>>>
>>>>>
>>>>>
>>>>>
>>>>>
>>>>
>>>>
>>>>
>>>
>>>
>>>
>>>
>>
>>
>>
>>
>
>




More information about the petsc-users mailing list