OpenCMISS-Iron Internal API Documentation
interface_conditions_routines.f90
Go to the documentation of this file.
1 
43 
45 MODULE interface_conditions_routines
46 
47  USE base_routines
48  USE basis_routines
49  USE field_routines
50  USE input_output
54  USE interface_matrices_routines
55  USE interface_operators_routines
57  USE kinds
58  USE matrix_vector
59  USE strings
60  USE timer
61  USE types
62 
63 #include "macros.h"
64 
65  IMPLICIT NONE
66 
67  !Module types
68 
69  !Module variables
70 
71  !Interfaces
72 
73  PUBLIC interface_condition_create_finish,interface_condition_create_start
74 
75  PUBLIC interface_condition_dependent_variable_add
76 
77  PUBLIC interface_condition_destroy
78 
79  PUBLIC interface_condition_equations_create_finish,interface_condition_equations_create_start
80 
81  PUBLIC interface_condition_equations_destroy
82 
83  PUBLIC interfacecondition_integrationtypeget,interfacecondition_integrationtypeset
84 
85  PUBLIC interfacecondition_lagrangefieldcreatefinish,interfacecondition_lagrangefieldcreatestart
86 
87  PUBLIC interface_condition_method_get,interface_condition_method_set
88 
89  PUBLIC interface_condition_operator_get,interface_condition_operator_set
90 
91  PUBLIC interfacecondition_penaltyfieldcreatefinish,interfacecondition_penaltyfieldcreatestart
92 
93  PUBLIC interface_condition_user_number_find
94 
95  PUBLIC interface_conditions_finalise,interface_conditions_initialise
96 
97 CONTAINS
98 
99  !
100  !================================================================================================================================
101  !
102 
104  SUBROUTINE interface_condition_assemble(INTERFACE_CONDITION,ERR,ERROR,*)
105 
106  !Argument variables
107  TYPE(interface_condition_type), POINTER :: interface_condition
108  INTEGER(INTG), INTENT(OUT) :: err
109  TYPE(varying_string), INTENT(OUT) :: error
110  !Local Variables
111  TYPE(interface_equations_type), POINTER :: interface_equations
112  TYPE(varying_string) :: local_error
113 
114  enters("INTERFACE_CONDITION_ASSEMBLE",err,error,*999)
115 
116 ! CALL WRITE_STRING(GENERAL_OUTPUT_TYPE,"********************Interface Condition Assemble******************",ERR,ERROR,*999)
117  IF(ASSOCIATED(interface_condition)) THEN
118  interface_equations=>interface_condition%INTERFACE_EQUATIONS
119  IF(ASSOCIATED(interface_equations)) THEN
120  IF(interface_equations%INTERFACE_EQUATIONS_FINISHED) THEN
121  SELECT CASE(interface_condition%METHOD)
123  CALL interface_condition_assemble_fem(interface_condition,err,error,*999)
125  CALL flagerror("Not implemented.",err,error,*999)
127  CALL flagerror("Not implemented.",err,error,*999)
128  CASE DEFAULT
129  local_error="The interface condition method of "// &
130  & trim(number_to_vstring(interface_condition%METHOD,"*",err,error))// &
131  & " is invalid."
132  CALL flagerror(local_error,err,error,*999)
133  END SELECT
134  ELSE
135  CALL flagerror("Interface equations have not been finished.",err,error,*999)
136  ENDIF
137  ELSE
138  CALL flagerror("Interface condition interface equations is not associated.",err,error,*999)
139  ENDIF
140  ELSE
141  CALL flagerror("Interface condition is not associated.",err,error,*999)
142  ENDIF
143 
144  exits("INTERFACE_CONDITION_ASSEMBLE")
145  RETURN
146 999 errorsexits("INTERFACE_CONDITION_ASSEMBLE",err,error)
147  RETURN 1
148  END SUBROUTINE interface_condition_assemble
149 
150  !
151  !================================================================================================================================
152  !
153 
155  SUBROUTINE interface_condition_assemble_fem(INTERFACE_CONDITION,ERR,ERROR,*)
156 
157  !Argument variables
158  TYPE(interface_condition_type), POINTER :: interface_condition
159  INTEGER(INTG), INTENT(OUT) :: err
160  TYPE(varying_string), INTENT(OUT) :: error
161  !Local Variables
162  INTEGER(INTG) :: element_idx,ne,number_of_times
163  REAL(SP) :: element_user_elapsed,element_system_elapsed,user_elapsed,user_time1(1),user_time2(1),user_time3(1),user_time4(1), &
164  & USER_TIME5(1),USER_TIME6(1),SYSTEM_ELAPSED,SYSTEM_TIME1(1),SYSTEM_TIME2(1),SYSTEM_TIME3(1),SYSTEM_TIME4(1), &
165  & SYSTEM_TIME5(1),SYSTEM_TIME6(1)
166  TYPE(domain_mapping_type), POINTER :: elements_mapping
167  TYPE(interface_equations_type), POINTER :: interface_equations
168  TYPE(interface_matrices_type), POINTER :: interface_matrices
169  TYPE(field_type), POINTER :: lagrange_field
170 
171 !#ifdef TAUPROF
172 ! CHARACTER(28) :: CVAR
173 ! INTEGER :: PHASE(2) = (/ 0, 0 /)
174 ! SAVE PHASE
175 !#endif
176 
177  enters("INTERFACE_CONDITION_ASSEMBLE_FEM",err,error,*999)
178 
179  IF(ASSOCIATED(interface_condition)) THEN
180  IF(ASSOCIATED(interface_condition%LAGRANGE)) THEN
181  lagrange_field=>interface_condition%LAGRANGE%LAGRANGE_FIELD
182  IF(ASSOCIATED(lagrange_field)) THEN
183  interface_equations=>interface_condition%INTERFACE_EQUATIONS
184  IF(ASSOCIATED(interface_equations)) THEN
185  interface_matrices=>interface_equations%INTERFACE_MATRICES
186  IF(ASSOCIATED(interface_matrices)) THEN
187  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
188  CALL cpu_timer(user_cpu,user_time1,err,error,*999)
189  CALL cpu_timer(system_cpu,system_time1,err,error,*999)
190  ENDIF
191  !Initialise the matrices and rhs vector
192 #ifdef TAUPROF
193  CALL tau_static_phase_start("INTERFACE_MATRICES_VALUES_INITIALISE()")
194 #endif
195  CALL interface_matrices_values_initialise(interface_matrices,0.0_dp,err,error,*999)
196 #ifdef TAUPROF
197  CALL tau_static_phase_stop("INTERFACE_MATRICES_VALUES_INITIALISE()")
198 #endif
199  !Assemble the elements
200  !Allocate the element matrices
201 #ifdef TAUPROF
202  CALL tau_static_phase_start("InterfaceMatrices_ElementInitialise()")
203 #endif
204  CALL interfacematrices_elementinitialise(interface_matrices,err,error,*999)
205  elements_mapping=>lagrange_field%DECOMPOSITION%DOMAIN(lagrange_field%DECOMPOSITION%MESH_COMPONENT_NUMBER)%PTR% &
206  & mappings%ELEMENTS
207 #ifdef TAUPROF
208  CALL tau_static_phase_stop("InterfaceMatrices_ElementInitialise()")
209 #endif
210  !Output timing information if required
211  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
212  CALL cpu_timer(user_cpu,user_time2,err,error,*999)
213  CALL cpu_timer(system_cpu,system_time2,err,error,*999)
214  user_elapsed=user_time2(1)-user_time1(1)
215  system_elapsed=system_time2(1)-system_time1(1)
216  CALL write_string_value(general_output_type,"User time for interface equations setup and initialisation = ", &
217  & user_elapsed,err,error,*999)
218  CALL write_string_value(general_output_type,"System time for interface equations setup and initialisation = ", &
219  & system_elapsed,err,error,*999)
220  element_user_elapsed=0.0_sp
221  element_system_elapsed=0.0_sp
222  ENDIF
223  number_of_times=0
224  !Loop over the internal elements
225 
226 #ifdef TAUPROF
227  CALL tau_static_phase_start("Internal Elements Loop")
228 #endif
229  DO element_idx=elements_mapping%INTERNAL_START,elements_mapping%INTERNAL_FINISH
230 !#ifdef TAUPROF
231 ! WRITE (CVAR,'(a23,i3)') 'Internal Elements Loop ',element_idx
232 ! CALL TAU_PHASE_CREATE_DYNAMIC(PHASE,CVAR)
233 ! CALL TAU_PHASE_START(PHASE)
234 !#endif
235  ne=elements_mapping%DOMAIN_LIST(element_idx)
236  number_of_times=number_of_times+1
237  CALL interfacematrices_elementcalculate(interface_matrices,ne,err,error,*999)
238  CALL interfacecondition_finiteelementcalculate(interface_condition,ne,err,error,*999)
239  CALL interface_matrices_element_add(interface_matrices,err,error,*999)
240 !#ifdef TAUPROF
241 ! CALL TAU_PHASE_STOP(PHASE)
242 !#endif
243  ENDDO !element_idx
244 #ifdef TAUPROF
245  CALL tau_static_phase_stop("Internal Elements Loop")
246 #endif
247 
248  !Output timing information if required
249  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
250  CALL cpu_timer(user_cpu,user_time3,err,error,*999)
251  CALL cpu_timer(system_cpu,system_time3,err,error,*999)
252  user_elapsed=user_time3(1)-user_time2(1)
253  system_elapsed=system_time3(1)-system_time2(1)
254  element_user_elapsed=user_elapsed
255  element_system_elapsed=system_elapsed
256  CALL write_string_value(general_output_type,"User time for internal interface equations assembly = ", &
257  & user_elapsed, err,error,*999)
258  CALL write_string_value(general_output_type,"System time for internal interface equations assembly = ", &
259  & system_elapsed,err,error,*999)
260  ENDIF
261  !Output timing information if required
262  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
263  CALL cpu_timer(user_cpu,user_time4,err,error,*999)
264  CALL cpu_timer(system_cpu,system_time4,err,error,*999)
265  user_elapsed=user_time4(1)-user_time3(1)
266  system_elapsed=system_time4(1)-system_time3(1)
267  CALL write_string_value(general_output_type,"User time for parameter transfer completion = ",user_elapsed, &
268  & err,error,*999)
269  CALL write_string_value(general_output_type,"System time for parameter transfer completion = ",system_elapsed, &
270  & err,error,*999)
271  ENDIF
272  !Loop over the boundary and ghost elements
273 #ifdef TAUPROF
274  CALL tau_static_phase_start("Boundary and Ghost Elements Loop")
275 #endif
276  DO element_idx=elements_mapping%BOUNDARY_START,elements_mapping%GHOST_FINISH
277  ne=elements_mapping%DOMAIN_LIST(element_idx)
278  number_of_times=number_of_times+1
279  CALL interfacematrices_elementcalculate(interface_matrices,ne,err,error,*999)
280  CALL interfacecondition_finiteelementcalculate(interface_condition,ne,err,error,*999)
281  CALL interface_matrices_element_add(interface_matrices,err,error,*999)
282  ENDDO !element_idx
283 #ifdef TAUPROF
284  CALL tau_static_phase_stop("Boundary and Ghost Elements Loop")
285 #endif
286  !Output timing information if required
287  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
288  CALL cpu_timer(user_cpu,user_time5,err,error,*999)
289  CALL cpu_timer(system_cpu,system_time5,err,error,*999)
290  user_elapsed=user_time5(1)-user_time4(1)
291  system_elapsed=system_time5(1)-system_time4(1)
292  element_user_elapsed=element_user_elapsed+user_elapsed
293  element_system_elapsed=element_system_elapsed+user_elapsed
294  CALL write_string_value(general_output_type,"User time for boundary+ghost equations assembly = ",user_elapsed, &
295  & err,error,*999)
296  CALL write_string_value(general_output_type,"System time for boundary+ghost equations assembly = ",system_elapsed, &
297  & err,error,*999)
298  IF(number_of_times>0) THEN
299  CALL write_string_value(general_output_type,"Average element user time for equations assembly = ", &
300  & element_user_elapsed/number_of_times,err,error,*999)
301  CALL write_string_value(general_output_type,"Average element system time for equations assembly = ", &
302  & element_system_elapsed/number_of_times,err,error,*999)
303  ENDIF
304  ENDIF
305  !Finalise the element matrices
306 #ifdef TAUPROF
307  CALL tau_static_phase_start("INTERFACE_MATRICES_ELEMENT_FINALISE()")
308 #endif
309  CALL interface_matrices_element_finalise(interface_matrices,err,error,*999)
310 #ifdef TAUPROF
311  CALL tau_static_phase_stop("INTERFACE_MATRICES_ELEMENT_FINALISE()")
312 #endif
313  !Output equations matrices and vector if required
314  IF(interface_equations%OUTPUT_TYPE>=interface_equations_matrix_output) THEN
315  CALL interface_matrices_output(general_output_type,interface_matrices,err,error,*999)
316  ENDIF
317  !Output timing information if required
318  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
319  CALL cpu_timer(user_cpu,user_time6,err,error,*999)
320  CALL cpu_timer(system_cpu,system_time6,err,error,*999)
321  user_elapsed=user_time6(1)-user_time1(1)
322  system_elapsed=system_time6(1)-system_time1(1)
323  CALL write_string(general_output_type,"***",err,error,*999)
324  CALL write_string_value(general_output_type,"Total user time for equations assembly = ",user_elapsed, &
325  & err,error,*999)
326  CALL write_string_value(general_output_type,"Total system time for equations assembly = ",system_elapsed, &
327  & err,error,*999)
328  CALL write_string(general_output_type,"***",err,error,*999)
329  ENDIF
330  ELSE
331  CALL flagerror("Interface matrices is not associated.",err,error,*999)
332  ENDIF
333  ELSE
334  CALL flagerror("Interface matrices is not associated.",err,error,*999)
335  ENDIF
336  ELSE
337  CALL flagerror("Lagrange field is not associated.",err,error,*999)
338  ENDIF
339  ELSE
340  CALL flagerror("Interface condition Lagrange is not associated.",err,error,*999)
341  ENDIF
342  ELSE
343  CALL flagerror("Interface condition is not associated",err,error,*999)
344  ENDIF
345 
346  exits("INTERFACE_CONDITION_ASSEMBLE_FEM")
347  RETURN
348 999 errorsexits("INTERFACE_CONDITION_ASSEMBLE_FEM",err,error)
349  RETURN 1
350  END SUBROUTINE interface_condition_assemble_fem
351 
352  !
353  !==================================================================================================================================
354  !
355 
357  SUBROUTINE interface_condition_create_finish(INTERFACE_CONDITION,ERR,ERROR,*)
358 
359  !Argument variables
360  TYPE(interface_condition_type), POINTER :: interface_condition
361  INTEGER(INTG), INTENT(OUT) :: err
362  TYPE(varying_string), INTENT(OUT) :: error
363  !Local Variables
364  INTEGER(INTG) :: mesh_idx,mesh_idx_count,number_of_components,variable_idx
365  INTEGER(INTG), POINTER :: new_variable_mesh_indices(:)
366  TYPE(field_variable_type), POINTER :: field_variable
367  TYPE(field_variable_ptr_type), POINTER :: new_field_variables(:)
368  TYPE(interface_type), POINTER :: interface
369  TYPE(interface_dependent_type), POINTER :: interface_dependent
370  TYPE(varying_string) :: local_error
371 
372  NULLIFY(new_field_variables)
373  NULLIFY(new_variable_mesh_indices)
374 
375  enters("INTERFACE_CONDITION_CREATE_FINISH",err,error,*999)
376 
377  IF(ASSOCIATED(interface_condition)) THEN
378  IF(interface_condition%INTERFACE_CONDITION_FINISHED) THEN
379  CALL flagerror("Interface condition has already been finished.",err,error,*999)
380  ELSE
381  interface=>interface_condition%INTERFACE
382  IF(ASSOCIATED(interface)) THEN
383  !Test various inputs have been set up.
384  SELECT CASE(interface_condition%METHOD)
386  interface_dependent=>interface_condition%DEPENDENT
387  IF(ASSOCIATED(interface_dependent)) THEN
388  !Check the dependent field variables have been set.
389  IF(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES<2) THEN
390  local_error="The number of added dependent variables of "// &
391  & trim(number_to_vstring(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES,"*",err,error))// &
392  & " is invalid. The number must be >= 2."
393  CALL flagerror(local_error,err,error,*999)
394  ENDIF
395 
396  !\todo check if interface mesh connectivity basis has same number of gauss points as interface geometric field IF(INTERFACE_CONDITION%INTERFACE%MESH_CONNECTIVITY%BASIS%QUADRATURE%NUMBER_OF_GAUSS_XI/=)
397 
398  !Note There is no need to check that the dependent variables have the same number of components.
399  !The user will need to set a fixed BC on the interface dof relating to the field components
400  !not present in each of the coupled bodies, eliminating this dof from the solver matrices
401  SELECT CASE(interface_condition%OPERATOR)
404  !Check that the dependent variables have the same number of components
405  field_variable=>interface_dependent%FIELD_VARIABLES(1)%PTR
406  IF(ASSOCIATED(field_variable)) THEN
407  number_of_components=field_variable%NUMBER_OF_COMPONENTS
408  DO variable_idx=2,interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES
409  field_variable=>interface_dependent%FIELD_VARIABLES(variable_idx)%PTR
410  IF(ASSOCIATED(field_variable)) THEN
411  !do nothing
412  ELSE
413  local_error="The interface condition field variables is not associated for variable index "// &
414  & trim(number_to_vstring(variable_idx,"*",err,error))
415  CALL flagerror(local_error,err,error,*999)
416  ENDIF
417  ENDDO !variable_idx
418  ELSE
419  CALL flagerror("Interface field variable is not associated.",err,error,*999)
420  ENDIF
422  CALL flagerror("Not implemented.",err,error,*999)
424  CALL flagerror("Not implemented.",err,error,*999)
425  CASE DEFAULT
426  local_error="The interface condition operator of "// &
427  & trim(number_to_vstring(interface_condition%OPERATOR,"*",err,error))//" is invalid."
428  CALL flagerror(local_error,err,error,*999)
429  END SELECT
430 
431  !Reorder the dependent variables based on mesh index order
432  ALLOCATE(new_field_variables(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES),stat=err)
433  IF(err/=0) CALL flagerror("Could not allocate new field variables.",err,error,*999)
434  ALLOCATE(new_variable_mesh_indices(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES),stat=err)
435  IF(err/=0) CALL flagerror("Could not allocate new variable mesh indices.",err,error,*999)
436  new_variable_mesh_indices=0
437  mesh_idx_count=0
438  DO mesh_idx=1,interface%NUMBER_OF_COUPLED_MESHES
439  DO variable_idx=1,interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES
440  IF(interface_dependent%VARIABLE_MESH_INDICES(variable_idx)==mesh_idx) THEN
441  mesh_idx_count=mesh_idx_count+1
442  new_field_variables(mesh_idx_count)%PTR=>interface_dependent%FIELD_VARIABLES(variable_idx)%PTR
443  new_variable_mesh_indices(mesh_idx_count)=mesh_idx
444  ENDIF
445  ENDDO !variable_idx
446  ENDDO !mesh_idx
447  IF(mesh_idx_count/=interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES) &
448  & CALL flagerror("Invalid dependent variable mesh index setup.",err,error,*999)
449  IF(ASSOCIATED(interface_dependent%FIELD_VARIABLES)) DEALLOCATE(interface_dependent%FIELD_VARIABLES)
450  IF(ASSOCIATED(interface_dependent%VARIABLE_MESH_INDICES)) DEALLOCATE(interface_dependent%VARIABLE_MESH_INDICES)
451  interface_dependent%FIELD_VARIABLES=>new_field_variables
452  interface_dependent%VARIABLE_MESH_INDICES=>new_variable_mesh_indices
453  ELSE
454  CALL flagerror("Interface condition dependent is not associated.",err,error,*999)
455  ENDIF
457  CALL flagerror("Not implemented.",err,error,*999)
459  CALL flagerror("Not implemented.",err,error,*999)
460  CASE DEFAULT
461  local_error="The interface condition method of "//trim(number_to_vstring(interface_condition%METHOD,"*",err,error))// &
462  & " is invalid."
463  CALL flagerror(local_error,err,error,*999)
464  END SELECT
465  !Finish the interface condition creation
466  interface_condition%INTERFACE_CONDITION_FINISHED=.true.
467  ELSE
468  CALL flagerror("Interface condition interface is not associated.",err,error,*999)
469  ENDIF
470  ENDIF
471  ELSE
472  CALL flagerror("Interface condition is not associated.",err,error,*999)
473  ENDIF
474 
475  exits("INTERFACE_CONDITION_CREATE_FINISH")
476  RETURN
477 999 IF(ASSOCIATED(new_field_variables)) DEALLOCATE(new_field_variables)
478  IF(ASSOCIATED(new_variable_mesh_indices)) DEALLOCATE(new_variable_mesh_indices)
479  errorsexits("INTERFACE_CONDITION_CREATE_FINISH",err,error)
480  RETURN 1
481 
482  END SUBROUTINE interface_condition_create_finish
483 
484  !
485  !================================================================================================================================
486  !
487 
489  SUBROUTINE interface_condition_create_start(USER_NUMBER,INTERFACE,GEOMETRIC_FIELD,INTERFACE_CONDITION,ERR,ERROR,*)
490 
491  !Argument variables
492  INTEGER(INTG), INTENT(IN) :: user_number
493  TYPE(interface_type), POINTER :: interface
494  TYPE(field_type), POINTER :: geometric_field
495  TYPE(interface_condition_type), POINTER :: interface_condition
496  INTEGER(INTG), INTENT(OUT) :: err
497  TYPE(varying_string), INTENT(OUT) :: error
498  !Local Variables
499  INTEGER(INTG) :: dummy_err,interface_conditions_idx
500  TYPE(interface_type), POINTER :: geometric_interface
501  TYPE(interface_condition_type), POINTER :: new_interface_condition
502  TYPE(interface_condition_ptr_type), POINTER :: new_interface_conditions(:)
503  TYPE(region_type), POINTER :: geometric_region,geometric_interface_parent_region,interface_parent_region
504  TYPE(varying_string) :: dummy_error,local_error
505 
506  NULLIFY(new_interface_condition)
507  NULLIFY(new_interface_conditions)
508 
509  enters("INTERFACE_CONDITION_CREATE_START",err,error,*997)
510 
511  IF(ASSOCIATED(interface)) THEN
512  IF(ASSOCIATED(interface%INTERFACE_CONDITIONS)) THEN
513  CALL interface_condition_user_number_find(user_number,interface,new_interface_condition,err,error,*997)
514  IF(ASSOCIATED(new_interface_condition)) THEN
515  local_error="Interface condition user number "//trim(number_to_vstring(user_number,"*",err,error))// &
516  & " has already been created on interface number "//trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//"."
517  CALL flagerror(local_error,err,error,*997)
518  ELSE
519  IF(ASSOCIATED(geometric_field)) THEN
520  IF(geometric_field%FIELD_FINISHED) THEN
521  !Check the geometric field is defined on the interface
522  geometric_interface=>geometric_field%INTERFACE
523  IF(ASSOCIATED(geometric_interface)) THEN
524  IF(ASSOCIATED(geometric_interface,interface)) THEN
525  NULLIFY(new_interface_condition)
526  !Initialise the new interface condition
527  CALL interface_condition_initialise(new_interface_condition,err,error,*999)
528  !Set default interface condition values
529  new_interface_condition%USER_NUMBER=user_number
530  new_interface_condition%GLOBAL_NUMBER=interface%INTERFACE_CONDITIONS%NUMBER_OF_INTERFACE_CONDITIONS+1
531  new_interface_condition%INTERFACE_CONDITIONS=>interface%INTERFACE_CONDITIONS
532  new_interface_condition%INTERFACE=>INTERFACE
533  !Default attributes
534  new_interface_condition%GEOMETRY%GEOMETRIC_FIELD=>geometric_field
535  new_interface_condition%METHOD=interface_condition_lagrange_multipliers_method
536  new_interface_condition%OPERATOR=interface_condition_field_continuity_operator
537  IF(ASSOCIATED(interface%pointsConnectivity)) THEN
538  new_interface_condition%integrationType=interface_condition_data_points_integration
539  ELSE
540  new_interface_condition%integrationType=interface_condition_gauss_integration
541  ENDIF
542  CALL interface_condition_dependent_initialise(new_interface_condition,err,error,*999)
543  !Add new interface condition into list of interface conditions in the interface
544  ALLOCATE(new_interface_conditions(interface%INTERFACE_CONDITIONS%NUMBER_OF_INTERFACE_CONDITIONS+1),stat=err)
545  IF(err/=0) CALL flagerror("Could not allocate new interface conditions.",err,error,*999)
546  DO interface_conditions_idx=1,interface%INTERFACE_CONDITIONS%NUMBER_OF_INTERFACE_CONDITIONS
547  new_interface_conditions(interface_conditions_idx)%PTR=>interface%INTERFACE_CONDITIONS% &
548  & interface_conditions(interface_conditions_idx)%PTR
549  ENDDO !interface_conditions_idx
550  new_interface_conditions(interface%INTERFACE_CONDITIONS%NUMBER_OF_INTERFACE_CONDITIONS+1)%PTR=> &
551  & new_interface_condition
552  IF(ASSOCIATED(interface%INTERFACE_CONDITIONS%INTERFACE_CONDITIONS)) DEALLOCATE(interface%INTERFACE_CONDITIONS% &
553  & interface_conditions)
554  interface%INTERFACE_CONDITIONS%INTERFACE_CONDITIONS=>new_interface_conditions
555  interface%INTERFACE_CONDITIONS%NUMBER_OF_INTERFACE_CONDITIONS=interface%INTERFACE_CONDITIONS% &
556  number_of_interface_conditions+1
557  !Return the pointer
558  interface_condition=>new_interface_condition
559  ELSE
560  interface_parent_region=>interface%PARENT_REGION
561  IF(ASSOCIATED(interface_parent_region)) THEN
562  geometric_interface_parent_region=>geometric_interface%PARENT_REGION
563  IF(ASSOCIATED(geometric_interface_parent_region)) THEN
564  local_error="Geometric field interface does not match specified interface. "// &
565  "The geometric field was created on interface number "// &
566  & trim(number_to_vstring(geometric_interface%USER_NUMBER,"*",err,error))// &
567  & " of parent region number "// &
568  & trim(number_to_vstring(geometric_interface_parent_region%USER_NUMBER,"*",err,error))// &
569  & " and the specified interface was created as number "// &
570  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//" on parent region number "// &
571  & trim(number_to_vstring(interface_parent_region%USER_NUMBER,"*",err,error))//"."
572  CALL flagerror(local_error,err,error,*999)
573  ELSE
574  CALL flagerror("Geometric interface parent region is not associated.",err,error,*999)
575  ENDIF
576  ELSE
577  CALL flagerror("Interface parent region is not associated.",err,error,*999)
578  ENDIF
579  ENDIF
580  ELSE
581  geometric_region=>geometric_field%REGION
582  IF(ASSOCIATED(geometric_region)) THEN
583  local_error="The geometric field was created on region number "// &
584  & trim(number_to_vstring(geometric_region%USER_NUMBER,"*",err,error))// &
585  & " and not on the specified interface."
586  CALL flagerror(local_error,err,error,*999)
587  ELSE
588  CALL flagerror("The geometric field does not have a region or interface created.",err,error,*999)
589  ENDIF
590  ENDIF
591  ELSE
592  CALL flagerror("Geometric field has not been finished.",err,error,*999)
593  ENDIF
594  ELSE
595  CALL flagerror("Geometric field is not finished.",err,error,*999)
596  ENDIF
597  ENDIF
598  ELSE
599  local_error="The interface conditions on interface number "// &
600  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//" are not associated."
601  CALL flagerror(local_error,err,error,*997)
602  ENDIF
603  ELSE
604  CALL flagerror("Interface is not associated.",err,error,*997)
605  ENDIF
606 
607  exits("INTERFACE_CONDITION_CREATE_START")
608  RETURN
609 999 IF(ASSOCIATED(new_interface_condition)) CALL interface_condition_finalise(new_interface_condition,dummy_err,dummy_error,*998)
610 998 IF(ASSOCIATED(new_interface_conditions)) DEALLOCATE(new_interface_conditions)
611 997 errorsexits("INTERFACE_CONDITION_CREATE_START",err,error)
612  RETURN 1
613  END SUBROUTINE interface_condition_create_start
614 
615  !
616  !================================================================================================================================
617  !
618 
620  SUBROUTINE interface_condition_dependent_finalise(INTERFACE_DEPENDENT,ERR,ERROR,*)
621 
622  !Argument variables
623  TYPE(interface_dependent_type), POINTER :: interface_dependent
624  INTEGER(INTG), INTENT(OUT) :: err
625  TYPE(varying_string), INTENT(OUT) :: error
626  !Local Variables
627 
628  enters("INTERFACE_CONDITION_DEPENDENT_FINALISE",err,error,*999)
629 
630  IF(ASSOCIATED(interface_dependent)) THEN
631  IF(ASSOCIATED(interface_dependent%EQUATIONS_SETS)) DEALLOCATE(interface_dependent%EQUATIONS_SETS)
632  IF(ASSOCIATED(interface_dependent%FIELD_VARIABLES)) DEALLOCATE(interface_dependent%FIELD_VARIABLES)
633  IF(ASSOCIATED(interface_dependent%VARIABLE_MESH_INDICES)) DEALLOCATE(interface_dependent%VARIABLE_MESH_INDICES)
634  DEALLOCATE(interface_dependent)
635  ENDIF
636 
637  exits("INTERFACE_CONDITION_DEPENDENT_FINALISE")
638  RETURN
639 999 errorsexits("INTERFACE_CONDITION_DEPENDENT_FINALISE",err,error)
640  RETURN 1
641  END SUBROUTINE interface_condition_dependent_finalise
642 
643  !
644  !================================================================================================================================
645  !
646 
648  SUBROUTINE interface_condition_dependent_initialise(INTERFACE_CONDITION,ERR,ERROR,*)
649 
650  !Argument variables
651  TYPE(interface_condition_type), POINTER :: interface_condition
652  INTEGER(INTG), INTENT(OUT) :: err
653  TYPE(varying_string), INTENT(OUT) :: error
654  !Local Variables
655  INTEGER(INTG) :: dummy_err
656  TYPE(varying_string) :: dummy_error
657 
658  enters("INTERFACE_CONDITION_DEPENDENT_INITIALISE",err,error,*998)
659 
660  IF(ASSOCIATED(interface_condition)) THEN
661  IF(ASSOCIATED(interface_condition%DEPENDENT)) THEN
662  CALL flagerror("Interface condition dependent is already associated.",err,error,*999)
663  ELSE
664  ALLOCATE(interface_condition%DEPENDENT,stat=err)
665  IF(err/=0) CALL flagerror("Could not allocate interface condition dependent.",err,error,*999)
666  interface_condition%DEPENDENT%INTERFACE_CONDITION=>interface_condition
667  interface_condition%DEPENDENT%NUMBER_OF_DEPENDENT_VARIABLES=0
668  NULLIFY(interface_condition%DEPENDENT%EQUATIONS_SETS)
669  NULLIFY(interface_condition%DEPENDENT%FIELD_VARIABLES)
670  NULLIFY(interface_condition%DEPENDENT%VARIABLE_MESH_INDICES)
671  ENDIF
672  ELSE
673  CALL flagerror("Interface condition is not associated.",err,error,*999)
674  ENDIF
675 
676  exits("INTERFACE_CONDITION_DEPENDENT_INITIALISE")
677  RETURN
678 999 CALL interface_condition_dependent_finalise(interface_condition%DEPENDENT,dummy_err,dummy_error,*998)
679 998 errorsexits("INTERFACE_CONDITION_DEPENDENT_INITIALISE",err,error)
680  RETURN 1
681  END SUBROUTINE interface_condition_dependent_initialise
682 
683  !
684  !================================================================================================================================
685  !
686 
688  SUBROUTINE interface_condition_dependent_variable_add(INTERFACE_CONDITION,MESH_INDEX,EQUATIONS_SET,VARIABLE_TYPE,ERR,ERROR,*)
689 
690  !Argument variables
691  TYPE(interface_condition_type), POINTER :: interface_condition
692  INTEGER(INTG), INTENT(IN) :: mesh_index
693  TYPE(equations_set_type), POINTER :: equations_set
694  INTEGER(INTG), INTENT(IN) :: variable_type
695  INTEGER(INTG), INTENT(OUT) :: err
696  TYPE(varying_string), INTENT(OUT) :: error
697  !Local Variables
698  INTEGER(INTG) :: variable_idx
699  INTEGER(INTG), POINTER :: new_variable_mesh_indices(:)
700  LOGICAL :: found_mesh_index
701  TYPE(decomposition_type), POINTER :: decomposition
702  TYPE(equations_set_ptr_type), POINTER :: new_equations_sets(:)
703  TYPE(field_type), POINTER :: dependent_field
704  TYPE(field_variable_type), POINTER :: field_variable,interface_variable
705  TYPE(field_variable_ptr_type), POINTER :: new_field_variables(:)
706  TYPE(interface_type), POINTER :: interface
707  TYPE(interface_dependent_type), POINTER :: interface_dependent
708  TYPE(mesh_type), POINTER :: dependent_mesh,interface_mesh
709  TYPE(varying_string) :: local_error
710 
711  enters("INTERFACE_CONDITION_DEPENDENT_VARIABLE_ADD",err,error,*999)
712 
713  IF(ASSOCIATED(interface_condition)) THEN
714  interface_dependent=>interface_condition%DEPENDENT
715  IF(ASSOCIATED(interface_dependent)) THEN
716  interface=>interface_condition%INTERFACE
717  IF(ASSOCIATED(interface)) THEN
718  IF(mesh_index>0.AND.mesh_index<=interface%NUMBER_OF_COUPLED_MESHES) THEN
719  IF(ASSOCIATED(equations_set)) THEN
720  dependent_field=>equations_set%DEPENDENT%DEPENDENT_FIELD
721  IF(ASSOCIATED(dependent_field)) THEN
722  IF(variable_type>=1.AND.variable_type<=field_number_of_variable_types) THEN
723  field_variable=>dependent_field%VARIABLE_TYPE_MAP(variable_type)%PTR
724  IF(ASSOCIATED(field_variable)) THEN
725  !Check that the field variable hasn't already been added.
726  variable_idx=1
727  NULLIFY(interface_variable)
728  DO WHILE(variable_idx<=interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES.AND. &
729  & .NOT.ASSOCIATED(interface_variable))
730  IF(ASSOCIATED(field_variable,interface_dependent%FIELD_VARIABLES(variable_idx)%PTR)) THEN
731  interface_variable=>interface_dependent%FIELD_VARIABLES(variable_idx)%PTR
732  ELSE
733  variable_idx=variable_idx+1
734  ENDIF
735  ENDDO
736  IF(ASSOCIATED(interface_variable)) THEN
737  !Check if we are dealing with the same mesh index.
738  IF(mesh_index/=interface_dependent%VARIABLE_MESH_INDICES(variable_idx)) THEN
739  local_error="The dependent variable has already been added to the interface condition at "// &
740  & "position index "//trim(number_to_vstring(variable_idx,"*",err,error))//"."
741  CALL flagerror(local_error,err,error,*999)
742  ENDIF
743  ELSE
744  !Check the dependent variable and the mesh index match.
745  interface_mesh=>interface%COUPLED_MESHES(mesh_index)%PTR
746  IF(ASSOCIATED(interface_mesh)) THEN
747  decomposition=>dependent_field%DECOMPOSITION
748  IF(ASSOCIATED(decomposition)) THEN
749  dependent_mesh=>decomposition%MESH
750  IF(ASSOCIATED(dependent_mesh)) THEN
751  IF(ASSOCIATED(interface_mesh,dependent_mesh)) THEN
752  !The meshes match. Check if the dependent variable has already been added for the mesh index.
753  found_mesh_index=.false.
754  DO variable_idx=1,interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES
755  IF(interface_dependent%VARIABLE_MESH_INDICES(variable_idx)==mesh_index) THEN
756  found_mesh_index=.true.
757  EXIT
758  ENDIF
759  ENDDO !variable_idx
760  IF(found_mesh_index) THEN
761  !The mesh index has already been added to replace the dependent variable with the specified variable
762  interface_dependent%FIELD_VARIABLES(variable_idx)%PTR=>dependent_field% &
763  & variable_type_map(variable_type)%PTR
764  ELSE
765  !The mesh index has not been found so add a new dependent variable.
766  ALLOCATE(new_equations_sets(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1),stat=err)
767  IF(err/=0) CALL flagerror("Could not allocate new equations sets.",err,error,*999)
768  ALLOCATE(new_field_variables(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1),stat=err)
769  IF(err/=0) CALL flagerror("Could not allocate new field variables.",err,error,*999)
770  ALLOCATE(new_variable_mesh_indices(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1),stat=err)
771  IF(err/=0) CALL flagerror("Could not allocate new variable mesh indices.",err,error,*999)
772  DO variable_idx=1,interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES
773  new_equations_sets(variable_idx)%PTR=>interface_dependent%EQUATIONS_SETS(variable_idx)%PTR
774  new_field_variables(variable_idx)%PTR=>interface_dependent%FIELD_VARIABLES(variable_idx)%PTR
775  new_variable_mesh_indices(variable_idx)=interface_dependent%VARIABLE_MESH_INDICES(variable_idx)
776  ENDDO !variable_idx
777  new_equations_sets(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1)%PTR=>equations_set
778  new_field_variables(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1)%PTR=>dependent_field% &
779  & variable_type_map(variable_type)%PTR
780  new_variable_mesh_indices(interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1)=mesh_index
781  IF(ASSOCIATED(interface_dependent%EQUATIONS_SETS)) DEALLOCATE(interface_dependent%EQUATIONS_SETS)
782  IF(ASSOCIATED(interface_dependent%FIELD_VARIABLES)) DEALLOCATE(interface_dependent%FIELD_VARIABLES)
783  IF(ASSOCIATED(interface_dependent%VARIABLE_MESH_INDICES)) &
784  & DEALLOCATE(interface_dependent%VARIABLE_MESH_INDICES)
785  interface_dependent%EQUATIONS_SETS=>new_equations_sets
786  interface_dependent%FIELD_VARIABLES=>new_field_variables
787  interface_dependent%VARIABLE_MESH_INDICES=>new_variable_mesh_indices
788  interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES= &
789  & interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1
790  ENDIF
791  ELSE
792  CALL flagerror("The dependent field mesh does not match the interface mesh.",err,error,*999)
793  ENDIF
794  ELSE
795  CALL flagerror("The dependent field decomposition mesh is not associated.",err,error,*999)
796  ENDIF
797  ELSE
798  CALL flagerror("The dependent field decomposition is not associated.",err,error,*999)
799  ENDIF
800  ELSE
801  local_error="The interface mesh for mesh index "//trim(number_to_vstring(mesh_index,"*",err,error))// &
802  & " is not associated."
803  CALL flagerror(local_error,err,error,*999)
804  ENDIF
805  ENDIF
806  ELSE
807  local_error="The field variable type of "//trim(number_to_vstring(variable_type,"*",err,error))// &
808  & " has not been created on field number "// &
809  & trim(number_to_vstring(dependent_field%USER_NUMBER,"*",err,error))//"."
810  CALL flagerror(local_error,err,error,*999)
811  ENDIF
812  ELSE
813  local_error="The field variable type of "//trim(number_to_vstring(variable_type,"*",err,error))// &
814  & " is invalid. The variable type must be between 1 and "// &
815  & trim(number_to_vstring(field_number_of_variable_types,"*",err,error))//"."
816  CALL flagerror(local_error,err,error,*999)
817  ENDIF
818  ELSE
819  CALL flagerror("Equations set dependent field is not associated.",err,error,*999)
820  ENDIF
821  ELSE
822  CALL flagerror("Equations set is not associated.",err,error,*999)
823  ENDIF
824  ELSE
825  local_error="The specificed mesh index of "//trim(number_to_vstring(mesh_index,"*",err,error))// &
826  & " is invalid. The mesh index must be > 0 and <= "// &
827  & trim(number_to_vstring(interface%NUMBER_OF_COUPLED_MESHES,"*",err,error))//"."
828  CALL flagerror(local_error,err,error,*999)
829  ENDIF
830  ELSE
831  CALL flagerror("Interface condition interface is not associated.",err,error,*999)
832  ENDIF
833  ELSE
834  CALL flagerror("Interface condition dependent is not associated.",err,error,*999)
835  ENDIF
836  ELSE
837  CALL flagerror("Interface conditions is not associated.",err,error,*999)
838  ENDIF
839 
840  exits("INTERFACE_CONDITION_DEPENDENT_VARIABLE_ADD")
841  RETURN
842 999 errorsexits("INTERFACE_CONDITION_DEPENDENT_VARIABLE_ADD",err,error)
843  RETURN 1
844  END SUBROUTINE interface_condition_dependent_variable_add
845 
846  !
847  !================================================================================================================================
848  !
849 
851  SUBROUTINE interface_condition_destroy(INTERFACE_CONDITION,ERR,ERROR,*)
852 
853  !Argument variables
854  TYPE(interface_condition_type), POINTER :: interface_condition
855  INTEGER(INTG), INTENT(OUT) :: err
856  TYPE(varying_string), INTENT(OUT) :: error
857  !Local Variables
858  INTEGER(INTG) :: interface_condition_idx,interface_condition_position
859  TYPE(interface_condition_ptr_type), POINTER :: new_interface_conditions(:)
860  TYPE(interface_conditions_type), POINTER :: interface_conditions
861 
862  NULLIFY(new_interface_conditions)
863 
864  enters("INTERFACE_CONDITION_DESTROY",err,error,*999)
865 
866  IF(ASSOCIATED(interface_condition)) THEN
867  interface_conditions=>interface_condition%INTERFACE_CONDITIONS
868  IF(ASSOCIATED(interface_conditions)) THEN
869  interface_condition_position=interface_condition%GLOBAL_NUMBER
870 
871  !Destroy all the interface condition components
872  CALL interface_condition_finalise(interface_condition,err,error,*999)
873 
874  !Remove the interface condition from the list of interface conditions
875  IF(interface_conditions%NUMBER_OF_INTERFACE_CONDITIONS>1) THEN
876  ALLOCATE(new_interface_conditions(interface_conditions%NUMBER_OF_INTERFACE_CONDITIONS-1),stat=err)
877  IF(err/=0) CALL flagerror("Could not allocate new interface conditions.",err,error,*999)
878  DO interface_condition_idx=1,interface_conditions%NUMBER_OF_INTERFACE_CONDITIONS
879  IF(interface_condition_idx<interface_condition_position) THEN
880  new_interface_conditions(interface_condition_idx)%PTR=>interface_conditions% &
881  & interface_conditions(interface_condition_idx)%PTR
882  ELSE IF(interface_condition_idx>interface_condition_position) THEN
883  interface_conditions%INTERFACE_CONDITIONS(interface_condition_idx)%PTR%GLOBAL_NUMBER=interface_conditions% &
884  & interface_conditions(interface_condition_idx)%PTR%GLOBAL_NUMBER-1
885  new_interface_conditions(interface_condition_idx-1)%PTR=>interface_conditions% &
886  & interface_conditions(interface_condition_idx)%PTR
887  ENDIF
888  ENDDO !interface_conditions_idx
889  IF(ASSOCIATED(interface_conditions%INTERFACE_CONDITIONS)) DEALLOCATE(interface_conditions%INTERFACE_CONDITIONS)
890  interface_conditions%INTERFACE_CONDITIONS=>new_interface_conditions
891  interface_conditions%NUMBER_OF_INTERFACE_CONDITIONS=interface_conditions%NUMBER_OF_INTERFACE_CONDITIONS-1
892  ELSE
893  DEALLOCATE(interface_conditions%INTERFACE_CONDITIONS)
894  interface_conditions%NUMBER_OF_INTERFACE_CONDITIONS=0
895  ENDIF
896 
897  ELSE
898  CALL flagerror("Interface conditions interface conditions is not associated.",err,error,*999)
899  ENDIF
900  ELSE
901  CALL flagerror("Interface conditions is not associated.",err,error,*998)
902  ENDIF
903 
904  exits("INTERFACE_CONDITIONS_DESTROY")
905  RETURN
906 999 IF(ASSOCIATED(new_interface_conditions)) DEALLOCATE(new_interface_conditions)
907 998 errorsexits("INTERFACE_CONDITION_DESTROY",err,error)
908  RETURN 1
909  END SUBROUTINE interface_condition_destroy
910 
911  !
912  !================================================================================================================================
913  !
914 
916  SUBROUTINE interface_condition_equations_create_finish(INTERFACE_CONDITION,ERR,ERROR,*)
917 
918  !Argument variables
919  TYPE(interface_condition_type), POINTER :: interface_condition
920  INTEGER(INTG), INTENT(OUT) :: err
921  TYPE(varying_string), INTENT(OUT) :: error
922  !Local Variables
923  INTEGER(INTG), ALLOCATABLE :: storage_type(:),structure_type(:)
924  LOGICAL, ALLOCATABLE :: matrices_transpose(:)
925  INTEGER(INTG) :: number_of_dependent_variables
926  TYPE(interface_dependent_type), POINTER :: interface_dependent
927  TYPE(interface_equations_type), POINTER :: interface_equations
928  TYPE(interface_mapping_type), POINTER :: interface_mapping
929  TYPE(interface_matrices_type), POINTER :: interface_matrices
930  TYPE(varying_string) :: local_error
931 
932  enters("INTERFACE_CONDITIONS_EQUATIONS_CREATE_FINISH",err,error,*999)
933 
934  IF(ASSOCIATED(interface_condition)) THEN
935  SELECT CASE(interface_condition%METHOD)
937  !Finish the interface equations creation
938  NULLIFY(interface_equations)
939  CALL interface_condition_equations_get(interface_condition,interface_equations,err,error,*999)
940  IF(interface_equations%INTERFACE_EQUATIONS_FINISHED) THEN
941  CALL flagerror("Interface condition equations have already been finished.",err,error,*999)
942  ELSE
943  CALL interface_equations_create_finish(interface_equations,err,error,*999)
944  interface_dependent=>interface_condition%DEPENDENT
945  IF(ASSOCIATED(interface_dependent)) THEN
946  !Create the interface mapping.
947  NULLIFY(interface_mapping)
948  CALL interface_mapping_create_start(interface_equations,interface_mapping,err,error,*999)
949  CALL interfacemapping_lagrangevariableset(interface_mapping,field_u_variable_type,err,error,*999)
950  SELECT CASE(interface_condition%METHOD)
952  number_of_dependent_variables=interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES
954  number_of_dependent_variables=interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES+1
955  ENDSELECT
956  CALL interface_mapping_matrices_number_set(interface_mapping,number_of_dependent_variables,err,error,*999)
957  ALLOCATE(matrices_transpose(number_of_dependent_variables),stat=err)
958  IF(err/=0) CALL flagerror("Could not allocate matrices transpose.",err,error,*999)
959  matrices_transpose=.true.
960  SELECT CASE(interface_condition%METHOD)
962  !Set the last interface matrix to have no transpose
963  matrices_transpose(number_of_dependent_variables)=.false.
964  ENDSELECT
965  CALL interface_mapping_matrices_transpose_set(interface_mapping,matrices_transpose,err,error,*999)
966  IF(ALLOCATED(matrices_transpose)) DEALLOCATE(matrices_transpose)
967  CALL interface_mapping_rhs_variable_type_set(interface_mapping,field_deludeln_variable_type,err,error,*999)
968  CALL interface_mapping_create_finish(interface_mapping,err,error,*999)
969  !Create the interface matrices
970  NULLIFY(interface_matrices)
971  CALL interface_matrices_create_start(interface_equations,interface_matrices,err,error,*999)
972  ALLOCATE(storage_type(interface_matrices%NUMBER_OF_INTERFACE_MATRICES),stat=err)
973  IF(err/=0) CALL flagerror("Could not allocate storage type.",err,error,*999)
974  SELECT CASE(interface_equations%SPARSITY_TYPE)
975  CASE(interface_matrices_full_matrices)
976  storage_type=matrix_block_storage_type
977  CALL interface_matrices_storage_type_set(interface_matrices,storage_type,err,error,*999)
978  CASE(interface_matrices_sparse_matrices)
979  ALLOCATE(structure_type(interface_matrices%NUMBER_OF_INTERFACE_MATRICES),stat=err)
980  IF(err/=0) CALL flagerror("Could not allocate structure type.",err,error,*999)
982  structure_type=interface_matrix_fem_structure
983  CALL interface_matrices_storage_type_set(interface_matrices,storage_type,err,error,*999)
984  CALL interface_matrices_structure_type_set(interface_matrices,structure_type,err,error,*999)
985  IF(ALLOCATED(structure_type)) DEALLOCATE(structure_type)
986  CASE DEFAULT
987  local_error="The interface equations sparsity type of "// &
988  & trim(number_to_vstring(interface_equations%SPARSITY_TYPE,"*",err,error))//" is invalid."
989  CALL flagerror(local_error,err,error,*999)
990  END SELECT
991  IF(ALLOCATED(storage_type)) DEALLOCATE(storage_type)
992  CALL interface_matrices_create_finish(interface_matrices,err,error,*999)
993  ELSE
994  CALL flagerror("Interface condition dependent is not associated.",err,error,*999)
995  ENDIF
996  ENDIF
998  CALL flagerror("Not implemented.",err,error,*999)
1000  CALL flagerror("Not implemented.",err,error,*999)
1001  CASE DEFAULT
1002  local_error="The interface condition method of "//trim(number_to_vstring(interface_condition%METHOD,"*",err,error))// &
1003  & " is invalid."
1004  CALL flagerror(local_error,err,error,*999)
1005  END SELECT
1006  ELSE
1007  CALL flagerror("Interface conditions is not associated.",err,error,*999)
1008  ENDIF
1009 
1010  exits("INTERFACE_CONDITION_EQUATIONS_CREATE_FINISH")
1011  RETURN
1012 999 IF(ALLOCATED(matrices_transpose)) DEALLOCATE(matrices_transpose)
1013  IF(ALLOCATED(storage_type)) DEALLOCATE(storage_type)
1014  IF(ALLOCATED(structure_type)) DEALLOCATE(structure_type)
1015  errorsexits("INTERFACE_CONDITION_EQUATIONS_CREATE_FINISH",err,error)
1016  RETURN 1
1017 
1018  END SUBROUTINE interface_condition_equations_create_finish
1019 
1020  !
1021  !================================================================================================================================
1022  !
1023 
1028  SUBROUTINE interface_condition_equations_create_start(INTERFACE_CONDITION,INTERFACE_EQUATIONS,ERR,ERROR,*)
1029 
1030  !Argument variables
1031  TYPE(interface_condition_type), POINTER :: interface_condition
1032  TYPE(interface_equations_type), POINTER :: interface_equations
1033  INTEGER(INTG), INTENT(OUT) :: err
1034  TYPE(varying_string), INTENT(OUT) :: error
1035  !Local Variables
1036  INTEGER(INTG) :: variable_idx
1037  TYPE(interface_dependent_type), POINTER :: interface_dependent
1038  TYPE(varying_string) :: local_error
1039 
1040  enters("INTERFACE_CONDITION_EQUATIONS_CREATE_START",err,error,*999)
1041 
1042  IF(ASSOCIATED(interface_condition)) THEN
1043  IF(ASSOCIATED(interface_equations)) THEN
1044  CALL flagerror("Interface equations is already associated.",err,error,*999)
1045  ELSE
1046  NULLIFY(interface_equations)
1047  SELECT CASE(interface_condition%METHOD)
1049  IF(ASSOCIATED(interface_condition%LAGRANGE)) THEN
1050  IF(interface_condition%LAGRANGE%LAGRANGE_FINISHED) THEN
1051  interface_dependent=>interface_condition%DEPENDENT
1052  IF(ASSOCIATED(interface_dependent)) THEN
1053  !Initialise the setup
1054  CALL interface_equations_create_start(interface_condition,interface_equations,err,error,*999)
1055  !Set the number of interpolation sets
1056  CALL interfaceequations_interfaceinterpsetsnumberset(interface_equations,1,1,1,err,error,*999)
1057  DO variable_idx=1,interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES
1058  CALL interfaceequations_variableinterpsetsnumberset(interface_equations,variable_idx,1,1,0, &
1059  & err,error,*999)
1060  ENDDO !variable_idx
1061  ELSE
1062  CALL flagerror("Interface condition dependent is not associated.",err,error,*999)
1063  ENDIF
1064  !Return the pointer
1065  interface_equations=>interface_condition%INTERFACE_EQUATIONS
1066  ELSE
1067  CALL flagerror("Interface condition Lagrange field has not been finished.",err,error,*999)
1068  ENDIF
1069  ELSE
1070  CALL flagerror("Interface condition Lagrange is not associated.",err,error,*999)
1071  ENDIF
1073  CALL flagerror("Not implemented.",err,error,*999)
1075  CALL flagerror("Not implemented.",err,error,*999)
1076  CASE DEFAULT
1077  local_error="The interface condition method of "//trim(number_to_vstring(interface_condition%METHOD,"*",err,error))// &
1078  & " is invalid."
1079  CALL flagerror(local_error,err,error,*999)
1080  END SELECT
1081  ENDIF
1082  ELSE
1083  CALL flagerror("Interface condition is not associated.",err,error,*999)
1084  ENDIF
1085 
1086  exits("INTERFACE_CONDITION_EQUATIONS_CREATE_START")
1087  RETURN
1088 999 errorsexits("INTERFACE_CONDITION_EQUATIONS_CREATE_START",err,error)
1089  RETURN 1
1090  END SUBROUTINE interface_condition_equations_create_start
1091 
1092  !
1093  !================================================================================================================================
1094  !
1095 
1097  SUBROUTINE interface_condition_equations_destroy(INTERFACE_CONDITION,ERR,ERROR,*)
1098 
1099  !Argument variables
1100  TYPE(interface_condition_type), POINTER :: interface_condition
1101  INTEGER(INTG), INTENT(OUT) :: err
1102  TYPE(varying_string), INTENT(OUT) :: error
1103  !Local Variables
1104 
1105  enters("INTERFACE_CONDITION_EQUATIONS_DESTROY",err,error,*999)
1106 
1107  IF(ASSOCIATED(interface_condition)) THEN
1108  IF(ASSOCIATED(interface_condition%INTERFACE_EQUATIONS)) THEN
1109  CALL interface_equations_destroy(interface_condition%INTERFACE_EQUATIONS,err,error,*999)
1110  ELSE
1111  CALL flagerror("Interface condition interface equations is not associated.",err,error,*999)
1112  ENDIF
1113  ELSE
1114  CALL flagerror("Interface condition is not associated.",err,error,*999)
1115  ENDIF
1116 
1117  exits("INTERFACE_CONDITION_EQUATIONS_DESTROY")
1118  RETURN
1119 999 errorsexits("INTERFACE_CONDITION_EQUATIONS_DESTROY",err,error)
1120  RETURN 1
1121  END SUBROUTINE interface_condition_equations_destroy
1122 
1123  !
1124  !================================================================================================================================
1125  !
1126 
1128  SUBROUTINE interface_condition_finalise(INTERFACE_CONDITION,ERR,ERROR,*)
1129 
1130  !Argument variables
1131  TYPE(interface_condition_type), POINTER :: interface_condition
1132  INTEGER(INTG), INTENT(OUT) :: err
1133  TYPE(varying_string), INTENT(OUT) :: error
1134  !Local Variables
1135 
1136  enters("INTERFACE_CONDITION_FINALISE",err,error,*999)
1137 
1138  IF(ASSOCIATED(interface_condition)) THEN
1139  CALL interface_condition_geometry_finalise(interface_condition%GEOMETRY,err,error,*999)
1140  CALL interface_condition_lagrange_finalise(interface_condition%LAGRANGE,err,error,*999)
1141  CALL interface_condition_penalty_finalise(interface_condition%PENALTY,err,error,*999)
1142  CALL interface_condition_dependent_finalise(interface_condition%DEPENDENT,err,error,*999)
1143  IF(ASSOCIATED(interface_condition%INTERFACE_EQUATIONS)) &
1144  & CALL interface_equations_destroy(interface_condition%INTERFACE_EQUATIONS,err,error,*999)
1145  DEALLOCATE(interface_condition)
1146  ENDIF
1147 
1148  exits("INTERFACE_CONDITION_FINALISE")
1149  RETURN
1150 999 errorsexits("INTERFACE_CONDITION_FINALISE",err,error)
1151  RETURN 1
1152  END SUBROUTINE interface_condition_finalise
1153 
1154 !
1155  !================================================================================================================================
1156  !
1157 
1159  SUBROUTINE interfacecondition_integrationtypeget(interfaceCondition,interfaceConditionIntegrationType,err,error,*)
1160 
1161  !Argument variables
1162  TYPE(interface_condition_type), POINTER :: interfacecondition
1163  INTEGER(INTG), INTENT(OUT) :: interfaceconditionintegrationtype
1164  INTEGER(INTG), INTENT(OUT) :: err
1165  TYPE(varying_string), INTENT(OUT) :: error
1166  !Local Variables
1167 
1168  enters("InterfaceCondition_IntegrationTypeGet",err,error,*999)
1169 
1170  IF(ASSOCIATED(interfacecondition)) THEN
1171  IF(interfacecondition%INTERFACE_CONDITION_FINISHED) THEN
1172  interfaceconditionintegrationtype=interfacecondition%integrationType
1173  ELSE
1174  CALL flagerror("Interface condition has not been finished.",err,error,*999)
1175  ENDIF
1176  ELSE
1177  CALL flagerror("Interface condition is not associated.",err,error,*999)
1178  ENDIF
1179 
1180  exits("InterfaceCondition_IntegrationTypeGet")
1181  RETURN
1182 999 errorsexits("InterfaceCondition_IntegrationTypeGet",err,error)
1183  RETURN 1
1184  END SUBROUTINE interfacecondition_integrationtypeget
1185 
1186  !
1187  !================================================================================================================================
1188  !
1189 
1191  SUBROUTINE interfacecondition_integrationtypeset(interfaceCondition,interfaceConditionIntegrationType,err,error,*)
1192 
1193  !Argument variables
1194  TYPE(interface_condition_type), POINTER :: interfacecondition
1195  INTEGER(INTG), INTENT(IN) :: interfaceconditionintegrationtype
1196  INTEGER(INTG), INTENT(OUT) :: err
1197  TYPE(varying_string), INTENT(OUT) :: error
1198  !Local Variables
1199  TYPE(varying_string) :: localerror
1200 
1201  enters("InterfaceCondition_IntegrationTypeSet",err,error,*999)
1202 
1203  IF(ASSOCIATED(interfacecondition)) THEN
1204  IF(interfacecondition%INTERFACE_CONDITION_FINISHED) THEN
1205  CALL flagerror("Interface condition has been finished.",err,error,*999)
1206  ELSE
1207  SELECT CASE(interfaceconditionintegrationtype)
1209  interfacecondition%integrationType=interface_condition_gauss_integration
1211  interfacecondition%integrationType=interface_condition_data_points_integration
1212  CASE DEFAULT
1213  localerror="The specified interface condition operator of "// &
1214  & trim(number_to_vstring(interfaceconditionintegrationtype,"*",err,error))//" is not valid."
1215  CALL flagerror(localerror,err,error,*999)
1216  END SELECT
1217  ENDIF
1218  ELSE
1219  CALL flagerror("Interface condition is not associated.",err,error,*999)
1220  ENDIF
1221 
1222  exits("InterfaceCondition_IntegrationTypeSet")
1223  RETURN
1224 999 errorsexits("InterfaceCondition_IntegrationTypeSet",err,error)
1225  RETURN 1
1226  END SUBROUTINE interfacecondition_integrationtypeset
1227 
1228 
1229  !
1230  !================================================================================================================================
1231  !
1232 
1234  SUBROUTINE interface_condition_geometry_finalise(INTERFACE_GEOMETRY,ERR,ERROR,*)
1235 
1236  !Argument variables
1237  TYPE(interface_geometry_type) :: interface_geometry
1238  INTEGER(INTG), INTENT(OUT) :: err
1239  TYPE(varying_string), INTENT(OUT) :: error
1240  !Local Variables
1241 
1242  enters("INTERFACE_CONDITION_GEOMETRY_FINALISE",err,error,*999)
1243 
1244  NULLIFY(interface_geometry%INTERFACE_CONDITION)
1245  NULLIFY(interface_geometry%GEOMETRIC_FIELD)
1246 
1247  exits("INTERFACE_CONDITION_GEOMETRY_FINALISE")
1248  RETURN
1249 999 errorsexits("INTERFACE_CONDITION_GEOMETRY_FINALISE",err,error)
1250  RETURN 1
1251  END SUBROUTINE interface_condition_geometry_finalise
1252 
1253  !
1254  !================================================================================================================================
1255  !
1256 
1258  SUBROUTINE interface_condition_geometry_initialise(INTERFACE_CONDITION,ERR,ERROR,*)
1259 
1260  !Argument variables
1261  TYPE(interface_condition_type), POINTER :: interface_condition
1262  INTEGER(INTG), INTENT(OUT) :: err
1263  TYPE(varying_string), INTENT(OUT) :: error
1264  !Local Variables
1265  INTEGER(INTG) :: dummy_err
1266  TYPE(varying_string) :: dummy_error
1267 
1268  enters("INTERFACE_CONDITION_GEOMETRY_INITIALISE",err,error,*998)
1269 
1270  IF(ASSOCIATED(interface_condition)) THEN
1271  interface_condition%GEOMETRY%INTERFACE_CONDITION=>interface_condition
1272  NULLIFY(interface_condition%GEOMETRY%GEOMETRIC_FIELD)
1273  ELSE
1274  CALL flagerror("Interface condition is not associated.",err,error,*999)
1275  ENDIF
1276 
1277  exits("INTERFACE_CONDITION_GEOMETRY_INITIALISE")
1278  RETURN
1279 999 CALL interface_condition_geometry_finalise(interface_condition%GEOMETRY,dummy_err,dummy_error,*998)
1280 998 errorsexits("INTERFACE_CONDITION_GEOMETRY_INITIALISE",err,error)
1281  RETURN 1
1282  END SUBROUTINE interface_condition_geometry_initialise
1283 
1284  !
1285  !================================================================================================================================
1286  !
1287 
1289  SUBROUTINE interface_condition_initialise(INTERFACE_CONDITION,ERR,ERROR,*)
1290 
1291  !Argument variables
1292  TYPE(interface_condition_type), POINTER :: interface_condition
1293  INTEGER(INTG), INTENT(OUT) :: err
1294  TYPE(varying_string), INTENT(OUT) :: error
1295  !Local Variables
1296  INTEGER(INTG) :: dummy_err
1297  TYPE(varying_string) :: dummy_error
1298 
1299  enters("INTERFACE_CONDITION_INITIALISE",err,error,*998)
1300 
1301  IF(ASSOCIATED(interface_condition)) THEN
1302  CALL flagerror("Interface condition is already associated.",err,error,*998)
1303  ELSE
1304  ALLOCATE(interface_condition,stat=err)
1305  IF(err/=0) CALL flagerror("Could not allocate interface condition.",err,error,*999)
1306  interface_condition%USER_NUMBER=0
1307  interface_condition%GLOBAL_NUMBER=0
1308  interface_condition%INTERFACE_CONDITION_FINISHED=.false.
1309  NULLIFY(interface_condition%INTERFACE_CONDITIONS)
1310  NULLIFY(interface_condition%INTERFACE)
1311  interface_condition%METHOD=0
1312  interface_condition%OPERATOR=0
1313  NULLIFY(interface_condition%LAGRANGE)
1314  NULLIFY(interface_condition%PENALTY)
1315  NULLIFY(interface_condition%DEPENDENT)
1316  NULLIFY(interface_condition%INTERFACE_EQUATIONS)
1317  CALL interface_condition_geometry_initialise(interface_condition,err,error,*999)
1318  NULLIFY(interface_condition%BOUNDARY_CONDITIONS)
1319  ENDIF
1320 
1321  exits("INTERFACE_CONDITION_INITIALISE")
1322  RETURN
1323 999 CALL interface_condition_finalise(interface_condition,dummy_err,dummy_error,*998)
1324 998 errorsexits("INTERFACE_CONDITION_INITIALISE",err,error)
1325  RETURN 1
1326  END SUBROUTINE interface_condition_initialise
1327 
1328  !
1329  !================================================================================================================================
1330  !
1331 
1333  SUBROUTINE interfacecondition_lagrangefieldcreatefinish(INTERFACE_CONDITION,ERR,ERROR,*)
1334 
1335  !Argument variables
1336  TYPE(interface_condition_type), POINTER :: interface_condition
1337  INTEGER(INTG), INTENT(OUT) :: err
1338  TYPE(varying_string), INTENT(OUT) :: error
1339  !Local Variables
1340  INTEGER(INTG) :: lagrangefielduvariablenumberofcomponents,lagrangefielddeludelnvariablenumberofcomponents
1341 
1342  enters("InterfaceCondition_LagrangeFieldCreateFinish",err,error,*999)
1343 
1344  IF(ASSOCIATED(interface_condition)) THEN
1345  IF(ASSOCIATED(interface_condition%LAGRANGE)) THEN
1346  IF(interface_condition%LAGRANGE%LAGRANGE_FINISHED) THEN
1347  CALL flagerror("Interface condition Lagrange field has already been finished.",err,error,*999)
1348  ELSE
1349  !Finish the Lagrange field creation
1350  IF(interface_condition%LAGRANGE%LAGRANGE_FIELD_AUTO_CREATED) THEN
1351  CALL field_create_finish(interface_condition%LAGRANGE%LAGRANGE_FIELD,err,error,*999)
1352  ENDIF
1353  interface_condition%LAGRANGE%LAGRANGE_FINISHED=.true.
1354  !\todo test following condition using some other method since FIELD_NUMBER_OF_COMPONENTS_GET requires the field to be finished which is what occurs above, but below condition needs to be checked before this.
1355  CALL field_number_of_components_get(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_u_variable_type, &
1356  & lagrangefielduvariablenumberofcomponents,err,error,*999)
1357  CALL field_number_of_components_get(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_deludeln_variable_type, &
1358  & lagrangefielddeludelnvariablenumberofcomponents,err,error,*999)
1359  IF (lagrangefielduvariablenumberofcomponents /= lagrangefielddeludelnvariablenumberofcomponents) THEN
1360  CALL flagerror("Interface Lagrange field U and DelUDelN variable components do not match.",err,error,*999)
1361  ENDIF
1362  ENDIF
1363  ELSE
1364  CALL flagerror("Interface condition Lagrange is not associated.",err,error,*999)
1365  ENDIF
1366  ELSE
1367  CALL flagerror("Interface condition is not associated.",err,error,*999)
1368  ENDIF
1369 
1370  exits("InterfaceCondition_LagrangeFieldCreateFinish")
1371  RETURN
1372 999 errors("InterfaceCondition_LagrangeFieldCreateFinish",err,error)
1373  exits("InterfaceCondition_LagrangeFieldCreateFinish")
1374  RETURN 1
1375 
1376  END SUBROUTINE interfacecondition_lagrangefieldcreatefinish
1377 
1378  !
1379  !================================================================================================================================
1380  !
1381 
1383  SUBROUTINE interfacecondition_lagrangefieldcreatestart(INTERFACE_CONDITION,LAGRANGE_FIELD_USER_NUMBER,LAGRANGE_FIELD, &
1384  & err,error,*)
1385 
1386  !Argument variables
1387  TYPE(interface_condition_type), POINTER :: interface_condition
1388  INTEGER(INTG), INTENT(IN) :: lagrange_field_user_number
1389  TYPE(field_type), POINTER :: lagrange_field
1390  INTEGER(INTG), INTENT(OUT) :: err
1391  TYPE(varying_string), INTENT(OUT) :: error
1392  !Local Variables
1393  INTEGER(INTG) :: component_idx,interpolation_type,geometric_scaling_type,dependent_variable_number
1394  TYPE(decomposition_type), POINTER :: geometric_decomposition
1395  TYPE(field_type), POINTER :: field
1396  TYPE(interface_type), POINTER :: interface
1397  TYPE(interface_dependent_type), POINTER :: interface_dependent
1398  TYPE(region_type), POINTER :: interface_region,lagrange_field_region
1399  TYPE(varying_string) :: local_error
1400 
1401  enters("InterfaceCondition_LagrangeFieldCreateStart",err,error,*999)
1402 
1403  IF(ASSOCIATED(interface_condition)) THEN
1404  IF(ASSOCIATED(interface_condition%LAGRANGE)) THEN
1405  CALL flagerror("Interface condition Lagrange is already associated.",err,error,*999)
1406  ELSE
1407  interface_dependent=>interface_condition%DEPENDENT
1408  IF(ASSOCIATED(interface_dependent)) THEN
1409  interface=>interface_condition%INTERFACE
1410  IF(ASSOCIATED(interface)) THEN
1411  interface_region=>interface%PARENT_REGION
1412  IF(ASSOCIATED(interface_region)) THEN
1413  IF(ASSOCIATED(lagrange_field)) THEN
1414  !Check the Lagrange field has been finished
1415  IF(lagrange_field%FIELD_FINISHED) THEN
1416  !Check the user numbers match
1417  IF(lagrange_field_user_number/=lagrange_field%USER_NUMBER) THEN
1418  local_error="The specified Lagrange field user number of "// &
1419  & trim(number_to_vstring(lagrange_field_user_number,"*",err,error))// &
1420  & " does not match the user number of the specified Lagrange field of "// &
1421  & trim(number_to_vstring(lagrange_field%USER_NUMBER,"*",err,error))//"."
1422  CALL flagerror(local_error,err,error,*999)
1423  ENDIF
1424  lagrange_field_region=>lagrange_field%REGION
1425  IF(ASSOCIATED(lagrange_field_region)) THEN
1426  !Check the field is defined on the same region as the interface
1427  IF(lagrange_field_region%USER_NUMBER/=interface_region%USER_NUMBER) THEN
1428  local_error="Invalid region setup. The specified Lagrange field has been created on interface number "// &
1429  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//" in parent region number "// &
1430  & trim(number_to_vstring(lagrange_field_region%USER_NUMBER,"*",err,error))// &
1431  & " and the specified interface has been created in parent region number "// &
1432  & trim(number_to_vstring(interface_region%USER_NUMBER,"*",err,error))//"."
1433  CALL flagerror(local_error,err,error,*999)
1434  ENDIF
1435  ELSE
1436  CALL flagerror("The Lagrange field region is not associated.",err,error,*999)
1437  ENDIF
1438  ELSE
1439  CALL flagerror("The specified Lagrange field has not been finished.",err,error,*999)
1440  ENDIF
1441  ELSE
1442  !Check the user number has not already been used for a field in this region.
1443  NULLIFY(field)
1444  CALL field_user_number_find(lagrange_field_user_number,interface,field,err,error,*999)
1445  IF(ASSOCIATED(field)) THEN
1446  local_error="The specified Lagrange field user number of "// &
1447  & trim(number_to_vstring(lagrange_field_user_number,"*",err,error))// &
1448  & " has already been used to create a field on interface number "// &
1449  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//"."
1450  CALL flagerror(local_error,err,error,*999)
1451  ENDIF
1452  ENDIF
1453  CALL interface_condition_lagrange_initialise(interface_condition,err,error,*999)
1454  IF(.NOT.ASSOCIATED(lagrange_field)) THEN
1455  !Create the Lagrange field
1456  interface_condition%LAGRANGE%LAGRANGE_FIELD_AUTO_CREATED=.true.
1457  CALL field_create_start(lagrange_field_user_number,interface_condition%INTERFACE,interface_condition%LAGRANGE% &
1458  & lagrange_field,err,error,*999)
1459  CALL field_label_set(interface_condition%LAGRANGE%LAGRANGE_FIELD,"Lagrange Multipliers Field",err,error,*999)
1460  CALL field_type_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_general_type,err,error,*999)
1461  CALL field_dependent_type_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_dependent_type, &
1462  & err,error,*999)
1463  NULLIFY(geometric_decomposition)
1464  CALL field_mesh_decomposition_get(interface_condition%GEOMETRY%GEOMETRIC_FIELD,geometric_decomposition, &
1465  & err,error,*999)
1466  CALL field_mesh_decomposition_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,geometric_decomposition, &
1467  & err,error,*999)
1468  CALL field_geometric_field_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,interface_condition%GEOMETRY% &
1469  & geometric_field,err,error,*999)
1470  CALL field_number_of_variables_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,2,err,error,*999)
1471  CALL field_variable_types_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,[field_u_variable_type, &
1472  & field_deludeln_variable_type],err,error,*999)
1473  CALL field_variable_label_set(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_u_variable_type,"Lambda", &
1474  & err,error,*999)
1475  CALL field_variable_label_set(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_deludeln_variable_type, &
1476  & "Lambda RHS",err,error,*999)
1477  CALL field_dimension_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_u_variable_type, &
1478  & field_vector_dimension_type,err,error,*999)
1479  CALL field_dimension_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_deludeln_variable_type, &
1480  & field_vector_dimension_type,err,error,*999)
1481  CALL field_data_type_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_u_variable_type, &
1482  & field_dp_type,err,error,*999)
1483  CALL field_data_type_set_and_lock(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_deludeln_variable_type, &
1484  & field_dp_type,err,error,*999)
1485  !Note that only components present in both the coupled meshes interface dependent fields can be coupled
1486  !Default the number of component to be the minimum number of components across all the coupled dependent variables
1487  !\todo Check ordering of variable components which are coupled and uncoupled are handled correctly to ensure that
1488  !coupled variable components don't have to always come before the uncoupled variable components
1489  interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS=0
1490  DO dependent_variable_number=1,interface_dependent%NUMBER_OF_DEPENDENT_VARIABLES
1491  IF (interface_dependent%FIELD_VARIABLES(dependent_variable_number)%PTR%NUMBER_OF_COMPONENTS< &
1492  & interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS) THEN
1493  interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS= &
1494  & interface_dependent%FIELD_VARIABLES(dependent_variable_number)%PTR%NUMBER_OF_COMPONENTS
1495  ELSEIF (interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS==0) THEN
1496  interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS= &
1497  & interface_dependent%FIELD_VARIABLES(dependent_variable_number)%PTR%NUMBER_OF_COMPONENTS
1498  ENDIF
1499  ENDDO
1500  ! Remove pressure component from number of coupled components
1501  ! INTERFACE_CONDITION_SOLID_FLUID_OPERATOR might not be used as it is equivalent to
1502  ! INTERFACE_CONDITION_FIELD_CONTINUITY_OPERATOR if set up correctly
1503  IF (interface_condition%OPERATOR==interface_condition_solid_fluid_operator) THEN
1504  interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS=interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS-1
1505  ENDIF
1506  CALL field_number_of_components_set(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_u_variable_type, &
1507  & interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS,err,error,*999)
1508  CALL field_number_of_components_set(interface_condition%LAGRANGE%LAGRANGE_FIELD,field_deludeln_variable_type, &
1509  & interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS,err,error,*999)
1510  DO component_idx=1,interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS
1511  CALL field_component_interpolation_get(interface_dependent%FIELD_VARIABLES(1)%PTR%FIELD,field_u_variable_type, &
1512  & component_idx,interpolation_type,err,error,*999)
1513  CALL field_component_interpolation_set(interface_condition%LAGRANGE%LAGRANGE_FIELD, &
1514  & field_u_variable_type,component_idx,interpolation_type,err,error,*999)
1515  CALL field_component_interpolation_set(interface_condition%LAGRANGE%LAGRANGE_FIELD, &
1516  & field_deludeln_variable_type,component_idx,interpolation_type,err,error,*999)
1517  ENDDO !component_idx
1518  CALL field_scaling_type_get(interface_condition%GEOMETRY%GEOMETRIC_FIELD,geometric_scaling_type, &
1519  & err,error,*999)
1520  CALL field_scaling_type_set(interface_condition%LAGRANGE%LAGRANGE_FIELD,geometric_scaling_type, &
1521  & err,error,*999)
1522  ELSE
1523  !Check the Lagrange field
1524  CALL flagerror("Not implemented.",err,error,*999)
1525  ENDIF
1526  !Set pointers
1527  IF(interface_condition%LAGRANGE%LAGRANGE_FIELD_AUTO_CREATED) THEN
1528  lagrange_field=>interface_condition%LAGRANGE%LAGRANGE_FIELD
1529  ELSE
1530  interface_condition%LAGRANGE%LAGRANGE_FIELD=>lagrange_field
1531  ENDIF
1532  ELSE
1533  CALL flagerror("The interface parent region is not associated.",err,error,*999)
1534  ENDIF
1535  ELSE
1536  CALL flagerror("The interface interface conditions is not associated.",err,error,*999)
1537  ENDIF
1538  ELSE
1539  CALL flagerror("Interface condition dependent is not associated.",err,error,*999)
1540  ENDIF
1541  ENDIF
1542  ELSE
1543  CALL flagerror("Interface conditions is not associated.",err,error,*999)
1544  ENDIF
1545 
1546  exits("InterfaceCondition_LagrangeFieldCreateStart")
1547  RETURN
1548 999 errorsexits("InterfaceCondition_LagrangeFieldCreateStart",err,error)
1549  RETURN 1
1550 
1551  END SUBROUTINE interfacecondition_lagrangefieldcreatestart
1552 
1553  !
1554  !================================================================================================================================
1555  !
1556 
1558  SUBROUTINE interface_condition_lagrange_finalise(INTERFACE_LAGRANGE,ERR,ERROR,*)
1559 
1560  !Argument variables
1561  TYPE(interface_lagrange_type), POINTER :: interface_lagrange
1562  INTEGER(INTG), INTENT(OUT) :: err
1563  TYPE(varying_string), INTENT(OUT) :: error
1564  !Local Variables
1565 
1566  enters("INTERFACE_CONDITION_LAGRANGE_FINALISE",err,error,*999)
1567 
1568  IF(ASSOCIATED(interface_lagrange)) THEN
1569  DEALLOCATE(interface_lagrange)
1570  ENDIF
1571 
1572  exits("INTERFACE_CONDITION_LAGRANGE_FINALISE")
1573  RETURN
1574 999 errorsexits("INTERFACE_CONDITION_LAGRANGE_FINALISE",err,error)
1575  RETURN 1
1576  END SUBROUTINE interface_condition_lagrange_finalise
1577 
1578  !
1579  !================================================================================================================================
1580  !
1581 
1583  SUBROUTINE interface_condition_lagrange_initialise(INTERFACE_CONDITION,ERR,ERROR,*)
1584 
1585  !Argument variables
1586  TYPE(interface_condition_type), POINTER :: interface_condition
1587  INTEGER(INTG), INTENT(OUT) :: err
1588  TYPE(varying_string), INTENT(OUT) :: error
1589  !Local Variables
1590  INTEGER(INTG) :: dummy_err
1591  TYPE(varying_string) :: dummy_error
1592 
1593  enters("INTERFACE_CONDITION_LAGRANGE_INITIALISE",err,error,*998)
1594 
1595  IF(ASSOCIATED(interface_condition)) THEN
1596  IF(ASSOCIATED(interface_condition%LAGRANGE)) THEN
1597  CALL flagerror("Interface condition Lagrange is already associated.",err,error,*999)
1598  ELSE
1599  ALLOCATE(interface_condition%LAGRANGE,stat=err)
1600  IF(err/=0) CALL flagerror("Could not allocate interface condition Lagrange.",err,error,*999)
1601  interface_condition%LAGRANGE%INTERFACE_CONDITION=>interface_condition
1602  interface_condition%LAGRANGE%LAGRANGE_FINISHED=.false.
1603  interface_condition%LAGRANGE%LAGRANGE_FIELD_AUTO_CREATED=.false.
1604  NULLIFY(interface_condition%LAGRANGE%LAGRANGE_FIELD)
1605  interface_condition%LAGRANGE%NUMBER_OF_COMPONENTS=0
1606  ENDIF
1607  ELSE
1608  CALL flagerror("Interface condition is not associated.",err,error,*999)
1609  ENDIF
1610 
1611  exits("INTERFACE_CONDITION_LAGRANGE_INITIALISE")
1612  RETURN
1613 999 CALL interface_condition_lagrange_finalise(interface_condition%LAGRANGE,dummy_err,dummy_error,*998)
1614 998 errorsexits("INTERFACE_CONDITION_LAGRANGE_INITIALISE",err,error)
1615  RETURN 1
1616  END SUBROUTINE interface_condition_lagrange_initialise
1617 
1618  !
1619  !================================================================================================================================
1620  !
1621 
1623  SUBROUTINE interfacecondition_penaltyfieldcreatefinish(INTERFACE_CONDITION,ERR,ERROR,*)
1624 
1625  !Argument variables
1626  TYPE(interface_condition_type), POINTER :: interface_condition
1627  INTEGER(INTG), INTENT(OUT) :: err
1628  TYPE(varying_string), INTENT(OUT) :: error
1629  !Local Variables
1630 
1631  enters("InterfaceCondition_PenaltyFieldCreateFinish",err,error,*999)
1632 
1633  IF(ASSOCIATED(interface_condition)) THEN
1634  IF(ASSOCIATED(interface_condition%PENALTY)) THEN
1635  IF(interface_condition%PENALTY%PENALTY_FINISHED) THEN
1636  CALL flagerror("Interface condition penalty field has already been finished.",err,error,*999)
1637  ELSE
1638  !Finish the penalty field creation
1639  IF(interface_condition%PENALTY%PENALTY_FIELD_AUTO_CREATED) THEN
1640  CALL field_create_finish(interface_condition%PENALTY%PENALTY_FIELD,err,error,*999)
1641  ENDIF
1642  interface_condition%PENALTY%PENALTY_FINISHED=.true.
1643  ENDIF
1644  ELSE
1645  CALL flagerror("Interface condition penalty is not associated.",err,error,*999)
1646  ENDIF
1647  ELSE
1648  CALL flagerror("Interface condition is not associated.",err,error,*999)
1649  ENDIF
1650 
1651  exits("InterfaceCondition_PenaltyFieldCreateFinish")
1652  RETURN
1653 999 errorsexits("InterfaceCondition_PenaltyFieldCreateFinish",err,error)
1654  RETURN 1
1655 
1656  END SUBROUTINE interfacecondition_penaltyfieldcreatefinish
1657 
1658  !
1659  !================================================================================================================================
1660  !
1661 
1663  SUBROUTINE interfacecondition_penaltyfieldcreatestart(INTERFACE_CONDITION,PENALTY_FIELD_USER_NUMBER,PENALTY_FIELD, &
1664  & err,error,*)
1665 
1666  !Argument variables
1667  TYPE(interface_condition_type), POINTER :: interface_condition
1668  INTEGER(INTG), INTENT(IN) :: penalty_field_user_number
1669  TYPE(field_type), POINTER :: penalty_field
1670  INTEGER(INTG), INTENT(OUT) :: err
1671  TYPE(varying_string), INTENT(OUT) :: error
1672  !Local Variables
1673  INTEGER(INTG) :: component_idx,geometric_scaling_type
1674  TYPE(decomposition_type), POINTER :: geometric_decomposition
1675  TYPE(field_type), POINTER :: field
1676  TYPE(interface_type), POINTER :: interface
1677  TYPE(interface_dependent_type), POINTER :: interface_dependent
1678  TYPE(region_type), POINTER :: interface_region,penalty_field_region
1679  TYPE(varying_string) :: local_error
1680 
1681  enters("InterfaceCondition_PenaltyFieldCreateStart",err,error,*999)
1682 
1683  IF(ASSOCIATED(interface_condition)) THEN
1684  IF(ASSOCIATED(interface_condition%PENALTY)) THEN
1685  CALL flagerror("Interface condition penalty is already associated.",err,error,*999)
1686  ELSE
1687  interface_dependent=>interface_condition%DEPENDENT
1688  IF(ASSOCIATED(interface_dependent)) THEN
1689  interface=>interface_condition%INTERFACE
1690  IF(ASSOCIATED(interface)) THEN
1691  interface_region=>interface%PARENT_REGION
1692  IF(ASSOCIATED(interface_region)) THEN
1693  IF(ASSOCIATED(penalty_field)) THEN
1694  !Check the penalty field has been finished
1695  IF(penalty_field%FIELD_FINISHED) THEN
1696  !Check the user numbers match
1697  IF(penalty_field_user_number/=penalty_field%USER_NUMBER) THEN
1698  local_error="The specified penalty field user number of "// &
1699  & trim(number_to_vstring(penalty_field_user_number,"*",err,error))// &
1700  & " does not match the user number of the specified penalty field of "// &
1701  & trim(number_to_vstring(penalty_field%USER_NUMBER,"*",err,error))//"."
1702  CALL flagerror(local_error,err,error,*999)
1703  ENDIF
1704  penalty_field_region=>penalty_field%REGION
1705  IF(ASSOCIATED(penalty_field_region)) THEN
1706  !Check the field is defined on the same region as the interface
1707  IF(penalty_field_region%USER_NUMBER/=interface_region%USER_NUMBER) THEN
1708  local_error="Invalid region setup. The specified penalty field has been created on interface number "// &
1709  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//" in parent region number "// &
1710  & trim(number_to_vstring(penalty_field_region%USER_NUMBER,"*",err,error))// &
1711  & " and the specified interface has been created in parent region number "// &
1712  & trim(number_to_vstring(interface_region%USER_NUMBER,"*",err,error))//"."
1713  CALL flagerror(local_error,err,error,*999)
1714  ENDIF
1715  ELSE
1716  CALL flagerror("The penalty field region is not associated.",err,error,*999)
1717  ENDIF
1718  ELSE
1719  CALL flagerror("The specified penalty field has not been finished.",err,error,*999)
1720  ENDIF
1721  ELSE
1722  !Check the user number has not already been used for a field in this region.
1723  NULLIFY(field)
1724  CALL field_user_number_find(penalty_field_user_number,interface,field,err,error,*999)
1725  IF(ASSOCIATED(field)) THEN
1726  local_error="The specified penalty field user number of "// &
1727  & trim(number_to_vstring(penalty_field_user_number,"*",err,error))// &
1728  & " has already been used to create a field on interface number "// &
1729  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//"."
1730  CALL flagerror(local_error,err,error,*999)
1731  ENDIF
1732  ENDIF
1733  CALL interface_condition_penalty_initialise(interface_condition,err,error,*999)
1734  IF(.NOT.ASSOCIATED(penalty_field)) THEN
1735  !Create the penalty field
1736  interface_condition%PENALTY%PENALTY_FIELD_AUTO_CREATED=.true.
1737  CALL field_create_start(penalty_field_user_number,interface_condition%INTERFACE,interface_condition%PENALTY% &
1738  & penalty_field,err,error,*999)
1739  CALL field_label_set(interface_condition%PENALTY%PENALTY_FIELD,"Penalty Field",err,error,*999)
1740  CALL field_type_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,field_general_type,err,error,*999)
1741  CALL field_dependent_type_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,field_dependent_type, &
1742  & err,error,*999)
1743  NULLIFY(geometric_decomposition)
1744  CALL field_mesh_decomposition_get(interface_condition%GEOMETRY%GEOMETRIC_FIELD,geometric_decomposition, &
1745  & err,error,*999)
1746  CALL field_mesh_decomposition_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,geometric_decomposition, &
1747  & err,error,*999)
1748  CALL field_geometric_field_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,interface_condition%GEOMETRY% &
1749  & geometric_field,err,error,*999)
1750  CALL field_number_of_variables_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,1,err,error,*999)
1751  CALL field_variable_types_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,[field_u_variable_type], &
1752  & err,error,*999)
1753  CALL field_variable_label_set(interface_condition%PENALTY%PENALTY_FIELD,field_u_variable_type,"Alpha", &
1754  & err,error,*999)
1755  CALL field_dimension_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,field_u_variable_type, &
1756  & field_vector_dimension_type,err,error,*999)
1757  CALL field_data_type_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD,field_u_variable_type, &
1758  & field_dp_type,err,error,*999)
1759  IF(interface_condition%OPERATOR==interface_condition_fls_contact_operator .OR. &
1760  & interface_condition%OPERATOR==interface_condition_fls_contact_reproject_operator) THEN
1761  !Default 1 component for the contact lagrange variable in a frictionless contact problem
1762  CALL field_number_of_components_set(interface_condition%PENALTY%PENALTY_FIELD,field_u_variable_type, &
1763  & 1,err,error,*999)
1764  CALL field_component_interpolation_set(interface_condition%PENALTY%PENALTY_FIELD, &
1765  & field_u_variable_type,1,field_constant_interpolation,err,error,*999)
1766  ELSE
1767  !Default the number of component to the first variable of the interface dependent field's number of components,
1768  CALL field_number_of_components_set(interface_condition%PENALTY%PENALTY_FIELD,field_u_variable_type, &
1769  & interface_dependent%FIELD_VARIABLES(1)%PTR%NUMBER_OF_COMPONENTS,err,error,*999)
1770  DO component_idx=1,interface_dependent%FIELD_VARIABLES(1)%PTR%NUMBER_OF_COMPONENTS
1771  CALL field_component_interpolation_set_and_lock(interface_condition%PENALTY%PENALTY_FIELD, &
1772  & field_u_variable_type,component_idx,field_constant_interpolation,err,error,*999)
1773  ENDDO !component_idx
1774  ENDIF
1775  CALL field_scaling_type_get(interface_condition%GEOMETRY%GEOMETRIC_FIELD,geometric_scaling_type, &
1776  & err,error,*999)
1777  CALL field_scaling_type_set(interface_condition%PENALTY%PENALTY_FIELD,geometric_scaling_type, &
1778  & err,error,*999)
1779  ELSE
1780  !Check the penalty field
1781  CALL flagerror("Not implemented.",err,error,*999)
1782  ENDIF
1783  !Set pointers
1784  IF(interface_condition%PENALTY%PENALTY_FIELD_AUTO_CREATED) THEN
1785  penalty_field=>interface_condition%PENALTY%PENALTY_FIELD
1786  ELSE
1787  interface_condition%PENALTY%PENALTY_FIELD=>penalty_field
1788  ENDIF
1789  ELSE
1790  CALL flagerror("The interface parent region is not associated.",err,error,*999)
1791  ENDIF
1792  ELSE
1793  CALL flagerror("The interface interface conditions is not associated.",err,error,*999)
1794  ENDIF
1795  ELSE
1796  CALL flagerror("Interface condition dependent is not associated.",err,error,*999)
1797  ENDIF
1798  ENDIF
1799  ELSE
1800  CALL flagerror("Interface conditions is not associated.",err,error,*999)
1801  ENDIF
1802 
1803  exits("InterfaceCondition_PenaltyFieldCreateStart")
1804  RETURN
1805 999 errorsexits("InterfaceCondition_PenaltyFieldCreateStart",err,error)
1806  RETURN 1
1807 
1808  END SUBROUTINE interfacecondition_penaltyfieldcreatestart
1809 
1810  !
1811  !================================================================================================================================
1812  !
1813 
1815  SUBROUTINE interface_condition_penalty_finalise(INTERFACE_PENALTY,ERR,ERROR,*)
1816 
1817  !Argument variables
1818  TYPE(interface_penalty_type), POINTER :: interface_penalty
1819  INTEGER(INTG), INTENT(OUT) :: err
1820  TYPE(varying_string), INTENT(OUT) :: error
1821  !Local Variables
1822 
1823  enters("INTERFACE_CONDITION_PENALTY_FINALISE",err,error,*999)
1824 
1825  IF(ASSOCIATED(interface_penalty)) THEN
1826  DEALLOCATE(interface_penalty)
1827  ENDIF
1828 
1829  exits("INTERFACE_CONDITION_PENALTY_FINALISE")
1830  RETURN
1831 999 errorsexits("INTERFACE_CONDITION_PENALTY_FINALISE",err,error)
1832  RETURN 1
1833  END SUBROUTINE interface_condition_penalty_finalise
1834 
1835  !
1836  !================================================================================================================================
1837  !
1838 
1840  SUBROUTINE interface_condition_penalty_initialise(INTERFACE_CONDITION,ERR,ERROR,*)
1841 
1842  !Argument variables
1843  TYPE(interface_condition_type), POINTER :: interface_condition
1844  INTEGER(INTG), INTENT(OUT) :: err
1845  TYPE(varying_string), INTENT(OUT) :: error
1846  !Local Variables
1847  INTEGER(INTG) :: dummy_err
1848  TYPE(varying_string) :: dummy_error
1849 
1850  enters("INTERFACE_CONDITION_PENALTY_INITIALISE",err,error,*998)
1851 
1852  IF(ASSOCIATED(interface_condition)) THEN
1853  IF(ASSOCIATED(interface_condition%PENALTY)) THEN
1854  CALL flagerror("Interface condition penalty is already associated.",err,error,*999)
1855  ELSE
1856  ALLOCATE(interface_condition%PENALTY,stat=err)
1857  IF(err/=0) CALL flagerror("Could not allocate interface condition penalty.",err,error,*999)
1858  interface_condition%PENALTY%INTERFACE_CONDITION=>interface_condition
1859  interface_condition%PENALTY%PENALTY_FINISHED=.false.
1860  interface_condition%PENALTY%PENALTY_FIELD_AUTO_CREATED=.false.
1861  NULLIFY(interface_condition%PENALTY%PENALTY_FIELD)
1862  ENDIF
1863  ELSE
1864  CALL flagerror("Interface condition is not associated.",err,error,*999)
1865  ENDIF
1866 
1867  exits("INTERFACE_CONDITION_PENALTY_INITIALISE")
1868  RETURN
1869 999 CALL interface_condition_penalty_finalise(interface_condition%PENALTY,dummy_err,dummy_error,*998)
1870 998 errorsexits("INTERFACE_CONDITION_PENALTY_INITIALISE",err,error)
1871  RETURN 1
1872  END SUBROUTINE interface_condition_penalty_initialise
1873 
1874  !
1875  !================================================================================================================================
1876  !
1877 
1879  SUBROUTINE interface_condition_method_get(INTERFACE_CONDITION,INTERFACE_CONDITION_METHOD,ERR,ERROR,*)
1880 
1881  !Argument variables
1882  TYPE(interface_condition_type), POINTER :: interface_condition
1883  INTEGER(INTG), INTENT(OUT) :: interface_condition_method
1884  INTEGER(INTG), INTENT(OUT) :: err
1885  TYPE(varying_string), INTENT(OUT) :: error
1886  !Local Variables
1887 
1888  enters("INTERFACE_CONDITION_METHOD_GET",err,error,*999)
1889 
1890  IF(ASSOCIATED(interface_condition)) THEN
1891  IF(interface_condition%INTERFACE_CONDITION_FINISHED) THEN
1892  interface_condition_method=interface_condition%METHOD
1893  ELSE
1894  CALL flagerror("Interface condition has not been finished.",err,error,*999)
1895  ENDIF
1896  ELSE
1897  CALL flagerror("Interface condition is not associated.",err,error,*999)
1898  ENDIF
1899 
1900  exits("INTERFACE_CONDITION_METHOD_GET")
1901  RETURN
1902 999 errorsexits("INTERFACE_CONDITION_METHOD_GET",err,error)
1903  RETURN 1
1904  END SUBROUTINE interface_condition_method_get
1905 
1906  !
1907  !================================================================================================================================
1908  !
1909 
1911  SUBROUTINE interface_condition_method_set(INTERFACE_CONDITION,INTERFACE_CONDITION_METHOD,ERR,ERROR,*)
1912 
1913  !Argument variables
1914  TYPE(interface_condition_type), POINTER :: interface_condition
1915  INTEGER(INTG), INTENT(IN) :: interface_condition_method
1916  INTEGER(INTG), INTENT(OUT) :: err
1917  TYPE(varying_string), INTENT(OUT) :: error
1918  !Local Variables
1919  TYPE(varying_string) :: local_error
1920 
1921  enters("INTERFACE_CONDITION_METHOD_SET",err,error,*999)
1922 
1923  IF(ASSOCIATED(interface_condition)) THEN
1924  IF(interface_condition%INTERFACE_CONDITION_FINISHED) THEN
1925  CALL flagerror("Interface condition has been finished.",err,error,*999)
1926  ELSE
1927  SELECT CASE(interface_condition_method)
1929  interface_condition%METHOD=interface_condition_point_to_point_method
1931  interface_condition%METHOD=interface_condition_lagrange_multipliers_method
1933  interface_condition%METHOD=interface_condition_augmented_lagrange_method
1935  interface_condition%METHOD=interface_condition_penalty_method
1936  CASE DEFAULT
1937  local_error="The specified interface condition method of "// &
1938  & trim(number_to_vstring(interface_condition_method,"*",err,error))//" is not valid."
1939  CALL flagerror(local_error,err,error,*999)
1940  END SELECT
1941  ENDIF
1942  ELSE
1943  CALL flagerror("Interface condition is not associated.",err,error,*999)
1944  ENDIF
1945 
1946  exits("INTERFACE_CONDITION_METHOD_SET")
1947  RETURN
1948 999 errorsexits("INTERFACE_CONDITION_METHOD_SET",err,error)
1949  RETURN 1
1950  END SUBROUTINE interface_condition_method_set
1951 
1952  !
1953  !================================================================================================================================
1954  !
1955 
1957  SUBROUTINE interface_condition_operator_get(INTERFACE_CONDITION,INTERFACE_CONDITION_OPERATOR,ERR,ERROR,*)
1958 
1959  !Argument variables
1960  TYPE(interface_condition_type), POINTER :: interface_condition
1961  INTEGER(INTG), INTENT(OUT) :: interface_condition_operator
1962  INTEGER(INTG), INTENT(OUT) :: err
1963  TYPE(varying_string), INTENT(OUT) :: error
1964  !Local Variables
1965 
1966  enters("INTERFACE_CONDITION_OPERATOR_GET",err,error,*999)
1967 
1968  IF(ASSOCIATED(interface_condition)) THEN
1969  IF(interface_condition%INTERFACE_CONDITION_FINISHED) THEN
1970  interface_condition_operator=interface_condition%OPERATOR
1971  ELSE
1972  CALL flagerror("Interface condition has not been finished.",err,error,*999)
1973  ENDIF
1974  ELSE
1975  CALL flagerror("Interface condition is not associated.",err,error,*999)
1976  ENDIF
1977 
1978  exits("INTERFACE_CONDITION_OPERATOR_GET")
1979  RETURN
1980 999 errorsexits("INTERFACE_CONDITION_OPERATOR_GET",err,error)
1981  RETURN 1
1982  END SUBROUTINE interface_condition_operator_get
1983 
1984  !
1985  !================================================================================================================================
1986  !
1987 
1989  SUBROUTINE interface_condition_operator_set(INTERFACE_CONDITION,INTERFACE_CONDITION_OPERATOR,ERR,ERROR,*)
1990 
1991  !Argument variables
1992  TYPE(interface_condition_type), POINTER :: interface_condition
1993  INTEGER(INTG), INTENT(IN) :: interface_condition_operator
1994  INTEGER(INTG), INTENT(OUT) :: err
1995  TYPE(varying_string), INTENT(OUT) :: error
1996  !Local Variables
1997  TYPE(varying_string) :: local_error
1998 
1999  enters("INTERFACE_CONDITION_OPERATOR_SET",err,error,*999)
2000 
2001  IF(ASSOCIATED(interface_condition)) THEN
2002  IF(interface_condition%INTERFACE_CONDITION_FINISHED) THEN
2003  CALL flagerror("Interface condition has been finished.",err,error,*999)
2004  ELSE
2005  SELECT CASE(interface_condition_operator)
2007  interface_condition%OPERATOR=interface_condition_field_continuity_operator
2009  interface_condition%OPERATOR=interface_condition_field_normal_continuity_operator
2011  interface_condition%OPERATOR=interface_condition_fls_contact_operator
2013  interface_condition%OPERATOR=interface_condition_fls_contact_reproject_operator
2015  interface_condition%OPERATOR=interface_condition_solid_fluid_operator
2017  interface_condition%OPERATOR=interface_condition_solid_fluid_normal_operator
2018  CASE DEFAULT
2019  local_error="The specified interface condition operator of "// &
2020  & trim(number_to_vstring(interface_condition_operator,"*",err,error))//" is not valid."
2021  CALL flagerror(local_error,err,error,*999)
2022  END SELECT
2023  ENDIF
2024  ELSE
2025  CALL flagerror("Interface condition is not associated.",err,error,*999)
2026  ENDIF
2027 
2028  exits("INTERFACE_CONDITION_OPERATOR_SET")
2029  RETURN
2030 999 errorsexits("INTERFACE_CONDITION_OPERATOR_SET",err,error)
2031  RETURN 1
2032  END SUBROUTINE interface_condition_operator_set
2033 
2034  !
2035  !================================================================================================================================
2036  !
2037 
2039  SUBROUTINE interface_condition_residual_evaluate(INTERFACE_CONDITION,ERR,ERROR,*)
2040 
2041  !Argument variables
2042  TYPE(interface_condition_type), POINTER :: interface_condition
2043  INTEGER(INTG), INTENT(OUT) :: err
2044  TYPE(varying_string), INTENT(OUT) :: error
2045  !Local Variables
2046  TYPE(interface_equations_type), POINTER :: interface_equations
2047  TYPE(varying_string) :: local_error
2048 
2049  enters("INTERFACE_CONDITION_RESIDUAL_EVALUATE",err,error,*999)
2050 
2051  IF(ASSOCIATED(interface_condition)) THEN
2052  interface_equations=>interface_condition%INTERFACE_EQUATIONS
2053  IF(ASSOCIATED(interface_equations)) THEN
2054  IF(interface_equations%INTERFACE_EQUATIONS_FINISHED) THEN
2055  SELECT CASE(interface_condition%METHOD)
2057  CALL interface_condition_residual_evaluate_fem(interface_condition,err,error,*999)
2059  CALL flagerror("Not implemented.",err,error,*999)
2061  CALL flagerror("Not implemented.",err,error,*999)
2062  CASE DEFAULT
2063  local_error="The interface condition method of "// &
2064  & trim(number_to_vstring(interface_condition%METHOD,"*",err,error))// &
2065  & " is invalid."
2066  CALL flagerror(local_error,err,error,*999)
2067  END SELECT
2068  ELSE
2069  CALL flagerror("Interface equations have not been finished.",err,error,*999)
2070  ENDIF
2071  ELSE
2072  CALL flagerror("Interface condition equations is not associated.",err,error,*999)
2073  ENDIF
2074  ELSE
2075  CALL flagerror("Interface condition is not associated.",err,error,*999)
2076  ENDIF
2077 
2078  exits("INTERFACE_CONDITION_RESIDUAL_EVALUATE")
2079  RETURN
2080 999 errorsexits("INTERFACE_CONDITION_RESIDUAL_EVALUATE",err,error)
2081  RETURN 1
2082 
2083  END SUBROUTINE interface_condition_residual_evaluate
2084 
2085  !
2086  !================================================================================================================================
2087  !
2088 
2090  SUBROUTINE interface_condition_residual_evaluate_fem(INTERFACE_CONDITION,ERR,ERROR,*)
2091 
2092  !Argument variables
2093  TYPE(interface_condition_type), POINTER :: interface_condition
2094  INTEGER(INTG), INTENT(OUT) :: err
2095  TYPE(varying_string), INTENT(OUT) :: error
2096  !Local Variables
2097  INTEGER(INTG) :: element_idx,ne,number_of_times
2098  REAL(SP) :: element_user_elapsed,element_system_elapsed,user_elapsed,user_time1(1),user_time2(1),user_time3(1),user_time4(1), &
2099  & USER_TIME5(1),USER_TIME6(1),SYSTEM_ELAPSED,SYSTEM_TIME1(1),SYSTEM_TIME2(1),SYSTEM_TIME3(1),SYSTEM_TIME4(1), &
2100  & SYSTEM_TIME5(1),SYSTEM_TIME6(1)
2101  TYPE(domain_mapping_type), POINTER :: elements_mapping
2102  TYPE(interface_equations_type), POINTER :: interface_equations
2103  TYPE(interface_lagrange_type), POINTER :: lagrange
2104  TYPE(interface_matrices_type), POINTER :: interface_matrices
2105  TYPE(field_type), POINTER :: lagrange_field
2106 
2107  enters("INTERFACE_CONDITION_RESIDUAL_EVALUATE_FEM",err,error,*999)
2108 
2109  IF(ASSOCIATED(interface_condition)) THEN
2110  lagrange=>interface_condition%LAGRANGE
2111  IF(ASSOCIATED(lagrange)) THEN
2112  lagrange_field=>interface_condition%LAGRANGE%LAGRANGE_FIELD
2113  IF(ASSOCIATED(lagrange_field)) THEN
2114  interface_equations=>interface_condition%INTERFACE_EQUATIONS
2115  IF(ASSOCIATED(interface_equations)) THEN
2116  interface_matrices=>interface_equations%INTERFACE_MATRICES
2117  IF(ASSOCIATED(interface_matrices)) THEN
2118  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
2119  CALL cpu_timer(user_cpu,user_time1,err,error,*999)
2120  CALL cpu_timer(system_cpu,system_time1,err,error,*999)
2121  ENDIF
2122 !!Do we need to transfer parameter sets???
2123  !Initialise the matrices and rhs vector
2124 #ifdef TAUPROF
2125  CALL tau_static_phase_start("INTERFACE_MATRICES_VALUES_INITIALISE()")
2126 #endif
2127  CALL interface_matrices_values_initialise(interface_matrices,0.0_dp,err,error,*999)
2128 #ifdef TAUPROF
2129  CALL tau_static_phase_stop("INTERFACE_MATRICES_VALUES_INITIALISE()")
2130 #endif
2131  !Assemble the elements
2132  !Allocate the element matrices
2133 #ifdef TAUPROF
2134  CALL tau_static_phase_start("InterfaceMatrices_ElementInitialise()")
2135 #endif
2136  CALL interfacematrices_elementinitialise(interface_matrices,err,error,*999)
2137  elements_mapping=>lagrange_field%DECOMPOSITION%DOMAIN(lagrange_field%DECOMPOSITION%MESH_COMPONENT_NUMBER)%PTR% &
2138  & mappings%ELEMENTS
2139 #ifdef TAUPROF
2140  CALL tau_static_phase_stop("InterfaceMatrices_ElementInitialise()")
2141 #endif
2142  !Output timing information if required
2143  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
2144  CALL cpu_timer(user_cpu,user_time2,err,error,*999)
2145  CALL cpu_timer(system_cpu,system_time2,err,error,*999)
2146  user_elapsed=user_time2(1)-user_time1(1)
2147  system_elapsed=system_time2(1)-system_time1(1)
2148  CALL write_string_value(general_output_type,"User time for interface setup and initialisation = ",user_elapsed, &
2149  & err,error,*999)
2150  CALL write_string_value(general_output_type,"System time for interface setup and initialisation = ", &
2151  & system_elapsed,err,error,*999)
2152  element_user_elapsed=0.0_sp
2153  element_system_elapsed=0.0_sp
2154  ENDIF
2155  number_of_times=0
2156  !Loop over the internal elements
2157 #ifdef TAUPROF
2158  CALL tau_static_phase_start("Internal Elements Loop")
2159 #endif
2160  DO element_idx=elements_mapping%INTERNAL_START,elements_mapping%INTERNAL_FINISH
2161  ne=elements_mapping%DOMAIN_LIST(element_idx)
2162  number_of_times=number_of_times+1
2163  CALL interfacematrices_elementcalculate(interface_matrices,ne,err,error,*999)
2164  CALL interfacecondition_finiteelementcalculate(interface_condition,ne,err,error,*999)
2165  CALL interface_matrices_element_add(interface_matrices,err,error,*999)
2166  ENDDO !element_idx
2167 #ifdef TAUPROF
2168  CALL tau_static_phase_stop("Internal Elements Loop")
2169 #endif
2170  !Output timing information if required
2171  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
2172  CALL cpu_timer(user_cpu,user_time3,err,error,*999)
2173  CALL cpu_timer(system_cpu,system_time3,err,error,*999)
2174  user_elapsed=user_time3(1)-user_time2(1)
2175  system_elapsed=system_time3(1)-system_time2(1)
2176  element_user_elapsed=user_elapsed
2177  element_system_elapsed=system_elapsed
2178  CALL write_string_value(general_output_type,"User time for internal interface assembly = ",user_elapsed, &
2179  & err,error,*999)
2180  CALL write_string_value(general_output_type,"System time for internal interface assembly = ",system_elapsed, &
2181  & err,error,*999)
2182  ENDIF
2183  !Output timing information if required
2184  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
2185  CALL cpu_timer(user_cpu,user_time4,err,error,*999)
2186  CALL cpu_timer(system_cpu,system_time4,err,error,*999)
2187  user_elapsed=user_time4(1)-user_time3(1)
2188  system_elapsed=system_time4(1)-system_time3(1)
2189  CALL write_string_value(general_output_type,"User time for parameter transfer completion = ",user_elapsed, &
2190  & err,error,*999)
2191  CALL write_string_value(general_output_type,"System time for parameter transfer completion = ",system_elapsed, &
2192  & err,error,*999)
2193  ENDIF
2194  !Loop over the boundary and ghost elements
2195 #ifdef TAUPROF
2196  CALL tau_static_phase_start("Boundary and Ghost Elements Loop")
2197 #endif
2198  DO element_idx=elements_mapping%BOUNDARY_START,elements_mapping%GHOST_FINISH
2199  ne=elements_mapping%DOMAIN_LIST(element_idx)
2200  number_of_times=number_of_times+1
2201  CALL interfacematrices_elementcalculate(interface_matrices,ne,err,error,*999)
2202  CALL interfacecondition_finiteelementcalculate(interface_condition,ne,err,error,*999)
2203  CALL interface_matrices_element_add(interface_matrices,err,error,*999)
2204  ENDDO !element_idx
2205 #ifdef TAUPROF
2206  CALL tau_static_phase_stop("Boundary and Ghost Elements Loop")
2207 #endif
2208  !Output timing information if required
2209  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
2210  CALL cpu_timer(user_cpu,user_time5,err,error,*999)
2211  CALL cpu_timer(system_cpu,system_time5,err,error,*999)
2212  user_elapsed=user_time5(1)-user_time4(1)
2213  system_elapsed=system_time5(1)-system_time4(1)
2214  element_user_elapsed=element_user_elapsed+user_elapsed
2215  element_system_elapsed=element_system_elapsed+user_elapsed
2216  CALL write_string_value(general_output_type,"User time for boundary+ghost interface assembly = ",user_elapsed, &
2217  & err,error,*999)
2218  CALL write_string_value(general_output_type,"System time for boundary+ghost interface assembly = ", &
2219  & system_elapsed,err,error,*999)
2220  IF(number_of_times>0) THEN
2221  CALL write_string_value(general_output_type,"Average element user time for interface assembly = ", &
2222  & element_user_elapsed/number_of_times,err,error,*999)
2223  CALL write_string_value(general_output_type,"Average element system time for interface assembly = ", &
2224  & element_system_elapsed/number_of_times,err,error,*999)
2225  ENDIF
2226  ENDIF
2227  !Finalise the element matrices
2228 #ifdef TAUPROF
2229  CALL tau_static_phase_start("INTERFACE_MATRICES_ELEMENT_FINALISE()")
2230 #endif
2231  CALL interface_matrices_element_finalise(interface_matrices,err,error,*999)
2232 #ifdef TAUPROF
2233  CALL tau_static_phase_stop("INTERFACE_MATRICES_ELEMENT_FINALISE()")
2234 #endif
2235  !Output equations matrices and RHS vector if required
2236  IF(interface_equations%OUTPUT_TYPE>=interface_equations_matrix_output) THEN
2237  CALL interface_matrices_output(general_output_type,interface_matrices,err,error,*999)
2238  ENDIF
2239  !Output timing information if required
2240  IF(interface_equations%OUTPUT_TYPE>=interface_equations_timing_output) THEN
2241  CALL cpu_timer(user_cpu,user_time6,err,error,*999)
2242  CALL cpu_timer(system_cpu,system_time6,err,error,*999)
2243  user_elapsed=user_time6(1)-user_time1(1)
2244  system_elapsed=system_time6(1)-system_time1(1)
2245  CALL write_string(general_output_type,"***",err,error,*999)
2246  CALL write_string_value(general_output_type,"Total user time for interface equations assembly = ",user_elapsed, &
2247  & err,error,*999)
2248  CALL write_string_value(general_output_type,"Total system time for interface equations assembly = ", &
2249  & system_elapsed,err,error,*999)
2250  CALL write_string(general_output_type,"***",err,error,*999)
2251  ENDIF
2252  ELSE
2253  CALL flagerror("Interface matrices is not associated.",err,error,*999)
2254  ENDIF
2255  ELSE
2256  CALL flagerror("Interface equations is not associated.",err,error,*999)
2257  ENDIF
2258  ELSE
2259  CALL flagerror("Interface condition Lagrange field is not associated.",err,error,*999)
2260  ENDIF
2261  ELSE
2262  CALL flagerror("Interface condition Lagrange is not associated",err,error,*999)
2263  ENDIF
2264  ELSE
2265  CALL flagerror("Interface condition is not associated.",err,error,*999)
2266  ENDIF
2267 
2268  exits("INTERFACE_CONDITION_RESIDUAL_EVALUATE_FEM")
2269  RETURN
2270 999 errorsexits("INTERFACE_CONDITION_RESIDUAL_EVALUATE_FEM",err,error)
2271  RETURN 1
2272 
2273  END SUBROUTINE interface_condition_residual_evaluate_fem
2274 
2275  !
2276  !================================================================================================================================
2277  !
2278 
2280  SUBROUTINE interfacecondition_finiteelementcalculate(interfaceCondition,interfaceElementNumber,err,error,*)
2281 
2282  !Argument variables
2283  TYPE(interface_condition_type), POINTER :: interfacecondition
2284  INTEGER(INTG), INTENT(IN) :: interfaceelementnumber
2285  INTEGER(INTG), INTENT(OUT) :: err
2286  TYPE(varying_string), INTENT(OUT) :: error
2287  !Local Variables
2288  TYPE(interface_equations_type), POINTER :: interfaceequations
2289  TYPE(interface_matrices_type), POINTER :: interfacematrices
2290  TYPE(element_matrix_type), POINTER :: elementmatrix
2291  INTEGER(INTG) :: interfacematrixidx
2292  TYPE(varying_string) :: localerror
2293 
2294 #ifdef TAUPROF
2295  CALL tau_static_phase_start("InterfaceCondition_FiniteElementCalculate")
2296 #endif
2297 
2298  enters("InterfaceCondition_FiniteElementCalculate",err,error,*999)
2299 
2300  IF(ASSOCIATED(interfacecondition)) THEN
2301  interfaceequations=>interfacecondition%INTERFACE_EQUATIONS
2302  IF(ASSOCIATED(interfaceequations)) THEN
2303  SELECT CASE(interfacecondition%OPERATOR)
2305  CALL fieldcontinuity_finiteelementcalculate(interfacecondition,interfaceelementnumber,err,error,*999)
2307  CALL flagerror("Not implemented!",err,error,*999)
2309  CALL frictionlesscontact_finiteelementcalculate(interfacecondition,interfaceelementnumber,err,error,*999)
2311  CALL solidfluidoperator_finiteelementcalculate(interfacecondition,interfaceelementnumber,err,error,*999)
2312  !CALL FlagError("Not implemented!",ERR,ERROR,*999)
2314  CALL flagerror("Not implemented!",err,error,*999)
2315  CASE DEFAULT
2316  localerror="The interface condition operator of "//trim(number_to_vstring(interfacecondition%OPERATOR,"*",err,error))// &
2317  & " is invalid."
2318  CALL flagerror(localerror,err,error,*999)
2319  END SELECT
2320 
2321  IF(interfaceequations%OUTPUT_TYPE>=interface_equations_element_matrix_output) THEN
2322  interfacematrices=>interfaceequations%INTERFACE_MATRICES
2323  IF(ASSOCIATED(interfacematrices)) THEN
2324  CALL write_string(general_output_type,"Finite element interface matrices:",err,error,*999)
2325  CALL write_string_value(general_output_type,"Element number = ",interfaceelementnumber,err,error,*999)
2326  CALL write_string_value(general_output_type,"Number of element matrices = ",interfacematrices% &
2327  & number_of_interface_matrices,err,error,*999)
2328  DO interfacematrixidx=1,interfacematrices%NUMBER_OF_INTERFACE_MATRICES
2329  CALL write_string_value(general_output_type,"Element matrix : ",interfacematrixidx,err,error,*999)
2330  CALL write_string_value(general_output_type," Update matrix = ",interfacematrices%MATRICES(interfacematrixidx)%PTR% &
2331  & update_matrix,err,error,*999)
2332  IF(interfacematrices%MATRICES(interfacematrixidx)%PTR%UPDATE_MATRIX) THEN
2333  elementmatrix=>interfacematrices%MATRICES(interfacematrixidx)%PTR%ELEMENT_MATRIX
2334  CALL write_string_value(general_output_type," Number of rows = ",elementmatrix%NUMBER_OF_ROWS,err,error,*999)
2335  CALL write_string_value(general_output_type," Number of columns = ",elementmatrix%NUMBER_OF_COLUMNS, &
2336  & err,error,*999)
2337  CALL write_string_value(general_output_type," Maximum number of rows = ",elementmatrix%MAX_NUMBER_OF_ROWS, &
2338  & err,error,*999)
2339  CALL write_string_value(general_output_type," Maximum number of columns = ",elementmatrix% &
2340  & max_number_of_columns,err,error,*999)
2341  CALL write_string_vector(general_output_type,1,1,elementmatrix%NUMBER_OF_ROWS,8,8,elementmatrix%ROW_DOFS, &
2342  & '(" Row dofs :",8(X,I13))','(16X,8(X,I13))',err,error,*999)
2343  CALL write_string_vector(general_output_type,1,1,elementmatrix%NUMBER_OF_COLUMNS,8,8,elementmatrix% &
2344  & column_dofs,'(" Column dofs :",8(X,I13))','(16X,8(X,I13))',err,error,*999)
2345  CALL write_string_matrix(general_output_type,1,1,elementmatrix%NUMBER_OF_ROWS,1,1,elementmatrix% &
2346  & number_of_columns,8,8,elementmatrix%MATRIX(1:elementmatrix%NUMBER_OF_ROWS,1:elementmatrix% &
2347  & number_of_columns),write_string_matrix_name_and_indices,'(" Matrix','(",I2,",:)',' :",8(X,E13.6))', &
2348  & '(16X,8(X,E13.6))',err,error,*999)
2349  ENDIF
2350  ENDDO !interfaceMatrixIdx
2351  ENDIF
2352  ENDIF
2353  ELSE
2354  CALL flagerror("Interface equations is not associated.",err,error,*999)
2355  ENDIF
2356  ELSE
2357  CALL flagerror("Interface condition is not associated.",err,error,*999)
2358  ENDIF
2359 
2360 #ifdef TAUPROF
2361  CALL tau_static_phase_stop("InterfaceCondition_FiniteElementCalculate")
2362 #endif
2363 
2364  exits("InterfaceCondition_FiniteElementCalculate")
2365  RETURN
2366 999 errorsexits("InterfaceCondition_FiniteElementCalculate",err,error)
2367  RETURN 1
2368 
2369  END SUBROUTINE interfacecondition_finiteelementcalculate
2370 
2371  !
2372  !================================================================================================================================
2373  !
2374 
2376  SUBROUTINE interface_condition_user_number_find(USER_NUMBER,INTERFACE,INTERFACE_CONDITION,ERR,ERROR,*)
2377 
2378  !Argument variables
2379  INTEGER(INTG), INTENT(IN) :: user_number
2380  TYPE(interface_type), POINTER :: interface
2381  TYPE(interface_condition_type), POINTER :: interface_condition
2382  INTEGER(INTG), INTENT(OUT) :: err
2383  TYPE(varying_string), INTENT(OUT) :: error
2384  !Local Variables
2385  INTEGER(INTG) :: interface_condition_idx
2386  TYPE(varying_string) :: local_error
2387 
2388  enters("INTERFACE_CONDITION_USER_NUMBER_FIND",err,error,*999)
2389 
2390  IF(ASSOCIATED(interface)) THEN
2391  IF(ASSOCIATED(interface_condition)) THEN
2392  CALL flagerror("Interface condition is already associated.",err,error,*999)
2393  ELSE
2394  NULLIFY(interface_condition)
2395  IF(ASSOCIATED(interface%INTERFACE_CONDITIONS)) THEN
2396  interface_condition_idx=1
2397  DO WHILE(interface_condition_idx<=interface%INTERFACE_CONDITIONS%NUMBER_OF_INTERFACE_CONDITIONS.AND. &
2398  & .NOT.ASSOCIATED(interface_condition))
2399  IF(interface%INTERFACE_CONDITIONS%INTERFACE_CONDITIONS(interface_condition_idx)%PTR%USER_NUMBER==user_number) THEN
2400  interface_condition=>interface%INTERFACE_CONDITIONS%INTERFACE_CONDITIONS(interface_condition_idx)%PTR
2401  ELSE
2402  interface_condition_idx=interface_condition_idx+1
2403  ENDIF
2404  ENDDO
2405  ELSE
2406  local_error="The interface conditions on interface number "// &
2407  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//" are not associated."
2408  CALL flagerror(local_error,err,error,*999)
2409  ENDIF
2410  ENDIF
2411  ELSE
2412  CALL flagerror("Interface is not associated.",err,error,*999)
2413  ENDIF
2414 
2415  exits("INTERFACE_CONDITION_USER_NUMBER_FIND")
2416  RETURN
2417 999 errorsexits("INTERFACE_CONDITION_USER_NUMBER_FIND",err,error)
2418  RETURN 1
2419  END SUBROUTINE interface_condition_user_number_find
2420 
2421  !
2422  !================================================================================================================================
2423  !
2424 
2426  SUBROUTINE interface_conditions_finalise(INTERFACE_CONDITIONS,ERR,ERROR,*)
2427 
2428  !Argument variables
2429  TYPE(interface_conditions_type), POINTER :: interface_conditions
2430  INTEGER(INTG), INTENT(OUT) :: err
2431  TYPE(varying_string), INTENT(OUT) :: error
2432  !Local Variables
2433  TYPE(interface_condition_type), POINTER :: interface_condition
2434 
2435  enters("INTERFACE_CONDITIONS_FINALISE",err,error,*999)
2436 
2437  IF(ASSOCIATED(interface_conditions)) THEN
2438  DO WHILE(interface_conditions%NUMBER_OF_INTERFACE_CONDITIONS>0)
2439  interface_condition=>interface_conditions%INTERFACE_CONDITIONS(1)%PTR
2440  CALL interface_condition_destroy(interface_condition,err,error,*999)
2441  ENDDO
2442  IF(ASSOCIATED(interface_conditions%INTERFACE_CONDITIONS)) DEALLOCATE(interface_conditions%INTERFACE_CONDITIONS)
2443  DEALLOCATE(interface_conditions)
2444  ENDIF
2445 
2446  exits("INTERFACE_CONDITIONS_FINALISE")
2447  RETURN
2448 999 errorsexits("INTERFACE_CONDITIONS_FINALISE",err,error)
2449  RETURN 1
2450  END SUBROUTINE interface_conditions_finalise
2451 
2452  !
2453  !================================================================================================================================
2454  !
2455 
2457  SUBROUTINE interface_conditions_initialise(INTERFACE,ERR,ERROR,*)
2458 
2459  !Argument variables
2460  TYPE(interface_type), POINTER :: interface
2461  INTEGER(INTG), INTENT(OUT) :: err
2462  TYPE(varying_string), INTENT(OUT) :: error
2463  !Local Variables
2464  INTEGER(INTG) :: dummy_err
2465  TYPE(varying_string) :: dummy_error,local_error
2466 
2467  enters("INTERFACE_CONDITIONS_INITIALISE",err,error,*998)
2468 
2469  IF(ASSOCIATED(interface)) THEN
2470  IF(ASSOCIATED(interface%INTERFACE_CONDITIONS)) THEN
2471  local_error="Interface conditions is already associated for interface number "// &
2472  & trim(number_to_vstring(interface%USER_NUMBER,"*",err,error))//"."
2473  CALL flagerror(local_error,err,error,*999)
2474  ELSE
2475  ALLOCATE(interface%INTERFACE_CONDITIONS,stat=err)
2476  IF(err/=0) CALL flagerror("Could not allocate interface interface conditions.",err,error,*999)
2477  interface%INTERFACE_CONDITIONS%INTERFACE=>INTERFACE
2478  interface%INTERFACE_CONDITIONS%NUMBER_OF_INTERFACE_CONDITIONS=0
2479  NULLIFY(interface%INTERFACE_CONDITIONS%INTERFACE_CONDITIONS)
2480  ENDIF
2481  ELSE
2482  CALL flagerror("Interface is not associated.",err,error,*998)
2483  ENDIF
2484 
2485  exits("INTERFACE_CONDITIONS_INITIALISE")
2486  RETURN
2487 999 CALL interface_conditions_finalise(interface%INTERFACE_CONDITIONS,dummy_err,dummy_error,*998)
2488 998 errorsexits("INTERFACE_CONDITIONS_INITIALISE",err,error)
2489  RETURN 1
2490  END SUBROUTINE interface_conditions_initialise
2491 
2492  !
2493  !================================================================================================================================
2494  !
2495 
2496 END MODULE interface_conditions_routines
This module contains all basis function routines.
subroutine, public enters(NAME, ERR, ERROR,)
Records the entry into the named procedure and initialises the error code.
Write a string followed by a value to a given output stream.
Contains information for a region.
Definition: types.f90:3252
Converts a number to its equivalent varying string representation.
Definition: strings.f90:161
Contains information on the mesh decomposition.
Definition: types.f90:1063
integer(intg), parameter interface_condition_fls_contact_reproject_operator
Frictionless contact operator, reproject at each newton iteration and has geometric linearisation ter...
integer(intg), parameter interface_condition_lagrange_multipliers_method
Lagrange multipliers interface condition method.
Contains information about the penalty field information for an interface condition.
Definition: types.f90:2129
Contains information on an equations set.
Definition: types.f90:1941
This module contains all string manipulation and transformation routines.
Definition: strings.f90:45
Contains information for the interface condition data.
Definition: types.f90:2155
This module contains routines for timing the program.
Definition: timer_f.f90:45
Write a string followed by a matrix to a specified output stream.
Contains information for a field defined on a region.
Definition: types.f90:1346
This module provides an iso_varying_string module, conformant to the API specified in ISO/IEC 1539-2:...
Contains information on an interface mapping. TODO: Generalise to non-Lagrange multipler mappings...
Definition: types.f90:2065
integer(intg), parameter interface_condition_augmented_lagrange_method
Augmented Lagrange multiplers interface condition method.
integer(intg), parameter interface_condition_gauss_integration
Gauss points integration type, i.e. Loop over element Gauss points and sum up their contribution...
This module contains all interface mapping routines.
integer(intg), parameter, public user_cpu
User CPU time type.
Definition: timer_f.f90:68
integer(intg), parameter interface_condition_field_normal_continuity_operator
Continuous field normal operator, i.e., lambda(u_1.n_1-u_2.n_2).
Contains information about the dependent field information for an interface condition.
Definition: types.f90:2146
A buffer type to allow for an array of pointers to a EQUATIONS_SET_TYPE.
Definition: types.f90:1962
integer(intg), parameter, public interface_equations_timing_output
Timing information output.
integer(intg), parameter interface_condition_solid_fluid_normal_operator
Solid fluid normal operator, i.e., lambda(v_f.n_f-du_s/dt.n_s).
subroutine, public exits(NAME)
Records the exit out of the named procedure.
This module contains all type definitions in order to avoid cyclic module references.
Definition: types.f90:70
integer(intg), parameter interface_condition_data_points_integration
Data points integration type i.e. Loop over data points and sum up their contribution.
integer(intg), parameter, public interface_equations_element_matrix_output
All below and element matrices output .
Write a string to a given output stream.
This module contains all the low-level base routines e.g., all debug, control, and low-level communic...
integer(intg), parameter, public general_output_type
General output type.
Contains information about the Lagrange field information for an interface condition.
Definition: types.f90:2137
integer(intg), parameter interface_condition_penalty_method
Penalty interface condition method.
A buffer type to allow for an array of pointers to a FIELD_VARIABLE_TYPE.
Definition: types.f90:1311
integer(intg), parameter, public system_cpu
System CPU time type.
Definition: timer_f.f90:69
integer(intg), parameter, public write_string_matrix_name_and_indices
Write the matrix name together with the matrix indices.
Contains information for an element matrix.
Definition: types.f90:1387
integer(intg), parameter, public matrix_compressed_row_storage_type
Matrix compressed row storage type.
Contains information on a mesh defined on a region.
Definition: types.f90:503
integer(intg), parameter interface_condition_fls_contact_operator
Frictionless contact operator, i.e., lambda.(x_1.n-x_2.n).
A buffer type to allow for an array of pointers to a INTERFACE_CONDITION_TYPE.
Definition: types.f90:2173
This module defines all constants shared across interface condition routines.
Contains information for interface region specific data that is not of &#39;region&#39; importance. <<>>
Definition: types.f90:2178
integer(intg), parameter interface_condition_point_to_point_method
Point to point interface condition method.
subroutine, public cpu_timer(TIME_TYPE, TIME, ERR, ERROR,)
CPU_TIMER returns the CPU time in TIME(1). TIME_TYPE indicates the type of time required.
Definition: timer_f.f90:99
This module handles all interface equations routines.
This module contains all routines dealing with (non-distributed) matrix and vectors types...
Write a string followed by a vector to a specified output stream.
integer(intg), parameter interface_condition_field_continuity_operator
Continuous field operator, i.e., lambda.(u_1-u_2).
Contains information for a field variable defined on a field.
Definition: types.f90:1289
Contains information on the domain mappings (i.e., local and global numberings).
Definition: types.f90:904
Contains information on the interface matrices.
Definition: types.f90:2012
subroutine, public errors(NAME, ERR, ERROR)
Records the exiting error of the subroutine.
Contains information for the interface data.
Definition: types.f90:2228
integer(intg), parameter, public matrix_block_storage_type
Matrix block storage type.
Contains information on the geometry for an interface condition.
Definition: types.f90:2123
integer(intg), parameter interface_condition_solid_fluid_operator
Solid fluid operator, i.e., lambda.(v_f-du_s/dt).
Flags an error condition.
integer(intg), parameter, public interface_equations_matrix_output
All below and equation matrices output.
This module contains all kind definitions.
Definition: kinds.f90:45
Contains information about the interface equations for an interface condition.
Definition: types.f90:2110
This module handles all formating and input and output.