← Back to team overview

dolfin team mailing list archive

Re: Additions to assembler interface

 

On 27 September 2011 10:36, Andre Massing <massing@xxxxxxxxx> wrote:
> -----BEGIN PGP SIGNED MESSAGE-----
> Hash: SHA1
>
> On 09/26/2011 09:21 PM, Garth N. Wells wrote:
>> On 26 September 2011 15:21, Andre Massing <massing@xxxxxxxxx>
>> wrote:
>>> On 09/26/2011 01:15 PM, Anders Logg wrote:
>>>> On Mon, Sep 26, 2011 at 12:06:55AM -0700, Johan Hake wrote:
>>>>> On Sunday September 25 2011 23:54:57 Anders Logg wrote:
>>>>>> As part of assembling systems on overlapping meshes (for
>>>>>> Nitsche type methods), my student Andre Massing is calling
>>>>>> the regular DOLFIN assemble function for part of the domain
>>>>>> and then calling a special purpose assembler (currently
>>>>>> residing in a branch called dolfin-olm) for terms that are
>>>>>> not supported by the regular assembler.
>>>>>>
>>>>>> For this to work out well, I'd like to propose the
>>>>>> following changes to the assembler:
>>>>>>
>>>>>> 1. Make the assemble_cells,
>>>>>> assemble_exterior/interior_facets public
>>>>>>
>>>>>> We need to access these separately.
>>>>>>
>>>>>> 2. Add options for controlling if and how we call
>>>>>> A.apply()
>>>>>
>>>>> Wouldn't this be unnesseary if 1 is implemented? Then you
>>>>> just impement your own more flexible assemble routine?
>>>>
>>>> Good point.
>>>>
>>>> Andre: Would that be enough? Note that apply() is not called
>>>> from the assemble_cells/interior/exterior_facets functions,
>>>> only from the main (currently public) assemble function.
>>>
>>> In principle that would be okay. But it enforces the user to
>>> reimplement (or basically to copy paste) the entire
>>>
>>> void Assembler::assemble(GenericTensor& A, const Form& a, const
>>> MeshFunction<uint>* cell_domains, const MeshFunction<uint>*
>>> exterior_facet_domains, const MeshFunction<uint>*
>>> interior_facet_domains, bool reset_sparsity, bool add_values)
>>>
>>> even if the user only wants to first do an regular assemble (via
>>> that function) and then add/insert manually new entries. Top copy
>>> paste the entire assemble function just to get rid of one apply
>>> call at the end seems not to be really user friendly.
>>>
>>> One nice example is the use of the ident_zeros function which as
>>> of today is "damm" slow :) if there are many zeros along the
>>> diagonal. That is because the Matrix will be called with an
>>> apply() at the end of the assemble routine and hence adding
>>> values afterwards via ident_zeros is very expensive.
>>
>> My understanding is that it does not add entries.
>
> Yes, your are right, I did not mean adding in the sense of altering
> the sparsity pattern, merely adding 1 to zero entries. Nevertheless it
> is slow after an apply call since Petsc then does compressing etc.
> which makes adding even existing elements very slow.
>
>> Someone did add some hacks which I never understood the need for.
>
> Which hack do you mean? The ident_zeros function? At least that one has
> proven rather useful in cases where you define the FEM on the whole
> mesh but only assemble part of the mesh. IIRC that was also needed in
> Kristoffer's FSI implementation and it is very useful in the
> overlapping mesh case where you have to constraint those dofs which
> did not contribute to the solution at all (since the the corresponding
> cells where completely overlapped).
>
> A good solution is to not
>> use Assembler but SystemAssembler. The latter should really be the
>> default in DOLFIN since it preserves symmetry.
>
> System assembler is a good solution to what exactly? You lost me :)
>

With  SystemAssembler, essential boundary conditions are applied
during the assembly process, therefore there is no call to zero out
rows in order to apply bcs.

>>
>> I haven't tested, but I don't know if, for PETSc, MatZeroRows can
>> be called post set/insert but pre apply.
>
> How does MatZerosRow exactly come into play?
>

