getfem-users
[Top][All Lists]
Advanced

[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]

Re: [Getfem-users] High level generic assembly procedures


From: Marco Pischedda
Subject: Re: [Getfem-users] High level generic assembly procedures
Date: Wed, 2 Apr 2014 17:21:41 +0200

Hi,

I have other questions:

- I tried to assemble separately the following terms:
  "Grad_u.Test_p" , "p.Grad_Test_u", "p.Test_u" and the code works correctly.
  Then I want to assemble the sum of this terms, i.e: "Grad_u.Test_p +
p.Grad_Test_u+p.Test_u"
  but I receive the following error:

   Addition or substraction of incompatible expressions or of different sizes
   terminate called after throwing an instance of 'gmm::gmm_error'
   what():  Error in getfem_generic_assembly.cc, line 3591 :
   Error in assembly string
   Aborted

- I'm working in 1d problem with vectors as unknowns. Grad_u is
therefore a vector or is a        tensor? When I do Grad_u.Test_p the
result is a scalar or a vector? Test_p is a vector or a scalar?

Thanks in Advance

Marco



2014-04-01 17:51 GMT+02:00 Marco Pischedda <address@hidden>:
> Ok thank you,
>
> I tried it and it works. I will let you know if there are other problems.
>
> Thanks
>
> Marco
>
> 2014-04-01 17:16 GMT+02:00 Yves Renard <address@hidden>:
>>
>> Dear Marco,
>>
>> All seems to me correct in your implementation. This is probably just
>> the test line 2074 of getfem_generic_assembly.cc
>> which is not correct for one-dimensionnal problems.
>>
>> I think the test should be
>>
>>  GA_DEBUG_ASSERT((qdim == 1 && t.sizes()[0] == N) ||
>>                       (t.sizes()[1] == N && t.sizes()[0] == qdim) ||
>>                       (N == 1 && t.sizes()[0] == qdim),
>>                       "dimensions mismatch");
>>
>> May be you can try this. I will validate it if it works.
>>
>> Regards,
>>
>> Yves.
>>
>>
>>
>> Le 01/04/2014 16:56, Marco Pischedda a écrit :
>>> Dear Yves,
>>>
>>> thank you for your fast answer.
>>>
>>> I have another question:
>>>
>>> - I have a monodimensional problem but the unknowns are vectorial. I
>>> set the vectorial
>>> dimension of the unknows with mf_u.set_qdim(3) and with mf_p.set_qdim(3).
>>> When I define the expression " Grad_u.Test_p "  i receive the following 
>>> error:
>>>
>>> terminate called after throwing an instance of 'gmm::gmm_error'
>>>   what():  Error in
>>> ../../getfem-svn/getfem/src/getfem_generic_assembly.cc, line 2074 :
>>> dimensions mismatch
>>>
>>> Grad_u should be a 2nd order tensor, while Test_p should be a vector.
>>> I suppose the
>>> problem is that Grad_u is interpreted as a vector while Test_p is
>>> interpreted as a scalar.
>>> How can I use the high level generic assembly procedures on vectorial 
>>> problems
>>> defined on a monodimensional computational domain?
>>>
>>> Thanks in advance
>>>
>>> Marco
>>>
>>>
>>>
>>>
>>>
>>> 2014-04-01 8:42 GMT+02:00 Yves Renard <address@hidden>:
>>>> Dear Marco,
>>>>
>>>> Yes, your code is correct and it should work correctly.
>>>> Note that you can also use the high level generic assembly with the
>>>> model object using a generic assembly brick.
>>>> Note also that high level assembly is only available in the svn
>>>> repository version of Getfem
>>>> (and a very recent version, at least r4570 because a bug have been
>>>> corrected on coupled problems recently).
>>>>
>>>> Yves.
>>>>
>>>>
>>>> Le 31/03/2014 14:00, Marco Pischedda a écrit :
>>>>> Dear all,
>>>>>
>>>>> I want to use the high level generic assembly procedures in order to
>>>>> use the automatic differentation for calculating the Jacobian of a
>>>>> non-linear problem.
>>>>>
>>>>> For example I want to calculate the Jacobian of two vectorial
>>>>> equations in the unknowns "u" and "p". The first vectorial equation is
>>>>> (u+p).Test_u while the second is (u+p).Test_p.
>>>>>
>>>>>  It is correct the following code?
>>>>>
>>>>>  gmm::sub_interval Iu(0, ndofs_u);
>>>>>  gmm::sub_interval Ip(ndofs_u,ndofs_p);
>>>>>  workspace.add_fem_variable("u", mf_u, Iu, U);
>>>>>  workspace.add_fem_variable("p", mf_p, Ip, P);
>>>>>  workspace.add_expression("(u+p).Test_u + (u+p).Test_p",mim);
>>>>>  getfem::model_real_sparse_matrix Jac(ndofs_u+ndofs_p, ndofs_u+ndofs_p);
>>>>>  workspace.set_assembled_matrix(Jac);
>>>>>  workspace.assembly(2);
>>>>>
>>>>> In the matrix Jac I have the Jacobian of the system?
>>>>>
>>>>>
>>>>> Thanks in advance
>>>>>
>>>>> Marco
>>>>>
>>>>> _______________________________________________
>>>>> Getfem-users mailing list
>>>>> address@hidden
>>>>> https://mail.gna.org/listinfo/getfem-users
>>>>
>>>> --
>>>>
>>>>   Yves Renard (address@hidden)       tel : (33) 04.72.43.87.08
>>>>   Pole de Mathematiques, INSA-Lyon             fax : (33) 04.72.43.85.29
>>>>   20, rue Albert Einstein
>>>>   69621 Villeurbanne Cedex, FRANCE
>>>>   http://math.univ-lyon1.fr/~renard
>>>>
>>>> ---------
>>>>
>>
>>
>> --
>>
>>   Yves Renard (address@hidden)       tel : (33) 04.72.43.87.08
>>   Pole de Mathematiques, INSA-Lyon             fax : (33) 04.72.43.85.29
>>   20, rue Albert Einstein
>>   69621 Villeurbanne Cedex, FRANCE
>>   http://math.univ-lyon1.fr/~renard
>>
>> ---------
>>



reply via email to

[Prev in Thread] Current Thread [Next in Thread]