[petsc-users] Create a DM given sets of IS's
Matthew Knepley
knepley at gmail.com
Wed Jan 2 16:38:54 CST 2019
On Wed, Jan 2, 2019 at 4:55 PM Justin Chang <jychang48 at gmail.com> wrote:
> Okay that's fine.
>
> Now consider one more case:
>
> Suppose I create a PetscSection with 9 fields:
>
> section = PETSc.Section().create()
> section.setNumFields(9)
> section.setFieldName(0,'u1')
> section.setFieldName(1,'c1')
> section.setFieldName(2,'t1')
> section.setFieldName(3,'u2')
> section.setFieldName(4,'c2')
> section.setFieldName(5,'t2')
> section.setFieldName(6,'u3')
> section.setFieldName(7,'c3')
> section.setFieldName(8,'t3')
> section.setFieldComponents(0,1)
> section.setFieldComponents(1,1)
> section.setFieldComponents(2,1)
> section.setFieldComponents(3,1)
> section.setFieldComponents(4,1)
> section.setFieldComponents(5,1)
> section.setFieldComponents(6,1)
> section.setFieldComponents(7,1)
> section.setFieldComponents(8,1)
>
> ....
>
> where u is potential, c is concentration, and t is temperature. The
> numbers refer to different domains of a battery (anode, electrolyte,
> cathode).
>
> I want three levels of splits:
> 1) Split by domain
> 2) Split potential and concentration from temperature.
> 3) Split potential and concentration for schur complementing purpose.
>
> Do the following PETSc command line options describe the split I mentioned
> above:
>
> -pc_type fieldsplit
> -pc_fieldsplit_0_fields 0,1,2
> -pc_fieldsplit_1_fields 3,4,5
> -pc_fieldsplit_2_fields 6,7,8
> -fieldsplit_0_pc_type fieldsplit
> -fieldsplit_1_pc_type fieldsplit
> -fieldsplit_2_pc_type fieldsplit
> -fieldsplit_0_pc_fieldsplit_0_fields 0,1
> -fieldsplit_0_pc_fieldsplit_1_fields 2
> -fieldsplit_1_pc_fieldsplit_0_fields 3,4
> -fieldsplit_1_pc_fieldsplit_1_fields 5
> -fieldsplit_2_pc_fieldsplit_0_fields 6,7
> -fieldsplit_2_pc_fieldsplit_1_fields 8
>
> I don't have an actual code to test these out with yet as I've never done
> anything beyond two levels of split. So I am wondering if the above will
> break anything as far as you can tell.
>
This should work as you expect. Send me the -ksp_view output when you run
it. It would be very cool to see 3-level in action.
Thanks,
Matt
> Thanks,
> Justin
>
> On Tue, Jan 1, 2019 at 7:10 AM Matthew Knepley <knepley at gmail.com> wrote:
>
>> On Tue, Jan 1, 2019 at 5:06 AM Justin Chang <jychang48 at gmail.com> wrote:
>>
>>> Okay so I managed to successfully do this by translating the IS's into a
>>> PetscSection, assigning it to an empty DMShell, and then assigning that to
>>> KSP.
>>>
>>> However, I seem to be unable to rename the outer fields that aggregate
>>> two or more inner fields.
>>>
>>
>> Right now you cannot do that because of the way that I check for the
>> option:
>>
>>
>> https://bitbucket.org/petsc/petsc/src/0a5f382fff62a328ec919e3e9fe959e6cbbcf413/src/ksp/pc/impls/fieldsplit/fieldsplit.c#lines-369
>>
>> I guess we could replace this by some search code that looked for any
>> option of that form and then read
>> out the splitname using scanf. Right now, I do not see a pattern search
>> function for the options database.
>>
>> Thanks,
>>
>> Matt
>>
>>
>>> Consider this snippet of a four-field dual porosity/permeability
>>> FEniCS/petsc4py code:
>>>
>>> ...
>>>
>>> ## Extract FEniCS dof layout, global indices ##
>>> dof_total = np.array(W.dofmap().dofs())
>>> dof_v1 = np.array(W.sub(0).dofmap().dofs())
>>> dof_p1 = np.array(W.sub(1).dofmap().dofs())
>>> dof_v2 = np.array(W.sub(2).dofmap().dofs())
>>> dof_p2 = np.array(W.sub(3).dofmap().dofs())
>>> offset = np.min(dof_total)
>>>
>>> ## Create PetscSection ##
>>> section = PETSc.Section().create()
>>> section.setNumFields(4)
>>> section.setFieldName(0,'v1')
>>> section.setFieldName(1,'p1')
>>> section.setFieldName(2,'v2')
>>> section.setFieldName(3,'p2')
>>> section.setFieldComponents(0,1)
>>> section.setFieldComponents(1,1)
>>> section.setFieldComponents(2,1)
>>> section.setFieldComponents(3,1)
>>> section.setChart(0,len(dof_total))
>>> for i in np.nditer(dof_v1):
>>> section.setDof(i-offset,1)
>>> section.setFieldDof(i-offset,0,1)
>>> for i in np.nditer(dof_p1):
>>> section.setDof(i-offset,1)
>>> section.setFieldDof(i-offset,1,1)
>>> for i in np.nditer(dof_v2):
>>> section.setDof(i-offset,1)
>>> section.setFieldDof(i-offset,2,1)
>>> for i in np.nditer(dof_p2):
>>> section.setDof(i-offset,1)
>>> section.setFieldDof(i-offset,3,1)
>>> section.setUp()
>>>
>>> ## Create DM and assign PetscSection ##
>>> dm = PETSc.DMShell().create()
>>> dm.setDefaultSection(section)
>>> dm.setUp()
>>>
>>> ## Create KSP and assign DM ##
>>> ksp = PETSc.KSP().create()
>>> ksp.setDM(dm)
>>> ksp.setDMActive(False)
>>>
>>> ### PETSc Command-line options ##
>>> PETScOptions.set('ksp_monitor_true_residual')
>>> PETScOptions.set('ksp_view')
>>> PETScOptions.set('ksp_type','gmres')
>>> PETScOptions.set('pc_type','fieldsplit')
>>> PETScOptions.set('pc_fieldsplit_0_fields','0,1')
>>> PETScOptions.set('pc_fieldsplit_1_fields','2,3')
>>> PETScOptions.set('pc_fieldsplit_type','additive')
>>> PETScOptions.set('fieldsplit_0_ksp_type','preonly')
>>> PETScOptions.set('fieldsplit_0_pc_type','fieldsplit')
>>> PETScOptions.set('fieldsplit_0_pc_fieldsplit_type','schur')
>>> PETScOptions.set('fieldsplit_0_pc_fieldsplit_schur_fact_type','full')
>>> PETScOptions.set('fieldsplit_0_pc_fieldsplit_schur_precondition','selfp')
>>> PETScOptions.set('fieldsplit_0_fieldsplit_v1_ksp_type','preonly')
>>> PETScOptions.set('fieldsplit_0_fieldsplit_v1_pc_type','bjacobi')
>>> PETScOptions.set('fieldsplit_0_fieldsplit_p1_ksp_type','preonly')
>>> PETScOptions.set('fieldsplit_0_fieldsplit_p1_pc_type','hypre')
>>> PETScOptions.set('fieldsplit_1_ksp_type','preonly')
>>> PETScOptions.set('fieldsplit_1_pc_type','fieldsplit')
>>> PETScOptions.set('fieldsplit_1_pc_fieldsplit_type','schur')
>>> PETScOptions.set('fieldsplit_1_pc_fieldsplit_schur_fact_type','full')
>>> PETScOptions.set('fieldsplit_1_pc_fieldsplit_schur_precondition','selfp')
>>> PETScOptions.set('fieldsplit_1_fieldsplit_v2_ksp_type','preonly')
>>> PETScOptions.set('fieldsplit_1_fieldsplit_v2_pc_type','bjacobi')
>>> PETScOptions.set('fieldsplit_1_fieldsplit_p2_ksp_type','preonly')
>>> PETScOptions.set('fieldsplit_1_fieldsplit_p2_pc_type','hypre')
>>>
>>> ...
>>>
>>> Is it possible to rename the outer splits (aka make
>>> '-pc_fieldsplit_0_fields 0,1' be something like
>>> '-pc_fieldsplit_micro_fields 0,1')?
>>>
>>> Thanks,
>>> Justin
>>>
>>>
>>> On Mon, Dec 31, 2018 at 7:40 AM Matthew Knepley <knepley at gmail.com>
>>> wrote:
>>>
>>>> On Mon, Dec 31, 2018 at 2:40 AM Justin Chang via petsc-users <
>>>> petsc-users at mcs.anl.gov> wrote:
>>>>
>>>>> Hi all,
>>>>>
>>>>> I am solving a six field battery problem (concentration and potential
>>>>> for each of the two solid and one electrolyte domains) and I want to
>>>>> experiment with nested/recursice fieldsplitting. I have the IS's and am
>>>>> able to use these to define my PCFieldSplitsSetIS(). However, I can imagine
>>>>> this getting really messy from a programming standpoint, especially once I
>>>>> need to add temperature into the mix, so it is my hope that I can translate
>>>>> these index sets and fields into a DM (maybe DMShell?) so that I can just
>>>>> rely on command line options to play around with various combinations of
>>>>> field assignments and splits (e.g. -pc_fieldsplit_X_fields)
>>>>>
>>>>> However, it doesn't seem clear to me how I would create a DM when you
>>>>> already know the IS's for each fields. If I understand functions like
>>>>> DMCreateFieldDecomposition() correctly, it seems that it returns to you the
>>>>> IS's and sub DM's associated with the original DM, whereas I want to do it
>>>>> the other way around. Perhaps the "reversal" of something like
>>>>> DMCreateFieldIS()
>>>>> <https://www.mcs.anl.gov/petsc/petsc-dev/src/dm/interface/dm.c.html#DMCreateFieldIS>,
>>>>> where you convert the IS into a PetscSection or is there an easier/better
>>>>> way?
>>>>>
>>>>> Any thoughts/help appreciated!
>>>>>
>>>>
>>>> Paul has recently done this for LibMesh. I believe that constructing a
>>>> PetscSection is enough to get you minimally started. That allows
>>>> DMCreateSubDM() to work by subsetting the Section, and that should
>>>> allow the command line to work. CreateFieldDecomposition() should
>>>> be removed I think.
>>>>
>>>> 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/>
>>
>
--
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/20190102/bbe646fa/attachment.html>
More information about the petsc-users
mailing list