MatZerosRow is the PETSc function that is eventually called when
zeroing rows and placing 1 on the diagonal.

Garth


> - --
> Andre
>
>>
>> The function Assembler::assemble referred to above could have
>> something like 'bool finalise_matrix=true' added to the function
>> signature.
>>
>> Garth
>>
>>> Naturally a user want to assemble,  ident_zeros and then call
>>> apply. Basically copy pasting the the whole assemble function
>>> only to omit the postpone the one apply call at the end of
>>> assemble function seems to be an overkill.
>>>
>>> So one thing is the possibility to omit the apply at the end the
>>> assemble. The other thing is to make assemble_cells, etc. public
>>> or protected to ease the implementation of new assemblers. As you
>>> said there are related, but I can think of simple scenarios where
>>> a completely new assembler is not needed but just a tiny function
>>> to add some elements. There, the missing option of postponing the
>>> apply makes it rather cumbersome to reuse the nice functionality
>>> which already is provided in DOLFIN.
>>>
>>> -- Andre
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>
>>>>
>>>> -- Anders
>>>>
>>>>
>>>>> Johan
>>>>>
>>>>>> A.apply() calls PETSc MatAssemblyBegin/End with
>>>>>> MAT_FINAL_ASSEMBLY. This makes it *very* expensive for us
>>>>>> since we need to switch between inserting and adding values
>>>>>> to the matrix. We need to instead use the option
>>>>>> MAT_FLUSH_ASSEMBLY. It would further be useful to have an
>>>>>> option of not calling A.apply() at all.
>>>>>>
>>>>>> I'm not sure how this should be handled in the
>>>>>> GenericTensor interface since Trilinos does not seem to
>>>>>> have a corresponding option.
>>>>>>
>>>>>> A first step would be to just add a boolean option to the
>>>>>> assembler:
>>>>>>
>>>>>> call_apply=false
>>>>>>
>>>>>> Then we could manually flush the values in between calls.
>>>>>>
>>>>>> Opinions?
>>>>>>
>>>>
>>>> _______________________________________________ Mailing list:
>>>> https://launchpad.net/~dolfin Post to     :
>>>> dolfin@xxxxxxxxxxxxxxxxxxx Unsubscribe :
>>>> https://launchpad.net/~dolfin More help   :
>>>> https://help.launchpad.net/ListHelp
>>>
>>>
>>> _______________________________________________ Mailing list:
>>> https://launchpad.net/~dolfin Post to     :
>>> dolfin@xxxxxxxxxxxxxxxxxxx Unsubscribe :
>>> https://launchpad.net/~dolfin More help   :
>>> https://help.launchpad.net/ListHelp
>>>
>
> -----BEGIN PGP SIGNATURE-----
> Version: GnuPG v2.0.18 (GNU/Linux)
> Comment: Using GnuPG with Mozilla - http://enigmail.mozdev.org/
>
> iQEcBAEBAgAGBQJOgZkIAAoJEA79ggnbq9dm9igH/2n2YHygDXYZCt9AI2qkyA2k
> xjylLmNrvoRzbXO+HB90Xf039xVcS91TnlDKK6Co49i+uGb4rmLaSK1SltmDIJQE
> 5xCix5PlE2kxf/ACPg9z31oaA8erEFDXeKIwXZMCl8L8S0sc15F74s3SEDzz5alF
> xTokgb/9Cnjg+r0etawBy+rWWWK3m2sc/vORzUcl5h9liGQ//4fUd1/3PZuqiU7r
> TLW6Xi1VhrtS/OOGUMF5HVo995XUrkVvblYxOlzp5YFdKgm/d9ho1u3J3NiA94dR
> S/nqdi5M38oCGXOGYDjMF/4NLNGOHIX6dCIZJeRBPJ3EoDdld8gjEBhhKDbJWx8=
> =DCLT
> -----END PGP SIGNATURE-----
>
> _______________________________________________
> Mailing list: https://launchpad.net/~dolfin
> Post to     : dolfin@xxxxxxxxxxxxxxxxxxx
> Unsubscribe : https://launchpad.net/~dolfin
> More help   : https://help.launchpad.net/ListHelp
>


References